27 #ifndef OOMPH_FOEPPLVONKARMAN_ELEMENTS_HEADER
28 #define OOMPH_FOEPPLVONKARMAN_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
39 #include "../generic/projection.h"
40 #include "../generic/nodes.h"
41 #include "../generic/Qelements.h"
42 #include "../generic/oomph_utilities.h"
86 const double&
eta()
const
105 if (data_pt->
nvalue() != 1)
107 throw OomphLibError(
"Data object that contains volume control pressure "
108 "should only contain a single value. ",
139 const unsigned n_plot = 5;
145 void output(std::ostream& outfile,
const unsigned& n_plot);
150 const unsigned n_plot = 5;
156 void output(FILE* file_pt,
const unsigned& n_plot);
160 const unsigned& n_plot,
167 std::ostream& outfile,
168 const unsigned& n_plot,
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);
251 double& airy_forcing)
const
261 (*Airy_forcing_fct_pt)(
x, airy_forcing);
270 const unsigned n_node =
nnode();
283 for (
unsigned j = 0;
j < 2;
j++)
289 for (
unsigned l = 0; l < n_node; l++)
292 for (
unsigned j = 0;
j < 2;
j++)
294 gradient[
j] += this->
nodal_value(l, w_nodal_index) * dpsidx(l,
j);
309 unsigned index = 0)
const
312 const unsigned n_node =
nnode();
324 double interpolated_w = 0.0;
327 for (
unsigned l = 0; l < n_node; l++)
329 interpolated_w += this->
nodal_value(l, w_nodal_index) * psi[l];
332 return (interpolated_w);
348 const unsigned n_node =
nnode();
359 double integral_w = 0;
362 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
374 double interpolated_w = 0;
375 double w_nodal_value;
378 for (
unsigned l = 0; l < n_node; l++)
383 interpolated_w += w_nodal_value * psi(l);
387 integral_w += interpolated_w *
W;
410 unsigned total_fvk_nodal_indicies = 8;
413 unsigned n_node =
nnode();
416 for (
unsigned index = first_fvk_nodal_index + 2;
417 index < first_fvk_nodal_index + total_fvk_nodal_indicies;
421 for (
unsigned inod = 0; inod < n_node; inod++)
437 DShape& dtestdx)
const = 0;
447 DShape& dtestdx)
const = 0;
489 template<
unsigned NNODE_1D>
530 void output(std::ostream& outfile,
const unsigned& n_plot)
546 void output(FILE* file_pt,
const unsigned& n_plot)
555 const unsigned& n_plot,
566 const unsigned& n_plot,
571 outfile, n_plot, time, exact_soln_pt);
604 template<
unsigned NNODE_1D>
613 const double J = this->dshape_eulerian(
s, psi, dpsidx);
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>
692 std::stringstream error_stream;
694 <<
"Foeppl von Karman elements only store eight fields so fld must be"
695 <<
"0 to 7 rather than " << fld << std::endl;
702 unsigned nnod = this->
nnode();
706 for (
unsigned j = 0;
j < nnod;
j++)
709 data_values[
j] = std::make_pair(this->
node_pt(
j), fld);
729 std::stringstream error_stream;
731 <<
"Foeppl von Karman elements only store eight fields so fld must be"
732 <<
"0 to 7 rather than " << fld << std::endl;
757 std::stringstream error_stream;
759 <<
"Foeppl von Karman elements only store eight fields so fld must be"
760 <<
"0 to 7 rather than " << fld << std::endl;
765 unsigned n_dim = this->
dim();
766 unsigned n_node = this->
nnode();
768 DShape dpsidx(n_node, n_dim), dtestdx(n_node, n_dim);
770 this->dshape_and_dtest_eulerian_fvk(
s, psi, dpsidx,
test, dtestdx);
784 std::stringstream error_stream;
786 <<
"Foeppl von Karman elements only store eight fields so fld must be"
787 <<
"0 to 7 rather than " << fld << std::endl;
793 unsigned w_nodal_index = this->nodal_index_fvk(fld);
796 unsigned n_node = this->
nnode();
803 double interpolated_w = 0.0;
806 for (
unsigned l = 0; l < n_node; l++)
808 interpolated_w += this->
nodal_value(t, l, w_nodal_index) * psi[l];
810 return interpolated_w;
820 std::stringstream error_stream;
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();
838 std::stringstream error_stream;
840 <<
"Foeppl von Karman elements only store eight fields so fld must be"
841 <<
"0 to 7 rather than " << fld << std::endl;
846 const unsigned w_nodal_index = this->nodal_index_fvk(fld);
855 template<
class ELEMENT>
868 template<
class ELEMENT>
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
RowVector3d w
Definition: Matrix_resize_int.cpp:3
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
Definition: nodes.h:483
unsigned ntstorage() const
Definition: nodes.cc:879
FaceGeometry()
Definition: foeppl_von_karman_elements.h:873
FaceGeometry()
Definition: foeppl_von_karman_elements.h:860
FaceGeometry()
Definition: foeppl_von_karman_elements.h:668
Definition: elements.h:4998
Definition: elements.h:1313
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual void shape(const Vector< double > &s, Shape &psi) const =0
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned dim() const
Definition: elements.h:2611
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3325
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2576
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Definition: elements.h:1765
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3298
Definition: foeppl_von_karman_elements.h:57
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: foeppl_von_karman_elements.cc:547
void output(FILE *file_pt)
C_style output with default number of plot points.
Definition: foeppl_von_karman_elements.h:148
double interpolated_w_fvk(const Vector< double > &s, unsigned index=0) const
Definition: foeppl_von_karman_elements.h:308
FoepplvonKarmanPressureFctPt pressure_fct_pt() const
Access function: Pointer to pressure function. Const version.
Definition: foeppl_von_karman_elements.h:208
void get_gradient_of_deflection(const Vector< double > &s, Vector< double > &gradient) const
Get gradient of deflection: gradient[i] = dw/dx_i.
Definition: foeppl_von_karman_elements.h:266
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.
Definition: foeppl_von_karman_elements.cc:496
FoepplvonKarmanPressureFctPt airy_forcing_fct_pt() const
Access function: Pointer to Airy forcing function. Const version.
Definition: foeppl_von_karman_elements.h:220
void(* FoepplvonKarmanPressureFctPt)(const Vector< double > &x, double &f)
Definition: foeppl_von_karman_elements.h:61
double *& eta_pt()
Pointer to eta.
Definition: foeppl_von_karman_elements.h:92
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: foeppl_von_karman_elements.h:166
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: foeppl_von_karman_elements.h:137
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals with this element's contribution.
Definition: foeppl_von_karman_elements.cc:50
FoepplvonKarmanPressureFctPt & pressure_fct_pt()
Access function: Pointer to pressure function.
Definition: foeppl_von_karman_elements.h:202
virtual void get_airy_forcing_fvk(const unsigned &ipt, const Vector< double > &x, double &airy_forcing) const
Definition: foeppl_von_karman_elements.h:249
FoepplvonKarmanEquations()
Definition: foeppl_von_karman_elements.h:67
virtual double dshape_and_dtest_eulerian_fvk(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
int Volume_constraint_pressure_external_data_index
Definition: foeppl_von_karman_elements.h:472
void use_linear_bending_model()
Definition: foeppl_von_karman_elements.h:400
virtual double get_bounded_volume() const
Definition: foeppl_von_karman_elements.h:345
bool Linear_bending_model
Definition: foeppl_von_karman_elements.h:466
virtual void get_pressure_fvk(const unsigned &ipt, const Vector< double > &x, double &pressure) const
Definition: foeppl_von_karman_elements.h:229
double * Eta_pt
Pointer to global eta.
Definition: foeppl_von_karman_elements.h:452
virtual unsigned nodal_index_fvk(const unsigned &i=0) const
Definition: foeppl_von_karman_elements.h:131
unsigned self_test()
Self-test: Return 0 for OK.
Definition: foeppl_von_karman_elements.cc:349
FoepplvonKarmanPressureFctPt & airy_forcing_fct_pt()
Access function: Pointer to Airy forcing function.
Definition: foeppl_von_karman_elements.h:214
virtual double dshape_and_dtest_eulerian_at_knot_fvk(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
const double & eta() const
Eta.
Definition: foeppl_von_karman_elements.h:86
FoepplvonKarmanPressureFctPt Pressure_fct_pt
Pointer to pressure function:
Definition: foeppl_von_karman_elements.h:455
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.
Definition: foeppl_von_karman_elements.h:188
static double Default_Physical_Constant_Value
Default value for physical constants.
Definition: foeppl_von_karman_elements.h:462
void interpolated_stress(const Vector< double > &s, double &sigma_xx, double &sigma_yy, double &sigma_xy)
Compute in-plane stresses.
Definition: foeppl_von_karman_elements.cc:374
FoepplvonKarmanPressureFctPt Airy_forcing_fct_pt
Definition: foeppl_von_karman_elements.h:458
void set_volume_constraint_pressure_data_as_external_data(Data *data_pt)
Definition: foeppl_von_karman_elements.h:102
void operator=(const FoepplvonKarmanEquations &)=delete
Broken assignment operator.
unsigned add_external_data(Data *const &data_pt, const bool &fd=true)
Definition: elements.cc:307
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.
Definition: nodes.h:1022
Definition: oomph_definitions.h:222
Definition: projection.h:183
Foeppl von Karman upgraded to become projectable.
Definition: foeppl_von_karman_elements.h:682
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: foeppl_von_karman_elements.h:750
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: foeppl_von_karman_elements.h:777
unsigned nhistory_values_for_coordinate_projection()
Definition: foeppl_von_karman_elements.h:743
unsigned nfields_for_projection()
Number of fields to be projected: Just two.
Definition: foeppl_von_karman_elements.h:717
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of field fld of node j.
Definition: foeppl_von_karman_elements.h:833
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
Definition: foeppl_von_karman_elements.h:815
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: foeppl_von_karman_elements.h:687
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: foeppl_von_karman_elements.h:724
Definition: Qelements.h:459
Definition: foeppl_von_karman_elements.h:492
void output(FILE *file_pt, const unsigned &n_plot)
Definition: foeppl_von_karman_elements.h:546
double dshape_and_dtest_eulerian_fvk(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: foeppl_von_karman_elements.h:605
void output(std::ostream &outfile)
Definition: foeppl_von_karman_elements.h:522
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
Definition: foeppl_von_karman_elements.h:632
static const unsigned Initial_Nvalue
Set the data for the number of Variables at each node - 8.
Definition: foeppl_von_karman_elements.h:496
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: foeppl_von_karman_elements.h:530
void output(FILE *file_pt)
Definition: foeppl_von_karman_elements.h:538
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: foeppl_von_karman_elements.h:554
unsigned required_nvalue(const unsigned &n) const
Definition: foeppl_von_karman_elements.h:515
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: foeppl_von_karman_elements.h:565
QFoepplvonKarmanElement()
Definition: foeppl_von_karman_elements.h:501
void operator=(const QFoepplvonKarmanElement< NNODE_1D > &)=delete
Broken assignment operator.
unsigned ntstorage() const
Definition: timesteppers.h:601
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237
RealScalar s
Definition: level1_cplx_impl.h:130
squared absolute sa ArrayBase::abs2 DOXCOMMA MatrixBase::cwiseAbs2 sa Eigen::abs2 DOXCOMMA Eigen::pow DOXCOMMA ArrayBase::square nearest sa Eigen::floor DOXCOMMA Eigen::ceil DOXCOMMA ArrayBase::round nearest integer not less than the given sa Eigen::floor DOXCOMMA ArrayBase::ceil not a number test
Definition: GlobalFunctions.h:109
int error
Definition: calibrate.py:297
@ W
Definition: quadtree.h:63
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
list x
Definition: plotDoE.py:28
t
Definition: plotPSD.py:36
Definition: indexed_view.cpp:20
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2