refineable_gen_axisym_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-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 file for elements that solve the advection diffusion equation
27 // and that can be refined.
28 
29 #ifndef OOMPH_REFINEABLE_GEN_AXISYM_ADVECTION_DIFFUSION_ELEMENTS_HEADER
30 #define OOMPH_REFINEABLE_GEN_AXISYM_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 GeneralisedAxisymAdvectionDiffusion
47  /// equations that can be
48  /// used with non-uniform mesh refinement. In essence, the class overloads
49  /// the fill_in_generic_residual_contribution_cons_axisym_adv_diff()
50  /// function so that contributions
51  /// from hanging nodes (or alternatively in-compatible function values)
52  /// are taken into account.
53  //======================================================================
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  // Commented out broken assignment operator because this can lead to a
76  // conflict warning when used in the virtual inheritence hierarchy.
77  // Essentially the compiler doesn't realise that two separate
78  // implementations of the broken function are the same and so, quite
79  // rightly, it shouts.
80  /*void operator=(const
81  RefineableGeneralisedAxisymAdvectionDiffusionEquations&) = delete;*/
82 
83  /// Number of 'flux' terms for Z2 error estimation
84  unsigned num_Z2_flux_terms()
85  {
86  return 2;
87  }
88 
89  /// Get 'flux' for Z2 error recovery:
90  /// Standard flux.from GeneralisedAxisymAdvectionDiffusion equations
92  {
93  this->get_flux(s, flux);
94  }
95 
96 
97  /// Get the function value u in Vector.
98  /// Note: Given the generality of the interface (this function
99  /// is usually called from black-box documentation or interpolation
100  /// routines), the values Vector sets its own size in here.
102  Vector<double>& values)
103  {
104  // Set size of Vector: u
105  values.resize(1);
106 
107  // Find number of nodes
108  const unsigned n_node = nnode();
109 
110  // Find the index at which the unknown is stored
111  const unsigned u_nodal_index = this->u_index_cons_axisym_adv_diff();
112 
113  // Local shape function
114  Shape psi(n_node);
115 
116  // Find values of shape function
117  shape(s, psi);
118 
119  // Initialise value of u
120  values[0] = 0.0;
121 
122  // Loop over the local nodes and sum
123  for (unsigned l = 0; l < n_node; l++)
124  {
125  values[0] += this->nodal_value(l, u_nodal_index) * psi[l];
126  }
127  }
128 
129  /// Get the function value u in Vector.
130  /// Note: Given the generality of the interface (this function
131  /// is usually called from black-box documentation or interpolation
132  /// routines), the values Vector sets its own size in here.
133  void get_interpolated_values(const unsigned& t,
134  const Vector<double>& s,
135  Vector<double>& values)
136  {
137  // Set size of Vector:
138  values.resize(1);
139 
140  // Find out how many nodes there are
141  const unsigned n_node = nnode();
142 
143  // Find the nodal index at which the unknown is stored
144  const unsigned u_nodal_index = this->u_index_cons_axisym_adv_diff();
145 
146  // Shape functions
147  Shape psi(n_node);
148 
149  // Find values of shape function
150  shape(s, psi);
151 
152  // Initialise the value of u
153  values[0] = 0.0;
154 
155  // Calculate value
156  for (unsigned l = 0; l < n_node; l++)
157  {
158  values[0] += this->nodal_value(t, l, u_nodal_index) * psi[l];
159  }
160  }
161 
162  /// Fill in the geometric Jacobian, which in this case is r
164  {
165  return x[0];
166  }
167 
168 
169  /// Further build: Copy source function pointer from father element
171  {
173  cast_father_element_pt =
175  this->father_element_pt());
176 
177  // Set the values of the pointers from the father
178  this->Source_fct_pt = cast_father_element_pt->source_fct_pt();
179  this->Wind_fct_pt = cast_father_element_pt->wind_fct_pt();
180  this->Conserved_wind_fct_pt =
181  cast_father_element_pt->conserved_wind_fct_pt();
182  this->Diff_fct_pt = cast_father_element_pt->diff_fct_pt();
183  this->Pe_pt = cast_father_element_pt->pe_pt();
184  this->PeSt_pt = cast_father_element_pt->pe_st_pt();
185 
186  // Set the ALE status
187  this->ALE_is_disabled = cast_father_element_pt->ALE_is_disabled;
188  }
189 
190  protected:
191  /// Add the element's contribution to the elemental residual vector
192  /// and/or Jacobian matrix
193  /// flag=1: compute both
194  /// flag=0: compute only residual vector
196  Vector<double>& residuals,
197  DenseMatrix<double>& jacobian,
198  DenseMatrix<double>& mass_matrix,
199  unsigned flag);
200  };
201 
202 
203  //======================================================================
204  /// Refineable version of QGeneralisedAxisymAdvectionDiffusionElement.
205  /// Inherit from the standard QGeneralisedAxisymAdvectionDiffusionElement
206  /// and the
207  /// appropriate refineable geometric element and the refineable equations.
208  //======================================================================
209  template<unsigned NNODE_1D>
213  public virtual RefineableQElement<2>
214  {
215  public:
216  /// Empty Constructor:
218  : RefineableElement(),
220  RefineableQElement<2>(),
222  {
223  }
224 
225 
226  /// Broken copy constructor
229  dummy) = delete;
230 
231  /// Broken assignment operator
232  /*void operator=(const
233  RefineableQGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D>&) =
234  delete;*/
235 
236  /// Number of continuously interpolated values: 1
237  unsigned ncont_interpolated_values() const
238  {
239  return 1;
240  }
241 
242  /// Number of vertex nodes in the element
243  unsigned nvertex_node() const
244  {
246  NNODE_1D>::nvertex_node();
247  }
248 
249  /// Pointer to the j-th vertex node in the element
250  Node* vertex_node_pt(const unsigned& j) const
251  {
253  NNODE_1D>::vertex_node_pt(j);
254  }
255 
256  /// Rebuild from sons: empty
257  void rebuild_from_sons(Mesh*& mesh_pt) {}
258 
259  /// Order of recovery shape functions for Z2 error estimation:
260  /// Same order as shape functions.
261  unsigned nrecovery_order()
262  {
263  return (NNODE_1D - 1);
264  }
265 
266  /// Perform additional hanging node procedures for variables
267  /// that are not interpolated by all nodes. Empty.
269  };
270 
271  /// /////////////////////////////////////////////////////////////////////
272  /// /////////////////////////////////////////////////////////////////////
273  /// /////////////////////////////////////////////////////////////////////
274 
275 
276  //=======================================================================
277  /// Face geometry for the
278  /// RefineableQuadGeneralisedAxisymAdvectionDiffusionElement elements:
279  /// The spatial
280  /// dimension of the face elements is one lower than that of the
281  /// bulk element but they have the same number of points
282  /// along their 1D edges.
283  //=======================================================================
284  template<unsigned NNODE_1D>
287  : public virtual QElement<1, NNODE_1D>
288  {
289  public:
290  /// Constructor: Call the constructor for the
291  /// appropriate lower-dimensional QElement
292  FaceGeometry() : QElement<1, NNODE_1D>() {}
293  };
294 
295 } // namespace oomph
296 
297 #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:4998
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:2593
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:2210
GeneralisedAxisymAdvectionDiffusionDiffFctPt Diff_fct_pt
Pointer to diffusivity funciton.
double *& pe_st_pt()
Pointer to Peclet number multipled by Strouha number.
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
GeneralisedAxisymAdvectionDiffusionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
GeneralisedAxisymAdvectionDiffusionWindFctPt & conserved_wind_fct_pt()
Access function: Pointer to additional (conservative) wind function.
GeneralisedAxisymAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
GeneralisedAxisymAdvectionDiffusionDiffFctPt & diff_fct_pt()
Access function: Pointer to diffusion function.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
GeneralisedAxisymAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
GeneralisedAxisymAdvectionDiffusionWindFctPt Conserved_wind_fct_pt
Pointer to additional (conservative) wind function:
GeneralisedAxisymAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source 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 GeneralisedAxisymAdvectionDiffusion equations that can be used with non-uniform mesh...
RefineableGeneralisedAxisymAdvectionDiffusionEquations(const RefineableGeneralisedAxisymAdvectionDiffusionEquations &dummy)=delete
Broken copy constructor.
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...
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...
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from GeneralisedAxisymAdvectionDiffusion equations.
double geometric_jacobian(const Vector< double > &x)
Fill in the geometric Jacobian, which in this case is r.
void further_build()
Further build: Copy source function pointer from father element.
void fill_in_generic_residual_contribution_cons_axisym_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...
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 QGeneralisedAxisymAdvectionDiffusionElement. Inherit from the standard QGeneral...
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
RefineableQGeneralisedAxisymAdvectionDiffusionElement(const RefineableQGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
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
//////////////////////////////////////////////////////////////////// ////////////////////////////////...