oomph::AxisymmetricPVDEquations Class Reference

#include <axisym_solid_elements.h>

+ Inheritance diagram for oomph::AxisymmetricPVDEquations:

Public Member Functions

 AxisymmetricPVDEquations ()
 Constructor. More...
 
ConstitutiveLaw *& constitutive_law_pt ()
 Return the constitutive law pointer. More...
 
void get_stress (const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma)
 Return the stress tensor, as calculated from the constitutive law. More...
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Fill in the residuals by calling the generic function. More...
 
void fill_in_contribution_to_residuals_axisym_pvd (Vector< double > &residuals)
 Return the residuals for the equations of solid mechanics. More...
 
double compute_physical_size () const
 Overload/implement the function to calculate the volume of the element. More...
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Assign the contribution to the residual using only finite differences. More...
 
void output (std::ostream &outfile)
 Overload the output function. More...
 
void output (std::ostream &outfile, const unsigned &n_plot)
 Output function. More...
 
void output (FILE *file_pt)
 Overload the output function. More...
 
void output (FILE *file_pt, const unsigned &n_plot)
 Output function. More...
 
- Public Member Functions inherited from oomph::SolidFiniteElement
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...
 
Datageom_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)
 
MacroElementundeformed_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
 
Nodeconstruct_node (const unsigned &n)
 Construct the local node n and return a pointer to it. More...
 
Nodeconstruct_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
Nodeconstruct_boundary_node (const unsigned &n)
 
Nodeconstruct_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 &current_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...
 
MultiplierFctPtmultiplier_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
 
- Public Member Functions inherited from oomph::FiniteElement
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)
 
MacroElementmacro_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 &current_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 Nodevertex_node_pt (const unsigned &j) const
 
int get_node_number (Node *const &node_pt) const
 
virtual Nodeget_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 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 (const unsigned &t, std::ostream &outfile, const unsigned &n_plot) const
 
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
 
- Public Member Functions inherited from oomph::GeneralisedElement
 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 &current_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 &parameter_pt, Vector< double > &dres_dparam)
 
virtual void get_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void get_djacobian_and_dmass_matrix_dparameter (double *const &parameter_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)
 
virtual unsigned ndof_types () const
 
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 
- Public Member Functions inherited from oomph::GeomObject
 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 ()
 
TimeSteppertime_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
 

Private Attributes

ConstitutiveLawConstitutive_law_pt
 Pointer to constitutive law. More...
 

Additional Inherited Members

- Public Types inherited from oomph::SolidFiniteElement
typedef double(* MultiplierFctPt) (const Vector< double > &xi)
 
- Public Types inherited from oomph::FiniteElement
typedef void(* SteadyExactSolutionFctPt) (const Vector< double > &, Vector< double > &)
 
typedef void(* UnsteadyExactSolutionFctPt) (const double &, const Vector< double > &, Vector< double > &)
 
- Static Public Attributes inherited from oomph::FiniteElement
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 Public Attributes inherited from oomph::GeneralisedElement
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
 
- Protected Member Functions inherited from oomph::SolidFiniteElement
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 &current_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)
 
- Protected Member Functions inherited from oomph::FiniteElement
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
 
- Protected Member Functions inherited from oomph::GeneralisedElement
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)
 
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 &parameter_pt, Vector< double > &dres_dparam)
 
virtual void fill_in_contribution_to_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void fill_in_contribution_to_djacobian_and_dmass_matrix_dparameter (double *const &parameter_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)
 
- Protected Attributes inherited from oomph::SolidFiniteElement
MacroElementUndeformed_macro_elem_pt
 Pointer to the element's "undeformed" macro element (NULL by default) More...
 
SolidInitialConditionSolid_ic_pt
 Pointer to object that specifies the initial condition. More...
 
bool Solve_for_consistent_newmark_accel_flag
 
- Protected Attributes inherited from oomph::FiniteElement
MacroElementMacro_elem_pt
 Pointer to the element's macro element (NULL by default) More...
 
