Tfourier_decomposed_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 TFourierDecomposedHelmholtz elements
31 #ifndef OOMPH_TFOURIER_DECOMPOSED_HELMHOLTZ_ELEMENTS_HEADER
32 #define OOMPH_TFOURIER_DECOMPOSED_HELMHOLTZ_ELEMENTS_HEADER
33 
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/nodes.h"
43 #include "../generic/oomph_utilities.h"
44 #include "../generic/Telements.h"
45 #include "../generic/error_estimator.h"
47 
48 namespace oomph
49 {
50 
51 /////////////////////////////////////////////////////////////////////////
52 /////////////////////////////////////////////////////////////////////////
53 // TFourierDecomposedHelmholtzElement
54 ////////////////////////////////////////////////////////////////////////
55 ////////////////////////////////////////////////////////////////////////
56 
57 
58 
59 //======================================================================
60 /// TFourierDecomposedHelmholtzElement<NNODE_1D> elements are
61 /// isoparametric triangular
62 /// FourierDecomposedHelmholtz elements with NNODE_1D nodal points along each
63 /// element edge. Inherits from TElement and FourierDecomposedHelmholtzEquations
64 //======================================================================
65 template <unsigned NNODE_1D>
66 class TFourierDecomposedHelmholtzElement : public TElement<2,NNODE_1D>,
68  public virtual ElementWithZ2ErrorEstimator
69 {
70 
71  public:
72 
73  ///\short Constructor: Call constructors for TElement and
74  /// FourierDecomposedHelmholtz equations
77  { }
78 
79 
80  /// Broken copy constructor
83  {
84  BrokenCopy::broken_copy("TFourierDecomposedHelmholtzElement");
85  }
86 
87  /// Broken assignment operator
88 //Commented out broken assignment operator because this can lead to a conflict warning
89 //when used in the virtual inheritence hierarchy. Essentially the compiler doesn't
90 //realise that two separate implementations of the broken function are the same and so,
91 //quite rightly, it shouts.
92  /*void operator=(const TFourierDecomposedHelmholtzElement<NNODE_1D>&)
93  {
94  BrokenCopy::broken_assign("TFourierDecomposedHelmholtzElement");
95  }*/
96 
97  /// \short Access function for Nvalue: # of `values' (pinned or dofs)
98  /// at node n (always returns the same value at every node, 2)
99  inline unsigned required_nvalue(const unsigned &n) const
100  {return Initial_Nvalue;}
101 
102  /// \short Output function:
103  /// r,z,u
104  void output(std::ostream &outfile)
105  {
107  }
108 
109  /// \short Output function:
110  /// r,z,u n_plot^2 plot points
111  void output(std::ostream &outfile, const unsigned &n_plot)
112  {
114  }
115 
116 
117  /// \short C-style output function:
118  /// r,z,u or x,y,z,u
119  void output(FILE* file_pt)
120  {
122  }
123 
124 
125  /// \short C-style output function:
126  /// r,z,u at n_plot^2 plot points
127  void output(FILE* file_pt, const unsigned &n_plot)
128  {
130  }
131 
132 
133  /// \short Output function for an exact solution:
134  /// r,z,u_exact
135  void output_fct(std::ostream &outfile, const unsigned &n_plot,
137  {
138  FourierDecomposedHelmholtzEquations::output_fct(outfile,n_plot,exact_soln_pt);
139  }
140 
141 
142  /// \short Output function for a time-dependent exact solution.
143  /// x,y,u_exact (calls the steady version)
144  void output_fct(std::ostream &outfile, const unsigned &n_plot,
145  const double& time,
147  {
149  outfile,n_plot,time,exact_soln_pt);
150  }
151 
152  protected:
153 
154 /// Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
156  const Vector<double> &s,
157  Shape &psi,
158  DShape &dpsidx,
159  Shape &test,
160  DShape &dtestdx) const;
161 
162 
163 /// Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
165  const unsigned &ipt,
166  Shape &psi,
167  DShape &dpsidx,
168  Shape &test,
169  DShape &dtestdx)
170  const;
171 
172 
173  /// \short Order of recovery shape functions for Z2 error estimation:
174  /// Same order as shape functions.
175  unsigned nrecovery_order() {return (NNODE_1D-1);}
176 
177  /// Number of 'flux' terms for Z2 error estimation
178  unsigned num_Z2_flux_terms() {return 2*2;}
179 
180  /// \short Get 'flux' for Z2 error recovery: Standard flux from
181  /// UnsteadyHeat equations
183  {
184  Vector<std::complex <double> > complex_flux(2);
185  this->get_flux(s,complex_flux);
186  unsigned count=0;
187  for (unsigned i=0;i<2;i++)
188  {
189  flux[count++]=complex_flux[i].real();
190  flux[count++]=complex_flux[i].imag();
191  }
192  }
193 
194  /// \short Number of vertex nodes in the element
195  unsigned nvertex_node() const
197 
198  /// \short Pointer to the j-th vertex node in the element
199  Node* vertex_node_pt(const unsigned& j) const
201 
202  private:
203 
204  /// Static unsigned that holds the (same) number of variables at every node
205  static const unsigned Initial_Nvalue;
206 
207 
208 };
209 
210 
211 
212 
213 //Inline functions:
214 
215 
216 //======================================================================
217 /// Define the shape functions and test functions and derivatives
218 /// w.r.t. global coordinates and return Jacobian of mapping.
219 ///
220 /// Galerkin: Test functions = shape functions
221 //======================================================================
222 template<unsigned NNODE_1D>
225  const Vector<double> &s,
226  Shape &psi,
227  DShape &dpsidx,
228  Shape &test,
229  DShape &dtestdx) const
230  {
231  unsigned n_node = this->nnode();
232 
233  //Call the geometrical shape functions and derivatives
234  double J = this->dshape_eulerian(s,psi,dpsidx);
235 
236  //Loop over the test functions and derivatives and set them equal to the
237  //shape functions
238  for(unsigned i=0;i<n_node;i++)
239  {
240  test[i] = psi[i];
241  dtestdx(i,0) = dpsidx(i,0);
242  dtestdx(i,1) = dpsidx(i,1);
243  }
244 
245  //Return the jacobian
246  return J;
247  }
248 
249 
250 
251 //======================================================================
252 /// Define the shape functions and test functions and derivatives
253 /// w.r.t. global coordinates and return Jacobian of mapping.
254 ///
255 /// Galerkin: Test functions = shape functions
256 //======================================================================
257 template<unsigned NNODE_1D>
260  const unsigned &ipt,
261  Shape &psi,
262  DShape &dpsidx,
263  Shape &test,
264  DShape &dtestdx) const
265  {
266 
267  //Call the geometrical shape functions and derivatives
268  double J = this->dshape_eulerian_at_knot(ipt,psi,dpsidx);
269 
270  //Set the pointers of the test functions
271  test = psi;
272  dtestdx = dpsidx;
273 
274  //Return the jacobian
275  return J;
276 
277  }
278 
279 
280 
281 //=======================================================================
282 /// Face geometry for the TFourierDecomposedHelmholtzElement elements:
283 /// The spatial dimension of the face elements is one lower than that of the
284 /// bulk element but they have the same number of points
285 /// along their 1D edges.
286 //=======================================================================
287 template<unsigned NNODE_1D>
289  public virtual TElement<1,NNODE_1D>
290  {
291 
292  public:
293 
294  /// \short Constructor: Call the constructor for the
295  /// appropriate lower-dimensional TElement
296  FaceGeometry() : TElement<1,NNODE_1D>() {}
297 
298 };
299 
300 
301 }
302 
303 #endif
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
virtual double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Return the geometric shape functions and also first derivatives w.r.t. global coordinates at the ipt-...
Definition: elements.cc:3254
cstr elem_len * i
Definition: cfortran.h:607
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
void output(std::ostream &outfile)
Output with default number of plot points.
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as ...
Definition: elements.h:1729
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:852
void output(FILE *file_pt)
C-style output function: r,z,u or x,y,z,u.
double dshape_and_dtest_eulerian_fourier_decomposed_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 nvertex_node() const
Number of vertex nodes in the element.
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(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as ...
Definition: elements.h:1723
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
Definition: elements.cc:3227
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement. ...
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: r,z,u at n_plot^2 plot points.
static char t char * s
Definition: cfortran.h:572
void output(std::ostream &outfile)
Output function: r,z,u.
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) ...
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
unsigned num_Z2_flux_terms()
Number of &#39;flux&#39; terms for Z2 error estimation.
double dshape_and_dtest_eulerian_at_knot_fourier_decomposed_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.
TFourierDecomposedHelmholtzElement()
Constructor: Call constructors for TElement and FourierDecomposedHelmholtz equations.
TFourierDecomposedHelmholtzElement(const TFourierDecomposedHelmholtzElement< NNODE_1D > &dummy)
Broken copy constructor.
static const unsigned Initial_Nvalue
Static unsigned that holds the (same) number of variables at every node.
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2146
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: r,z,u n_plot^2 plot points.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get &#39;flux&#39; for Z2 error recovery: Standard flux from UnsteadyHeat equations.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: r,z,u_exact.