31 #ifndef OOMPH_TFOEPPLVONKARMAN_DISPLACEMENT_ELEMENTS_HEADER 32 #define OOMPH_TFOEPPLVONKARMAN_DISPLACEMENT_ELEMENTS_HEADER 37 #include <oomph-lib-config.h> 42 #include "../generic/nodes.h" 43 #include "../generic/oomph_utilities.h" 44 #include "../generic/Telements.h" 45 #include "../generic/error_estimator.h" 71 template <
unsigned NNODE_1D>
132 std::list<std::pair<unsigned long,unsigned> >& dof_lookup_list)
const 135 const unsigned n_node = this->
nnode();
138 std::pair<unsigned,unsigned> dof_lookup;
141 for (
unsigned n = 0; n < n_node; n++)
151 if (local_eqn_number >= 0)
155 dof_lookup.first = this->
eqn_number(local_eqn_number);
158 dof_lookup.second = 2;
161 dof_lookup_list.push_front(dof_lookup);
172 if (local_eqn_number >= 0)
176 dof_lookup.first = this->
eqn_number(local_eqn_number);
179 if (
node_pt(n)->is_on_boundary(0) ||
node_pt(n)->is_on_boundary(1))
181 dof_lookup.second = 1;
186 dof_lookup.second = 0;
190 dof_lookup_list.push_front(dof_lookup);
201 if (local_eqn_number >= 0)
205 dof_lookup.first = this->
eqn_number(local_eqn_number);
208 dof_lookup.second = 3;
211 dof_lookup_list.push_front(dof_lookup);
222 if (local_eqn_number >= 0)
226 dof_lookup.first = this->
eqn_number(local_eqn_number);
229 dof_lookup.second = 4;
232 dof_lookup_list.push_front(dof_lookup);
248 void output(std::ostream &outfile,
const unsigned &n_plot)
264 void output(FILE* file_pt,
const unsigned &n_plot)
272 void output_fct(std::ostream &outfile,
const unsigned &n_plot,
281 void output_fct(std::ostream &outfile,
const unsigned &n_plot,
344 template<
unsigned NNODE_1D>
352 unsigned n_node = this->
nnode();
359 for(
unsigned i=0;
i<n_node;
i++)
362 dtestdx(
i,0) = dpsidx(
i,0);
363 dtestdx(
i,1) = dpsidx(
i,1);
378 template<
unsigned NNODE_1D>
407 template<
unsigned NNODE_1D>
TDisplacementBasedFoepplvonKarmanElement()
Constructor: Call constructors for TElement and Foeppl von Karman equations.
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
void operator=(const TDisplacementBasedFoepplvonKarmanElement< NNODE_1D > &)
Broken assignment operator.
void output(std::ostream &outfile)
Output with default number of plot points.
unsigned nvertex_node() const
Number of vertex nodes in the element.
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.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
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-...
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from FvK equations.
unsigned ndof_types() const
The number of dof types that degrees of freedom in this element are sub-divided into.
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 function for an exact solution: x,y,w_exact.
void output(FILE *file_pt)
C-style output function: x,y,w.
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as ...
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: x,y,w at n_plot^2 plot points.
TDisplacementBasedFoepplvonKarmanElement(const TDisplacementBasedFoepplvonKarmanElement< NNODE_1D > &dummy)
Broken copy constructor.
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...
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.
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...
void output(std::ostream &outfile)
Output function: x,y,w.
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned > > &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number...
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
static const unsigned Initial_Nvalue
Static unsigned that holds the (same) number of variables at every node.
void broken_assign(const std::string &class_name)
Issue error message and terminate execution.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional TElement. ...
int local_eqn_number(const unsigned long &ieqn_global) const
Return the local equation number corresponding to the ieqn_global-th global equation number...
unsigned nnode() const
Return the number of nodes.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: x,y,w at n_plot^2 plot points.
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) ...
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...