- Protected Attributes inherited from oomph::GeomObject
unsigned NLagrangian
 Number of Lagrangian (intrinsic) coordinates. More...
 
unsigned Ndim
 Number of Eulerian coordinates. More...
 
TimeStepperGeom_object_time_stepper_pt
 
- Static Protected Attributes inherited from oomph::FiniteElement
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 Protected Attributes inherited from oomph::GeneralisedElement
static DenseMatrix< doubleDummy_matrix
 
static std::deque< double * > Dof_pt_deque
 

Detailed Description

A class for elements that solve the equations of solid mechanics, based on the principle of virtual displacements in an axisymmetric formulation. In this case x[0] is the component of displacement in the radial direction and x[1] is that in the theta direction.

Constructor & Destructor Documentation

◆ AxisymmetricPVDEquations()

oomph::AxisymmetricPVDEquations::AxisymmetricPVDEquations ( )
inline

Constructor.

57 : Constitutive_law_pt(0) {}
ConstitutiveLaw * Constitutive_law_pt
Pointer to constitutive law.
Definition: axisym_solid_elements.h:53

Member Function Documentation

◆ compute_physical_size()

double oomph::AxisymmetricPVDEquations::compute_physical_size ( ) const
inlinevirtual

Overload/implement the function to calculate the volume of the element.

Reimplemented from oomph::FiniteElement.

