refineable_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-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 refineable unsteady heat elements
27 
28 #ifndef OOMPH_REFINEABLE_UNSTEADY_HEAT_ELEMENTS_HEADER
29 #define OOMPH_REFINEABLE_UNSTEADY_HEAT_ELEMENTS_HEADER
30 
31 // Config header generated by autoconfig
32 #ifdef HAVE_CONFIG_H
33 #include <oomph-lib-config.h>
34 #endif
35 
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"
41 #include "unsteady_heat_elements.h"
42 
43 
44 namespace oomph
45 {
46  /// ////////////////////////////////////////////////////////////////////////
47  /// ////////////////////////////////////////////////////////////////////////
48 
49 
50  //======================================================================
51  /// Refineable version of Unsteady HEat equations
52  ///
53  ///
54  //======================================================================
55  template<unsigned DIM>
57  : public virtual UnsteadyHeatEquations<DIM>,
58  public virtual RefineableElement,
59  public virtual ElementWithZ2ErrorEstimator
60  {
61  public:
62  /// Constructor
64  : UnsteadyHeatEquations<DIM>(),
67  {
68  }
69 
70 
71  /// Broken copy constructor
73  const RefineableUnsteadyHeatEquations<DIM>& dummy) = delete;
74 
75  /// Broken assignment operator
76  // Commented out broken assignment operator because this can lead to a
77  // conflict warning when used in the virtual inheritence hierarchy.
78  // Essentially the compiler doesn't realise that two separate
79  // implementations of the broken function are the same and so, quite
80  // rightly, it shouts.
81  /*void operator=(const RefineableUnsteadyHeatEquations<DIM>&) = delete;*/
82 
83  /// Number of 'flux' terms for Z2 error estimation
84  unsigned num_Z2_flux_terms()
85  {
86  return DIM;
87  }
88 
89  /// Get 'flux' for Z2 error recovery:
90  /// Standard flux.from UnsteadyHeat equations
92  {
93  this->get_flux(s, flux);
94  }
95 
96  /// Get the function value u in Vector.
97  /// Note: Given the generality of the interface (this function
98  /// is usually called from black-box documentation or interpolation
99  /// routines), the values Vector sets its own size in here.
101  Vector<double>& values)
102  {
103  // Set size of Vector: u
104  values.resize(1);
105 
106  // Find number of nodes
107  unsigned n_node = nnode();
108 
109  // Find the nodal index at which the unknown is stored
110  unsigned u_nodal_index = this->u_index_ust_heat();
111 
112  // Local shape function
113  Shape psi(n_node);
114 
115  // Find values of shape function
116  shape(s, psi);
117 
118  // Initialise value of u
119  values[0] = 0.0;
120 
121  // Loop over the local nodes and sum
122  for (unsigned l = 0; l < n_node; l++)
123  {
124  values[0] += this->nodal_value(l, u_nodal_index) * psi[l];
125  }
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  // Initialise
141  values[0] = 0.0;
142 
143  // Find out how many nodes there are
144  unsigned n_node = nnode();
145 
146  // Find the nodal index at which the unknown is stored
147  unsigned u_nodal_index = this->u_index_ust_heat();
148 
149  // Shape functions
150  Shape psi(n_node);
151  shape(s, psi);
152 
153  // Calculate value
154  for (unsigned l = 0; l < n_node; l++)
155  {
156  values[0] += this->nodal_value(t, l, u_nodal_index) * psi[l];
157  }
158  }
159 
160 
161  /// Further build: Copy source function pointer from father element
163  {
164  // Get pointer to the father
165  RefineableUnsteadyHeatEquations<DIM>* cast_father_element_pt =
167  this->father_element_pt());
168 
169  // Set the source function from the parent
170  this->Source_fct_pt = cast_father_element_pt->source_fct_pt();
171 
172  // Set the ALE status from the father
173  this->ALE_is_disabled = cast_father_element_pt->ALE_is_disabled;
174  }
175 
176 
177  private:
178  /// Add element's contribution to elemental residual vector and/or
179  /// Jacobian matrix
180  /// flag=1: compute both
181  /// flag=0: compute only residual vector
183  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag);
184  };
185 
186 
187  //======================================================================
188  /// Refineable version of 2D QUnsteadyHeatElement elements
189  ///
190  ///
191  //======================================================================
192  template<unsigned DIM, unsigned NNODE_1D>
194  : public QUnsteadyHeatElement<DIM, NNODE_1D>,
195  public virtual RefineableUnsteadyHeatEquations<DIM>,
196  public virtual RefineableQElement<DIM>
197  {
198  public:
199  /// Constructor
201  : RefineableElement(),
203  RefineableQElement<DIM>(),
204  QUnsteadyHeatElement<DIM, NNODE_1D>()
205  {
206  }
207 
208 
209  /// Broken copy constructor
211  const RefineableQUnsteadyHeatElement<DIM, NNODE_1D>& dummy) = delete;
212 
213  /// Broken assignment operator
214  /*void operator=(const RefineableQUnsteadyHeatElement<DIM,NNODE_1D>&) =
215  * delete;*/
216 
217  /// Number of continuously interpolated values: 1
218  unsigned ncont_interpolated_values() const
219  {
220  return 1;
221  }
222 
223  /// Number of vertex nodes in the element
224  unsigned nvertex_node() const
225  {
227  }
228 
229  /// Pointer to the j-th vertex node in the element
230  Node* vertex_node_pt(const unsigned& j) const
231  {
233  }
234 
235  /// Rebuild from sons: empty
236  void rebuild_from_sons(Mesh*& mesh_pt) {}
237 
238  /// Order of recovery shape functions for Z2 error estimation:
239  /// Same order as shape functions.
240  unsigned nrecovery_order()
241  {
242  return (NNODE_1D - 1);
243  }
244 
245  /// Perform additional hanging node procedures for variables
246  /// that are not interpolated by all nodes. Empty.
248  };
249  /// /////////////////////////////////////////////////////////////////////
250  /// /////////////////////////////////////////////////////////////////////
251  /// /////////////////////////////////////////////////////////////////////
252 
253 
254  //=======================================================================
255  /// Face geometry for the RefineableQuadUnsteadyHeatElement elements:
256  /// The spatial
257  /// dimension of the face elements is one lower than that of the
258  /// bulk element but they have the same number of points
259  /// along their 1D edges.
260  //=======================================================================
261  template<unsigned DIM, unsigned NNODE_1D>
263  : public virtual QElement<DIM - 1, NNODE_1D>
264  {
265  public:
266  /// Constructor: Call the constructor for the
267  /// appropriate lower-dimensional QElement
268  FaceGeometry() : QElement<DIM - 1, NNODE_1D>() {}
269  };
270 
271 } // namespace oomph
272 
273 #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 unsigned nvertex_node() const
Return the number of vertex nodes in this element. Broken virtual function in "pure" finite elements.
Definition: elements.h:2491
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
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:2500
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 QUnsteadyHeatElement elements.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
unsigned ncont_interpolated_values() const
Broken assignment operator.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
unsigned nvertex_node() const
Number of vertex nodes in the element.
RefineableQUnsteadyHeatElement(const RefineableQUnsteadyHeatElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
////////////////////////////////////////////////////////////////////////
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 fill_in_generic_residual_contribution_ust_heat(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Add element's contribution to elemental residual vector and/or Jacobian matrix flag=1: compute both f...
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 UnsteadyHeat equations.
unsigned num_Z2_flux_terms()
Broken assignment operator.
void further_build()
Further build: Copy source function pointer from father element.
RefineableUnsteadyHeatEquations(const RefineableUnsteadyHeatEquations< DIM > &dummy)=delete
Broken copy constructor.
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 UnsteadyHeat equations.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
UnsteadyHeatSourceFctPt Source_fct_pt
Pointer to source function:
UnsteadyHeatSourceFctPt & 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....
virtual unsigned u_index_ust_heat() const
Broken assignment operator.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...