shell_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 shell elements
27#ifndef OOMPH_KIRCHHOFF_LOVE_SHELL_ELEMENTS_HEADER
28#define OOMPH_KIRCHHOFF_LOVE_SHELL_ELEMENTS_HEADER
29
30// Config header
31#ifdef HAVE_CONFIG_H
32#include <oomph-lib-config.h>
33#endif
34
35// OOMPH-LIB header
38#include "generic/fsi.h"
40#include "generic/matrices.h"
41
42namespace oomph
43{
44 //========================================================================
45 /// A class for elements that solves the equations of Kirchhoff Love shell
46 /// thin-shell theory
47 //========================================================================
49 {
50 private:
51 /// Static default value for the Poisson's ratio
52 static double Default_nu_value;
53
54 /// Static default value for the timescale ratio (1.0 for natural scaling)
56
57 /// Static default value for the thickness ratio
58 static double Default_h_value;
59
60 /// Boolean flag to ignore membrane terms
62
63 /// Pointer to Poisson's ratio
64 double* Nu_pt;
65
66 /// Pointer to dimensionless wall thickness
67 double* H_pt;
68
69 /// Pointer to timescale ratio (non-dimensional density)
70 double* Lambda_sq_pt;
71
72 /// Pointer to membrane pre-stress terms -- should
73 /// probably generalise this to function pointers at some point
75
76 /// Static default for prestress (set to zero)
77 static double Zero_prestress;
78
79 protected:
80 /// Invert a DIM by DIM matrix
81 inline double calculate_contravariant(double A[2][2], double Aup[2][2]);
82
83 /// Default load function (zero traction)
84 static void Zero_traction_fct(const Vector<double>& xi,
85 const Vector<double>& x,
86 const Vector<double>& N,
88
89 /// Pointer to load vector function: Its arguments are:
90 /// Lagrangian coordinate, Eulerian coordinate, normal vector and
91 /// load vector itself (not all of the input arguments will be
92 /// required for all specific load functions but the list should
93 /// cover all cases)
95 const Vector<double>& x,
96 const Vector<double>& N,
98
99
100 /// Pointer to the GeomObject that specifies the undeformed midplane
101 /// of the shell
103
104
105 /// Helper function to Return the residuals for
106 /// the equations of KL shell theory. This is used to prevent
107 /// a call of fill_in_contribution_to_residuals in
108 /// the function fill_in_contribution_to_jacobian that can
109 /// lead to virtual inheritance woes if this element is ever
110 /// used as part of a multi-physics element.
112
113 /// Get the load vector for the computation of the rate of work
114 /// done by the load. Here we simply forward this to
115 /// load_vector(...) but allow it to be overloaded in derived classes
116 /// to allow the computation of the rate of work due to constituent
117 /// parts of the load vector (e.g. the fluid traction in an FSI
118 /// problem). Pass number of integration point (dummy),
119 /// Lagr. coordinate and normal vector and return the load vector
120 /// (not all of the input arguments will be
121 /// required for all specific load functions but the list should
122 /// cover all cases).
124 const unsigned& intpt,
125 const Vector<double>& xi,
126 const Vector<double>& x,
127 const Vector<double>& N,
129 {
130 load_vector(intpt, xi, x, N, load);
131 }
132
133
134 public:
135 /// Constructor: Initialise physical parameter values to defaults
137 {
138 // Set the default values of physical parameters
142
143 // Don't ignore membrane terms
144 Ignore_membrane_terms = false;
145
146 // Default load is zero traction
148
149 // Default prestress is zero
150 Prestress_pt.resize(2, 2);
155 }
156
157 /// Access to the load vector function pointer
159 const Vector<double>& x,
160 const Vector<double>& N,
162 {
163 return Load_vector_fct_pt;
164 }
165
166 /// Get the load vector: Pass number of integration point (dummy),
167 /// Lagr. coordinate and normal vector and return the load vector
168 /// (not all of the input arguments will be
169 /// required for all specific load functions but the list should
170 /// cover all cases). This function is virtual so it can be
171 /// overloaded for FSI.
172 virtual void load_vector(const unsigned& intpt,
173 const Vector<double>& xi,
174 const Vector<double>& x,
175 const Vector<double>& N,
177 {
178 Load_vector_fct_pt(xi, x, N, load);
179 }
180
181
182 /// Set pointer to (i,j)-th component of second Piola Kirchhoff
183 /// membrane prestress to specified value (automatically imposes
184 /// symmetry for off-diagonal entries)
185 void set_prestress_pt(const unsigned& i,
186 const unsigned& j,
187 double* value_pt)
188 {
189 Prestress_pt(i, j) = value_pt;
190 Prestress_pt(j, i) = value_pt;
191 }
192
193 /// Return (i,j)-th component of second Piola Kirchhoff membrane
194 /// prestress
195 double prestress(const unsigned& i, const unsigned& j)
196 {
197 return *Prestress_pt(i, j);
198 }
199
200 /// Set to disable the calculation of membrane terms
202 {
204 }
205
206 /// Set to renable the calculation of membrane terms (default)
208 {
209 Ignore_membrane_terms = false;
210 }
211
212 /// Return the Poisson's ratio
213 const double& nu() const
214 {
215 return *Nu_pt;
216 }
217
218 /// Return the wall thickness to undeformed radius ratio
219 const double& h() const
220 {
221 return *H_pt;
222 }
223
224 /// Return the timescale ratio (non-dimensional density)
225 const double& lambda_sq() const
226 {
227 return *Lambda_sq_pt;
228 }
229
230 /// Return a pointer to the Poisson ratio
231 double*& nu_pt()
232 {
233 return Nu_pt;
234 }
235
236 /// Return a pointer to the non-dim wall thickness
237 double*& h_pt()
238 {
239 return H_pt;
240 }
241
242 /// Return a pointer to timescale ratio (nondim density)
243 double*& lambda_sq_pt()
244 {
245 return Lambda_sq_pt;
246 }
247
248 /// Return a reference to geometric object that specifies the shell's
249 /// undeformed geometry
254
255 /// Get normal vector
256 void get_normal(const Vector<double>& s, Vector<double>& N);
257
258 /// Overload the standard fill in residuals contribution
264
265 /// Return the jacobian is calculated by finite differences by default,
267 DenseMatrix<double>& jacobian);
268
269 /// Get potential (strain) and kinetic energy of the element
270 void get_energy(double& pot_en, double& kin_en);
271
272
273 /// Get strain and bending tensors; returns pair comprising the
274 /// determinant of the undeformed (*.first) and deformed (*.second)
275 /// midplane metric tensor.
276 std::pair<double, double> get_strain_and_bend(const Vector<double>& s,
279
280
281 /// Get integral of instantaneous rate of work done on
282 /// the wall due to the load returned by the virtual
283 /// function load_vector_for_rate_of_work_computation(...).
284 /// In the current class
285 /// the latter function simply de-references the external
286 /// load but this can be overloaded in derived classes
287 /// (e.g. in FSI elements) to determine the rate of work done
288 /// by individual constituents of this load (e.g. the fluid load
289 /// in an FSI problem).
290 double load_rate_of_work();
291
292 /// Generic FiniteElement output function
293 void output(std::ostream& outfile)
294 {
296 }
297
298 /// Generic FiniteElement output function
299 void output(std::ostream& outfile, const unsigned& n_plot)
300 {
302 }
303
304 /// Generic FiniteElement output function
309
310 /// Generic FiniteElement output function
311 void output(FILE* file_pt, const unsigned& n_plot)
312 {
314 }
315 };
316
317
318 //==================================================================
319 /// Matrix inversion for 2 dimensions
320 //==================================================================
322 double Aup[2][2])
323 {
324 // Calculate determinant
325 double det = A[0][0] * A[1][1] - A[0][1] * A[1][0];
326
327 // Calculate entries of the inverse
328 Aup[0][0] = A[1][1] / det;
329 Aup[0][1] = -A[0][1] / det;
330 Aup[1][0] = -A[1][0] / det;
331 Aup[1][1] = A[0][0] / det;
332
333 // Return determinant
334 return (det);
335 }
336
337
338 //=======================================================================
339 /// An element that solves the Kirchhoff-Love shell theory equations
340 /// using Hermite interpolation (displacements
341 /// and slopes are interpolated separately. The local and global
342 /// (Lagrangian) coordinates are not assumed to be aligned.
343 /// N.B. It will be DOG SLOW.
344 //=======================================================================
345 class HermiteShellElement : public virtual SolidQHermiteElement<2>,
347 {
348 public:
349 /// Constructor, there are no internal data points
352 {
353 // Set the number of dimensions at each node (3D nodes on 2D surface)
355 }
356
357 /// Output position veloc and accel.
358 void output_with_time_dep_quantities(std::ostream& outfile,
359 const unsigned& n_plot);
360
361 /// Overload the output function
362 void output(std::ostream& outfile)
363 {
365 }
366
367 /// Output function
368 void output(std::ostream& outfile, const unsigned& n_plot);
369
370 /// Overload the output function
375
376 /// Output function
377 void output(FILE* file_pt, const unsigned& n_plot);
378 };
379
380
381 //=========================================================================
382 /// An element that solves the Kirchhoff-Love shell theory equations
383 /// using Hermite interpolation (displacements
384 /// and slopes are interpolated separately. The local and global
385 /// (Lagrangian) coordinates are assumed to be aligned so that the
386 /// Jacobian of the mapping between these coordinates is diagonal.
387 /// This significantly simplifies (and speeds up) the computation of the
388 /// derivatives of the shape functions.
389 //=========================================================================
392 {
393 public:
394 /// Constructor, there are no internal data points
399 };
400
401
402 ////////////////////////////////////////////////////////////////////////
403 ////////////////////////////////////////////////////////////////////////
404 ////////////////////////////////////////////////////////////////////////
405
406
407 //=======================================================================
408 /// Face geometry for the HermiteShell elements: 1D SolidQHermiteElement
409 //=======================================================================
410 template<>
412 : public virtual SolidQHermiteElement<1>
413 {
414 public:
415 /// Constructor [this was only required explicitly
416 /// from gcc 4.5.2 onwards...]
418 };
419
420
421 ////////////////////////////////////////////////////////////////////////
422 ////////////////////////////////////////////////////////////////////////
423 ////////////////////////////////////////////////////////////////////////
424
425
426 //=========================================================================
427 /// Diag Hermite Kirchhoff Love shell "upgraded" to a FSIWallElement
428 /// (and thus, by inheritance, a GeomObject), so it can be used in FSI.
429 //=========================================================================
431 public virtual FSIWallElement
432 {
433 private:
434 /// Get the load vector for the computation of the rate of work
435 /// done by the load. Can switch between full load and fluid load only.
436 /// Overloads the version in the shell element base class.
437 /// Pass number of integration point (dummy)
438 /// Lagr. coordinate and normal vector and return the load vector
439 /// (not all of the input arguments will be
440 /// required for all specific load functions but the list should
441 /// cover all cases).
443 const unsigned& intpt,
444 const Vector<double>& xi,
445 const Vector<double>& x,
446 const Vector<double>& N,
448 {
449 /// Get fluid-only load vector
451 {
453 }
454 // Get full load vector as per default
455 else
456 {
457 load_vector(intpt, xi, x, N, load);
458 }
459 }
460
461 /// Boolean flag to indicate whether the normal is directed into the fluid
463
464 /// Boolean flag to indicate if rate-of-work by load is to be
465 /// based on the fluid traction only
467
468 public:
469 /// Constructor: Create shell element as FSIWallElement (and thus,
470 /// by inheritance, a GeomObject) with two Lagrangian coordinates
471 /// and 3 Eulerian coordinates. By default, we assume that the
472 /// normal vector computed by KirchhoffLoveShellEquations::get_normal(...)
473 /// points into the fluid.
474 /// If this is not the case, call the access function
475 /// FSIDiagHermiteShellElement::set_normal_pointing_out_of_fluid()
485
486 /// Destructor: empty
488
489 /// Set the normal computed by
490 /// KirchhoffLoveShellEquations::get_normal(...) to point into the fluid
495
496 /// Set the normal computed by
497 /// KirchhoffLoveShellEquations::get_normal(...) to point out of the fluid
502
503
504 /// Derivative of position vector w.r.t. the SolidFiniteElement's
505 /// Lagrangian coordinates; evaluated at current time.
508
509
510 /// Get integral of instantaneous rate of work done on
511 /// the wall due to the fluid load returned by the function
512 /// fluid_load_vector(...).
520
521
522 /// Get the load vector: Pass number of the integration point,
523 /// Lagr. coordinate, Eulerian coordinate and normal vector
524 /// and return the load vector. (Not all of the input arguments will be
525 /// required for all specific load functions but the list should
526 /// cover all cases). We first evaluate the load function defined via
527 /// KirchhoffLoveShellEquations::load_vector_fct_pt() -- this
528 /// represents the non-FSI load on the shell, e.g. an external
529 /// pressure load. Then we add to this the FSI load due to
530 /// the traction exerted by the adjacent FSIFluidElements, taking
531 /// the sign of the normal into account.
532 void load_vector(const unsigned& intpt,
533 const Vector<double>& xi,
534 const Vector<double>& x,
535 const Vector<double>& N,
537 {
538 // Initially call the standard Load_vector_fct_pt
539 Load_vector_fct_pt(xi, x, N, load);
540
541 // Memory for the FSI load
543
544 // Get the fluid load on the wall stress scale
546
547 // If the normal is outer to the fluid switch the direction
548 double sign = 1.0;
550 {
551 sign = -1.0;
552 }
553
554 // Add the FSI load to the load vector
555 for (unsigned i = 0; i < 3; i++)
556 {
557 load[i] += sign * fsi_load[i];
558 }
559 }
560
561 /// Get the Jacobian and residuals. Wrapper to generic FSI version;
562 /// that catches the case when we replace the Jacobian by the
563 /// mass matrix (for the consistent assignment of initial conditions).
565 DenseMatrix<double>& jacobian)
566 {
567 // Call the basic shell jacobian
569 jacobian);
570 // Fill in the external interaction by finite differences
572 }
573
574
575 /// The number of "DOF types" that degrees of freedom in this element
576 /// are sub-divided into: Just the solid degrees of freedom themselves.
577 unsigned ndof_types() const
578 {
579 return 1;
580 }
581
582 /// Create a list of pairs for all unknowns in this element,
583 /// so that the first entry in each pair contains the global equation
584 /// number of the unknown, while the second one contains the number
585 /// of the "DOF types" that this unknown is associated with.
586 /// (Function can obviously only be called if the equation numbering
587 /// scheme has been set up.)
589 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list) const;
590 };
591
592
593 ////////////////////////////////////////////////////////////////////////
594 ////////////////////////////////////////////////////////////////////////
595 ////////////////////////////////////////////////////////////////////////
596
597
598 //======================================================================
599 /// Element that allows the imposition of boundary
600 /// conditions for a shell that is clamped to a 2D plane
601 /// that is specified by its normal. Constraint is applied by
602 /// a Lagrange multiplier.
603 /// \n
604 /// \b Note \b 1: Note that the introduction of the Lagrange
605 /// multiplier adds two additional values (relative to the number
606 /// of values before the addition of the FaceElement) to the nodes.
607 /// This ensures that nodes that are shared by adjacent FaceElements
608 /// are not resized repeatedly but also means that this won't work
609 /// if two "edges" of the shell (that share a node) are subject
610 /// to different constraints, each applied with its own
611 /// independent Lagrange multiplier. In such cases a modified
612 /// version of this class must be written.
613 /// \n
614 /// \b Note \b 2: The FaceGeometry for a HermiteShellElement is
615 /// the 1D two-node element SolidQHermiteElement<1> which has
616 /// four shape functions (two nodes, two types -- representing the
617 /// shape functions that interpolate the value and the derivative).
618 /// These are the "correct" shape functions for the interpolation
619 /// of the Lagrange multiplier and the isoparametric representation
620 /// of the geometry. However, when applying the contribution from the
621 /// constraint equation to the bulk equations, we have to take
622 /// all four types of dof into account so the element has to
623 /// reset the number of positional dofs to four. To avoid any clashes
624 /// we overload (the relevant subset of) the access functions to the
625 /// shape functions and their derivatives and set the shape functions
626 /// associated with the spurious positional dofs to zero. This is a bit
627 /// hacky but the only way (?) this can be done...
628 //======================================================================
630 : public virtual FaceGeometry<HermiteShellElement>,
631 public virtual SolidFaceElement
632 {
633 public:
634 /// Constructor, takes the pointer to the "bulk" element and the
635 /// face index.
637 FiniteElement* const& bulk_el_pt, const int& face_index);
638
639
640 /// Broken empty constructor
642 {
643 throw OomphLibError("Don't call empty constructor for "
644 "ClampedHermiteShellBoundaryConditionElement",
647 }
648
649 /// Broken copy constructor
652
653 /// Broken assignment operator
655
656 /// Set normal vector to clamping plane
663
664 /// Fill in the element's contribution to its residual vector
666
667
668 //////////////////////////////////////////////////////////////////
669 // Note: We should also overload all other versions of shape(...)
670 // and dshape(...) but these are the only ones that are used.
671 //////////////////////////////////////////////////////////////////
672
673
674 /// Calculate the geometric shape functions
675 /// at local coordinate s. Set any "superfluous" shape functions to zero.
676 void shape(const Vector<double>& s, Shape& psi) const
677 {
678 // Initialise all of them to zero
679 unsigned n = psi.nindex1();
680 unsigned m = psi.nindex2();
681 for (unsigned i = 0; i < n; i++)
682 {
683 for (unsigned j = 0; j < m; j++)
684 {
685 psi(i, j) = 0.0;
686 }
687 }
689 }
690
691
692 /// Calculate the geometric shape functions
693 /// at local coordinate s. Set any "superfluous" shape functions to zero.
695 {
696 // Initialise all of them to zero
697 unsigned n = psi.nindex1();
698 unsigned m = psi.nindex2();
699 for (unsigned i = 0; i < n; i++)
700 {
701 for (unsigned j = 0; j < m; j++)
702 {
703 psi(i, j) = 0.0;
704 }
705 }
706 unsigned n1 = dpsids.nindex1();
707 unsigned n2 = dpsids.nindex2();
708 unsigned n3 = dpsids.nindex3();
709 for (unsigned i = 0; i < n1; i++)
710 {
711 for (unsigned j = 0; j < n2; j++)
712 {
713 for (unsigned k = 0; k < n3; k++)
714 {
715 dpsids(i, j, k) = 0.0;
716 }
717 }
718 }
720 }
721
722 /// Output function -- forward to broken version in FiniteElement
723 /// until somebody decides what exactly they want to plot here...
724 void output(std::ostream& outfile)
725 {
727 }
728
729 /// Output function
730 void output(std::ostream& outfile, const unsigned& n_plot);
731
732 /// C-style output function -- forward to broken version in FiniteElement
733 /// until somebody decides what exactly they want to plot here...
738
739 /// C-style output function -- forward to broken version in
740 /// FiniteElement until somebody decides what exactly they want to plot
741 /// here...
742 void output(FILE* file_pt, const unsigned& n_plot)
743 {
745 }
746
747 /// The number of "DOF types" that degrees of freedom in this element
748 /// are sub-divided into: Just the solid degrees of freedom themselves.
749 unsigned ndof_types() const
750 {
751 return 1;
752 }
753
754 /// Create a list of pairs for all unknowns in this element,
755 /// so that the first entry in each pair contains the global equation
756 /// number of the unknown, while the second one contains the number
757 /// of the "DOF types" that this unknown is associated with.
758 /// (Function can obviously only be called if the equation numbering
759 /// scheme has been set up.)
761 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list) const;
762
763
764 private:
765 /// Normal vector to the clamping plane
767 };
768
769
770} // namespace oomph
771
772#endif
static char t char * s
Definition cfortran.h:568
cstr elem_len * i
Definition cfortran.h:603
Element that allows the imposition of boundary conditions for a shell that is clamped to a 2D plane t...
void output(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: Just the soli...
ClampedHermiteShellBoundaryConditionElement()
Broken empty constructor.
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...
void set_symmetry_line(const Vector< double > &normal_to_clamping_plane)
Set normal vector to clamping plane.
void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Calculate the geometric shape functions at local coordinate s. Set any "superfluous" shape functions ...
void operator=(const ClampedHermiteShellBoundaryConditionElement &)=delete
Broken assignment operator.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the element's contribution to its residual vector.
ClampedHermiteShellBoundaryConditionElement(const ClampedHermiteShellBoundaryConditionElement &dummy)=delete
Broken copy constructor.
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 > Normal_to_clamping_plane
Normal vector to the clamping plane.
void shape(const Vector< double > &s, Shape &psi) const
Calculate the geometric shape functions at local coordinate s. Set any "superfluous" shape functions ...
void output(FILE *file_pt)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition shape.h:278
Class of matrices containing doubles, and stored as a DenseMatrix<double>, but with solving functiona...
Definition matrices.h:1271
void resize(const unsigned long &n)
Resize to a square nxn matrix; any values already present will be transfered.
Definition matrices.h:498
An element that solves the Kirchhoff-Love shell theory equations using Hermite interpolation (displac...
DiagHermiteShellElement()
Constructor, there are no internal data points.
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...
Diag Hermite Kirchhoff Love shell "upgraded" to a FSIWallElement (and thus, by inheritance,...
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 ...
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...
void set_normal_pointing_out_of_fluid()
Set the normal computed by KirchhoffLoveShellEquations::get_normal(...) to point out of the fluid.
bool Compute_rate_of_work_by_load_with_fluid_load_only
Boolean flag to indicate if rate-of-work by load is to be based on the fluid traction only.
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...
FSIDiagHermiteShellElement()
Constructor: Create shell element as FSIWallElement (and thus, by inheritance, a GeomObject) with two...
~FSIDiagHermiteShellElement()
Destructor: empty.
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...
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: Just the soli...
bool Normal_points_into_fluid
Boolean flag to indicate whether the normal is directed into the fluid.
void set_normal_pointing_into_fluid()
Set the normal computed by KirchhoffLoveShellEquations::get_normal(...) to point into the fluid.
double fluid_load_rate_of_work()
Get integral of instantaneous rate of work done on the wall due to the fluid load returned by the fun...
virtual void load_vector_for_rate_of_work_computation(const unsigned &intpt, const Vector< double > &xi, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Get the load vector for the computation of the rate of work done by the load. Can switch between full...
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
virtual void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the elemental contribution to the jacobian matrix. and the residuals vector. Note that this funct...
Definition elements.h:499
A geometric object is an object that provides a parametrised description of its shape via the functio...
An element that solves the Kirchhoff-Love shell theory equations using Hermite interpolation (displac...
void output(FILE *file_pt)
Overload the output function.
void output_with_time_dep_quantities(std::ostream &outfile, const unsigned &n_plot)
Output position veloc and accel.
void output(std::ostream &outfile)
Overload the output function.
HermiteShellElement()
Constructor, there are no internal data points.
A class for elements that solves the equations of Kirchhoff Love shell thin-shell theory.
void fill_in_contribution_to_residuals_shell(Vector< double > &residuals)
Helper function to Return the residuals for the equations of KL shell theory. This is used to prevent...
const double & lambda_sq() const
Return the timescale ratio (non-dimensional density)
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)
void output(FILE *file_pt, const unsigned &n_plot)
Generic FiniteElement output function.
bool Ignore_membrane_terms
Boolean flag to ignore membrane terms.
double *& nu_pt()
Return a pointer to the Poisson ratio.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Return the jacobian is calculated by finite differences by default,.
double * H_pt
Pointer to dimensionless wall thickness.
double *& h_pt()
Return a pointer to the non-dim wall thickness.
double *& lambda_sq_pt()
Return a pointer to timescale ratio (nondim density)
KirchhoffLoveShellEquations()
Constructor: Initialise physical parameter values to defaults.
void output(std::ostream &outfile, const unsigned &n_plot)
Generic FiniteElement output function.
const double & h() const
Return the wall thickness to undeformed radius ratio.
double prestress(const unsigned &i, const unsigned &j)
Return (i,j)-th component of second Piola Kirchhoff membrane prestress.
static double Default_lambda_sq_value
Static default value for the timescale ratio (1.0 for natural scaling)
double calculate_contravariant(double A[2][2], double Aup[2][2])
Invert a DIM by DIM matrix.
DenseMatrix< double * > Prestress_pt
Pointer to membrane pre-stress terms – should probably generalise this to function pointers at some p...
double * Lambda_sq_pt
Pointer to timescale ratio (non-dimensional density)
static double Zero_prestress
Static default for prestress (set to zero)
void get_energy(double &pot_en, double &kin_en)
Get potential (strain) and kinetic energy of the element.
void output(std::ostream &outfile)
Generic FiniteElement output function.
void get_normal(const Vector< double > &s, Vector< double > &N)
Get normal vector.
double * Nu_pt
Pointer to Poisson's ratio.
static double Default_h_value
Static default value for the thickness ratio.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Overload the standard fill in residuals contribution.
double load_rate_of_work()
Get integral of instantaneous rate of work done on the wall due to the load returned by the virtual f...
void(*&)(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &N, Vector< double > &load) load_vector_fct_pt()
Access to the load vector 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,...
static double Default_nu_value
Static default value for the Poisson's ratio.
void set_prestress_pt(const unsigned &i, const unsigned &j, double *value_pt)
Set pointer to (i,j)-th component of second Piola Kirchhoff membrane prestress to specified value (au...
GeomObject *& undeformed_midplane_pt()
Return a reference to geometric object that specifies the shell's undeformed geometry.
const double & nu() const
Return the Poisson's ratio.
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. coordinate and normal vector and...
std::pair< double, double > get_strain_and_bend(const Vector< double > &s, DenseDoubleMatrix &strain, DenseDoubleMatrix &bend)
Get strain and bending tensors; returns pair comprising the determinant of the undeformed (*....
void enable_membrane_terms()
Set to renable the calculation of membrane terms (default)
void output(FILE *file_pt)
Generic FiniteElement output function.
void disable_membrane_terms()
Set to disable the calculation of membrane terms.
virtual void load_vector_for_rate_of_work_computation(const unsigned &intpt, const Vector< double > &xi, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Get the load vector for the computation of the rate of work done by the load. Here we simply forward ...
GeomObject * Undeformed_midplane_pt
Pointer to the GeomObject that specifies the undeformed midplane of the shell.
An OomphLibError object which should be thrown when an run-time error is encountered....
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition shape.h:76
SolidQHermiteElements in which we assume the local and global coordinates to be aligned so that the J...
SolidFaceElements combine FaceElements and SolidFiniteElements and overload various functions so they...
Definition elements.h:4918
SolidFiniteElement class.
Definition elements.h:3565
SolidQHermiteElement elements are Hermite elements whose Jacobian matrices include derivatives w....
void output(std::ostream &outfile)
Overload the output function.
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).