![]() |
|
#include <solid_elements.h>
Public Types | |
typedef void(* | IsotropicGrowthFctPt) (const Vector< double > &xi, double &gamma) |
typedef double(* | PrestressFctPt) (const unsigned &i, const unsigned &j, const Vector< double > &xi) |
typedef void(* | BodyForceFctPt) (const double &t, const Vector< double > &xi, Vector< double > &b) |
![]() | |
typedef double(* | MultiplierFctPt) (const Vector< double > &xi) |
![]() | |
typedef void(* | SteadyExactSolutionFctPt) (const Vector< double > &, Vector< double > &) |
typedef void(* | UnsteadyExactSolutionFctPt) (const double &, const Vector< double > &, Vector< double > &) |
Public Member Functions | |
PVDEquationsBase () | |
ConstitutiveLaw *& | constitutive_law_pt () |
Return the constitutive law pointer. More... | |
const double & | lambda_sq () const |
Access function for timescale ratio (nondim density) More... | |
double *& | lambda_sq_pt () |
Access function for pointer to timescale ratio (nondim density) More... | |
IsotropicGrowthFctPt & | isotropic_growth_fct_pt () |
Access function: Pointer to isotropic growth function. More... | |
PrestressFctPt & | prestress_fct_pt () |
Access function: Pointer to pre-stress function. More... | |
IsotropicGrowthFctPt | isotropic_growth_fct_pt () const |
Access function: Pointer to isotropic growth function (const version) More... | |
BodyForceFctPt & | body_force_fct_pt () |
Access function: Pointer to body force function. More... | |
BodyForceFctPt | body_force_fct_pt () const |
Access function: Pointer to body force function (const version) More... | |
void | enable_inertia () |
Switch on solid inertia. More... | |
void | disable_inertia () |
Switch off solid inertia. More... | |
bool | is_inertia_enabled () const |
Access function to flag that switches inertia on/off (const version) More... | |
virtual unsigned | npres_solid () const |
virtual int | solid_p_local_eqn (const unsigned &i) const |
virtual int | solid_p_nodal_index () const |
virtual void | unpin_elemental_solid_pressure_dofs ()=0 |
Unpin all solid pressure dofs in the element. More... | |
virtual void | pin_elemental_redundant_nodal_solid_pressures () |
Pin the element's redundant solid pressures (needed for refinement) More... | |
virtual void | get_stress (const Vector< double > &s, DenseMatrix< double > &sigma)=0 |
void | get_strain (const Vector< double > &s, DenseMatrix< double > &strain) const |
Return the strain tensor. More... | |
void | get_energy (double &pot_en, double &kin_en) |
Get potential (strain) and kinetic energy. More... | |
void | get_deformed_covariant_basis_vectors (const Vector< double > &s, DenseMatrix< double > &def_covariant_basis) |
void | get_principal_stress (const Vector< double > &s, DenseMatrix< double > &principal_stress_vector, Vector< double > &principal_stress) |
virtual void | get_isotropic_growth (const unsigned &ipt, const Vector< double > &s, const Vector< double > &xi, double &gamma) const |
void | body_force (const Vector< double > &xi, Vector< double > &b) const |
unsigned | ndof_types () const |
returns the number of DOF types associated with this element. More... | |
void | get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const |
void | enable_evaluate_jacobian_by_fd () |
Set Jacobian to be evaluated by FD? Else: Analytically. More... | |
void | disable_evaluate_jacobian_by_fd () |
Set Jacobian to be evaluated analytically Else: by FD. More... | |
bool | is_jacobian_evaluated_by_fd () const |
Return the flag indicating whether the jacobian is evaluated by fd. More... | |
double | prestress (const unsigned &i, const unsigned &j, const Vector< double > xi) |
![]() | |
void | set_lagrangian_dimension (const unsigned &lagrangian_dimension) |
virtual bool | has_internal_solid_data () |
SolidFiniteElement () | |
Constructor: Set defaults. More... | |
virtual | ~SolidFiniteElement () |
Destructor to clean up any allocated memory. More... | |
SolidFiniteElement (const SolidFiniteElement &)=delete | |
Broken copy constructor. More... | |
unsigned | ngeom_data () const |
Broken assignment operator. More... | |
Data * | geom_data_pt (const unsigned &j) |
void | identify_geometric_data (std::set< Data * > &geometric_data_pt) |
double | zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const |
virtual void | get_x_and_xi (const Vector< double > &s, Vector< double > &x_fe, Vector< double > &x, Vector< double > &xi_fe, Vector< double > &xi) const |
virtual void | set_macro_elem_pt (MacroElement *macro_elem_pt) |
virtual void | set_macro_elem_pt (MacroElement *macro_elem_pt, MacroElement *undeformed_macro_elem_pt) |
void | set_undeformed_macro_elem_pt (MacroElement *undeformed_macro_elem_pt) |
MacroElement * | undeformed_macro_elem_pt () |
Access function to pointer to "undeformed" macro element. More... | |
double | dshape_lagrangian (const Vector< double > &s, Shape &psi, DShape &dpsidxi) const |
virtual double | dshape_lagrangian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidxi) const |
double | d2shape_lagrangian (const Vector< double > &s, Shape &psi, DShape &dpsidxi, DShape &d2psidxi) const |
virtual double | d2shape_lagrangian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidxi, DShape &d2psidxi) const |
unsigned | lagrangian_dimension () const |
unsigned | nnodal_lagrangian_type () const |
Node * | construct_node (const unsigned &n) |
Construct the local node n and return a pointer to it. More... | |
Node * | construct_node (const unsigned &n, TimeStepper *const &time_stepper_pt) |
Node * | construct_boundary_node (const unsigned &n) |
Node * | construct_boundary_node (const unsigned &n, TimeStepper *const &time_stepper_pt) |
virtual void | assign_all_generic_local_eqn_numbers (const bool &store_local_dof_pt) |
void | describe_local_dofs (std::ostream &out, const std::string ¤t_string) const |
double | raw_lagrangian_position (const unsigned &n, const unsigned &i) const |
double | raw_lagrangian_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const |
double | lagrangian_position (const unsigned &n, const unsigned &i) const |
Return i-th Lagrangian coordinate at local node n. More... | |
double | lagrangian_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const |
virtual double | interpolated_xi (const Vector< double > &s, const unsigned &i) const |
virtual void | interpolated_xi (const Vector< double > &s, Vector< double > &xi) const |
virtual void | interpolated_dxids (const Vector< double > &s, DenseMatrix< double > &dxids) const |
virtual void | J_lagrangian (const Vector< double > &s) const |
virtual double | J_lagrangian_at_knot (const unsigned &ipt) const |
SolidInitialCondition *& | solid_ic_pt () |
Pointer to object that describes the initial condition. More... | |
void | enable_solve_for_consistent_newmark_accel () |
void | disable_solve_for_consistent_newmark_accel () |
Set to reset the problem being solved to be the standard problem. More... | |
MultiplierFctPt & | multiplier_fct_pt () |
MultiplierFctPt | multiplier_fct_pt () const |
virtual void | get_residuals_for_solid_ic (Vector< double > &residuals) |
void | fill_in_residuals_for_solid_ic (Vector< double > &residuals) |
void | fill_in_jacobian_for_solid_ic (Vector< double > &residuals, DenseMatrix< double > &jacobian) |
void | fill_in_jacobian_for_newmark_accel (DenseMatrix< double > &jacobian) |
void | compute_norm (double &el_norm) |
int | position_local_eqn (const unsigned &n, const unsigned &k, const unsigned &j) const |
![]() | |
void | set_dimension (const unsigned &dim) |
void | set_nodal_dimension (const unsigned &nodal_dim) |
void | set_nnodal_position_type (const unsigned &nposition_type) |
Set the number of types required to interpolate the coordinate. More... | |
void | set_n_node (const unsigned &n) |
int | nodal_local_eqn (const unsigned &n, const unsigned &i) const |
double | dJ_eulerian_at_knot (const unsigned &ipt, Shape &psi, DenseMatrix< double > &djacobian_dX) const |
FiniteElement () | |
Constructor. More... | |
virtual | ~FiniteElement () |
FiniteElement (const FiniteElement &)=delete | |
Broken copy constructor. More... | |
virtual bool | local_coord_is_valid (const Vector< double > &s) |
Broken assignment operator. More... | |
virtual void | move_local_coord_back_into_element (Vector< double > &s) const |
void | get_centre_of_gravity_and_max_radius_in_terms_of_zeta (Vector< double > &cog, double &max_radius) const |
virtual void | local_coordinate_of_node (const unsigned &j, Vector< double > &s) const |
virtual void | local_fraction_of_node (const unsigned &j, Vector< double > &s_fraction) |
virtual double | local_one_d_fraction_of_node (const unsigned &n1d, const unsigned &i) |
MacroElement * | macro_elem_pt () |
Access function to pointer to macro element. More... | |
void | get_x (const Vector< double > &s, Vector< double > &x) const |
void | get_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) |
virtual void | get_x_from_macro_element (const Vector< double > &s, Vector< double > &x) const |
virtual void | get_x_from_macro_element (const unsigned &t, const Vector< double > &s, Vector< double > &x) |
virtual void | set_integration_scheme (Integral *const &integral_pt) |
Set the spatial integration scheme. More... | |
Integral *const & | integral_pt () const |
Return the pointer to the integration scheme (const version) More... | |
virtual void | shape (const Vector< double > &s, Shape &psi) const =0 |
virtual void | shape_at_knot (const unsigned &ipt, Shape &psi) const |
virtual void | dshape_local (const Vector< double > &s, Shape &psi, DShape &dpsids) const |
virtual void | dshape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids) const |
virtual void | d2shape_local (const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const |
virtual void | d2shape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids, DShape &d2psids) const |
virtual double | J_eulerian (const Vector< double > &s) const |
virtual double | J_eulerian_at_knot (const unsigned &ipt) const |
void | check_J_eulerian_at_knots (bool &passed) const |
void | check_jacobian (const double &jacobian) const |
double | dshape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx) const |
virtual double | dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx) const |
virtual double | dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsi, DenseMatrix< double > &djacobian_dX, RankFourTensor< double > &d_dpsidx_dX) const |
double | d2shape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx, DShape &d2psidx) const |
virtual double | d2shape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx, DShape &d2psidx) const |
virtual void | assign_nodal_local_eqn_numbers (const bool &store_local_dof_pt) |
virtual void | describe_nodal_local_dofs (std::ostream &out, const std::string ¤t_string) const |
Node *& | node_pt (const unsigned &n) |
Return a pointer to the local node n. More... | |
Node *const & | node_pt (const unsigned &n) const |
Return a pointer to the local node n (const version) More... | |
unsigned | nnode () const |
Return the number of nodes. More... | |
virtual unsigned | nnode_1d () const |
double | raw_nodal_position (const unsigned &n, const unsigned &i) const |
double | raw_nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const |
double | raw_dnodal_position_dt (const unsigned &n, const unsigned &i) const |
double | raw_dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const |
double | raw_nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const |
double | raw_nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const |
double | raw_dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const |
double | raw_dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const |
double | nodal_position (const unsigned &n, const unsigned &i) const |
double | nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const |
double | dnodal_position_dt (const unsigned &n, const unsigned &i) const |
Return the i-th component of nodal velocity: dx/dt at local node n. More... | |
double | dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const |
double | nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const |
double | nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const |
double | dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const |
double | dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const |
virtual void | get_dresidual_dnodal_coordinates (RankThreeTensor< double > &dresidual_dnodal_coordinates) |
virtual void | disable_ALE () |
virtual void | enable_ALE () |
virtual unsigned | required_nvalue (const unsigned &n) const |
unsigned | nnodal_position_type () const |
bool | has_hanging_nodes () const |
unsigned | nodal_dimension () const |
Return the required Eulerian dimension of the nodes in this element. More... | |
virtual unsigned | nvertex_node () const |
virtual Node * | vertex_node_pt (const unsigned &j) const |
int | get_node_number (Node *const &node_pt) const |
virtual Node * | get_node_at_local_coordinate (const Vector< double > &s) const |
double | raw_nodal_value (const unsigned &n, const unsigned &i) const |
double | raw_nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const |
double | nodal_value (const unsigned &n, const unsigned &i) const |
double | nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const |
unsigned | dim () const |
virtual ElementGeometry::ElementGeometry | element_geometry () const |
Return the geometry type of the element (either Q or T usually). More... | |
virtual double | interpolated_x (const Vector< double > &s, const unsigned &i) const |
Return FE interpolated coordinate x[i] at local coordinate s. More... | |
virtual double | interpolated_x (const unsigned &t, const Vector< double > &s, const unsigned &i) const |
virtual void | interpolated_x (const Vector< double > &s, Vector< double > &x) const |
Return FE interpolated position x[] at local coordinate s as Vector. More... | |
virtual void | interpolated_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) const |
virtual double | interpolated_dxdt (const Vector< double > &s, const unsigned &i, const unsigned &t) |
virtual void | interpolated_dxdt (const Vector< double > &s, const unsigned &t, Vector< double > &dxdt) |
void | position (const Vector< double > &zeta, Vector< double > &r) const |
void | position (const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const |
void | dposition_dt (const Vector< double > &zeta, const unsigned &t, Vector< double > &drdt) |
void | interpolated_zeta (const Vector< double > &s, Vector< double > &zeta) const |
void | locate_zeta (const Vector< double > &zeta, GeomObject *&geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false) |
virtual void | node_update () |
virtual void | identify_field_data_for_interactions (std::set< std::pair< Data *, unsigned >> &paired_field_data) |
virtual double | s_min () const |
Min value of local coordinate. More... | |
virtual double | s_max () const |
Max. value of local coordinate. More... | |
double | size () const |
virtual double | compute_physical_size () const |
virtual void | point_output_data (const Vector< double > &s, Vector< double > &data) |
void | point_output (std::ostream &outfile, const Vector< double > &s) |
virtual unsigned | nplot_points_paraview (const unsigned &nplot) const |
virtual unsigned | nsub_elements_paraview (const unsigned &nplot) const |
void | output_paraview (std::ofstream &file_out, const unsigned &nplot) const |
virtual void | write_paraview_output_offset_information (std::ofstream &file_out, const unsigned &nplot, unsigned &counter) const |
virtual void | write_paraview_type (std::ofstream &file_out, const unsigned &nplot) const |
virtual void | write_paraview_offsets (std::ofstream &file_out, const unsigned &nplot, unsigned &offset_sum) const |
virtual unsigned | nscalar_paraview () const |
virtual void | scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const |
virtual void | scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const |
virtual void | scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const |
virtual std::string | scalar_name_paraview (const unsigned &i) const |
virtual void | output (std::ostream &outfile) |
virtual void | output (std::ostream &outfile, const unsigned &n_plot) |
virtual void | output (const unsigned &t, std::ostream &outfile, const unsigned &n_plot) const |
virtual void | output (FILE *file_pt) |
virtual void | output (FILE *file_pt, const unsigned &n_plot) |
virtual void | output_fct (std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) |
Output an exact solution over the element. More... | |
virtual void | output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) |
Output a time-dependent exact solution over the element. More... | |
virtual void | output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, const SolutionFunctorBase &exact_soln) const |
Output a time-dependent exact solution over the element. More... | |
virtual void | get_s_plot (const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const |
virtual std::string | tecplot_zone_string (const unsigned &nplot) const |
virtual void | write_tecplot_zone_footer (std::ostream &outfile, const unsigned &nplot) const |
virtual void | write_tecplot_zone_footer (FILE *file_pt, const unsigned &nplot) const |
virtual unsigned | nplot_points (const unsigned &nplot) const |
virtual void | compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm) |
Calculate the norm of the error and that of the exact solution. More... | |
virtual void | compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm) |
Calculate the norm of the error and that of the exact solution. More... | |
virtual void | compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm) |
virtual void | compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm) |
virtual void | compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm) |
virtual void | compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm) |
virtual void | compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm) |
virtual void | compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm) |
virtual void | compute_abs_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error) |
void | integrate_fct (FiniteElement::SteadyExactSolutionFctPt integrand_fct_pt, Vector< double > &integral) |
Integrate Vector-valued function over element. More... | |
void | integrate_fct (FiniteElement::UnsteadyExactSolutionFctPt integrand_fct_pt, const double &time, Vector< double > &integral) |
Integrate Vector-valued time-dep function over element. More... | |
virtual void | build_face_element (const int &face_index, FaceElement *face_element_pt) |
virtual unsigned | self_test () |
virtual unsigned | get_bulk_node_number (const int &face_index, const unsigned &i) const |
virtual int | face_outer_unit_normal_sign (const int &face_index) const |
Get the sign of the outer unit normal on the face given by face_index. More... | |
virtual unsigned | nnode_on_face () const |
void | face_node_number_error_check (const unsigned &i) const |
Range check for face node numbers. More... | |
virtual CoordinateMappingFctPt | face_to_bulk_coordinate_fct_pt (const int &face_index) const |
virtual BulkCoordinateDerivativesFctPt | bulk_coordinate_derivatives_fct_pt (const int &face_index) const |
![]() | |
GeneralisedElement () | |
Constructor: Initialise all pointers and all values to zero. More... | |
virtual | ~GeneralisedElement () |
Virtual destructor to clean up any memory allocated by the object. More... | |
GeneralisedElement (const GeneralisedElement &)=delete | |
Broken copy constructor. More... | |
void | operator= (const GeneralisedElement &)=delete |
Broken assignment operator. More... | |
Data *& | internal_data_pt (const unsigned &i) |
Return a pointer to i-th internal data object. More... | |
Data *const & | internal_data_pt (const unsigned &i) const |
Return a pointer to i-th internal data object (const version) More... | |
Data *& | external_data_pt (const unsigned &i) |
Return a pointer to i-th external data object. More... | |
Data *const & | external_data_pt (const unsigned &i) const |
Return a pointer to i-th external data object (const version) More... | |
unsigned long | eqn_number (const unsigned &ieqn_local) const |
int | local_eqn_number (const unsigned long &ieqn_global) const |
unsigned | add_external_data (Data *const &data_pt, const bool &fd=true) |
bool | external_data_fd (const unsigned &i) const |
void | exclude_external_data_fd (const unsigned &i) |
void | include_external_data_fd (const unsigned &i) |
void | flush_external_data () |
Flush all external data. More... | |
void | flush_external_data (Data *const &data_pt) |
Flush the object addressed by data_pt from the external data array. More... | |
unsigned | ninternal_data () const |
Return the number of internal data objects. More... | |
unsigned | nexternal_data () const |
Return the number of external data objects. More... | |
unsigned | ndof () const |
Return the number of equations/dofs in the element. More... | |
void | dof_vector (const unsigned &t, Vector< double > &dof) |
Return the vector of dof values at time level t. More... | |
void | dof_pt_vector (Vector< double * > &dof_pt) |
Return the vector of pointers to dof values. More... | |
void | set_internal_data_time_stepper (const unsigned &i, TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data) |
void | assign_internal_eqn_numbers (unsigned long &global_number, Vector< double * > &Dof_pt) |
void | describe_dofs (std::ostream &out, const std::string ¤t_string) const |
void | add_internal_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt) |
virtual void | assign_local_eqn_numbers (const bool &store_local_dof_pt) |
virtual void | complete_setup_of_dependencies () |
virtual void | get_residuals (Vector< double > &residuals) |
virtual void | get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian) |
virtual void | get_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix) |
virtual void | get_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix) |
virtual void | get_dresiduals_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam) |
virtual void | get_djacobian_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam) |
virtual void | get_djacobian_and_dmass_matrix_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam) |
virtual void | get_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product) |
virtual void | get_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product) |
virtual void | get_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector) |
virtual void | compute_norm (Vector< double > &norm) |
![]() | |
GeomObject () | |
Default constructor. More... | |
GeomObject (const unsigned &ndim) | |
GeomObject (const unsigned &nlagrangian, const unsigned &ndim) | |
GeomObject (const unsigned &nlagrangian, const unsigned &ndim, TimeStepper *time_stepper_pt) | |
GeomObject (const GeomObject &dummy)=delete | |
Broken copy constructor. More... | |
void | operator= (const GeomObject &)=delete |
Broken assignment operator. More... | |
virtual | ~GeomObject () |
(Empty) destructor More... | |
unsigned | nlagrangian () const |
Access function to # of Lagrangian coordinates. More... | |
unsigned | ndim () const |
Access function to # of Eulerian coordinates. More... | |
void | set_nlagrangian_and_ndim (const unsigned &n_lagrangian, const unsigned &n_dim) |
Set # of Lagrangian and Eulerian coordinates. More... | |
TimeStepper *& | time_stepper_pt () |
TimeStepper * | time_stepper_pt () const |
virtual void | position (const double &t, const Vector< double > &zeta, Vector< double > &r) const |
virtual void | dposition (const Vector< double > &zeta, DenseMatrix< double > &drdzeta) const |
virtual void | d2position (const Vector< double > &zeta, RankThreeTensor< double > &ddrdzeta) const |
virtual void | d2position (const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const |
Static Public Member Functions | |
static void | pin_redundant_nodal_solid_pressures (const Vector< GeneralisedElement * > &element_pt) |
static void | unpin_all_solid_pressure_dofs (const Vector< GeneralisedElement * > &element_pt) |
Unpin all pressure dofs in elements listed in vector. More... | |
Protected Attributes | |
IsotropicGrowthFctPt | Isotropic_growth_fct_pt |
Pointer to isotropic growth function. More... | |
PrestressFctPt | Prestress_fct_pt |
Pointer to prestress function. More... | |
ConstitutiveLaw * | Constitutive_law_pt |
Pointer to the constitutive law. More... | |
double * | Lambda_sq_pt |
Timescale ratio (non-dim. density) More... | |
bool | Unsteady |
Flag that switches inertia on/off. More... | |
BodyForceFctPt | Body_force_fct_pt |
Pointer to body force function. More... | |
bool | Evaluate_jacobian_by_fd |
Use FD to evaluate Jacobian. More... | |
![]() | |
MacroElement * | Undeformed_macro_elem_pt |
Pointer to the element's "undeformed" macro element (NULL by default) More... | |
SolidInitialCondition * | Solid_ic_pt |
Pointer to object that specifies the initial condition. More... | |
bool | Solve_for_consistent_newmark_accel_flag |
![]() | |
MacroElement * | Macro_elem_pt |
Pointer to the element's macro element (NULL by default) More... | |
![]() | |
unsigned | NLagrangian |
Number of Lagrangian (intrinsic) coordinates. More... | |
unsigned | Ndim |
Number of Eulerian coordinates. More... | |
TimeStepper * | Geom_object_time_stepper_pt |
Static Protected Attributes | |
static double | Default_lambda_sq_value = 1.0 |
Static default value for timescale ratio (1.0 – for natural scaling) More... | |
![]() | |
static const unsigned | Default_Initial_Nvalue = 0 |
Default value for the number of values at a node. More... | |
static const double | Node_location_tolerance = 1.0e-14 |
static const unsigned | N2deriv [] = {0, 1, 3, 6} |
![]() | |
static DenseMatrix< double > | Dummy_matrix |
static std::deque< double * > | Dof_pt_deque |
Static Private Attributes | |
static int | Solid_pressure_not_stored_at_node = -100 |
Additional Inherited Members | |
![]() | |
static double | Tolerance_for_singular_jacobian = 1.0e-16 |
Tolerance below which the jacobian is considered singular. More... | |
static bool | Accept_negative_jacobian = false |
static bool | Suppress_output_while_checking_for_inverted_elements |
![]() | |
static bool | Suppress_warning_about_repeated_internal_data |
static bool | Suppress_warning_about_repeated_external_data = true |
static double | Default_fd_jacobian_step = 1.0e-8 |
![]() | |
void | fill_in_generic_jacobian_for_solid_ic (Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag) |
void | set_nnodal_lagrangian_type (const unsigned &nlagrangian_type) |
virtual double | local_to_lagrangian_mapping (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const |
double | local_to_lagrangian_mapping (const DShape &dpsids, DenseMatrix< double > &inverse_jacobian) const |
virtual double | local_to_lagrangian_mapping_diagonal (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const |
virtual void | assign_solid_local_eqn_numbers (const bool &store_local_dof) |
Assign local equation numbers for the solid equations in the element. More... | |
void | describe_solid_local_dofs (std::ostream &out, const std::string ¤t_string) const |
Classifies dofs locally for solid specific aspects. More... | |
void | fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian) |
virtual void | fill_in_jacobian_from_solid_position_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian) |
void | fill_in_jacobian_from_solid_position_by_fd (DenseMatrix< double > &jacobian) |
virtual void | update_before_solid_position_fd () |
virtual void | reset_after_solid_position_fd () |
virtual void | update_in_solid_position_fd (const unsigned &i) |
virtual void | reset_in_solid_position_fd (const unsigned &i) |
![]() | |
virtual void | assemble_local_to_eulerian_jacobian (const DShape &dpsids, DenseMatrix< double > &jacobian) const |
virtual void | assemble_local_to_eulerian_jacobian2 (const DShape &d2psids, DenseMatrix< double > &jacobian2) const |
virtual void | assemble_eulerian_base_vectors (const DShape &dpsids, DenseMatrix< double > &interpolated_G) const |
template<unsigned DIM> | |
double | invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const |
virtual double | invert_jacobian_mapping (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const |
virtual double | local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const |
double | local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &inverse_jacobian) const |
virtual double | local_to_eulerian_mapping_diagonal (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const |
virtual void | dJ_eulerian_dnodal_coordinates (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const |
template<unsigned DIM> | |
void | dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const |
virtual void | d_dshape_eulerian_dnodal_coordinates (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const |
template<unsigned DIM> | |
void | d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const |
virtual void | transform_derivatives (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const |
void | transform_derivatives_diagonal (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const |
virtual void | transform_second_derivatives (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const |
template<unsigned DIM> | |
void | transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const |
template<unsigned DIM> | |
void | transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const |
virtual void | fill_in_jacobian_from_nodal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian) |
void | fill_in_jacobian_from_nodal_by_fd (DenseMatrix< double > &jacobian) |
virtual void | update_before_nodal_fd () |
virtual void | reset_after_nodal_fd () |
virtual void | update_in_nodal_fd (const unsigned &i) |
virtual void | reset_in_nodal_fd (const unsigned &i) |
template<> | |
double | invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const |
Zero-d specialisation of function to calculate inverse of jacobian mapping. More... | |
template<> | |
double | invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const |
One-d specialisation of function to calculate inverse of jacobian mapping. More... | |
template<> | |
double | invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const |
Two-d specialisation of function to calculate inverse of jacobian mapping. More... | |
template<> | |
double | invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const |
template<> | |
void | dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const |
template<> | |
void | dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const |
template<> | |
void | dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const |
template<> | |
void | dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const |
template<> | |
void | d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const |
template<> | |
void | d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const |
template<> | |
void | d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const |
template<> | |
void | d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const |
template<> | |
void | transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const |
template<> | |
void | transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const |
template<> | |
void | transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const |
template<> | |
void | transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const |
![]() | |
unsigned | add_internal_data (Data *const &data_pt, const bool &fd=true) |
bool | internal_data_fd (const unsigned &i) const |
void | exclude_internal_data_fd (const unsigned &i) |
void | include_internal_data_fd (const unsigned &i) |
void | clear_global_eqn_numbers () |
void | add_global_eqn_numbers (std::deque< unsigned long > const &global_eqn_numbers, std::deque< double * > const &global_dof_pt) |
virtual void | assign_internal_and_external_local_eqn_numbers (const bool &store_local_dof_pt) |
virtual void | assign_additional_local_eqn_numbers () |
int | internal_local_eqn (const unsigned &i, const unsigned &j) const |
int | external_local_eqn (const unsigned &i, const unsigned &j) |
virtual void | fill_in_contribution_to_residuals (Vector< double > &residuals) |
void | fill_in_jacobian_from_internal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false) |
void | fill_in_jacobian_from_internal_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false) |
void | fill_in_jacobian_from_external_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false) |
void | fill_in_jacobian_from_external_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false) |
virtual void | update_before_internal_fd () |
virtual void | reset_after_internal_fd () |
virtual void | update_in_internal_fd (const unsigned &i) |
virtual void | reset_in_internal_fd (const unsigned &i) |
virtual void | update_before_external_fd () |
virtual void | reset_after_external_fd () |
virtual void | update_in_external_fd (const unsigned &i) |
virtual void | reset_in_external_fd (const unsigned &i) |
virtual void | fill_in_contribution_to_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix) |
virtual void | fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix) |
virtual void | fill_in_contribution_to_dresiduals_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam) |
virtual void | fill_in_contribution_to_djacobian_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam) |
virtual void | fill_in_contribution_to_djacobian_and_dmass_matrix_dparameter (double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam) |
virtual void | fill_in_contribution_to_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product) |
virtual void | fill_in_contribution_to_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product) |
virtual void | fill_in_contribution_to_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector) |
A base class for elements that solve the equations of solid mechanics, based on the principle of virtual displacements in Cartesian coordinates. Combines a few generic functions that are shared by PVDEquations and PVDEquationsWithPressure.
typedef void(* oomph::PVDEquationsBase< DIM >::BodyForceFctPt) (const double &t, const Vector< double > &xi, Vector< double > &b) |
Function pointer to function that specifies the body force as a function of the Lagrangian coordinates and time FCT(t,xi,b) – xi and b are Vectors!
typedef void(* oomph::PVDEquationsBase< DIM >::IsotropicGrowthFctPt) (const Vector< double > &xi, double &gamma) |
Function pointer to function that specifies the isotropic growth as a function of the Lagrangian coordinates FCT(xi,gamma(xi)) – xi is a Vector!
typedef double(* oomph::PVDEquationsBase< DIM >::PrestressFctPt) (const unsigned &i, const unsigned &j, const Vector< double > &xi) |
Function pointer to function that specifies the pre-stress sigma_0(i,j) as a function of the Lagrangian coordinates FCT(i,j,xi) – xi is a Vector!
|
inline |
Constructor: Set null pointers for constitutive law and for isotropic growth function. Set physical parameter values to default values, enable inertia and set body force to zero. Default evaluation of Jacobian: analytically rather than by FD.
|
inline |
Evaluate body force at Lagrangian coordinate xi at present time (returns zero vector if no body force function pointer has been set)
References b, oomph::PVDEquationsBase< DIM >::Body_force_fct_pt, oomph::FiniteElement::dim(), i, n, oomph::FiniteElement::node_pt(), oomph::Time::time(), oomph::TimeStepper::time_pt(), and oomph::Data::time_stepper_pt().
Referenced by oomph::RefineableQDPVDElement< DIM, NNODE_1D >::fill_in_generic_contribution_to_residuals_pvd().
|
inline |
Access function: Pointer to body force function.
References oomph::PVDEquationsBase< DIM >::Body_force_fct_pt.
Referenced by oomph::RefineablePVDEquations< DIM >::further_build(), and oomph::RefineablePVDEquationsWithPressure< DIM >::further_build().
|
inline |
Access function: Pointer to body force function (const version)
References oomph::PVDEquationsBase< DIM >::Body_force_fct_pt.
|
inline |
Return the constitutive law pointer.
References oomph::PVDEquationsBase< DIM >::Constitutive_law_pt.
Referenced by oomph::RefineablePVDEquations< DIM >::further_build(), and oomph::RefineablePVDEquationsWithPressure< DIM >::further_build().
|
inline |
Set Jacobian to be evaluated analytically Else: by FD.
References oomph::PVDEquationsBase< DIM >::Evaluate_jacobian_by_fd.
|
inline |
Switch off solid inertia.
References oomph::PVDEquationsBase< DIM >::Unsteady.
|
inline |
Set Jacobian to be evaluated by FD? Else: Analytically.
References oomph::PVDEquationsBase< DIM >::Evaluate_jacobian_by_fd.
|
inline |
Switch on solid inertia.
References oomph::PVDEquationsBase< DIM >::Unsteady.
void oomph::PVDEquationsBase< DIM >::get_deformed_covariant_basis_vectors | ( | const Vector< double > & | s, |
DenseMatrix< double > & | def_covariant_basis | ||
) |
Return the deformed covariant basis vectors at specified local coordinate: def_covariant_basis(i,j)
is the j-th component of the i-th basis vector.
|
inlinevirtual |
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contains the global equation number of the unknown, while the second one contains the number of the "DOF" that this unknown is associated with. (Function can obviously only be called if the equation numbering scheme has been set up.) E.g. in a 3D problem there are 3 types of DOF: 0 - x displacement 1 - y displacement 2 - z displacement
Reimplemented from oomph::GeneralisedElement.
Reimplemented in oomph::PVDEquationsWithPressure< DIM >, oomph::PseudoSolidNodeUpdateElement< TCrouzeixRaviartElement< 2 >, TPVDBubbleEnrichedElement< 2, 3 > >, oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >, oomph::PseudoSolidNodeUpdateElement< GeneralisedNewtonianAxisymmetricTTaylorHoodElement, TPVDElement< 2, 3 > >, oomph::PseudoSolidNodeUpdateElement< AxisymmetricTTaylorHoodElement, TPVDElement< 2, 3 > >, and oomph::PseudoSolidNodeUpdateElement< GeneralisedNewtonianTTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.
References oomph::GeneralisedElement::eqn_number(), i, k, n, oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_dimension(), and oomph::SolidFiniteElement::position_local_eqn().
void oomph::PVDEquationsBase< DIM >::get_energy | ( | double & | pot_en, |
double & | kin_en | ||
) |
Get potential (strain) and kinetic energy.
References DIM, mathsFunc::gamma(), i, J, j, k, s, calibrate::sigma, w, and oomph::QuadTreeNames::W.
|
inlinevirtual |
Evaluate isotropic growth function at Lagrangian coordinate xi and/or local coordinate s. (returns 1, i.e. no growth, if no function pointer has been set) This function is virtual to allow overloading in multi-physics problems where the growth function might be determined by another system of equations
Reimplemented in QThermalPVDElement< DIM >.
References mathsFunc::gamma(), and oomph::PVDEquationsBase< DIM >::Isotropic_growth_fct_pt.
Referenced by oomph::RefineableQDPVDElement< DIM, NNODE_1D >::fill_in_generic_contribution_to_residuals_pvd().
void oomph::PVDEquationsBase< DIM >::get_principal_stress | ( | const Vector< double > & | s, |
DenseMatrix< double > & | principal_stress_vector, | ||
Vector< double > & | principal_stress | ||
) |
Compute principal stress vectors and (scalar) principal stresses at specified local coordinate. principal_stress_vector(i,j)
is the j-th component of the i-th principal stress vector.
/////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////// Compute principal stress vectors and (scalar) principal stresses at specified local coordinate: principal_stress_vector(i,j)
is the j-th component of the i-th principal stress vector.
References DIM, i, j, k, Eigen::bfloat16_impl::pow(), s, calibrate::sigma, oomph::DoubleMatrixBase::solve(), and sqrt().
void oomph::PVDEquationsBase< DIM >::get_strain | ( | const Vector< double > & | s, |
DenseMatrix< double > & | strain | ||
) | const |
Return the strain tensor.
Compute the strain tensor at local coordinate s.
/////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////
References DIM, G, mathsFunc::gamma(), i, j, k, oomph::DenseMatrix< T >::ncol(), oomph::DenseMatrix< T >::nrow(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, Eigen::bfloat16_impl::pow(), and s.
Referenced by oomph::RefineablePVDEquations< DIM >::get_Z2_flux(), oomph::RefineablePVDEquationsWithPressure< DIM >::get_Z2_flux(), oomph::TPVDElement< DIM, NNODE_1D >::get_Z2_flux(), oomph::TPVDBubbleEnrichedElement< DIM, NNODE_1D >::get_Z2_flux(), and oomph::TPVDElementWithContinuousPressure< DIM >::get_Z2_flux().
|
pure virtual |
Return the 2nd Piola Kirchoff stress tensor, as calculated from the constitutive law at specified local coordinate (needed by get_principal_stress(...)
, so I'm afraid I will have to insist that you implement it...
Implemented in oomph::PVDEquationsWithPressure< DIM >, and oomph::PVDEquations< DIM >.
|
inline |
Access function to flag that switches inertia on/off (const version)
References oomph::PVDEquationsBase< DIM >::Unsteady.
Referenced by oomph::RefineablePVDEquations< DIM >::further_build(), and oomph::RefineablePVDEquationsWithPressure< DIM >::further_build().
|
inline |
Return the flag indicating whether the jacobian is evaluated by fd.
References oomph::PVDEquationsBase< DIM >::Evaluate_jacobian_by_fd.
Referenced by oomph::RefineablePVDEquations< DIM >::further_build(), and oomph::RefineablePVDEquationsWithPressure< DIM >::further_build().
|
inline |
Access function: Pointer to isotropic growth function.
References oomph::PVDEquationsBase< DIM >::Isotropic_growth_fct_pt.
Referenced by oomph::RefineablePVDEquations< DIM >::further_build(), and oomph::RefineablePVDEquationsWithPressure< DIM >::further_build().
|
inline |
Access function: Pointer to isotropic growth function (const version)
References oomph::PVDEquationsBase< DIM >::Isotropic_growth_fct_pt.
|
inline |
Access function for timescale ratio (nondim density)
References oomph::PVDEquationsBase< DIM >::Lambda_sq_pt.
Referenced by oomph::RefineableQDPVDElement< DIM, NNODE_1D >::fill_in_generic_contribution_to_residuals_pvd().
|
inline |
Access function for pointer to timescale ratio (nondim density)
References oomph::PVDEquationsBase< DIM >::Lambda_sq_pt.
Referenced by oomph::RefineablePVDEquations< DIM >::further_build(), and oomph::RefineablePVDEquationsWithPressure< DIM >::further_build().
|
inlinevirtual |
returns the number of DOF types associated with this element.
Reimplemented from oomph::GeneralisedElement.
Reimplemented in oomph::PVDEquationsWithPressure< DIM >, oomph::PseudoSolidNodeUpdateElement< TCrouzeixRaviartElement< 2 >, TPVDBubbleEnrichedElement< 2, 3 > >, oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >, oomph::PseudoSolidNodeUpdateElement< GeneralisedNewtonianAxisymmetricTTaylorHoodElement, TPVDElement< 2, 3 > >, oomph::PseudoSolidNodeUpdateElement< AxisymmetricTTaylorHoodElement, TPVDElement< 2, 3 > >, and oomph::PseudoSolidNodeUpdateElement< GeneralisedNewtonianTTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.
References DIM.
|
inlinevirtual |
Return the number of solid pressure degrees of freedom Default is that there are no solid pressures
Reimplemented in oomph::TPVDElementWithContinuousPressure< DIM >, oomph::TPVDElementWithContinuousPressure< DIM >, oomph::TPVDElementWithContinuousPressure< DIM >, oomph::QPVDElementWithContinuousPressure< DIM >, and oomph::QPVDElementWithPressure< DIM >.
Referenced by oomph::PVDEquationsWithPressure< DIM >::get_dof_numbers_for_unknowns(), and oomph::PVDEquationsWithPressure< DIM >::interpolated_solid_p().
|
inlinevirtual |
Pin the element's redundant solid pressures (needed for refinement)
Reimplemented in oomph::RefineableQPVDElementWithContinuousPressure< DIM >.
Referenced by oomph::PVDEquationsBase< DIM >::pin_redundant_nodal_solid_pressures().
|
inlinestatic |
Loop over all elements in Vector (which typically contains all the elements in a refineable solid mesh) and pin the nodal solid pressure degrees of freedom that are not being used. Function uses the member function
PVDEquationsBase<DIM>::
pin_elemental_redundant_nodal_pressure_dofs()which is empty by default and should be implemented for elements with nodal solid pressure degrees of freedom (e.g. solid elements with continuous pressure interpolation.)
References e(), and oomph::PVDEquationsBase< DIM >::pin_elemental_redundant_nodal_solid_pressures().
Referenced by SolidProblem< ELEMENT_TYPE >::prepareForSolve().
|
inline |
Return (i,j)-th component of second Piola Kirchhoff membrane prestress at Lagrangian coordinate xi
References i, j, and oomph::PVDEquationsBase< DIM >::Prestress_fct_pt.
Referenced by oomph::RefineableQDPVDElement< DIM, NNODE_1D >::fill_in_generic_contribution_to_residuals_pvd().
|
inline |
Access function: Pointer to pre-stress function.
References oomph::PVDEquationsBase< DIM >::Prestress_fct_pt.
|
inlinevirtual |
Return the local degree of freedom associated with the i-th solid pressure. Default is that there are none.
Reimplemented in oomph::TPVDElementWithContinuousPressure< DIM >, oomph::QPVDElementWithContinuousPressure< DIM >, and oomph::QPVDElementWithPressure< DIM >.
Referenced by oomph::PVDEquationsWithPressure< DIM >::get_dof_numbers_for_unknowns().
|
inlinevirtual |
Return the index at which the solid pressure is stored if it is stored at the nodes. If not stored at the nodes this will return a negative number.
Reimplemented in oomph::TPVDElementWithContinuousPressure< DIM >, oomph::QPVDElementWithContinuousPressure< DIM >, oomph::PseudoSolidNodeUpdateElement< TCrouzeixRaviartElement< 2 >, TPVDBubbleEnrichedElement< 2, 3 > >, oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >, oomph::PseudoSolidNodeUpdateElement< GeneralisedNewtonianAxisymmetricTTaylorHoodElement, TPVDElement< 2, 3 > >, oomph::PseudoSolidNodeUpdateElement< AxisymmetricTTaylorHoodElement, TPVDElement< 2, 3 > >, and oomph::PseudoSolidNodeUpdateElement< GeneralisedNewtonianTTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.
References oomph::PVDEquationsBase< DIM >::Solid_pressure_not_stored_at_node.
|
inlinestatic |
Unpin all pressure dofs in elements listed in vector.
References e(), and oomph::PVDEquationsBase< DIM >::unpin_elemental_solid_pressure_dofs().
|
pure virtual |
Unpin all solid pressure dofs in the element.
Implemented in oomph::TPVDElementWithContinuousPressure< DIM >, oomph::QPVDElementWithContinuousPressure< DIM >, oomph::QPVDElementWithPressure< DIM >, oomph::PVDEquations< DIM >, oomph::RefineableQPVDElementWithContinuousPressure< DIM >, and oomph::RefineableQPVDElementWithPressure< DIM >.
Referenced by oomph::PVDEquationsBase< DIM >::unpin_all_solid_pressure_dofs().
|
protected |
Pointer to body force function.
Referenced by oomph::PVDEquationsBase< DIM >::body_force(), oomph::PVDEquationsBase< DIM >::body_force_fct_pt(), oomph::RefineablePVDEquations< DIM >::further_build(), and oomph::RefineablePVDEquationsWithPressure< DIM >::further_build().
|
protected |
Pointer to the constitutive law.
Referenced by oomph::PVDEquationsBase< DIM >::constitutive_law_pt(), oomph::RefineableQDPVDElement< DIM, NNODE_1D >::fill_in_generic_contribution_to_residuals_pvd(), oomph::RefineablePVDEquations< DIM >::further_build(), oomph::RefineablePVDEquationsWithPressure< DIM >::further_build(), oomph::PVDEquationsWithPressure< DIM >::get_d_stress_dG_upper(), oomph::PVDEquations< DIM >::get_d_stress_dG_upper(), oomph::PVDEquations< DIM >::get_stress(), and oomph::PVDEquationsWithPressure< DIM >::get_stress().
|
staticprotected |
Static default value for timescale ratio (1.0 – for natural scaling)
|
protected |
Use FD to evaluate Jacobian.
Referenced by oomph::PVDEquationsBase< DIM >::disable_evaluate_jacobian_by_fd(), oomph::PVDEquationsBase< DIM >::enable_evaluate_jacobian_by_fd(), oomph::PVDEquations< DIM >::fill_in_contribution_to_jacobian(), oomph::PVDEquationsWithPressure< DIM >::fill_in_contribution_to_jacobian(), oomph::PVDEquationsWithPressure< DIM >::fill_in_contribution_to_jacobian_and_mass_matrix(), oomph::RefineablePVDEquations< DIM >::further_build(), oomph::RefineablePVDEquationsWithPressure< DIM >::further_build(), and oomph::PVDEquationsBase< DIM >::is_jacobian_evaluated_by_fd().
|
protected |
Pointer to isotropic growth function.
Referenced by oomph::RefineablePVDEquations< DIM >::further_build(), oomph::RefineablePVDEquationsWithPressure< DIM >::further_build(), oomph::PVDEquationsBase< DIM >::get_isotropic_growth(), and oomph::PVDEquationsBase< DIM >::isotropic_growth_fct_pt().
|
protected |
Timescale ratio (non-dim. density)
Referenced by oomph::RefineablePVDEquations< DIM >::further_build(), oomph::RefineablePVDEquationsWithPressure< DIM >::further_build(), oomph::PVDEquationsBase< DIM >::lambda_sq(), and oomph::PVDEquationsBase< DIM >::lambda_sq_pt().
|
protected |
Pointer to prestress function.
Referenced by oomph::PVDEquationsBase< DIM >::prestress(), and oomph::PVDEquationsBase< DIM >::prestress_fct_pt().
|
staticprivate |
Static "magic" number that indicates that the solid pressure is not stored at a node
/////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////// "Magic" number that indicates that the solid pressure is not stored at a node. It is a negative number that cannot be -1 because that is used to represent the positional hanging scheme in Hanging_pt objects
Referenced by oomph::PVDEquationsBase< DIM >::solid_p_nodal_index().
|
protected |
Flag that switches inertia on/off.
Referenced by oomph::PVDEquationsBase< DIM >::disable_inertia(), oomph::PVDEquationsBase< DIM >::enable_inertia(), oomph::RefineableQDPVDElement< DIM, NNODE_1D >::fill_in_generic_contribution_to_residuals_pvd(), oomph::RefineablePVDEquations< DIM >::further_build(), oomph::RefineablePVDEquationsWithPressure< DIM >::further_build(), and oomph::PVDEquationsBase< DIM >::is_inertia_enabled().