27#ifndef OOMPH_FOEPPLVONKARMAN_ELEMENTS_HEADER
28#define OOMPH_FOEPPLVONKARMAN_ELEMENTS_HEADER
33#include <oomph-lib-config.h>
86 const double&
eta()
const
107 throw OomphLibError(
"Data object that contains volume control pressure "
108 "should only contain a single value. ",
139 const unsigned n_plot = 5;
150 const unsigned n_plot = 5;
173 "There is no time-dependent output_fct() for Foeppl von Karman"
195 "There is no time-dependent compute_error() for Foeppl von Karman"
231 double& pressure)
const
241 (*Pressure_fct_pt)(x, pressure);
283 for (
unsigned j = 0;
j < 2;
j++)
292 for (
unsigned j = 0;
j < 2;
j++)
309 unsigned index = 0)
const
489 template<
unsigned NNODE_1D>
604 template<
unsigned NNODE_1D>
630 template<
unsigned NNODE_1D>
632 NNODE_1D>::dshape_and_dtest_eulerian_at_knot_fvk(
const unsigned&
ipt,
639 const double J = this->dshape_eulerian_at_knot(
ipt,
psi,
dpsidx);
661 template<
unsigned NNODE_1D>
663 :
public virtual QElement<1, NNODE_1D>
679 template<
class FVK_ELEMENT>
694 <<
"Foeppl von Karman elements only store eight fields so fld must be"
695 <<
"0 to 7 rather than " <<
fld << std::endl;
706 for (
unsigned j = 0;
j <
nnod;
j++)
731 <<
"Foeppl von Karman elements only store eight fields so fld must be"
732 <<
"0 to 7 rather than " <<
fld << std::endl;
759 <<
"Foeppl von Karman elements only store eight fields so fld must be"
760 <<
"0 to 7 rather than " <<
fld << std::endl;
786 <<
"Foeppl von Karman elements only store eight fields so fld must be"
787 <<
"0 to 7 rather than " <<
fld << std::endl;
822 <<
"Foeppl von Karman elements only store eight fields so fld must be"
823 <<
"0 to 7 rather than " <<
fld << std::endl;
828 return this->
nnode();
840 <<
"Foeppl von Karman elements only store eight fields so fld must be"
841 <<
"0 to 7 rather than " <<
fld << std::endl;
855 template<
class ELEMENT>
868 template<
class ELEMENT>
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
A class that represents a collection of data; each Data object may contain many different individual ...
void pin(const unsigned &i)
Pin the i-th stored variable.
unsigned ntstorage() const
Return total number of doubles stored per value to record time history of each value (one for steady ...
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 ...
A general Finite Element class.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
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...
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-...
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.
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
unsigned nnode() const
Return the number of nodes.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
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...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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.
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
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.
double *& eta_pt()
Pointer to eta.
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.
const double & eta() const
Eta.
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....
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.
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.
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!...
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...
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).