axisym_displ_based_fvk_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-2024 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 axisymmetric FoepplvonKarman elements
27 #ifndef OOMPH_AXISYM_DISPL_BASED_FOEPPLVONKARMAN_ELEMENTS_HEADER
28 #define OOMPH_AXISYM_DISPL_BASED_FOEPPLVONKARMAN_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 // OOMPH-LIB headers
37 #include "../generic/nodes.h"
38 #include "../generic/Qelements.h"
39 #include "../generic/oomph_utilities.h"
40 
41 namespace oomph
42 {
43  //=============================================================
44  /// A class for all isoparametric elements that solve the
45  /// axisYm Foeppl von Karman equations in a displacement based formulation.
46  ///
47  /// This contains the generic maths. Shape functions, geometric
48  /// mapping etc. must get implemented in derived class.
49  //=============================================================
51  {
52  public:
53  /// Function pointer to pressure function fct(r,f(r)) --
54  /// r is a Vector!
55  typedef void (*AxisymFoepplvonKarmanPressureFctPt)(const double& r,
56  double& f);
57 
58  /// Constructor (must initialise the Pressure_fct_pt). Also
59  /// set physical parameters to their default values.
61 
62  /// Broken copy constructor
64  const AxisymFoepplvonKarmanEquations& dummy) = delete;
65 
66  /// Broken assignment operator
68 
69  /// Poisson's ratio
70  const double& nu() const
71  {
72 #ifdef PARANOID
73  if (Nu_pt == 0)
74  {
75  std::stringstream error_stream;
76  error_stream << "Nu has not yet been set!" << std::endl;
77  throw OomphLibError(
78  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
79  }
80 #endif
81  return *Nu_pt;
82  }
83 
84  /// Pointer to Poisson's ratio
85  double*& nu_pt()
86  {
87  return Nu_pt;
88  }
89 
90  /// FvK parameter
91  const double& eta() const
92  {
93  return *Eta_pt;
94  }
95 
96  /// Pointer to FvK parameter
97  double*& eta_pt()
98  {
99  return Eta_pt;
100  }
101 
102  /// Return the index at which the i-th unknown value
103  /// is stored. The default value, i, is appropriate for single-physics
104  /// problems. By default, these are:
105  /// 0: transverse displacement w
106  /// 1: laplacian w
107  /// 2: radial displacement u
108  /// In derived multi-physics elements, this function should be overloaded
109  /// to reflect the chosen storage scheme. Note that these equations require
110  /// that the unknown is always stored at the same index at each node.
111  virtual inline unsigned nodal_index_fvk(const unsigned& i = 0) const
112  {
113  return i;
114  }
115 
116  /// Output with default number of plot points
117  void output(std::ostream& outfile)
118  {
119  const unsigned n_plot = 5;
120  output(outfile, n_plot);
121  }
122 
123  /// Output FE representation of soln: r,w,u,sigma_r_r,sigma_phi_phi
124  /// at n_plot plot points
125  void output(std::ostream& outfile, const unsigned& n_plot);
126 
127  /// C_style output with default number of plot points
128  void output(FILE* file_pt)
129  {
130  const unsigned n_plot = 5;
131  output(file_pt, n_plot);
132  }
133 
134  /// C-style output FE representation of soln: r,w at
135  /// n_plot plot points
136  void output(FILE* file_pt, const unsigned& n_plot);
137 
138  /// Output exact soln: r,w_exact at n_plot plot points
139  void output_fct(std::ostream& outfile,
140  const unsigned& n_plot,
142 
143  /// Output exact soln: r,w_exact at
144  /// n_plot plot points (dummy time-dependent version to
145  /// keep intel compiler happy)
146  virtual void output_fct(
147  std::ostream& outfile,
148  const unsigned& n_plot,
149  const double& time,
151  {
152  throw OomphLibError(
153  "There is no time-dependent output_fct() for Foeppl von Karman"
154  "elements ",
155  OOMPH_CURRENT_FUNCTION,
156  OOMPH_EXCEPTION_LOCATION);
157  }
158 
159  /// Get error against and norm of exact solution
160  void compute_error(std::ostream& outfile,
162  double& error,
163  double& norm);
164 
165 
166  /// Dummy, time dependent error checker
167  void compute_error(std::ostream& outfile,
169  const double& time,
170  double& error,
171  double& norm)
172  {
173  throw OomphLibError(
174  "There is no time-dependent compute_error() for Foeppl von Karman"
175  "elements",
176  OOMPH_CURRENT_FUNCTION,
177  OOMPH_EXCEPTION_LOCATION);
178  }
179 
180  /// Access function: Pointer to pressure function
182  {
183  return Pressure_fct_pt;
184  }
185 
186  /// Access function: Pointer to pressure function. Const version
188  {
189  return Pressure_fct_pt;
190  }
191 
192  /// Get pressure term at (Eulerian) position r. This function is
193  /// virtual to allow overloading in multi-physics problems where
194  /// the strength of the pressure function might be determined by
195  /// another system of equations.
196  inline virtual void get_pressure_fvk(const unsigned& ipt,
197  const double& r,
198  double& pressure) const
199  {
200  // If no pressure function has been set, return zero
201  if (Pressure_fct_pt == 0)
202  {
203  pressure = 0.0;
204  }
205  else
206  {
207  // Get pressure strength
208  (*Pressure_fct_pt)(r, pressure);
209  }
210  }
211 
212  /// Get gradient of deflection: gradient[i] = dw/dr_i */
214  Vector<double>& gradient) const
215  {
216  // Find out how many nodes there are in the element
217  const unsigned n_node = nnode();
218 
219  // Get the index at which the unknown is stored
220  const unsigned w_nodal_index = nodal_index_fvk(0);
221 
222  // Set up memory for the shape and test functions
223  Shape psi(n_node);
224  DShape dpsidr(n_node, 1);
225 
226  // Call the derivatives of the shape and test functions
227  dshape_eulerian(s, psi, dpsidr);
228 
229  // Initialise to zero
230  gradient[0] = 0.0;
231 
232  // Loop over nodes
233  for (unsigned l = 0; l < n_node; l++)
234  {
235  gradient[0] += this->nodal_value(l, w_nodal_index) * dpsidr(l, 0);
236  }
237  }
238 
239  /// Fill in the residuals with this element's contribution
241 
242 
243  // hierher Jacobian not yet implemented
244  // void fill_in_contribution_to_jacobian(Vector<double> &residuals,
245  // DenseMatrix<double> &jacobian);
246 
247 
248  /// Return FE representation of transverse displacement
249  inline double interpolated_w_fvk(const Vector<double>& s) const
250  {
251  // Find number of nodes
252  const unsigned n_node = nnode();
253 
254  // Get the index at which the transverse displacement unknown is stored
255  const unsigned w_nodal_index = nodal_index_fvk(0);
256 
257  // Local shape function
258  Shape psi(n_node);
259 
260  // Find values of shape function
261  shape(s, psi);
262 
263  // Initialise value of u
264  double interpolated_w = 0.0;
265 
266  // Loop over the local nodes and sum
267  for (unsigned l = 0; l < n_node; l++)
268  {
269  interpolated_w += this->nodal_value(l, w_nodal_index) * psi[l];
270  }
271 
272  return (interpolated_w);
273  }
274 
275 
276  /// Return FE representation of radial displacement
277  inline double interpolated_u_fvk(const Vector<double>& s) const
278  {
279  // Find number of nodes
280  const unsigned n_node = nnode();
281 
282  // Get the index at which the radial displacement unknown is stored
283  const unsigned u_nodal_index = nodal_index_fvk(2);
284 
285  // Local shape function
286  Shape psi(n_node);
287 
288  // Find values of shape function
289  shape(s, psi);
290 
291  // Initialise value of u
292  double interpolated_u = 0.0;
293 
294  // Loop over the local nodes and sum
295  for (unsigned l = 0; l < n_node; l++)
296  {
297  interpolated_u += this->nodal_value(l, u_nodal_index) * psi[l];
298  }
299 
300  return (interpolated_u);
301  }
302 
303 
304  /// Compute in-plane stresses
306  double& sigma_r_r,
307  double& sigma_phi_phi) const;
308 
309 
310  /// Self-test: Return 0 for OK
311  unsigned self_test();
312 
313 
314  // switch back on and test!
315 
316 
317  /// Sets a flag to signify that we are solving the linear,
318  /// pure bending equations, and pin all the nodal values that will
319  /// not be used in this case
321  {
322  // Set the boolean flag
323  Linear_bending_model = true;
324 
325  // Get the index of the first FvK nodal value
326  unsigned first_fvk_nodal_index = nodal_index_fvk();
327 
328  // Get the total number of FvK nodal values (assuming they are stored
329  // contiguously) at node 0 (it's the same at all nodes anyway)
330  unsigned total_fvk_nodal_indices = 3;
331 
332  // Get the number of nodes in this element
333  unsigned n_node = nnode();
334 
335  // Loop over the appropriate nodal indices
336  for (unsigned index = first_fvk_nodal_index + 2;
337  index < first_fvk_nodal_index + total_fvk_nodal_indices;
338  index++)
339  {
340  // Loop over the nodes in the element
341  for (unsigned inod = 0; inod < n_node; inod++)
342  {
343  // Pin the nodal value at the current index
344  node_pt(inod)->pin(index);
345  }
346  }
347  }
348 
349 
350  protected:
351  /// Shape/test functions and derivs w.r.t. to global coords at
352  /// local coord. s; return Jacobian of mapping
354  const Vector<double>& s,
355  Shape& psi,
356  DShape& dpsidr,
357  Shape& test,
358  DShape& dtestdr) const = 0;
359 
360 
361  /// Shape/test functions and derivs w.r.t. to global coords at
362  /// integration point ipt; return Jacobian of mapping
364  const unsigned& ipt,
365  Shape& psi,
366  DShape& dpsidr,
367  Shape& test,
368  DShape& dtestdr) const = 0;
369 
370  /// Pointer to FvK parameter
371  double* Eta_pt;
372 
373  /// Pointer to pressure function:
375 
376  /// Pointer to Poisson's ratio
377  double* Nu_pt;
378 
379  /// Flag which stores whether we are using a linear,
380  /// pure bending model instead of the full non-linear Foeppl-von Karman
382 
383  private:
384  /// Default value for physical constants
386  };
387 
388 
389  /// ////////////////////////////////////////////////////////////////////////
390  /// ////////////////////////////////////////////////////////////////////////
391  /// ////////////////////////////////////////////////////////////////////////
392 
393 
394  //======================================================================
395  /// Axisym FoepplvonKarmanElement elements are 1D
396  /// Foeppl von Karman elements with isoparametric interpolation for the
397  /// function.
398  //======================================================================
399  template<unsigned NNODE_1D>
401  : public virtual QElement<1, NNODE_1D>,
402  public virtual AxisymFoepplvonKarmanEquations
403  {
404  private:
405  /// Static int that holds the number of variables at
406  /// nodes: always the same
407  static const unsigned Initial_Nvalue;
408 
409  public:
410  /// Constructor: Call constructors for QElement and
411  /// AxisymFoepplvonKarmanEquations
413  : QElement<1, NNODE_1D>(), AxisymFoepplvonKarmanEquations()
414  {
415  }
416 
417  /// Broken copy constructor
419  const AxisymFoepplvonKarmanElement<NNODE_1D>& dummy) = delete;
420 
421  /// Broken assignment operator
423 
424  /// Required # of `values' (pinned or dofs)
425  /// at node n
426  inline unsigned required_nvalue(const unsigned& n) const
427  {
428  return Initial_Nvalue;
429  }
430 
431 
432  /// Output function:
433  /// r,w,u,sigma_r_r,sigma_phi_phi
434  void output(std::ostream& outfile)
435  {
437  }
438 
439  /// Output function:
440  /// r,w,u,sigma_r_r,sigma_phi_phi at n_plot plot points
441  void output(std::ostream& outfile, const unsigned& n_plot)
442  {
444  }
445 
446  /// C-style output function:
447  /// r,w
448  void output(FILE* file_pt)
449  {
451  }
452 
453  /// C-style output function:
454  /// r,w at n_plot plot points
455  void output(FILE* file_pt, const unsigned& n_plot)
456  {
458  }
459 
460  /// Output function for an exact solution:
461  /// r,w_exact at n_plot plot points
462  void output_fct(std::ostream& outfile,
463  const unsigned& n_plot,
465  {
467  outfile, n_plot, exact_soln_pt);
468  }
469 
470  /// Output function for a time-dependent exact solution.
471  /// r,w_exact at n_plot plot points
472  /// (Calls the steady version)
473  void output_fct(std::ostream& outfile,
474  const unsigned& n_plot,
475  const double& time,
477  {
479  outfile, n_plot, time, exact_soln_pt);
480  }
481 
482 
483  protected:
484  /// Shape, test functions & derivs. w.r.t. to global coords.
485  /// Return Jacobian.
487  Shape& psi,
488  DShape& dpsidr,
489  Shape& test,
490  DShape& dtestdr) const;
491 
492  /// Shape, test functions & derivs. w.r.t. to global coords. at
493  /// integration point ipt. Return Jacobian.
495  const unsigned& ipt,
496  Shape& psi,
497  DShape& dpsidr,
498  Shape& test,
499  DShape& dtestdr) const;
500  };
501 
502 
503  // Inline functions:
504 
505  //======================================================================
506  /// Define the shape functions and test functions and derivatives
507  /// w.r.t. global coordinates and return Jacobian of mapping.
508  ///
509  /// Galerkin: Test functions = shape functions
510  //======================================================================
511  template<unsigned NNODE_1D>
512  double AxisymFoepplvonKarmanElement<
513  NNODE_1D>::dshape_and_dtest_eulerian_axisym_fvk(const Vector<double>& s,
514  Shape& psi,
515  DShape& dpsidr,
516  Shape& test,
517  DShape& dtestdr) const
518 
519  {
520  // Call the geometrical shape functions and derivatives
521  const double J = this->dshape_eulerian(s, psi, dpsidr);
522 
523  // Set the test functions equal to the shape functions
524  test = psi;
525  dtestdr = dpsidr;
526 
527  // Return the jacobian
528  return J;
529  }
530 
531 
532  //======================================================================
533  /// Define the shape functions and test functions and derivatives
534  /// w.r.t. global coordinates and return Jacobian of mapping.
535  ///
536  /// Galerkin: Test functions = shape functions
537  //======================================================================
538  template<unsigned NNODE_1D>
541  Shape& psi,
542  DShape& dpsidr,
543  Shape& test,
544  DShape& dtestdr) const
545  {
546  // Call the geometrical shape functions and derivatives
547  const double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidr);
548 
549  // Set the pointers of the test functions
550  test = psi;
551  dtestdr = dpsidr;
552 
553  // Return the jacobian
554  return J;
555  }
556 
557 } // namespace oomph
558 
559 #endif
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: r,w at n_plot plot points.
double dshape_and_dtest_eulerian_axisym_fvk(const Vector< double > &s, Shape &psi, DShape &dpsidr, Shape &test, DShape &dtestdr) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
void output(FILE *file_pt)
C-style output function: r,w.
double dshape_and_dtest_eulerian_at_knot_axisym_fvk(const unsigned &ipt, Shape &psi, DShape &dpsidr, Shape &test, DShape &dtestdr) const
Shape, test functions & derivs. w.r.t. to global coords. at integration point ipt....
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. r,w_exact at n_plot plot points (Calls the stead...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: r,w,u,sigma_r_r,sigma_phi_phi at n_plot plot points.
AxisymFoepplvonKarmanElement(const AxisymFoepplvonKarmanElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile)
Output function: r,w,u,sigma_r_r,sigma_phi_phi.
static const unsigned Initial_Nvalue
Static int that holds the number of variables at nodes: always the same.
unsigned required_nvalue(const unsigned &n) const
Required # of ‘values’ (pinned or dofs) at node n.
void operator=(const AxisymFoepplvonKarmanElement< NNODE_1D > &)=delete
Broken assignment operator.
AxisymFoepplvonKarmanElement()
Constructor: Call constructors for QElement and AxisymFoepplvonKarmanEquations.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: r,w_exact at n_plot plot points.
A class for all isoparametric elements that solve the axisYm Foeppl von Karman equations in a displac...
AxisymFoepplvonKarmanPressureFctPt pressure_fct_pt() const
Access function: Pointer to pressure function. Const version.
AxisymFoepplvonKarmanPressureFctPt Pressure_fct_pt
Pointer to pressure function:
static double Default_Physical_Constant_Value
Default value for physical constants.
void operator=(const AxisymFoepplvonKarmanEquations &)=delete
Broken assignment operator.
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,w_exact at n_plot plot points (dummy time-dependent version to keep intel compil...
void output(std::ostream &outfile)
Output with default number of plot points.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals with this element's contribution.
void output(FILE *file_pt)
C_style output with default number of plot points.
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
unsigned self_test()
Self-test: Return 0 for OK.
virtual double dshape_and_dtest_eulerian_at_knot_axisym_fvk(const unsigned &ipt, Shape &psi, DShape &dpsidr, Shape &test, DShape &dtestdr) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
virtual unsigned nodal_index_fvk(const unsigned &i=0) const
Return the index at which the i-th unknown value is stored. The default value, i, is appropriate for ...
void(* AxisymFoepplvonKarmanPressureFctPt)(const double &r, double &f)
Function pointer to pressure function fct(r,f(r)) – r is a Vector!
void get_gradient_of_deflection(const Vector< double > &s, Vector< double > &gradient) const
Get gradient of deflection: gradient[i] = dw/dr_i *‍/.
void interpolated_stress(const Vector< double > &s, double &sigma_r_r, double &sigma_phi_phi) const
Compute in-plane stresses.
AxisymFoepplvonKarmanPressureFctPt & pressure_fct_pt()
Access function: Pointer to pressure function.
AxisymFoepplvonKarmanEquations()
Constructor (must initialise the Pressure_fct_pt). Also set physical parameters to their default valu...
virtual double dshape_and_dtest_eulerian_axisym_fvk(const Vector< double > &s, Shape &psi, DShape &dpsidr, Shape &test, DShape &dtestdr) const =0
Shape/test functions and derivs w.r.t. to global coords at local coord. s; return Jacobian of mapping...
bool Linear_bending_model
Flag which stores whether we are using a linear, pure bending model instead of the full non-linear Fo...
const double & eta() const
FvK parameter.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,w_exact at n_plot plot points.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
void use_linear_bending_model()
Sets a flag to signify that we are solving the linear, pure bending equations, and pin all the nodal ...
double interpolated_u_fvk(const Vector< double > &s) const
Return FE representation of radial displacement.
double interpolated_w_fvk(const Vector< double > &s) const
Return FE representation of transverse displacement.
virtual void get_pressure_fvk(const unsigned &ipt, const double &r, double &pressure) const
Get pressure term at (Eulerian) position r. This function is virtual to allow overloading in multi-ph...
double *& eta_pt()
Pointer to FvK parameter.
AxisymFoepplvonKarmanEquations(const AxisymFoepplvonKarmanEquations &dummy)=delete
Broken copy constructor.
const double & nu() const
Poisson's ratio.
double *& nu_pt()
Pointer to Poisson's ratio.
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition: shape.h:278
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
A general Finite Element class.
Definition: elements.h:1317
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2179
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:2597
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...
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2214
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition: elements.h:1763
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:3328
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
Definition: elements.h:1769
An OomphLibError object which should be thrown when an run-time error is encountered....
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: Qelements.h:459
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
//////////////////////////////////////////////////////////////////// ////////////////////////////////...