refineable_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// 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_AXISYM_ADVECTION_DIFFUSION_ELEMENTS_HEADER
34 #define OOMPH_REFINEABLE_AXISYM_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 Advection Diffusion in axisym
52 /// coordinates equations that can be
53 /// used with non-uniform mesh refinement. In essence, the class overloads
54 /// the fill_in_generic_residual_contribution_axisym_adv_diff()
55 /// function so that contributions
56 /// from hanging nodes (or alternatively in-compatible function values)
57 /// are taken into account.
58 //======================================================================
60  public virtual AxisymAdvectionDiffusionEquations,
61  public virtual RefineableElement,
62  public virtual ElementWithZ2ErrorEstimator
63 {
64  public:
65 
66  /// \short Empty Constructor
71  {}
72 
73 
74  /// Broken copy constructor
77  {
78  BrokenCopy::broken_copy("RefineableAxisymAdvectionDiffusionEquations");
79  }
80 
81  /// Broken assignment operator
82 //Commented out broken assignment operator because this can lead to a conflict warning
83 //when used in the virtual inheritence hierarchy. Essentially the compiler doesn't
84 //realise that two separate implementations of the broken function are the same and so,
85 //quite rightly, it shouts.
86  /*void operator=(const RefineableAxisymAdvectionDiffusionEquations&)
87  {
88  BrokenCopy::broken_assign("RefineableAxisymAdvectionDiffusionEquations");
89  }*/
90 
91  /// Number of 'flux' terms for Z2 error estimation
92  unsigned num_Z2_flux_terms() {return 2;}
93 
94  /// \short Get 'flux' for Z2 error recovery:
95  /// Standard flux.from AdvectionDiffusion equations
97  {this->get_flux(s,flux);}
98 
99 
100  /// \short Get the function value u in Vector.
101  /// Note: Given the generality of the interface (this function
102  /// is usually called from black-box documentation or interpolation routines),
103  /// the values Vector sets its own size in here.
105  {
106  // Set size of Vector: u
107  values.resize(1);
108 
109  //Find number of nodes
110  unsigned n_node = nnode();
111 
112  //Find the index at which the unknown is stored
113  unsigned u_nodal_index = this->u_index_axi_adv_diff();
114 
115  //Local shape function
116  Shape psi(n_node);
117 
118  //Find values of shape function
119  shape(s,psi);
120 
121  //Initialise value of u
122  values[0] = 0.0;
123 
124  //Loop over the local nodes and sum
125  for(unsigned l=0;l<n_node;l++)
126  {
127  values[0] += this->nodal_value(l,u_nodal_index)*psi[l];
128  }
129  }
130 
131  /// \short Get the function value u in Vector.
132  /// Note: Given the generality of the interface (this function
133  /// is usually called from black-box documentation or interpolation routines),
134  /// the values Vector sets its own size in here.
135  void get_interpolated_values(const unsigned& t, const Vector<double>&s,
136  Vector<double>& values)
137  {
138  // Set size of Vector: u
139  values.resize(1);
140 
141  //Find number of nodes
142  const unsigned n_node = nnode();
143 
144  //Find the index at which the unknown is stored
145  const unsigned u_nodal_index = this->u_index_axi_adv_diff();
146 
147  //Local shape function
148  Shape psi(n_node);
149 
150  //Find values of shape function
151  shape(s,psi);
152 
153  //Initialise value of u
154  values[0] = 0.0;
155 
156  //Loop over the local nodes and sum
157  for(unsigned l=0;l<n_node;l++)
158  {
159  values[0] += this->nodal_value(t,l,u_nodal_index)*psi[l];
160  }
161  // }
162 
163  if (t!=0)
164  {
165  std::string error_message =
166  "Time-dependent version of get_interpolated_values() ";
167  error_message += "not implemented for this element \n";
168  throw
170  error_message,
171  "RefineableAxisymAdvectionDiffusionEquations::get_interpolated_values()",
172  OOMPH_EXCEPTION_LOCATION);
173  }
174  else
175  {
176  //Make sure that we call the appropriate steady version
177  //(the entire function might be overloaded lower down)
179  get_interpolated_values(s,values);
180  }
181  }
182 
183  /// Fill in the geometric Jacobian, which in this case is r
185  {return x[0];}
186 
187  /// Further build: Copy source function pointer from father element
189  {
190  RefineableAxisymAdvectionDiffusionEquations* cast_father_element_pt
192  this->father_element_pt());
193 
194  //Set the values of the pointers from the father
195  this->Source_fct_pt = cast_father_element_pt->source_fct_pt();
196  this->Wind_fct_pt = cast_father_element_pt->wind_fct_pt();
197  this->Pe_pt = cast_father_element_pt->pe_pt();
198  this->PeSt_pt = cast_father_element_pt->pe_st_pt();
199 
200  //Set the ALE status
201  //this->ALE_is_disabled = cast_father_element_pt->ALE_is_disabled;
202  }
203 
204  /// \short Compute the derivatives of the i-th component of
205  /// velocity at point s with respect
206  /// to all data that can affect its value. In addition, return the global
207  /// equation numbers corresponding to the data.
208  /// Overload the non-refineable version to take account of hanging node
209  /// information
211  Vector<double> &du_ddata,
212  Vector<unsigned> &global_eqn_number)
213  {
214  //Find number of nodes
215  unsigned n_node = this->nnode();
216  //Local shape function
217  Shape psi(n_node);
218  //Find values of shape function at the given local coordinate
219  this->shape(s,psi);
220 
221  //Find the index at which the velocity component is stored
222  const unsigned u_nodal_index = this->u_index_axi_adv_diff();
223 
224  //Storage for hang info pointer
225  HangInfo* hang_info_pt=0;
226  //Storage for global equation
227  int global_eqn = 0;
228 
229  //Find the number of dofs associated with interpolated u
230  unsigned n_u_dof=0;
231  for(unsigned l=0;l<n_node;l++)
232  {
233  unsigned n_master = 1;
234 
235  //Local bool (is the node hanging)
236  bool is_node_hanging = this->node_pt(l)->is_hanging();
237 
238  //If the node is hanging, get the number of master nodes
239  if(is_node_hanging)
240  {
241  hang_info_pt = this->node_pt(l)->hanging_pt();
242  n_master = hang_info_pt->nmaster();
243  }
244  //Otherwise there is just one master node, the node itself
245  else
246  {
247  n_master = 1;
248  }
249 
250  //Loop over the master nodes
251  for(unsigned m=0;m<n_master;m++)
252  {
253  //Get the equation number
254  if(is_node_hanging)
255  {
256  //Get the equation number from the master node
257  global_eqn = hang_info_pt->master_node_pt(m)->
258  eqn_number(u_nodal_index);
259  }
260  else
261  {
262  // Global equation number
263  global_eqn = this->node_pt(l)->eqn_number(u_nodal_index);
264  }
265 
266  //If it's positive add to the count
267  if (global_eqn >= 0) {++n_u_dof;}
268  }
269  }
270 
271  //Now resize the storage schemes
272  du_ddata.resize(n_u_dof,0.0);
273  global_eqn_number.resize(n_u_dof,0);
274 
275  //Loop over th nodes again and set the derivatives
276  unsigned count=0;
277  //Loop over the local nodes and sum
278  for(unsigned l=0;l<n_node;l++)
279  {
280  unsigned n_master = 1;
281  double hang_weight = 1.0;
282 
283  //Local bool (is the node hanging)
284  bool is_node_hanging = this->node_pt(l)->is_hanging();
285 
286  //If the node is hanging, get the number of master nodes
287  if(is_node_hanging)
288  {
289  hang_info_pt = this->node_pt(l)->hanging_pt();
290  n_master = hang_info_pt->nmaster();
291  }
292  //Otherwise there is just one master node, the node itself
293  else
294  {
295  n_master = 1;
296  }
297 
298  //Loop over the master nodes
299  for(unsigned m=0;m<n_master;m++)
300  {
301  //If the node is hanging get weight from master node
302  if(is_node_hanging)
303  {
304  //Get the hang weight from the master node
305  hang_weight = hang_info_pt->master_weight(m);
306  }
307  else
308  {
309  // Node contributes with full weight
310  hang_weight = 1.0;
311  }
312 
313  //Get the equation number
314  if(is_node_hanging)
315  {
316  //Get the equation number from the master node
317  global_eqn = hang_info_pt->master_node_pt(m)->
318  eqn_number(u_nodal_index);
319  }
320  else
321  {
322  // Global equation number
323  global_eqn = this->node_pt(l)->eqn_number(u_nodal_index);
324  }
325 
326  if (global_eqn >= 0)
327  {
328  //Set the global equation number
329  global_eqn_number[count] = global_eqn;
330  //Set the derivative with respect to the unknown
331  du_ddata[count] = psi[l]*hang_weight;
332  //Increase the counter
333  ++count;
334  }
335  }
336  }
337  }
338 
339 
340  protected:
341 
342 /// \short Add the element's contribution to the elemental residual vector
343 /// and/or Jacobian matrix
344 /// flag=1: compute both
345 /// flag=0: compute only residual vector
347  Vector<double> &residuals, DenseMatrix<double> &jacobian,
348  DenseMatrix<double> &mass_matrix, unsigned flag);
349 
350 };
351 
352 
353 //======================================================================
354 /// \short Refineable version of QAxisymAdvectionDiffusionElement.
355 /// Inherit from the standard QAxisymAdvectionDiffusionElement and the
356 /// appropriate refineable geometric element and the refineable equations.
357 //======================================================================
358 template <unsigned NNODE_1D>
360 public QAxisymAdvectionDiffusionElement<NNODE_1D>,
362  public virtual RefineableQElement<2>
363 {
364  public:
365 
366  /// \short Empty Constructor:
370  RefineableQElement<2>(),
372  {}
373 
374 
375  /// Broken copy constructor
378  dummy)
379  {
380  BrokenCopy::broken_copy("RefineableQAxisymAdvectionDiffusionElement");
381  }
382 
383  /// Broken assignment operator
384  /*void operator=(const
385  RefineableQAxisymAdvectionDiffusionElement<NNODE_1D>&)
386  {
387  BrokenCopy::broken_assign("RefineableQAxisymAdvectionDiffusionElement");
388  }*/
389 
390  /// Number of continuously interpolated values: 1
391  unsigned ncont_interpolated_values() const {return 1;}
392 
393  /// \short Number of vertex nodes in the element
394  unsigned nvertex_node() const
396 
397  /// \short Pointer to the j-th vertex node in the element
398  Node* vertex_node_pt(const unsigned& j) const
400 
401  /// Rebuild from sons: empty
402  void rebuild_from_sons(Mesh* &mesh_pt) {}
403 
404  /// \short Order of recovery shape functions for Z2 error estimation:
405  /// Same order as shape functions.
406  unsigned nrecovery_order() {return (NNODE_1D-1);}
407 
408  /// \short Perform additional hanging node procedures for variables
409  /// that are not interpolated by all nodes. Empty.
411 
412 };
413 
414 ////////////////////////////////////////////////////////////////////////
415 ////////////////////////////////////////////////////////////////////////
416 ////////////////////////////////////////////////////////////////////////
417 
418 
419 
420 //=======================================================================
421 /// Face geometry for the RefineableQAxisymAdvectionDiffusionElement
422 /// elements: The spatial
423 /// dimension of the face elements is one lower than that of the
424 /// bulk element but they have the same number of points
425 /// along their 1D edges.
426 //=======================================================================
427 template<unsigned NNODE_1D>
429  public virtual QElement<1,NNODE_1D>
430 {
431 
432  public:
433 
434  /// \short Constructor: Call the constructor for the
435  /// appropriate lower-dimensional QElement
436  FaceGeometry() : QElement<1,NNODE_1D>() {}
437 
438 };
439 
440 }
441 
442 #endif
443 
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: Qelements.h:910
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes...
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.
RefineableQAxisymAdvectionDiffusionElement(const RefineableQAxisymAdvectionDiffusionElement< NNODE_1D > &dummy)
Broken copy constructor.
HangInfo *const & hanging_pt() const
Return pointer to hanging node data (this refers to the geometric hanging node status) (const version...
Definition: nodes.h:1148
RefineableAxisymAdvectionDiffusionEquations(const RefineableAxisymAdvectionDiffusionEquations &dummy)
Broken copy constructor.
double *& pe_st_pt()
Pointer to Peclet number multipled by Strouha number.
char t
Definition: cfortran.h:572
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: Qelements.h:906
double const & master_weight(const unsigned &i) const
Return weight for dofs on i-th master node.
Definition: nodes.h:753
unsigned nvertex_node() const
Number of vertex nodes in the element.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:852
unsigned nmaster() const
Return the number of master nodes.
Definition: nodes.h:733
bool is_hanging() const
Test whether the node is geometrically hanging.
Definition: nodes.h:1207
long & eqn_number(const unsigned &i)
Return the equation number of the i-th stored variable.
Definition: nodes.h:365
double * Pe_pt
Pointer to global Peclet number.
AxisymAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
AxisymAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
QAxisymAdvectionDiffusionElement elements are linear/quadrilateral/brick-shaped Axisymmetric Advectio...
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement. ...
static char t char * s
Definition: cfortran.h:572
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux:
Refineable version of QAxisymAdvectionDiffusionElement. Inherit from the standard QAxisymAdvectionDif...
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number...
Definition: elements.h:709
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_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get &#39;flux&#39; for Z2 error recovery: Standard flux.from AdvectionDiffusion equations.
void further_build()
Further build: Copy source function pointer from father element.
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
Definition: nodes.h:736
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2109
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 dinterpolated_u_adv_diff_ddata(const Vector< double > &s, Vector< double > &du_ddata, Vector< unsigned > &global_eqn_number)
Compute the derivatives of the i-th component of velocity at point s with respect to all data that ca...
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
Class that contains data for hanging nodes.
Definition: nodes.h:684
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
AxisymAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn&#39;t been defined.
AxisymAdvectionDiffusionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
double geometric_jacobian(const Vector< double > &x)
Fill in the geometric Jacobian, which in this case is r.
void fill_in_generic_residual_contribution_axi_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...
virtual unsigned u_index_axi_adv_diff() const
Broken assignment operator.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
A version of the Advection Diffusion in axisym coordinates equations that can be used with non-unifor...
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2146
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
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
A class for all elements that solve the Advection Diffusion equations in a cylindrical polar coordina...