foeppl_von_karman_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 FoepplvonKarman elements
27#ifndef OOMPH_FOEPPLVONKARMAN_ELEMENTS_HEADER
28#define OOMPH_FOEPPLVONKARMAN_ELEMENTS_HEADER
29
30
31// Config header
32#ifdef HAVE_CONFIG_H
33#include <oomph-lib-config.h>
34#endif
35
36#include <sstream>
37
38// OOMPH-LIB headers
39#include "generic/projection.h"
40#include "generic/nodes.h"
41#include "generic/Qelements.h"
43
44
45namespace oomph
46{
47 //=============================================================
48 /// A class for all isoparametric elements that solve the
49 /// Foeppl von Karman equations.
50 /// \f[ \nabla^4 w - \eta Diamond^4(w,\phi) = p(x,y) \f]
51 /// and
52 /// \f[ \nabla^4 \phi + \frac{1}{2} Diamond^4(w,w) = 0 \f]
53 /// This contains the generic maths. Shape functions, geometric
54 /// mapping etc. must get implemented in derived class.
55 //=============================================================
57 {
58 public:
59 /// Function pointer to pressure function fct(x,f(x)) --
60 /// x is a Vector!
62 double& f);
63
64 /// Constructor (must initialise the Pressure_fct_pt and
65 /// Airy_forcing_fct_pt to null). Also set physical parameters to their
66 /// default values. No volume constraint applied by default.
68 {
69 // Set all the physical constants to the default value (zero)
72
73 // No volume constraint
75 }
76
77 /// Broken copy constructor
79
80 /// Broken assignment operator
81 void operator=(const FoepplvonKarmanEquations&) = delete;
82
83 // Access functions for the physical constants
84
85 /// Eta
86 const double& eta() const
87 {
88 return *Eta_pt;
89 }
90
91 /// Pointer to eta
92 double*& eta_pt()
93 {
94 return Eta_pt;
95 }
96
97
98 /// Set Data value containing a single value which represents the
99 /// volume control pressure as external data for this element.
100 /// Only used for volume controlled problems in conjunction with
101 /// FoepplvonKarmanVolumeConstraintElement.
103 {
104#ifdef PARANOID
105 if (data_pt->nvalue() != 1)
106 {
107 throw OomphLibError("Data object that contains volume control pressure "
108 "should only contain a single value. ",
111 }
112#endif
113
114 // Add as external data and remember the index in the element's storage
115 // scheme
118 }
119
120 /// Return the index at which the i-th unknown value
121 /// is stored. The default value, i, is appropriate for single-physics
122 /// problems. By default, these are:
123 /// 0: w
124 /// 1: laplacian w
125 /// 2: phi
126 /// 3: laplacian phi
127 /// 4-8: smooth first derivatives
128 /// In derived multi-physics elements, this function should be overloaded
129 /// to reflect the chosen storage scheme. Note that these equations require
130 /// that the unknown is always stored at the same index at each node.
131 virtual inline unsigned nodal_index_fvk(const unsigned& i = 0) const
132 {
133 return i;
134 }
135
136 /// Output with default number of plot points
137 void output(std::ostream& outfile)
138 {
139 const unsigned n_plot = 5;
141 }
142
143 /// Output FE representation of soln: x,y,w at
144 /// n_plot^DIM plot points
145 void output(std::ostream& outfile, const unsigned& n_plot);
146
147 /// C_style output with default number of plot points
149 {
150 const unsigned n_plot = 5;
152 }
153
154 /// C-style output FE representation of soln: x,y,w at
155 /// n_plot^DIM plot points
156 void output(FILE* file_pt, const unsigned& n_plot);
157
158 /// Output exact soln: x,y,w_exact at n_plot^DIM plot points
159 void output_fct(std::ostream& outfile,
160 const unsigned& n_plot,
162
163 /// Output exact soln: x,y,w_exact at
164 /// n_plot^DIM plot points (dummy time-dependent version to
165 /// keep intel compiler happy)
166 virtual void output_fct(
167 std::ostream& outfile,
168 const unsigned& n_plot,
169 const double& time,
171 {
172 throw OomphLibError(
173 "There is no time-dependent output_fct() for Foeppl von Karman"
174 "elements ",
177 }
178
179
180 /// Get error against and norm of exact solution
181 void compute_error(std::ostream& outfile,
183 double& error,
184 double& norm);
185
186
187 /// Dummy, time dependent error checker
188 void compute_error(std::ostream& outfile,
190 const double& time,
191 double& error,
192 double& norm)
193 {
194 throw OomphLibError(
195 "There is no time-dependent compute_error() for Foeppl von Karman"
196 "elements",
199 }
200
201 /// Access function: Pointer to pressure function
206
207 /// Access function: Pointer to pressure function. Const version
212
213 /// Access function: Pointer to Airy forcing function
218
219 /// Access function: Pointer to Airy forcing function. Const version
224
225 /// Get pressure term at (Eulerian) position x. This function is
226 /// virtual to allow overloading in multi-physics problems where
227 /// the strength of the pressure function might be determined by
228 /// another system of equations.
229 inline virtual void get_pressure_fvk(const unsigned& ipt,
230 const Vector<double>& x,
231 double& pressure) const
232 {
233 // If no pressure function has been set, return zero
234 if (Pressure_fct_pt == 0)
235 {
236 pressure = 0.0;
237 }
238 else
239 {
240 // Get pressure strength
241 (*Pressure_fct_pt)(x, pressure);
242 }
243 }
244
245 /// Get Airy forcing term at (Eulerian) position x. This function is
246 /// virtual to allow overloading in multi-physics problems where
247 /// the strength of the pressure function might be determined by
248 /// another system of equations.
249 inline virtual void get_airy_forcing_fvk(const unsigned& ipt,
250 const Vector<double>& x,
251 double& airy_forcing) const
252 {
253 // If no pressure function has been set, return zero
254 if (Airy_forcing_fct_pt == 0)
255 {
256 airy_forcing = 0.0;
257 }
258 else
259 {
260 // Get pressure strength
261 (*Airy_forcing_fct_pt)(x, airy_forcing);
262 }
263 }
264
265 /// Get gradient of deflection: gradient[i] = dw/dx_i
268 {
269 // Find out how many nodes there are in the element
270 const unsigned n_node = nnode();
271
272 // Get the index at which the unknown is stored
273 const unsigned w_nodal_index = nodal_index_fvk(0);
274
275 // Set up memory for the shape and test functions
277 DShape dpsidx(n_node, 2);
278
279 // Call the derivatives of the shape and test functions
281
282 // Initialise to zero
283 for (unsigned j = 0; j < 2; j++)
284 {
285 gradient[j] = 0.0;
286 }
287
288 // Loop over nodes
289 for (unsigned l = 0; l < n_node; l++)
290 {
291 // Loop over derivative directions
292 for (unsigned j = 0; j < 2; j++)
293 {
294 gradient[j] += this->nodal_value(l, w_nodal_index) * dpsidx(l, j);
295 }
296 }
297 }
298
299 /// Fill in the residuals with this element's contribution
301
302 // void fill_in_contribution_to_jacobian(Vector<double> &residuals,
303 // DenseMatrix<double> &jacobian);
304
305 /// Return FE representation of function value w_fvk(s)
306 /// at local coordinate s (by default - if index > 0, returns
307 /// FE representation of valued stored at index^th nodal index
308 inline double interpolated_w_fvk(const Vector<double>& s,
309 unsigned index = 0) const
310 {
311 // Find number of nodes
312 const unsigned n_node = nnode();
313
314 // Get the index at which the poisson unknown is stored
315 const unsigned w_nodal_index = nodal_index_fvk(index);
316
317 // Local shape function
319
320 // Find values of shape function
321 shape(s, psi);
322
323 // Initialise value of u
324 double interpolated_w = 0.0;
325
326 // Loop over the local nodes and sum
327 for (unsigned l = 0; l < n_node; l++)
328 {
330 }
331
332 return (interpolated_w);
333 }
334
335 /// Compute in-plane stresses
337 double& sigma_xx,
338 double& sigma_yy,
339 double& sigma_xy);
340
341 /// Return the integral of the displacement over the current
342 /// element, effectively calculating its contribution to the volume under
343 /// the membrane. Virtual so it can be overloaded in multi-physics
344 /// where the volume may incorporate an offset, say.
345 virtual double get_bounded_volume() const
346 {
347 // Number of nodes and integration points for the current element
348 const unsigned n_node = nnode();
349 const unsigned n_intpt = integral_pt()->nweight();
350
351 // Shape functions and their derivatives
353 DShape dpsidx(n_node, 2);
354
355 // The nodal index at which the displacement is stored
356 const unsigned w_nodal_index = nodal_index_fvk();
357
358 // Initalise the integral variable
359 double integral_w = 0;
360
361 // Loop over the integration points
362 for (unsigned ipt = 0; ipt < n_intpt; ipt++)
363 {
364 // Get the integral weight
365 double w = integral_pt()->weight(ipt);
366
367 // Get determinant of the Jacobian of the mapping
369
370 // Premultiply the weight and Jacobian
371 double W = w * J;
372
373 // Initialise storage for the w value and nodal value
374 double interpolated_w = 0;
375 double w_nodal_value;
376
377 // Loop over the shape functions/nodes
378 for (unsigned l = 0; l < n_node; l++)
379 {
380 // Get the current nodal value
382 // Add the contribution to interpolated w
384 }
385
386 // Add the contribution from the current integration point
388 }
389
390 // Return the calculated integral
391 return integral_w;
392 }
393
394 /// Self-test: Return 0 for OK
395 unsigned self_test();
396
397 /// Sets a flag to signify that we are solving the linear, pure
398 /// bending equations, and pin all the nodal values that will not be used in
399 /// this case
401 {
402 // Set the boolean flag
404
405 // Get the index of the first FvK nodal value
407
408 // Get the total number of FvK nodal values (assuming they are stored
409 // contiguously) at node 0 (it's the same at all nodes anyway)
410 unsigned total_fvk_nodal_indicies = 8;
411
412 // Get the number of nodes in this element
413 unsigned n_node = nnode();
414
415 // Loop over the appropriate nodal indices
416 for (unsigned index = first_fvk_nodal_index + 2;
418 index++)
419 {
420 // Loop over the nodes in the element
421 for (unsigned inod = 0; inod < n_node; inod++)
422 {
423 // Pin the nodal value at the current index
424 node_pt(inod)->pin(index);
425 }
426 }
427 }
428
429
430 protected:
431 /// Shape/test functions and derivs w.r.t. to global coords at
432 /// local coord. s; return Jacobian of mapping
434 Shape& psi,
435 DShape& dpsidx,
436 Shape& test,
437 DShape& dtestdx) const = 0;
438
439
440 /// Shape/test functions and derivs w.r.t. to global coords at
441 /// integration point ipt; return Jacobian of mapping
443 const unsigned& ipt,
444 Shape& psi,
445 DShape& dpsidx,
446 Shape& test,
447 DShape& dtestdx) const = 0;
448
449 // Pointers to global physical constants
450
451 /// Pointer to global eta
452 double* Eta_pt;
453
454 /// Pointer to pressure function:
456
457 // mjr Ridicu-hack for made-up second pressure, q
459
460 private:
461 /// Default value for physical constants
463
464 /// Flag which stores whether we are using a linear, pure bending
465 /// model instead of the full non-linear Foeppl-von Karman
467
468 /// Index of the external Data object that represents the volume
469 /// constraint pressure (initialised to -1 indicating no such constraint)
470 /// Gets overwritten when calling
471 /// set_volume_constraint_pressure_data_as_external_data(...)
473 };
474
475
476 ///////////////////////////////////////////////////////////////////////////
477 ///////////////////////////////////////////////////////////////////////////
478 ///////////////////////////////////////////////////////////////////////////
479
480
481 // mjr QFoepplvonKarmanElements are untested!
482
483 //======================================================================
484 /// QFoepplvonKarmanElement elements are linear/quadrilateral/brick-shaped
485 /// Foeppl von Karman elements with isoparametric interpolation for the
486 /// function.
487 //======================================================================
488
489 template<unsigned NNODE_1D>
490 class QFoepplvonKarmanElement : public virtual QElement<2, NNODE_1D>,
491 public virtual FoepplvonKarmanEquations
492 {
493 private:
494 /// Static int that holds the number of variables at
495 /// nodes: always the same
496 static const unsigned Initial_Nvalue;
497
498 public:
499 /// Constructor: Call constructors for QElement and
500 /// FoepplvonKarmanEquations
505
506 /// Broken copy constructor
508 delete;
509
510 /// Broken assignment operator
512
513 /// Required # of `values' (pinned or dofs)
514 /// at node n
515 inline unsigned required_nvalue(const unsigned& n) const
516 {
517 return Initial_Nvalue;
518 }
519
520 /// Output function:
521 /// x,y,w
522 void output(std::ostream& outfile)
523 {
525 }
526
527
528 /// Output function:
529 /// x,y,w at n_plot^DIM plot points
530 void output(std::ostream& outfile, const unsigned& n_plot)
531 {
533 }
534
535
536 /// C-style output function:
537 /// x,y,w
542
543
544 /// C-style output function:
545 /// x,y,w at n_plot^DIM plot points
546 void output(FILE* file_pt, const unsigned& n_plot)
547 {
549 }
550
551
552 /// Output function for an exact solution:
553 /// x,y,w_exact at n_plot^DIM plot points
560
561
562 /// Output function for a time-dependent exact solution.
563 /// x,y,w_exact at n_plot^DIM plot points
564 /// (Calls the steady version)
565 void output_fct(std::ostream& outfile,
566 const unsigned& n_plot,
567 const double& time,
569 {
572 }
573
574
575 protected:
576 /// Shape, test functions & derivs. w.r.t. to global coords. Return
577 /// Jacobian.
579 Shape& psi,
580 DShape& dpsidx,
581 Shape& test,
582 DShape& dtestdx) const;
583
584
585 /// Shape, test functions & derivs. w.r.t. to global coords. at
586 /// integration point ipt. Return Jacobian.
587 inline double dshape_and_dtest_eulerian_at_knot_fvk(const unsigned& ipt,
588 Shape& psi,
589 DShape& dpsidx,
590 Shape& test,
591 DShape& dtestdx) const;
592 };
593
594
595 // Inline functions:
596
597
598 //======================================================================
599 /// Define the shape functions and test functions and derivatives
600 /// w.r.t. global coordinates and return Jacobian of mapping.
601 ///
602 /// Galerkin: Test functions = shape functions
603 //======================================================================
604 template<unsigned NNODE_1D>
606 const Vector<double>& s,
607 Shape& psi,
608 DShape& dpsidx,
609 Shape& test,
610 DShape& dtestdx) const
611 {
612 // Call the geometrical shape functions and derivatives
613 const double J = this->dshape_eulerian(s, psi, dpsidx);
614
615 // Set the test functions equal to the shape functions
616 test = psi;
617 dtestdx = dpsidx;
618
619 // Return the jacobian
620 return J;
621 }
622
623
624 //======================================================================
625 /// Define the shape functions and test functions and derivatives
626 /// w.r.t. global coordinates and return Jacobian of mapping.
627 ///
628 /// Galerkin: Test functions = shape functions
629 //======================================================================
630 template<unsigned NNODE_1D>
632 NNODE_1D>::dshape_and_dtest_eulerian_at_knot_fvk(const unsigned& ipt,
633 Shape& psi,
634 DShape& dpsidx,
635 Shape& test,
636 DShape& dtestdx) const
637 {
638 // Call the geometrical shape functions and derivatives
639 const double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
640
641 // Set the pointers of the test functions
642 test = psi;
643 dtestdx = dpsidx;
644
645 // Return the jacobian
646 return J;
647 }
648
649
650 ////////////////////////////////////////////////////////////////////////
651 ////////////////////////////////////////////////////////////////////////
652 ////////////////////////////////////////////////////////////////////////
653
654
655 //=======================================================================
656 /// Face geometry for the QFoepplvonKarmanElement elements: The spatial
657 /// dimension of the face elements is one lower than that of the
658 /// bulk element but they have the same number of points
659 /// along their 1D edges.
660 //=======================================================================
661 template<unsigned NNODE_1D>
663 : public virtual QElement<1, NNODE_1D>
664 {
665 public:
666 /// Constructor: Call the constructor for the
667 /// appropriate lower-dimensional QElement
669 };
670
671 ////////////////////////////////////////////////////////////////////////
672 ////////////////////////////////////////////////////////////////////////
673 ////////////////////////////////////////////////////////////////////////
674
675
676 //==========================================================
677 /// Foeppl von Karman upgraded to become projectable
678 //==========================================================
679 template<class FVK_ELEMENT>
681 : public virtual ProjectableElement<FVK_ELEMENT>
682 {
683 public:
684 /// Specify the values associated with field fld.
685 /// The information is returned in a vector of pairs which comprise
686 /// the Data object and the value within it, that correspond to field fld.
688 {
689#ifdef PARANOID
690 if (fld > 7)
691 {
692 std::stringstream error_stream;
694 << "Foeppl von Karman elements only store eight fields so fld must be"
695 << "0 to 7 rather than " << fld << std::endl;
696 throw OomphLibError(
698 }
699#endif
700
701 // Create the vector
702 unsigned nnod = this->nnode();
704
705 // Loop over all nodes
706 for (unsigned j = 0; j < nnod; j++)
707 {
708 // Add the data value associated field: The node itself
709 data_values[j] = std::make_pair(this->node_pt(j), fld);
710 }
711
712 // Return the vector
713 return data_values;
714 }
715
716 /// Number of fields to be projected: Just two
718 {
719 return 8;
720 }
721
722 /// Number of history values to be stored for fld-th field.
723 /// (Note: count includes current value!)
724 unsigned nhistory_values_for_projection(const unsigned& fld)
725 {
726#ifdef PARANOID
727 if (fld > 7)
728 {
729 std::stringstream error_stream;
731 << "Foeppl von Karman elements only store eight fields so fld must be"
732 << "0 to 7 rather than " << fld << std::endl;
733 throw OomphLibError(
735 }
736#endif
737 return this->node_pt(0)->ntstorage();
738 }
739
740
741 /// Number of positional history values
742 /// (Note: count includes current value!)
747
748 /// Return Jacobian of mapping and shape functions of field fld
749 /// at local coordinate s
750 double jacobian_and_shape_of_field(const unsigned& fld,
751 const Vector<double>& s,
752 Shape& psi)
753 {
754#ifdef PARANOID
755 if (fld > 7)
756 {
757 std::stringstream error_stream;
759 << "Foeppl von Karman elements only store eight fields so fld must be"
760 << "0 to 7 rather than " << fld << std::endl;
761 throw OomphLibError(
763 }
764#endif
765 unsigned n_dim = this->dim();
766 unsigned n_node = this->nnode();
769 double J =
770 this->dshape_and_dtest_eulerian_fvk(s, psi, dpsidx, test, dtestdx);
771 return J;
772 }
773
774
775 /// Return interpolated field fld at local coordinate s, at time
776 /// level t (t=0: present; t>0: history values)
777 double get_field(const unsigned& t,
778 const unsigned& fld,
779 const Vector<double>& s)
780 {
781#ifdef PARANOID
782 if (fld > 7)
783 {
784 std::stringstream error_stream;
786 << "Foeppl von Karman elements only store eight fields so fld must be"
787 << "0 to 7 rather than " << fld << std::endl;
788 throw OomphLibError(
790 }
791#endif
792 // Find the index at which the variable is stored
793 unsigned w_nodal_index = this->nodal_index_fvk(fld);
794
795 // Local shape function
796 unsigned n_node = this->nnode();
798
799 // Find values of shape function
800 this->shape(s, psi);
801
802 // Initialise value of u
803 double interpolated_w = 0.0;
804
805 // Sum over the local nodes
806 for (unsigned l = 0; l < n_node; l++)
807 {
809 }
810 return interpolated_w;
811 }
812
813
814 /// Return number of values in field fld: One per node
815 unsigned nvalue_of_field(const unsigned& fld)
816 {
817#ifdef PARANOID
818 if (fld > 7)
819 {
820 std::stringstream error_stream;
822 << "Foeppl von Karman elements only store eight fields so fld must be"
823 << "0 to 7 rather than " << fld << std::endl;
824 throw OomphLibError(
826 }
827#endif
828 return this->nnode();
829 }
830
831
832 /// Return local equation number of field fld of node j.
833 int local_equation(const unsigned& fld, const unsigned& j)
834 {
835#ifdef PARANOID
836 if (fld > 7)
837 {
838 std::stringstream error_stream;
840 << "Foeppl von Karman elements only store eight fields so fld must be"
841 << "0 to 7 rather than " << fld << std::endl;
842 throw OomphLibError(
844 }
845#endif
846 const unsigned w_nodal_index = this->nodal_index_fvk(fld);
847 return this->nodal_local_eqn(j, w_nodal_index);
848 }
849 };
850
851 //=======================================================================
852 /// Face geometry for element is the same as that for the underlying
853 /// wrapped element
854 //=======================================================================
855 template<class ELEMENT>
857 : public virtual FaceGeometry<ELEMENT>
858 {
859 public:
860 FaceGeometry() : FaceGeometry<ELEMENT>() {}
861 };
862
863
864 //=======================================================================
865 /// Face geometry of the Face Geometry for element is the same as
866 /// that for the underlying wrapped element
867 //=======================================================================
868 template<class ELEMENT>
870 : public virtual FaceGeometry<FaceGeometry<ELEMENT>>
871 {
872 public:
874 };
875
876} // namespace oomph
877
878#endif
static char t char * s
Definition cfortran.h:568
cstr elem_len * i
Definition cfortran.h:603
char t
Definition cfortran.h:568
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition shape.h:278
A class that represents a collection of data; each Data object may contain many different individual ...
Definition nodes.h:86
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition nodes.h:385
unsigned ntstorage() const
Return total number of doubles stored per value to record time history of each value (one for steady ...
Definition nodes.cc:879
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
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
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition elements.h:1967
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
Definition elements.h:2597
virtual 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-...
Definition elements.cc:3355
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
Definition elements.h:1436
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
Definition elements.h:2615
unsigned nnode() const
Return the number of nodes.
Definition elements.h:2214
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition elements.h:1763
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
Definition elements.cc:3328
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition elements.h:2179
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n but do NOT take hanging nodes into account.
Definition elements.h:2580
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
Definition elements.h:1769
A class for all isoparametric elements that solve the Foeppl von Karman equations.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
FoepplvonKarmanPressureFctPt & pressure_fct_pt()
Access function: Pointer to pressure function.
void output(FILE *file_pt)
C_style output with default number of plot points.
double interpolated_w_fvk(const Vector< double > &s, unsigned index=0) const
Return FE representation of function value w_fvk(s) at local coordinate s (by default - if index > 0,...
FoepplvonKarmanPressureFctPt pressure_fct_pt() const
Access function: Pointer to pressure function. Const version.
void get_gradient_of_deflection(const Vector< double > &s, Vector< double > &gradient) const
Get gradient of deflection: gradient[i] = dw/dx_i.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: x,y,w_exact at n_plot^DIM plot points.
FoepplvonKarmanPressureFctPt & airy_forcing_fct_pt()
Access function: Pointer to Airy forcing function.
FoepplvonKarmanPressureFctPt airy_forcing_fct_pt() const
Access function: Pointer to Airy forcing function. Const version.
void(* FoepplvonKarmanPressureFctPt)(const Vector< double > &x, double &f)
Function pointer to pressure function fct(x,f(x)) – x is a Vector!
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: x,y,w_exact at n_plot^DIM plot points (dummy time-dependent version to keep intel ...
void output(std::ostream &outfile)
Output with default number of plot points.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals with this element's contribution.
virtual void get_airy_forcing_fvk(const unsigned &ipt, const Vector< double > &x, double &airy_forcing) const
Get Airy forcing term at (Eulerian) position x. This function is virtual to allow overloading in mult...
FoepplvonKarmanEquations()
Constructor (must initialise the Pressure_fct_pt and Airy_forcing_fct_pt to null)....
virtual double dshape_and_dtest_eulerian_fvk(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at local coord. s; return Jacobian of mapping...
int Volume_constraint_pressure_external_data_index
Index of the external Data object that represents the volume constraint pressure (initialised to -1 i...
void use_linear_bending_model()
Sets a flag to signify that we are solving the linear, pure bending equations, and pin all the nodal ...
virtual double get_bounded_volume() const
Return the integral of the displacement over the current element, effectively calculating its contrib...
bool Linear_bending_model
Flag which stores whether we are using a linear, pure bending model instead of the full non-linear Fo...
virtual void get_pressure_fvk(const unsigned &ipt, const Vector< double > &x, double &pressure) const
Get pressure term at (Eulerian) position x. This function is virtual to allow overloading in multi-ph...
double * Eta_pt
Pointer to global eta.
virtual unsigned nodal_index_fvk(const unsigned &i=0) const
Return the index at which the i-th unknown value is stored. The default value, i, is appropriate for ...
unsigned self_test()
Self-test: Return 0 for OK.
virtual double dshape_and_dtest_eulerian_at_knot_fvk(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
FoepplvonKarmanPressureFctPt Pressure_fct_pt
Pointer to pressure function:
FoepplvonKarmanEquations(const FoepplvonKarmanEquations &dummy)=delete
Broken copy constructor.
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
static double Default_Physical_Constant_Value
Default value for physical constants.
void interpolated_stress(const Vector< double > &s, double &sigma_xx, double &sigma_yy, double &sigma_xy)
Compute in-plane stresses.
FoepplvonKarmanPressureFctPt Airy_forcing_fct_pt
void set_volume_constraint_pressure_data_as_external_data(Data *data_pt)
Set Data value containing a single value which represents the volume control pressure as external dat...
void operator=(const FoepplvonKarmanEquations &)=delete
Broken assignment operator.
unsigned add_external_data(Data *const &data_pt, const bool &fd=true)
Add a (pointer to an) external data object to the element and return its index (i....
Definition elements.cc:312
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
Definition nodes.h:1022
An OomphLibError object which should be thrown when an run-time error is encountered....
Wrapper class for projectable elements. Adds "projectability" to the underlying ELEMENT.
Definition projection.h:183
Foeppl von Karman upgraded to become projectable.
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Return Jacobian of mapping and shape functions of field fld at local coordinate s.
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Return interpolated field fld at local coordinate s, at time level t (t=0: present; t>0: history valu...
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (Note: count includes current value!)
unsigned nfields_for_projection()
Number of fields to be projected: Just two.
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of field fld of node j.
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Specify the values associated with field fld. The information is returned in a vector of pairs which ...
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. (Note: count includes current value!...
General QElement class.
Definition Qelements.h:459
QFoepplvonKarmanElement elements are linear/quadrilateral/brick-shaped Foeppl von Karman elements wit...
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: x,y,w at n_plot^DIM plot points.
double dshape_and_dtest_eulerian_fvk(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
void output(std::ostream &outfile)
Output function: x,y,w.
QFoepplvonKarmanElement(const QFoepplvonKarmanElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
double dshape_and_dtest_eulerian_at_knot_fvk(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape, test functions & derivs. w.r.t. to global coords. at integration point ipt....
static const unsigned Initial_Nvalue
Static int that holds the number of variables at nodes: always the same.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: x,y,w at n_plot^DIM plot points.
void output(FILE *file_pt)
C-style output function: x,y,w.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: x,y,w_exact at n_plot^DIM plot points.
unsigned required_nvalue(const unsigned &n) const
Required # of ‘values’ (pinned or dofs) at node n.
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output function for a time-dependent exact solution. x,y,w_exact at n_plot^DIM plot points (Calls the...
QFoepplvonKarmanElement()
Constructor: Call constructors for QElement and FoepplvonKarmanEquations.
void operator=(const QFoepplvonKarmanElement< NNODE_1D > &)=delete
Broken assignment operator.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition shape.h:76
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
TAdvectionDiffusionReactionElement()
Constructor: Call constructors for TElement and AdvectionDiffusionReaction equations.
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).