dg_elements.h
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2023 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 // Header functions for classes that define Discontinuous Galerkin elements
27 
28 // Include guards to prevent multiple inclusion of the header
29 #ifndef OOMPH_DG_ELEMENT_HEADER
30 #define OOMPH_DG_ELEMENT_HEADER
31 
32 // Config header generated by autoconfig
33 #ifdef HAVE_CONFIG_H
34 #include <oomph-lib-config.h>
35 #endif
36 
37 // Oomph-lib header
38 #include "elements.h"
39 #include "mesh.h"
40 
41 namespace oomph
42 {
43  //=============================================================
44  /// Base class for Discontinuous Galerkin Faces.
45  /// These are responsible for calculating the normal fluxes
46  /// that provide the communication between the discontinuous
47  /// elements.
48  //===============================================================
49  class DGFaceElement : public virtual FaceElement
50  {
51  /// Vector of neighbouring face elements at the integration points
53 
54  /// Vector of neighbouring local coordinates at the integration points
56 
57  /// Vector of the vectors that will store the number of the
58  /// bulk external data that correspond to the dofs in the neighbouring face
59  /// This is only used if we are using implict timestepping and wish
60  /// to assemble the jacobian matrix. In order that the data be correctly
61  /// set up DGMesh::setup_face_neighbour_info() must be called with the
62  /// boolean flag set to true.
64 
65  protected:
66  /// Return the index at which the i-th unknown flux is stored.
67  // The default return is suitable for single-physics problem
68  virtual inline unsigned flux_index(const unsigned& i) const
69  {
70  return i;
71  }
72 
73  /// Set the number of flux components
74  virtual unsigned required_nflux()
75  {
76  return 0;
77  }
78 
79  public:
80  /// Empty Constructor
82 
83  /// Empty Destructor
84  virtual ~DGFaceElement() {}
85 
86  /// Access function for neighbouring face information
87  FaceElement* neighbour_face_pt(const unsigned& i)
88  {
89  return Neighbour_face_pt[i];
90  }
91 
92  /// Setup information from the neighbouring face, return a set
93  /// of nodes (as data) in the neighour. The boolean flag
94  /// is used to determine whether the data in the neighbour are
95  /// added as external data to the bulk element --- required when
96  /// computing the jacobian of the system
97  void setup_neighbour_info(const bool& add_neighbour_data_to_bulk);
98 
99  /// Output information about the present element and its neighbour
100  void report_info();
101 
102  // Get the value of the unknowns
103  virtual void interpolated_u(const Vector<double>& s, Vector<double>& f);
104 
105  /// Get the data that are used to interpolate the unkowns
106  /// in the element. These must be returned in order.
107  virtual void get_interpolation_data(Vector<Data*>& interpolation_data);
108 
109 
110  /// Calculate the normal numerical flux at the integration point.
111  /// This is the most general interface that can be overloaded if desired
112  virtual void numerical_flux_at_knot(const unsigned& ipt,
113  const Shape& psi,
114  Vector<double>& flux,
115  DenseMatrix<double>& dflux_du_int,
116  DenseMatrix<double>& dflux_du_ext,
117  unsigned flag);
118 
119  /// Calculate the normal flux, which is the dot product of our
120  /// approximation to the flux with the outer unit normal
121  virtual void numerical_flux(const Vector<double>& n_out,
122  const Vector<double>& u_int,
123  const Vector<double>& u_ext,
124  Vector<double>& flux)
125  {
126  std::ostringstream error_stream;
127  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";
131  throw OomphLibError(
132  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
133  }
134 
135  /// Calculate the derivative of the
136  /// normal flux, which is the dot product of our
137  /// approximation to the flux with the outer unit normal,
138  /// with respect to the interior and exterior variables
139  /// Default is to use finite differences
140  virtual void dnumerical_flux_du(const Vector<double>& n_out,
141  const Vector<double>& u_int,
142  const Vector<double>& u_ext,
143  DenseMatrix<double>& dflux_du_int,
144  DenseMatrix<double>& dflux_du_ext);
145 
146 
147  /// Add the contribution from integrating the numerical flux
148  // over the face to the residuals
149  void add_flux_contributions(Vector<double>& residuals,
150  DenseMatrix<double>& jacobian,
151  unsigned flag);
152  };
153 
154  class DGMesh;
155  class SlopeLimiter;
156 
157  //==================================================================
158  /// A Base class for DGElements
159  //=================================================================
160  class DGElement : public virtual FiniteElement
161  {
162  protected:
163  /// The DGFaceElement requires access to the nodal local equation
164  /// information, so it's a friend
165  friend class DGFaceElement;
166 
167  /// Vector of pointers to faces of the element
169 
170  /// Pointer to Mesh, which will be responsible for the neighbour finding
172 
173  /// Pointer to storage for a mass matrix that can be recycled if
174  /// desired
176 
177  /// Pointer to storage for the average values of the of the
178  /// variables over the element
179  double* Average_value;
180 
181  /// Boolean flag to indicate whether to reuse the mass matrix
183 
184  /// Boolean flag to indicate whether the mass matrix has been
185  /// computed
187 
188  /// Boolean flag to indicate whether the mass matrix can be
189  /// deleted (i.e. was it created by this element)
191 
192  /// Set the number of flux components
193  virtual unsigned required_nflux()
194  {
195  return 0;
196  }
197 
198  public:
199  /// Constructor, initialise the pointers to zero
201  : DG_mesh_pt(0),
202  M_pt(0),
203  Average_value(0),
207  {
208  }
209 
210  /// Virtual destructor, destroy the mass matrix, if we created it
211  /// Clean-up storage associated with average values
212  virtual ~DGElement()
213  {
214  if ((M_pt != 0) && Can_delete_mass_matrix)
215  {
216  delete M_pt;
217  }
218  if (this->Average_value != 0)
219  {
220  delete[] Average_value;
221  Average_value = 0;
222  }
223  }
224 
225  /// Access function for the boolean to indicate whether the
226  /// mass matrix has been computed
228  {
230  }
231 
232  /// Function that allows the reuse of the mass matrix
234  {
237  }
238 
239  /// Function that disables the reuse of the mass matrix
241  {
242  // If we are using another element's mass matrix then reset our pointer
243  // to zero
245  {
246  M_pt = 0;
247  }
248  // Otherwise we do not reuse the mass matrix
250  // Recalculate the mass matrix
252  }
253 
254 
255  /// Set the mass matrix to point to one in another element
256  virtual void set_mass_matrix_from_element(DGElement* const& element_pt)
257  {
258  // If the element's mass matrix has not been computed, compute it!
259  if (!element_pt->mass_matrix_has_been_computed())
260  {
261  element_pt->pre_compute_mass_matrix();
262  }
263 
264  // Now set the mass matrix in this element to address that
265  // of element_pt
266  this->M_pt = element_pt->M_pt;
267  // We must reuse the mass matrix, or there will be trouble
268  // Because we will recalculate it in the original element
271  // We cannot delete the mass matrix
272  Can_delete_mass_matrix = false;
273  }
274 
275  /// Function that computes and stores the (inverse) mass matrix
277 
278  // Function that is used to construct all the faces of the DGElement
279  virtual void build_all_faces() = 0;
280 
281  /// Function that returns the current value of the residuals
282  /// multiplied by the inverse mass matrix (virtual so that it can be
283  /// overloaded specific elements in which time saving tricks can be applied)
285  Vector<double>& minv_res);
286 
287  /// Construct all nodes and faces of the element.
288  /// The vector of booleans boundary should be the same size
289  /// as the number of nodes and if any entries are true
290  /// that node will be constructed as a boundary node.
292  std::vector<bool>& boundary_flag,
294  {
295  // Construct the nodes (This should not be used in a base class)
296  const unsigned n_node = this->nnode();
297 #ifdef PARANOID
298  if (boundary_flag.size() != n_node)
299  {
300  std::ostringstream error_stream;
301  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 "
304  << n_node << "\n";
305 
306  throw OomphLibError(
307  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
308  }
309 #endif
310 
311  // Loop over the nodes
312  for (unsigned n = 0; n < n_node; n++)
313  {
314  // If the node is marked on a boundary construct a boundary node
315  if (boundary_flag[n])
316  {
318  }
319  // Otherwise, construct a normal node
320  else
321  {
322  (void)this->construct_node(n, time_stepper_pt);
323  }
324  }
325 
326  // Make the faces
327  this->build_all_faces();
328 
329  // Set the Mesh pointer
330  DG_mesh_pt = mesh_pt;
331  }
332 
333 
334  /// Construct the nodes and faces of the element
335  void construct_nodes_and_faces(DGMesh* const& mesh_pt,
337  {
338  // Loop over the nodes
339  const unsigned n_node = this->nnode();
340  for (unsigned n = 0; n < n_node; n++)
341  {
342  // Construct the node and ignore the return value
343  (void)this->construct_node(n, time_stepper_pt);
344  }
345 
346  // Make the faces
347  this->build_all_faces();
348 
349  // Set the Mesh pointer
350  DG_mesh_pt = mesh_pt;
351  }
352 
353  // Set the mesh pointer of the element
354  void set_mesh_pt(DGMesh*& mesh_pt)
355  {
356  DG_mesh_pt = mesh_pt;
357  }
358 
359  /// Return the number of faces
360  unsigned nface() const
361  {
362  return Face_element_pt.size();
363  }
364 
365  /// Access function for the faces
366  DGFaceElement* face_element_pt(const unsigned& i)
367  {
368  return dynamic_cast<DGFaceElement*>(Face_element_pt[i]);
369  }
370 
371  /// Output the faces of the element
372  void output_faces(std::ostream& outfile)
373  {
374  // Loop over the faces
375  unsigned n_face = nface();
376  for (unsigned f = 0; f < n_face; f++)
377  {
378  face_element_pt(f)->output(outfile);
379  }
380  }
381 
382  /// Return the neighbour info
384  const int& face_index,
385  const Vector<double>& s,
387  Vector<double>& s_face);
388 
389  // Setup the face information
390  /// The boolean flag determines whether the data from the neighbouring
391  /// elements is added as external data to the element (required for
392  /// correct computation of the jacobian)
393  void setup_face_neighbour_info(const bool& add_face_data_as_external)
394  {
395  unsigned n_face = this->nface();
396  for (unsigned f = 0; f < n_face; f++)
397  {
398  face_element_pt(f)->setup_neighbour_info(add_face_data_as_external);
400  }
401  }
402 
403 
404  /// Loop over all faces and add their integrated numerical fluxes
405  /// to the residuals
407  DenseMatrix<double>& jacobian,
408  unsigned flag)
409  {
410  // Add up the contributions from each face
411  unsigned n_face = this->nface();
412  for (unsigned f = 0; f < n_face; f++)
413  {
414  face_element_pt(f)->add_flux_contributions(residuals, jacobian, flag);
415  }
416  }
417 
418  /// Limit the slope within the element
419  void slope_limit(SlopeLimiter* const& slope_limiter_pt);
420 
421  /// Calculate the averages in the element
422  virtual void calculate_element_averages(double*& average_values)
423  {
424  throw OomphLibError("Default (empty) version called",
425  OOMPH_CURRENT_FUNCTION,
426  OOMPH_EXCEPTION_LOCATION);
427  }
428 
429  /// Calculate the elemental averages
431  {
432  this->calculate_element_averages(this->Average_value);
433  }
434 
435  /// Return the average values
436  double& average_value(const unsigned& i)
437  {
438  if (Average_value == 0)
439  {
440  throw OomphLibError("Averages not calculated yet",
441  OOMPH_CURRENT_FUNCTION,
442  OOMPH_EXCEPTION_LOCATION);
443  }
444  return Average_value[i];
445  }
446 
447 
448  /// Return the average values
449  const double& average_value(const unsigned& i) const
450  {
451  if (Average_value == 0)
452  {
453  throw OomphLibError("Averages not calculated yet",
454  OOMPH_CURRENT_FUNCTION,
455  OOMPH_EXCEPTION_LOCATION);
456  }
457  return Average_value[i];
458  }
459  };
460 
461  class DGMesh : public Mesh
462  {
463  public:
464  static double FaceTolerance;
465 
466  DGMesh() : Mesh() {}
467 
468  virtual ~DGMesh() {}
469 
470  virtual void neighbour_finder(FiniteElement* const& bulk_element_pt,
471  const int& face_index,
472  const Vector<double>& s_bulk,
473  FaceElement*& face_element_pt,
474  Vector<double>& s_face)
475  {
476  std::string error_message = "Empty neighbour_finder() has been called.\n";
477  error_message +=
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";
480 
481  throw OomphLibError(
482  error_message, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
483  }
484 
485 
486  // Setup the face information for all the elements
487  /// If the boolean flag is set to true then the data from neighbouring
488  /// faces will be added as external data to the bulk element.
490  const bool& add_face_data_as_external = false)
491  {
492  // Loop over all the elements and setup their face neighbour information
493  const unsigned n_element = this->nelement();
494  for (unsigned e = 0; e < n_element; e++)
495  {
496  dynamic_cast<DGElement*>(this->element_pt(e))
497  ->setup_face_neighbour_info(add_face_data_as_external);
498  }
499  }
500 
501  // Limit the slopes on the entire mesh
502  void limit_slopes(SlopeLimiter* const& slope_limiter_pt)
503  {
504  // Loop over all the elements and calculate the averages
505  const unsigned n_element = this->nelement();
506  for (unsigned e = 0; e < n_element; e++)
507  {
508  dynamic_cast<DGElement*>(this->element_pt(e))->calculate_averages();
509  }
510 
511  // Now loop over again and limit the values
512  for (unsigned e = 0; e < n_element; e++)
513  {
514  dynamic_cast<DGElement*>(this->element_pt(e))
515  ->slope_limit(slope_limiter_pt);
516  }
517  }
518  };
519 
520 
521  //======================================================
522  /// Base class for slope limiters
523  //=====================================================
525  {
526  public:
527  /// Empty constructor
529 
530  /// virtual destructor
531  virtual ~SlopeLimiter() {}
532 
533  /// Basic function
534  virtual void limit(const unsigned& i,
535  const Vector<DGElement*>& required_element_pt)
536  {
537  throw OomphLibError("Calling default empty limiter\n",
538  OOMPH_CURRENT_FUNCTION,
539  OOMPH_EXCEPTION_LOCATION);
540  }
541  };
542 
544  {
545  /// Constant that is used in the modified min-mod function to
546  /// provide better properties at extrema
547  double M;
548 
549  /// Boolean flag to indicate a MUSCL or straight min-mod limiter
550  bool MUSCL;
551 
552  public:
553  /// Constructor takes a value for the modification parameter M
554  /// (default to zero --- classic min mod) and a flag to indicate whether
555  /// we use MUSCL limiting or not --- default false
556  MinModLimiter(const double& m = 0.0, const bool& muscl = false)
557  : SlopeLimiter(), M(m), MUSCL(muscl)
558  {
559  }
560 
561  /// Empty destructor
562  virtual ~MinModLimiter() {}
563 
564  /// The basic minmod function
565  double minmod(Vector<double>& args);
566 
567 
568  /// The modified minmod function
569  double minmodB(Vector<double>& args, const double& h);
570 
571  /// The limit function
572  void limit(const unsigned& i,
573  const Vector<DGElement*>& required_element_pt);
574  };
575 
576 
577 } // namespace oomph
578 #endif
e
Definition: cfortran.h:571
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
A Base class for DGElements.
Definition: dg_elements.h:161
bool Can_delete_mass_matrix
Boolean flag to indicate whether the mass matrix can be deleted (i.e. was it created by this element)
Definition: dg_elements.h:190
void disable_mass_matrix_reuse()
Function that disables the reuse of the mass matrix.
Definition: dg_elements.h:240
void calculate_averages()
Calculate the elemental averages.
Definition: dg_elements.h:430
DGMesh * DG_mesh_pt
Pointer to Mesh, which will be responsible for the neighbour finding.
Definition: dg_elements.h:171
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...
Definition: dg_elements.h:393
virtual void calculate_element_averages(double *&average_values)
Calculate the averages in the element.
Definition: dg_elements.h:422
virtual void build_all_faces()=0
void pre_compute_mass_matrix()
Function that computes and stores the (inverse) mass matrix.
Definition: dg_elements.cc:578
bool mass_matrix_has_been_computed()
Access function for the boolean to indicate whether the mass matrix has been computed.
Definition: dg_elements.h:227
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.
Definition: dg_elements.h:406
void enable_mass_matrix_reuse()
Function that allows the reuse of the mass matrix.
Definition: dg_elements.h:233
double * Average_value
Pointer to storage for the average values of the of the variables over the element.
Definition: dg_elements.h:179
void slope_limit(SlopeLimiter *const &slope_limiter_pt)
Limit the slope within the element.
Definition: dg_elements.cc:686
DGElement()
Constructor, initialise the pointers to zero.
Definition: dg_elements.h:200
void construct_nodes_and_faces(DGMesh *const &mesh_pt, TimeStepper *const &time_stepper_pt)
Construct the nodes and faces of the element.
Definition: dg_elements.h:335
Vector< FaceElement * > Face_element_pt
Vector of pointers to faces of the element.
Definition: dg_elements.h:168
unsigned nface() const
Return the number of faces.
Definition: dg_elements.h:360
virtual void set_mass_matrix_from_element(DGElement *const &element_pt)
Set the mass matrix to point to one in another element.
Definition: dg_elements.h:256
virtual unsigned required_nflux()
Set the number of flux components.
Definition: dg_elements.h:193
void set_mesh_pt(DGMesh *&mesh_pt)
Definition: dg_elements.h:354
bool Mass_matrix_has_been_computed
Boolean flag to indicate whether the mass matrix has been computed.
Definition: dg_elements.h:186
double & average_value(const unsigned &i)
Return the average values.
Definition: dg_elements.h:436
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...
Definition: dg_elements.h:291
bool Mass_matrix_reuse_is_enabled
Boolean flag to indicate whether to reuse the mass matrix.
Definition: dg_elements.h:182
void output_faces(std::ostream &outfile)
Output the faces of the element.
Definition: dg_elements.h:372
DGFaceElement * face_element_pt(const unsigned &i)
Access function for the faces.
Definition: dg_elements.h:366
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.
Definition: dg_elements.cc:675
virtual ~DGElement()
Virtual destructor, destroy the mass matrix, if we created it Clean-up storage associated with averag...
Definition: dg_elements.h:212
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...
Definition: dg_elements.cc:614
const double & average_value(const unsigned &i) const
Return the average values.
Definition: dg_elements.h:449
DenseDoubleMatrix * M_pt
Pointer to storage for a mass matrix that can be recycled if desired.
Definition: dg_elements.h:175
Base class for Discontinuous Galerkin Faces. These are responsible for calculating the normal fluxes ...
Definition: dg_elements.h:50
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...
Definition: dg_elements.h:121
virtual ~DGFaceElement()
Empty Destructor.
Definition: dg_elements.h:84
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...
Definition: dg_elements.cc:241
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....
Definition: dg_elements.cc:42
void report_info()
Output information about the present element and its neighbour.
Definition: dg_elements.cc:123
Vector< FaceElement * > Neighbour_face_pt
Vector of neighbouring face elements at the integration points.
Definition: dg_elements.h:52
void add_flux_contributions(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Add the contribution from integrating the numerical flux.
Definition: dg_elements.cc:416
DGFaceElement()
Empty Constructor.
Definition: dg_elements.h:81
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...
Definition: dg_elements.h:63
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...
Definition: dg_elements.cc:338
virtual void interpolated_u(const Vector< double > &s, Vector< double > &f)
Return the interpolated values of the unknown fluxes.
Definition: dg_elements.cc:196
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...
Definition: dg_elements.cc:266
Vector< Vector< double > > Neighbour_local_coordinate
Vector of neighbouring local coordinates at the integration points.
Definition: dg_elements.h:55
virtual unsigned flux_index(const unsigned &i) const
Return the index at which the i-th unknown flux is stored.
Definition: dg_elements.h:68
FaceElement * neighbour_face_pt(const unsigned &i)
Access function for neighbouring face information.
Definition: dg_elements.h:87
virtual unsigned required_nflux()
Set the number of flux components.
Definition: dg_elements.h:74
void limit_slopes(SlopeLimiter *const &slope_limiter_pt)
Definition: dg_elements.h:502
virtual ~DGMesh()
Definition: dg_elements.h:468
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)
Definition: dg_elements.h:470
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...
Definition: dg_elements.h:489
static double FaceTolerance
Definition: dg_elements.h:464
Class of matrices containing doubles, and stored as a DenseMatrix<double>, but with solving functiona...
Definition: matrices.h:1271
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
Definition: elements.h:4338
A general Finite Element class.
Definition: elements.h:1313
virtual Node * construct_node(const unsigned &n)
Construct the local node n and return a pointer to the newly created node object.
Definition: elements.h:2509
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
Definition: elements.h:3050
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
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...
Definition: elements.h:2538
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
Definition: geom_objects.h:192
A general mesh class.
Definition: mesh.h:67
const Vector< GeneralisedElement * > & element_pt() const
Return reference to the Vector of elements.
Definition: mesh.h:460
unsigned long nelement() const
Return number of elements in the mesh.
Definition: mesh.h:590
double minmod(Vector< double > &args)
The basic minmod function.
Definition: dg_elements.cc:740
bool MUSCL
Boolean flag to indicate a MUSCL or straight min-mod limiter.
Definition: dg_elements.h:550
double minmodB(Vector< double > &args, const double &h)
The modified minmod function.
Definition: dg_elements.cc:800
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 ...
Definition: dg_elements.h:556
virtual ~MinModLimiter()
Empty destructor.
Definition: dg_elements.h:562
void limit(const unsigned &i, const Vector< DGElement * > &required_element_pt)
The limit function.
Definition: dg_elements.cc:822
double M
Constant that is used in the modified min-mod function to provide better properties at extrema.
Definition: dg_elements.h:547
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...
Definition: shape.h:76
Base class for slope limiters.
Definition: dg_elements.h:525
virtual ~SlopeLimiter()
virtual destructor
Definition: dg_elements.h:531
SlopeLimiter()
Empty constructor.
Definition: dg_elements.h:528
virtual void limit(const unsigned &i, const Vector< DGElement * > &required_element_pt)
Basic function.
Definition: dg_elements.h:534
////////////////////////////////////////////////////////////////////// //////////////////////////////...
Definition: timesteppers.h:231
A slight extension to the standard template vector class so that we can include "graceful" array rang...
Definition: Vector.h:58
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...