refineable_gen_advection_diffusion_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-2024 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 file for elements that solve the advection diffusion equation
27 // and that can be refined.
28 
29 #ifndef OOMPH_REFINEABLE_GEN_ADVECTION_DIFFUSION_ELEMENTS_HEADER
30 #define OOMPH_REFINEABLE_GEN_ADVECTION_DIFFUSION_ELEMENTS_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 headers
38 #include "../generic/refineable_quad_element.h"
39 #include "../generic/refineable_brick_element.h"
40 #include "../generic/error_estimator.h"
42 
43 namespace oomph
44 {
45  //======================================================================
46  /// A version of the GeneralisedAdvection Diffusion equations that can
47  /// be used with non-uniform mesh refinement. In essence, the class overloads
48  /// the fill_in_generic_residual_contribution_cons_adv_diff()
49  /// function so that contributions
50  /// from hanging nodes (or alternatively in-compatible function values)
51  /// are taken into account.
52  //======================================================================
53  template<unsigned DIM>
55  : public virtual GeneralisedAdvectionDiffusionEquations<DIM>,
56  public virtual RefineableElement,
57  public virtual ElementWithZ2ErrorEstimator
58  {
59  public:
60  /// Empty Constructor
65  {
66  }
67 
68 
69  /// Broken copy constructor
72  delete;
73 
74  /// Broken assignment operator
75  void operator=(
77 
78  /// Number of 'flux' terms for Z2 error estimation
79  unsigned num_Z2_flux_terms()
80  {
81  return DIM;
82  }
83 
84  /// Get 'flux' for Z2 error recovery:
85  /// Standard flux.from GeneralisedAdvectionDiffusion equations
87  {
88  this->get_flux(s, flux);
89  }
90 
91 
92  /// Get the function value u in Vector.
93  /// Note: Given the generality of the interface (this function
94  /// is usually called from black-box documentation or interpolation
95  /// routines), the values Vector sets its own size in here.
97  Vector<double>& values)
98  {
99  // Set size of Vector: u
100  values.resize(1);
101 
102  // Find number of nodes
103  const unsigned n_node = nnode();
104 
105  // Find the index at which the unknown is stored
106  const unsigned u_nodal_index = this->u_index_cons_adv_diff();
107 
108  // Local shape function
109  Shape psi(n_node);
110 
111  // Find values of shape function
112  shape(s, psi);
113 
114  // Initialise value of u
115  values[0] = 0.0;
116 
117  // Loop over the local nodes and sum
118  for (unsigned l = 0; l < n_node; l++)
119  {
120  values[0] += this->nodal_value(l, u_nodal_index) * psi[l];
121  }
122  }
123 
124  /// Get the function value u in Vector.
125  /// Note: Given the generality of the interface (this function
126  /// is usually called from black-box documentation or interpolation
127  /// routines), the values Vector sets its own size in here.
128  void get_interpolated_values(const unsigned& t,
129  const Vector<double>& s,
130  Vector<double>& values)
131  {
132  // Set size of Vector:
133  values.resize(1);
134 
135  // Find out how many nodes there are
136  const unsigned n_node = nnode();
137 
138  // Find the nodal index at which the unknown is stored
139  const unsigned u_nodal_index = this->u_index_cons_adv_diff();
140 
141  // Shape functions
142  Shape psi(n_node);
143 
144  // Find values of shape function
145  shape(s, psi);
146 
147  // Initialise the value of u
148  values[0] = 0.0;
149 
150  // Calculate value
151  for (unsigned l = 0; l < n_node; l++)
152  {
153  values[0] += this->nodal_value(t, l, u_nodal_index) * psi[l];
154  }
155  }
156 
157 
158  /// Further build: Copy source function pointer from father element
160  {
162  cast_father_element_pt =
164  this->father_element_pt());
165 
166  // Set the values of the pointers from the father
167  this->Source_fct_pt = cast_father_element_pt->source_fct_pt();
168  this->Wind_fct_pt = cast_father_element_pt->wind_fct_pt();
169  this->Conserved_wind_fct_pt =
170  cast_father_element_pt->conserved_wind_fct_pt();
171  this->Diff_fct_pt = cast_father_element_pt->diff_fct_pt();
172  this->Pe_pt = cast_father_element_pt->pe_pt();
173  this->PeSt_pt = cast_father_element_pt->pe_st_pt();
174 
175  // Set the ALE status
176  this->ALE_is_disabled = cast_father_element_pt->ALE_is_disabled;
177  }
178 
179  protected:
180  /// Add the element's contribution to the elemental residual vector
181  /// and/or Jacobian matrix
182  /// flag=1: compute both
183  /// flag=0: compute only residual vector
185  Vector<double>& residuals,
186  DenseMatrix<double>& jacobian,
187  DenseMatrix<double>& mass_matrix,
188  unsigned flag);
189  };
190 
191 
192  //======================================================================
193  /// Refineable version of QGeneralisedAdvectionDiffusionElement.
194  /// Inherit from the standard QGeneralisedAdvectionDiffusionElement and the
195  /// appropriate refineable geometric element and the refineable equations.
196  //======================================================================
197  template<unsigned DIM, unsigned NNODE_1D>
199  : public QGeneralisedAdvectionDiffusionElement<DIM, NNODE_1D>,
201  public virtual RefineableQElement<DIM>
202  {
203  public:
204  /// Empty Constructor:
206  : RefineableElement(),
208  RefineableQElement<DIM>(),
210  {
211  }
212 
213 
214  /// Broken copy constructor
217  dummy) = delete;
218 
219  /// Broken assignment operator
220  void operator=(
222  delete;
223 
224  /// Number of continuously interpolated values: 1
225  unsigned ncont_interpolated_values() const
226  {
227  return 1;
228  }
229 
230  /// Number of vertex nodes in the element
231  unsigned nvertex_node() const
232  {
234  NNODE_1D>::nvertex_node();
235  }
236 
237  /// Pointer to the j-th vertex node in the element
238  Node* vertex_node_pt(const unsigned& j) const
239  {
241  NNODE_1D>::vertex_node_pt(j);
242  }
243 
244  /// Rebuild from sons: empty
245  void rebuild_from_sons(Mesh*& mesh_pt) {}
246 
247  /// Order of recovery shape functions for Z2 error estimation:
248  /// Same order as shape functions.
249  unsigned nrecovery_order()
250  {
251  return (NNODE_1D - 1);
252  }
253 
254  /// Perform additional hanging node procedures for variables
255  /// that are not interpolated by all nodes. Empty.
257  };
258 
259  /// /////////////////////////////////////////////////////////////////////
260  /// /////////////////////////////////////////////////////////////////////
261  /// /////////////////////////////////////////////////////////////////////
262 
263 
264  //=======================================================================
265  /// Face geometry for the
266  /// RefineableQuadGeneralisedAdvectionDiffusionElement elements: The spatial
267  /// dimension of the face elements is one lower than that of the
268  /// bulk element but they have the same number of points
269  /// along their 1D edges.
270  //=======================================================================
271  template<unsigned DIM, unsigned NNODE_1D>
274  : public virtual QElement<DIM - 1, NNODE_1D>
275  {
276  public:
277  /// Constructor: Call the constructor for the
278  /// appropriate lower-dimensional QElement
279  FaceGeometry() : QElement<DIM - 1, NNODE_1D>() {}
280  };
281 
282 } // namespace oomph
283 
284 #endif
static char t char * s
Definition: cfortran.h:568
char t
Definition: cfortran.h:568
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: elements.h:5002
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
Definition: elements.h:2597
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2214
A class for all elements that solve the Advection Diffusion equations in conservative form using isop...
GeneralisedAdvectionDiffusionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
GeneralisedAdvectionDiffusionDiffFctPt Diff_fct_pt
Pointer to diffusivity funciton.
double *& pe_st_pt()
Pointer to Peclet number multipled by Strouha number.
GeneralisedAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
GeneralisedAdvectionDiffusionWindFctPt Conserved_wind_fct_pt
Pointer to additional (conservative) wind function:
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
GeneralisedAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
virtual unsigned u_index_cons_adv_diff() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
GeneralisedAdvectionDiffusionDiffFctPt & diff_fct_pt()
Access function: Pointer to diffusion function.
GeneralisedAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
GeneralisedAdvectionDiffusionWindFctPt & conserved_wind_fct_pt()
Access function: Pointer to additional (conservative) wind function.
A general mesh class.
Definition: mesh.h:67
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:906
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: Qelements.h:459
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
A version of the GeneralisedAdvection Diffusion equations that can be used with non-uniform mesh refi...
void further_build()
Further build: Copy source function pointer from father element.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from GeneralisedAdvectionDiffusion equations.
void operator=(const RefineableGeneralisedAdvectionDiffusionEquations< DIM > &)=delete
Broken assignment operator.
RefineableGeneralisedAdvectionDiffusionEquations(const RefineableGeneralisedAdvectionDiffusionEquations< DIM > &dummy)=delete
Broken copy constructor.
void fill_in_generic_residual_contribution_cons_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Add the element's contribution to the elemental residual vector and/or Jacobian matrix flag=1: comput...
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Get the function value u in Vector. Note: Given the generality of the interface (this function is usu...
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Get the function value u in Vector. Note: Given the generality of the interface (this function is usu...
A class that is used to template the refineable Q elements by dimension. It's really nothing more tha...
Definition: Qelements.h:2259
Refineable version of QGeneralisedAdvectionDiffusionElement. Inherit from the standard QGeneralisedAd...
void operator=(const RefineableQGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
RefineableQGeneralisedAdvectionDiffusionElement(const RefineableQGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
unsigned nvertex_node() const
Number of vertex nodes in the element.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
//////////////////////////////////////////////////////////////////// ////////////////////////////////...