29 #ifndef OOMPH_DG_ELEMENT_HEADER
30 #define OOMPH_DG_ELEMENT_HEADER
34 #include <oomph-lib-config.h>
126 std::ostringstream error_stream;
128 <<
"Empty numerical flux function called\n"
129 <<
"This function should be overloaded with a specific flux\n"
130 <<
"that is appropriate to the equations being solved.\n";
132 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
218 if (this->Average_value != 0)
266 this->M_pt = element_pt->
M_pt;
292 std::vector<bool>& boundary_flag,
296 const unsigned n_node = this->
nnode();
298 if (boundary_flag.size() != n_node)
300 std::ostringstream error_stream;
302 <<
"Size of boundary_flag vector is " << boundary_flag.size() <<
"\n"
303 <<
"It must be the same as the number of nodes in the element "
307 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
312 for (
unsigned n = 0; n < n_node; n++)
315 if (boundary_flag[n])
330 DG_mesh_pt = mesh_pt;
339 const unsigned n_node = this->
nnode();
340 for (
unsigned n = 0; n < n_node; n++)
350 DG_mesh_pt = mesh_pt;
375 unsigned n_face =
nface();
376 for (
unsigned f = 0; f < n_face; f++)
384 const int& face_index,
395 unsigned n_face = this->
nface();
396 for (
unsigned f = 0; f < n_face; f++)
411 unsigned n_face = this->
nface();
412 for (
unsigned f = 0; f < n_face; f++)
425 OOMPH_CURRENT_FUNCTION,
426 OOMPH_EXCEPTION_LOCATION);
441 OOMPH_CURRENT_FUNCTION,
442 OOMPH_EXCEPTION_LOCATION);
454 OOMPH_CURRENT_FUNCTION,
455 OOMPH_EXCEPTION_LOCATION);
471 const int& face_index,
476 std::string error_message =
"Empty neighbour_finder() has been called.\n";
478 "This function is implemented in the base class of a DGMesh.\n";
479 error_message +=
"It must be overloaded in a specific DGMesh\n";
482 error_message, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
490 const bool& add_face_data_as_external =
false)
493 const unsigned n_element = this->
nelement();
494 for (
unsigned e = 0;
e < n_element;
e++)
505 const unsigned n_element = this->
nelement();
506 for (
unsigned e = 0;
e < n_element;
e++)
512 for (
unsigned e = 0;
e < n_element;
e++)
515 ->slope_limit(slope_limiter_pt);
538 OOMPH_CURRENT_FUNCTION,
539 OOMPH_EXCEPTION_LOCATION);
572 void limit(
const unsigned&
i,
A Base class for DGElements.
bool Can_delete_mass_matrix
Boolean flag to indicate whether the mass matrix can be deleted (i.e. was it created by this element)
void disable_mass_matrix_reuse()
Function that disables the reuse of the mass matrix.
void calculate_averages()
Calculate the elemental averages.
DGMesh * DG_mesh_pt
Pointer to Mesh, which will be responsible for the neighbour finding.
void setup_face_neighbour_info(const bool &add_face_data_as_external)
The boolean flag determines whether the data from the neighbouring elements is added as external data...
virtual void calculate_element_averages(double *&average_values)
Calculate the averages in the element.
virtual void build_all_faces()=0
void pre_compute_mass_matrix()
Function that computes and stores the (inverse) mass matrix.
bool mass_matrix_has_been_computed()
Access function for the boolean to indicate whether the mass matrix has been computed.
void add_flux_contributions_to_residuals(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Loop over all faces and add their integrated numerical fluxes to the residuals.
void enable_mass_matrix_reuse()
Function that allows the reuse of the mass matrix.
double * Average_value
Pointer to storage for the average values of the of the variables over the element.
void slope_limit(SlopeLimiter *const &slope_limiter_pt)
Limit the slope within the element.
DGElement()
Constructor, initialise the pointers to zero.
void construct_nodes_and_faces(DGMesh *const &mesh_pt, TimeStepper *const &time_stepper_pt)
Construct the nodes and faces of the element.
Vector< FaceElement * > Face_element_pt
Vector of pointers to faces of the element.
unsigned nface() const
Return the number of faces.
virtual void set_mass_matrix_from_element(DGElement *const &element_pt)
Set the mass matrix to point to one in another element.
virtual unsigned required_nflux()
Set the number of flux components.
void set_mesh_pt(DGMesh *&mesh_pt)
bool Mass_matrix_has_been_computed
Boolean flag to indicate whether the mass matrix has been computed.
double & average_value(const unsigned &i)
Return the average values.
void construct_boundary_nodes_and_faces(DGMesh *const &mesh_pt, std::vector< bool > &boundary_flag, TimeStepper *const &time_stepper_pt)
Construct all nodes and faces of the element. The vector of booleans boundary should be the same size...
bool Mass_matrix_reuse_is_enabled
Boolean flag to indicate whether to reuse the mass matrix.
void output_faces(std::ostream &outfile)
Output the faces of the element.
DGFaceElement * face_element_pt(const unsigned &i)
Access function for the faces.
void get_neighbouring_face_and_local_coordinate(const int &face_index, const Vector< double > &s, FaceElement *&face_element_pt, Vector< double > &s_face)
Return the neighbour info.
virtual ~DGElement()
Virtual destructor, destroy the mass matrix, if we created it Clean-up storage associated with averag...
virtual void get_inverse_mass_matrix_times_residuals(Vector< double > &minv_res)
Function that returns the current value of the residuals multiplied by the inverse mass matrix (virtu...
const double & average_value(const unsigned &i) const
Return the average values.
DenseDoubleMatrix * M_pt
Pointer to storage for a mass matrix that can be recycled if desired.
Base class for Discontinuous Galerkin Faces. These are responsible for calculating the normal fluxes ...
virtual void numerical_flux(const Vector< double > &n_out, const Vector< double > &u_int, const Vector< double > &u_ext, Vector< double > &flux)
Calculate the normal flux, which is the dot product of our approximation to the flux with the outer u...
virtual ~DGFaceElement()
Empty Destructor.
virtual void get_interpolation_data(Vector< Data * > &interpolation_data)
Get the data that are used to interpolate the unkowns in the element. These must be returned in order...
void setup_neighbour_info(const bool &add_neighbour_data_to_bulk)
Setup information from the neighbouring face, return a set of nodes (as data) in the neighour....
void report_info()
Output information about the present element and its neighbour.
Vector< FaceElement * > Neighbour_face_pt
Vector of neighbouring face elements at the integration points.
void add_flux_contributions(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Add the contribution from integrating the numerical flux.
DGFaceElement()
Empty Constructor.
Vector< Vector< unsigned > > Neighbour_external_data
Vector of the vectors that will store the number of the bulk external data that correspond to the dof...
virtual void dnumerical_flux_du(const Vector< double > &n_out, const Vector< double > &u_int, const Vector< double > &u_ext, DenseMatrix< double > &dflux_du_int, DenseMatrix< double > &dflux_du_ext)
Calculate the derivative of the normal flux, which is the dot product of our approximation to the flu...
virtual void interpolated_u(const Vector< double > &s, Vector< double > &f)
Return the interpolated values of the unknown fluxes.
virtual void numerical_flux_at_knot(const unsigned &ipt, const Shape &psi, Vector< double > &flux, DenseMatrix< double > &dflux_du_int, DenseMatrix< double > &dflux_du_ext, unsigned flag)
Calculate the normal numerical flux at the integration point. This is the most general interface that...
Vector< Vector< double > > Neighbour_local_coordinate
Vector of neighbouring local coordinates at the integration points.
virtual unsigned flux_index(const unsigned &i) const
Return the index at which the i-th unknown flux is stored.
FaceElement * neighbour_face_pt(const unsigned &i)
Access function for neighbouring face information.
virtual unsigned required_nflux()
Set the number of flux components.
void limit_slopes(SlopeLimiter *const &slope_limiter_pt)
virtual void neighbour_finder(FiniteElement *const &bulk_element_pt, const int &face_index, const Vector< double > &s_bulk, FaceElement *&face_element_pt, Vector< double > &s_face)
void setup_face_neighbour_info(const bool &add_face_data_as_external=false)
If the boolean flag is set to true then the data from neighbouring faces will be added as external da...
static double FaceTolerance
Class of matrices containing doubles, and stored as a DenseMatrix<double>, but with solving functiona...
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
A general Finite Element class.
virtual Node * construct_node(const unsigned &n)
Construct the local node n and return a pointer to the newly created node object.
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
unsigned nnode() const
Return the number of nodes.
virtual Node * construct_boundary_node(const unsigned &n)
Construct the local node n as a boundary node; that is a node that MAY be placed on a mesh boundary a...
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
const Vector< GeneralisedElement * > & element_pt() const
Return reference to the Vector of elements.
unsigned long nelement() const
Return number of elements in the mesh.
double minmod(Vector< double > &args)
The basic minmod function.
bool MUSCL
Boolean flag to indicate a MUSCL or straight min-mod limiter.
double minmodB(Vector< double > &args, const double &h)
The modified minmod function.
MinModLimiter(const double &m=0.0, const bool &muscl=false)
Constructor takes a value for the modification parameter M (default to zero — classic min mod) and a ...
virtual ~MinModLimiter()
Empty destructor.
void limit(const unsigned &i, const Vector< DGElement * > &required_element_pt)
The limit function.
double M
Constant that is used in the modified min-mod function to provide better properties at extrema.
An OomphLibError object which should be thrown when an run-time error is encountered....
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Base class for slope limiters.
virtual ~SlopeLimiter()
virtual destructor
SlopeLimiter()
Empty constructor.
virtual void limit(const unsigned &i, const Vector< DGElement * > &required_element_pt)
Basic function.
////////////////////////////////////////////////////////////////////// //////////////////////////////...
A slight extension to the standard template vector class so that we can include "graceful" array rang...
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...