refineable_space_time_unsteady_heat_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 refineable unsteady heat elements
27 #ifndef OOMPH_REFINEABLE_SPACE_TIME_UNSTEADY_HEAT_ELEMENTS_HEADER
28 #define OOMPH_REFINEABLE_SPACE_TIME_UNSTEADY_HEAT_ELEMENTS_HEADER
29 
30 // Config header generated by autoconfig
31 #ifdef HAVE_CONFIG_H
32 #include <oomph-lib-config.h>
33 #endif
34 
35 // Oomph-lib headers
40 
41 /// //////////////////////////////////////////////////////////////////////
42 /// //////////////////////////////////////////////////////////////////////
43 /// //////////////////////////////////////////////////////////////////////
44 
45 namespace oomph
46 {
47  //======================================================================
48  /// Refineable version of Unsteady Heat equations
49  //======================================================================
50  template<unsigned SPATIAL_DIM>
52  : public virtual SpaceTimeUnsteadyHeatEquations<SPATIAL_DIM>,
53  public virtual RefineableElement,
54  public virtual ElementWithZ2ErrorEstimator
55  {
56  public:
57  /// Constructor
59  : SpaceTimeUnsteadyHeatEquations<SPATIAL_DIM>(),
62  {
63  }
64 
65 
66  /// Broken copy constructor
69  delete;
70 
71  /// Number of 'flux' terms for Z2 error estimation
72  unsigned num_Z2_flux_terms()
73  {
74  // The flux terms are associated with spatial AND temporal derivatives
75  return SPATIAL_DIM + 1;
76  } // End of num_Z2_flux_terms
77 
78 
79  /// Get 'flux' for Z2 error recovery:
80  /// Different to the get_flux function in the base class as we also
81  /// have to include du/dt as we're doing temporal adaptivity too
83  {
84  // Find out how many nodes there are in the element
85  unsigned n_node = nnode();
86 
87  // Find the index at which the variable is stored
88  unsigned u_nodal_index = this->u_index_ust_heat();
89 
90  // Set up memory for the shape and test functions
91  Shape psi(n_node);
92 
93  // Set up memory for the derivatives of the shape and test functions
94  DShape dpsidx(n_node, SPATIAL_DIM + 1);
95 
96  // Call the derivatives of the shape and test functions
97  dshape_eulerian(s, psi, dpsidx);
98 
99  // Loop over the entries of the flux vector
100  for (unsigned j = 0; j < SPATIAL_DIM + 1; j++)
101  {
102  // Initialise j-th flux entry to zero
103  flux[j] = 0.0;
104  }
105 
106  // Loop over nodes
107  for (unsigned l = 0; l < n_node; l++)
108  {
109  // Loop over derivative directions
110  for (unsigned j = 0; j < SPATIAL_DIM + 1; j++)
111  {
112  // Update the flux value
113  flux[j] += this->nodal_value(l, u_nodal_index) * dpsidx(l, j);
114  }
115  } // for (unsigned l=0;l<n_node;l++)
116  } // End of get_Z2_flux
117 
118 
119  /// Get the function value u in Vector.
120  /// Note: Given the generality of the interface (this function is usually
121  /// called from black-box documentation or interpolation routines), the
122  /// values Vector sets its own size in here.
124  Vector<double>& values)
125  {
126  // Set the size of the vector u
127  values.resize(1);
128 
129  // Find the number of nodes
130  unsigned n_node = nnode();
131 
132  // Find the nodal index at which the unknown is stored
133  unsigned u_nodal_index = this->u_index_ust_heat();
134 
135  // Local shape function
136  Shape psi(n_node);
137 
138  // Find values of shape function
139  shape(s, psi);
140 
141  // Initialise the value of u
142  values[0] = 0.0;
143 
144  // Loop over the local nodes and sum
145  for (unsigned l = 0; l < n_node; l++)
146  {
147  // Update the solution value
148  values[0] += this->nodal_value(l, u_nodal_index) * psi[l];
149  }
150  } // End of get_interpolated_values
151 
152 
153  /// Get the function value u in Vector.
154  /// Note: Given the generality of the interface (this function is usually
155  /// called from black-box documentation or interpolation routines), the
156  /// values Vector sets its own size in here.
157  void get_interpolated_values(const unsigned& t,
158  const Vector<double>& s,
159  Vector<double>& values)
160  {
161  // Set the size of the vector u
162  values.resize(1);
163 
164  // Find the number of nodes
165  unsigned n_node = nnode();
166 
167  // Find the nodal index at which the unknown is stored
168  unsigned u_nodal_index = this->u_index_ust_heat();
169 
170  // Local shape function
171  Shape psi(n_node);
172 
173  // Find values of shape function
174  shape(s, psi);
175 
176  // Initialise the value of u
177  values[0] = 0.0;
178 
179  // Loop over the local nodes and sum
180  for (unsigned l = 0; l < n_node; l++)
181  {
182  // Update the solution value
183  values[0] += this->nodal_value(t, l, u_nodal_index) * psi[l];
184  }
185  } // End of get_interpolated_values
186 
187 
188  /// Further build: Copy source function pointer from father element
190  {
191  // Get pointer to the father
193  cast_father_element_pt =
195  this->father_element_pt());
196 
197  // Get the source function from the parent and store it
198  this->Source_fct_pt = cast_father_element_pt->source_fct_pt();
199 
200  // Set the ALE status from the father
201  this->ALE_is_disabled = cast_father_element_pt->ALE_is_disabled;
202  } // End of further_build
203 
204  private:
205  /// Add element's contribution to elemental residual vector and/or
206  /// Jacobian matrix
207  /// flag=0: compute residual vector only
208  /// flag=1: compute both
210  Vector<double>& residuals,
211  DenseMatrix<double>& jacobian,
212  const unsigned& flag);
213  }; // End of RefineableSpaceTimeUnsteadyHeatEquations class
214 
215 
216  //======================================================================
217  /// Refineable version of 2D QUnsteadyHeatSpaceTimeElement elements
218  //======================================================================
219  template<unsigned SPATIAL_DIM, unsigned NNODE_1D>
221  : public QUnsteadyHeatSpaceTimeElement<SPATIAL_DIM, NNODE_1D>,
222  public virtual RefineableSpaceTimeUnsteadyHeatEquations<SPATIAL_DIM>,
223  public virtual RefineableQElement<SPATIAL_DIM + 1>
224  {
225  public:
226  /// Constructor
228  : RefineableElement(),
230  RefineableQElement<SPATIAL_DIM + 1>(),
231  QUnsteadyHeatSpaceTimeElement<SPATIAL_DIM, NNODE_1D>()
232  {
233  }
234 
235 
236  /// Broken copy constructor
239  dummy) = delete;
240 
241  /// Rebuild from sons (empty)
242  void rebuild_from_sons(Mesh*& mesh_pt) {}
243 
244 
245  /// Perform additional hanging node procedures for variables
246  /// that are not interpolated by all nodes (empty).
248 
249 
250  /// Number of continuously interpolated values: 1
251  unsigned ncont_interpolated_values() const
252  {
253  // Return the appropriate value
254  return 1;
255  } // End of ncont_interpolated_values
256 
257 
258  /// Number of vertex nodes in the element
259  unsigned nvertex_node() const
260  {
261  // Call the base class function
262  return QUnsteadyHeatSpaceTimeElement<SPATIAL_DIM,
263  NNODE_1D>::nvertex_node();
264  } // End of nvertex_node
265 
266 
267  /// Pointer to the j-th vertex node in the element
268  Node* vertex_node_pt(const unsigned& j) const
269  {
270  // Call the base class function
271  return QUnsteadyHeatSpaceTimeElement<SPATIAL_DIM,
272  NNODE_1D>::vertex_node_pt(j);
273  } // End of vertex_node_pt
274 
275 
276  /// Order of recovery shape functions for Z2 error estimation:
277  /// Same order as shape functions.
278  unsigned nrecovery_order()
279  {
280  // Return the approriate value
281  return (NNODE_1D - 1);
282  } // End of nrecovery_order
283  }; // End of RefineableQUnsteadyHeatSpaceTimeElement class
284 
285 
286  /// /////////////////////////////////////////////////////////////////////
287  /// /////////////////////////////////////////////////////////////////////
288  /// /////////////////////////////////////////////////////////////////////
289 
290 
291  //=======================================================================
292  /// Face geometry for the RefineableQuadUnsteadyHeatSpaceTimeElement elements:
293  /// The spatial
294  /// dimension of the face elements is one lower than that of the
295  /// bulk element but they have the same number of points
296  /// along their 1D edges.
297  //=======================================================================
298  template<unsigned SPATIAL_DIM, unsigned NNODE_1D>
300  RefineableQUnsteadyHeatSpaceTimeElement<SPATIAL_DIM, NNODE_1D>>
301  : public virtual QElement<SPATIAL_DIM, NNODE_1D>
302  {
303  public:
304  /// Constructor: Call the constructor for the
305  /// appropriate lower-dimensional QElement
306  FaceGeometry() : QElement<SPATIAL_DIM, NNODE_1D>() {}
307  }; // End of FaceGeometry<RefineableQUnsteadyHeatSpaceTimeElement... class
308 } // End of namespace oomph
309 
310 #endif
static char t char * s
Definition: cfortran.h:568
char t
Definition: cfortran.h:568
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition: shape.h:278
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
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
Definition: elements.cc:3328
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 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 2D QUnsteadyHeatSpaceTimeElement elements.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
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 further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes (empt...
unsigned nvertex_node() const
Number of vertex nodes in the element.
RefineableQUnsteadyHeatSpaceTimeElement(const RefineableQUnsteadyHeatSpaceTimeElement< SPATIAL_DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
virtual void fill_in_generic_residual_contribution_ust_heat(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Add element's contribution to elemental residual vector and/or Jacobian matrix flag=0: compute residu...
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...
RefineableSpaceTimeUnsteadyHeatEquations(const RefineableSpaceTimeUnsteadyHeatEquations< SPATIAL_DIM > &dummy)=delete
Broken copy constructor.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
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: Different to the get_flux function in the base class as we also hav...
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...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
A class for all isoparametric elements that solve the SpaceTimeUnsteadyHeat equations.
virtual unsigned u_index_ust_heat() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
SpaceTimeUnsteadyHeatSourceFctPt Source_fct_pt
Pointer to source function:
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
SpaceTimeUnsteadyHeatSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...