27 #ifndef OOMPH_FOEPPLVONKARMAN_DISPLACEMENT_ELEMENTS_HEADER 
   28 #define OOMPH_FOEPPLVONKARMAN_DISPLACEMENT_ELEMENTS_HEADER 
   32 #include <oomph-lib-config.h> 
   38 #include "../generic/projection.h" 
   39 #include "../generic/nodes.h" 
   40 #include "../generic/oomph_utilities.h" 
   91     const double& 
nu()
 const 
  132       const unsigned n_plot = 5;
 
  138     void output(std::ostream& outfile, 
const unsigned& n_plot);
 
  143       const unsigned n_plot = 5;
 
  149     void output(FILE* file_pt, 
const unsigned& n_plot);
 
  153                     const unsigned& n_plot,
 
  160       std::ostream& outfile,
 
  161       const unsigned& n_plot,
 
  166         "There is no time-dependent output_fct() for Foeppl von Karman" 
  168         OOMPH_CURRENT_FUNCTION,
 
  169         OOMPH_EXCEPTION_LOCATION);
 
  188         "There is no time-dependent compute_error() for Foeppl von Karman" 
  190         OOMPH_CURRENT_FUNCTION,
 
  191         OOMPH_EXCEPTION_LOCATION);
 
  224                                          double& pressure)
 const 
  234         (*Pressure_fct_pt)(x, pressure);
 
  251         (*Traction_fct_pt)(x, traction);
 
  260       const unsigned n_node = 
nnode();
 
  278       for (
unsigned j = 0; j < 2; j++)
 
  284       for (
unsigned l = 0; l < n_node; l++)
 
  287         for (
unsigned j = 0; j < 2; j++)
 
  289           gradient[j] += this->
nodal_value(l, w_nodal_index) * dpsidx(l, j);
 
  302                                const unsigned& index)
 const 
  305       const unsigned n_node = 
nnode();
 
  323       for (
unsigned j = 0; j < 2; j++)
 
  329       for (
unsigned l = 0; l < n_node; l++)
 
  332         for (
unsigned j = 0; j < 2; j++)
 
  334           gradient[j] += this->
nodal_value(l, w_nodal_index) * dpsidx(l, j);
 
  347       double e_xx = interpolated_duxdx[0];
 
  348       double e_yy = interpolated_duydx[1];
 
  349       double e_xy = 0.5 * (interpolated_duxdx[1] + interpolated_duydx[0]);
 
  350       e_xx += 0.5 * interpolated_dwdx[0] * interpolated_dwdx[0];
 
  351       e_yy += 0.5 * interpolated_dwdx[1] * interpolated_dwdx[1];
 
  352       e_xy += 0.5 * interpolated_dwdx[0] * interpolated_dwdx[1];
 
  355       const double _nu = 
nu();
 
  357       double inv_denom = 1.0 / (1.0 - _nu * _nu);
 
  361       sigma(0, 0) = (e_xx + _nu * e_yy) * inv_denom;
 
  364       sigma(1, 1) = (e_yy + _nu * e_xx) * inv_denom;
 
  367       sigma(0, 1) = sigma(1, 0) = e_xy / (1.0 + _nu);
 
  376       const unsigned n_node = 
nnode();
 
  406       for (
unsigned l = 0; l < n_node; l++)
 
  416         for (
unsigned j = 0; j < 2; j++)
 
  418           interpolated_dwdx[j] += 
nodal_value[0] * dpsidx(l, j);
 
  419           interpolated_duxdx[j] += 
nodal_value[2] * dpsidx(l, j);
 
  420           interpolated_duydx[j] += 
nodal_value[3] * dpsidx(l, j);
 
  429         sigma, interpolated_dwdx, interpolated_duxdx, interpolated_duydx);
 
  434       strain(0, 0) = interpolated_duxdx[0];
 
  435       strain(0, 0) += 0.5 * interpolated_dwdx[0] * interpolated_dwdx[0];
 
  438       strain(1, 1) = interpolated_duydx[1];
 
  439       strain(1, 1) += 0.5 * interpolated_dwdx[1] * interpolated_dwdx[1];
 
  442       strain(0, 1) = 0.5 * (interpolated_duxdx[1] + interpolated_duydx[0]);
 
  443       strain(0, 1) += 0.5 * interpolated_dwdx[0] * interpolated_dwdx[1];
 
  446       strain(1, 0) = strain(0, 1);
 
  462       for (
unsigned i = 0; 
i < 
ndof; 
i++)
 
  464         mass_matrix(
i, 
i) += 1.0;
 
  473       const unsigned n_node = 
nnode();
 
  476       Shape psi(n_node), test(n_node);
 
  477       DShape dpsidx(n_node, 2), dtestdx(n_node, 2);
 
  492       for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
 
  499           ipt, psi, dpsidx, test, dtestdx);
 
  508         double interpolated_laplacian_w = 0;
 
  519         for (
unsigned l = 0; l < n_node; l++)
 
  528           interpolated_laplacian_w += 
