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// Version 1.0; svn revision $LastChangedRevision$
7 //LIC//
8 //LIC// $LastChangedDate$
9 //LIC//
10 //LIC// Copyright (C) 2006-2016 Matthias Heil and Andrew Hazel
11 //LIC//
12 //LIC// This library is free software; you can redistribute it and/or
13 //LIC// modify it under the terms of the GNU Lesser General Public
14 //LIC// License as published by the Free Software Foundation; either
15 //LIC// version 2.1 of the License, or (at your option) any later version.
16 //LIC//
17 //LIC// This library is distributed in the hope that it will be useful,
18 //LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
19 //LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 //LIC// Lesser General Public License for more details.
21 //LIC//
22 //LIC// You should have received a copy of the GNU Lesser General Public
23 //LIC// License along with this library; if not, write to the Free Software
24 //LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
25 //LIC// 02110-1301 USA.
26 //LIC//
27 //LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
28 //LIC//
29 //LIC//====================================================================
30 //Header file for elements that solve the advection diffusion equation
31 //and that can be refined.
32 
33 #ifndef OOMPH_REFINEABLE_GEN_ADVECTION_DIFFUSION_ELEMENTS_HEADER
34 #define OOMPH_REFINEABLE_GEN_ADVECTION_DIFFUSION_ELEMENTS_HEADER
35 
36 // Config header generated by autoconfig
37 #ifdef HAVE_CONFIG_H
38  #include <oomph-lib-config.h>
39 #endif
40 
41 //oomph-lib headers
42 #include "../generic/refineable_quad_element.h"
43 #include "../generic/refineable_brick_element.h"
44 #include "../generic/error_estimator.h"
46 
47 namespace oomph
48 {
49 
50 //======================================================================
51 /// \short A version of the GeneralisedAdvection Diffusion equations that can be
52 /// used with non-uniform mesh refinement. In essence, the class overloads
53 /// the fill_in_generic_residual_contribution_cons_adv_diff()
54 /// function so that contributions
55 /// from hanging nodes (or alternatively in-compatible function values)
56 /// are taken into account.
57 //======================================================================
58 template <unsigned DIM>
60  public virtual GeneralisedAdvectionDiffusionEquations<DIM>,
61  public virtual RefineableElement,
62  public virtual ElementWithZ2ErrorEstimator
63 {
64  public:
65 
66  /// \short Empty Constructor
71  {}
72 
73 
74  /// Broken copy constructor
77  {
79  "RefineableGeneralisedAdvectionDiffusionEquations");
80  }
81 
82  /// Broken assignment operator
84  {
86  "RefineableGeneralisedAdvectionDiffusionEquations");
87  }
88 
89  /// Number of 'flux' terms for Z2 error estimation
90  unsigned num_Z2_flux_terms() {return DIM;}
91 
92  /// \short Get 'flux' for Z2 error recovery:
93  /// Standard flux.from GeneralisedAdvectionDiffusion equations
95  {this->get_flux(s,flux);}
96 
97 
98  /// \short Get the function value u in Vector.
99  /// Note: Given the generality of the interface (this function
100  /// is usually called from black-box documentation or interpolation routines),
101  /// the values Vector sets its own size in here.
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_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  /// \short 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 routines),
132  /// the values Vector sets its own size in here.
133  void get_interpolated_values(const unsigned& t, const Vector<double>&s,
134  Vector<double>& values)
135  {
136  // Set size of Vector:
137  values.resize(1);
138 
139  //Find out how many nodes there are
140  const unsigned n_node = nnode();
141 
142  //Find the nodal index at which the unknown is stored
143  const unsigned u_nodal_index = this->u_index_cons_adv_diff();
144 
145  // Shape functions
146  Shape psi(n_node);
147 
148  //Find values of shape function
149  shape(s,psi);
150 
151  //Initialise the value of u
152  values[0] = 0.0;
153 
154  //Calculate value
155  for(unsigned l=0;l<n_node;l++)
156  {
157  values[0] += this->nodal_value(t,l,u_nodal_index)*psi[l];
158  }
159  }
160 
161 
162  /// Further build: Copy source function pointer from father element
164  {
166  cast_father_element_pt
168  this->father_element_pt());
169 
170  //Set the values of the pointers from the father
171  this->Source_fct_pt = cast_father_element_pt->source_fct_pt();
172  this->Wind_fct_pt = cast_father_element_pt->wind_fct_pt();
173  this->Conserved_wind_fct_pt = cast_father_element_pt->
175  this->Diff_fct_pt = cast_father_element_pt->diff_fct_pt();
176  this->Pe_pt = cast_father_element_pt->pe_pt();
177  this->PeSt_pt = cast_father_element_pt->pe_st_pt();
178 
179  //Set the ALE status
180  this->ALE_is_disabled = cast_father_element_pt->ALE_is_disabled;
181  }
182 
183  protected:
184 
185 /// \short Add the element's contribution to the elemental residual vector
186 /// and/or Jacobian matrix
187 /// flag=1: compute both
188 /// flag=0: compute only residual vector
190  Vector<double> &residuals, DenseMatrix<double> &jacobian,
191  DenseMatrix<double> &mass_matrix, unsigned flag);
192 
193 };
194 
195 
196 //======================================================================
197 /// \short Refineable version of QGeneralisedAdvectionDiffusionElement.
198 /// Inherit from the standard QGeneralisedAdvectionDiffusionElement and the
199 /// appropriate refineable geometric element and the refineable equations.
200 //======================================================================
201 template <unsigned DIM, unsigned NNODE_1D>
203 public QGeneralisedAdvectionDiffusionElement<DIM,NNODE_1D>,
205  public virtual RefineableQElement<DIM>
206 {
207  public:
208 
209  /// \short Empty Constructor:
213  RefineableQElement<DIM>(),
215  {}
216 
217 
218  /// Broken copy constructor
221  dummy)
222  {
224  "RefineableQuadGeneralisedAdvectionDiffusionElement");
225  }
226 
227  /// Broken assignment operator
229  {
231  "RefineableQuadGeneralisedAdvectionDiffusionElement");
232  }
233 
234  /// Number of continuously interpolated values: 1
235  unsigned ncont_interpolated_values() const {return 1;}
236 
237  /// \short Number of vertex nodes in the element
238  unsigned nvertex_node() const
239  {return
241 
242  /// \short Pointer to the j-th vertex node in the element
243  Node* vertex_node_pt(const unsigned& j) const
244  {return
246 
247  /// Rebuild from sons: empty
248  void rebuild_from_sons(Mesh* &mesh_pt) {}
249 
250  /// \short Order of recovery shape functions for Z2 error estimation:
251  /// Same order as shape functions.
252  unsigned nrecovery_order() {return (NNODE_1D-1);}
253 
254  /// \short 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 
266 //=======================================================================
267 /// Face geometry for the
268 /// RefineableQuadGeneralisedAdvectionDiffusionElement elements: The spatial
269 /// dimension of the face elements is one lower than that of the
270 /// bulk element but they have the same number of points
271 /// along their 1D edges.
272 //=======================================================================
273 template<unsigned DIM, unsigned NNODE_1D>
275  public virtual QElement<DIM-1,NNODE_1D>
276 {
277 
278  public:
279 
280  /// \short Constructor: Call the constructor for the
281  /// appropriate lower-dimensional QElement
282  FaceGeometry() : QElement<DIM-1,NNODE_1D>() {}
283 
284 };
285 
286 }
287 
288 #endif
289 
Refineable version of QGeneralisedAdvectionDiffusionElement. Inherit from the standard QGeneralisedAd...
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
unsigned num_Z2_flux_terms()
Number of &#39;flux&#39; terms for Z2 error estimation.
void further_build()
Further build: Copy source function pointer from father element.
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
char t
Definition: cfortran.h:572
GeneralisedAdvectionDiffusionDiffFctPt & diff_fct_pt()
Access function: Pointer to diffusion function.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:852
unsigned nvertex_node() const
Number of vertex nodes in the element.
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.
void fill_in_generic_residual_contribution_cons_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Add the element&#39;s contribution to the elemental residual vector and/or Jacobian matrix flag=1: comput...
A class for all elements that solve the Advection Diffusion equations in conservative form using isop...
A version of the GeneralisedAdvection Diffusion equations that can be used with non-uniform mesh refi...
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...
GeneralisedAdvectionDiffusionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
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 operator=(const RefineableGeneralisedAdvectionDiffusionEquations< DIM > &)
Broken assignment operator.
RefineableGeneralisedAdvectionDiffusionEquations(const RefineableGeneralisedAdvectionDiffusionEquations< DIM > &dummy)
Broken copy constructor.
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed...
virtual unsigned nvertex_node() const
Definition: elements.h:2374
GeneralisedAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
static char t char * s
Definition: cfortran.h:572
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
GeneralisedAdvectionDiffusionWindFctPt Conserved_wind_fct_pt
Pointer to additional (conservative) wind function:
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:2470
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes...
void operator=(const RefineableQGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D > &)
Broken assignment operator.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
QGeneralisedAdvectionDiffusionElement elements are linear/quadrilateral/brick-shaped Advection Diffus...
void broken_assign(const std::string &class_name)
Issue error message and terminate execution.
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
Pointer to diffusivity funciton.
RefineableQGeneralisedAdvectionDiffusionElement(const RefineableQGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D > &dummy)
Broken copy constructor.
GeneralisedAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
virtual Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element. Broken virtual function in "pure" finite elements...
Definition: elements.h:2383
double *& pe_st_pt()
Pointer to Peclet number multipled by Strouha number.
GeneralisedAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get &#39;flux&#39; for Z2 error recovery: Standard flux.from GeneralisedAdvectionDiffusion equations...
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2146
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement. ...
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...
A general mesh class.
Definition: mesh.h:74
GeneralisedAdvectionDiffusionWindFctPt & conserved_wind_fct_pt()
Access function: Pointer to additional (conservative) wind function.