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// Copyright (C) 2006-2023 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 TPoisson elements
27 #ifndef OOMPH_TPOISSON_ELEMENTS_HEADER
28 #define OOMPH_TPOISSON_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 "poisson_elements.h"
43 
44 namespace oomph
45 {
46  /// //////////////////////////////////////////////////////////////////////
47  /// //////////////////////////////////////////////////////////////////////
48  // TPoissonElement
49  /// /////////////////////////////////////////////////////////////////////
50  /// /////////////////////////////////////////////////////////////////////
51 
52 
53  //======================================================================
54  /// TPoissonElement<DIM,NNODE_1D> elements are isoparametric triangular
55  /// DIM-dimensional Poisson elements with NNODE_1D nodal points along each
56  /// element edge. Inherits from TElement and PoissonEquations
57  //======================================================================
58  template<unsigned DIM, unsigned NNODE_1D>
59  class TPoissonElement : public virtual TElement<DIM, NNODE_1D>,
60  public virtual PoissonEquations<DIM>,
61  public virtual ElementWithZ2ErrorEstimator
62  {
63  public:
64  /// Constructor: Call constructors for TElement and
65  /// Poisson equations
66  TPoissonElement() : TElement<DIM, NNODE_1D>(), PoissonEquations<DIM>() {}
67 
68 
69  /// Broken copy constructor
71 
72  /// Broken assignment operator
74 
75  /// Access function for Nvalue: # of `values' (pinned or dofs)
76  /// at node n (always returns the same value at every node, 1)
77  inline unsigned required_nvalue(const unsigned& n) const
78  {
79  return Initial_Nvalue;
80  }
81 
82  /// Output function:
83  /// x,y,u or x,y,z,u
84  void output(std::ostream& outfile)
85  {
87  }
88 
89  /// Output function:
90  /// x,y,u or x,y,z,u at n_plot^DIM plot points
91  void output(std::ostream& outfile, const unsigned& n_plot)
92  {
93  PoissonEquations<DIM>::output(outfile, n_plot);
94  }
95 
96 
97  /// C-style output function:
98  /// x,y,u or x,y,z,u
99  void output(FILE* file_pt)
100  {
102  }
103 
104 
105  /// C-style output function:
106  /// x,y,u or x,y,z,u at n_plot^DIM plot points
107  void output(FILE* file_pt, const unsigned& n_plot)
108  {
109  PoissonEquations<DIM>::output(file_pt, n_plot);
110  }
111 
112 
113  /// Output function for an exact solution:
114  /// x,y,u_exact
115  void output_fct(std::ostream& outfile,
116  const unsigned& n_plot,
118  {
119  PoissonEquations<DIM>::output_fct(outfile, n_plot, exact_soln_pt);
120  }
121 
122 
123  /// Output function for a time-dependent exact solution.
124  /// x,y,u_exact (calls the steady version)
125  void output_fct(std::ostream& outfile,
126  const unsigned& n_plot,
127  const double& time,
129  {
130  PoissonEquations<DIM>::output_fct(outfile, n_plot, time, exact_soln_pt);
131  }
132 
133  protected:
134  /// Shape, test functions & derivs. w.r.t. to global coords. Return
135  /// Jacobian.
137  Shape& psi,
138  DShape& dpsidx,
139  Shape& test,
140  DShape& dtestdx) const;
141 
142 
143  /// Shape, test functions & derivs. w.r.t. to global coords. Return
144  /// Jacobian.
146  const unsigned& ipt,
147  Shape& psi,
148  DShape& dpsidx,
149  Shape& test,
150  DShape& dtestdx) const;
151 
152  /// Shape/test functions and derivs w.r.t. to global coords at
153  /// integration point ipt; return Jacobian of mapping (J). Also compute
154  /// derivatives of dpsidx, dtestdx and J w.r.t. nodal coordinates.
156  const unsigned& ipt,
157  Shape& psi,
158  DShape& dpsidx,
159  RankFourTensor<double>& d_dpsidx_dX,
160  Shape& test,
161  DShape& dtestdx,
162  RankFourTensor<double>& d_dtestdx_dX,
163  DenseMatrix<double>& djacobian_dX) const;
164 
165  /// Order of recovery shape functions for Z2 error estimation:
166  /// Same order as shape functions.
167  unsigned nrecovery_order()
168  {
169  return (NNODE_1D - 1);
170  }
171 
172  /// Number of 'flux' terms for Z2 error estimation
173  unsigned num_Z2_flux_terms()
174  {
175  return DIM;
176  }
177 
178  /// Get 'flux' for Z2 error recovery: Standard flux.from Poisson equations
180  {
181  this->get_flux(s, flux);
182  }
183 
184  /// Number of vertex nodes in the element
185  unsigned nvertex_node() const
186  {
188  }
189 
190  /// Pointer to the j-th vertex node in the element
191  Node* vertex_node_pt(const unsigned& j) const
192  {
194  }
195 
196  private:
197  /// Static unsigned that holds the (same) number of variables at every node
198  static const unsigned Initial_Nvalue;
199  };
200 
201 
202  // Inline functions:
203 
204 
205  //======================================================================
206  /// Define the shape functions and test functions and derivatives
207  /// w.r.t. global coordinates and return Jacobian of mapping.
208  ///
209  /// Galerkin: Test functions = shape functions
210  //======================================================================
211  template<unsigned DIM, unsigned NNODE_1D>
213  const Vector<double>& s,
214  Shape& psi,
215  DShape& dpsidx,
216  Shape& test,
217  DShape& dtestdx) const
218  {
219  unsigned n_node = this->nnode();
220 
221  // Call the geometrical shape functions and derivatives
222  double J = this->dshape_eulerian(s, psi, dpsidx);
223 
224  // Loop over the test functions and derivatives and set them equal to the
225  // shape functions
226  for (unsigned i = 0; i < n_node; i++)
227  {
228  test[i] = psi[i];
229  dtestdx(i, 0) = dpsidx(i, 0);
230  dtestdx(i, 1) = dpsidx(i, 1);
231  }
232 
233  // Return the jacobian
234  return J;
235  }
236 
237 
238  //======================================================================
239  /// Define the shape functions and test functions and derivatives
240  /// w.r.t. global coordinates and return Jacobian of mapping.
241  ///
242  /// Galerkin: Test functions = shape functions
243  //======================================================================
244  template<unsigned DIM, unsigned NNODE_1D>
247  Shape& psi,
248  DShape& dpsidx,
249  Shape& test,
250  DShape& dtestdx) const
251  {
252  // Call the geometrical shape functions and derivatives
253  double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
254 
255  // Set the pointers of the test functions
256  test = psi;
257  dtestdx = dpsidx;
258 
259  // Return the jacobian
260  return J;
261  }
262 
263 
264  //======================================================================
265  /// Define the shape functions (psi) and test functions (test) and
266  /// their derivatives w.r.t. global coordinates (dpsidx and dtestdx)
267  /// and return Jacobian of mapping (J). Additionally compute the
268  /// derivatives of dpsidx, dtestdx and J w.r.t. nodal coordinates.
269  ///
270  /// Galerkin: Test functions = shape functions
271  //======================================================================
272  template<unsigned DIM, unsigned NNODE_1D>
275  const unsigned& ipt,
276  Shape& psi,
277  DShape& dpsidx,
278  RankFourTensor<double>& d_dpsidx_dX,
279  Shape& test,
280  DShape& dtestdx,
281  RankFourTensor<double>& d_dtestdx_dX,
282  DenseMatrix<double>& djacobian_dX) const
283  {
284  // Call the geometrical shape functions and derivatives
285  const double J = this->dshape_eulerian_at_knot(
286  ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
287 
288  // Set the pointers of the test functions
289  test = psi;
290  dtestdx = dpsidx;
291  d_dtestdx_dX = d_dpsidx_dX;
292 
293  // Return the jacobian
294  return J;
295  }
296 
297 
298  //=======================================================================
299  /// Face geometry for the TPoissonElement elements: The spatial
300  /// dimension of the face elements is one lower than that of the
301  /// bulk element but they have the same number of points
302  /// along their 1D edges.
303  //=======================================================================
304  template<unsigned DIM, unsigned NNODE_1D>
305  class FaceGeometry<TPoissonElement<DIM, NNODE_1D>>
306  : public virtual TElement<DIM - 1, NNODE_1D>
307  {
308  public:
309  /// Constructor: Call the constructor for the
310  /// appropriate lower-dimensional TElement
311  FaceGeometry() : TElement<DIM - 1, NNODE_1D>() {}
312  };
313 
314  //=======================================================================
315  /// Face geometry for the 1D TPoissonElement elements: Point elements
316  //=======================================================================
317  template<unsigned NNODE_1D>
318  class FaceGeometry<TPoissonElement<1, NNODE_1D>> : public virtual PointElement
319  {
320  public:
321  /// Constructor: Call the constructor for the
322  /// appropriate lower-dimensional TElement
324  };
325 
326 
327 } // namespace oomph
328 
329 #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:4998
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition: elements.h:1759
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
Definition: elements.h:1765
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:3439
A class for all isoparametric elements that solve the Poisson equations.
void output(std::ostream &outfile)
Output with default number of plot points.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
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.
////////////////////////////////////////////////////////////////// //////////////////////////////////...
Definition: matrices.h:1701
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)
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.
TPoissonElement()
Constructor: Call constructors for TElement and Poisson equations.
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.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
unsigned nvertex_node() const
Number of vertex nodes in the element.
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.
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.
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.
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.
unsigned required_nvalue(const unsigned &n) const
Access function for Nvalue: # of ‘values’ (pinned or dofs) at node n (always returns the same value a...
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
void operator=(const TPoissonElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
TPoissonElement(const TPoissonElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from Poisson equations.
static const unsigned Initial_Nvalue
Static unsigned that holds the (same) number of variables at every node.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...