nodal_value[1] * psi(l);
 
  531           for (
unsigned j = 0; j < 2; j++)
 
  534             interpolated_dwdx[j] += 
nodal_value[0] * dpsidx(l, j);
 
  535             interpolated_dlaplacian_wdx[j] += 
nodal_value[1] * dpsidx(l, j);
 
  536             interpolated_duxdx[j] += 
nodal_value[2] * dpsidx(l, j);
 
  537             interpolated_duydx[j] += 
nodal_value[3] * dpsidx(l, j);
 
  552             sigma, interpolated_dwdx, interpolated_duxdx, interpolated_duydx);
 
  557         double pressure = 0.0;
 
  564         for (
unsigned l = 0; l < n_node; l++)
 
  572             residuals[local_eqn] += pressure * test(l) * 
W;
 
  575             for (
unsigned k = 0; k < 2; k++)
 
  577               residuals[local_eqn] +=
 
  578                 interpolated_dlaplacian_wdx[k] * dtestdx(l, k) * 
W;
 
  585               for (
unsigned alpha = 0; alpha < 2; alpha++)
 
  588                 for (
unsigned beta = 0; beta < 2; beta++)
 
  590                   residuals[local_eqn] -= 
eta() * sigma(alpha, beta) *
 
  591                                           interpolated_dwdx[alpha] *
 
  592                                           dtestdx(l, beta) * 
W;
 
  605             residuals[local_eqn] += interpolated_laplacian_w * test(l) * 
W;
 
  607             for (
unsigned k = 0; k < 2; k++)
 
  609               residuals[local_eqn] += interpolated_dwdx[k] * dtestdx(l, k) * 
W;
 
  624             residuals[local_eqn] += traction[0] * test(l) * 
W;
 
  627             for (
unsigned beta = 0; beta < 2; beta++)
 
  629               residuals[local_eqn] += sigma(0, beta) * dtestdx(l, beta) * 
W;
 
  640             residuals[local_eqn] += traction[1] * test(l) * 
W;
 
  643             for (
unsigned beta = 0; beta < 2; beta++)
 
  645               residuals[local_eqn] += sigma(1, beta) * dtestdx(l, beta) * 
W;
 
  659                                      unsigned index = 0)
 const 
  662       const unsigned n_node = 
nnode();
 
  674       double interpolated_w = 0.0;
 
  677       for (
unsigned l = 0; l < n_node; l++)
 
  679         interpolated_w += this->
nodal_value(l, w_nodal_index) * psi[l];
 
  682       return (interpolated_w);
 
  701       unsigned total_fvk_nodal_indices = 4;
 
  704       unsigned n_node = 
nnode();
 
  707       for (
unsigned index = first_fvk_nodal_index + 2; 
 
  709            index < first_fvk_nodal_index + total_fvk_nodal_indices;
 
  713         for (
unsigned inod = 0; inod < n_node; inod++)
 
  729                                                  DShape& dtestdx) 
const = 0;
 
  739       DShape& dtestdx) 
const = 0;
 
  770   template<
class FVK_ELEMENT>
 
  784         std::stringstream error_stream;
 
  786           << 
"Foeppl von Karman elements only store four fields so fld must be" 
  787           << 
"0 to 3 rather than " << fld << std::endl;
 
  789           error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
 
  794       unsigned nnod = this->
nnode();
 
  798       for (
unsigned j = 0; j < nnod; j++)
 
  801         data_values[j] = std::make_pair(this->
node_pt(j), fld);
 
  821         std::stringstream error_stream;
 
  823           << 
"Foeppl von Karman elements only store four fields so fld must be" 
  824           << 
"0 to 3 rather than " << fld << std::endl;
 
  826           error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
 
  849         std::stringstream error_stream;
 
  851           << 
"Foeppl von Karman elements only store four fields so fld must be" 
  852           << 
"0 to 3 rather than " << fld << std::endl;
 
  854           error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
 
  857       unsigned n_dim = this->
dim();
 
  858       unsigned n_node = this->
nnode();
 
  860       DShape dpsidx(n_node, n_dim), dtestdx(n_node, n_dim);
 
  862         this->dshape_and_dtest_eulerian_fvk(
s, psi, dpsidx, test, dtestdx);
 
  876         std::stringstream error_stream;
 
  878           << 
"Foeppl von Karman elements only store four fields so fld must be" 
  879           << 
"0 to 3 rather than " << fld << std::endl;
 
  881           error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
 
  885       unsigned w_nodal_index = this->nodal_index_fvk(fld);
 
  888       unsigned n_node = this->
nnode();
 
  895       double interpolated_w = 0.0;
 
  898       for (
unsigned l = 0; l < n_node; l++)
 
  900         interpolated_w += this->
