27 #ifndef OOMPH_AXISYM_FOEPPLVONKARMAN_ELEMENTS_HEADER
28 #define OOMPH_AXISYM_FOEPPLVONKARMAN_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
36 #include "../generic/nodes.h"
37 #include "../generic/Qelements.h"
38 #include "../generic/oomph_utilities.h"
39 #include "../generic/element_with_external_element.h"
50 class AxisymFoepplvonKarmanEquations :
public virtual FiniteElement
77 const double&
eta()
const
107 const unsigned n_plot = 5;
113 void output(std::ostream& outfile,
const unsigned& n_plot);
118 const unsigned n_plot = 5;
124 void output(FILE* file_pt,
const unsigned& n_plot);
128 const unsigned& n_plot,
135 std::ostream& outfile,
136 const unsigned& n_plot,
141 "There is no time-dependent output_fct() for Foeppl von Karman"
143 OOMPH_CURRENT_FUNCTION,
144 OOMPH_EXCEPTION_LOCATION);
162 "There is no time-dependent compute_error() for Foeppl von Karman"
164 OOMPH_CURRENT_FUNCTION,
165 OOMPH_EXCEPTION_LOCATION);
198 double& pressure)
const
208 (*Pressure_fct_pt)(r, pressure);
218 double& airy_forcing)
const
228 (*Airy_forcing_fct_pt)(r, airy_forcing);
237 const unsigned n_node =
nnode();
253 for (
unsigned l = 0; l < n_node; l++)
255 gradient[0] += this->
nodal_value(l, w_nodal_index) * dpsidr(l, 0);
272 const unsigned n_node =
nnode();
284 double interpolated_w = 0.0;
287 for (
unsigned l = 0; l < n_node; l++)
289 interpolated_w += this->
nodal_value(l, w_nodal_index) * psi[l];
292 return (interpolated_w);
299 double& sigma_phi_phi);
318 unsigned total_fvk_nodal_indices = 6;
321 unsigned n_node =
nnode();
324 for (
unsigned index = first_fvk_nodal_index + 2;
325 index < first_fvk_nodal_index + total_fvk_nodal_indices;
329 for (
unsigned inod = 0; inod < n_node; inod++)
346 DShape& dtestdr)
const = 0;
356 DShape& dtestdr)
const = 0;
386 template<
unsigned NNODE_1D>
388 :
public virtual QElement<1, NNODE_1D>,
428 void output(std::ostream& outfile,
const unsigned& n_plot)
442 void output(FILE* file_pt,
const unsigned& n_plot)
450 const unsigned& n_plot,
454 outfile, n_plot, exact_soln_pt);
461 const unsigned& n_plot,
466 outfile, n_plot, time, exact_soln_pt);
498 template<
unsigned NNODE_1D>
500 NNODE_1D>::dshape_and_dtest_eulerian_axisym_fvk(
const Vector<double>&
s,
508 const double J = this->dshape_eulerian(
s, psi, dpsidr);
525 template<
unsigned NNODE_1D>
531 DShape& dtestdr)
const
534 const double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidr);
556 template<
unsigned NNODE_1D,
class FLUID_ELEMENT>
576 const double&
q()
const
592 return this->
nnode();
607 const unsigned t = 0;
619 const unsigned n_node = this->
nnode();
628 this->
shape(zeta, psi);
631 double interpolated_w = 0.0;
632 double interpolated_r = 0.0;
635 for (
unsigned l = 0; l < n_node; l++)
637 interpolated_w += this->
nodal_value(t, l, w_nodal_index) * psi[l];
638 interpolated_r += this->
node_pt(l)->
x(
t, 0) * psi[l];
642 r[0] = interpolated_r;
643 r[1] = interpolated_w;
653 const unsigned n_node = this->
nnode();
662 this->
shape(zeta, psi);
665 double interpolated_dwdt = 0.0;
666 double interpolated_drdt = 0.0;
669 for (
unsigned l = 0; l < n_node; l++)
681 for (
unsigned t = 0;
t < n_time;
t++)
692 drdt[0] = interpolated_drdt;
693 drdt[1] = interpolated_dwdt;
703 double& pressure)
const
711 const double q_value =
q();
714 FLUID_ELEMENT* ext_el_pt =
726 ext_el_pt->traction(s_ext, normal, traction);
729 pressure -= q_value * traction[1];
737 unsigned nnod = this->
nnode();
745 outfile <<
"ZONE I=" << n_intpt << std::endl;
746 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
754 double interpolated_w = 0.0;
755 double interpolated_r = 0.0;
758 for (
unsigned l = 0; l < nnod; l++)
760 interpolated_w += this->
nodal_value(l, w_nodal_index) * psi[l];
761 interpolated_r += this->
node_pt(l)->
x(0) * psi[l];
765 FLUID_ELEMENT* ext_el_pt =
771 ext_el_pt->interpolated_u_axi_nst(s_ext, veloc);
773 ext_el_pt->interpolated_x(s_ext, x);
775 outfile << interpolated_r <<
" " << interpolated_w <<
" " << veloc[0]
776 <<
" " << veloc[1] <<
" " << x[0] <<
" " << x[1] <<
" "
784 const unsigned& nplot)
788 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
791 FLUID_ELEMENT* ext_el_pt =
795 ext_el_pt->output(outfile, nplot);
835 std::map<FLUID_ELEMENT*, bool> done;
841 for (
unsigned iint = 0; iint < n_intpt; iint++)
844 FLUID_ELEMENT* el_f_pt =
854 el_f_pt->node_update();
855 done[el_f_pt] =
true;
864 void output(std::ostream& outfile,
const unsigned& n_plot)
874 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
884 double sigma_r_r = 0.0;
885 double sigma_phi_phi = 0.0;
891 << drdt[1] <<
" " << sigma_r_r <<
" " << sigma_phi_phi
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: r,w at n_plot plot points.
double dshape_and_dtest_eulerian_axisym_fvk(const Vector< double > &s, Shape &psi, DShape &dpsidr, Shape &test, DShape &dtestdr) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
void output(FILE *file_pt)
C-style output function: r,w.
double dshape_and_dtest_eulerian_at_knot_axisym_fvk(const unsigned &ipt, Shape &psi, DShape &dpsidr, Shape &test, DShape &dtestdr) const
Shape, test functions & derivs. w.r.t. to global coords. at integration point ipt....
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. r,w_exact at n_plot plot points (Calls the stead...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: r,w,sigma_r_r,sigma_phi_phi at n_plot plot points.
AxisymFoepplvonKarmanElement(const AxisymFoepplvonKarmanElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile)
Output function: r,w,sigma_r_r,sigma_phi_phi.
static const unsigned Initial_Nvalue
Static int that holds the number of variables at nodes: always the same.
unsigned required_nvalue(const unsigned &n) const
Required # of ‘values’ (pinned or dofs) at node n.
void operator=(const AxisymFoepplvonKarmanElement< NNODE_1D > &)=delete
Broken assignment operator.
AxisymFoepplvonKarmanElement()
Constructor: Call constructors for QElement and AxisymFoepplvonKarmanEquations.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: r,w_exact at n_plot plot points.
A class for all isoparametric elements that solve the axisYm Foeppl von Karman equations in a displac...
AxisymFoepplvonKarmanPressureFctPt pressure_fct_pt() const
Access function: Pointer to pressure function. Const version.
AxisymFoepplvonKarmanPressureFctPt Pressure_fct_pt
Pointer to pressure function:
static double Default_Physical_Constant_Value
Default value for physical constants.
void operator=(const AxisymFoepplvonKarmanEquations &)=delete
Broken assignment operator.
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,w_exact at n_plot plot points (dummy time-dependent version to keep intel compil...
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.
bool interpolated_stress(const Vector< double > &s, double &sigma_r_r, double &sigma_phi_phi) const
Compute in-plane stresses. Return boolean to indicate success (false if attempt to evaluate stresses ...
void output(FILE *file_pt)
C_style output with default number of plot points.
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
unsigned self_test()
Self-test: Return 0 for OK.
AxisymFoepplvonKarmanPressureFctPt Airy_forcing_fct_pt
Pointer to Airy forcing function.
double * Eta_pt
Pointer to FvK parameter.
virtual double dshape_and_dtest_eulerian_at_knot_axisym_fvk(const unsigned &ipt, Shape &psi, DShape &dpsidr, Shape &test, DShape &dtestdr) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
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 ...
void(* AxisymFoepplvonKarmanPressureFctPt)(const double &r, double &f)
Function pointer to pressure function fct(r,f(r)) – r is a Vector!
void output(FILE *file_pt, const unsigned &n_plot)
C-style output FE representation of soln: r,w at n_plot plot points.
void get_gradient_of_deflection(const Vector< double > &s, Vector< double > &gradient) const
Get gradient of deflection: gradient[i] = dw/dr_i */.
AxisymFoepplvonKarmanPressureFctPt & airy_forcing_fct_pt()
Access function: Pointer to Airy forcing function.
AxisymFoepplvonKarmanPressureFctPt & pressure_fct_pt()
Access function: Pointer to pressure function.
void output(std::ostream &outfile, const unsigned &n_plot)
Output FE representation of soln: r,w,sigma_r_r,sigma_phi_phi at n_plot plot points.
AxisymFoepplvonKarmanEquations()
Constructor (must initialise the Pressure_fct_pt and Airy_forcing_fct_pt to null)....
virtual double dshape_and_dtest_eulerian_axisym_fvk(const Vector< double > &s, Shape &psi, DShape &dpsidr, Shape &test, DShape &dtestdr) const =0
Shape/test functions and derivs w.r.t. to global coords at local coord. s; return Jacobian of mapping...
bool Linear_bending_model
Flag which stores whether we are using a linear, pure bending model instead of the full non-linear Fo...
const double & eta() const
FvK parameter.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,w_exact at n_plot plot points.
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 ...
double interpolated_w_fvk(const Vector< double > &s) const
Return FE representation of vertical displacement, w_fvk(s) at local coordinate s.
virtual void get_pressure_fvk(const unsigned &ipt, const double &r, double &pressure) const
Get pressure term at (Eulerian) position r. This function is virtual to allow overloading in multi-ph...
double *& eta_pt()
Pointer to FvK parameter.
AxisymFoepplvonKarmanEquations(const AxisymFoepplvonKarmanEquations &dummy)=delete
Broken copy constructor.
virtual void get_airy_forcing_fvk(const unsigned &ipt, const double &r, double &airy_forcing) const
Get Airy forcing term at (Eulerian) position r. This function is virtual to allow overloading in mult...
AxisymFoepplvonKarmanPressureFctPt airy_forcing_fct_pt() const
Access function: Pointer to Airy forcing function. Const version.
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.
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
This is a base class for all elements that require external sources (e.g. FSI, multi-domain problems ...
Vector< double > & external_element_local_coord(const unsigned &interaction_index, const unsigned &ipt)
Access function to get source element's local coords for specified interaction index at specified int...
void set_ninteraction(const unsigned &n_interaction)
Set the number of interactions in the element This function is usually called in the specific element...
FiniteElement *& external_element_pt(const unsigned &interaction_index, const unsigned &ipt)
Access function to source element for specified interaction index at specified integration point.
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
void output(std::ostream &outfile, const unsigned &n_plot)
Output FE representation of soln: r,w,dwdt,sigma_r_r,sigma_phi_phi at n_plot plot points.
void reset_after_external_interaction_geometric_fd()
Perform any auxiliary node update fcts of the adjacent fluid elements.
void node_update_adjacent_fluid_elements()
Update the nodal positions in all fluid elements that affect the traction on this element.
virtual ~FSIAxisymFoepplvonKarmanElement()
Empty virtual destructor.
double *& q_pt()
Return a pointer the ratio of stress scales used to non-dimensionalise the fluid and solid equations.
virtual Data * geom_data_pt(const unsigned &j)
Return pointer to the j-th Data item that the object's shape depends on.
void update_before_external_interaction_geometric_fd()
Perform any auxiliary node update fcts of the adjacent fluid elements.
double * Q_pt
Pointer to the ratio, , of the stress used to non-dimensionalise the fluid stresses to the stress us...
void output_adjacent_fluid_elements(std::ostream &outfile, const unsigned &nplot)
Output adjacent fluid elements (for checking of fsi setup)
void dposition_dt(const Vector< double > &zeta, const unsigned &j, Vector< double > &drdt)
j-th time-derivative on object at current time: .
void output_integration_points(std::ostream &outfile)
Output integration points (for checking of fsi setup)
virtual void get_pressure_fvk(const unsigned &ipt, const double &r, double &pressure) const
Overload pressure term at (Eulerian) position r. Adds fluid traction to pressure imposed by "pressure...
void position(const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const
Overloaded position function: Return 2D position vector: (r(zeta),z(zeta)) of material point whose "L...
void position(const Vector< double > &zeta, Vector< double > &r) const
Overloaded position function: Return 2D position vector: (r(zeta),z(zeta)) of material point whose "L...
const double & q() const
Return the ratio of the stress scales used to non-dimensionalise the fluid and elasticity equations.
void update_in_external_interaction_geometric_fd(const unsigned &i)
Perform any auxiliary node update fcts of the adjacent fluid elements.
virtual unsigned ngeom_data() const
How many items of Data does the shape of the object depend on? All nodal data.
FSIAxisymFoepplvonKarmanElement()
Constructor.
void reset_in_external_interaction_geometric_fd(const unsigned &i)
Perform any auxiliary node update fcts of the adjacent fluid elements.
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.
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)
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
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 .
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
double & x(const unsigned &i)
Return the i-th nodal coordinate.
An OomphLibError object which should be thrown when an run-time error is encountered....
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
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)
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
bool is_steady() const
Flag to indicate if a timestepper has been made steady (possibly temporarily to switch off time-depen...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...