26#ifndef OOMPH_AXISYM_POROELASTICITY_ELEMENTS_HEADER
27#define OOMPH_AXISYM_POROELASTICITY_ELEMENTS_HEADER
31#include <oomph-lib-config.h>
96 const double&
nu()
const
101 std::ostringstream error_message;
102 error_message <<
"No pointer to Poisson's ratio set. Please set one!\n";
241 for (
unsigned i = 0;
i <
n;
i++)
248 (*Solid_body_force_fct_pt)(time, x, b);
263 for (
unsigned i = 0;
i <
n;
i++)
270 (*Fluid_body_force_fct_pt)(time, x, b);
287 (*Mass_source_fct_pt)(time, x, b);
321 virtual double q_edge(
const unsigned&
j)
const = 0;
325 virtual double q_edge(
const unsigned&
t,
const unsigned&
j)
const = 0;
329 const unsigned&
j)
const = 0;
344 virtual double q_internal(
const unsigned&
t,
const unsigned&
j)
const = 0;
347 virtual void set_q_edge(
const unsigned&
j,
const double& value) = 0;
356 const double& value) = 0;
362 const double& value) = 0;
403 const unsigned&
edge,
const unsigned&
j)
const = 0;
412 const double& value) = 0;
421 virtual double p_value(
const unsigned&
j)
const = 0;
496 for (
unsigned i = 0;
i < 2;
i++)
544 for (
unsigned i = 0;
i < 2;
i++)
587 for (
unsigned i = 0;
i < 2;
i++)
635 const unsigned&
i)
const
675 for (
unsigned i = 0;
i < 2;
i++)
695 for (
unsigned i = 0;
i < 2;
i++)
719 for (
unsigned i = 0;
i < 2;
i++)
759 const unsigned i)
const
882 double du_dt(
const unsigned&
n,
const unsigned&
i)
const
912 double d2u_dt2(
const unsigned&
n,
const unsigned&
i)
const
1022 for (
unsigned j = 0;
j <
np;
j++)
1030 for (
unsigned j = 0;
j <
nq;
j++)
1037 for (
unsigned j = 0;
j <
nq;
j++)
1061 const unsigned&
nplot)
const
1109 <<
"Axisymmetric poroelasticity elements only store 6 fields "
1126 return "radial displacement";
1130 return "axial displacement";
1134 return "radial flux";
1138 return "axial flux";
1142 return "divergence of flux";
1146 return "pore pressure";
1150 return "radial skeleton velocity";
1154 return "axial skeleton velocity";
1161 <<
"AxisymmetricPoroelasticityEquations only store 8 fields "
1177 for (
unsigned i = 0;
i < 2;
i++)
1183 for (
unsigned i = 0;
i < 2;
i++)
1189 for (
unsigned i = 0;
i < 2;
i++)
1203 data.push_back(
du_dt[0]);
1204 data.push_back(
du_dt[1]);
1234 const unsigned&
nplot,
1240 const unsigned&
nplot,
1246 const unsigned&
nplot,
1300 const unsigned&
ipt,
1395 template<
class AXISYMMETRIC_POROELASTICITY_ELEMENT>
1415 <<
"AxisymmetricPoroelasticity elements only store 4 fields so fld = "
1416 <<
fld <<
" is illegal \n";
1431 unsigned nvalue =
dat_pt->nvalue();
1432 for (
unsigned i = 0;
i < nvalue;
i++)
1443 for (
unsigned j = 0;
j <
n;
j++)
1445 unsigned nvalue = this->nedge_flux_interpolation_point();
1446 for (
unsigned i = 0;
i < nvalue;
i++)
1452 if (this->nq_basis_internal() > 0)
1456 for (
unsigned i = 0;
i < nvalue;
i++)
1467 for (
unsigned j = 0;
j <
nnod;
j++)
1472 this->
node_pt(
j), this->u_index_axisym_poroelasticity(
fld)));
1496 <<
"AxisymmetricPoroelasticity elements only store 4 fields so fld = "
1497 <<
fld <<
" is illegal\n";
1514 return return_value;
1535 <<
"AxisymmetricPoroelasticity elements only store 4 fields so fld = "
1536 <<
fld <<
" is illegal.\n";
1545 const unsigned n_dim = this->
dim();
1547 const unsigned n_q_basis = this->nq_basis();
1548 const unsigned n_p_basis = this->np_basis();
1560 double J = this->shape_basis_test_local(
s,
1577 for (
unsigned i = 0;
i <
n;
i++)
1587 for (
unsigned i = 0;
i <
n;
i++)
1589 for (
unsigned j = 0;
j <
m;
j++)
1611 const unsigned&
fld,
1619 <<
"AxisymmetricPoroelasticity elements only store 4 fields so fld = "
1620 <<
fld <<
" is illegal\n";
1636 "Pressure in AxisymmetricPoroelasticity has no time-dep.!",
1647 "Don't call this function for AxisymmetricPoroelasticity!",
1668 <<
"AxisymmetricPoroelasticity elements only store 4 fields so fld = "
1669 <<
fld <<
" is illegal\n";
1692 return return_value;
1704 <<
"AxisymmetricPoroelasticity elements only store 4 fields so fld = "
1705 <<
fld <<
" is illegal\n";
1721 unsigned nedge = this->nq_basis_edge();
1745 AXISYMMETRIC_POROELASTICITY_ELEMENT::output(
outfile,
nplot);
1754 const unsigned&
flag)
1783 for (
unsigned i = 0;
i <
n_dim;
i++)
1810 throw OomphLibError(
"Trying to project Lagrangian coordinates in "
1811 "non-SolidElement\n",
1935 "positions not in SolidElement\n",
2063 for (
unsigned i = 0;
i <
n_dim;
i++)
2089 "Wrong field specified. This should never happen\n",
2098 throw OomphLibError(
"Wrong type specificied in Projection_type. "
2099 "This should never happen\n",
2115 template<
class ELEMENT>
Class implementing the generic maths of the axisym poroelasticity equations: axisym linear elasticity...
static double Default_lambda_sq_value
Static default value for timescale ratio.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Z2 flux (use Darcy flux)
void mass_source(const double &time, const Vector< double > &x, double &b) const
Indirect access to the mass source function - returns 0 if no mass source function has been set.
static double Default_youngs_modulus_value
Static default value for Young's modulus (1.0 – for natural scaling, i.e. all stresses have been non-...
virtual unsigned q_internal_index() const =0
Return the index of the internal data where the q internal degrees of freedom are stored.
const double & porosity() const
Access function for the porosity.
double transform_basis(const Vector< double > &s, const Shape &q_basis_local, Shape &psi, Shape &q_basis) const
Performs a div-conserving transformation of the vector basis functions from the reference element to ...
void scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specif...
double * Lambda_sq_pt
Timescale ratio (non-dim. density)
virtual unsigned nedge_flux_interpolation_point() const =0
Returns the number of flux_interpolation points along each edge of the element.
double dq_internal_dt(const unsigned &n) const
dq_internal/dt for the n-th internal degree of freedom
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
Compute the error between the FE solution and the exact solution using the H(div) norm for q and L^2 ...
void point_output_data(const Vector< double > &s, Vector< double > &data)
Output solution in data vector at local cordinates s: r,z,u_r,u_z,q_r,q_z,div_q,p,...
unsigned num_Z2_flux_terms()
Number off flux terms for Z2 error estimator (use Darcy flux)
void output_with_projected_flux(std::ostream &outfile, const unsigned &nplot, const Vector< double > unit_normal)
Output incl. projection of fluxes into direction of the specified unit vector.
virtual void scale_basis(Shape &basis) const =0
Scale the edge basis to allow arbitrary edge mappings.
virtual unsigned u_index_axisym_poroelasticity(const unsigned &j) const =0
Return the nodal index of the j-th solid displacement unknown.
unsigned self_test()
Self test.
double *& porosity_pt()
Access function for pointer to the porosity.
unsigned nscalar_paraview() const
Number of scalars/fields output by this element. Reimplements broken virtual function in base class.
const double & nu() const
Access function for Poisson's ratio.
virtual int q_edge_local_eqn(const unsigned &j) const =0
Return the equation number of the j-th edge (flux) degree of freedom.
MassSourceFctPt & mass_source_fct_pt()
Access function: Pointer to mass source function.
void interpolated_u(const Vector< double > &s, Vector< double > &disp) const
Calculate the FE representation of u.
virtual unsigned q_edge_node_number(const unsigned &j) const =0
Return the number of the node where the jth edge unknown is stored.
virtual void set_q_edge(const unsigned &t, const unsigned &j, const double &value)=0
Set the values of the j-th edge (flux) degree of freedom at time history level t.
virtual int q_internal_local_eqn(const unsigned &j) const =0
Return the equation number of the j-th internal degree of freedom.
double interpolated_div_q(const Vector< double > &s) const
Calculate the FE representation of div q and return it.
const double & alpha() const
Access function for alpha, the Biot parameter.
virtual Data * p_data_pt() const =0
Return pointer to the Data object that stores the pressure values.
double *& lambda_sq_pt()
Access function for pointer to timescale ratio (nondim density)
double interpolated_u(const Vector< double > &s, const unsigned &i) const
Calculate the FE representation of the i-th component of u.
virtual double shape_basis_test_local_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsi, Shape &u_basis, Shape &u_test, DShape &du_basis_dx, DShape &du_test_dx, Shape &q_basis, Shape &q_test, Shape &p_basis, Shape &p_test, Shape &div_q_basis_ds, Shape &div_q_test_ds) const =0
Compute the geometric basis, and the q, p and divergence basis functions and test functions at integr...
virtual void set_q_edge(const unsigned &j, const double &value)=0
Set the values of the j-th edge (flux) degree of freedom.
void interpolated_q(const Vector< double > &s, Vector< double > &q) const
Calculate the FE representation of q.
const double & permeability() const
Access function for the nondim permeability.
double interpolated_div_u(const Vector< double > &s, Vector< double > &div_u_components) const
Calculate the FE representation of the divergence of the skeleton displ, div(u), and its components: ...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in contribution to residuals for the Darcy equations.
AxisymmetricPoroelasticityEquations()
Constructor.
static double Default_porosity_value
Static default value for the porosity.
void(* MassSourceFctPt)(const double &time, const Vector< double > &x, double &f)
Mass source function pointer typedef.
double * Alpha_pt
Alpha – the biot parameter.
double * Nu_pt
Pointer to Poisson's ratio.
virtual void set_q_internal(const unsigned &j, const double &value)=0
Set the values of the j-th internal degree of freedom.
virtual double shape_basis_test_local(const Vector< double > &s, Shape &psi, DShape &dpsi, Shape &u_basis, Shape &u_test, DShape &du_basis_dx, DShape &du_test_dx, Shape &q_basis, Shape &q_test, Shape &p_basis, Shape &p_test, Shape &div_q_basis_ds, Shape &div_q_test_ds) const =0
Compute the geometric basis, and the q, p and divergence basis functions and test functions at local ...
double du_dt(const unsigned &n, const unsigned &i) const
du/dt at local node n
double *& nu_pt()
Access function for pointer to Poisson's ratio.
virtual Vector< Data * > q_edge_data_pt() const =0
Return vector of pointers to the Data objects that store the edge flux values.
double * Porosity_pt
Porosity.
virtual unsigned q_edge_index(const unsigned &j) const =0
Return the nodal index at which the jth edge unknown is stored.
SourceFctPt Fluid_body_force_fct_pt
Pointer to fluid source function.
virtual void face_local_coordinate_of_flux_interpolation_point(const unsigned &edge, const unsigned &n, Vector< double > &s) const =0
Compute the face element coordinates of the nth flux interpolation point along an edge.
double *& alpha_pt()
Access function for pointer to alpha, the Biot parameter.
MassSourceFctPt mass_source_fct_pt() const
Access function: Pointer to mass source function (const version)
virtual unsigned face_index_of_q_edge_basis_fct(const unsigned &j) const =0
Return the face index associated with j-th edge flux degree of freedom.
double interpolated_q(const unsigned &t, const Vector< double > &s, const unsigned i) const
Calculate the FE representation of the i-th component of q at time level t (t=0: current)
static double Default_permeability_value
Static default value for the permeability (1.0 for natural scaling; i.e. timescale is given by L^2/(k...
void get_q_basis(const Vector< double > &s, Shape &q_basis) const
Compute the transformed basis at local coordinate s.
virtual unsigned face_index_of_edge(const unsigned &j) const =0
Return the face index associated with specified edge.
SourceFctPt & fluid_body_force_fct_pt()
Access function: Pointer to fluid force function.
static double Default_permeability_ratio_value
Static default value for the ratio of the material's actual permeability to the permeability used to ...
void(* SourceFctPt)(const double &time, const Vector< double > &x, Vector< double > &f)
Source function pointer typedef.
virtual void pin_q_internal_value(const unsigned &j, const double &value)=0
Pin the jth internal q value and set it to specified value.
void interpolated_q(const unsigned &t, const Vector< double > &s, Vector< double > &q) const
Calculate the FE representation of q at time level t (t=0: current)
MassSourceFctPt Mass_source_fct_pt
Pointer to the mass source function.
void interpolated_div_q(const Vector< double > &s, double &div_q) const
Calculate the FE representation of div u.
double transform_basis(const Vector< double > &s, const Shape &q_basis_local, Shape &psi, DShape &dpsi, Shape &q_basis) const
Performs a div-conserving transformation of the vector basis functions from the reference element to ...
virtual unsigned np_basis() const =0
Return the total number of pressure basis functions.
virtual unsigned nq_basis() const
Return the total number of computational basis functions for q.
virtual double q_edge(const unsigned &t, const unsigned &j) const =0
Return the values of the j-th edge (flux) degree of freedom at time history level t.
virtual unsigned required_nvalue(const unsigned &n) const =0
Number of values required at node n.
virtual void edge_flux_interpolation_point_global(const unsigned &edge, const unsigned &j, Vector< double > &x) const =0
Compute the global coordinates of the jth flux_interpolation point along an edge.
double *& permeability_pt()
Access function for pointer to the nondim permeability.
double * Permeability_pt
permeability
virtual void set_p_value(const unsigned &j, const double &value)=0
Set the jth pressure value.
double * Youngs_modulus_pt
Pointer to the nondim Young's modulus.
void switch_off_darcy()
Switch off Darcy flow.
virtual void pin_p_value(const unsigned &j, const double &p)=0
Pin the jth pressure value and set it to p.
virtual unsigned nq_basis_edge() const =0
Return the number of edge basis functions for q.
double *& permeability_ratio_pt()
Access function for pointer to ratio of the material's actual permeability to the permeability used i...
virtual Data * q_internal_data_pt() const =0
Return pointer to the Data object that stores the internal flux values.
virtual double q_internal(const unsigned &j) const =0
Return the values of the j-th internal degree of freedom.
void interpolated_p(const Vector< double > &s, double &p) const
Calculate the FE representation of p.
const double & youngs_modulus() const
Access function to non-dim Young's modulus (ratio of actual Young's modulus to reference stress used ...
void set_q_internal_timestepper(TimeStepper *const time_stepper_pt)
Set the timestepper of the q internal data object.
void interpolated_du_dt(const Vector< double > &s, Vector< double > &du_dt) const
Calculate the FE representation of du_dt.
bool darcy_is_switched_off()
Is Darcy flow switched off?
virtual double p_value(const unsigned &j) const =0
Return the jth pressure value.
std::string scalar_name_paraview(const unsigned &i) const
Name of the i-th scalar field. Default implementation returns V1 for the first one,...
virtual Vector< double > edge_flux_interpolation_point(const unsigned &edge, const unsigned &j) const =0
Returns the local coordinate of the jth flux_interpolation point along the specified edge.
double *& youngs_modulus_pt()
Pointer to non-dim Young's modulus (ratio of actual Young's modulus to reference stress used to non-d...
double interpolated_div_du_dt(const Vector< double > &s, Vector< double > &div_dudt_components) const
Calculate the FE representation of the divergence of the skeleton velocity, div(du/dt),...
SourceFctPt solid_body_force_fct_pt() const
Access function: Pointer to solid body force function (const version)
virtual void fill_in_generic_residual_contribution(Vector< double > &residuals, DenseMatrix< double > &jacobian, bool flag)
Fill in residuals and, if flag==true, jacobian.
double interpolated_p(const Vector< double > &s) const
Calculate the FE representation of p and return it.
virtual void get_div_q_basis_local(const Vector< double > &s, Shape &div_q_basis_ds) const =0
Compute the local form of the q basis and dbasis/ds at local coordinate s.
SourceFctPt fluid_body_force_fct_pt() const
Access function: Pointer to fluid force function (const version)
SourceFctPt & solid_body_force_fct_pt()
Access function: Pointer to solid body force function.
double interpolated_q(const Vector< double > &s, const unsigned i) const
Calculate the FE representation of the i-th component of q.
double *& density_ratio_pt()
Access function for pointer to the density ratio (fluid to solid)
double d2u_dt2(const unsigned &n, const unsigned &i) const
d^2u/dt^2 at local node n
static double Default_density_ratio_value
Static default value for the density ratio.
const double & lambda_sq() const
Access function for timescale ratio (nondim density)
virtual double q_internal(const unsigned &t, const unsigned &j) const =0
Return the values of the j-th internal degree of freedom at time history level t.
bool Darcy_is_switched_off
Boolean to record that darcy has been switched off.
virtual void pin_q_edge_value(const unsigned &j, const double &value)=0
Pin the j-th edge (flux) degree of freedom and set it to specified value.
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output FE representation of exact soln: x,y,u1,u2,div_q,p at Nplot^2 plot points.
SourceFctPt Solid_body_force_fct_pt
Pointer to solid body force function.
void fluid_body_force(const double &time, const Vector< double > &x, Vector< double > &b) const
Indirect access to the fluid body force function - returns 0 if no forcing function has been set.
static double Default_alpha_value
Static default value for alpha, the biot parameter.
double * Permeability_ratio_pt
Ratio of the material's actual permeability to the permeability used in the non-dimensionalisation of...
double interpolated_u(const unsigned &t, const Vector< double > &s, const unsigned &i) const
Calculate the FE representation of the i-th component of u at time level t (t=0: current)
void solid_body_force(const double &time, const Vector< double > &x, Vector< double > &b) const
Indirect access to the solid body force function - returns 0 if no forcing function has been set.
virtual double q_edge(const unsigned &j) const =0
Return the values of the j-th edge (flux) degree of freedom.
virtual void get_p_basis(const Vector< double > &s, Shape &p_basis) const =0
Compute the pressure basis.
const double & density_ratio() const
Access function for the density ratio (fluid to solid)
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in the Jacobian matrix for the Newton method.
double dq_edge_dt(const unsigned &n) const
dq_edge/dt for the n-th edge degree of freedom
double * Density_ratio_pt
Density ratio.
const double & permeability_ratio() const
Access function for the ratio of the material's actual permeability to the permeability used in the n...
virtual void set_q_internal(const unsigned &t, const unsigned &j, const double &value)=0
Set the values of the j-th internal degree of freedom at time history level t.
virtual unsigned nq_basis_internal() const =0
Return the number of internal basis functions for q.
void output(std::ostream &outfile)
Output with default number of plot points.
virtual int p_local_eqn(const unsigned &j) const =0
Return the equation number of the j-th pressure degree of freedom.
virtual void get_q_basis_local(const Vector< double > &s, Shape &q_basis) const =0
Comute the local form of the q basis at local coordinate s.
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 ...
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
unsigned ntstorage() const
Return total number of doubles stored per value to record time history of each value (one for steady ...
void set_time_stepper(TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
Set a new timestepper by resizing the appropriate storage. If already assigned the equation numbering...
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...
FiniteElement *& external_element_pt(const unsigned &interaction_index, const unsigned &ipt)
Access function to source element for specified interaction index at specified integration point.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
A general Finite Element class.
virtual unsigned nplot_points_paraview(const unsigned &nplot) const
Return the number of actual plot points for paraview plot with parameter nplot. Broken virtual; can b...
virtual double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
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...
unsigned nnodal_position_type() const
Return the number of coordinate types that the element requires to interpolate the geometry between t...
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
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 .
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 double local_to_eulerian_mapping(const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
Calculate the mapping from local to Eulerian coordinates, given the derivatives of the shape function...
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...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
virtual void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Function to compute the geometric shape functions and derivatives w.r.t. local 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 .
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
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.
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....
Axisymmetric poro elasticity upgraded to become projectable.
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (Note: count includes current value!)
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
void residual_for_projection(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Residual for the projection step. Flag indicates if we want the Jacobian (1) or not (0)....
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 ...
ProjectableAxisymmetricPoroelasticityElement()
Constructor [this was only required explicitly from gcc 4.5.2 onwards...].
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.
unsigned nfields_for_projection()
Number of fields to be projected: 4 (two displacement components, pressure, Darcy flux)
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
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...
void output(std::ostream &outfile, const unsigned &nplot)
Output FE representation of soln as in underlying element.
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. (Note: count includes current value!...
Template-free Base class for projectable elements.
unsigned Projected_lagrangian
When projecting the Lagrangain coordinates indicate which coordiante is to be projected.
Projection_Type Projection_type
Variable to indicate if we're projecting the history values of the nodal coordinates (Coordinate) the...
unsigned Time_level_for_projection
Time level we are projecting (0=current values; >0: history values)
unsigned Projected_field
Field that is currently being projected.
unsigned Projected_coordinate
When projecting the history values of the nodal coordinates, this is the coordinate we're projecting.
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...
SolidFiniteElement class.
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
Base class for time-stepping schemes. Timestepper provides an approximation of the temporal derivativ...
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...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).