27 #ifndef OOMPH_TFOEPPLVONKARMAN_ELEMENTS_HEADER
28 #define OOMPH_TFOEPPLVONKARMAN_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
38 #include "../generic/nodes.h"
39 #include "../generic/oomph_utilities.h"
40 #include "../generic/Telements.h"
41 #include "../generic/error_estimator.h"
60 template<
unsigned NNODE_1D>
97 void output(std::ostream& outfile,
const unsigned& n_plot)
113 void output(FILE* file_pt,
const unsigned& n_plot)
122 const unsigned& n_plot,
132 const unsigned& n_plot,
137 outfile, n_plot, time, exact_soln_pt);
162 return (NNODE_1D - 1);
204 template<
unsigned NNODE_1D>
212 unsigned n_node = this->nnode();
215 double J = this->dshape_eulerian(
s, psi, dpsidx);
219 for (
unsigned i = 0;
i < n_node;
i++)
222 dtestdx(
i, 0) = dpsidx(
i, 0);
223 dtestdx(
i, 1) = dpsidx(
i, 1);
237 template<
unsigned NNODE_1D>
239 NNODE_1D>::dshape_and_dtest_eulerian_at_knot_fvk(
const unsigned& ipt,
246 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
263 template<
unsigned NNODE_1D>
265 :
public virtual TElement<1, NNODE_1D>
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
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 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.
void output(std::ostream &outfile)
Output with default number of plot points.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
//////////////////////////////////////////////////////////////////////
TFoepplvonKarmanElement()
Constructor: Call constructors for TElement and Foeppl von Karman equations.
unsigned required_nvalue(const unsigned &n) const
Access function for Nvalue: # of ‘values’ (pinned or dofs) at node n (always returns the same value a...
static const unsigned Initial_Nvalue
Static unsigned that holds the (same) number of variables at every node.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
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. Return Jacobian.
TFoepplvonKarmanElement(const TFoepplvonKarmanElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
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 (calls the steady version)
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
void operator=(const TFoepplvonKarmanElement< NNODE_1D > &)=delete
Broken assignment operator.
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.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from FvK equations.
void output(std::ostream &outfile)
Output function: x,y,w.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: x,y,w at n_plot^2 plot points.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: x,y,w at n_plot^2 plot points.
void output(FILE *file_pt)
C-style output function: x,y,w.
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.
unsigned nvertex_node() const
Number of vertex nodes in the element.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...