oomph::TimeHarmonicLinearElasticityTractionElement< ELEMENT > Class Template Reference

#include <time_harmonic_linear_elasticity_traction_elements.h>

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

Public Member Functions

 TimeHarmonicLinearElasticityTractionElement (FiniteElement *const &element_pt, const int &face_index)
 
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 &nplot)
 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...
 
void traction (const Vector< double > &s, Vector< std::complex< double >> &traction)
 
- 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 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 Vector< double > &x, const Vector< double > &n, Vector< std::complex< double >> &tractiontraction_fct_pt ()
 Reference to the traction function pointer. More...
 

Protected Member Functions

virtual void get_traction (const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, Vector< std::complex< double >> &traction)
 
void fill_in_contribution_to_residuals_time_harmonic_linear_elasticity_traction (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)
 
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

Vector< std::complex< unsigned > > U_index_time_harmonic_linear_elasticity_traction
 Index at which the i-th displacement component is stored. More...
 
void(* Traction_fct_pt )(const Vector< double > &x, const Vector< double > &n, Vector< std::complex< 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::TimeHarmonicLinearElasticityTractionElement< ELEMENT >

A class for elements that allow the imposition of an applied traction in the equations of time-harmonic linear elasticity. The geometrical information can be read from the FaceGeometry<ELEMENT> class and and thus, we can be generic enough without the need to have a separate equations class.

Constructor & Destructor Documentation

◆ TimeHarmonicLinearElasticityTractionElement()

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

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

120  : FaceGeometry<ELEMENT>(), FaceElement()
121  {
122  // Attach the geometrical information to the element. N.B. This function
123  // also assigns nbulk_value from the required_nvalue of the bulk element
124  element_pt->build_face_element(face_index, this);
125 
126  // Find the dimension of the problem
127  unsigned n_dim = element_pt->nodal_dimension();
128 
129  // Find the index at which the displacement unknowns are stored
130  ELEMENT* cast_element_pt = dynamic_cast<ELEMENT*>(element_pt);
132  for (unsigned i = 0; i < n_dim; i++)
133  {
135  cast_element_pt->u_index_time_harmonic_linear_elasticity(i);
136  }
137 
138  // Zero traction
141 
142 #ifdef PARANOID
143  {
144  // Check that the element is not a refineable 3d element
145  ELEMENT* elem_pt = dynamic_cast<ELEMENT*>(element_pt);
146  // If it's three-d
147  if (elem_pt->dim() == 3)
148  {
149  // Is it refineable
150  RefineableElement* ref_el_pt =
151  dynamic_cast<RefineableElement*>(elem_pt);
152  if (ref_el_pt != 0)
153  {
154  if (this->has_hanging_nodes())
155  {
156  throw OomphLibError("This flux element will not work correctly "
157  "if nodes are hanging\n",
160  }
161  }
162  }
163  }
164 #endif
165  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
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(* Traction_fct_pt)(const Vector< double > &x, const Vector< double > &n, Vector< std::complex< double >> &result)
Definition: time_harmonic_linear_elasticity_traction_elements.h:88
Vector< std::complex< unsigned > > U_index_time_harmonic_linear_elasticity_traction
Index at which the i-th displacement component is stored.
Definition: time_harmonic_linear_elasticity_traction_elements.h:81
void Zero_traction_fct(const Vector< double > &x, const Vector< double > &N, Vector< std::complex< double >> &load)
Default load function (zero traction)
Definition: time_harmonic_linear_elasticity_traction_elements.h:52
#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(), i, oomph::FiniteElement::nodal_dimension(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::TimeHarmonicLinearElasticityTractionElement< ELEMENT >::Traction_fct_pt, oomph::TimeHarmonicLinearElasticityTractionElement< ELEMENT >::U_index_time_harmonic_linear_elasticity_traction, and oomph::TimeHarmonicLinearElasticityTractionElementHelper::Zero_traction_fct().

Member Function Documentation

◆ fill_in_contribution_to_jacobian()

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

Fill in contribution from Jacobian.

Reimplemented from oomph::FiniteElement.

188  {
189  // Call the residuals
191  residuals);
192  }
void fill_in_contribution_to_residuals_time_harmonic_linear_elasticity_traction(Vector< double > &residuals)
Helper function that actually calculates the residuals.
Definition: time_harmonic_linear_elasticity_traction_elements.h:317

References oomph::TimeHarmonicLinearElasticityTractionElement< ELEMENT >::fill_in_contribution_to_residuals_time_harmonic_linear_elasticity_traction().

◆ fill_in_contribution_to_residuals()

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

◆ fill_in_contribution_to_residuals_time_harmonic_linear_elasticity_traction()

template<class ELEMENT >
void oomph::TimeHarmonicLinearElasticityTractionElement< ELEMENT >::fill_in_contribution_to_residuals_time_harmonic_linear_elasticity_traction ( Vector< double > &  residuals)
protected

Helper function that actually calculates the residuals.

Return the residuals for the TimeHarmonicLinearElasticityTractionElement equations

319  {
320  // Find out how many nodes there are
321  unsigned n_node = nnode();
322 
323 #ifdef PARANOID
324  // Find out how many positional dofs there are
325  unsigned n_position_type = this->nnodal_position_type();
326  if (n_position_type != 1)
327  {
328  throw OomphLibError("TimeHarmonicLinearElasticity is not yet implemented "
329  "for more than one position type",
332  }
333 #endif
334 
335  // Find out the dimension of the node
336  const unsigned n_dim = this->nodal_dimension();
337 
338  // Cache the nodal indices at which the displacement components are stored
339  std::vector<std::complex<unsigned>> u_nodal_index(n_dim);
340  for (unsigned i = 0; i < n_dim; i++)
341  {
342  // u_nodal_index[i].real() =
343  // this->U_index_time_harmonic_linear_elasticity_traction[i].real();
344  //
345  // u_nodal_index[i].imag() =
346  // this->U_index_time_harmonic_linear_elasticity_traction[i].imag();
347 
348  u_nodal_index[i] =
350  }
351 
352  // Integer to hold the local equation number
353  int local_eqn = 0;
354 
355  // Set up memory for the shape functions
356  // Note that in this case, the number of lagrangian coordinates is always
357  // equal to the dimension of the nodes
358  Shape psi(n_node);
359  DShape dpsids(n_node, n_dim - 1);
360 
361  // Set the value of n_intpt
362  unsigned n_intpt = integral_pt()->nweight();
363 
364  // Loop over the integration points
365  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
366  {
367  // Get the integral weight
368  double w = integral_pt()->weight(ipt);
369 
370  // Only need to call the local derivatives
371  dshape_local_at_knot(ipt, psi, dpsids);
372 
373  // Calculate the Eulerian and Lagrangian coordinates
374  Vector<double> interpolated_x(n_dim, 0.0);
375 
376  // Also calculate the surface Vectors (derivatives wrt local coordinates)
377  DenseMatrix<double> interpolated_A(n_dim - 1, n_dim, 0.0);
378 
379  // Calculate displacements and derivatives
380  for (unsigned l = 0; l < n_node; l++)
381  {
382  // Loop over directions
383  for (unsigned i = 0; i < n_dim; i++)
384  {
385  // Calculate the Eulerian coords
386  const double x_local = nodal_position(l, i);
387  interpolated_x[i] += x_local * psi(l);
388 
389  // Loop over LOCAL derivative directions, to calculate the tangent(s)
390  for (unsigned j = 0; j < n_dim - 1; j++)
391  {
392  interpolated_A(j, i) += x_local * dpsids(l, j);
393  }
394  }
395  }
396 
397  // Now find the local metric tensor from the tangent Vectors
398  DenseMatrix<double> A(n_dim - 1);
399  for (unsigned i = 0; i < n_dim - 1; i++)
400  {
401  for (unsigned j = 0; j < n_dim - 1; j++)
402  {
403  // Initialise surface metric tensor to zero
404  A(i, j) = 0.0;
405 
406  // Take the dot product
407  for (unsigned k = 0; k < n_dim; k++)
408  {
409  A(i, j) += interpolated_A(i, k) * interpolated_A(j, k);
410  }
411  }
412  }
413 
414  // Get the outer unit normal
415  Vector<double> interpolated_normal(n_dim);
416  outer_unit_normal(ipt, interpolated_normal);
417 
418  // Find the determinant of the metric tensor
419  double Adet = 0.0;
420  switch (n_dim)
421  {
422  case 2:
423  Adet = A(0, 0);
424  break;
425  case 3:
426  Adet = A(0, 0) * A(1, 1) - A(0, 1) * A(1, 0);
427  break;
428  default:
429  throw OomphLibError(
430  "Wrong dimension in TimeHarmonicLinearElasticityTractionElement",
431  "TimeHarmonicLinearElasticityTractionElement::fill_in_contribution_"
432  "to_residuals()",
434  }
435 
436  // Premultiply the weights and the square-root of the determinant of
437  // the metric tensor
438  double W = w * sqrt(Adet);
439 
440  // Now calculate the load
441  Vector<std::complex<double>> traction(n_dim);
442  get_traction(ipt, interpolated_x, interpolated_normal, traction);
443 
444  // Loop over the test functions, nodes of the element
445  for (unsigned l = 0; l < n_node; l++)
446  {
447  // Loop over the displacement components
448  for (unsigned i = 0; i < n_dim; i++)
449  {
450  // Real eqn
451  local_eqn = this->nodal_local_eqn(l, u_nodal_index[i].real());
452  /*IF it's not a boundary condition*/
453  if (local_eqn >= 0)
454  {
455  // Add the loading terms to the residuals
456  residuals[local_eqn] -= traction[i].real() * psi(l) * W;
457  }
458 
459  // Imag eqn
460  local_eqn = this->nodal_local_eqn(l, u_nodal_index[i].imag());
461  /*IF it's not a boundary condition*/
462  if (local_eqn >= 0)
463  {
464  // Add the loading terms to the residuals
465  residuals[local_eqn] -= traction[i].imag() * psi(l) * W;
466  }
467  }
468  } // End of loop over shape functions
469  } // End of loop over integration points
470  }
AnnoyingScalar imag(const AnnoyingScalar &)
Definition: AnnoyingScalar.h:132
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
RowVector3d w
Definition: Matrix_resize_int.cpp:3
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
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
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Definition: elements.h:4528
unsigned nnodal_position_type() const
Definition: elements.h:2463
virtual void dshape_local_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsids) const
Definition: elements.cc:3239
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
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
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
Definition: elements.h:2484
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.
virtual void get_traction(const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, Vector< std::complex< double >> &traction)
Definition: time_harmonic_linear_elasticity_traction_elements.h:98
void traction(const Vector< double > &s, Vector< std::complex< double >> &traction)
Definition: time_harmonic_linear_elasticity_traction_elements.h:290
float real
Definition: datatypes.h:10
char char char int int * k
Definition: level2_impl.h:374
@ W
Definition: quadtree.h:63
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

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

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

◆ get_traction()

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

Get the traction vector: Pass number of integration point (dummy), Eulerian coordinate and normal vector and return the load vector (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.

102  {
104  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
list x
Definition: plotDoE.py:28

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

◆ output() [1/4]

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

C_style output function.

Reimplemented from oomph::FiniteElement.

262  {
263  FiniteElement::output(file_pt);
264  }
virtual void output(std::ostream &outfile)
Definition: elements.h:3050

References oomph::FiniteElement::output().

◆ output() [2/4]

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

C-style output function.

Reimplemented from oomph::FiniteElement.

268  {
269  FiniteElement::output(file_pt, n_plot);
270  }

References oomph::FiniteElement::output().

◆ output() [3/4]

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

Output function.

Reimplemented from oomph::FiniteElement.

208  {
209  unsigned nplot = 5;
210  output(outfile, nplot);
211  }
void output(std::ostream &outfile)
Output function.
Definition: time_harmonic_linear_elasticity_traction_elements.h:207

◆ output() [4/4]

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

Output function.

Reimplemented from oomph::FiniteElement.

215  {
216  unsigned ndim = dim();
217  Vector<double> s(ndim);
218  Vector<double> x(ndim + 1);
219  Vector<std::complex<double>> traction(ndim + 1);
220 
221  // Tecplot header info
222  outfile << this->tecplot_zone_string(nplot);
223 
224  // Loop over plot points
225  unsigned num_plot_points = this->nplot_points(nplot);
226  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
227  {
228  // Get local coordinates of plot point
229  this->get_s_plot(iplot, nplot, s);
230 
231  // Get Eulerian coordinates and displacements
232  this->interpolated_x(s, x);
233  this->traction(s, traction);
234 
235  // Output the x,y,..
236  for (unsigned i = 0; i < ndim + 1; i++)
237  {
238  outfile << x[i] << " ";
239  }
240 
241  // Output u,v,..
242  for (unsigned i = 0; i < ndim + 1; i++)
243  {
244  outfile << traction[i].real() << " ";
245  }
246 
247  // Output u,v,..
248  for (unsigned i = 0; i < ndim + 1; i++)
249  {
250  outfile << traction[i].imag() << " ";
251  }
252 
253  outfile << std::endl;
254  }
255 
256  // Write tecplot footer (e.g. FE connectivity lists)
257  this->write_tecplot_zone_footer(outfile, nplot);
258  }
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
unsigned dim() const
Definition: elements.h:2611
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
unsigned ndim() const
Access function to # of Eulerian coordinates.
Definition: geom_objects.h:177
RealScalar s
Definition: level1_cplx_impl.h:130

References oomph::FiniteElement::dim(), oomph::FiniteElement::get_s_plot(), i, oomph::FaceElement::interpolated_x(), oomph::GeomObject::ndim(), oomph::FiniteElement::nplot_points(), s, oomph::FiniteElement::tecplot_zone_string(), oomph::TimeHarmonicLinearElasticityTractionElement< ELEMENT >::traction(), oomph::FiniteElement::write_tecplot_zone_footer(), and plotDoE::x.

◆ traction()

template<class ELEMENT >
void oomph::TimeHarmonicLinearElasticityTractionElement< ELEMENT >::traction ( const Vector< double > &  s,
Vector< std::complex< 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!

292  {
293  unsigned n_dim = this->nodal_dimension();
294 
295  // Position vector
296  Vector<double> x(n_dim);
297  interpolated_x(s, x);
298 
299  // Outer unit normal
300  Vector<double> unit_normal(n_dim);
301  outer_unit_normal(s, unit_normal);
302 
303  // Dummy
304  unsigned ipt = 0;
305 
306  // Traction vector
307  get_traction(ipt, x, unit_normal, traction);
308  }

References s, and plotDoE::x.

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

◆ zeta_nodal()

template<class ELEMENT >
double oomph::TimeHarmonicLinearElasticityTractionElement< 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::FaceElement.

202  {
203  return FaceElement::zeta_nodal(n, k, i);
204  }
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

◆ traction_fct_pt

◆ Traction_fct_pt

template<class ELEMENT >
void(* oomph::TimeHarmonicLinearElasticityTractionElement< ELEMENT >::Traction_fct_pt) (const Vector< double > &x, const Vector< double > &n, Vector< std::complex< 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::TimeHarmonicLinearElasticityTractionElement< ELEMENT >::get_traction(), and oomph::TimeHarmonicLinearElasticityTractionElement< ELEMENT >::TimeHarmonicLinearElasticityTractionElement().

◆ U_index_time_harmonic_linear_elasticity_traction

template<class ELEMENT >
Vector<std::complex<unsigned> > oomph::TimeHarmonicLinearElasticityTractionElement< ELEMENT >::U_index_time_harmonic_linear_elasticity_traction
protected

Index at which the i-th displacement component is stored.

Referenced by oomph::TimeHarmonicLinearElasticityTractionElement< ELEMENT >::TimeHarmonicLinearElasticityTractionElement().


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