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// Version 1.0; svn revision $LastChangedRevision$
7 //LIC//
8 //LIC// $LastChangedDate$
9 //LIC//
10 //LIC// Copyright (C) 2006-2016 Matthias Heil and Andrew Hazel
11 //LIC//
12 //LIC// This library is free software; you can redistribute it and/or
13 //LIC// modify it under the terms of the GNU Lesser General Public
14 //LIC// License as published by the Free Software Foundation; either
15 //LIC// version 2.1 of the License, or (at your option) any later version.
16 //LIC//
17 //LIC// This library is distributed in the hope that it will be useful,
18 //LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
19 //LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 //LIC// Lesser General Public License for more details.
21 //LIC//
22 //LIC// You should have received a copy of the GNU Lesser General Public
23 //LIC// License along with this library; if not, write to the Free Software
24 //LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
25 //LIC// 02110-1301 USA.
26 //LIC//
27 //LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
28 //LIC//
29 //LIC//====================================================================
30 //Header file for refineable QHelmholtzElement elements
31 
32 #ifndef OOMPH_REFINEABLE_HELMHOLTZ_ELEMENTS_HEADER
33 #define OOMPH_REFINEABLE_HELMHOLTZ_ELEMENTS_HEADER
34 
35 // Config header generated by autoconfig
36 #ifdef HAVE_CONFIG_H
37  #include <oomph-lib-config.h>
38 #endif
39 
40 
41 //oomph-lib headers
42 #include "../generic/refineable_quad_element.h"
43 #include "../generic/refineable_brick_element.h"
44 #include "../generic/error_estimator.h"
45 #include "helmholtz_elements.h"
46 #include "math.h"
47 #include <complex>
48 
49 namespace oomph
50 {
51 
52 ///////////////////////////////////////////////////////////////////////////
53 ///////////////////////////////////////////////////////////////////////////
54 
55 
56 
57 //======================================================================
58 /// Refineable version of Helmholtz equations
59 ///
60 ///
61 //======================================================================
62 template <unsigned DIM>
64  public virtual RefineableElement,
65  public virtual ElementWithZ2ErrorEstimator
66  {
67  public:
68 
69  /// \short Constructor, simply call other constructors
72 
73  /// Broken copy constructor
75  {
76  BrokenCopy::broken_copy("RefineableHelmholtzEquations");
77  }
78 
79  /// Broken assignment operator
80 //Commented out broken assignment operator because this can lead to a conflict warning
81 //when used in the virtual inheritence hierarchy. Essentially the compiler doesn't
82 //realise that two separate implementations of the broken function are the same and so,
83 //quite rightly, it shouts.
84  /*void operator=(const RefineableHelmholtzEquations<DIM>&)
85  {
86  BrokenCopy::broken_assign("RefineableHelmholtzEquations");
87  }*/
88 
89  /// Number of 'flux' terms for Z2 error estimation
90  unsigned num_Z2_flux_terms() {return 2*DIM;}
91 
92  /// Get 'flux' for Z2 error recovery: Standard flux.from Helmholtz equations
94  {
95  Vector<std::complex<double> > actual_flux(DIM);
96  this->get_flux(s,actual_flux);
97  unsigned count=0;
98  for (unsigned i=0;i<DIM;i++)
99  {
100  flux[count] =actual_flux[i].real();
101  flux[count+1]=actual_flux[i].imag();
102  count+=2;
103  }
104  }
105 
106 /// \short Get the function value u in Vector.
107 /// Note: Given the generality of the interface (this function
108 /// is usually called from black-box documentation or interpolation routines),
109 /// the values Vector sets its own size in here.
111  {
112 
113  //Resize nitialise
114  values.resize(2,0.0);
115 
116  //Find number of nodes
117  unsigned n_node = nnode();
118 
119  //Local shape function
120  Shape psi(n_node);
121 
122  //Find values of shape function
123  shape(s,psi);
124 
125  //Loop over the local nodes and sum up the values
126  for(unsigned l=0;l<n_node;l++)
127  {
128  values[0]+=this->nodal_value(l,this->u_index_helmholtz().real())*psi[l];
129  values[1]+=this->nodal_value(l,this->u_index_helmholtz().imag())*psi[l];
130  }
131  }
132 
133 
134  /// \short 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 routines),
137  /// the values Vector sets its own size in here.
138  void get_interpolated_values(const unsigned& t, const Vector<double>&s,
139  Vector<double>& values)
140  {
141  if (t!=0)
142  {
143  std::string error_message =
144  "Time-dependent version of get_interpolated_values() ";
145  error_message += "not implemented for this element \n";
146  throw
147  OomphLibError(error_message,
148  "RefineableHelmholtzEquations::get_interpolated_values()",
149  OOMPH_EXCEPTION_LOCATION);
150  }
151  else
152  {
153  //Make sure that we call this particular object's steady
154  //get_interpolated_values (it could get overloaded lower down)
156  }
157  }
158 
159 
160  /// Further build: Copy source function pointer from father element
162  {
164  this->father_element_pt())->source_fct_pt();
165 
166  this->K_squared_pt=dynamic_cast<RefineableHelmholtzEquations<DIM>*>(
167  this->father_element_pt())->k_squared_pt();
168  }
169 
170 
171 
172 
173  private:
174 
175 
176  /// \short Add element's contribution to elemental residual vector and/or
177  /// Jacobian matrix
178  /// flag=1: compute both
179  /// flag=0: compute only residual vector
181  Vector<double> &residuals, DenseMatrix<double> &jacobian,
182  const unsigned& flag);
183 
184 
185  };
186 
187 
188 //======================================================================
189 /// Refineable version of 2D QHelmholtzElement elements
190 ///
191 ///
192 //======================================================================
193 template <unsigned DIM, unsigned NNODE_1D>
195  public QHelmholtzElement<DIM,NNODE_1D>,
196  public virtual RefineableHelmholtzEquations<DIM>,
197  public virtual RefineableQElement<DIM>
198  {
199  public:
200 
201  /// \short Constructor, simply call the other constructors
205  RefineableQElement<DIM>(),
206  QHelmholtzElement<DIM,NNODE_1D>()
207  {}
208 
209 
210  /// Broken copy constructor
212  dummy)
213  {
214  BrokenCopy::broken_copy("RefineableQuadHelmholtzElement");
215  }
216 
217  /// Broken assignment operator
218  /*void operator=(const RefineableQHelmholtzElement<DIM,NNODE_1D>&)
219  {
220  BrokenCopy::broken_assign("RefineableQuadHelmholtzElement");
221  }*/
222 
223  /// Number of continuously interpolated values: 2
224  unsigned ncont_interpolated_values() const {return 2;}
225 
226  /// \short Number of vertex nodes in the element
227  unsigned nvertex_node() const
229 
230  /// \short Pointer to the j-th vertex node in the element
231  Node* vertex_node_pt(const unsigned& j) const
233 
234  /// Rebuild from sons: empty
235  void rebuild_from_sons(Mesh* &mesh_pt) {}
236 
237  /// \short Order of recovery shape functions for Z2 error estimation:
238  /// Same order as shape functions.
239  unsigned nrecovery_order() {return (NNODE_1D-1);}
240 
241  /// \short Perform additional hanging node procedures for variables
242  /// that are not interpolated by all nodes. Empty.
244 
245  };
246 
247 ////////////////////////////////////////////////////////////////////////
248 ////////////////////////////////////////////////////////////////////////
249 ////////////////////////////////////////////////////////////////////////
250 
251 
252 
253 //=======================================================================
254 /// Face geometry for the RefineableQuadHelmholtzElement elements: The spatial
255 /// dimension of the face elements is one lower than that of the
256 /// bulk element but they have the same number of points
257 /// along their 1D edges.
258 //=======================================================================
259 template<unsigned DIM, unsigned NNODE_1D>
261  public virtual QElement<DIM-1,NNODE_1D>
262  {
263 
264  public:
265 
266  /// \short Constructor: Call the constructor for the
267  /// appropriate lower-dimensional QElement
268  FaceGeometry() : QElement<DIM-1,NNODE_1D>() {}
269 
270  };
271 
272 }
273 
274 #endif
275 
HelmholtzSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
RefineableQHelmholtzElement(const RefineableQHelmholtzElement< DIM, NNODE_1D > &dummy)
Broken copy constructor.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get &#39;flux&#39; for Z2 error recovery: Standard flux.from Helmholtz equations.
HelmholtzSourceFctPt Source_fct_pt
Pointer to source function:
cstr elem_len * i
Definition: cfortran.h:607
unsigned nvertex_node() const
Number of vertex nodes in the element.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes...
char t
Definition: cfortran.h:572
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.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:852
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...
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
RefineableHelmholtzEquations()
Constructor, simply call other constructors.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
virtual unsigned nvertex_node() const
Definition: elements.h:2374
static char t char * s
Definition: cfortran.h:572
double *& k_squared_pt()
Get pointer to square of wavenumber.
double * K_squared_pt
Pointer to square of wavenumber.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement. ...
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:2470
RefineableHelmholtzEquations(const RefineableHelmholtzEquations< DIM > &dummy)
Broken copy constructor.
virtual std::complex< unsigned > u_index_helmholtz() const
Broken assignment operator.
unsigned num_Z2_flux_terms()
Broken assignment operator.
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.
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn&#39;t been defined.
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:2383
void fill_in_generic_residual_contribution_helmholtz(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Add element&#39;s contribution to elemental residual vector and/or Jacobian matrix flag=1: compute both f...
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2146
RefineableQHelmholtzElement()
Constructor, simply call the other constructors.
unsigned ncont_interpolated_values() const
Broken assignment operator.
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...
A general mesh class.
Definition: mesh.h:74