spectral_poisson_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 Spectral Poisson elements
27 #ifndef OOMPH_SPECTRAL_POISSON_ELEMENTS_HEADER
28 #define OOMPH_SPECTRAL_POISSON_ELEMENTS_HEADER
29 
30 // Config header generated by autoconfig
31 #ifdef HAVE_CONFIG_H
32 #include <oomph-lib-config.h>
33 #endif
34 
35 // OOMPH-LIB headers
36 #include "poisson_elements.h"
37 #include "../generic/Qspectral_elements.h"
38 
39 namespace oomph
40 {
41  //======================================================================
42  /// QSpectralPoissonElement elements are linear/quadrilateral/brick-shaped
43  /// Poisson elements with isoparametric spectral interpolation for the
44  /// function. Note that the implementation is PoissonEquations<DIM> does
45  /// not use sum factorisation for the evaluation of the residuals and is,
46  /// therefore, not optimal for higher dimensions.
47  //======================================================================
48  template<unsigned DIM, unsigned NNODE_1D>
50  : public virtual QSpectralElement<DIM, NNODE_1D>,
51  public virtual PoissonEquations<DIM>
52  {
53  private:
54  /// Static array of ints to hold number of variables at
55  /// nodes: Initial_Nvalue[n]
56  static const unsigned Initial_Nvalue;
57 
58  public:
59  /// Constructor: Call constructors for QSpectralElement and
60  /// Poisson equations
62  : QSpectralElement<DIM, NNODE_1D>(), PoissonEquations<DIM>()
63  {
64  }
65 
66  /// Broken copy constructor
68  const QSpectralPoissonElement<DIM, NNODE_1D>& dummy) = delete;
69 
70  /// Broken assignment operator
71  // Commented out broken assignment operator because this can lead to a
72  // conflict warning when used in the virtual inheritence hierarchy.
73  // Essentially the compiler doesn't realise that two separate
74  // implementations of the broken function are the same and so, quite
75  // rightly, it shouts.
76  /*void operator=(const QSpectralPoissonElement<DIM,NNODE_1D>&) =
77  delete;*/
78 
79  /// Required # of `values' (pinned or dofs)
80  /// at node n
81  inline unsigned required_nvalue(const unsigned& n) const
82  {
83  return Initial_Nvalue;
84  }
85 
86  /// Output function:
87  /// x,y,u or x,y,z,u
88  void output(std::ostream& outfile)
89  {
91  }
92 
93  /// Output function:
94  /// x,y,u or x,y,z,u at n_plot^DIM plot points
95  void output(std::ostream& outfile, const unsigned& n_plot)
96  {
97  PoissonEquations<DIM>::output(outfile, n_plot);
98  }
99 
100 
101  /// C-style output function:
102  /// x,y,u or x,y,z,u
103  void output(FILE* file_pt)
104  {
106  }
107 
108 
109  /// C-style output function:
110  /// x,y,u or x,y,z,u at n_plot^DIM plot points
111  void output(FILE* file_pt, const unsigned& n_plot)
112  {
113  PoissonEquations<DIM>::output(file_pt, n_plot);
114  }
115 
116 
117  /// Output function for an exact solution:
118  /// x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points
119  void output_fct(std::ostream& outfile,
120  const unsigned& n_plot,
122  {
123  PoissonEquations<DIM>::output_fct(outfile, n_plot, exact_soln_pt);
124  }
125 
126 
127  /// Output function for a time-dependent exact solution.
128  /// x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points
129  /// (Calls the steady version)
130  void output_fct(std::ostream& outfile,
131  const unsigned& n_plot,
132  const double& time,
134  {
135  PoissonEquations<DIM>::output_fct(outfile, n_plot, time, exact_soln_pt);
136  }
137 
138 
139  protected:
140  /// Shape, test functions & derivs. w.r.t. to global coords. Return
141  /// Jacobian.
143  Shape& psi,
144  DShape& dpsidx,
145  Shape& test,
146  DShape& dtestdx) const;
147 
148 
149  /// Shape, test functions & derivs. w.r.t. to global coords. at
150  /// integration point ipt. Return Jacobian.
152  const unsigned& ipt,
153  Shape& psi,
154  DShape& dpsidx,
155  Shape& test,
156  DShape& dtestdx) const;
157 
158  /// Shape/test functions and derivs w.r.t. to global coords at
159  /// integration point ipt; return Jacobian of mapping (J). Also compute
160  /// derivatives of dpsidx, dtestdx and J w.r.t. nodal coordinates.
162  const unsigned& ipt,
163  Shape& psi,
164  DShape& dpsidx,
165  RankFourTensor<double>& d_dpsidx_dX,
166  Shape& test,
167  DShape& dtestdx,
168  RankFourTensor<double>& d_dtestdx_dX,
169  DenseMatrix<double>& djacobian_dX) const;
170  };
171 
172 
173  // Inline functions:
174 
175 
176  //======================================================================
177  /// Define the shape functions and test functions and derivatives
178  /// w.r.t. global coordinates and return Jacobian of mapping.
179  ///
180  /// Galerkin: Test functions = shape functions
181  //======================================================================
182  template<unsigned DIM, unsigned NNODE_1D>
185  Shape& psi,
186  DShape& dpsidx,
187  Shape& test,
188  DShape& dtestdx) const
189  {
190  // Call the geometrical shape functions and derivatives
191  double J = this->dshape_eulerian(s, psi, dpsidx);
192 
193  // Loop over the test functions and derivatives and set them equal to the
194  // shape functions
195  unsigned nnod = this->nnode();
196  for (unsigned i = 0; i < nnod; i++)
197  {
198  test[i] = psi[i];
199  for (unsigned j = 0; j < DIM; j++)
200  {
201  dtestdx(i, j) = dpsidx(i, j);
202  }
203  }
204 
205  // Return the jacobian
206  return J;
207  }
208 
209  //======================================================================
210  /// Define the shape functions and test functions and derivatives
211  /// w.r.t. global coordinates and return Jacobian of mapping.
212  ///
213  /// Galerkin: Test functions = shape functions
214  //======================================================================
215  template<unsigned DIM, unsigned NNODE_1D>
218  Shape& psi,
219  DShape& dpsidx,
220  Shape& test,
221  DShape& dtestdx) const
222  {
223  // Call the geometrical shape functions and derivatives
224  double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
225 
226  // Set the pointers of the test functions
227  test = psi;
228  dtestdx = dpsidx;
229 
230  // Return the jacobian
231  return J;
232  }
233 
234  //======================================================================
235  /// Define the shape functions (psi) and test functions (test) and
236  /// their derivatives w.r.t. global coordinates (dpsidx and dtestdx)
237  /// and return Jacobian of mapping (J). Additionally compute the
238  /// derivatives of dpsidx, dtestdx and J w.r.t. nodal coordinates.
239  ///
240  /// Galerkin: Test functions = shape functions
241  //======================================================================
242  template<unsigned DIM, unsigned NNODE_1D>
245  const unsigned& ipt,
246  Shape& psi,
247  DShape& dpsidx,
248  RankFourTensor<double>& d_dpsidx_dX,
249  Shape& test,
250  DShape& dtestdx,
251  RankFourTensor<double>& d_dtestdx_dX,
252  DenseMatrix<double>& djacobian_dX) const
253  {
254  // Call the geometrical shape functions and derivatives
255  const double J = this->dshape_eulerian_at_knot(
256  ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
257 
258  // Set the pointers of the test functions
259  test = psi;
260  dtestdx = dpsidx;
261  d_dtestdx_dX = d_dpsidx_dX;
262 
263  // Return the jacobian
264  return J;
265  }
266 
267  /// /////////////////////////////////////////////////////////////////////
268  /// /////////////////////////////////////////////////////////////////////
269  /// /////////////////////////////////////////////////////////////////////
270 
271 
272  //=======================================================================
273  /// Face geometry for the QSpectralPoissonElement elements: The spatial
274  /// dimension of the face elements is one lower than that of the
275  /// bulk element but they have the same number of points
276  /// along their 1D edges.
277  //=======================================================================
278  template<unsigned DIM, unsigned NNODE_1D>
279  class FaceGeometry<QSpectralPoissonElement<DIM, NNODE_1D>>
280  : public virtual QSpectralElement<DIM - 1, NNODE_1D>
281  {
282  public:
283  /// Constructor: Call the constructor for the
284  /// appropriate lower-dimensional QElement
285  FaceGeometry() : QSpectralElement<DIM - 1, NNODE_1D>() {}
286  };
287 
288 
289  //=======================================================================
290  /// Face geometry for the 1D QPoissonElement elements: Point elements
291  //=======================================================================
292  template<unsigned NNODE_1D>
294  : public virtual PointElement
295  {
296  public:
297  /// Constructor: Call the constructor for the
298  /// appropriate lower-dimensional QElement
300  };
301 
302 } // namespace oomph
303 
304 #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
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
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
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
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 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.
General QLegendreElement class.
QSpectralPoissonElement elements are linear/quadrilateral/brick-shaped Poisson elements with isoparam...
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_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points.
QSpectralPoissonElement()
Constructor: Call constructors for QSpectralElement and Poisson equations.
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. at integration point ipt....
QSpectralPoissonElement(const QSpectralPoissonElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
void output(std::ostream &outfile)
Output function: x,y,u or x,y,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 or x,y,z,u_exact at n_plot^DIM plot ...
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.
static const unsigned Initial_Nvalue
Static array of ints to hold number of variables at nodes: Initial_Nvalue[n].
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(FILE *file_pt)
C-style output function: x,y,u or x,y,z,u.
////////////////////////////////////////////////////////////////// //////////////////////////////////...
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
//////////////////////////////////////////////////////////////////// ////////////////////////////////...