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