oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT > Class Template Reference

#include <axisym_poroelasticity_face_elements.h>

+ Inheritance diagram for oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >:

Public Member Functions

 AxisymmetricPoroelasticityTractionElement (FiniteElement *const &element_pt, const int &face_index)
 
 AxisymmetricPoroelasticityTractionElement ()
 Default constructor. More...
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Return the residuals. More...
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Fill in contribution from Jacobian. More...
 
double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 
void output (std::ostream &outfile)
 Output function. More...
 
void output (std::ostream &outfile, const unsigned &n_plot)
 Output function. More...
 
void output (FILE *file_pt)
 C_style output function. More...
 
void output (FILE *file_pt, const unsigned &n_plot)
 C-style output function. More...
 
double lagrangian_eulerian_translation_factor (const Vector< double > &s)
 
void traction (const double &time, const Vector< double > &s, Vector< double > &traction)
 
void pressure (const double &time, const Vector< double > &s, double &pressure)
 
void contribution_to_total_porous_flux (double &skeleton_flux_contrib, double &seepage_flux_contrib)
 
- Public Member Functions inherited from oomph::FaceElement
 FaceElement ()
 Constructor: Initialise all appropriate member data. More...
 
virtual ~FaceElement ()
 Empty virtual destructor. More...
 
 FaceElement (const FaceElement &)=delete
 Broken copy constructor. More...
 
const unsignedboundary_number_in_bulk_mesh () const
 Broken assignment operator. More...
 
void set_boundary_number_in_bulk_mesh (const unsigned &b)
 Set function for the boundary number in bulk mesh. More...
 
double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double J_eulerian (const Vector< double > &s) const
 
double J_eulerian_at_knot (const unsigned &ipt) const
 
void check_J_eulerian_at_knots (bool &passed) const
 
double interpolated_x (const Vector< double > &s, const unsigned &i) const
 
double interpolated_x (const unsigned &t, const Vector< double > &s, const unsigned &i) const
 
void interpolated_x (const Vector< double > &s, Vector< double > &x) const
 
void interpolated_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) const
 
double interpolated_dxdt (const Vector< double > &s, const unsigned &i, const unsigned &t)
 
void interpolated_dxdt (const Vector< double > &s, const unsigned &t, Vector< double > &dxdt)
 
intnormal_sign ()
 
int normal_sign () const
 
intface_index ()
 
int face_index () const
 
const Vector< double > * tangent_direction_pt () const
 Public access function for the tangent direction pointer. More...
 
void set_tangent_direction (Vector< double > *tangent_direction_pt)
 Set the tangent direction vector. More...
 
void turn_on_warning_for_discontinuous_tangent ()
 
void turn_off_warning_for_discontinuous_tangent ()
 
void continuous_tangent_and_outer_unit_normal (const Vector< double > &s, Vector< Vector< double >> &tang_vec, Vector< double > &unit_normal) const
 
void continuous_tangent_and_outer_unit_normal (const unsigned &ipt, Vector< Vector< double >> &tang_vec, Vector< double > &unit_normal) const
 
void outer_unit_normal (const Vector< double > &s, Vector< double > &unit_normal) const
 Compute outer unit normal at the specified local coordinate. More...
 
void outer_unit_normal (const unsigned &ipt, Vector< double > &unit_normal) const
 Compute outer unit normal at ipt-th integration point. More...
 
FiniteElement *& bulk_element_pt ()
 Pointer to higher-dimensional "bulk" element. More...
 
FiniteElementbulk_element_pt () const
 Pointer to higher-dimensional "bulk" element (const version) More...
 
CoordinateMappingFctPtface_to_bulk_coordinate_fct_pt ()
 
CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt () const
 
BulkCoordinateDerivativesFctPtbulk_coordinate_derivatives_fct_pt ()
 
BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt () const
 
Vector< doublelocal_coordinate_in_bulk (const Vector< double > &s) const
 
void get_local_coordinate_in_bulk (const Vector< double > &s, Vector< double > &s_bulk) const
 
void get_ds_bulk_ds_face (const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction) const
 
unsignedbulk_position_type (const unsigned &i)
 
const unsignedbulk_position_type (const unsigned &i) const
 
void bulk_node_number_resize (const unsigned &i)
 Resize the storage for the bulk node numbers. More...
 
unsignedbulk_node_number (const unsigned &n)
 
const unsignedbulk_node_number (const unsigned &n) const
 
void bulk_position_type_resize (const unsigned &i)
 Resize the storage for bulk_position_type to i entries. More...
 
unsignednbulk_value (const unsigned &n)
 
unsigned nbulk_value (const unsigned &n) const
 
void nbulk_value_resize (const unsigned &i)
 
void resize_nodes (Vector< unsigned > &nadditional_data_values)
 
void output_zeta (std::ostream &outfile, const unsigned &nplot)
 Output boundary coordinate zeta. More...
 
- 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)
 
virtual void set_macro_elem_pt (MacroElement *macro_elem_pt)
 
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
 
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_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void describe_nodal_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void assign_all_generic_local_eqn_numbers (const bool &store_local_dof_pt)
 
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
 
virtual Nodeconstruct_node (const unsigned &n)
 