nodal_value(t, l, w_nodal_index) * psi[l];
 
  902       return interpolated_w;
 
  912         std::stringstream error_stream;
 
  914           << 
"Foeppl von Karman elements only store four fields so fld must be" 
  915           << 
"0 to 3 rather than " << fld << std::endl;
 
  917           error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
 
  920       return this->
nnode();
 
  930         std::stringstream error_stream;
 
  932           << 
"Foeppl von Karman elements only store four fields so fld must be" 
  933           << 
"0 to 3 rather than " << fld << std::endl;
 
  935           error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
 
  938       const unsigned w_nodal_index = this->nodal_index_fvk(fld);
 
  947   template<
class ELEMENT>
 
  961   template<
class ELEMENT>
 
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
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 ...
unsigned long nrow() const
Return the number of rows of the matrix.
A class for all isoparametric elements that solve the Foeppl von Karman equations.
bool Linear_bending_model
Flag which stores whether we are using a linear, pure bending model instead of the full non-linear Fo...
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
void use_linear_bending_model()
Sets a flag to signify that we are solving the linear, pure bending equations, and pin all the nodal ...
DisplacementBasedFoepplvonKarmanEquations(const DisplacementBasedFoepplvonKarmanEquations &dummy)=delete
Broken copy constructor.
unsigned self_test()
Self-test: Return 0 for OK.
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.
void get_stress_and_strain_for_output(const Vector< double > &s, DenseMatrix< double > &sigma, DenseMatrix< double > &strain)
const double & nu() const
Poisson's ratio.
FoepplvonKarmanPressureFctPt & pressure_fct_pt()
Access function: Pointer to pressure function.
double * Nu_pt
Pointer to global Poisson's ratio.
FoepplvonKarmanTractionFctPt traction_fct_pt() const
Access function: Pointer to in-plane traction function. Const version.
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 ...
DisplacementBasedFoepplvonKarmanEquations()
Constructor (must initialise the Pressure_fct_pt and the Traction_fct_pt. Also set physical parameter...
FoepplvonKarmanPressureFctPt Pressure_fct_pt
Pointer to pressure function:
void output(std::ostream &outfile)
Output with default number of plot points.
void get_gradient_of_field(const Vector< double > &s, Vector< double > &gradient, const unsigned &index) const
Get gradient of field: gradient[i] = d[.]/dx_i,.
FoepplvonKarmanPressureFctPt pressure_fct_pt() const
Access function: Pointer to pressure function. Const version.
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 ...
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...
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 ...
void operator=(const DisplacementBasedFoepplvonKarmanEquations &)=delete
Broken assignment operator.
void get_gradient_of_deflection(const Vector< double > &s, Vector< double > &gradient) const
Get gradient of deflection: gradient[i] = dw/dx_i.
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
hierher dummy
void(* FoepplvonKarmanPressureFctPt)(const Vector< double > &x, double &f)
Function pointer to pressure function fct(x,f(x)) – x is a Vector!
double *& nu_pt()
Pointer to Poisson's ratio.
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
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,...
void(* FoepplvonKarmanTractionFctPt)(const Vector< double > &x, Vector< double > &f)
Function pointer to traction function fct(x,f(x)) – x is a Vector!
double * Eta_pt
Pointer to global eta.
void output(FILE *file_pt)
C_style output with default number of plot points.
void get_sigma(DenseMatrix< double > &sigma, const Vector< double > &interpolated_dwdx, const Vector< double > &interpolated_duxdx, const Vector< double > &interpolated_duydx)
Get the in-plane stress (sigma) as a fct of the pre=computed displcement derivatives.
const double & eta() const
Eta.
FoepplvonKarmanTractionFctPt & traction_fct_pt()
Access function: Pointer to in-plane traction function.
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...
static double Default_Nu_Value
Default value for Poisson's ratio.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals with this element's contribution.
FoepplvonKarmanTractionFctPt Traction_fct_pt
Pointer to traction function:
double *& eta_pt()
Pointer to eta.
static double Default_Physical_Constant_Value
Default value for physical constants.
virtual void get_traction_fvk(Vector< double > &x, Vector< double > &traction) const
Get in-plane traction term at (Eulerian) position x.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
A general Finite Element class.
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...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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 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...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
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 .
Integral *const  & integral_pt() const
Return the pointer to the integration scheme (const version)
double nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. If the node is hanging, the appropriate interpolation is ...
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...
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function  as .
unsigned ndof() const
Return the number of equations/dofs in the element.
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....
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.
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of field fld of node j.
unsigned nfields_for_projection()
Number of fields to be projected: Just two.
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (Note: count includes current value!)
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. (Note: count includes current value!...
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 nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
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 ...
Wrapper class for projectable elements. Adds "projectability" to the underlying ELEMENT.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
//////////////////////////////////////////////////////////////////// ////////////////////////////////...