refineable_linear_wave_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 QPoissonElement elements
27 #ifndef OOMPH_REFINEABLE_LINEAR_WAVE_ELEMENTS_HEADER
28 #define OOMPH_REFINEABLE_LINEAR_WAVE_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
36 #include "../generic/refineable_quad_element.h"
37 #include "../generic/refineable_brick_element.h"
38 #include "../generic/error_estimator.h"
39 #include "linear_wave_elements.h"
40 
41 namespace oomph
42 {
43  //======================================================================
44  /// Refineable version of LinearWave equations.
45  //======================================================================
46  template<unsigned DIM>
48  : public virtual LinearWaveEquations<DIM>,
49  public virtual RefineableElement,
50  public virtual ElementWithZ2ErrorEstimator
51  {
52  public:
53  /// Constructor
55  : LinearWaveEquations<DIM>(),
58  {
59  }
60 
61 
62  /// Broken copy constructor
64  const RefineableLinearWaveEquations<DIM>& dummy) = delete;
65 
66  /// Broken assignment operator
68 
69  /// Number of 'flux' terms for Z2 error estimation
70  unsigned num_Z2_flux_terms()
71  {
72  return DIM;
73  }
74 
75  /// Get 'flux' for Z2 error recovery: Standard flux.from LinearWave
76  /// equations
78  {
79  this->get_flux(s, flux);
80  }
81 
82 
83  /// Get the function value u in Vector.
84  /// Note: Given the generality of the interface (this function
85  /// is usually called from black-box documentation or interpolation
86  /// routines), the values Vector sets its own size in here.
88  Vector<double>& values)
89  {
90  // Set size of Vector: u
91  values.resize(1);
92 
93  // Find number of nodes
94  unsigned n_node = nnode();
95 
96  // Find the index at which the unknown is stored
97  unsigned u_nodal_index = this->u_index_lin_wave();
98 
99  // Local shape function
100  Shape psi(n_node);
101 
102  // Find values of shape function
103  shape(s, psi);
104 
105  // Initialise value of u
106  values[0] = 0.0;
107 
108  // Loop over the local nodes and sum
109  for (unsigned l = 0; l < n_node; l++)
110  {
111  values[0] += this->nodal_value(l, u_nodal_index) * psi[l];
112  }
113  }
114 
115 
116  /// Get the function value u in Vector.
117  /// Note: Given the generality of the interface (this function
118  /// is usually called from black-box documentation or interpolation
119  /// routines), the values Vector sets its own size in here.
120  void get_interpolated_values(const unsigned& t,
121  const Vector<double>& s,
122  Vector<double>& values)
123  {
124  // Set size of Vector:
125  values.resize(1);
126 
127  // Initialise
128  values[0] = 0.0;
129 
130  // Find out how many nodes there are
131  unsigned n_node = nnode();
132 
133  // Find the index at which the unknown is stored
134  unsigned u_nodal_index = this->u_index_lin_wave();
135 
136  // Shape functions
137  Shape psi(n_node);
138  shape(s, psi);
139 
140  // Calculate value
141  for (unsigned l = 0; l < n_node; l++)
142  {
143  values[0] += this->nodal_value(t, l, u_nodal_index) * psi[l];
144  }
145  }
146 
147 
148  /// Further build: Copy source function pointer from father element
150  {
151  this->Source_fct_pt = dynamic_cast<RefineableLinearWaveEquations<DIM>*>(
152  this->father_element_pt())
153  ->source_fct_pt();
154  }
155 
156 
157  private:
158  /// Add element's contribution to elemental residual vector and/or
159  /// Jacobian matrix
160  /// flag=1: compute both
161  /// flag=0: compute only residual vector
163  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag);
164  };
165 
166 
167  //======================================================================
168  /// Refineable version of 2D QLinearWaveElement elements
169  ///
170  ///
171  //======================================================================
172  template<unsigned DIM, unsigned NNODE_1D>
174  : public QLinearWaveElement<DIM, NNODE_1D>,
175  public virtual RefineableLinearWaveEquations<DIM>,
176  public virtual RefineableQElement<DIM>
177  {
178  public:
179  /// Constructor
181  : RefineableElement(),
183  RefineableQElement<DIM>(),
184  QLinearWaveElement<DIM, NNODE_1D>()
185  {
186  }
187 
188 
189  /// Broken copy constructor
191  const RefineableQLinearWaveElement<DIM, NNODE_1D>& dummy) = delete;
192 
193  /// Broken assignment operator
195 
196  /// Number of continuously interpolated values: 1
197  unsigned ncont_interpolated_values() const
198  {
199  return 1;
200  }
201 
202  /// Number of vertex nodes in the element
203  unsigned nvertex_node() const
204  {
206  }
207 
208  /// Pointer to the j-th vertex node in the element
209  Node* vertex_node_pt(const unsigned& j) const
210  {
212  }
213 
214  /// Rebuild from sons: empty
215  void rebuild_from_sons(Mesh*& mesh_pt) {}
216 
217  /// Order of recovery shape functions for Z2 error estimation:
218  /// Same order as shape functions.
219  unsigned nrecovery_order()
220  {
221  return (NNODE_1D - 1);
222  }
223 
224  /// Perform additional hanging node procedures for variables
225  /// that are not interpolated by all nodes. Empty.
227  };
228  /// /////////////////////////////////////////////////////////////////////
229  /// /////////////////////////////////////////////////////////////////////
230  /// /////////////////////////////////////////////////////////////////////
231 
232 
233  //=======================================================================
234  /// Face geometry for the RefineableQuadLinearWaveElement elements: The
235  /// spatial dimension of the face elements is one lower than that of the bulk
236  /// element but they have the same number of points along their 1D edges.
237  //=======================================================================
238  template<unsigned DIM, unsigned NNODE_1D>
240  : public virtual QElement<DIM - 1, NNODE_1D>
241  {
242  public:
243  /// Constructor: Call the constructor for the
244  /// appropriate lower-dimensional QElement
245  FaceGeometry() : QElement<DIM - 1, NNODE_1D>() {}
246  };
247 
248 } // namespace oomph
249 
250 #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 class for all isoparametric elements that solve the LinearWave equations.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
virtual unsigned u_index_lin_wave() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
LinearWaveSourceFctPt Source_fct_pt
Pointer to source function:
LinearWaveSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
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.
Refineable version of LinearWave equations.
void fill_in_generic_residual_contribution_lin_wave(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 Vector< double > &s, Vector< double > &values)
Get the function value u in Vector. Note: Given the generality of the interface (this function is usu...
RefineableLinearWaveEquations(const RefineableLinearWaveEquations< DIM > &dummy)=delete
Broken copy constructor.
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...
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 operator=(const RefineableLinearWaveEquations< DIM > &)=delete
Broken assignment operator.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from LinearWave equations.
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 QLinearWaveElement elements.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
unsigned nvertex_node() const
Number of vertex nodes in the element.
void operator=(const RefineableQLinearWaveElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
RefineableQLinearWaveElement(const RefineableQLinearWaveElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
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.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
//////////////////////////////////////////////////////////////////// ////////////////////////////////...