Thelmholtz_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 THelmholtz elements
27 #ifndef OOMPH_THELMHOLTZ_ELEMENTS_HEADER
28 #define OOMPH_THELMHOLTZ_ELEMENTS_HEADER
29 
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/nodes.h"
39 #include "../generic/oomph_utilities.h"
40 #include "../generic/Telements.h"
41 #include "../generic/error_estimator.h"
42 #include "helmholtz_elements.h"
43 
44 namespace oomph
45 {
46  /// //////////////////////////////////////////////////////////////////////
47  /// //////////////////////////////////////////////////////////////////////
48  // THelmholtzElement
49  /// /////////////////////////////////////////////////////////////////////
50  /// /////////////////////////////////////////////////////////////////////
51 
52 
53  //======================================================================
54  /// THelmholtzElement<DIM,NNODE_1D> elements are isoparametric triangular
55  /// DIM-dimensional Helmholtz elements with NNODE_1D nodal points along each
56  /// element edge. Inherits from TElement and HelmholtzEquations
57  //======================================================================
58  template<unsigned DIM, unsigned NNODE_1D>
59  class THelmholtzElement : public TElement<DIM, NNODE_1D>,
60  public HelmholtzEquations<DIM>,
61  public virtual ElementWithZ2ErrorEstimator
62  {
63  public:
64  /// Constructor: Call constructors for TElement and
65  /// Helmholtz equations
66  THelmholtzElement() : TElement<DIM, NNODE_1D>(), HelmholtzEquations<DIM>()
67  {
68  }
69 
70 
71  /// Broken copy constructor
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 THelmholtzElement<DIM,NNODE_1D>&) = delete;*/
81 
82  /// Access function for Nvalue: # of `values' (pinned or dofs)
83  /// at node n (always returns the same value at every node, 1)
84  inline unsigned required_nvalue(const unsigned& n) const
85  {
86  return Initial_Nvalue;
87  }
88 
89  /// Output function:
90  /// x,y,u or x,y,z,u
91  void output(std::ostream& outfile)
92  {
94  }
95 
96  /// Output function:
97  /// x,y,u or x,y,z,u at n_plot^DIM plot points
98  void output(std::ostream& outfile, const unsigned& n_plot)
99  {
100  HelmholtzEquations<DIM>::output(outfile, n_plot);
101  }
102 
103 
104  /// C-style output function:
105  /// x,y,u or x,y,z,u
106  void output(FILE* file_pt)
107  {
109  }
110 
111 
112  /// C-style output function:
113  /// x,y,u or x,y,z,u at n_plot^DIM plot points
114  void output(FILE* file_pt, const unsigned& n_plot)
115  {
116  HelmholtzEquations<DIM>::output(file_pt, n_plot);
117  }
118 
119 
120  /// Output function for an exact solution:
121  /// x,y,u_exact
122  void output_fct(std::ostream& outfile,
123  const unsigned& n_plot,
125  {
126  HelmholtzEquations<DIM>::output_fct(outfile, n_plot, exact_soln_pt);
127  }
128 
129 
130  /// Output function for a time-dependent exact solution.
131  /// x,y,u_exact (calls the steady version)
132  void output_fct(std::ostream& outfile,
133  const unsigned& n_plot,
134  const double& time,
136  {
137  HelmholtzEquations<DIM>::output_fct(outfile, n_plot, time, exact_soln_pt);
138  }
139 
140  protected:
141  /// Shape, test functions & derivs. w.r.t. to global coords. Return
142  /// Jacobian.
144  Shape& psi,
145  DShape& dpsidx,
146  Shape& test,
147  DShape& dtestdx) const;
148 
149 
150  /// Shape, test functions & derivs. w.r.t. to global coords. Return
151  /// Jacobian.
153  const unsigned& ipt,
154  Shape& psi,
155  DShape& dpsidx,
156  Shape& test,
157  DShape& dtestdx) const;
158 
159 
160  /// Order of recovery shape functions for Z2 error estimation:
161  /// Same order as shape functions.
162  unsigned nrecovery_order()
163  {
164  return (NNODE_1D - 1);
165  }
166 
167  /// Number of 'flux' terms for Z2 error estimation
168  unsigned num_Z2_flux_terms()
169  {
170  return 2 * DIM;
171  }
172 
173  /// Get 'flux' for Z2 error recovery: Standard flux from
174  /// UnsteadyHeat equations
176  {
177  Vector<std::complex<double>> complex_flux(DIM);
178  this->get_flux(s, complex_flux);
179  unsigned count = 0;
180  for (unsigned i = 0; i < DIM; i++)
181  {
182  flux[count++] = complex_flux[i].real();
183  flux[count++] = complex_flux[i].imag();
184  }
185  }
186 
187  /// Number of vertex nodes in the element
188  unsigned nvertex_node() const
189  {
191  }
192 
193  /// Pointer to the j-th vertex node in the element
194  Node* vertex_node_pt(const unsigned& j) const
195  {
197  }
198 
199  private:
200  /// Static unsigned that holds the (same) number of variables at every node
201  static const unsigned Initial_Nvalue;
202  };
203 
204 
205  // Inline functions:
206 
207 
208  //======================================================================
209  /// Define the shape functions and test functions and derivatives
210  /// w.r.t. global coordinates and return Jacobian of mapping.
211  ///
212  /// Galerkin: Test functions = shape functions
213  //======================================================================
214  template<unsigned DIM, unsigned NNODE_1D>
216  const Vector<double>& s,
217  Shape& psi,
218  DShape& dpsidx,
219  Shape& test,
220  DShape& dtestdx) const
221  {
222  unsigned n_node = this->nnode();
223 
224  // Call the geometrical shape functions and derivatives
225  double J = this->dshape_eulerian(s, psi, dpsidx);
226 
227  // Loop over the test functions and derivatives and set them equal to the
228  // shape functions
229  for (unsigned i = 0; i < n_node; i++)
230  {
231  test[i] = psi[i];
232  dtestdx(i, 0) = dpsidx(i, 0);
233  dtestdx(i, 1) = dpsidx(i, 1);
234  }
235 
236  // Return the jacobian
237  return J;
238  }
239 
240 
241  //======================================================================
242  /// Define the shape functions and test functions and derivatives
243  /// w.r.t. global coordinates and return Jacobian of mapping.
244  ///
245  /// Galerkin: Test functions = shape functions
246  //======================================================================
247  template<unsigned DIM, unsigned NNODE_1D>
250  Shape& psi,
251  DShape& dpsidx,
252  Shape& test,
253  DShape& dtestdx) const
254  {
255  // Call the geometrical shape functions and derivatives
256  double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
257 
258  // Set the pointers of the test functions
259  test = psi;
260  dtestdx = dpsidx;
261 
262  // Return the jacobian
263  return J;
264  }
265 
266 
267  //=======================================================================
268  /// Face geometry for the THelmholtzElement elements: The spatial
269  /// dimension of the face elements is one lower than that of the
270  /// bulk element but they have the same number of points
271  /// along their 1D edges.
272  //=======================================================================
273  template<unsigned DIM, unsigned NNODE_1D>
274  class FaceGeometry<THelmholtzElement<DIM, NNODE_1D>>
275  : public virtual TElement<DIM - 1, NNODE_1D>
276  {
277  public:
278  /// Constructor: Call the constructor for the
279  /// appropriate lower-dimensional TElement
280  FaceGeometry() : TElement<DIM - 1, NNODE_1D>() {}
281  };
282 
283  //=======================================================================
284  /// Face geometry for the 1D THelmholtzElement elements: Point elements
285  //=======================================================================
286  template<unsigned NNODE_1D>
287  class FaceGeometry<THelmholtzElement<1, NNODE_1D>>
288  : public virtual PointElement
289  {
290  public:
291  /// Constructor: Call the constructor for the
292  /// appropriate lower-dimensional TElement
294  };
295 
296 
297 } // namespace oomph
298 
299 #endif
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition: shape.h:278
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 TElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: elements.h:5002
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition: elements.h:1763
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
Definition: elements.h:1769
A class for all isoparametric elements that solve the Helmholtz equations.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: x,y,u_re_exact,u_im_exact or x,y,z,u_re_exact,u_im_exact at n_plot^DIM plot points...
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.
void output(std::ostream &outfile)
Output with default number of plot points.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:906
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: elements.h:3443
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
General TElement class.
Definition: Telements.h:1208
//////////////////////////////////////////////////////////////////////
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output function for a time-dependent exact solution. x,y,u_exact (calls the steady version)
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: x,y,u or x,y,z,u at n_plot^DIM plot points.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
static const unsigned Initial_Nvalue
Static unsigned that holds the (same) number of variables at every node.
THelmholtzElement(const THelmholtzElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: x,y,u_exact.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux from UnsteadyHeat equations.
void output(FILE *file_pt)
C-style output function: x,y,u or x,y,z,u.
void output(std::ostream &outfile)
Output function: x,y,u or x,y,z,u.
double dshape_and_dtest_eulerian_helmholtz(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
double dshape_and_dtest_eulerian_at_knot_helmholtz(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: x,y,u or x,y,z,u at n_plot^DIM plot points.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
THelmholtzElement()
Constructor: Call constructors for TElement and Helmholtz equations.
unsigned nvertex_node() const
Number of vertex nodes in the element.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...