Tpoisson_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 TPoisson elements
31 #ifndef OOMPH_TPOISSON_ELEMENTS_HEADER
32 #define OOMPH_TPOISSON_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"
46 #include "poisson_elements.h"
47 
48 namespace oomph
49 {
50 
51 /////////////////////////////////////////////////////////////////////////
52 /////////////////////////////////////////////////////////////////////////
53 // TPoissonElement
54 ////////////////////////////////////////////////////////////////////////
55 ////////////////////////////////////////////////////////////////////////
56 
57 
58 
59 //======================================================================
60 /// TPoissonElement<DIM,NNODE_1D> elements are isoparametric triangular
61 /// DIM-dimensional Poisson elements with NNODE_1D nodal points along each
62 /// element edge. Inherits from TElement and PoissonEquations
63 //======================================================================
64 template <unsigned DIM, unsigned NNODE_1D>
65 class TPoissonElement : public virtual TElement<DIM,NNODE_1D>,
66  public virtual PoissonEquations<DIM>,
67  public virtual ElementWithZ2ErrorEstimator
68 {
69 
70  public:
71 
72  ///\short Constructor: Call constructors for TElement and
73  /// Poisson equations
74  TPoissonElement() : TElement<DIM,NNODE_1D>(), PoissonEquations<DIM>()
75  { }
76 
77 
78  /// Broken copy constructor
80  {
81  BrokenCopy::broken_copy("TPoissonElement");
82  }
83 
84  /// Broken assignment operator
86  {
87  BrokenCopy::broken_assign("TPoissonElement");
88  }
89 
90  /// \short Access function for Nvalue: # of `values' (pinned or dofs)
91  /// at node n (always returns the same value at every node, 1)
92  inline unsigned required_nvalue(const unsigned &n) const
93  {return Initial_Nvalue;}
94 
95  /// \short Output function:
96  /// x,y,u or x,y,z,u
97  void output(std::ostream &outfile)
98  {
100  }
101 
102  /// \short Output function:
103  /// x,y,u or x,y,z,u at n_plot^DIM plot points
104  void output(std::ostream &outfile, const unsigned &n_plot)
105  {
106  PoissonEquations<DIM>::output(outfile,n_plot);
107  }
108 
109 
110  /// \short C-style output function:
111  /// x,y,u or x,y,z,u
112  void output(FILE* file_pt)
113  {
115  }
116 
117 
118  /// \short C-style output function:
119  /// x,y,u or x,y,z,u at n_plot^DIM plot points
120  void output(FILE* file_pt, const unsigned &n_plot)
121  {
122  PoissonEquations<DIM>::output(file_pt,n_plot);
123  }
124 
125 
126  /// \short Output function for an exact solution:
127  /// x,y,u_exact
128  void output_fct(std::ostream &outfile, const unsigned &n_plot,
130  {
131  PoissonEquations<DIM>::output_fct(outfile,n_plot,exact_soln_pt);
132  }
133 
134 
135  /// \short Output function for a time-dependent exact solution.
136  /// x,y,u_exact (calls the steady version)
137  void output_fct(std::ostream &outfile, const unsigned &n_plot,
138  const double& time,
140  {
141  PoissonEquations<DIM>::output_fct(outfile,n_plot,time,exact_soln_pt);
142  }
143 
144 protected:
145 
146 /// Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
147  inline double dshape_and_dtest_eulerian_poisson(const Vector<double> &s,
148  Shape &psi,
149  DShape &dpsidx,
150  Shape &test,
151  DShape &dtestdx) const;
152 
153 
154 /// Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
155  inline double dshape_and_dtest_eulerian_at_knot_poisson(const unsigned &ipt,
156  Shape &psi,
157  DShape &dpsidx,
158  Shape &test,
159  DShape &dtestdx)
160  const;
161 
162  /// \short Shape/test functions and derivs w.r.t. to global coords at
163  /// integration point ipt; return Jacobian of mapping (J). Also compute
164  /// derivatives of dpsidx, dtestdx and J w.r.t. nodal coordinates.
166  const unsigned &ipt,
167  Shape &psi,
168  DShape &dpsidx,
169  RankFourTensor<double> &d_dpsidx_dX,
170  Shape &test,
171  DShape &dtestdx,
172  RankFourTensor<double> &d_dtestdx_dX,
173  DenseMatrix<double> &djacobian_dX) const;
174 
175  /// \short Order of recovery shape functions for Z2 error estimation:
176  /// Same order as shape functions.
177  unsigned nrecovery_order() {return (NNODE_1D-1);}
178 
179  /// Number of 'flux' terms for Z2 error estimation
180  unsigned num_Z2_flux_terms() {return DIM;}
181 
182  /// Get 'flux' for Z2 error recovery: Standard flux.from Poisson equations
184  {this->get_flux(s,flux);}
185 
186  /// \short Number of vertex nodes in the element
187  unsigned nvertex_node() const
189 
190  /// \short Pointer to the j-th vertex node in the element
191  Node* vertex_node_pt(const unsigned& j) const
193 
194 private:
195 
196  /// Static unsigned that holds the (same) number of variables at every node
197  static const unsigned Initial_Nvalue;
198 
199 
200 };
201 
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 //======================================================================
243 /// Define the shape functions and test functions and derivatives
244 /// w.r.t. global coordinates and return Jacobian of mapping.
245 ///
246 /// Galerkin: Test functions = shape functions
247 //======================================================================
248 template<unsigned DIM, unsigned NNODE_1D>
251  const unsigned &ipt,
252  Shape &psi,
253  DShape &dpsidx,
254  Shape &test,
255  DShape &dtestdx) const
256 {
257 
258  //Call the geometrical shape functions and derivatives
259  double J = this->dshape_eulerian_at_knot(ipt,psi,dpsidx);
260 
261  //Set the pointers of the test functions
262  test = psi;
263  dtestdx = dpsidx;
264 
265  //Return the jacobian
266  return J;
267 
268 }
269 
270 
271 
272 //======================================================================
273 /// Define the shape functions (psi) and test functions (test) and
274 /// their derivatives w.r.t. global coordinates (dpsidx and dtestdx)
275 /// and return Jacobian of mapping (J). Additionally compute the
276 /// derivatives of dpsidx, dtestdx and J w.r.t. nodal coordinates.
277 ///
278 /// Galerkin: Test functions = shape functions
279 //======================================================================
280 template<unsigned DIM, unsigned NNODE_1D>
283  const unsigned &ipt,
284  Shape &psi,
285  DShape &dpsidx,
286  RankFourTensor<double> &d_dpsidx_dX,
287  Shape &test,
288  DShape &dtestdx,
289  RankFourTensor<double> &d_dtestdx_dX,
290  DenseMatrix<double> &djacobian_dX) const
291  {
292  // Call the geometrical shape functions and derivatives
293  const double J = this->dshape_eulerian_at_knot(ipt,psi,dpsidx,
294  djacobian_dX,d_dpsidx_dX);
295 
296  // Set the pointers of the test functions
297  test = psi;
298  dtestdx = dpsidx;
299  d_dtestdx_dX = d_dpsidx_dX;
300 
301  //Return the jacobian
302  return J;
303 }
304 
305 
306 
307 //=======================================================================
308 /// Face geometry for the TPoissonElement elements: The spatial
309 /// dimension of the face elements is one lower than that of the
310 /// bulk element but they have the same number of points
311 /// along their 1D edges.
312 //=======================================================================
313 template<unsigned DIM, unsigned NNODE_1D>
314 class FaceGeometry<TPoissonElement<DIM,NNODE_1D> >:
315  public virtual TElement<DIM-1,NNODE_1D>
316 {
317 
318  public:
319 
320  /// \short Constructor: Call the constructor for the
321  /// appropriate lower-dimensional TElement
322  FaceGeometry() : TElement<DIM-1,NNODE_1D>() {}
323 
324 };
325 
326 //=======================================================================
327 /// Face geometry for the 1D TPoissonElement elements: Point elements
328 //=======================================================================
329 template<unsigned NNODE_1D>
330 class FaceGeometry<TPoissonElement<1,NNODE_1D> >:
331  public virtual PointElement
332 {
333 
334  public:
335 
336  /// \short Constructor: Call the constructor for the
337  /// appropriate lower-dimensional TElement
339 
340 };
341 
342 
343 
344 
345 
346 
347 
348 
349 
350 
351 
352 }
353 
354 #endif
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
unsigned nvertex_node() const
Number of vertex nodes in the element.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
unsigned required_nvalue(const unsigned &n) const
Access function for Nvalue: # of `values&#39; (pinned or dofs) at node n (always returns the same value a...
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
unsigned num_Z2_flux_terms()
Number of &#39;flux&#39; terms for Z2 error estimation.
cstr elem_len * i
Definition: cfortran.h:607
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
double dshape_and_dtest_eulerian_at_knot_poisson(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
TPoissonElement()
Constructor: Call constructors for TElement and Poisson equations.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:852
double dshape_and_dtest_eulerian_poisson(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.
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) ...
A Rank 4 Tensor class.
Definition: matrices.h:1625
static const unsigned Initial_Nvalue
Static unsigned that holds the (same) number of variables at every node.
void output(std::ostream &outfile)
Output with default number of plot points.
void output(FILE *file_pt)
C-style output function: x,y,u or x,y,z,u.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement. ...
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. ...
static char t char * s
Definition: cfortran.h:572
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get &#39;flux&#39; for Z2 error recovery: Standard flux.from Poisson equations.
TPoissonElement(const TPoissonElement< DIM, NNODE_1D > &dummy)
Broken copy constructor.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
void broken_assign(const std::string &class_name)
Issue error message and terminate execution.
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 output(std::ostream &outfile, const unsigned &n_plot)
Output function: x,y,u or x,y,z,u at n_plot^DIM plot points.
void output(std::ostream &outfile)
Output function: x,y,u or x,y,z,u.
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.
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2146
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points.
void operator=(const TPoissonElement< DIM, NNODE_1D > &)
Broken assignment operator.