285  {
286  unsigned n_node = nnode();
287  unsigned n_position_type = 1;
288 
289  // Set up memory for the shape functions
290  Shape psi(n_node, n_position_type);
291  DShape dpsidxi(n_node, n_position_type, 2);
292 
293  // Set sum to zero
294  double sum = 0.0;
295  // Set the value of n_intpt
296  unsigned n_intpt = integral_pt()->nweight();
297 
298  // Loop over the integration points
299  // Loop in s1 direction*
300  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
301  {
302  // Get the integral weight
303  double w = integral_pt()->weight(ipt);
304  // Call the derivatives of the shape function wrt lagrangian coordinates
305  double J = dshape_lagrangian_at_knot(ipt, psi, dpsidxi);
306  // Premultiply the weights and the Jacobian
307  double W = w * J;
308 
309  // Calculate the local Lagrangian coordinates, position components
310  // and the derivatives of global position components
311  // wrt lagrangian coordinates
312  double interpolated_xi[2] = {0.0, 0.0};
313  double interpolated_X[2] = {0.0, 0.0};
314  double interpolated_dXdxi[2][2];
315 
316  // Initialise interpolated_dXdxi to zero
317  for (unsigned i = 0; i < 2; i++)
318  {
319  for (unsigned j = 0; j < 2; j++)
320  {
321  interpolated_dXdxi[i][j] = 0.0;
322  }
323  }
324 
325  // Calculate displacements and derivatives
326  for (unsigned l = 0; l < n_node; l++)
327  {
328  // Loop over positional dofs
329  for (unsigned k = 0; k < n_position_type; k++)
330  {
331  // Loop over displacement components (deformed position)
332  for (unsigned i = 0; i < 2; i++)
333  {
334  // Set the value of the lagrangian coordinate
335  interpolated_xi[i] +=
336  lagrangian_position_gen(l, k, i) * psi(l, k);
337  // Set the value of the position component
338  interpolated_X[i] += nodal_position_gen(l, k, i) * psi(l, k);
339  // Loop over Lagrangian derivative directions
340  for (unsigned j = 0; j < 2; j++)
341  {
342  // Calculate dX[i]/dxi_{j}
343  interpolated_dXdxi[i][j] +=
344  nodal_position_gen(l, k, i) * dpsidxi(l, k, j);
345  }
346  }
347  }
348  }
349 
350  // Now calculate the deformed metric tensor
352  // r row
353  G(0, 0) = interpolated_dXdxi[0][0] * interpolated_dXdxi[0][0] +
354  interpolated_dXdxi[1][0] * interpolated_dXdxi[1][0];
355  G(0, 1) = interpolated_dXdxi[0][0] *
356  (interpolated_dXdxi[0][1] - interpolated_X[1]) +
357  interpolated_dXdxi[1][0] *
358  (interpolated_dXdxi[1][1] + interpolated_X[0]);
359  G(0, 2) = 0.0;
360  // theta row
361  G(1, 0) = G(0, 1);
362  G(1, 1) = (interpolated_dXdxi[0][1] - interpolated_X[1]) *
363  (interpolated_dXdxi[0][1] - interpolated_X[1]) +
364  (interpolated_dXdxi[1][1] + interpolated_X[0]) *
365  (interpolated_dXdxi[1][1] + interpolated_X[0]);
366  G(1, 2) = 0.0;
367  // phi row
368  G(2, 0) = 0.0;
369  G(2, 1) = 0.0;
370  G(2, 2) = (interpolated_X[0] * sin(interpolated_xi[1]) +
371  interpolated_X[1] * cos(interpolated_xi[1])) *
372  (interpolated_X[0] * sin(interpolated_xi[1]) +
373  interpolated_X[1] * cos(interpolated_xi[1]));
374 
375  // Calculate the determinant of the metric tensor
376  double detG = G(0, 0) * G(1, 1) * G(2, 2) - G(0, 1) * G(1, 0) * G(2, 2);
377 
378  // Add the appropriate weight to the integral, i.e. sqrt of metric
379  // tensor
380  sum += W * sqrt(detG);
381  }
382 
383  // Return the volume, need to multiply by 2pi
384  return (2.0 * MathematicalConstants::Pi * sum);
385  }
AnnoyingScalar cos(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:136
AnnoyingScalar sin(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:137
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
JacobiRotation< float > G
Definition: Jacobi_makeGivens.cpp:2
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
RowVector3d w
Definition: Matrix_resize_int.cpp:3
double nodal_position_gen(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: elements.h:2349
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
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.
double lagrangian_position_gen(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: elements.h:3912
virtual double interpolated_xi(const Vector< double > &s, const unsigned &i) const
Definition: elements.cc:7104
virtual double dshape_lagrangian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidxi) const
Definition: elements.cc:6737
char char char int int * k
Definition: level2_impl.h:374
const double Pi
50 digits from maple
Definition: oomph_utilities.h:157
@ W
Definition: quadtree.h:63
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References cos(), oomph::SolidFiniteElement::dshape_lagrangian_at_knot(), G, i, oomph::FiniteElement::integral_pt(), oomph::SolidFiniteElement::interpolated_xi(), J, j, k, oomph::SolidFiniteElement::lagrangian_position_gen(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_position_gen(), oomph::Integral::nweight(), oomph::MathematicalConstants::Pi, sin(), sqrt(), w, oomph::QuadTreeNames::W, and oomph::Integral::weight().

◆ constitutive_law_pt()

ConstitutiveLaw*& oomph::AxisymmetricPVDEquations::constitutive_law_pt ( )
inline

Return the constitutive law pointer.

61  {
62  return Constitutive_law_pt;
63  }

References Constitutive_law_pt.

◆ fill_in_contribution_to_jacobian()

void oomph::AxisymmetricPVDEquations::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
inlinevirtual

Assign the contribution to the residual using only finite differences.

Reimplemented from oomph::GeneralisedElement.

390  {
391  // Add the solid contribution to the residuals
393 
394  // Solve for the consistent acceleration in Newmark scheme?
396  {
398  return;
399  }
400 
401  // Get the solid entries in the jacobian using finite differences
403  }
void fill_in_contribution_to_residuals_axisym_pvd(Vector< double > &residuals)
Return the residuals for the equations of solid mechanics.
Definition: axisym_solid_elements.h:94
bool Solve_for_consistent_newmark_accel_flag
Definition: elements.h:4302
virtual void fill_in_jacobian_from_solid_position_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: elements.cc:6985
void fill_in_jacobian_for_newmark_accel(DenseMatrix< double > &jacobian)
Definition: elements.cc:7227

References fill_in_contribution_to_residuals_axisym_pvd(), oomph::SolidFiniteElement::fill_in_jacobian_for_newmark_accel(), oomph::SolidFiniteElement::fill_in_jacobian_from_solid_position_by_fd(), and oomph::SolidFiniteElement::Solve_for_consistent_newmark_accel_flag.

◆ fill_in_contribution_to_residuals()

void oomph::AxisymmetricPVDEquations::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inlinevirtual

Fill in the residuals by calling the generic function.

Reimplemented from oomph::GeneralisedElement.

89  {
91  }

References fill_in_contribution_to_residuals_axisym_pvd().

◆ fill_in_contribution_to_residuals_axisym_pvd()

void oomph::AxisymmetricPVDEquations::fill_in_contribution_to_residuals_axisym_pvd ( Vector< double > &  residuals)
inline

Return the residuals for the equations of solid mechanics.

96  {
97  // Set the number of Lagrangian coordinates
98  unsigned n_lagrangian = 2;
99  // Find out how many nodes there are
100  unsigned n_node = nnode();
101  // Find out how many positional dofs there are
102  unsigned n_position_type = nnodal_position_type();
103 
104  // Integer to store local equation number
105  int local_eqn = 0;
106 
107  // Set up memory for the shape functions
108  Shape psi(n_node, n_position_type);
109  DShape dpsidxi(n_node, n_position_type, n_lagrangian);
110 
111  // Set the value of n_intpt
112  unsigned n_intpt = integral_pt()->nweight();
113 
114  // Loop over the integration points
115  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
116  {
117  // Get the integral weight
118  double w = integral_pt()->weight(ipt);
119  // Call the derivatives of the shape functions
120  double J = dshape_lagrangian_at_knot(ipt, psi, dpsidxi);
121  // Premultiply the weights and the Jacobian
122  double W = w * J;
123 
124  // Calculate the local Lagrangian coordinates, position components
125  // and the derivatives of global position components
126  // wrt lagrangian coordinates
127  double interpolated_xi[2] = {0.0, 0.0};
128  double interpolated_X[2] = {0.0, 0.0};
129  double interpolated_dXdxi[2][2];
130 
131  // Initialise interpolated_dXdxi to zero
132  for (unsigned i = 0; i < 2; i++)
133  {
134  for (unsigned j = 0; j < 2; j++)
135  {
136  interpolated_dXdxi[i][j] = 0.0;
137  }
138  }
139 
140  // Calculate displacements and derivatives
141  for (unsigned l = 0; l < n_node; l++)
142  {
143  // Loop over positional dofs
144  for (unsigned k = 0; k < n_position_type; k++)
145  {
146  // Loop over displacement components (deformed position)
147  for (unsigned i = 0; i < 2; i++)
148  {
149  // Set the value of the lagrangian coordinate
150  interpolated_xi[i] +=
151  lagrangian_position_gen(l, k, i) * psi(l, k);
152  // Set the value of the position component
153  interpolated_X[i] += nodal_position_gen(l, k, i) * psi(l, k);
154  // Loop over Lagrangian derivative directions
155  for (unsigned j = 0; j < 2; j++)
156  {
157  // Calculate dX[i]/dxi_{j}
158  interpolated_dXdxi[i][j] +=
159  nodal_position_gen(l, k, i) * dpsidxi(l, k, j);
160  }
161  }
162  }
163  }
164 
165  // We are now in a position to calculate the undeformed metric tensor
166  DenseMatrix<double> g(3);
167  // r row
168  g(0, 0) = 1.0;
169  g(0, 1) = 0.0;
170  g(0, 2) = 0.0;
171  // theta row
172  g(1, 0) = 0.0;
173  g(1, 1) = interpolated_xi[0] * interpolated_xi[0];
174  g(1, 2) = 0.0;
175  // phi row
176  g(2, 0) = 0.0;
177  g(2, 1) = 0.0;
178  g(2, 2) = interpolated_xi[0] * interpolated_xi[0] *
180 
181  // Now multiply the weight by the square-root of the undeformed metric
182  // tensor r^2 sin(theta)
183  W *= sqrt(g(0, 0) * g(1, 1) * g(2, 2));
184 
185  // Now calculate the deformed metric tensor
187  // r row
188  G(0, 0) = interpolated_dXdxi[0][0] * interpolated_dXdxi[0][0] +
189  interpolated_dXdxi[1][0] * interpolated_dXdxi[1][0];
190  G(0, 1) = interpolated_dXdxi[0][0] *
191  (interpolated_dXdxi[0][1] - interpolated_X[1]) +
192  interpolated_dXdxi[1][0] *
193  (interpolated_dXdxi[1][1] + interpolated_X[0]);
194  G(0, 2) = 0.0;
195  // theta row
196  G(1, 0) = G(0, 1);
197  G(1, 1) = (interpolated_dXdxi[0][1] - interpolated_X[1]) *
198  (interpolated_dXdxi[0][1] - interpolated_X[1]) +
199  (interpolated_dXdxi[1][1] + interpolated_X[0]) *
200  (interpolated_dXdxi[1][1] + interpolated_X[0]);
201  G(1, 2) = 0.0;
202  // phi row
203  G(2, 0) = 0.0;
204  G(2, 1) = 0.0;
205  G(2, 2) = (interpolated_X[0] * sin(interpolated_xi[1]) +
206  interpolated_X[1] * cos(interpolated_xi[1])) *
207  (interpolated_X[0] * sin(interpolated_xi[1]) +
208  interpolated_X[1] * cos(interpolated_xi[1]));
209 
210 
211  // Now calculate the stress tensor from the constitutive law
213  get_stress(g, G, sigma);
214 
215  //=====EQUATIONS OF ELASTICITY FROM PRINCIPLE OF VIRTUAL
216  // DISPLACEMENTS========
217 
218  // Loop over the test functions, nodes of the element
219  for (unsigned l = 0; l < n_node; l++)
220  {
221  // Loop of types of dofs
222  for (unsigned k = 0; k < n_position_type; k++)
223  {
224  // Radial displacemenet component
225  unsigned i = 0;
226  local_eqn = position_local_eqn(l, k, i);
227  /*IF it's not a boundary condition*/
228  if (local_eqn >= 0)
229  {
230  // Add the term for variations in radial position
231  residuals[local_eqn] +=
232  (sigma(0, 1) * interpolated_dXdxi[1][0] +
233  sigma(1, 1) * (interpolated_dXdxi[1][1] + interpolated_X[0]) +
234  sigma(2, 2) * sin(interpolated_xi[1]) *
235  (interpolated_X[0] * sin(interpolated_xi[1]) +
236  interpolated_X[1] * cos(interpolated_xi[1]))) *
237  psi(l, k) * W;
238 
239  // Add the terms for the variations in dX_{r}/dxi_{j}
240  for (unsigned j = 0; j < 2; j++)
241  {
242  residuals[local_eqn] +=
243  (sigma(j, 0) * interpolated_dXdxi[0][0] +
244  sigma(j, 1) *
245  (interpolated_dXdxi[0][1] - interpolated_X[1])) *
246  dpsidxi(l, k, j) * W;
247  }
248  }
249 
250  // Theta displacement component
251  i = 1;
252  local_eqn = position_local_eqn(l, k, i);
253  /*IF it's not a boundary condition*/
254  if (local_eqn >= 0)
255  {
256  // Add the term for variations in azimuthal position
257  residuals[local_eqn] +=
258  (-sigma(0, 1) * interpolated_dXdxi[0][0] -
259  sigma(1, 1) * (interpolated_dXdxi[0][1] - interpolated_X[1]) +
260  sigma(2, 2) * cos(interpolated_xi[1]) *
261  (interpolated_X[0] * sin(interpolated_xi[1]) +
262  interpolated_X[1] * cos(interpolated_xi[1]))) *
263  psi(l, k) * W;
264 
265  // Add the terms for the variations in dX_{theta}/dxi_{j}
266  for (unsigned j = 0; j < 2; j++)
267  {
268  residuals[local_eqn] +=
269  (sigma(j, 0) * interpolated_dXdxi[1][0] +
270  sigma(j, 1) *
271  (interpolated_dXdxi[1][1] + interpolated_X[0])) *
272  dpsidxi(l, k, j) * W;
273  }
274  }
275  } // End of loop over type of dof
276  } // End of loop over shape functions
277  } // End of loop over integration points
278  }
void get_stress(const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma)
Return the stress tensor, as calculated from the constitutive law.
Definition: axisym_solid_elements.h:66
unsigned nnodal_position_type() const
Definition: elements.h:2463
int position_local_eqn(const unsigned &n, const unsigned &k, const unsigned &j) const
Definition: elements.h:4137
int sigma
Definition: calibrate.py:179

References cos(), oomph::SolidFiniteElement::dshape_lagrangian_at_knot(), G, get_stress(), i, oomph::FiniteElement::integral_pt(), oomph::SolidFiniteElement::interpolated_xi(), J, j, k, oomph::SolidFiniteElement::lagrangian_position_gen(), oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_position_gen(), oomph::Integral::nweight(), oomph::SolidFiniteElement::position_local_eqn(), calibrate::sigma, sin(), sqrt(), w, oomph::QuadTreeNames::W, and oomph::Integral::weight().

Referenced by fill_in_contribution_to_jacobian(), and fill_in_contribution_to_residuals().

◆ get_stress()

void oomph::AxisymmetricPVDEquations::get_stress ( const DenseMatrix< double > &  g,
const DenseMatrix< double > &  G,
DenseMatrix< double > &  sigma 
)
inline

Return the stress tensor, as calculated from the constitutive law.

69  {
70 #ifdef PARANOID
71  // If the pointer to the constitutive law hasn't been set, issue an error
72  if (Constitutive_law_pt == 0)
73  {
74  std::string error_message =
75  "Elements derived from AxisymmetricPVDEquations";
76  error_message += " must have a constitutive law :\n ";
77  error_message +=
78  "set one using the constitutive_law_pt() member function\n";
79 
80  throw OomphLibError(
82  }
83 #endif
85  }
virtual void calculate_second_piola_kirchhoff_stress(const DenseMatrix< double > &g, const DenseMatrix< double > &G, DenseMatrix< double > &sigma)=0
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References oomph::ConstitutiveLaw::calculate_second_piola_kirchhoff_stress(), Constitutive_law_pt, G, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, calibrate::sigma, and oomph::Global_string_for_annotation::string().

Referenced by fill_in_contribution_to_residuals_axisym_pvd().

◆ output() [1/4]

void oomph::AxisymmetricPVDEquations::output ( FILE *  file_pt)
inlinevirtual

Overload the output function.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::AxisymDiagHermitePVDElement, and oomph::AxisymQPVDElement.

449  {
450  FiniteElement::output(file_pt);
451  }
virtual void output(std::ostream &outfile)
Definition: elements.h:3050

References oomph::FiniteElement::output().

◆ output() [2/4]

void oomph::AxisymmetricPVDEquations::output ( FILE *  file_pt,
const unsigned n_plot 
)
inlinevirtual

Output function.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::AxisymDiagHermitePVDElement, and oomph::AxisymQPVDElement.

456  {
457  // Set the output Vector
458  Vector<double> s(2);
459 
460  // Tecplot header info
461  fprintf(file_pt, "ZONE I=%i, J=%i\n", n_plot, n_plot);
462 
463  // Loop over element nodes
464  for (unsigned l2 = 0; l2 < n_plot; l2++)
465  {
466  s[1] = -1.0 + l2 * 2.0 / (n_plot - 1);
467  for (unsigned l1 = 0; l1 < n_plot; l1++)
468  {
469  s[0] = -1.0 + l1 * 2.0 / (n_plot - 1);
470 
471  double x_r = interpolated_x(s, 0), x_theta = interpolated_x(s, 1);
472  double theta = interpolated_xi(s, 1);
473  // Output the x,y,u,v
474  // First output x and y assuming phi = 0
475  // outfile << x_r*sin(theta) + x_theta*cos(theta) << " " <<
476  // x_r*cos(theta) - x_theta*sin(theta) << " ";
477  fprintf(file_pt,
478  "%g %g ",
479  x_r * sin(theta) + x_theta * cos(theta),
480  x_r * cos(theta) - x_theta * sin(theta));
481 
482  // Now output the true variables
483  for (unsigned i = 0; i < 2; i++)
484  fprintf(file_pt, "%g ", interpolated_x(s, i));
485  for (unsigned i = 0; i < 2; i++)
486  fprintf(file_pt, "%g ", interpolated_xi(s, i));
487  fprintf(file_pt, "\n");
488  }
489  }
490  fprintf(file_pt, "\n");
491  }
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
RealScalar s
Definition: level1_cplx_impl.h:130
double theta
Definition: two_d_biharmonic.cc:236

References cos(), i, oomph::FiniteElement::interpolated_x(), oomph::SolidFiniteElement::interpolated_xi(), s, sin(), and BiharmonicTestFunctions2::theta.

◆ output() [3/4]

void oomph::AxisymmetricPVDEquations::output ( std::ostream &  outfile)
inlinevirtual

Overload the output function.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::AxisymDiagHermitePVDElement, and oomph::AxisymQPVDElement.

408  {
409  FiniteElement::output(outfile);
410  }

References oomph::FiniteElement::output().

Referenced by oomph::AxisymQPVDElement::output().

◆ output() [4/4]

void oomph::AxisymmetricPVDEquations::output ( std::ostream &  outfile,
const unsigned n_plot 
)
inlinevirtual

Output function.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::AxisymDiagHermitePVDElement, and oomph::AxisymQPVDElement.

414  {
415  // Set the output Vector
416  Vector<double> s(2);
417 
418  // Tecplot header info
419  outfile << "ZONE I=" << n_plot << ", J=" << n_plot << std::endl;
420 
421  // Loop over element nodes
422  for (unsigned l2 = 0; l2 < n_plot; l2++)
423  {
424  s[1] = -1.0 + l2 * 2.0 / (n_plot - 1);
425  for (unsigned l1 = 0; l1 < n_plot; l1++)
426  {
427  s[0] = -1.0 + l1 * 2.0 / (n_plot - 1);
428 
429  double x_r = interpolated_x(s, 0), x_theta = interpolated_x(s, 1);
430  double theta = interpolated_xi(s, 1);
431  // Output the x,y,u,v
432  // First output x and y assuming phi = 0
433  outfile << x_r * sin(theta) + x_theta * cos(theta) << " "
434  << x_r * cos(theta) - x_theta * sin(theta) << " ";
435  // Now output the true variables
436  for (unsigned i = 0; i < 2; i++)
437  outfile << interpolated_x(s, i) << " ";
438  for (unsigned i = 0; i < 2; i++)
439  outfile << interpolated_xi(s, i) << " ";
440  outfile << std::endl;
441  }
442  }
443  outfile << std::endl;
444  }

References cos(), i, oomph::FiniteElement::interpolated_x(), oomph::SolidFiniteElement::interpolated_xi(), s, sin(), and BiharmonicTestFunctions2::theta.

Member Data Documentation

◆ Constitutive_law_pt

ConstitutiveLaw* oomph::AxisymmetricPVDEquations::Constitutive_law_pt
private

Pointer to constitutive law.

Referenced by constitutive_law_pt(), and get_stress().


The documentation for this class was generated from the following file: