refineable_pml_helmholtz_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 QPMLHelmholtzElement elements
27 
28 #ifndef OOMPH_REFINEABLE_PML_HELMHOLTZ_ELEMENTS_HEADER
29 #define OOMPH_REFINEABLE_PML_HELMHOLTZ_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/hp_refineable_elements.h"
41 #include "../generic/error_estimator.h"
42 #include "pml_helmholtz_elements.h"
43 
44 namespace oomph
45 {
46  /// ////////////////////////////////////////////////////////////////////////
47  /// ////////////////////////////////////////////////////////////////////////
48 
49 
50  //======================================================================
51  /// Refineable version of PMLHelmholtz equations
52  ///
53  ///
54  //======================================================================
55  template<unsigned DIM>
57  : public virtual PMLHelmholtzEquations<DIM>,
58  public virtual RefineableElement,
59  public virtual ElementWithZ2ErrorEstimator
60  {
61  public:
62  /// Constructor, simply call other constructors
64  : PMLHelmholtzEquations<DIM>(),
67  {
68  }
69 
70  /// Broken copy constructor
72  const RefineablePMLHelmholtzEquations<DIM>& dummy) = delete;
73 
74  /// Broken assignment operator
75  // Commented out broken assignment operator because this can lead to a
76  // conflict warning when used in the virtual inheritence hierarchy.
77  // Essentially the compiler doesn't realise that two separate
78  // implementations of the broken function are the same and so, quite
79  // rightly, it shouts.
80  /*void operator=(const RefineablePMLHelmholtzEquations<DIM>&) = delete;*/
81 
82  /// Number of 'flux' terms for Z2 error estimation
83  unsigned num_Z2_flux_terms()
84  {
85  return 2 * DIM;
86  }
87 
88  /// Get 'flux' for Z2 error recovery: Complex flux from
89  /// PMLHelmholtz equations, strung together
91  {
92  Vector<std::complex<double>> complex_flux(DIM);
93  this->get_flux(s, complex_flux);
94  unsigned count = 0;
95  for (unsigned i = 0; i < DIM; i++)
96  {
97  flux[count] = complex_flux[i].real();
98  count++;
99  flux[count] = complex_flux[i].imag();
100  count++;
101  }
102  }
103 
104 
105  /// Get the function value u in Vector.
106  /// Note: Given the generality of the interface (this function
107  /// is usually called from black-box documentation or interpolation
108  /// routines), the values Vector sets its own size in here.
110  Vector<double>& values)
111  {
112  // Set size of Vector: u
113  values.resize(2);
114 
115  // Find number of nodes
116  unsigned n_node = nnode();
117 
118  // Local shape function
119  Shape psi(n_node);
120 
121  // Find values of shape function
122  shape(s, psi);
123 
124  // Initialise value of u
125  values[0] = 0.0;
126  values[1] = 0.0;
127 
128  // Find the index at which the pml_helmholtz unknown is stored
129  std::complex<unsigned> u_nodal_index = this->u_index_helmholtz();
130 
131  // Loop over the local nodes and sum up the values
132  for (unsigned l = 0; l < n_node; l++)
133  {
134  values[0] += this->nodal_value(l, u_nodal_index.real()) * psi[l];
135  values[1] += this->nodal_value(l, u_nodal_index.imag()) * psi[l];
136  }
137  }
138 
139 
140  /// Get the function value u in Vector.
141  /// Note: Given the generality of the interface (this function
142  /// is usually called from black-box documentation or interpolation
143  /// routines), the values Vector sets its own size in here.
144  void get_interpolated_values(const unsigned& t,
145  const Vector<double>& s,
146  Vector<double>& values)
147  {
148  if (t != 0)
149  {
150  std::string error_message =
151  "Time-dependent version of get_interpolated_values() ";
152  error_message += "not implemented for this element \n";
153  throw OomphLibError(
154  error_message,
155  "RefineablePMLHelmholtzEquations::get_interpolated_values()",
156  OOMPH_EXCEPTION_LOCATION);
157  }
158  else
159  {
160  // Make sure that we call this particular object's steady
161  // get_interpolated_values (it could get overloaded lower down)
163  values);
164  }
165  }
166 
167 
168  /// Further build: Copy source function pointer from father element
170  {
172  this->father_element_pt())
173  ->source_fct_pt();
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,
184  DenseMatrix<double>& jacobian,
185  const unsigned& flag);
186  };
187 
188 
189  //======================================================================
190  /// Refineable version of QPMLHelmholtzElement elements
191  ///
192  ///
193  //======================================================================
194  template<unsigned DIM, unsigned NNODE_1D>
196  : public QPMLHelmholtzElement<DIM, NNODE_1D>,
197  public virtual RefineablePMLHelmholtzEquations<DIM>,
198  public virtual RefineableQElement<DIM>
199  {
200  public:
201  /// Constructor, simply call the other constructors
203  : RefineableElement(),
205  RefineableQElement<DIM>(),
206  QPMLHelmholtzElement<DIM, NNODE_1D>()
207  {
208  }
209 
210 
211  /// Broken copy constructor
213  const RefineableQPMLHelmholtzElement<DIM, NNODE_1D>& dummy) = delete;
214 
215  /// Broken assignment operator
216  /*void operator=(const RefineableQPMLHelmholtzElement<DIM,NNODE_1D>&) =
217  * delete;*/
218 
219  /// Number of continuously interpolated values: 2
220  unsigned ncont_interpolated_values() const
221  {
222  return 2;
223  }
224 
225  /// Number of vertex nodes in the element
226  unsigned nvertex_node() const
227  {
229  }
230 
231  /// Pointer to the j-th vertex node in the element
232  Node* vertex_node_pt(const unsigned& j) const
233  {
235  }
236 
237  /// Rebuild from sons: empty
238  void rebuild_from_sons(Mesh*& mesh_pt) {}
239 
240  /// Order of recovery shape functions for Z2 error estimation:
241  /// Same order as shape functions.
242  unsigned nrecovery_order()
243  {
244  return (NNODE_1D - 1);
245  }
246 
247  /// Perform additional hanging node procedures for variables
248  /// that are not interpolated by all nodes. Empty.
250  };
251 
252 
253  /// /////////////////////////////////////////////////////////////////////
254  /// /////////////////////////////////////////////////////////////////////
255  /// /////////////////////////////////////////////////////////////////////
256 
257 
258  //=======================================================================
259  /// Face geometry for the RefineableQuadPMLHelmholtzElement elements:
260  /// The spatial
261  /// dimension of the face elements is one lower than that of the
262  /// bulk element but they have the same number of points
263  /// along their 1D edges.
264  //=======================================================================
265  template<unsigned DIM, unsigned NNODE_1D>
267  : public virtual QElement<DIM - 1, NNODE_1D>
268  {
269  public:
270  /// Constructor: Call the constructor for the
271  /// appropriate lower-dimensional QElement
272  FaceGeometry() : QElement<DIM - 1, NNODE_1D>() {}
273  };
274 
275  //=======================================================================
276  /// Policy class defining the elements to be used in the actual
277  /// PML layers. Same!
278  //=======================================================================
279  template<unsigned NNODE_1D>
281  : public virtual RefineableQPMLHelmholtzElement<2, NNODE_1D>
282  {
283  public:
284  /// Constructor: Call the constructor for the
285  /// appropriate QElement
287  };
288 
289 } // namespace oomph
290 
291 #endif
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
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: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 unsigned nvertex_node() const
Return the number of vertex nodes in this element. Broken virtual function in "pure" finite elements.
Definition: elements.h:2495
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
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:2504
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
An OomphLibError object which should be thrown when an run-time error is encountered....
A class for all isoparametric elements that solve the Helmholtz equations with pml capabilities....
PMLHelmholtzSourceFctPt Source_fct_pt
Pointer to source function:
void get_flux(const Vector< double > &s, Vector< std::complex< double >> &flux) const
Get flux: flux[i] = du/dx_i for real and imag part.
PMLHelmholtzSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
virtual std::complex< unsigned > u_index_helmholtz() const
Broken assignment operator.
PMLLayerElement()
Constructor: Call the constructor for the appropriate QElement.
General definition of policy class defining the elements to be used in the actual PML layers....
Definition: pml_meshes.h:48
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
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.
////////////////////////////////////////////////////////////////////////
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 further_build()
Further build: Copy source function pointer from father element.
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...
unsigned num_Z2_flux_terms()
Broken assignment operator.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Complex flux from PMLHelmholtz equations, strung together.
void fill_in_generic_residual_contribution_helmholtz(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Add element's contribution to elemental residual vector and/or Jacobian matrix flag=1: compute both f...
RefineablePMLHelmholtzEquations()
Constructor, simply call other constructors.
RefineablePMLHelmholtzEquations(const RefineablePMLHelmholtzEquations< DIM > &dummy)=delete
Broken copy constructor.
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 QPMLHelmholtzElement elements.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
unsigned ncont_interpolated_values() const
Broken assignment operator.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
unsigned nvertex_node() const
Number of vertex nodes in the element.
RefineableQPMLHelmholtzElement()
Constructor, simply call the other constructors.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
RefineableQPMLHelmholtzElement(const RefineableQPMLHelmholtzElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...