axisym_linear_elasticity_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 
27 // Include guards to prevent multiple inclusion of the header
28 #ifndef OOMPH_AXISYMMETRIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
29 #define OOMPH_AXISYMMETRIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
30 
31 // Config header generated by autoconfig
32 #ifdef HAVE_CONFIG_H
33 #include <oomph-lib-config.h>
34 #endif
35 
36 
37 #ifdef OOMPH_HAS_MPI
38 #include "mpi.h"
39 #endif
40 
41 // OOMPH-LIB headers
42 #include "src/generic/Qelements.h"
43 #include "src/generic/Telements.h"
44 #include "src/generic/projection.h"
45 
46 
47 namespace oomph
48 {
49  //=======================================================================
50  /// A base class for elements that solve the axisymmetric (in
51  /// cylindrical polars) equations of linear elasticity.
52  //=======================================================================
54  {
55  public:
56  /// Return the index at which the i-th (i=0: r, i=1: z; i=2: theta)
57  /// unknown displacement component is stored at the nodes. The default
58  /// assignment here (u_r, u_z, u_theta) is appropriate for single-physics
59  /// problems.
60  virtual inline unsigned u_index_axisymmetric_linear_elasticity(
61  const unsigned& i) const
62  {
63  return i;
64  }
65 
66  /// d^2u/dt^2 at local node n
67  double d2u_dt2_axisymmetric_linear_elasticity(const unsigned& n,
68  const unsigned& i) const
69  {
70  // Get the timestepper
72 
73  // Storage for the derivative - initialise to 0
74  double d2u_dt2 = 0.0;
75 
76  // If we are doing an unsteady solve then calculate the derivative
77  if (!time_stepper_pt->is_steady())
78  {
79  // Get the nodal index
80  const unsigned u_nodal_index =
82 
83  // Get the number of values required to represent history
84  const unsigned n_time = time_stepper_pt->ntstorage();
85 
86  // Loop over history values
87  for (unsigned t = 0; t < n_time; t++)
88  {
89  // Add the contribution to the derivative
90  d2u_dt2 +=
91  time_stepper_pt->weight(2, t) * nodal_value(t, n, u_nodal_index);
92  }
93  }
94 
95  return d2u_dt2;
96  }
97 
98 
99  /// du/dt at local node n
100  double du_dt_axisymmetric_linear_elasticity(const unsigned& n,
101  const unsigned& i) const
102  {
103  // Get the timestepper
105 
106  // Storage for the derivative - initialise to 0
107  double du_dt = 0.0;
108 
109  // If we are doing an unsteady solve then calculate the derivative
110  if (!time_stepper_pt->is_steady())
111  {
112  // Get the nodal index
113  const unsigned u_nodal_index =
115 
116  // Get the number of values required to represent history
117  const unsigned n_time = time_stepper_pt->ntstorage();
118 
119  // Loop over history values
120  for (unsigned t = 0; t < n_time; t++)
121  {
122  // Add the contribution to the derivative
123  du_dt +=
124  time_stepper_pt->weight(1, t) * nodal_value(t, n, u_nodal_index);
125  }
126  }
127  return du_dt;
128  }
129 
130  /// Compute vector of FE interpolated displacement u at local coordinate s
132  const Vector<double>& s, Vector<double>& disp) const
133  {
134  // Find number of nodes
135  unsigned n_node = nnode();
136 
137  // Local shape function
138  Shape psi(n_node);
139 
140  // Find values of shape function
141  shape(s, psi);
142 
143  for (unsigned i = 0; i < 3; i++)
144  {
145  // Index at which the nodal value is stored
146  unsigned u_nodal_index = u_index_axisymmetric_linear_elasticity(i);
147 
148  // Initialise value of u
149  disp[i] = 0.0;
150 
151  // Loop over the local nodes and sum
152  for (unsigned l = 0; l < n_node; l++)
153  {
154  const double u_value = nodal_value(l, u_nodal_index);
155 
156  disp[i] += u_value * psi[l];
157  }
158  }
159  }
160 
161  /// Return FE interpolated displacement u[i] (i=0: r, i=1: z; i=2:
162  /// theta) at local coordinate s
164  const Vector<double>& s, const unsigned& i) const
165  {
166  // Find number of nodes
167  unsigned n_node = nnode();
168 
169  // Local shape function
170  Shape psi(n_node);
171 
172  // Find values of shape function
173  shape(s, psi);
174 
175  // Get nodal index at which i-th velocity is stored
176  unsigned u_nodal_index = u_index_axisymmetric_linear_elasticity(i);
177 
178  // Initialise value of u
179  double interpolated_u = 0.0;
180 
181  // Loop over the local nodes and sum
182  for (unsigned l = 0; l < n_node; l++)
183  {
184  const double u_value = nodal_value(l, u_nodal_index);
185 
186  interpolated_u += u_value * psi[l];
187  }
188 
189  return (interpolated_u);
190  }
191 
192 
193  /// Compute vector of FE interpolated velocity du/dt at local coordinate s
195  const Vector<double>& s, Vector<double>& du_dt) const
196  {
197  // Find number of nodes
198  unsigned n_node = nnode();
199 
200  // Local shape function
201  Shape psi(n_node);
202 
203  // Find values of shape function
204  shape(s, psi);
205 
206  // Loop over directions
207  for (unsigned i = 0; i < 3; i++)
208  {
209  // Initialise value of u
210  du_dt[i] = 0.0;
211 
212  // Loop over the local nodes and sum
213  for (unsigned l = 0; l < n_node; l++)
214  {
215  du_dt[i] += du_dt_axisymmetric_linear_elasticity(l, i) * psi[l];
216  }
217  }
218  }
219 
220  /// Compute vector of FE interpolated accel d2u/dt2 at local coordinate s
222  const Vector<double>& s, Vector<double>& d2u_dt2) const
223  {
224  // Find number of nodes
225  unsigned n_node = nnode();
226 
227  // Local shape function
228  Shape psi(n_node);
229 
230  // Find values of shape function
231  shape(s, psi);
232 
233  // Loop over directions
234  for (unsigned i = 0; i < 3; i++)
235  {
236  // Initialise value of u
237  d2u_dt2[i] = 0.0;
238 
239  // Loop over the local nodes and sum
240  for (unsigned l = 0; l < n_node; l++)
241  {
242  d2u_dt2[i] += d2u_dt2_axisymmetric_linear_elasticity(l, i) * psi[l];
243  }
244  }
245  }
246 
247  /// Function pointer to function that specifies the body force
248  /// as a function of the Cartesian coordinates and time FCT(x,b) --
249  /// x and b are Vectors!
250  typedef void (*BodyForceFctPt)(const double& time,
251  const Vector<double>& x,
252  Vector<double>& b);
253 
254  /// Constructor: Set null pointers for constitutive law.
255  /// Set physical parameter values to
256  /// default values, and set body force to zero.
259  Nu_pt(0),
262  {
263  }
264 
265  /// Return the pointer to Young's modulus
266  double*& youngs_modulus_pt()
267  {
268  return Youngs_modulus_pt;
269  }
270 
271  /// Access function to Young's modulus
272  inline double youngs_modulus() const
273  {
274  return (*Youngs_modulus_pt);
275  }
276 
277  /// Access function for Poisson's ratio
278  double& nu() const
279  {
280 #ifdef PARANOID
281  if (Nu_pt == 0)
282  {
283  std::ostringstream error_message;
284  error_message << "No pointer to Poisson's ratio set. Please set one!\n";
285  throw OomphLibError(error_message.str(),
286  OOMPH_CURRENT_FUNCTION,
287  OOMPH_EXCEPTION_LOCATION);
288  }
289 #endif
290  return *Nu_pt;
291  }
292 
293  /// Access function for pointer to Poisson's ratio
294  double*& nu_pt()
295  {
296  return Nu_pt;
297  }
298 
299  /// Access function for pointer to timescale ratio (nondim density)
300  double*& lambda_sq_pt()
301  {
302  return Lambda_sq_pt;
303  }
304 
305  /// Access function for timescale ratio (nondim density)
306  const double& lambda_sq() const
307  {
308  return *Lambda_sq_pt;
309  }
310 
311  /// Access function: Pointer to body force function
313  {
314  return Body_force_fct_pt;
315  }
316 
317  /// Access function: Pointer to body force function (const version)
319  {
320  return Body_force_fct_pt;
321  }
322 
323  /// Evaluate body force at Eulerian coordinate x at present time
324  /// (returns zero vector if no body force function pointer has been set)
325  inline void body_force(const double& time,
326  const Vector<double>& x,
327  Vector<double>& b) const
328  {
329  // If no function has been set, return zero vector
330  if (Body_force_fct_pt == 0)
331  {
332  // Get spatial dimension of element
333  unsigned n = dim();
334  for (unsigned i = 0; i < n; i++)
335  {
336  b[i] = 0.0;
337  }
338  }
339  else
340  {
341  (*Body_force_fct_pt)(time, x, b);
342  }
343  }
344 
345  /// The number of "DOF types" that degrees of freedom in this element
346  /// are sub-divided into: for now lump them all into one DOF type.
347  /// Can be adjusted later
348  unsigned ndof_types() const
349  {
350  return 1;
351  }
352 
353  /// Create a list of pairs for all unknowns in this element,
354  /// so that the first entry in each pair contains the global equation
355  /// number of the unknown, while the second one contains the number
356  /// of the "DOF type" that this unknown is associated with.
357  /// (Function can obviously only be called if the equation numbering
358  /// scheme has been set up.)
360  std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list) const
361  {
362  // temporary pair (used to store DOF lookup prior to being added
363  // to list)
364  std::pair<unsigned long, unsigned> dof_lookup;
365 
366  // number of nodes
367  const unsigned n_node = this->nnode();
368 
369  // Integer storage for local unknown
370  int local_unknown = 0;
371 
372  // Loop over the nodes
373  for (unsigned n = 0; n < n_node; n++)
374  {
375  // Loop over dimension
376  for (unsigned i = 0; i < 3; i++)
377  {
378  // If the variable is free
379  local_unknown = nodal_local_eqn(n, i);
380 
381  // ignore pinned values
382  if (local_unknown >= 0)
383  {
384  // store DOF type lookup in temporary pair: First entry in pair
385  // is global equation number; second entry is DOF type
386  dof_lookup.first = this->eqn_number(local_unknown);
387  dof_lookup.second = 0;
388 
389  // add to list
390  dof_lookup_list.push_front(dof_lookup);
391  }
392  }
393  }
394  }
395 
396 
397  protected:
398  /// Pointer to the Young's modulus
400 
401  /// Pointer to Poisson's ratio
402  double* Nu_pt;
403 
404  /// Timescale ratio (non-dim. density)
405  double* Lambda_sq_pt;
406 
407  /// Pointer to body force function
409 
410  /// Static default value for Young's modulus (1.0 -- for natural
411  /// scaling, i.e. all stresses have been non-dimensionalised by
412  /// the same reference Young's modulus. Setting the "non-dimensional"
413  /// Young's modulus (obtained by de-referencing Youngs_modulus_pt)
414  /// to a number larger than one means that the material is stiffer
415  /// than assumed in the non-dimensionalisation.
417 
418  /// Static default value for timescale ratio (1.0 for natural scaling)
419  static double Default_lambda_sq_value;
420  };
421 
422 
423  /// ////////////////////////////////////////////////////////////////////
424  /// ////////////////////////////////////////////////////////////////////
425  /// ////////////////////////////////////////////////////////////////////
426 
427 
428  //=======================================================================
429  /// A class for elements that solve the axisymmetric (in cylindrical
430  /// polars) equations of linear elasticity
431  //=======================================================================
434  {
435  public:
436  /// Constructor
438 
439  /// Number of values required at node n.
440  unsigned required_nvalue(const unsigned& n) const
441  {
442  return 3;
443  }
444 
445  /// Return the residuals for the equations (the discretised
446  /// principle of virtual displacements)
448  {
450  residuals, GeneralisedElement::Dummy_matrix, 0);
451  }
452 
453 
454  /// The jacobian is calculated by finite differences by default,
455  /// We need only to take finite differences w.r.t. positional variables
456  /// For this element
458  DenseMatrix<double>& jacobian)
459  {
460  // Add the contribution to the residuals
461  this
462  ->fill_in_generic_contribution_to_residuals_axisymmetric_linear_elasticity(
463  residuals, jacobian, 1);
464  }
465 
466 
467  /// Get strain (3x3 entries; r, z, phi)
468  void get_strain(const Vector<double>& s, DenseMatrix<double>& strain);
469 
470  /// Output exact solution: r,z, u_r, u_z, u_theta
471  void output_fct(std::ostream& outfile,
472  const unsigned& nplot,
474 
475  /// Output exact solution: r,z, u_r, u_z, u_theta
476  /// Time dependent version
477  void output_fct(std::ostream& outfile,
478  const unsigned& nplot,
479  const double& time,
481 
482  /// Output: r,z, u_r, u_z, u_theta
483  void output(std::ostream& outfile)
484  {
485  unsigned n_plot = 5;
486  output(outfile, n_plot);
487  }
488 
489  /// Output: r,z, u_r, u_z, u_theta
490  void output(std::ostream& outfile, const unsigned& n_plot);
491 
492  /// C-style output: r,z, u_r, u_z, u_theta
493  void output(FILE* file_pt)
494  {
495  unsigned n_plot = 5;
496  output(file_pt, n_plot);
497  }
498 
499  /// Output: r,z, u_r, u_z, u_theta
500  void output(FILE* file_pt, const unsigned& n_plot);
501 
502  /// Validate against exact solution.
503  /// Solution is provided via function pointer.
504  /// Plot at a given number of plot points and compute L2 error
505  /// and L2 norm of displacement solution over element
506  void compute_error(std::ostream& outfile,
508  double& error,
509  double& norm);
510 
511  /// Validate against exact solution.
512  /// Time-dependent version
513  void compute_error(std::ostream& outfile,
515  const double& time,
516  double& error,
517  double& norm);
518 
519 
520  protected:
521  /// Private helper function to compute residuals and (if requested
522  /// via flag) also the Jacobian matrix.
524  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag);
525  };
526 
527 
528  /// /////////////////////////////////////////////////////////////////////
529  /// /////////////////////////////////////////////////////////////////////
530  /// /////////////////////////////////////////////////////////////////////
531 
532 
533  //===========================================================================
534  /// An Element that solves the equations of axisymmetric (in cylindrical
535  /// polars) linear elasticity, using QElements for the geometry.
536  //============================================================================
537  template<unsigned NNODE_1D>
539  : public virtual QElement<2, NNODE_1D>,
541  {
542  public:
543  /// Constructor
546  {
547  }
548 
549  /// Output function
550  void output(std::ostream& outfile)
551  {
553  }
554 
555  /// Output function
556  void output(std::ostream& outfile, const unsigned& n_plot)
557  {
559  }
560 
561 
562  /// C-style output function
563  void output(FILE* file_pt)
564  {
566  }
567 
568  /// C-style output function
569  void output(FILE* file_pt, const unsigned& n_plot)
570  {
572  }
573  };
574 
575 
576  //============================================================================
577  /// FaceGeometry of a linear
578  /// QAxisymmetricLinearElasticityElement element
579  //============================================================================
580  template<unsigned NNODE_1D>
582  : public virtual QElement<1, NNODE_1D>
583  {
584  public:
585  /// Constructor must call the constructor of the underlying element
586  FaceGeometry() : QElement<1, NNODE_1D>() {}
587  };
588 
589 
590  /* //////////////////////////////////////////////////////////////////////// */
591  /* //////////////////////////////////////////////////////////////////////// */
592  /* //////////////////////////////////////////////////////////////////////// */
593 
594 
595  /* //===========================================================================
596  */
597  /* /// An Element that solves the equations of axisymmetric (in cylindrical */
598  /* /// polars) linear elasticity, using TElements for the geometry. */
599  /* //============================================================================
600  */
601  /* template<unsigned NNODE_1D> */
602  /* class TAxisymmetricLinearElasticityElement : */
603  /* public virtual TElement<2,NNODE_1D>, */
604  /* public virtual AxisymmetricLinearElasticityEquations */
605  /* { */
606  /* public: */
607 
608  /* /// Constructor */
609  /* TAxisymmetricLinearElasticityElement() : */
610  /* TElement<2,NNODE_1D>(), */
611  /* AxisymmetricLinearElasticityEquations() { } */
612 
613  /* /// Output function */
614  /* void output(std::ostream &outfile) */
615  /* {AxisymmetricLinearElasticityEquations::output(outfile);} */
616 
617  /* /// Output function */
618  /* void output(std::ostream &outfile, const unsigned &n_plot) */
619  /* {AxisymmetricLinearElasticityEquations:: */
620  /* output(outfile,n_plot);} */
621 
622  /* /// C-style output function */
623  /* void output(FILE* file_pt) */
624  /* {AxisymmetricLinearElasticityEquations::output(file_pt);} */
625 
626  /* /// C-style output function */
627  /* void output(FILE* file_pt, const unsigned &n_plot) */
628  /* {AxisymmetricLinearElasticityEquations:: */
629  /* output(file_pt,n_plot);} */
630 
631  /* }; */
632 
633 
634  /* //============================================================================
635  */
636  /* /// FaceGeometry of a linear */
637  /* /// TAxisymmetricLinearElasticityElement element */
638  /* //============================================================================
639  */
640  /* template<unsigned NNODE_1D> */
641  /* class FaceGeometry<TAxisymmetricLinearElasticityElement<NNODE_1D> > : */
642  /* public virtual TElement<1,NNODE_1D> */
643  /* { */
644  /* public: */
645  /* /// Constructor must call the constructor of the underlying element */
646  /* FaceGeometry() : TElement<1,NNODE_1D>() {} */
647  /* }; */
648 
649 
650  /// /////////////////////////////////////////////////////////////////
651  /// /////////////////////////////////////////////////////////////////
652  /// /////////////////////////////////////////////////////////////////
653 
654 
655  //==========================================================
656  /// Axisym linear elasticity upgraded to become projectable
657  //==========================================================
658  template<class AXISYM_LINEAR_ELAST_ELEMENT>
660  : public virtual ProjectableElement<AXISYM_LINEAR_ELAST_ELEMENT>
661  {
662  public:
663  /// Constructor [this was only required explicitly
664  /// from gcc 4.5.2 onwards...]
666 
667 
668  /// Specify the values associated with field fld.
669  /// The information is returned in a vector of pairs which comprise
670  /// the Data object and the value within it, that correspond to field fld.
671  /// In the underlying linear elasticity elements the
672  /// the displacements are stored at the nodal values
674  {
675  // Create the vector
677 
678  // Loop over all nodes and extract the fld-th nodal value
679  unsigned nnod = this->nnode();
680  for (unsigned j = 0; j < nnod; j++)
681  {
682  // Add the data value associated with the displacement components
683  data_values.push_back(std::make_pair(this->node_pt(j), fld));
684  }
685 
686  // Return the vector
687  return data_values;
688  }
689 
690  /// Number of fields to be projected: 3, corresponding to
691  /// the displacement components
693  {
694  return 3;
695  }
696 
697  /// Number of history values to be stored for fld-th field.
698  /// (includes present value!)
699  unsigned nhistory_values_for_projection(const unsigned& fld)
700  {
701 #ifdef PARANOID
702  if (fld > 2)
703  {
704  std::stringstream error_stream;
705  error_stream << "Elements only store two fields so fld can't be"
706  << " " << fld << std::endl;
707  throw OomphLibError(
708  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
709  }
710 #endif
711  return this->node_pt(0)->ntstorage();
712  }
713 
714  /// Number of positional history values: Read out from
715  /// positional timestepper (Note: count includes current value!)
717  {
718  return this->node_pt(0)->position_time_stepper_pt()->ntstorage();
719  }
720 
721  /// Return Jacobian of mapping and shape functions of field fld
722  /// at local coordinate s
723  double jacobian_and_shape_of_field(const unsigned& fld,
724  const Vector<double>& s,
725  Shape& psi)
726  {
727  unsigned n_dim = this->dim();
728  unsigned n_node = this->nnode();
729  DShape dpsidx(n_node, n_dim);
730 
731  // Call the derivatives of the shape functions and return
732  // the Jacobian
733  return this->dshape_eulerian(s, psi, dpsidx);
734  }
735 
736 
737  /// Return interpolated field fld at local coordinate s, at time
738  /// level t (t=0: present; t>0: history values)
739  double get_field(const unsigned& t,
740  const unsigned& fld,
741  const Vector<double>& s)
742  {
743  unsigned n_node = this->nnode();
744 
745  // Local shape function
746  Shape psi(n_node);
747 
748  // Find values of shape function
749  this->shape(s, psi);
750 
751  // Initialise value of u
752  double interpolated_u = 0.0;
753 
754  // Sum over the local nodes at that time
755  for (unsigned l = 0; l < n_node; l++)
756  {
757  interpolated_u += this->nodal_value(t, l, fld) * psi[l];
758  }
759  return interpolated_u;
760  }
761 
762 
763  /// Return number of values in field fld
764  unsigned nvalue_of_field(const unsigned& fld)
765  {
766  return this->nnode();
767  }
768 
769 
770  /// Return local equation number of value j in field fld.
771  int local_equation(const unsigned& fld, const unsigned& j)
772  {
773  return this->nodal_local_eqn(j, fld);
774  }
775  };
776 
777 
778  //=======================================================================
779  /// Face geometry for element is the same as that for the underlying
780  /// wrapped element
781  //=======================================================================
782  template<class ELEMENT>
784  : public virtual FaceGeometry<ELEMENT>
785  {
786  public:
787  FaceGeometry() : FaceGeometry<ELEMENT>() {}
788  };
789 
790 
791  //=======================================================================
792  /// Face geometry of the Face Geometry for element is the same as
793  /// that for the underlying wrapped element
794  //=======================================================================
795  template<class ELEMENT>
798  : public virtual FaceGeometry<FaceGeometry<ELEMENT>>
799  {
800  public:
802  };
803 
804 
805 } // namespace oomph
806 
807 
808 #endif
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
char t
Definition: cfortran.h:568
A base class for elements that solve the axisymmetric (in cylindrical polars) equations of linear ela...
double d2u_dt2_axisymmetric_linear_elasticity(const unsigned &n, const unsigned &i) const
d^2u/dt^2 at local node n
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: for now lump ...
void interpolated_d2u_dt2_axisymmetric_linear_elasticity(const Vector< double > &s, Vector< double > &d2u_dt2) const
Compute vector of FE interpolated accel d2u/dt2 at local coordinate s.
double & nu() const
Access function for Poisson's ratio.
const double & lambda_sq() const
Access function for timescale ratio (nondim density)
virtual unsigned u_index_axisymmetric_linear_elasticity(const unsigned &i) const
Return the index at which the i-th (i=0: r, i=1: z; i=2: theta) unknown displacement component is sto...
double *& nu_pt()
Access function for pointer to Poisson's ratio.
double youngs_modulus() const
Access function to Young's modulus.
static double Default_lambda_sq_value
Static default value for timescale ratio (1.0 for natural scaling)
void interpolated_du_dt_axisymmetric_linear_elasticity(const Vector< double > &s, Vector< double > &du_dt) const
Compute vector of FE interpolated velocity du/dt at local coordinate s.
double du_dt_axisymmetric_linear_elasticity(const unsigned &n, const unsigned &i) const
du/dt at local node n
AxisymmetricLinearElasticityEquationsBase()
Constructor: Set null pointers for constitutive law. Set physical parameter values to default values,...
double *& youngs_modulus_pt()
Return the pointer to Young's modulus.
double *& lambda_sq_pt()
Access function for pointer to timescale ratio (nondim density)
static double Default_youngs_modulus_value
Static default value for Young's modulus (1.0 – for natural scaling, i.e. all stresses have been non-...
void body_force(const double &time, const Vector< double > &x, Vector< double > &b) const
Evaluate body force at Eulerian coordinate x at present time (returns zero vector if no body force fu...
BodyForceFctPt body_force_fct_pt() const
Access function: Pointer to body force function (const version)
BodyForceFctPt Body_force_fct_pt
Pointer to body force function.
BodyForceFctPt & body_force_fct_pt()
Access function: Pointer to body force function.
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
double interpolated_u_axisymmetric_linear_elasticity(const Vector< double > &s, const unsigned &i) const
Return FE interpolated displacement u[i] (i=0: r, i=1: z; i=2: theta) at local coordinate s.
void(* BodyForceFctPt)(const double &time, const Vector< double > &x, Vector< double > &b)
Function pointer to function that specifies the body force as a function of the Cartesian coordinates...
void interpolated_u_axisymmetric_linear_elasticity(const Vector< double > &s, Vector< double > &disp) const
Compute vector of FE interpolated displacement u at local coordinate s.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void get_strain(const Vector< double > &s, DenseMatrix< double > &strain)
Get strain (3x3 entries; r, z, phi)
virtual void fill_in_generic_contribution_to_residuals_axisymmetric_linear_elasticity(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Private helper function to compute residuals and (if requested via flag) also the Jacobian matrix.
void output(std::ostream &outfile)
Output: r,z, u_r, u_z, u_theta.
void output(FILE *file_pt)
C-style output: r,z, u_r, u_z, u_theta.
unsigned required_nvalue(const unsigned &n) const
Number of values required at node n.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
The jacobian is calculated by finite differences by default, We need only to take finite differences ...
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Validate against exact solution. Solution is provided via function pointer. Plot at a given number of...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals for the equations (the discretised principle of virtual displacements)
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact solution: r,z, u_r, u_z, u_theta.
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition: shape.h:278
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
unsigned ntstorage() const
Return total number of doubles stored per value to record time history of each value (one for steady ...
Definition: nodes.cc:879
FaceGeometry()
Constructor must call the constructor of the underlying element.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: elements.h:4998
A general Finite Element class.
Definition: elements.h:1313
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
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:2593
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...
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
Definition: elements.h:1432
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
Definition: elements.h:2611
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition: elements.h:1759
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:3298
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
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number.
Definition: elements.h:704
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
Definition: elements.h:227
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
Definition: geom_objects.h:192
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
Definition: nodes.h:1022
An OomphLibError object which should be thrown when an run-time error is encountered....
///////////////////////////////////////////////////////////////// ///////////////////////////////////...
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Return interpolated field fld at local coordinate s, at time level t (t=0: present; t>0: history valu...
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values: Read out from positional timestepper (Note: count includes curre...
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Specify the values associated with field fld. The information is returned in a vector of pairs which ...
ProjectableAxisymLinearElasticityElement()
Constructor [this was only required explicitly from gcc 4.5.2 onwards...].
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Return Jacobian of mapping and shape functions of field fld at local coordinate s.
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. (includes present value!)
unsigned nfields_for_projection()
Number of fields to be projected: 3, corresponding to the displacement components.
Wrapper class for projectable elements. Adds "projectability" to the underlying ELEMENT.
Definition: projection.h:183
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
void output(FILE *file_pt)
C-style output function.
void output(std::ostream &outfile)
Output function.
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
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
////////////////////////////////////////////////////////////////////// //////////////////////////////...
Definition: timesteppers.h:231
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
Definition: timesteppers.h:601
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
Definition: timesteppers.h:594
bool is_steady() const
Flag to indicate if a timestepper has been made steady (possibly temporarily to switch off time-depen...
Definition: timesteppers.h:389
//////////////////////////////////////////////////////////////////// ////////////////////////////////...