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-2022 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 {
48 
49 
50  //======================================================================
54  //======================================================================
55  template<unsigned DIM>
57  : public virtual UnsteadyHeatEquations<DIM>,
58  public virtual RefineableElement,
59  public virtual ElementWithZ2ErrorEstimator
60  {
61  public:
67  {
68  }
69 
70 
73  const RefineableUnsteadyHeatEquations<DIM>& dummy) = delete;
74 
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 
84  unsigned num_Z2_flux_terms()
85  {
86  return DIM;
87  }
88 
92  {
93  this->get_flux(s, flux);
94  }
95 
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 
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 
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:
183  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag);
184  };
185 
186 
187  //======================================================================
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:
201  : RefineableElement(),
204  QUnsteadyHeatElement<DIM, NNODE_1D>()
205  {
206  }
207 
208 
211  const RefineableQUnsteadyHeatElement<DIM, NNODE_1D>& dummy) = delete;
212 
214  /*void operator=(const RefineableQUnsteadyHeatElement<DIM,NNODE_1D>&) =
215  * delete;*/
216 
218  unsigned ncont_interpolated_values() const
219  {
220  return 1;
221  }
222 
224  unsigned nvertex_node() const
225  {
227  }
228 
230  Node* vertex_node_pt(const unsigned& j) const
231  {
233  }
234 
236  void rebuild_from_sons(Mesh*& mesh_pt) {}
237 
240  unsigned nrecovery_order()
241  {
242  return (NNODE_1D - 1);
243  }
244 
248  };
252 
253 
254  //=======================================================================
260  //=======================================================================
261  template<unsigned DIM, unsigned NNODE_1D>
263  : public virtual QElement<DIM - 1, NNODE_1D>
264  {
265  public:
268  FaceGeometry() : QElement<DIM - 1, NNODE_1D>() {}
269  };
270 
271 } // namespace oomph
272 
273 #endif
Definition: error_estimator.h:79
FaceGeometry()
Definition: refineable_unsteady_heat_elements.h:268
Definition: elements.h:4998
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual unsigned nvertex_node() const
Definition: elements.h:2491
virtual void shape(const Vector< double > &s, Shape &psi) const =0
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
virtual Node * vertex_node_pt(const unsigned &j) const
Definition: elements.h:2500
Definition: mesh.h:67
Definition: nodes.h:906
Definition: Qelements.h:459
Definition: unsteady_heat_elements.h:482
Definition: refineable_elements.h:97
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
Definition: refineable_elements.h:539
Definition: Qelements.h:2259
Definition: refineable_unsteady_heat_elements.h:197
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_unsteady_heat_elements.h:230
unsigned ncont_interpolated_values() const
Broken assignment operator.
Definition: refineable_unsteady_heat_elements.h:218
void further_setup_hanging_nodes()
Definition: refineable_unsteady_heat_elements.h:247
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_unsteady_heat_elements.h:224
RefineableQUnsteadyHeatElement(const RefineableQUnsteadyHeatElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Definition: refineable_unsteady_heat_elements.h:236
RefineableQUnsteadyHeatElement()
Constructor.
Definition: refineable_unsteady_heat_elements.h:200
unsigned nrecovery_order()
Definition: refineable_unsteady_heat_elements.h:240
Definition: refineable_unsteady_heat_elements.h:60
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_unsteady_heat_elements.h:100
void fill_in_generic_residual_contribution_ust_heat(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: refineable_unsteady_heat_elements.cc:39
RefineableUnsteadyHeatEquations()
Constructor.
Definition: refineable_unsteady_heat_elements.h:63
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_unsteady_heat_elements.h:133
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: refineable_unsteady_heat_elements.h:91
unsigned num_Z2_flux_terms()
Broken assignment operator.
Definition: refineable_unsteady_heat_elements.h:84
void further_build()
Further build: Copy source function pointer from father element.
Definition: refineable_unsteady_heat_elements.h:162
RefineableUnsteadyHeatEquations(const RefineableUnsteadyHeatEquations< DIM > &dummy)=delete
Broken copy constructor.
Definition: shape.h:76
Definition: unsteady_heat_elements.h:72
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
Definition: unsteady_heat_elements.h:280
UnsteadyHeatSourceFctPt Source_fct_pt
Pointer to source function:
Definition: unsteady_heat_elements.h:446
UnsteadyHeatSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: unsteady_heat_elements.h:222
bool ALE_is_disabled
Definition: unsteady_heat_elements.h:451
virtual unsigned u_index_ust_heat() const
Broken assignment operator.
Definition: unsteady_heat_elements.h:112
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
t
Definition: plotPSD.py:36
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2