beam_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-2025 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 KL beam elements
27#ifndef OOMPH_KIRCHHOFF_LOVE_BEAM_ELEMENTS_HEADER
28#define OOMPH_KIRCHHOFF_LOVE_BEAM_ELEMENTS_HEADER
29
30// Config header
31#ifdef HAVE_CONFIG_H
32#include <oomph-lib-config.h>
33#endif
34
35// OOMPH-LIB header files
38#include "generic/fsi.h"
40
41namespace oomph
42{
43 //=======================================================================
44 /// A class for elements that solve the equations of Kirchhoff-Love
45 /// large-displacement (but linearly-elastic) thin-beam theory.
46 ///
47 /// The variational principle has the form
48 /// \f[ \int_0^{L} \left[ (\sigma_0 + \gamma) \ \delta \gamma + \frac{1}{12} \left(\frac{h}{R_0}\right)^2 \kappa \ \delta \kappa - \left( \left(\frac{R_0}{h}\right) {\bf f} - \Lambda^2 \frac{\partial^2 {\bf R}_w}{\partial t^2} \right) \cdot \delta {\bf R}_w \right] \ d\xi = 0, \f]
49 /// where all lengths have been non-dimensionalised w.r.t. \f$ R_0 \f$.
50 /// The strain and and bending "tensors" \f$\gamma\f$ and \f$\kappa\f$
51 /// are computed relative to the shape of the beam's undeformed shape
52 /// which is specified as a GeomObject.
53 ///
54 /// Time is scaled on the timescale \f$T\f$ and
55 /// \f[ \Lambda = \frac{a}{T} \sqrt{\frac{\rho}{E_{eff}}}, \f]
56 /// the ratio of the timescale used in the non-dimensionalisation of the
57 /// equations to the natural timescale of the wall oscillations (in the
58 /// wall's in-plane mode). \f$ \Lambda^2 \f$ can be interpreted as
59 /// the non-dimensional wall density, therefore \f$ \Lambda=0\f$
60 /// corresponds to the case without wall inertia.
61 ///
62 ///
63 /// Note that:
64 /// - the load vector \f$ {\bf f} \f$ is scaled
65 /// on the effective elastic modulus \f$ E_{eff}=E/(1-\nu^2)\f$
66 /// (rather than the
67 /// bending stiffness). Rescale the result yourself if you prefer
68 /// another non-dimensionalisation (the current version yields the
69 /// the most compact maths).
70 /// - Poisson's ratio does not appear explicitly since it only occurs
71 /// in combination with Young's modulus \f$E\f$.
72 ///
73 /// Default values:
74 /// - the 2nd Piola Kirchhoff pre-stress \f$ \sigma_0 \f$ is zero.
75 /// - the wall thickness \f$ h/R_0\f$ is 1/20.
76 /// - the timescale ratio \f$ \Lambda^2\f$ is 1.
77 /// - the traction vector \f$ f \f$ evaluates to zero.
78 ///
79 /// Need to specify:
80 /// - the undeformed wall shape (as a GeomObject).
81 ///
82 /// The governing equations can be switched from the principle of
83 /// virtual displacements to a system of equations that forces the
84 /// beam to deform into a shape specified by a SolidInitialCondition object.
85 /// If \c SolidFiniteElement::solid_ic_pt()!=0 we solve the
86 /// the equations
87 /// \f[ \int_0^{L} \left( \frac{\partial^i {\bf R}_{IC}}{\partial t^i} - {\bf R}_w \right) \psi_{jk} \ d\xi = 0, \f]
88 /// where \f$ \partial^i {\bf R}_{IC}/\partial t^i\f$ is
89 /// implemented by the SolidInitialCondition object, pointed to by
90 /// \c SolidFiniteElement::shell_ic_pt().
91 ///
92 //=======================================================================
94 {
95 private:
96 /// Static default value for 2nd Piola Kirchhoff prestress
97 static double Default_sigma0_value;
98
99 /// Static default value for timescale ratio (1.0 -- for natural scaling)
101
102 /// Static default value for non-dim wall thickness
103 // i.e. The reference value 'h_0'
104 static double Default_h_value;
105
106 /// Pointer to axial prestress
107 double* Sigma0_pt;
108
109 /// Pointer to wall thickness
110 // i.e. The reference value 'h_0'
111 double* H_pt;
112
113 /// Pointer to Timescale ratio (non-dim. density)
115
116 protected:
117 /// Default load function (zero traction)
118 static void Zero_traction_fct(const Vector<double>& xi,
119 const Vector<double>& x,
120 const Vector<double>& N,
122
123 /// Pointer to load vector function: Its arguments are:
124 /// Lagrangian coordinate, Eulerian coordinate, normal vector and
125 /// load vector itself (not all of the input arguments will be
126 /// required for all specific load functions but the list should
127 /// cover all cases)
129 const Vector<double>& x,
130 const Vector<double>& N,
132
133 /// Default profile function (constant thickness 'h_0')
134 static void Unit_profile_fct(const Vector<double>& xi,
135 const Vector<double>& x,
136 double& h_ratio);
137
138 /// Pointer to wall profile function: Its arguments are:
139 /// Lagrangian coordinate, Eulerian coordinate, and
140 /// profile itself (not all of the input arguments will be
141 /// required for all specific profile functions but the list should
142 /// cover all cases)
144 const Vector<double>& x,
145 double& h_ratio);
146
147 /// Pointer to the GeomObject that specifies the beam's
148 /// undeformed midplane
150
151 public:
152 /// Constructor. Set default values for all physical parameters
153 /// and zero traction.
155 {
156 // Set physical parameter pointers to the default values
159 // The reference thickness 'h_0'
161 // Zero traction
163 // Unit thickness profile
165 }
166
167
168 /// Reference to the load vector function pointer
170 const Vector<double>& x,
171 const Vector<double>& N,
173 {
174 return Load_vector_fct_pt;
175 }
176
177
178 /// Get the load vector: Pass number of integration point (dummy),
179 /// Lagr. and Eulerian coordinate and normal vector and return the load
180 /// vector (not all of the input arguments will be required for all specific
181 /// load functions but the list should cover all cases). This function is
182 /// virtual so it can be overloaded for FSI.
183 virtual void load_vector(const unsigned& intpt,
184 const Vector<double>& xi,
185 const Vector<double>& x,
186 const Vector<double>& N,
188 {
189 Load_vector_fct_pt(xi, x, N, load);
190 }
191
192 /// Reference to the wall thickness ratio profile function pointer
194 const Vector<double>& x,
195 double& h_ratio)
196 {
197 return Wall_profile_fct_pt;
198 }
199
200
201 /// Get the wall profile: Pass Lagrangian & Eulerian coordinate
202 /// and return the wall profile (not all of the input arguments will be
203 /// required for all specific thickness functions but the list should cover
204 /// all cases).
206 const Vector<double>& x,
207 double& h_ratio)
208 {
210 }
211
212
213 /// Return the non-dimensional wall thickness
214 // i.e. the reference value 'h_0'
215 const double& h() const
216 {
217 return *H_pt;
218 }
219
220 /// Return the timescale ratio (non-dimensional density)
221 const double& lambda_sq() const
222 {
223 return *Lambda_sq_pt;
224 }
225
226 /// Return the axial prestress
227 const double& sigma0() const
228 {
229 return *Sigma0_pt;
230 }
231
232 /// Return a pointer to axial prestress
233 double*& sigma0_pt()
234 {
235 return Sigma0_pt;
236 }
237
238 /// Return a pointer to non-dim. wall thickness
239 // i.e. the reference value 'h_0'
240 double*& h_pt()
241 {
242 return H_pt;
243 }
244
245 /// Return a pointer to timescale ratio (nondim density)
246 double*& lambda_sq_pt()
247 {
248 return Lambda_sq_pt;
249 }
250
251 /// Return a Pointer to geometric object that specifies the beam's
252 /// undeformed geometry
257
258 /// Get normal vector on wall
260 {
261 Vector<double> r(2);
262 get_normal(s, r, N);
263 }
264
265
266 /// Get position vector to and normal vector on wall
267 void get_normal(const Vector<double>& s,
269 Vector<double>& N);
270
271 /// Get position vector to and non-unit tangent vector on wall:
272 /// dr/ds
276
277 /// Return the residuals for the equations of Kirchhoff-Love beam
278 /// theory with linear constitutive equations; if Solid_ic_pt!=0, we
279 /// assign residuals which force the assignement of an initial shape/
280 /// veloc/accel to the dofs. This overloads the standard interface.
285
286
287 /// Return the residuals for the equations of Kirchhoff-Love beam
288 /// theory with linear constitutive equations; if Solid_ic_pt!=0, we
289 /// assign residuals which force the assignement of an initial shape/
290 /// veloc/accel to the dofs.
292
293
294 /// Get FE jacobian and residuals (Jacobian done by finite differences)
297
298 /// Get potential (strain) and kinetic energy of the element
299 void get_energy(double& pot_en, double& kin_en);
300
301 /// Get the potential energy due to stretching and bending and the
302 /// kinetic energy of the element
303 void get_energy(double& stretch, double& bend, double& kin_en);
304 };
305
306
307 //=========================================================================
308 /// Hermite Kirchhoff Love beam. Implements KirchhoffLoveBeamEquations
309 /// using 2-node Hermite elements as the underlying geometrical elements.
310 //=========================================================================
311 class HermiteBeamElement : public virtual SolidQHermiteElement<1>,
313 {
314 public:
315 /// Constructor (empty)
318 {
319 // Set the number of dimensions at each node (2D node on 1D surface)
321 }
322
323 /// Output function
324 void output(std::ostream& outfile);
325
326 /// Output function with specified number of plot points
327 void output(std::ostream& outfile, const unsigned& n_plot);
328
329 /// Output at previous time (t=0: present; t>0: previous)
330 /// with specified number of plot points
331 void output(const unsigned& t,
332 std::ostream& outfile,
333 const unsigned& n_plot) const;
334
335 /// C-style output function
336 void output(FILE* file_pt);
337
338 /// C-style output function with specified number of plot points
339 void output(FILE* file_pt, const unsigned& n_plot);
340
341 /// C-style output at previous time (t=0: present; t>0: previous)
342 /// with specified number of plot points
343 void output(const unsigned& t, FILE* file_pt, const unsigned& n_plot) const;
344 };
345
346 //=========================================================================
347 /// Hermite Kirchhoff Love beam "upgraded" to a FSIWallElement (and thus,
348 /// by inheritance, a GeomObject), so it can be used in FSI.
349 //=========================================================================
351 public virtual FSIWallElement
352 {
353 private:
354 // Boolean flag to indicate whether the normal is directed into the fluid
356
357 public:
358 /// Constructor: Create beam element as FSIWallElement (and thus,
359 /// by inheritance, a GeomObject). By default, we assume that the
360 /// normal vector computed by KirchhoffLoveBeamEquations::get_normal(...)
361 /// points into the fluid. If this is not the case, overwrite this
362 /// with the access function
363 /// FSIHermiteBeamElement::set_normal_pointing_out_of_fluid()
366 {
367 unsigned n_lagr = 1;
368 unsigned n_dim = 2;
370 }
371
372 /// Destructor: empty
374
375 /// Set the normal computed by
376 /// KirchhoffLoveBeamEquations::get_normal(...) to point into the fluid
381
382 /// Set the normal computed by
383 /// KirchhoffLoveBeamEquations::get_normal(...) to point out of the fluid
388
389
390 /// Derivative of position vector w.r.t. the SolidFiniteElement's
391 /// Lagrangian coordinates; evaluated at current time.
394
395 /// Get the load vector: Pass number of the integration point,
396 /// Lagr. coordinate, Eulerian coordinate and normal vector
397 /// and return the load vector. (Not all of the input arguments will be
398 /// required for all specific load functions but the list should
399 /// cover all cases). We first evaluate the load function defined via
400 /// KirchhoffLoveBeamEquations::load_vector_fct_pt() -- this
401 /// represents the non-FSI load on the beam, e.g. an external
402 /// pressure load. Then we add to this the FSI load due to
403 /// the traction exerted by the adjacent FSIFluidElements, taking
404 /// the sign of the normal into account.
405 void load_vector(const unsigned& intpt,
406 const Vector<double>& xi,
407 const Vector<double>& x,
408 const Vector<double>& N,
410 {
411 // Initially call the standard Load_vector_fct_pt
412 Load_vector_fct_pt(xi, x, N, load);
413
414 // Memory for the FSI load
416
417 // Get the fluid load on the wall stress scale
419
420 // If the normal is outer to the fluid switch the direction
421 double sign = 1.0;
423 {
424 sign = -1.0;
425 }
426
427 // Add the FSI load to the load vector
428 for (unsigned i = 0; i < 2; i++)
429 {
430 load[i] += sign * fsi_load[i];
431 }
432 }
433
434 /// Get the Jacobian and residuals. Wrapper to generic FSI version;
435 /// that catches the case when we replace the Jacobian by the
436 /// mass matrix (for the consistent assignment of initial conditions).
438 DenseMatrix<double>& jacobian)
439 {
440 // Call the standard beam element's jacobian function
442 // Now add the external interaction data by finite differences
444 }
445
446 /// Find the local coordinate s in this element
447 /// that corresponds to the global "intrinsic" coordinate \f$ \zeta \f$
448 /// (here identical to the Lagrangian coordinate \f$ \xi \f$).
449 /// If the coordinate is contained within this element, the
450 /// geom_object_pt points to "this" element; if the zeta coordinate
451 /// is not contained in this element geom_object_pt=NULL.
452 /// By default don't use any value passed in to the local coordinate s
453 /// as the initial guess in the Newton method
454 void locate_zeta(const Vector<double>& zeta,
455 GeomObject*& geom_object_pt,
457 const bool& use_coordinate_as_initial_guess = false);
458
459
460 /// The number of "DOF types" that degrees of freedom in this element
461 /// are sub-divided into: Just the solid degrees of freedom themselves.
462 unsigned ndof_types() const
463 {
464 return 1;
465 }
466
467 /// Create a list of pairs for all unknowns in this element,
468 /// so that the first entry in each pair contains the global equation
469 /// number of the unknown, while the second one contains the number
470 /// of the "DOF type" that this unknown is associated with.
471 /// (Function can obviously only be called if the equation numbering
472 /// scheme has been set up.)
474 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list) const;
475 };
476
477
478 ///////////////////////////////////////////////////////////////////////
479 ///////////////////////////////////////////////////////////////////////
480 ///////////////////////////////////////////////////////////////////////
481
482
483 //=======================================================================
484 /// Face geometry for the HermiteBeam elements: Solid point element
485 //=======================================================================
486 template<>
488 {
489 public:
490 /// Constructor [this was only required explicitly
491 /// from gcc 4.5.2 onwards...]
493 };
494
495
496 ///////////////////////////////////////////////////////////////////////
497 ///////////////////////////////////////////////////////////////////////
498 ///////////////////////////////////////////////////////////////////////
499
500
501 //======================================================================
502 /// Element that allows the imposition of boundary
503 /// conditions for a beam that is clamped but can slide
504 /// along a line which is specified by a position vector
505 /// to that line and the normal vector to it. The endpoint
506 /// of the beam is forced to stay on that line and meet
507 /// it at a right angle. This is achieved with Lagrange multipliers.
508 //======================================================================
510 : public virtual FaceGeometry<HermiteBeamElement>,
511 public virtual SolidFaceElement
512 {
513 public:
514 /// Constructor, takes the pointer to the "bulk" element, the
515 /// index of the fixed local coordinate and its value represented
516 /// by an integer (+/- 1), indicating that the face is located
517 /// at the max. or min. value of the "fixed" local coordinate
518 /// in the bulk element.
520 FiniteElement* const& bulk_el_pt, const int& face_index);
521
522 /// Broken empty constructor
524 {
525 throw OomphLibError("Don't call empty constructor for "
526 "ClampedSlidingHermiteBeamBoundaryConditionElement ",
529 }
530
531
532 /// Broken copy constructor
535
536 /// Broken assignment operator
537 // Commented out broken assignment operator because this can lead to a
538 // conflict warning when used in the virtual inheritence hierarchy.
539 // Essentially the compiler doesn't realise that two separate
540 // implementations of the broken function are the same and so, quite
541 // rightly, it shouts.
542 /*void operator=(const ClampedSlidingHermiteBeamBoundaryConditionElement&) =
543 delete;*/
544
545
546 /// Set vectors to some point on the symmetry line, and
547 /// normal to that line along which the end of the beam is sliding.
556
557
558 /// Fill in the element's contribution to its residual vector
560
561
562 /// Output function -- forward to broken version in FiniteElement
563 /// until somebody decides what exactly they want to plot here...
564 void output(std::ostream& outfile)
565 {
567 }
568
569 /// Output function -- forward to broken version in FiniteElement
570 /// until somebody decides what exactly they want to plot here...
571 void output(std::ostream& outfile, const unsigned& n_plot)
572 {
574 }
575
576 /// C-style output function -- forward to broken version in FiniteElement
577 /// until somebody decides what exactly they want to plot here...
582
583 /// C-style output function -- forward to broken version in
584 /// FiniteElement until somebody decides what exactly they want to plot
585 /// here...
586 void output(FILE* file_pt, const unsigned& n_plot)
587 {
589 }
590
591 private:
592 /// Vector to some point on the symmetry line along which the
593 /// end of the beam is sliding
595
596 /// Normal vector to the symmetry line along which the
597 /// end of the beam is sliding
599 };
600
601
602} // namespace oomph
603
604#endif
static char t char * s
Definition cfortran.h:568
cstr elem_len * i
Definition cfortran.h:603
char t
Definition cfortran.h:568
Element that allows the imposition of boundary conditions for a beam that is clamped but can slide al...
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
Vector< double > Vector_to_symmetry_line
Vector to some point on the symmetry line along which the end of the beam is sliding.
void set_symmetry_line(const Vector< double > &vector_to_symmetry_line, const Vector< double > &normal_to_symmetry_line)
Broken assignment operator.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
ClampedSlidingHermiteBeamBoundaryConditionElement(const ClampedSlidingHermiteBeamBoundaryConditionElement &dummy)=delete
Broken copy constructor.
void output(FILE *file_pt)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
void output(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
Vector< double > Normal_to_symmetry_line
Normal vector to the symmetry line along which the end of the beam is sliding.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the element's contribution to its residual vector.
ClampedSlidingHermiteBeamBoundaryConditionElement()
Broken empty constructor.
void fill_in_jacobian_from_external_interaction_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the contributions to the jacobian from all external interaction degrees of freedom (geometr...
Hermite Kirchhoff Love beam "upgraded" to a FSIWallElement (and thus, by inheritance,...
void dposition_dlagrangian_at_local_coordinate(const Vector< double > &s, DenseMatrix< double > &drdxi) const
Derivative of position vector w.r.t. the SolidFiniteElement's Lagrangian coordinates; evaluated at cu...
FSIHermiteBeamElement()
Constructor: Create beam element as FSIWallElement (and thus, by inheritance, a GeomObject)....
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: Just the soli...
void set_normal_pointing_into_fluid()
Set the normal computed by KirchhoffLoveBeamEquations::get_normal(...) to point into the fluid.
~FSIHermiteBeamElement()
Destructor: empty.
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...
virtual void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Get the Jacobian and residuals. Wrapper to generic FSI version; that catches the case when we replace...
void set_normal_pointing_out_of_fluid()
Set the normal computed by KirchhoffLoveBeamEquations::get_normal(...) to point out of the fluid.
void locate_zeta(const Vector< double > &zeta, GeomObject *&geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
Find the local coordinate s in this element that corresponds to the global "intrinsic" coordinate (h...
void load_vector(const unsigned &intpt, const Vector< double > &xi, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Get the load vector: Pass number of the integration point, Lagr. coordinate, Eulerian coordinate and ...
This is a base class for all SolidFiniteElements that participate in FSI computations....
Definition fsi.h:194
void setup_fsi_wall_element(const unsigned &nlagr_solid, const unsigned &ndim_fluid)
Setup: Assign storage – pass the Eulerian dimension of the "adjacent" fluid elements and the number o...
Definition fsi.cc:121
void fluid_load_vector(const unsigned &intpt, const Vector< double > &N, Vector< double > &load)
Get FE Jacobian by systematic finite differencing w.r.t. nodal positition Data, internal and external...
Definition fsi.cc:162
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
Definition elements.h:4630
FaceGeometry()
Constructor [this was only required explicitly from gcc 4.5.2 onwards...].
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
Definition elements.h:5002
A general Finite Element class.
Definition elements.h:1317
void set_nodal_dimension(const unsigned &nodal_dim)
Set the dimension of the nodes in the element. This will typically only be required when constructing...
Definition elements.h:1394
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
Definition elements.h:3054
A geometric object is an object that provides a parametrised description of its shape via the functio...
Hermite Kirchhoff Love beam. Implements KirchhoffLoveBeamEquations using 2-node Hermite elements as t...
void output(std::ostream &outfile)
Output function.
HermiteBeamElement()
Constructor (empty)
A class for elements that solve the equations of Kirchhoff-Love large-displacement (but linearly-elas...
KirchhoffLoveBeamEquations()
Constructor. Set default values for all physical parameters and zero traction.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals for the equations of Kirchhoff-Love beam theory with linear constitutive equatio...
static double Default_sigma0_value
Static default value for 2nd Piola Kirchhoff prestress.
static void Zero_traction_fct(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Default load function (zero traction)
GeomObject * Undeformed_beam_pt
Pointer to the GeomObject that specifies the beam's undeformed midplane.
const double & h() const
Return the non-dimensional wall thickness.
double *& h_pt()
Return a pointer to non-dim. wall thickness.
static void Unit_profile_fct(const Vector< double > &xi, const Vector< double > &x, double &h_ratio)
Default profile function (constant thickness 'h_0')
void fill_in_contribution_to_residuals_beam(Vector< double > &residuals)
Return the residuals for the equations of Kirchhoff-Love beam theory with linear constitutive equatio...
void(*&)(const Vector< double > &xi, const Vector< double > &x, double &h_ratio) wall_profile_fct_pt()
Reference to the wall thickness ratio profile function pointer.
void(* Load_vector_fct_pt)(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Pointer to load vector function: Its arguments are: Lagrangian coordinate, Eulerian coordinate,...
void(* Wall_profile_fct_pt)(const Vector< double > &xi, const Vector< double > &x, double &h_ratio)
Pointer to wall profile function: Its arguments are: Lagrangian coordinate, Eulerian coordinate,...
double * Lambda_sq_pt
Pointer to Timescale ratio (non-dim. density)
void get_non_unit_tangent(const Vector< double > &s, Vector< double > &r, Vector< double > &drds)
Get position vector to and non-unit tangent vector on wall: dr/ds.
const double & lambda_sq() const
Return the timescale ratio (non-dimensional density)
static double Default_h_value
Static default value for non-dim wall thickness.
double * Sigma0_pt
Pointer to axial prestress.
void get_normal(const Vector< double > &s, Vector< double > &N)
Get normal vector on wall.
static double Default_lambda_sq_value
Static default value for timescale ratio (1.0 – for natural scaling)
void wall_profile(const Vector< double > &xi, const Vector< double > &x, double &h_ratio)
Get the wall profile: Pass Lagrangian & Eulerian coordinate and return the wall profile (not all of t...
const double & sigma0() const
Return the axial prestress.
void get_energy(double &pot_en, double &kin_en)
Get potential (strain) and kinetic energy of the element.
double *& sigma0_pt()
Return a pointer to axial prestress.
double * H_pt
Pointer to wall thickness.
void(*&)(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &N, Vector< double > &load) load_vector_fct_pt()
Reference to the load vector function pointer.
virtual void load_vector(const unsigned &intpt, const Vector< double > &xi, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Get the load vector: Pass number of integration point (dummy), Lagr. and Eulerian coordinate and norm...
GeomObject *& undeformed_beam_pt()
Return a Pointer to geometric object that specifies the beam's undeformed geometry.
virtual void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Get FE jacobian and residuals (Jacobian done by finite differences)
double *& lambda_sq_pt()
Return a pointer to timescale ratio (nondim density)
An OomphLibError object which should be thrown when an run-time error is encountered....
SolidFaceElements combine FaceElements and SolidFiniteElements and overload various functions so they...
Definition elements.h:4918
SolidFiniteElement class.
Definition elements.h:3565
Solid point element.
Definition elements.h:4984
SolidQHermiteElement elements are Hermite elements whose Jacobian matrices include derivatives w....
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).