stored_shape_function_elements.h
Go to the documentation of this file.
1// LIC// ====================================================================
2// LIC// This file forms part of oomph-lib, the object-oriented,
3// LIC// multi-physics finite-element library, available
4// LIC// at http://www.oomph-lib.org.
5// LIC//
6// LIC// Copyright (C) 2006-2024 Matthias Heil and Andrew Hazel
7// LIC//
8// LIC// This library is free software; you can redistribute it and/or
9// LIC// modify it under the terms of the GNU Lesser General Public
10// LIC// License as published by the Free Software Foundation; either
11// LIC// version 2.1 of the License, or (at your option) any later version.
12// LIC//
13// LIC// This library is distributed in the hope that it will be useful,
14// LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15// LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16// LIC// Lesser General Public License for more details.
17// LIC//
18// LIC// You should have received a copy of the GNU Lesser General Public
19// LIC// License along with this library; if not, write to the Free Software
20// LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21// LIC// 02110-1301 USA.
22// LIC//
23// LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24// LIC//
25// LIC//====================================================================
26// Header file for classes that overload elements and implement
27// stored shape functions.
28
29// Include guard to prevent multiple inclusions of the header
30#ifndef OOMPH_STORED_SHAPE_FUNCTION_ELEMENTS_HEADER
31#define OOMPH_STORED_SHAPE_FUNCTION_ELEMENTS_HEADER
32
33
34// Config header generated by autoconfig
35#ifdef HAVE_CONFIG_H
36#include <oomph-lib-config.h>
37#endif
38
39#include "elements.h"
40
41
42namespace oomph
43{
44// Debugging flag
45#define OOMPH_STORED_SHAPE_FUNCTIONS_VERBOSE
46#undef OOMPH_STORED_SHAPE_FUNCTIONS_VERBOSE
47
48
49 //==========================================================================
50 /// Base class for elements that allow storage of precomputed shape functions
51 /// and their derivatives w.r.t to the local and global (Eulerian)
52 /// coordinates at the element's integration points.
53 //==========================================================================
55 {
56 private:
57 /// Pointer to storage for the pointers to the nodal shape functions
58 /// at the integration points (knots)
59 // ALH: Note that the vector must be Shape* because we do not know,
60 // a priori, how much storage to allocate in each Shape object.
61 // N.B. This could in 99% of cases be static, but that would not
62 // permit individual elements to have different integration schemes.
63 // If this proves to be a problem, one can use the the function
64 // set_local_shape_stored_from_element(), which sets this pointer to
65 // point to the pointer allocated by another element.
67
68 /// Pointer to storage for the pointers to the derivatives of the
69 /// nodal shape functions w.r.t. the local coordinates at integration points
71
72 /// Pointer to storage for the pointers to the second derivatives of
73 /// the nodal shape functions w.r.t. the local coordinates at integration
74 /// points
76
77 /// Boolean to determine whether the element can delete the stored
78 /// local shape functions
80
81 /// Pointer to storage for the derivatives of the
82 /// shape functions w.r.t. global coordinates at integration points
84
85 /// Pointer to storage for the second derivatives of the
86 /// shape functions w.r.t. global coordinates at integration points
88
89 /// Pointer to storage for the Jacobian of the element w.r.t
90 /// global coordinates
92
93 /// Boolean to determine whether the element can delete the stored
94 /// derivatives of shape functions w.r.t. global coordinates
96
97 public:
98 /// Constructor, set most storage pointers to NULL.
99 // By default the element can delete its own stored shape functions
112
113 /// The destructor cleans up the static memory allocated
114 /// for shape function storage. Internal and external data get
115 /// wiped by the GeneralisedElement destructor; nodes get
116 /// killed in mesh destructor.
118
119 /// Broken copy constructor
121
122 /// Broken assignment operator
123 void operator=(const StorableShapeElementBase&) = delete;
124
125 /// Delete all the objects stored in the [...]_local_stored_pt
126 /// vectors and delete the vectors themselves
128
129 /// Delete stored shape functions
131
132 /// Delete stored derivatives of shape functions w.r.t. to
133 /// local coordinates
135
136 /// Delete stored 2nd derivatives of shape functions w.r.t. to
137 /// local coordinates
139
140 /// Delete all storage related to deriatives of shape fcts
141 /// w.r.t. to global Eulerian coords
143
144 /// Delete stored deriatives of shape fcts w.r.t. to global
145 /// Eulerian coords
147
148 /// Delete stored second derivatives of shape functions w.r.t. to
149 /// global Eulerian coordinates
151
152 /// Delete stored Jacobian of mapping between local and global
153 /// (Eulerian) coordinates
155
156 /// Set the spatial integration scheme -- overloaded from the
157 /// finite element base class since a change in the integration scheme
158 /// forces a recomputation of the shape fcts at the integration points.
159 virtual void set_integration_scheme(Integral* const& integral_pt);
160
161 /// Return a pointer to the vector of pointers to the
162 /// stored shape functions
164 {
165 return Shape_stored_pt;
166 }
167
168 /// Return a pointer to the vector of pointers to the
169 /// stored shape functions (const version)
170 inline Vector<Shape*>* const& shape_stored_pt() const
171 {
172 return Shape_stored_pt;
173 }
174
175 /// Return a pointer to the stored shape function at the ipt-th
176 /// integration point
177 inline Shape*& shape_stored_pt(const unsigned& ipt)
178 {
179 return (*Shape_stored_pt)[ipt];
180 }
181
182 /// Return a pointer to the stored shape function at the ipt-th
183 /// integration point (const version)
184 inline Shape* const& shape_stored_pt(const unsigned& ipt) const
185 {
186 return (*Shape_stored_pt)[ipt];
187 }
188
189 /// Return a pointer to the vector of pointers to the stored first
190 /// derivatives of the shape functions w.r.t the local coordinates
195
196 /// Return a pointer to the vector of pointers to the stored first
197 /// derivatives of the shape functions w.r.t the local coordinates
198 /// (const version)
200 {
202 }
203
204 /// Return a pointer to the vector of pointers to the stored second
205 /// derivatives of the shape functions w.r.t the local coordinates
210
211 /// Return a pointer to the vector of pointers to the stored second
212 /// derivatives of the shape functions w.r.t the local coordinates
213 /// (const version)
215 {
217 }
218
219 /// Return a pointer to the vector of pointers to the stored first
220 /// derivatives of the shape functions w.r.t the global (eulerian)
221 /// coordinates
226
227 /// Return a pointer to the vector of pointers to the stored first
228 /// derivatives of the shape functions w.r.t the global (eulerian)
229 /// coordinates (const version)
231 {
233 }
234
235 /// Return a pointer to the vector of pointers to the stored second
236 /// derivatives of the shape functions w.r.t the global (eulerian)
237 /// coordinates
242
243 /// Return a pointer to the vector of pointers to the stored second
244 /// derivatives of the shape functions w.r.t the global (eulerian)
245 /// coordinates (const version)
247 {
249 }
250
251 /// Return a pointer to the vector of Jacobians of
252 /// the mapping between the local and global (eulerian) coordinates
257
258 /// Return a pointer to the vector of Jacobians of
259 /// the mapping between the local and global (eulerian) coordinates
260 /// (const version)
262 {
264 }
265
266 /// Set the shape functions pointed to internally to be
267 /// those pointed to by the FiniteElement element_pt (In most
268 /// cases all elements of the same type have the same number of
269 /// integration points so the shape function values and their
270 /// local derivatives are the same --> They only need to be
271 /// stored by one element). Calling this function deletes the locally
272 /// created storage and re-directs the pointers to the stored
273 /// shape function of the specified element.
275 StorableShapeElementBase* const& element_pt);
276
277 /// Set the derivatives of stored shape functions with respect
278 /// to the global coordinates to be the same as
279 /// those pointed to by the FiniteElement element_pt. Note that this
280 /// function also calls set_shape_local_stored_from_element(element_pt).
281 /// so that the local derivatives are also stored.
282 /// Calling this function only makes sense for uniformly-spaced meshes with
283 /// elements of equal sizes.
285 StorableShapeElementBase* const& element_pt);
286
287 /// Calculate the shape functions at the integration points
288 /// and store the results internally
290
291 /// Return the geometric shape function at the ipt-th integration
292 /// point
293 void shape_at_knot(const unsigned& ipt, Shape& psi) const;
294
295 /// Calculate the shape functions and first derivatives w.r.t. local
296 /// coordinatess at the integration points and store the results internally
298
299 /// Return the geometric shape function and its derivative w.r.t.
300 /// the local coordinates at the ipt-th integration point. If pre-computed
301 /// values have been stored, they will be used.
302 void dshape_local_at_knot(const unsigned& ipt,
303 Shape& psi,
304 DShape& dpsids) const;
305
306 /// Calculate the second derivatives of the shape functions
307 /// w.r.t. local coordinates at the integration points and store the
308 /// results internally
310
311 /// Return the geometric shape function and its first and
312 /// second derivatives w.r.t.
313 /// the local coordinates at the ipt-th integration point. If pre-computed
314 /// values have been stored, they will be used.
315 void d2shape_local_at_knot(const unsigned& ipt,
316 Shape& psi,
317 DShape& dpsids,
318 DShape& d2psids) const;
319
320 /// Calculate the Jacobian of the mapping from local to global
321 /// coordinates at the integration points and store the results
322 /// internally.
324
325 /// Return the Jacobian of the mapping from local to global
326 /// coordinates at the ipt-th integration point
327 double J_eulerian_at_knot(const unsigned& ipt) const;
328
329 /// Calculate the first derivatives of the shape functions
330 /// w.r.t the global coordinates at the integration points and store
331 /// the results internally
333
334 /// Return the geometric shape functions and also first
335 /// derivatives w.r.t. global coordinates at the ipt-th integration point.
336 /// If the values of the shape functions and derivatives have been
337 /// pre-computed, these will be used
338 double dshape_eulerian_at_knot(const unsigned& ipt,
339 Shape& psi,
340 DShape& dpsidx) const;
341
342
343 /// Calculate the first and second derivatives of the shape
344 /// functions w.r.t global coordinates at the integration points and
345 /// store the results internally.
347
348 /// Return the geometric shape functions and also first
349 /// and second derivatives w.r.t. global coordinates at ipt-th integration
350 /// point. If the values of the shape functions and derivatives have been
351 /// pre-computred, these will be used
352 double d2shape_eulerian_at_knot(const unsigned& ipt,
353 Shape& psi,
354 DShape& dpsidx,
355 DShape& d2psidx) const;
356
357
358 /* /// Diagnostic */
359 /* void tell_me() */
360 /* { */
361 /* oomph_info << "Diagnostic" << std::endl; */
362 /* oomph_info << Shape_stored_pt << " "; */
363 /* oomph_info << DShape_local_stored_pt << " "; */
364 /* oomph_info << D2Shape_local_stored_pt << " "; */
365 /* oomph_info << DShape_eulerian_stored_pt << " "; */
366 /* oomph_info << D2Shape_eulerian_stored_pt << " "; */
367 /* oomph_info << Jacobian_eulerian_stored_pt << std::endl; */
368 /* } */
369 };
370
371
372 //============================================================================
373 /// Base class for solid elements that allow storage of precomputed
374 /// shape functions and their derivatives w.r.t to the local and global
375 /// (Lagrangian) coordinates at the element's integration points.
376 //============================================================================
378 public virtual SolidFiniteElement
379 {
380 private:
381 /// Pointer to storage for the pointers to the derivatives of the
382 /// shape functions w.r.t. Lagrangian coordinates at integration points
384
385 /// Pointer to storage for the pointers to the second derivatives of
386 /// the shape functions w.r.t. Lagrangian coordinates at integration points
388
389 /// Pointer to storage for the Jacobian of the mapping between
390 /// the local and the global Lagrangian coordinates
392
393 /// Boolean to determine whether the element can delete the stored
394 /// shape function derivatives w.r.t. the Lagrangian coordinate
396
397 public:
398 /// Constructor: Set defaults: Nothing is stored
408
409 /// Destructor to clean up any allocated memory
414
415 /// Broken copy constructor
417 delete;
418
419 /// Broken assignment operator
421
422 /// Delete all the objects stored in the [...]_lagrangian_stored_pt
423 /// vectors and delete the vectors themselves
425
426 /// Delete all the objects stored in the
427 /// Lagrangian_stored vectors
429
430 /// Delete stored second derivatives of shape functions w.r.t.
431 /// Lagrangian coordinates
433
434 /// Delete stored Jaocbian of mapping between local and Lagrangian
435 /// coordinates
437
438 /// Overload the set_integration_scheme to recompute any stored
439 /// derivatives w.r.t. Lagrangian coordinates
441 {
443
444 // If we are storing Lagrangian first and second derivatives, recompute
445 // them
447 {
449 }
450 // If we are storing Lagrangian first derivatives, recompute them
451 else if (DShape_lagrangian_stored_pt != 0)
452 {
454 }
455 }
456
457 /// Return a pointer to the vector of pointers to the stored first
458 /// derivatives of the shape functions w.r.t the global (eulerian)
459 /// coordinates
464
465 /// Return a pointer to the vector of pointers to the stored first
466 /// derivatives of the shape functions w.r.t the global (eulerian)
467 /// coordinates (const version)
469 {
471 }
472
473 /// Return a pointer to the vector of pointers to the stored second
474 /// derivatives of the shape functions w.r.t the global (eulerian)
475 /// coordinates
480
481 /// Return a pointer to the vector of pointers to the stored second
482 /// derivatives of the shape functions w.r.t the global (eulerian)
483 /// coordinates (const version)
485 {
487 }
488
489 /// Return a pointer to the vector of Jacobians of
490 /// the mapping between the local and global (eulerian) coordinates
495
496 /// Return a pointer to the vector of Jacobians of
497 /// the mapping between the local and global (eulerian) coordinates
498 /// (const version)
500 {
502 }
503
504
505 /// Calculate the first derivatives of the shape functions w.r.t
506 /// Lagrangian coordinates at the integration points and store the
507 /// results internally.
509
510 /// Return the geometric shape functions and also first
511 /// derivatives w.r.t. Lagrangian coordinates at ipt-th integration point.
512 /// If the values of the shape function and derivatives have been
513 /// pre-computed, they will be used.
514 double dshape_lagrangian_at_knot(const unsigned& ipt,
515 Shape& psi,
516 DShape& dpsidxi) const;
517
518 /// Calculate the first and second derivatives of the
519 /// shape functions w.r.t
520 /// Lagrangian coordinates at the integration points and store the
521 /// results internally
523
524 /// Return the geometric shape functions and also first
525 /// and second derivatives w.r.t. Lagrangian coordinates at
526 /// the ipt-th integration point. If the values have been pre-computed,
527 /// they will be used.
528 /// Returns Jacobian of mapping from Lagrangian to local coordinates.
529 double d2shape_lagrangian_at_knot(const unsigned& ipt,
530 Shape& psi,
532 DShape& d2psidxi) const;
533
534
535 /// Set the derivatives of stored shape functions with respect
536 /// to the lagrangian coordinates to be the same as
537 /// those pointed to by the FiniteElement element_pt. Note that this
538 /// function also calls set_shape_local_stored_from_element(element_pt).
539 /// so that the local derivatives are also stored.
540 /// Calling this function only makes sense for uniformly-spaced meshes with
541 /// elements of equal sizes.
543 StorableShapeSolidElementBase* const& element_pt);
544 };
545
546
547 //========================================================================
548 /// Templated wrapper that attaches the ability to store the shape
549 /// functions and their derivatives w.r.t. to the local and global
550 /// (Eulerian) coordinates at the integration points to the
551 /// element specified by the template parameter.
552 //========================================================================
553 template<class ELEMENT>
555 public virtual ELEMENT
556 {
557 public:
558 /// Constructor, set most storage pointers to zero
559 // By default the element can delete its own stored shape functions
561 {
562 // Reset the integration scheme and force a (re)compute of the
563 // the stored shape functions and their derivatives w.r.t. to the
564 // local coordinates.
565 this->set_integration_scheme(this->integral_pt());
566 }
567
568 /// Empty virtual destructor
570
571 /// Broken copy constructor
573
574 /// Broken assignment operator
575 void operator=(const StorableShapeElement&) = delete;
576 };
577
578
579 /// //////////////////////////////////////////////////////////////////////
580 /// //////////////////////////////////////////////////////////////////////
581 /// //////////////////////////////////////////////////////////////////////
582
583
584 //============================================================================
585 /// Templated wrapper that attaches the ability to store the shape
586 /// functions and their derivatives w.r.t. to the local and global
587 /// (Eulerian) coordinates at the integration points to the
588 /// SolidFiniteElement specified by the template parameter.
589 //============================================================================
590 template<class ELEMENT>
592 : public virtual StorableShapeSolidElementBase,
593 public virtual ELEMENT
594 {
595 public:
596 /// Constructor: Set defaults
598 {
599 // Reset the integration scheme and force a (re-)computation
600 // of the shape functions and their derivatives w.r.t. to the
601 // local coordinates at the integration points.
602 this->set_integration_scheme(this->integral_pt());
603 }
604
605 /// Destructor to clean up any allocated memory
607
608
609 /// Broken copy constructor
611
612 /// Broken assignment operator
614 };
615
616} // namespace oomph
617
618#endif
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition shape.h:278
A general Finite Element class.
Definition elements.h:1317
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition elements.h:1967
Generic class for numerical integration schemes:
Definition integral.h:49
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition shape.h:76
////////////////////////////////////////////////////////////////////// //////////////////////////////...
Definition elements.h:3565
Base class for elements that allow storage of precomputed shape functions and their derivatives w....
Vector< DShape * > *& dshape_eulerian_stored_pt()
Return a pointer to the vector of pointers to the stored first derivatives of the shape functions w....
Vector< Shape * > *const & shape_stored_pt() const
Return a pointer to the vector of pointers to the stored shape functions (const version)
void shape_at_knot(const unsigned &ipt, Shape &psi) const
Return the geometric shape function at the ipt-th integration point.
virtual void set_integration_scheme(Integral *const &integral_pt)
Set the spatial integration scheme – overloaded from the finite element base class since a change in ...
Shape *& shape_stored_pt(const unsigned &ipt)
Return a pointer to the stored shape function at the ipt-th integration point.
void delete_all_dshape_eulerian_stored()
Delete all storage related to deriatives of shape fcts w.r.t. to global Eulerian coords.
Vector< DShape * > *& d2shape_eulerian_stored_pt()
Return a pointer to the vector of pointers to the stored second derivatives of the shape functions w....
Vector< DShape * > * DShape_local_stored_pt
Pointer to storage for the pointers to the derivatives of the nodal shape functions w....
Vector< DShape * > * D2Shape_local_stored_pt
Pointer to storage for the pointers to the second derivatives of the nodal shape functions w....
virtual ~StorableShapeElementBase()
The destructor cleans up the static memory allocated for shape function storage. Internal and externa...
Vector< Shape * > *& shape_stored_pt()
Return a pointer to the vector of pointers to the stored shape functions.
Vector< DShape * > *const & d2shape_local_stored_pt() const
Return a pointer to the vector of pointers to the stored second derivatives of the shape functions w....
Vector< DShape * > * D2Shape_eulerian_stored_pt
Pointer to storage for the second derivatives of the shape functions w.r.t. global coordinates at int...
void delete_dshape_local_stored()
Delete stored derivatives of shape functions w.r.t. to local coordinates.
void delete_shape_local_stored()
Delete stored shape functions.
Vector< DShape * > *const & dshape_eulerian_stored_pt() const
Return a pointer to the vector of pointers to the stored first derivatives of the shape functions w....
void delete_d2shape_local_stored()
Delete stored 2nd derivatives of shape functions w.r.t. to local coordinates.
StorableShapeElementBase(const StorableShapeElementBase &)=delete
Broken copy constructor.
double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Return the geometric shape functions and also first derivatives w.r.t. global coordinates at the ipt-...
Vector< Shape * > * Shape_stored_pt
Pointer to storage for the pointers to the nodal shape functions at the integration points (knots)
Vector< double > *& jacobian_eulerian_stored_pt()
Return a pointer to the vector of Jacobians of the mapping between the local and global (eulerian) co...
void pre_compute_d2shape_eulerian_at_knots()
Calculate the first and second derivatives of the shape functions w.r.t global coordinates at the int...
Vector< double > * Jacobian_eulerian_stored_pt
Pointer to storage for the Jacobian of the element w.r.t global coordinates.
void delete_d2shape_eulerian_stored()
Delete stored second derivatives of shape functions w.r.t. to global Eulerian coordinates.
void dshape_local_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsids) const
Return the geometric shape function and its derivative w.r.t. the local coordinates at the ipt-th int...
Vector< DShape * > *const & dshape_local_stored_pt() const
Return a pointer to the vector of pointers to the stored first derivatives of the shape functions w....
void pre_compute_dshape_eulerian_at_knots()
Calculate the first derivatives of the shape functions w.r.t the global coordinates at the integratio...
void set_shape_local_stored_from_element(StorableShapeElementBase *const &element_pt)
Set the shape functions pointed to internally to be those pointed to by the FiniteElement element_pt ...
Vector< DShape * > *const & d2shape_eulerian_stored_pt() const
Return a pointer to the vector of pointers to the stored second derivatives of the shape functions w....
void pre_compute_J_eulerian_at_knots()
Calculate the Jacobian of the mapping from local to global coordinates at the integration points and ...
double d2shape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
Return the geometric shape functions and also first and second derivatives w.r.t. global coordinates ...
void set_dshape_eulerian_stored_from_element(StorableShapeElementBase *const &element_pt)
Set the derivatives of stored shape functions with respect to the global coordinates to be the same a...
Vector< double > *const & jacobian_eulerian_stored_pt() const
Return a pointer to the vector of Jacobians of the mapping between the local and global (eulerian) co...
Vector< DShape * > *& d2shape_local_stored_pt()
Return a pointer to the vector of pointers to the stored second derivatives of the shape functions w....
void delete_J_eulerian_stored()
Delete stored Jacobian of mapping between local and global (Eulerian) coordinates.
void pre_compute_dshape_local_at_knots()
Calculate the shape functions and first derivatives w.r.t. local coordinatess at the integration poin...
void delete_dshape_eulerian_stored()
Delete stored deriatives of shape fcts w.r.t. to global Eulerian coords.
bool Can_delete_dshape_eulerian_stored
Boolean to determine whether the element can delete the stored derivatives of shape functions w....
StorableShapeElementBase()
Constructor, set most storage pointers to NULL.
bool Can_delete_shape_local_stored
Boolean to determine whether the element can delete the stored local shape functions.
void pre_compute_shape_at_knots()
Calculate the shape functions at the integration points and store the results internally.
double J_eulerian_at_knot(const unsigned &ipt) const
Return the Jacobian of the mapping from local to global coordinates at the ipt-th integration point.
Vector< DShape * > *& dshape_local_stored_pt()
Return a pointer to the vector of pointers to the stored first derivatives of the shape functions w....
Vector< DShape * > * DShape_eulerian_stored_pt
Pointer to storage for the derivatives of the shape functions w.r.t. global coordinates at integratio...
void delete_all_shape_local_stored()
Delete all the objects stored in the [...]_local_stored_pt vectors and delete the vectors themselves.
void pre_compute_d2shape_local_at_knots()
Calculate the second derivatives of the shape functions w.r.t. local coordinates at the integration p...
Shape *const & shape_stored_pt(const unsigned &ipt) const
Return a pointer to the stored shape function at the ipt-th integration point (const version)
void d2shape_local_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsids, DShape &d2psids) const
Return the geometric shape function and its first and second derivatives w.r.t. the local coordinates...
void operator=(const StorableShapeElementBase &)=delete
Broken assignment operator.
Templated wrapper that attaches the ability to store the shape functions and their derivatives w....
StorableShapeElement(const StorableShapeElement &)=delete
Broken copy constructor.
void operator=(const StorableShapeElement &)=delete
Broken assignment operator.
StorableShapeElement()
Constructor, set most storage pointers to zero.
virtual ~StorableShapeElement()
Empty virtual destructor.
Base class for solid elements that allow storage of precomputed shape functions and their derivatives...
Vector< double > * Jacobian_lagrangian_stored_pt
Pointer to storage for the Jacobian of the mapping between the local and the global Lagrangian coordi...
Vector< DShape * > *const & d2shape_lagrangian_stored_pt() const
Return a pointer to the vector of pointers to the stored second derivatives of the shape functions w....
virtual ~StorableShapeSolidElementBase()
Destructor to clean up any allocated memory.
Vector< DShape * > *& dshape_lagrangian_stored_pt()
Return a pointer to the vector of pointers to the stored first derivatives of the shape functions w....
void pre_compute_d2shape_lagrangian_at_knots()
Calculate the first and second derivatives of the shape functions w.r.t Lagrangian coordinates at the...
Vector< DShape * > *const & dshape_lagrangian_stored_pt() const
Return a pointer to the vector of pointers to the stored first derivatives of the shape functions w....
void pre_compute_dshape_lagrangian_at_knots()
Calculate the first derivatives of the shape functions w.r.t Lagrangian coordinates at the integratio...
double d2shape_lagrangian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidxi, DShape &d2psidxi) const
Return the geometric shape functions and also first and second derivatives w.r.t. Lagrangian coordina...
StorableShapeSolidElementBase(const StorableShapeSolidElementBase &)=delete
Broken copy constructor.
Vector< DShape * > * D2Shape_lagrangian_stored_pt
Pointer to storage for the pointers to the second derivatives of the shape functions w....
double dshape_lagrangian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidxi) const
Return the geometric shape functions and also first derivatives w.r.t. Lagrangian coordinates at ipt-...
Vector< double > *const & jacobian_lagrangian_stored_pt() const
Return a pointer to the vector of Jacobians of the mapping between the local and global (eulerian) co...
void delete_J_lagrangian_stored()
Delete stored Jaocbian of mapping between local and Lagrangian coordinates.
StorableShapeSolidElementBase()
Constructor: Set defaults: Nothing is stored.
Vector< DShape * > *& d2shape_lagrangian_stored_pt()
Return a pointer to the vector of pointers to the stored second derivatives of the shape functions w....
void delete_d2shape_lagrangian_stored()
Delete stored second derivatives of shape functions w.r.t. Lagrangian coordinates.
bool Can_delete_dshape_lagrangian_stored
Boolean to determine whether the element can delete the stored shape function derivatives w....
void set_integration_scheme(Integral *const &integral_pt)
Overload the set_integration_scheme to recompute any stored derivatives w.r.t. Lagrangian coordinates...
void delete_all_dshape_lagrangian_stored()
Delete all the objects stored in the [...]_lagrangian_stored_pt vectors and delete the vectors themse...
void operator=(const StorableShapeSolidElementBase &)=delete
Broken assignment operator.
void delete_dshape_lagrangian_stored()
Delete all the objects stored in the Lagrangian_stored vectors.
Vector< DShape * > * DShape_lagrangian_stored_pt
Pointer to storage for the pointers to the derivatives of the shape functions w.r....
Vector< double > *& jacobian_lagrangian_stored_pt()
Return a pointer to the vector of Jacobians of the mapping between the local and global (eulerian) co...
void set_dshape_lagrangian_stored_from_element(StorableShapeSolidElementBase *const &element_pt)
Set the derivatives of stored shape functions with respect to the lagrangian coordinates to be the sa...
////////////////////////////////////////////////////////////////////// //////////////////////////////...
void operator=(const StorableShapeSolidElement &)=delete
Broken assignment operator.
StorableShapeSolidElement(const StorableShapeSolidElement &)=delete
Broken copy constructor.
StorableShapeSolidElement()
Constructor: Set defaults.
virtual ~StorableShapeSolidElement()
Destructor to clean up any allocated memory.
//////////////////////////////////////////////////////////////////////
//////////////////////////////////////////////////////////////////// ////////////////////////////////...