virtual Nodeconstruct_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
virtual Nodeconstruct_boundary_node (const unsigned &n)
 
virtual Nodeconstruct_boundary_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
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...
 
unsigned ngeom_data () const
 
Datageom_data_pt (const unsigned &j)
 
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 void identify_geometric_data (std::set< Data * > &geometric_data_pt)
 
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 (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 void compute_norm (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
 

Public Attributes

void(*&)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &tractiontraction_fct_pt ()
 Reference to the traction function pointer. More...
 
void(*&)(const double &time, const Vector< double > &x, const Vector< double > &n, double &pressurepressure_fct_pt ()
 Reference to the pressure function pointer. More...
 

Protected Member Functions

virtual void get_traction (const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
 
virtual void get_pressure (const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, double &pressure)
 
void fill_in_contribution_to_residuals_axisym_poroelasticity_face (Vector< double > &residuals)
 Helper function that actually calculates the residuals. More...
 
- Protected Member Functions inherited from oomph::FaceElement
void add_additional_values (const Vector< unsigned > &nadditional_values, const unsigned &id)
 
- 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)
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
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

void(* Traction_fct_pt )(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
 
void(* Pressure_fct_pt )(const double &time, const Vector< double > &x, const Vector< double > &n, double &result)
 
- Protected Attributes inherited from oomph::FaceElement
unsigned Boundary_number_in_bulk_mesh
 The boundary number in the bulk mesh to which this element is attached. More...
 
FiniteElementBulk_element_pt
 Pointer to the associated higher-dimensional "bulk" element. More...
 
Vector< unsignedBulk_node_number
 
Vector< unsignedNbulk_value
 
Vector< double > * Tangent_direction_pt
 
- 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
 

Additional Inherited Members

- 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
 
- 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

template<class ELEMENT>
class oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >

A class for elements that allow the imposition of an applied combined traction and pore fluid pressure in the axisym poroelasticity equations. The geometrical information can be read from the FaceGeometry<ELEMENT> class and thus, we can be generic enough without the need to have a separate equations class.

Constructor & Destructor Documentation

◆ AxisymmetricPoroelasticityTractionElement() [1/2]

template<class ELEMENT >
oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::AxisymmetricPoroelasticityTractionElement ( FiniteElement *const &  element_pt,
const int face_index 
)
inline

Constructor, which takes a "bulk" element and the value of the index and its limit

158  : FaceGeometry<ELEMENT>(), FaceElement()
159  {
160 #ifdef PARANOID
161  {
162  // Check that the element is not a refineable 3d element
163  ELEMENT* elem_pt = dynamic_cast<ELEMENT*>(element_pt);
164  // If it's three-d
165  if (elem_pt->dim() == 3)
166  {
167  // Is it refineable
168  RefineableElement* ref_el_pt =
169  dynamic_cast<RefineableElement*>(elem_pt);
170  if (ref_el_pt != 0)
171  {
172  if (this->has_hanging_nodes())
173  {
174  throw OomphLibError("This flux element will not work correctly "
175  "if nodes are hanging\n",
178  }
179  }
180  }
181  }
182 #endif
183 
184  // Attach the geometrical information to the element. N.B. This function
185  // also assigns nbulk_value from the required_nvalue of the bulk element
186  element_pt->build_face_element(face_index, this);
187 
188  // Zero traction
191 
192  // Zero pressure
195  }
void(* Traction_fct_pt)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Definition: axisym_poroelasticity_face_elements.h:101
void(* Pressure_fct_pt)(const double &time, const Vector< double > &x, const Vector< double > &n, double &result)
Definition: axisym_poroelasticity_face_elements.h:110
int & face_index()
Definition: elements.h:4626
FaceElement()
Constructor: Initialise all appropriate member data.
Definition: elements.h:4446
bool has_hanging_nodes() const
Definition: elements.h:2470
void Zero_pressure_fct(const double &time, const Vector< double > &x, const Vector< double > &N, double &load)
Default load function (zero pressure)
Definition: axisym_poroelasticity_face_elements.h:68
void Zero_traction_fct(const double &time, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Default load function (zero traction)
Definition: axisym_poroelasticity_face_elements.h:53
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References oomph::FiniteElement::build_face_element(), oomph::FaceElement::face_index(), oomph::FiniteElement::has_hanging_nodes(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::Pressure_fct_pt, oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::Traction_fct_pt, oomph::AxisymmetricPoroelasticityTractionElementHelper::Zero_pressure_fct(), and oomph::AxisymmetricPoroelasticityTractionElementHelper::Zero_traction_fct().

◆ AxisymmetricPoroelasticityTractionElement() [2/2]

Default constructor.

199 {}

Member Function Documentation

◆ contribution_to_total_porous_flux()

template<class ELEMENT >
void oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::contribution_to_total_porous_flux ( double skeleton_flux_contrib,
double seepage_flux_contrib 
)
inline

Compute contributions to integrated porous flux over boundary: q_skeleton = \int \partial u_displ / \partial t \cdot n ds q_seepage = \int k q \cdot n ds

Calculate the FE representation of q

532  {
533  // Get pointer to bulk element
534  ELEMENT* bulk_el_pt = dynamic_cast<ELEMENT*>(bulk_element_pt());
535 
536  // Get permeability from the bulk poroelasticity element
537  const double permeability = bulk_el_pt->permeability();
538 
539  // Find out how many nodes there are
540  const unsigned n_node = this->nnode();
541 
542  // Set up memeory for the shape functions
543  Shape psi(n_node);
544  DShape dpsids(n_node, 1);
545 
546  // Get the value of Nintpt
547  const unsigned n_intpt = integral_pt()->nweight();
548 
549  // Set the Vector to hold local coordinates
550  Vector<double> s(1);
551  Vector<double> s_bulk(2);
552 
553  // Loop over the integration points
554  skeleton_flux_contrib = 0.0;
555  seepage_flux_contrib = 0.0;
556  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
557  {
558  // Assign values of s
559  s[0] = integral_pt()->knot(ipt, 0);
560 
561  // Get the integral weight
562  double W = this->integral_pt()->weight(ipt);
563 
564  // Call the derivatives of the shape function at the knot point
565  this->dshape_local_at_knot(ipt, psi, dpsids);
566 
567  // Get position and tangent vector
568  Vector<double> interpolated_t1(2, 0.0);
569  Vector<double> interpolated_x(2, 0.0);
570  for (unsigned l = 0; l < n_node; l++)
571  {
572  // Loop over directional components
573  for (unsigned i = 0; i < 2; i++)
574  {
575  interpolated_x[i] += this->nodal_position(l, i) * psi(l);
576  interpolated_t1[i] += this->nodal_position(l, i) * dpsids(l, 0);
577  }
578  }
579 
580  // Calculate the length of the tangent Vector
581  double tlength = interpolated_t1[0] * interpolated_t1[0] +
582  interpolated_t1[1] * interpolated_t1[1];
583 
584  // Set the Jacobian of the line element
585  double J = sqrt(tlength) * interpolated_x[0];
586 
587  // Get the outer unit normal
588  Vector<double> interpolated_normal(2);
589  outer_unit_normal(s, interpolated_normal);
590 
591  // Get coordinate in bulk element
592  s_bulk = this->local_coordinate_in_bulk(s);
593 
595  Vector<double> q(2);
596  bulk_el_pt->interpolated_q(s_bulk, q);
597 
598  // Skeleton velocity from bulk
599  Vector<double> du_dt(2);
600  bulk_el_pt->interpolated_du_dt(s_bulk, du_dt);
601 
602 #ifdef PARANOID
603  Vector<double> x_bulk(2);
604  x_bulk[0] = bulk_el_pt->interpolated_x(s_bulk, 0);
605  x_bulk[1] = bulk_el_pt->interpolated_x(s_bulk, 1);
606  double error = sqrt(
607  (interpolated_x[0] - x_bulk[0]) * (interpolated_x[0] - x_bulk[0]) +
608  (interpolated_x[1] - x_bulk[1]) * (interpolated_x[1] - x_bulk[1]));
609  double tol = 1.0e-10;
610  if (error > tol)
611  {
612  std::stringstream junk;
613  junk << "Gap between bulk and face element coordinate\n"
614  << "is suspiciously large: " << error
615  << "\nBulk at: " << x_bulk[0] << " " << x_bulk[1] << "\n"
616  << "Face at: " << interpolated_x[0] << " " << interpolated_x[1]
617  << "\n";
618  OomphLibWarning(
620  }
621 #endif
622 
623  // Get net flux through boundary
624  double q_flux = 0.0;
625  double dudt_flux = 0.0;
626  for (unsigned i = 0; i < 2; i++)
627  {
628  q_flux += permeability * q[i] * interpolated_normal[i];
629  dudt_flux += du_dt[i] * interpolated_normal[i];
630  }
631 
632  // Add
633  seepage_flux_contrib +=
634  2.0 * MathematicalConstants::Pi * q_flux * W * J;
635  skeleton_flux_contrib +=
636  2.0 * MathematicalConstants::Pi * dudt_flux * W * J;
637  }
638  }
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
Definition: elements.cc:6006
Vector< double > local_coordinate_in_bulk(const Vector< double > &s) const
Definition: elements.cc:6353
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
Definition: elements.h:4735
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Definition: elements.h:4528
virtual void dshape_local_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsids) const
Definition: elements.cc:3239
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
double nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.h:2317
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
RealScalar s
Definition: level1_cplx_impl.h:130
EIGEN_DEVICE_FUNC const Scalar & q
Definition: SpecialFunctionsImpl.h:2019
int error
Definition: calibrate.py:297
const double Pi
50 digits from maple
Definition: oomph_utilities.h:157
@ W
Definition: quadtree.h:63

References oomph::FaceElement::bulk_element_pt(), oomph::FiniteElement::dshape_local_at_knot(), calibrate::error, i, oomph::FiniteElement::integral_pt(), oomph::FaceElement::interpolated_x(), J, oomph::Integral::knot(), oomph::FaceElement::local_coordinate_in_bulk(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_position(), oomph::Integral::nweight(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::FaceElement::outer_unit_normal(), oomph::MathematicalConstants::Pi, Eigen::numext::q, s, sqrt(), oomph::QuadTreeNames::W, and oomph::Integral::weight().

◆ fill_in_contribution_to_jacobian()

template<class ELEMENT >
void oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
inlinevirtual

Fill in contribution from Jacobian.

Reimplemented from oomph::GeneralisedElement.

231  {
232  // Call the residuals (element makes no contribution to Jacobian)
234  }
void fill_in_contribution_to_residuals_axisym_poroelasticity_face(Vector< double > &residuals)
Helper function that actually calculates the residuals.
Definition: axisym_poroelasticity_face_elements.h:704

References oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::fill_in_contribution_to_residuals_axisym_poroelasticity_face().

◆ fill_in_contribution_to_residuals()

template<class ELEMENT >
void oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inlinevirtual

◆ fill_in_contribution_to_residuals_axisym_poroelasticity_face()

template<class ELEMENT >
void oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::fill_in_contribution_to_residuals_axisym_poroelasticity_face ( Vector< double > &  residuals)
protected

Helper function that actually calculates the residuals.

Return the residuals for the AxisymmetricPoroelasticityTractionElement equations

706  {
707  // Find out how many nodes there are
708  unsigned n_node = nnode();
709 
710  // Get continuous time from timestepper of first node
711  double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
712 
713 #ifdef PARANOID
714  // Find out how many positional dofs there are
715  unsigned n_position_type = this->nnodal_position_type();
716  if (n_position_type != 1)
717  {
718  throw OomphLibError("Poroelasticity equations are not yet implemented "
719  "for more than one position type",
722  }
723 #endif
724 
725  // Find out the dimension of the node
726  unsigned n_dim = this->nodal_dimension();
727 
728 
729  // Get bulk element
730  ELEMENT* bulk_el_pt = dynamic_cast<ELEMENT*>(bulk_element_pt());
731 
732  unsigned n_q_basis = bulk_el_pt->nq_basis();
733  unsigned n_q_basis_edge = bulk_el_pt->nq_basis_edge();
734 
735  // Integer to hold the local equation number
736  int local_eqn = 0;
737 
738  // Set up memory for the shape functions
739  // Note that in this case, the number of lagrangian coordinates is always
740  // equal to the dimension of the nodes
741  Shape psi(n_node);
742  DShape dpsids(n_node, n_dim - 1);
743  Shape q_basis(n_q_basis, n_dim);
744 
745  // Set the value of n_intpt
746  unsigned n_intpt = integral_pt()->nweight();
747 
748  // Storage for the local coordinates
749  Vector<double> s_face(n_dim - 1, 0.0), s_bulk(n_dim, 0.0);
750 
751  // Loop over the integration points
752  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
753  {
754  // Get the integral weight
755  double w = integral_pt()->weight(ipt);
756 
757  // Only need to call the local derivatives
758  dshape_local_at_knot(ipt, psi, dpsids);
759 
760  // Assign values of s in FaceElement and local coordinates in bulk element
761  for (unsigned i = 0; i < n_dim - 1; i++)
762  {
763  s_face[i] = integral_pt()->knot(ipt, i);
764  }
765  s_bulk = local_coordinate_in_bulk(s_face);
766 
767  // Get bulk element
768  ELEMENT* bulk_el_pt = dynamic_cast<ELEMENT*>(bulk_element_pt());
769 
770  // Get the q basis at bulk local coordinate s_bulk,
771  // corresponding to face local
772  // coordinate s_face
773  bulk_el_pt->get_q_basis(s_bulk, q_basis);
774 
775  // Calculate the Eulerian and Lagrangian coordinates
776  Vector<double> interpolated_x(n_dim, 0.0);
777 
778  // Also calculate the surface Vectors (derivatives wrt local coordinates)
779  DenseMatrix<double> interpolated_A(n_dim - 1, n_dim, 0.0);
780 
781  // Calculate displacements and derivatives
782  for (unsigned l = 0; l < n_node; l++)
783  {
784  // Loop over directions
785  for (unsigned i = 0; i < n_dim; i++)
786  {
787  // Calculate the Eulerian coords
788  const double x_local = nodal_position(l, i);
789  interpolated_x[i] += x_local * psi(l);
790 
791  // Loop over LOCAL derivative directions, to calculate the tangent(s)
792  for (unsigned j = 0; j < n_dim - 1; j++)
793  {
794  interpolated_A(j, i) += x_local * dpsids(l, j);
795  }
796  }
797  }
798 
799  // Now find the local metric tensor from the tangent vectors
800  DenseMatrix<double> A(n_dim - 1);
801  for (unsigned i = 0; i < n_dim - 1; i++)
802  {
803  for (unsigned j = 0; j < n_dim - 1; j++)
804  {
805  // Initialise surface metric tensor to zero
806  A(i, j) = 0.0;
807 
808  // Take the dot product
809  for (unsigned k = 0; k < n_dim; k++)
810  {
811  A(i, j) += interpolated_A(i, k) * interpolated_A(j, k);
812  }
813  }
814  }
815 
816  // Get the outer unit normal
817  Vector<double> interpolated_normal(n_dim);
818  outer_unit_normal(ipt, interpolated_normal);
819 
820  // Find the determinant of the metric tensor
821  double Adet = 0.0;
822  switch (n_dim)
823  {
824  case 2:
825  Adet = A(0, 0);
826  break;
827  case 3:
828  Adet = A(0, 0) * A(1, 1) - A(0, 1) * A(1, 0);
829  break;
830  default:
831  throw OomphLibError(
832  "Wrong dimension in AxisymmetricPoroelasticityTractionElement",
835  }
836 
837  // Premultiply the weights and the square-root of the determinant of
838  // the metric tensor
839  double W = w * sqrt(Adet);
840 
841  // Now calculate the traction load
842  Vector<double> traction(n_dim);
843  get_traction(time, ipt, interpolated_x, interpolated_normal, traction);
844 
845  // Now calculate the load
846  double pressure;
847  get_pressure(time, ipt, interpolated_x, interpolated_normal, pressure);
848 
849  // Loop over the test functions, nodes of the element
850  for (unsigned l = 0; l < n_node; l++)
851  {
852  // Loop over the displacement components
853  for (unsigned i = 0; i < n_dim; i++)
854  {
855  local_eqn = this->nodal_local_eqn(
856  l, bulk_el_pt->u_index_axisym_poroelasticity(i));
857  /*IF it's not a boundary condition*/
858  if (local_eqn >= 0)
859  {
860  // Add the traction loading terms to the residuals
861  residuals[local_eqn] -=
862  traction[i] * psi(l) * interpolated_x[0] * W;
863  } // End of if not boundary condition
864  }
865  } // End of loop over shape functions
866 
867  // Loop over the q edge test functions only (internal basis functions
868  // have zero normal component on the boundary)
869  for (unsigned l = 0; l < n_q_basis_edge; l++)
870  {
871  local_eqn = this->nodal_local_eqn(1, bulk_el_pt->q_edge_index(l));
872 
873  /*IF it's not a boundary condition*/
874  if (local_eqn >= 0)
875  {
876  // Loop over the displacement components
877  for (unsigned i = 0; i < n_dim; i++)
878  {
879  // Add the loading terms to the residuals
880  residuals[local_eqn] += pressure * q_basis(l, i) *
881  interpolated_normal[i] * interpolated_x[0] *
882  W;
883  }
884  } // End of if not boundary condition
885  } // End of loop over shape functions
886  } // End of loop over integration points
887  }
RowVector3d w
Definition: Matrix_resize_int.cpp:3
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
virtual void get_traction(const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Definition: axisym_poroelasticity_face_elements.h:121
virtual void get_pressure(const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, double &pressure)
Definition: axisym_poroelasticity_face_elements.h:135
void pressure(const double &time, const Vector< double > &s, double &pressure)
Definition: axisym_poroelasticity_face_elements.h:677
void traction(const double &time, const Vector< double > &s, Vector< double > &traction)
Definition: axisym_poroelasticity_face_elements.h:651
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
unsigned nnodal_position_type() const
Definition: elements.h:2463
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
Definition: elements.h:2484
Time *const & time_pt() const
Access function for the pointer to time (const version)
Definition: timesteppers.h:572
double & time()
Return the current value of the continuous time.
Definition: timesteppers.h:123
char char char int int * k
Definition: level2_impl.h:374
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References GlobalParameters::get_pressure(), i, j, k, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, sqrt(), w, and oomph::QuadTreeNames::W.

Referenced by oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::fill_in_contribution_to_jacobian(), and oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::fill_in_contribution_to_residuals().

◆ get_pressure()

template<class ELEMENT >
virtual void oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::get_pressure ( const double time,
const unsigned intpt,
const Vector< double > &  x,
const Vector< double > &  n,
double pressure 
)
inlineprotectedvirtual

Get the pressure value: Pass number of integration point (dummy), Eulerrian coordinate and normal vector and return the pressure (not all of the input arguments will be required for all specific load functions but the list should cover all cases). This function is virtual so it can be overloaded for FSI.

Reimplemented in oomph::FSILinearisedAxisymPoroelasticTractionElement< POROELASTICITY_BULK_ELEMENT, NAVIER_STOKES_BULK_ELEMENT >.

140  {
141  Pressure_fct_pt(time, x, n, pressure);
142  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
list x
Definition: plotDoE.py:28

References n, oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::pressure(), oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::Pressure_fct_pt, and plotDoE::x.

Referenced by oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::output().

◆ get_traction()

template<class ELEMENT >
virtual void oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::get_traction ( const double time,
const unsigned intpt,
const Vector< double > &  x,
const Vector< double > &  n,
Vector< double > &  traction 
)
inlineprotectedvirtual

Get the traction vector: Pass number of integration point (dummy), Eulerrian coordinate and normal vector and return the pressure (not all of the input arguments will be required for all specific load functions but the list should cover all cases). This function is virtual so it can be overloaded for FSI.

Reimplemented in oomph::FSILinearisedAxisymPoroelasticTractionElement< POROELASTICITY_BULK_ELEMENT, NAVIER_STOKES_BULK_ELEMENT >.

126  {
127  Traction_fct_pt(time, x, n, traction);
128  }

References n, oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::traction(), oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::Traction_fct_pt, and plotDoE::x.

Referenced by oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::output().

◆ lagrangian_eulerian_translation_factor()

template<class ELEMENT >
double oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::lagrangian_eulerian_translation_factor ( const Vector< double > &  s)
inline

Ratio of lengths of line elements (or annular surface areas) in the undeformed and deformed configuration (needed to translate normal flux correctly between small displacement formulation used here and the possibly large displacement NSt formulation used in FSI problems. Maths is as follows: Surface area: A = 2 \pi \int r \sqrt((dr/ds)^2+(dz/ds)^2) ds = 2 \pi \int J ds This function returns the ratio of J in the undeformed configuration (used here) to that in the deformed configuration where (r,z) = (r,z)_undef + (u_r,u_z).

Calculate the FE representation of u – the skeleton displacement

448  {
449  // Get continuous time from timestepper of first node
450  unsigned n_dim = this->nodal_dimension();
451 
452  // Find out how many nodes there are
453  unsigned n_node = nnode();
454 
455  Vector<double> x(n_dim);
456  Vector<double> s_bulk(n_dim);
457  Vector<double> disp(n_dim);
458  Shape psi(n_node);
459  DShape dpsids(n_node, n_dim - 1);
460 
461  // Call the derivatives of the shape function at the knot point
462  this->dshape_local(s, psi, dpsids);
463 
464  // Get pointer to bulk element
465  ELEMENT* bulk_pt = dynamic_cast<ELEMENT*>(bulk_element_pt());
466  s_bulk = local_coordinate_in_bulk(s);
467 
468  // Get Eulerian coordinates
469  this->interpolated_x(s, x);
470 
471  // Outer unit normal
472  Vector<double> unit_normal(n_dim);
473  outer_unit_normal(s, unit_normal);
474 
476  bulk_pt->interpolated_u(s_bulk, disp);
477 
478  // Get deformed and undeformed tangent vectors
479  Vector<double> interpolated_t1(2, 0.0);
480  Vector<double> interpolated_T1(2, 0.0);
481  for (unsigned l = 0; l < n_node; l++)
482  {
483  // Loop over directional components
484  for (unsigned i = 0; i < 2; i++)
485  {
486  // Index at which the nodal value is stored
487  unsigned u_nodal_index = bulk_pt->u_index_axisym_poroelasticity(i);
488 
489  interpolated_t1[i] += this->nodal_position(l, i) * dpsids(l, 0);
490  interpolated_T1[i] +=
491  (this->nodal_position(l, i) + nodal_value(l, u_nodal_index)) *
492  dpsids(l, 0);
493  }
494  }
495 
496  // Set the Jacobian of the undeformed line element
497  double J_undef = sqrt(interpolated_t1[0] * interpolated_t1[0] +
498  interpolated_t1[1] * interpolated_t1[1]) *
499  x[0];
500 
501 
502  // Set the Jacobian of the deformed line element
503  double J_def = sqrt(interpolated_T1[0] * interpolated_T1[0] +
504  interpolated_T1[1] * interpolated_T1[1]) *
505  (x[0] + disp[0]);
506 
507  double return_val = 1.0;
508  if (J_def != 0.0) return_val = J_undef / J_def;
509  return return_val;
510  }
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Definition: elements.h:1981

References oomph::FaceElement::bulk_element_pt(), oomph::FiniteElement::dshape_local(), i, oomph::FaceElement::interpolated_x(), oomph::FaceElement::local_coordinate_in_bulk(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_dimension(), oomph::FiniteElement::nodal_position(), oomph::FiniteElement::nodal_value(), oomph::FaceElement::outer_unit_normal(), s, sqrt(), and plotDoE::x.

Referenced by oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::output().

◆ output() [1/4]

template<class ELEMENT >
void oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::output ( FILE *  file_pt)
inlinevirtual

C_style output function.

Reimplemented from oomph::FiniteElement.

426  {
428  }
void output(std::ostream &outfile, const unsigned &nplot)
Overload output function.
Definition: overloaded_element_body.h:490

References output().

◆ output() [2/4]

template<class ELEMENT >
void oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::output ( FILE *  file_pt,
const unsigned n_plot 
)
inlinevirtual

C-style output function.

Reimplemented from oomph::FiniteElement.

432  {
433  FaceGeometry<ELEMENT>::output(file_pt, n_plot);
434  }

References output().

◆ output() [3/4]

template<class ELEMENT >
void oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::output ( std::ostream &  outfile)
inlinevirtual

Output function.

Reimplemented from oomph::FiniteElement.

250  {
251  unsigned n_plot = 5;
252  output(outfile, n_plot);
253  }
void output(std::ostream &outfile)
Output function.
Definition: axisym_poroelasticity_face_elements.h:249

◆ output() [4/4]

template<class ELEMENT >
void oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::output ( std::ostream &  outfile,
const unsigned n_plot 
)
inlinevirtual

Output function.

Calculate the FE representation of u – the skeleton displacement

Skeleton velocity

Reimplemented from oomph::FiniteElement.

258  {
259  // Get continuous time from timestepper of first node
260  double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
261  unsigned n_dim = this->nodal_dimension();
262 
263  // Find out how many nodes there are
264  unsigned n_node = nnode();
265 
266  Vector<double> x(n_dim);
267  Vector<double> s(n_dim - 1);
268  Vector<double> s_bulk(n_dim);
269  Vector<double> disp(n_dim);
270  Shape psi(n_node);
271  DShape dpsids(n_node, n_dim - 1);
272 
273  // Tecplot header info
274  outfile << this->tecplot_zone_string(n_plot);
275 
276  // Loop over plot points
277  unsigned num_plot_points = this->nplot_points(n_plot);
278  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
279  {
280  // Get local coordinates of plot point
281  this->get_s_plot(iplot, n_plot, s);
282 
283  // Call the derivatives of the shape function at the knot point
284  this->dshape_local(s, psi, dpsids);
285 
286  // Get pointer to bulk element
287  ELEMENT* bulk_pt = dynamic_cast<ELEMENT*>(bulk_element_pt());
288  s_bulk = local_coordinate_in_bulk(s);
289 
290 
291  // Get Eulerian coordinates
292  this->interpolated_x(s, x);
293 
294  // Outer unit normal
295  Vector<double> unit_normal(n_dim);
296  outer_unit_normal(s, unit_normal);
297 
299  bulk_pt->interpolated_u(s_bulk, disp);
300 
302  Vector<double> du_dt(2);
303  bulk_pt->interpolated_du_dt(s_bulk, du_dt);
304 
305  // Porous seepage flux
306  Vector<double> q(2);
307  bulk_pt->interpolated_q(s_bulk, q);
308 
309  // Get permeability from the bulk poroelasticity element
310  const double permeability = bulk_pt->permeability();
311 
312 
313  // Surface area: S = 2 \pi \int r \sqrt((dr/ds)^2+(dz/ds)^2) ds
314  // = 2 \pi \int r \sqrt( 1 + ( (dr/ds)/(dz/ds) )^2 )
315  // dz/ds ds = 2 \pi \int J dz
316  // where J is an objective measure of the length of the line element
317  // (indep of local coordinates) so should be the same from fluid
318  // and solid.
319 
320  // Get deformed and undeformed tangent vectors
321  Vector<double> interpolated_t1(2, 0.0);
322  Vector<double> interpolated_T1(2, 0.0);
323  for (unsigned l = 0; l < n_node; l++)
324  {
325  // Loop over directional components
326  for (unsigned i = 0; i < 2; i++)
327  {
328  // Index at which the nodal value is stored
329  unsigned u_nodal_index = bulk_pt->u_index_axisym_poroelasticity(i);
330 
331  interpolated_t1[i] += this->nodal_position(l, i) * dpsids(l, 0);
332  interpolated_T1[i] +=
333  (this->nodal_position(l, i) + nodal_value(l, u_nodal_index)) *
334  dpsids(l, 0);
335  }
336  }
337 
338  // Set the Jacobian of the undeformed line element
339  double J_undef =
340  sqrt(1.0 + (interpolated_t1[0] * interpolated_t1[0]) /
341  (interpolated_t1[1] * interpolated_t1[1])) *
342  x[0];
343 
344 
345  // Set the Jacobian of the deformed line element
346  double J_def = sqrt(1.0 + (interpolated_T1[0] * interpolated_T1[0]) /
347  (interpolated_T1[1] * interpolated_T1[1])) *
348  (x[0] + disp[0]);
349 
350  // Dummy
351  unsigned ipt = 0;
352 
353  // Now calculate the traction load
354  Vector<double> traction(n_dim);
355  get_traction(time, ipt, x, unit_normal, traction);
356 
357  // Now calculate the load
358  double pressure;
359  get_pressure(time, ipt, x, unit_normal, pressure);
360 
361  // Get correction factor for geometry
362  double lagr_euler_translation_factor =
364 
365  // Output the x,y,..
366  for (unsigned i = 0; i < n_dim; i++)
367  {
368  outfile << x[i] << " ";
369  } // column 1,2
370 
371  // Output displacement
372  for (unsigned i = 0; i < n_dim; i++)
373  {
374  outfile << disp[i] << " "; // column 3,4
375  }
376 
377  // Output imposed traction
378  for (unsigned i = 0; i < n_dim; i++)
379  {
380  outfile << traction[i] << " "; // column 5,6
381  }
382 
383  // Output imposed pressure
384  outfile << pressure << " "; // column 7
385 
386  // Output seepage flux
387  outfile << permeability * q[0] << " " // column 8
388  << permeability * q[1] << " "; // column 9
389 
390  // Output skeleton velocity
391  outfile << du_dt[0] << " " // column 10
392  << du_dt[1] << " "; // column 11
393 
394  // Total veloc
395  outfile << du_dt[0] + permeability * q[0] << " " // column 12
396  << du_dt[1] + permeability * q[1] << " "; // column 13
397 
398  // Outer unit normal
399  outfile << unit_normal[0] << " " // column 14
400  << unit_normal[1] << " "; // column 15
401 
402  // Output FE representation of div u at s_bulk
403  outfile << bulk_pt->interpolated_div_q(s_bulk) << " "; // column 16
404 
405  // Output FE representation of p at s_bulk
406  outfile << bulk_pt->interpolated_p(s_bulk) << " "; // column 17
407 
408  // Output undeformed jacobian
409  outfile << J_undef << " "; // column 18
410 
411  // Output deformed jacobian
412  outfile << J_def << " "; // column 19
413 
414  // Lagranian/Eulerian translation factor
415  outfile << lagr_euler_translation_factor << " "; // column 20
416 
417  outfile << std::endl;
418  }
419 
420  // Write tecplot footer (e.g. FE connectivity lists)
421  this->write_tecplot_zone_footer(outfile, n_plot);
422  }
double lagrangian_eulerian_translation_factor(const Vector< double > &s)
Definition: axisym_poroelasticity_face_elements.h:447
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Definition: elements.h:3174

References oomph::FaceElement::bulk_element_pt(), oomph::FiniteElement::dshape_local(), oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::get_pressure(), oomph::FiniteElement::get_s_plot(), oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::get_traction(), i, oomph::FaceElement::interpolated_x(), oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::lagrangian_eulerian_translation_factor(), oomph::FaceElement::local_coordinate_in_bulk(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_dimension(), oomph::FiniteElement::nodal_position(), oomph::FiniteElement::nodal_value(), oomph::FiniteElement::node_pt(), oomph::FiniteElement::nplot_points(), oomph::FaceElement::outer_unit_normal(), oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::pressure(), Eigen::numext::q, s, sqrt(), oomph::FiniteElement::tecplot_zone_string(), oomph::Time::time(), oomph::TimeStepper::time_pt(), oomph::Data::time_stepper_pt(), oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::traction(), oomph::FiniteElement::write_tecplot_zone_footer(), and plotDoE::x.

◆ pressure()

template<class ELEMENT >
void oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::pressure ( const double time,
const Vector< double > &  s,
double pressure 
)

Compute pressure value at specified local coordinate Should only be used for post-processing; ignores dependence on integration point!

679  {
680  unsigned n_dim = this->nodal_dimension();
681 
682  // Position vector
683  Vector<double> x(n_dim);
684  interpolated_x(s, x);
685 
686  // Outer unit normal
687  Vector<double> unit_normal(n_dim);
688  outer_unit_normal(s, unit_normal);
689 
690  // Dummy
691  unsigned ipt = 0;
692 
693  // Pressure value
694  get_pressure(time, ipt, x, unit_normal, pressure);
695  }

References GlobalParameters::get_pressure(), s, and plotDoE::x.

Referenced by oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::get_pressure(), and oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::output().

◆ traction()

template<class ELEMENT >
void oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::traction ( const double time,
const Vector< double > &  s,
Vector< double > &  traction 
)

Compute traction vector at specified local coordinate Should only be used for post-processing; ignores dependence on integration point!

//////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// Compute traction vector at specified local coordinate Should only be used for post-processing; ignores dependence on integration point!

653  {
654  unsigned n_dim = this->nodal_dimension();
655 
656  // Position vector
657  Vector<double> x(n_dim);
658  interpolated_x(s, x);
659 
660  // Outer unit normal
661  Vector<double> unit_normal(n_dim);
662  outer_unit_normal(s, unit_normal);
663 
664  // Dummy
665  unsigned ipt = 0;
666 
667  // Pressure value
668  get_traction(time, ipt, x, unit_normal, traction);
669  }

References s, and plotDoE::x.

Referenced by oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::get_traction(), and oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::output().

◆ zeta_nodal()

template<class ELEMENT >
double oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::zeta_nodal ( const unsigned n,
const unsigned k,
const unsigned i 
) const
inlinevirtual

Specify the value of nodal zeta from the face geometry The "global" intrinsic coordinate of the element when viewed as part of a geometric object should be given by the FaceElement representation, by default (needed to break indeterminacy if bulk element is SolidElement)

Reimplemented from oomph::FiniteElement.

244  {
245  return FaceElement::zeta_nodal(n, k, i);
246  }
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: elements.h:4497

References i, k, n, and oomph::FaceElement::zeta_nodal().

Member Data Documentation

◆ pressure_fct_pt

template<class ELEMENT >
void(*&)(const double& time, const Vector<double>& x, const Vector<double>& n, double& pressure) oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::pressure_fct_pt()
inline

Reference to the pressure function pointer.

216  {
217  return Pressure_fct_pt;
218  }

Referenced by AxiPoroProblem< ELEMENT, TIMESTEPPER >::create_pressure_elements().

◆ Pressure_fct_pt

template<class ELEMENT >
void(* oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::Pressure_fct_pt) (const double &time, const Vector< double > &x, const Vector< double > &n, double &result)
protected

Pointer to an imposed pressure function. Arguments: Eulerian coordinate; outer unit normal; applied pressure. (Not all of the input arguments will be required for all specific load functions but the list should cover all cases)

Referenced by oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::AxisymmetricPoroelasticityTractionElement(), and oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::get_pressure().

◆ traction_fct_pt

template<class ELEMENT >
void(*&)(const double& time, const Vector<double>& x, const Vector<double>& n, Vector<double>& traction) oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::traction_fct_pt()
inline

Reference to the traction function pointer.

206  {
207  return Traction_fct_pt;
208  }

Referenced by AxiPoroProblem< ELEMENT, TIMESTEPPER >::create_pressure_elements().

◆ Traction_fct_pt

template<class ELEMENT >
void(* oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::Traction_fct_pt) (const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
protected

Pointer to an imposed traction function. Arguments: Eulerian coordinate; outer unit normal; applied traction. (Not all of the input arguments will be required for all specific load functions but the list should cover all cases)

Referenced by oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::AxisymmetricPoroelasticityTractionElement(), and oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::get_traction().


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