refineable_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-2025 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 QHelmholtzElement elements
27
28#ifndef OOMPH_REFINEABLE_HELMHOLTZ_ELEMENTS_HEADER
29#define OOMPH_REFINEABLE_HELMHOLTZ_ELEMENTS_HEADER
30
31// Config header
32#ifdef HAVE_CONFIG_H
33#include <oomph-lib-config.h>
34#endif
35
36
37// oomph-lib headers
41#include "helmholtz_elements.h"
42#include "math.h"
43#include <complex>
44
45namespace oomph
46{
47 ///////////////////////////////////////////////////////////////////////////
48 ///////////////////////////////////////////////////////////////////////////
49
50
51 //======================================================================
52 /// Refineable version of Helmholtz equations
53 ///
54 ///
55 //======================================================================
56 template<unsigned DIM>
58 : public virtual HelmholtzEquations<DIM>,
59 public virtual RefineableElement,
60 public virtual ElementWithZ2ErrorEstimator
61 {
62 public:
63 /// Constructor, simply call other constructors
70
71 /// Broken copy constructor
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 RefineableHelmholtzEquations<DIM>&) = delete;*/
82
83 /// Number of 'flux' terms for Z2 error estimation
85 {
86 return 2 * DIM;
87 }
88
89 /// Get 'flux' for Z2 error recovery: Standard flux.from Helmholtz
90 /// equations
92 {
94 this->get_flux(s, actual_flux);
95 unsigned count = 0;
96 for (unsigned i = 0; i < DIM; i++)
97 {
98 flux[count] = actual_flux[i].real();
99 flux[count + 1] = actual_flux[i].imag();
100 count += 2;
101 }
102 }
103
104 /// Get the function value u in Vector.
105 /// Note: Given the generality of the interface (this function
106 /// is usually called from black-box documentation or interpolation
107 /// routines), the values Vector sets its own size in here.
109 Vector<double>& values)
110 {
111 // Resize nitialise
112 values.resize(2, 0.0);
113
114 // Find number of nodes
115 unsigned n_node = nnode();
116
117 // Local shape function
119
120 // Find values of shape function
121 shape(s, psi);
122
123 // Loop over the local nodes and sum up the values
124 for (unsigned l = 0; l < n_node; l++)
125 {
126 values[0] +=
127 this->nodal_value(l, this->u_index_helmholtz().real()) * psi[l];
128 values[1] +=
129 this->nodal_value(l, this->u_index_helmholtz().imag()) * psi[l];
130 }
131 }
132
133
134 /// 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
137 /// routines), the values Vector sets its own size in here.
138 void get_interpolated_values(const unsigned& t,
139 const Vector<double>& s,
140 Vector<double>& values)
141 {
142 if (t != 0)
143 {
144 std::string error_message =
145 "Time-dependent version of get_interpolated_values() ";
146 error_message += "not implemented for this element \n";
147 throw OomphLibError(
148 error_message,
149 "RefineableHelmholtzEquations::get_interpolated_values()",
151 }
152 else
153 {
154 // Make sure that we call this particular object's steady
155 // get_interpolated_values (it could get overloaded lower down)
157 }
158 }
159
160
161 /// Further build: Copy source function pointer from father element
163 {
165 this->father_element_pt())
166 ->source_fct_pt();
167
168 this->K_squared_pt = dynamic_cast<RefineableHelmholtzEquations<DIM>*>(
169 this->father_element_pt())
170 ->k_squared_pt();
171 }
172
173
174 private:
175 /// Add element's contribution to elemental residual vector and/or
176 /// Jacobian matrix
177 /// flag=1: compute both
178 /// flag=0: compute only residual vector
181 DenseMatrix<double>& jacobian,
182 const unsigned& flag);
183 };
184
185
186 //======================================================================
187 /// Refineable version of 2D QHelmholtzElement elements
188 ///
189 ///
190 //======================================================================
191 template<unsigned DIM, unsigned NNODE_1D>
193 : public QHelmholtzElement<DIM, NNODE_1D>,
194 public virtual RefineableHelmholtzEquations<DIM>,
195 public virtual RefineableQElement<DIM>
196 {
197 public:
198 /// Constructor, simply call the other constructors
206
207
208 /// Broken copy constructor
211
212 /// Broken assignment operator
213 /*void operator=(const RefineableQHelmholtzElement<DIM,NNODE_1D>&) =
214 * delete;*/
215
216 /// Number of continuously interpolated values: 2
218 {
219 return 2;
220 }
221
222 /// Number of vertex nodes in the element
223 unsigned nvertex_node() const
224 {
226 }
227
228 /// Pointer to the j-th vertex node in the element
229 Node* vertex_node_pt(const unsigned& j) const
230 {
232 }
233
234 /// Rebuild from sons: empty
235 void rebuild_from_sons(Mesh*& mesh_pt) {}
236
237 /// Order of recovery shape functions for Z2 error estimation:
238 /// Same order as shape functions.
240 {
241 return (NNODE_1D - 1);
242 }
243
244 /// Perform additional hanging node procedures for variables
245 /// that are not interpolated by all nodes. Empty.
247 };
248
249 ////////////////////////////////////////////////////////////////////////
250 ////////////////////////////////////////////////////////////////////////
251 ////////////////////////////////////////////////////////////////////////
252
253
254 //=======================================================================
255 /// Face geometry for the RefineableQuadHelmholtzElement elements: The spatial
256 /// dimension of the face elements is one lower than that of the
257 /// bulk element but they have the same number of points
258 /// along their 1D edges.
259 //=======================================================================
260 template<unsigned DIM, unsigned NNODE_1D>
262 : public virtual QElement<DIM - 1, NNODE_1D>
263 {
264 public:
265 /// Constructor: Call the constructor for the
266 /// appropriate lower-dimensional QElement
268 };
269
270} // namespace oomph
271
272#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.
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
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 class for all isoparametric elements that solve the Helmholtz equations.
double * K_squared_pt
Pointer to square of wavenumber.
double *& k_squared_pt()
Get pointer to square of wavenumber.
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.
HelmholtzSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
HelmholtzSourceFctPt Source_fct_pt
Pointer to source function:
virtual std::complex< unsigned > u_index_helmholtz() const
Broken assignment operator.
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....
General QElement class.
Definition Qelements.h:459
QHelmholtzElement elements are linear/quadrilateral/brick-shaped Helmholtz elements with isoparametri...
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 Helmholtz equations.
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...
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.
RefineableHelmholtzEquations()
Constructor, simply call other constructors.
RefineableHelmholtzEquations(const RefineableHelmholtzEquations< 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...
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from Helmholtz 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 QHelmholtzElement elements.
unsigned ncont_interpolated_values() const
Broken assignment operator.
RefineableQHelmholtzElement()
Constructor, simply call the other constructors.
RefineableQHelmholtzElement(const RefineableQHelmholtzElement< 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.
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 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.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition shape.h:76
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).