oomph::ElementWithMovingNodes Class Reference

#include <element_with_moving_nodes.h>

+ Inheritance diagram for oomph::ElementWithMovingNodes:

Public Types

enum  { Shape_derivs_by_chain_rule , Shape_derivs_by_direct_fd , Shape_derivs_by_fastest_method }
 Public enumeration to choose method for computing shape derivatives. More...
 
- Public Types inherited from oomph::FiniteElement
typedef void(* SteadyExactSolutionFctPt) (const Vector< double > &, Vector< double > &)
 
typedef void(* UnsteadyExactSolutionFctPt) (const double &, const Vector< double > &, Vector< double > &)
 

Public Member Functions

 ElementWithMovingNodes ()
 Constructor. More...
 
 ElementWithMovingNodes (const ElementWithMovingNodes &)=delete
 Broken copy constructor. More...
 
void operator= (const ElementWithMovingNodes &)=delete
 Broken assignment operator. More...
 
void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual ~ElementWithMovingNodes ()
 Virtual destructor (clean up and allocated memory) More...
 
unsigned ngeom_dof () const
 Number of geometric dofs. More...
 
int geometric_data_local_eqn (const unsigned &n, const unsigned &i)
 
void assemble_set_of_all_geometric_data (std::set< Data * > &unique_geom_data_pt)
 Return a set of all geometric data associated with the element. More...
 
void identify_geometric_data (std::set< Data * > &geometric_data_pt)
 
bool are_dresidual_dnodal_coordinates_always_evaluated_by_fd () const
 Return whether shape derivatives are evaluated by fd. More...
 
void enable_always_evaluate_dresidual_dnodal_coordinates_by_fd ()
 
void disable_always_evaluate_dresidual_dnodal_coordinates_by_fd ()
 
void evaluate_shape_derivs_by_direct_fd ()
 Evaluate shape derivatives by direct finite differencing. More...
 
void evaluate_shape_derivs_by_chain_rule (const bool &i_know_what_i_am_doing=false)
 
void evaluate_shape_derivs_by_fastest_method (const bool &i_know_what_i_am_doing=false)
 
intmethod_for_shape_derivs ()
 Access to method (enumerated flag) for determination of shape derivs. More...
 
void enable_bypass_fill_in_jacobian_from_geometric_data ()
 Bypass the call to fill_in_jacobian_from_geometric_data. More...
 
void disable_bypass_fill_in_jacobian_from_geometric_data ()
 Do not bypass the call to fill_in_jacobian_from_geometric_data. More...
 
bool is_fill_in_jacobian_from_geometric_data_bypassed () const
 
unsigned ngeom_data () const
 
- Public Member Functions inherited from oomph::FiniteElement
void set_dimension (const unsigned &dim)
 
void set_nodal_dimension (const unsigned &nodal_dim)
 
void set_nnodal_position_type (const unsigned &nposition_type)
 Set the number of types required to interpolate the coordinate. More...
 
void set_n_node (const unsigned &n)
 
int nodal_local_eqn (const unsigned &n, const unsigned &i) const
 
double dJ_eulerian_at_knot (const unsigned &ipt, Shape &psi, DenseMatrix< double > &djacobian_dX) const
 
 FiniteElement ()
 Constructor. More...
 
virtual ~FiniteElement ()
 
 FiniteElement (const FiniteElement &)=delete
 Broken copy constructor. More...
 
virtual bool local_coord_is_valid (const Vector< double > &s)
 Broken assignment operator. More...
 
virtual void move_local_coord_back_into_element (Vector< double > &s) const
 
void get_centre_of_gravity_and_max_radius_in_terms_of_zeta (Vector< double > &cog, double &max_radius) const
 
virtual void local_coordinate_of_node (const unsigned &j, Vector< double > &s) const
 
virtual void local_fraction_of_node (const unsigned &j, Vector< double > &s_fraction)
 
virtual double local_one_d_fraction_of_node (const unsigned &n1d, const unsigned &i)
 
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
 
virtual double J_eulerian (const Vector< double > &s) const
 
virtual double J_eulerian_at_knot (const unsigned &ipt) const
 
void check_J_eulerian_at_knots (bool &passed) const
 
void check_jacobian (const double &jacobian) const
 
double dshape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsi, DenseMatrix< double > &djacobian_dX, RankFourTensor< double > &d_dpsidx_dX) const
 
double d2shape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual double d2shape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual void assign_nodal_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void describe_nodal_local_dofs (std::ostream &out, const std::string &current_string) const
 
Node *& node_pt (const unsigned &n)
 Return a pointer to the local node n. More...
 
Node *const & node_pt (const unsigned &n) const
 Return a pointer to the local node n (const version) More...
 
unsigned nnode () const
 Return the number of nodes. More...
 
virtual unsigned nnode_1d () const
 
double raw_nodal_position (const unsigned &n, const unsigned &i) const
 
double raw_nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position (const unsigned &n, const unsigned &i) const
 
double nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double dnodal_position_dt (const unsigned &n, const unsigned &i) const
 Return the i-th component of nodal velocity: dx/dt at local node n. More...
 
double dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
virtual void get_dresidual_dnodal_coordinates (RankThreeTensor< double > &dresidual_dnodal_coordinates)
 
virtual void disable_ALE ()
 
virtual void enable_ALE ()
 
virtual unsigned required_nvalue (const unsigned &n) const
 
unsigned nnodal_position_type () const
 
bool has_hanging_nodes () const
 
unsigned nodal_dimension () const
 Return the required Eulerian dimension of the nodes in this element. More...
 
virtual unsigned nvertex_node () const
 
virtual Nodevertex_node_pt (const unsigned &j) const
 
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...
 
virtual double interpolated_x (const Vector< double > &s, const unsigned &i) const
 Return FE interpolated coordinate x[i] at local coordinate s. More...
 
virtual double interpolated_x (const unsigned &t, const Vector< double > &s, const unsigned &i) const
 
virtual void interpolated_x (const Vector< double > &s, Vector< double > &x) const
 Return FE interpolated position x[] at local coordinate s as Vector. More...
 
virtual void interpolated_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) const
 
virtual double interpolated_dxdt (const Vector< double > &s, const unsigned &i, const unsigned &t)
 
virtual void interpolated_dxdt (const Vector< double > &s, const unsigned &t, Vector< double > &dxdt)
 
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)
 
virtual double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 
void interpolated_zeta (const Vector< double > &s, Vector< double > &zeta) const
 
void locate_zeta (const Vector< double > &zeta, GeomObject *&geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
 
virtual void node_update ()
 
virtual void identify_field_data_for_interactions (std::set< std::pair< Data *, unsigned >> &paired_field_data)
 
virtual double s_min () const
 Min value of local coordinate. More...
 
virtual double s_max () const
 Max. value of local coordinate. More...
 
double size () const
 
virtual double compute_physical_size () const
 
virtual void point_output_data (const Vector< double > &s, Vector< double > &data)
 
void point_output (std::ostream &outfile, const Vector< double > &s)
 
virtual unsigned nplot_points_paraview (const unsigned &nplot) const
 
virtual unsigned nsub_elements_paraview (const unsigned &nplot) const
 
void output_paraview (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_output_offset_information (std::ofstream &file_out, const unsigned &nplot, unsigned &counter) const
 
virtual void write_paraview_type (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_offsets (std::ofstream &file_out, const unsigned &nplot, unsigned &offset_sum) const
 
virtual unsigned nscalar_paraview () const
 
virtual void scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
 
virtual std::string scalar_name_paraview (const unsigned &i) const
 
virtual void output (std::ostream &outfile)
 
virtual void output (std::ostream &outfile, const unsigned &n_plot)
 
virtual void output (const unsigned &t, std::ostream &outfile, const unsigned &n_plot) const
 
virtual void output (FILE *file_pt)
 
virtual void output (FILE *file_pt, const unsigned &n_plot)
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 Output an exact solution over the element. More...
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
 Output a time-dependent exact solution over the element. More...
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, const SolutionFunctorBase &exact_soln) const
 Output a time-dependent exact solution over the element. More...
 
virtual void get_s_plot (const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
 
virtual std::string tecplot_zone_string (const unsigned &nplot) const
 
virtual void write_tecplot_zone_footer (std::ostream &outfile, const unsigned &nplot) const
 
virtual void write_tecplot_zone_footer (FILE *file_pt, const unsigned &nplot) const
 
virtual unsigned nplot_points (const unsigned &nplot) const
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 Calculate the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Calculate the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_abs_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error)
 
void integrate_fct (FiniteElement::SteadyExactSolutionFctPt integrand_fct_pt, Vector< double > &integral)
 Integrate Vector-valued function over element. More...
 
void integrate_fct (FiniteElement::UnsteadyExactSolutionFctPt integrand_fct_pt, const double &time, Vector< double > &integral)
 Integrate Vector-valued time-dep function over element. More...
 
virtual void build_face_element (const int &face_index, FaceElement *face_element_pt)
 
virtual unsigned self_test ()
 
virtual unsigned get_bulk_node_number (const int &face_index, const unsigned &i) const
 
virtual int face_outer_unit_normal_sign (const int &face_index) const
 Get the sign of the outer unit normal on the face given by face_index. More...
 
virtual unsigned nnode_on_face () const
 
void face_node_number_error_check (const unsigned &i) const
 Range check for face node numbers. More...
 
virtual CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt (const int &face_index) const
 
virtual BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt (const int &face_index) const
 
- 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 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
 

Protected Member Functions

virtual void get_dnodal_coordinates_dgeom_dofs (RankThreeTensor< double > &dnodal_coordinates_dgeom_dofs)
 
void complete_setup_of_dependencies ()
 Construct the vector of (unique) geometric data. More...
 
virtual void assign_all_generic_local_eqn_numbers (const bool &store_local_dof_pt)
 
void fill_in_jacobian_from_geometric_data (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
void fill_in_jacobian_from_geometric_data (DenseMatrix< double > &jacobian)
 
- 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)
 
virtual void fill_in_contribution_to_residuals (Vector< double > &residuals)
 
void fill_in_jacobian_from_internal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_internal_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
virtual void update_before_internal_fd ()
 
virtual void reset_after_internal_fd ()
 
virtual void update_in_internal_fd (const unsigned &i)
 
virtual void reset_in_internal_fd (const unsigned &i)
 
virtual void update_before_external_fd ()
 
virtual void reset_after_external_fd ()
 
virtual void update_in_external_fd (const unsigned &i)
 
virtual void reset_in_external_fd (const unsigned &i)
 
virtual void fill_in_contribution_to_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_dresiduals_dparameter (double *const &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< Data * > Geom_data_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
 

Private Attributes

int ** Geometric_data_local_eqn
 
unsigned Ngeom_dof
 
bool Bypass_fill_in_jacobian_from_geometric_data
 
bool Evaluate_dresidual_dnodal_coordinates_by_fd
 
int Method_for_shape_derivs
 

Additional Inherited Members

- 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

/////////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////////// A policy class that serves to establish the common interfaces for elements that contain moving nodes. This class provides storage for the geometric data that affect the update of all the nodes of the element, i.e. USUALLY all data that are using during a call to the Element's node_update() function. In some cases (e.g. FluidInterfaceEdge elements), node_update() is overloaded to perform an update of the bulk element, in which case the additional bulk geometric data become external data of the element and the function GeneralisedElement::update_in_external_fd(i) is overloaded to also perform the bulk node update. The storage is populated during the assignment of the equation numbers via the complete_setup_of_dependencies() function and then local equations numbers are assigned to these data, accessible via geometric_data_local_eqn(n,i). Finally, a function is provided that calculates the terms in the jacobian matrix by due to these geometric data by finite differences.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum

Public enumeration to choose method for computing shape derivatives.

Enumerator
Shape_derivs_by_chain_rule 
Shape_derivs_by_direct_fd 
Shape_derivs_by_fastest_method 
69  {
73  };
@ Shape_derivs_by_chain_rule
Definition: element_with_moving_nodes.h:70
@ Shape_derivs_by_fastest_method
Definition: element_with_moving_nodes.h:72
@ Shape_derivs_by_direct_fd
Definition: element_with_moving_nodes.h:71

Constructor & Destructor Documentation

◆ ElementWithMovingNodes() [1/2]

oomph::ElementWithMovingNodes::ElementWithMovingNodes ( )
inline

Constructor.

81  // hierher: Anything other than the fd-based method is currently broken;
82  // at least for refineable elements -- this all needs to be checked
83  // VERY carefully again (see instructions in commit log). Until this
84  // has all been fixed/re-validated, we've frozen the default assignment
85  // (by breaking the functions that change the setting). Ultimately,
86  // we should use:
87  // Method_for_shape_derivs(Shape_derivs_by_fastest_method)
88  {
89  }
bool Evaluate_dresidual_dnodal_coordinates_by_fd
Definition: element_with_moving_nodes.h:367
bool Bypass_fill_in_jacobian_from_geometric_data
Definition: element_with_moving_nodes.h:361
int Method_for_shape_derivs
Definition: element_with_moving_nodes.h:371
int ** Geometric_data_local_eqn
Definition: element_with_moving_nodes.h:353

◆ ElementWithMovingNodes() [2/2]

oomph::ElementWithMovingNodes::ElementWithMovingNodes ( const ElementWithMovingNodes )
delete

Broken copy constructor.

◆ ~ElementWithMovingNodes()

virtual oomph::ElementWithMovingNodes::~ElementWithMovingNodes ( )
inlinevirtual

Virtual destructor (clean up and allocated memory)

110  {
112  {
113  delete[] Geometric_data_local_eqn[0];
114  delete[] Geometric_data_local_eqn;
115  }
116  }

References Geometric_data_local_eqn.

Member Function Documentation

◆ are_dresidual_dnodal_coordinates_always_evaluated_by_fd()

bool oomph::ElementWithMovingNodes::are_dresidual_dnodal_coordinates_always_evaluated_by_fd ( ) const
inline

Return whether shape derivatives are evaluated by fd.

187  {
189  }

References Evaluate_dresidual_dnodal_coordinates_by_fd.

◆ assemble_set_of_all_geometric_data()

void oomph::ElementWithMovingNodes::assemble_set_of_all_geometric_data ( std::set< Data * > &  unique_geom_data_pt)

Return a set of all geometric data associated with the element.

////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////// Return a set of all geometric data associated with the element's node update function

45  {
46  // First clear the set (just in case)
47  unique_geom_data_pt.clear();
48 
49  // Get number of nodes
50  const unsigned n_node = this->nnode();
51 
52  // Loop over all nodes
53  for (unsigned n = 0; n < n_node; n++)
54  {
55  // Cache pointer to the Node
56  Node* const nod_pt = this->node_pt(n);
57 
58  // Is the node hanging
59  const bool node_is_hanging = nod_pt->is_hanging();
60 
61  // Default number of master nodes
62  unsigned nmaster = 1;
63 
64  // Default: Node isn't hanging so it's its own master node
65  Node* master_node_pt = nod_pt;
66 
67  // Cache the hanging point
68  HangInfo* hang_info_pt = 0;
69 
70  // Find the number of master nodes if the node is hanging
71  if (node_is_hanging)
72  {
73  hang_info_pt = nod_pt->hanging_pt();
74  nmaster = hang_info_pt->nmaster();
75  }
76 
77  // Loop over all master nodes
78  for (unsigned imaster = 0; imaster < nmaster; imaster++)
79  {
80  // Get the master node
81  if (node_is_hanging)
82  {
83  master_node_pt = hang_info_pt->master_node_pt(imaster);
84  }
85 
86 
87  // Find the number of data
88  const unsigned n_geom_data = master_node_pt->ngeom_data();
89  // If there are geometric data add them to the set
90  if (n_geom_data > 0)
91  {
92  // Get vector of geometric data involved in the geometric
93  // change of this node
94  Data** node_geom_data_pt = master_node_pt->all_geom_data_pt();
95 
96  for (unsigned i = 0; i < n_geom_data; i++)
97  {
98  unique_geom_data_pt.insert(node_geom_data_pt[i]);
99  }
100  }
101 
102  // Find the number of geometric objects
103  unsigned n_geom_obj = master_node_pt->ngeom_object();
104 
105  // If there are geometric objects, add them to the set
106  if (n_geom_obj > 0)
107  {
108  // Get vector of geometric objects involved in the default
109  // update function for this (master) node.
110  // Vector is constructed by copy operation.
111  GeomObject** geom_object_pt = master_node_pt->all_geom_object_pt();
112 
113  // Loop over the geometric objects
114  for (unsigned i = 0; i < n_geom_obj; i++)
115  {
116  // Get the next geometric object
117  GeomObject* geom_obj_pt = geom_object_pt[i];
118 
119  // Number of items of geometric data that the geometric
120  // object depends on
121  unsigned n_geom_data = geom_obj_pt->ngeom_data();
122 
123  // Loop over geometric data and add to set (use set to ensure
124  // that each one is only counted once)
125  for (unsigned idata = 0; idata < n_geom_data; idata++)
126  {
127  unique_geom_data_pt.insert(geom_obj_pt->geom_data_pt(idata));
128  }
129  }
130  } // End of add geom object loop
131  }
132  }
133  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
GeomObject()
Default constructor.
Definition: geom_objects.h:104

References oomph::Node::all_geom_data_pt(), oomph::Node::all_geom_object_pt(), oomph::GeomObject::geom_data_pt(), oomph::Node::hanging_pt(), i, oomph::Node::is_hanging(), oomph::HangInfo::master_node_pt(), n, oomph::GeomObject::ngeom_data(), oomph::Node::ngeom_data(), oomph::Node::ngeom_object(), oomph::HangInfo::nmaster(), oomph::FiniteElement::nnode(), and oomph::FiniteElement::node_pt().

Referenced by complete_setup_of_dependencies(), oomph::FixedVolumeSpineLineMarangoniFluidInterfaceElement< ELEMENT >::make_bounding_element(), and oomph::SpineUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::make_bounding_element().

◆ assign_all_generic_local_eqn_numbers()

void oomph::ElementWithMovingNodes::assign_all_generic_local_eqn_numbers ( const bool store_local_dof_pt)
protectedvirtual

Assign local equation numbers for the geometric Data in the element If the boolean argument is true then the degrees of freedom are stored in Dof_pt

Assign local equation numbers for the geometric data associated with the element.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::ElementWithSpecificMovingNodes< ELEMENT, NODE_TYPE >, oomph::ElementWithSpecificMovingNodes< ELEMENT, SpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, PerturbedSpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, AlgebraicNode >, and oomph::ElementWithSpecificMovingNodes< ELEMENT, MacroElementNodeUpdateNode >.

199  {
200  // Get local number of dofs so far
201  unsigned local_eqn_number = this->ndof();
202 
203  // Set the number of data
204  const unsigned n_geom_data = ngeom_data();
205 
206  // Reset number of geometric dofs
207  Ngeom_dof = 0;
208 
209  // If we have any geometric data
210  if (n_geom_data > 0)
211  {
212  // Work out total number of values involved
213  // Initialise from the first object
214  unsigned n_total_values = Geom_data_pt[0]->nvalue();
215 
216  // Add the values from the other data
217  for (unsigned i = 1; i < n_geom_data; i++)
218  {
219  n_total_values += Geom_data_pt[i]->nvalue();
220  }
221 
222  // If allocated delete the old storage
224  {
225  delete[] Geometric_data_local_eqn[0];
226  delete[] Geometric_data_local_eqn;
227  }
228 
229  // If there are no values, we are done, null out the storage and
230  // return
231  if (n_total_values == 0)
232  {
234  return;
235  }
236 
237  // Resize the storage for the geometric data local equation numbers
238  // Firstly allocate the row pointers
239  Geometric_data_local_eqn = new int*[n_geom_data];
240 
241  // Now allocate storage for all the equation numbers
242  Geometric_data_local_eqn[0] = new int[n_total_values];
243 
244  // Initially all local equations are unclassified
245  for (unsigned i = 0; i < n_total_values; i++)
246  {
248  }
249 
250  // Loop over the remaining rows and set their pointers
251  for (unsigned i = 1; i < n_geom_data; ++i)
252  {
253  // Initially set the pointer to the i-th row to the pointer
254  // to the i-1th row
256 
257  // Now increase the row pointer by the number of values
258  // stored at the i-1th geometric data
259  Geometric_data_local_eqn[i] += Geom_data_pt[i - 1]->nvalue();
260  }
261 
262  // A local queue to store the global equation numbers
263  std::deque<unsigned long> global_eqn_number_queue;
264 
265  // Loop over the node update data
266  for (unsigned i = 0; i < n_geom_data; i++)
267  {
268  // Pointer to geometric Data
269  Data* data_pt = Geom_data_pt[i];
270 
271  // Loop over values at this Data item
272  unsigned n_value = data_pt->nvalue();
273  for (unsigned j = 0; j < n_value; j++)
274  {
275  // Get global equation number
276  long eqn_number = data_pt->eqn_number(j);
277 
278  // If equation number positive
279  if (eqn_number >= 0)
280  {
281  // Add the global equation number to our queue
282  global_eqn_number_queue.push_back(eqn_number);
283  // Add pointer to the dof to the queue if required
284  if (store_local_dof_pt)
285  {
286  GeneralisedElement::Dof_pt_deque.push_back(data_pt->value_pt(j));
287  }
288 
289  // Add to local value
292 
293  // Bump up number of geometric dofs
294  Ngeom_dof++;
295  }
296  else
297  {
298  // Set the local scheme to be pinned
300  }
301  }
302  }
303 
304  // Now add our global equations numbers to the internal element storage
305  this->add_global_eqn_numbers(global_eqn_number_queue,
307  // Clear the memory used in the deque
308  if (store_local_dof_pt)
309  {
310  std::deque<double*>().swap(GeneralisedElement::Dof_pt_deque);
311  }
312  }
313  }
static long Is_pinned
Static "Magic number" to indicate pinned values.
Definition: nodes.h:183
static long Is_unclassified
Definition: nodes.h:192
Vector< Data * > Geom_data_pt
Definition: element_with_moving_nodes.h:340
unsigned Ngeom_dof
Definition: element_with_moving_nodes.h:357
unsigned ngeom_data() const
Definition: element_with_moving_nodes.h:345
unsigned ndof() const
Return the number of equations/dofs in the element.
Definition: elements.h:835
unsigned long eqn_number(const unsigned &ieqn_local) const
Definition: elements.h:704
static std::deque< double * > Dof_pt_deque
Definition: elements.h:231
int local_eqn_number(const unsigned long &ieqn_global) const
Definition: elements.h:726
void add_global_eqn_numbers(std::deque< unsigned long > const &global_eqn_numbers, std::deque< double * > const &global_dof_pt)
Definition: elements.cc:156
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::GeneralisedElement::add_global_eqn_numbers(), oomph::GeneralisedElement::Dof_pt_deque, oomph::Data::eqn_number(), oomph::GeneralisedElement::eqn_number(), Geom_data_pt, Geometric_data_local_eqn, i, oomph::Data::Is_pinned, oomph::Data::Is_unclassified, j, oomph::GeneralisedElement::local_eqn_number(), oomph::GeneralisedElement::ndof(), ngeom_data(), Ngeom_dof, oomph::Data::nvalue(), and oomph::Data::value_pt().

Referenced by oomph::ElementWithSpecificMovingNodes< ELEMENT, NODE_TYPE >::assign_all_generic_local_eqn_numbers().

◆ complete_setup_of_dependencies()

void oomph::ElementWithMovingNodes::complete_setup_of_dependencies ( )
protectedvirtual

Construct the vector of (unique) geometric data.

Reimplemented from oomph::GeneralisedElement.

Reimplemented in oomph::SpineElement< ELEMENT >, oomph::SpineElement< PointElement >, oomph::SpineElement< FaceGeometry< ELEMENT > >, oomph::ElementWithSpecificMovingNodes< ELEMENT, NODE_TYPE >, oomph::ElementWithSpecificMovingNodes< ELEMENT, SpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, PerturbedSpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, AlgebraicNode >, and oomph::ElementWithSpecificMovingNodes< ELEMENT, MacroElementNodeUpdateNode >.

139  {
140  // This set will hold the pointers to all the unique (geometric) Data that
141  // affects the shape of this element
142  std::set<Data*> unique_geom_data_pt;
143  // Assemble that data
144  this->assemble_set_of_all_geometric_data(unique_geom_data_pt);
145 
146  // Resize storage for the pointers to the Data items that are
147  // involved in the element's node update operation.
148  Geom_data_pt.clear();
149 
150  // Loop over all the unique remaining Data items involved in the
151  // node update operations
152  typedef std::set<Data*>::iterator IT;
153  for (IT it = unique_geom_data_pt.begin(); it != unique_geom_data_pt.end();
154  it++)
155  {
156  Geom_data_pt.push_back(*it);
157  }
158  }
void assemble_set_of_all_geometric_data(std::set< Data * > &unique_geom_data_pt)
Return a set of all geometric data associated with the element.
Definition: element_with_moving_nodes.cc:43

References assemble_set_of_all_geometric_data(), and Geom_data_pt.

Referenced by oomph::ElementWithSpecificMovingNodes< ELEMENT, NODE_TYPE >::complete_setup_of_dependencies().

◆ describe_local_dofs()

void oomph::ElementWithMovingNodes::describe_local_dofs ( std::ostream &  out,
const std::string &  current_string 
) const
virtual

Function to describe the local dofs of the element. The ostream specifies the output stream to which the description is written; the string stores the currently assembled output that is ultimately written to the output stream by Data::describe_dofs(...); it is typically built up incrementally as we descend through the call hierarchy of this function when called from Problem::describe_dofs(...)

Function to describe the local dofs of the element[s]. The ostream specifies the output stream to which the description is written; the string stores the currently assembled output that is ultimately written to the output stream by Data::describe_dofs(...); it is typically built up incrementally as we descend through the call hierarchy of this function when called from Problem::describe_dofs(...)

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::ElementWithSpecificMovingNodes< ELEMENT, NODE_TYPE >, oomph::ElementWithSpecificMovingNodes< ELEMENT, SpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, PerturbedSpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, AlgebraicNode >, and oomph::ElementWithSpecificMovingNodes< ELEMENT, MacroElementNodeUpdateNode >.

172  {
173  // Call the standard finite element classification.
174  FiniteElement::describe_local_dofs(out, current_string);
175 
176  // Set the number of data
177  const unsigned n_geom_data = ngeom_data();
178 
179  // Loop over the node update data
180  for (unsigned i = 0; i < n_geom_data; i++)
181  {
182  // Pointer to geometric Data
183  Data* data_pt = Geom_data_pt[i];
184 
185  std::stringstream conversion;
186  conversion << " of Geometric Data " << i << current_string;
187  std::string in(conversion.str());
188  data_pt->describe_dofs(out, in);
189  }
190  }
virtual void describe_local_dofs(std::ostream &out, const std::string &current_string) const
Definition: elements.cc:1709
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
std::ofstream out("Result.txt")

References oomph::Data::describe_dofs(), oomph::FiniteElement::describe_local_dofs(), Geom_data_pt, i, ngeom_data(), out(), and oomph::Global_string_for_annotation::string().

Referenced by oomph::ElementWithSpecificMovingNodes< ELEMENT, NODE_TYPE >::describe_local_dofs().

◆ disable_always_evaluate_dresidual_dnodal_coordinates_by_fd()

void oomph::ElementWithMovingNodes::disable_always_evaluate_dresidual_dnodal_coordinates_by_fd ( )
inline

Insist that shape derivatives are always evaluated using the overloaded version of this function that may have been implemented in a derived class. (The default behaviour will still be finite differences unless the function has actually been overloaded

205  {
207  }

References Evaluate_dresidual_dnodal_coordinates_by_fd.

Referenced by FSIRingProblem::FSIRingProblem().

◆ disable_bypass_fill_in_jacobian_from_geometric_data()

void oomph::ElementWithMovingNodes::disable_bypass_fill_in_jacobian_from_geometric_data ( )
inline

Do not bypass the call to fill_in_jacobian_from_geometric_data.

281  {
283  }

References Bypass_fill_in_jacobian_from_geometric_data.

◆ enable_always_evaluate_dresidual_dnodal_coordinates_by_fd()

void oomph::ElementWithMovingNodes::enable_always_evaluate_dresidual_dnodal_coordinates_by_fd ( )
inline

Insist that shape derivatives are always evaluated by fd (using FiniteElement::get_dresidual_dnodal_coordinates())

195  {
197  }

References Evaluate_dresidual_dnodal_coordinates_by_fd.

Referenced by oomph::RefineableQElement< 3 >::build(), oomph::RefineableQElement< 1 >::build(), oomph::RefineableQElement< 2 >::build(), and FSIRingProblem::FSIRingProblem().

◆ enable_bypass_fill_in_jacobian_from_geometric_data()

void oomph::ElementWithMovingNodes::enable_bypass_fill_in_jacobian_from_geometric_data ( )
inline

◆ evaluate_shape_derivs_by_chain_rule()

void oomph::ElementWithMovingNodes::evaluate_shape_derivs_by_chain_rule ( const bool i_know_what_i_am_doing = false)
inline

Evaluate shape derivatives by chain rule. Currently disabled by default because it's broken; can re-enable use by setting optional boolean to true.

222  {
223  if (!i_know_what_i_am_doing)
224  {
225  std::ostringstream error_message;
226  error_message
227  << "Evaluation of shape derivatives by chain rule is currently \n"
228  << "disabled because it's broken, at least for refineable \n"
229  << "elements. This all needs to be checked again very carefully\n"
230  << "following the instructions in the commit log\n"
231  << "If you know what you're doing and want to force this "
232  "methodology\n"
233  << "call this function with the optional boolean set to true.\n";
234  throw OomphLibError(
235  error_message.str(),
236  "ElementWithMovingNodes::evaluate_shape_derivs_by_chain_rule()",
238  }
240  }
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61

References Method_for_shape_derivs, OOMPH_EXCEPTION_LOCATION, and Shape_derivs_by_chain_rule.

Referenced by FSIRingProblem::FSIRingProblem().

◆ evaluate_shape_derivs_by_direct_fd()

void oomph::ElementWithMovingNodes::evaluate_shape_derivs_by_direct_fd ( )
inline

◆ evaluate_shape_derivs_by_fastest_method()

void oomph::ElementWithMovingNodes::evaluate_shape_derivs_by_fastest_method ( const bool i_know_what_i_am_doing = false)
inline

Evaluate shape derivatives by (anticipated) fastest method. Currently disabled by default because it's broken; can re-enable use by setting optional boolean to true.

247  {
248  if (!i_know_what_i_am_doing)
249  {
250  std::ostringstream error_message;
251  error_message
252  << "Evaluation of shape derivatives by fastest method is currently \n"
253  << "disabled because it's broken, at least for refineable \n"
254  << "elements. This all needs to be checked again very carefully\n"
255  << "following the instructions in the commit log\n"
256  << "If you know what you're doing and want to force this "
257  "methodology\n"
258  << "call this function with the optional boolean set to true.\n";
259  throw OomphLibError(
260  error_message.str(),
261  "ElementWithMovingNodes::evaluate_shape_derivs_by_fastest_method()",
263  }
265  }

References Method_for_shape_derivs, OOMPH_EXCEPTION_LOCATION, and Shape_derivs_by_fastest_method.

Referenced by FSIRingProblem::FSIRingProblem().

◆ fill_in_jacobian_from_geometric_data() [1/2]

void oomph::ElementWithMovingNodes::fill_in_jacobian_from_geometric_data ( DenseMatrix< double > &  jacobian)
inlineprotected

Calculate the contributions to the Jacobian matrix from the geometric data. This version computes the residuals vector before calculating the Jacobian terms

321  {
323  {
324  // Allocate storage for the residuals
325  const unsigned n_dof = ndof();
326  Vector<double> residuals(n_dof);
327 
328  // Get the residuals for the entire element
329  get_residuals(residuals);
330 
331  // Call the jacobian calculation
332  fill_in_jacobian_from_geometric_data(residuals, jacobian);
333  }
334  }
void fill_in_jacobian_from_geometric_data(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: element_with_moving_nodes.cc:323
virtual void get_residuals(Vector< double > &residuals)
Definition: elements.h:980

References Bypass_fill_in_jacobian_from_geometric_data, fill_in_jacobian_from_geometric_data(), oomph::GeneralisedElement::get_residuals(), and oomph::GeneralisedElement::ndof().

◆ fill_in_jacobian_from_geometric_data() [2/2]

void oomph::ElementWithMovingNodes::fill_in_jacobian_from_geometric_data ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
protected

Calculate the contributions to the Jacobian matrix from the geometric data. This version assumes that the (full) residuals vector has already been calculated and is passed in as the first argument – needed in case the derivatives are computed by FD.

Calculate the node-update–related entries in the Jacobian. The vector passed in residuals has to contain the nonlinear residuals, evaluated for the current values of the unknowns, in case FDing is used to computed the derivatives.

325  {
327  {
328  // Get number of Data items involved in node update operations
329  const unsigned n_geometric_data = ngeom_data();
330 
331  // If there is nothing to be done, then leave
332  if (n_geometric_data == 0) return;
333 
334  // Number of dofs
335  const unsigned n_dof = this->ndof();
336 
337  // Number of nodes
338  unsigned n_nod = this->nnode();
339 
340  // If there are no dofs, return
341  if (n_nod == 0) return;
342 
343  // Get nodal dimension from first node
344  const unsigned dim_nod = node_pt(0)->ndim();
345 
346  // Number of shape controlling nodes for nonrefineable elements
347  unsigned n_shape_controlling_node = nnode();
348 
349  // Are we dealing with a refineable element?
350  RefineableElement* ref_el_pt = dynamic_cast<RefineableElement*>(this);
351  if (ref_el_pt != 0)
352  {
353  // Adjust number of shape controlling nodes
354  n_shape_controlling_node = ref_el_pt->nshape_controlling_nodes();
355  }
356 
357  // How are we going to evaluate the shape derivs?
358  unsigned method = 0;
360  {
361  method = 0;
362  }
364  {
365  method = 1;
366  }
368  {
369  // Direct FD-ing of residuals w.r.t. geometric dofs is likely to be
370  // faster if there are fewer geometric dofs than total nodal coordinates
371  // (nodes x dim) in element:
372  if (Ngeom_dof < (n_shape_controlling_node * dim_nod))
373  {
374  method = 0;
375  }
376  else
377  {
378  method = 1;
379  }
380  }
381 
382  // Choose method
383  //===============
384  switch (method)
385  {
386  // Direct FD:
387  //-----------
388  case 0:
389 
390  {
391  // Create newres vector
392  Vector<double> newres(n_dof);
393 
394  // Use the default finite difference step
395  const double fd_step = GeneralisedElement::Default_fd_jacobian_step;
396 
397  // Integer storage for the local unknown
398  int local_unknown = 0;
399 
400  // Loop over the Data items that affect the node update operations
401  for (unsigned i = 0; i < n_geometric_data; i++)
402  {
403  // Loop over values
404  unsigned n_value = Geom_data_pt[i]->nvalue();
405  for (unsigned j = 0; j < n_value; j++)
406  {
407  local_unknown = geometric_data_local_eqn(i, j);
408 
409  // If the value is free
410  if (local_unknown >= 0)
411  {
412  // Get a pointer to the geometric data value
413  double* value_pt = Geom_data_pt[i]->value_pt(j);
414 
415  // Save the old value
416  double old_var = *value_pt;
417 
418  // Increment the variable
419  *value_pt += fd_step;
420 
421  // Update the whole element (Bit inefficient)
422  this->node_update();
423 
424  // Calculate the new residuals
425  this->get_residuals(newres);
426 
427  // Now do finite differences
428  for (unsigned m = 0; m < n_dof; m++)
429  {
430  // Stick the entry into the Jacobian matrix
431  jacobian(m, local_unknown) =
432  (newres[m] - residuals[m]) / fd_step;
433  }
434 
435  // Reset the variable
436  *value_pt = old_var;
437 
438  // We're relying on the total node update in the next loop
439  }
440  }
441  }
442 
443  // Node update the element one final time to get things back to
444  // the original state
445  this->node_update();
446  }
447 
448  break;
449 
450  // Chain rule
451  //-----------
452  case 1:
453 
454  {
455  // Get derivatives of residuals w.r.t. all nodal coordinates
456  RankThreeTensor<double> dresidual_dnodal_coordinates(
457  n_dof, dim_nod, n_shape_controlling_node, 0.0);
458 
459  // Use FD-version in base class?
461  {
462  if (ref_el_pt != 0)
463  {
464  ref_el_pt->RefineableElement::get_dresidual_dnodal_coordinates(
465  dresidual_dnodal_coordinates);
466  }
467  else
468  {
470  dresidual_dnodal_coordinates);
471  }
472  }
473  // Otherwise use the overloaded analytical version in derived
474  // class (if it exists -- if it doesn't this just drops through
475  // to the default implementation in FiniteElement).
476  else
477  {
479  dresidual_dnodal_coordinates);
480  }
481 
482  // Get derivatives of nodal coordinates w.r.t. geometric dofs
483  RankThreeTensor<double> dnodal_coordinates_dgeom_dofs(
484  n_dof, dim_nod, n_shape_controlling_node, 0.0);
485 
486  get_dnodal_coordinates_dgeom_dofs(dnodal_coordinates_dgeom_dofs);
487 
488  // Assemble Jacobian via chain rule
489  for (unsigned l = 0; l < n_dof; l++)
490  {
491  // Loop over the Data items that affect the node update operations
492  for (unsigned i_data = 0; i_data < n_geometric_data; i_data++)
493  {
494  // Loop over values
495  unsigned n_value = Geom_data_pt[i_data]->nvalue();
496  for (unsigned j_val = 0; j_val < n_value; j_val++)
497  {
498  int k = geometric_data_local_eqn(i_data, j_val);
499 
500  // If the value is free
501  if (k >= 0)
502  {
503  jacobian(l, k) = 0.0;
504  for (unsigned i = 0; i < dim_nod; i++)
505  {
506  for (unsigned j = 0; j < n_shape_controlling_node; j++)
507  {
508  jacobian(l, k) += dresidual_dnodal_coordinates(l, i, j) *
509  dnodal_coordinates_dgeom_dofs(k, i, j);
510  }
511  }
512  }
513  }
514  }
515  }
516  }
517 
518  break;
519 
520  default:
521 
522  std::ostringstream error_message;
523  error_message << "Never get here: method " << method;
524  throw OomphLibError(error_message.str(),
527  }
528  }
529  }
int geometric_data_local_eqn(const unsigned &n, const unsigned &i)
Definition: element_with_moving_nodes.h:126
virtual void get_dnodal_coordinates_dgeom_dofs(RankThreeTensor< double > &dnodal_coordinates_dgeom_dofs)
Definition: element_with_moving_nodes.cc:538
virtual void node_update()
Definition: elements.cc:5072
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Definition: elements.cc:3744
static double Default_fd_jacobian_step
Definition: elements.h:1198
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
int * m
Definition: level2_cplx_impl.h:294
char char char int int * k
Definition: level2_impl.h:374
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References Bypass_fill_in_jacobian_from_geometric_data, oomph::GeneralisedElement::Default_fd_jacobian_step, Evaluate_dresidual_dnodal_coordinates_by_fd, Geom_data_pt, geometric_data_local_eqn(), get_dnodal_coordinates_dgeom_dofs(), oomph::FiniteElement::get_dresidual_dnodal_coordinates(), oomph::GeneralisedElement::get_residuals(), i, j, k, m, Method_for_shape_derivs, oomph::Node::ndim(), oomph::GeneralisedElement::ndof(), ngeom_data(), Ngeom_dof, oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), oomph::FiniteElement::node_update(), oomph::RefineableElement::nshape_controlling_nodes(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, Shape_derivs_by_chain_rule, Shape_derivs_by_direct_fd, and Shape_derivs_by_fastest_method.

Referenced by oomph::PerturbedSpineLinearisedAxisymmetricFluidInterfaceElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::SpineLineVolumeConstraintBoundingElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::SpineAxisymmetricVolumeConstraintBoundingElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::SpineSurfaceVolumeConstraintBoundingElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::SpineUpdateFluidInterfaceElement< EQUATION_CLASS, DERIVATIVE_CLASS, ELEMENT >::fill_in_contribution_to_jacobian(), oomph::SpinePointFluidInterfaceBoundingElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::SpineLineFluidInterfaceBoundingElement< ELEMENT >::fill_in_contribution_to_jacobian(), fill_in_jacobian_from_geometric_data(), oomph::ElementWithSpecificMovingNodes< ELEMENT, NODE_TYPE >::get_jacobian(), and oomph::ElementWithSpecificMovingNodes< ELEMENT, NODE_TYPE >::get_jacobian_and_mass_matrix().

◆ geometric_data_local_eqn()

int oomph::ElementWithMovingNodes::geometric_data_local_eqn ( const unsigned n,
const unsigned i 
)
inline

Return the local equation number corresponding to the i-th value at the n-th geometric data object.

127  {
128 #ifdef RANGE_CHECKING
129  unsigned n_data = Geom_data_pt.size();
130  if (n >= n_data)
131  {
132  std::ostringstream error_message;
133  error_message << "Range Error: Data number " << n
134  << " is not in the range (0," << n_data - 1 << ")";
135  throw OomphLibError(error_message.str(),
138  }
139  else
140  {
141  unsigned n_value = Geom_data_pt[n]->nvalue();
142  if (i >= n_value)
143  {
144  std::ostringstream error_message;
145  error_message << "Range Error: value " << i << " at data " << n
146  << " is not in the range (0," << n_value - 1 << ")";
147  throw OomphLibError(error_message.str(),
150  }
151  }
152 #endif
153 #ifdef PARANOID
154  // Check that the equations have been allocated
155  if (Geometric_data_local_eqn == 0)
156  {
157  throw OomphLibError(
158  "Geometric data local equation numbers have not been allocated",
161  }
162 #endif
163  return Geometric_data_local_eqn[n][i];
164  }

References Geom_data_pt, Geometric_data_local_eqn, i, n, OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by fill_in_jacobian_from_geometric_data(), get_dnodal_coordinates_dgeom_dofs(), oomph::SpineElement< ELEMENT >::spine_local_eqn(), and oomph::PerturbedSpineElement< ELEMENT >::spine_local_eqn().

◆ get_dnodal_coordinates_dgeom_dofs()

void oomph::ElementWithMovingNodes::get_dnodal_coordinates_dgeom_dofs ( RankThreeTensor< double > &  dnodal_coordinates_dgeom_dofs)
protectedvirtual

Compute derivatives of the nodal coordinates w.r.t. to the geometric dofs. Default implementation by FD can be overwritten for specific elements. dnodal_coordinates_dgeom_dofs(l,i,j) = dX_{ij} / d s_l

540  {
541  // Get number of Data items involved in node update operations
542  const unsigned n_geometric_data = ngeom_data();
543 
544  // If there is nothing to be done, then leave
545  if (n_geometric_data == 0)
546  {
547  return;
548  }
549 
550  // Number of nodes
551  const unsigned n_nod = nnode();
552 
553  // If the element has no nodes (why??!!) return straightaway
554  if (n_nod == 0) return;
555 
556  // Get dimension from first node
557  unsigned dim_nod = node_pt(0)->ndim();
558 
559  // Number of shape controlling nodes for nonrefineable elements
560  unsigned n_shape_controlling_node = n_nod;
561 
562  // Are we dealing with a refineable element?
563  RefineableElement* ref_el_pt = dynamic_cast<RefineableElement*>(this);
564  if (ref_el_pt != 0)
565  {
566  // Adjust number of shape controlling nodes
567  n_shape_controlling_node = ref_el_pt->nshape_controlling_nodes();
568  }
569 
570  // Current and advanced nodal positions
571  DenseMatrix<double> pos(dim_nod, n_shape_controlling_node);
572 
573  // Shape controlling nodes
574  std::map<Node*, unsigned> local_shape_controlling_node_lookup;
575 
576  // Refineable element:
577  if (ref_el_pt != 0)
578  {
579  local_shape_controlling_node_lookup =
580  ref_el_pt->shape_controlling_node_lookup();
581  }
582  // Non-refineable element: the nodes themselves
583  else
584  {
585  unsigned count = 0;
586  for (unsigned j = 0; j < n_nod; j++)
587  {
588  local_shape_controlling_node_lookup[node_pt(j)] = count;
589  count++;
590  }
591  }
592 
593  // Loop over all shape-controlling nodes to backup their original position
594  for (std::map<Node*, unsigned>::iterator it =
595  local_shape_controlling_node_lookup.begin();
596  it != local_shape_controlling_node_lookup.end();
597  it++)
598  {
599  // Get node
600  Node* nod_pt = it->first;
601 
602  // Get its number
603  unsigned node_number = it->second;
604 
605  // Backup
606  for (unsigned i = 0; i < dim_nod; i++)
607  {
608  pos(i, node_number) = nod_pt->position(i);
609  }
610  }
611 
612 
613  // Integer storage for the local unknown
614  int local_unknown = 0;
615 
616  // Use the default finite difference step
617  const double fd_step = GeneralisedElement::Default_fd_jacobian_step;
618 
619  // Loop over the Data items that affect the node update operations
620  for (unsigned i = 0; i < n_geometric_data; i++)
621  {
622  // Loop over values
623  unsigned n_value = Geom_data_pt[i]->nvalue();
624  for (unsigned j = 0; j < n_value; j++)
625  {
626  local_unknown = geometric_data_local_eqn(i, j);
627 
628  // If the value is free
629  if (local_unknown >= 0)
630  {
631  // Get a pointer to the geometric data value
632  double* value_pt = Geom_data_pt[i]->value_pt(j);
633 
634  // Save the old value
635  double old_var = *value_pt;
636 
637  // Increment the variable
638  *value_pt += fd_step;
639 
640  // Update the whole element
641  this->node_update();
642 
643  // Loop over all shape-controlling nodes
644  for (std::map<Node*, unsigned>::iterator it =
645  local_shape_controlling_node_lookup.begin();
646  it != local_shape_controlling_node_lookup.end();
647  it++)
648  {
649  // Get node
650  Node* nod_pt = it->first;
651 
652  // Get its number
653  unsigned node_number = it->second;
654 
655  // Get advanced position and FD
656  for (unsigned ii = 0; ii < dim_nod; ii++)
657  {
658  dnodal_coordinates_dgeom_dofs(local_unknown, ii, node_number) =
659  (nod_pt->position(ii) - pos(ii, node_number)) / fd_step;
660  }
661  }
662 
663  // Reset the variable
664  *value_pt = old_var;
665 
666  // We're relying on the total node update in the next loop
667  }
668  }
669  }
670  // Node update the element one final time to get things back to
671  // the original state
672  this->node_update();
673  }

References oomph::GeneralisedElement::Default_fd_jacobian_step, Geom_data_pt, geometric_data_local_eqn(), i, j, oomph::Node::ndim(), ngeom_data(), oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), oomph::FiniteElement::node_update(), oomph::RefineableElement::nshape_controlling_nodes(), oomph::Node::position(), and oomph::RefineableElement::shape_controlling_node_lookup().

Referenced by fill_in_jacobian_from_geometric_data().

◆ identify_geometric_data()

void oomph::ElementWithMovingNodes::identify_geometric_data ( std::set< Data * > &  geometric_data_pt)
inlinevirtual

Specify Data that affects the geometry of the element by adding the element's geometric Data to the set that's passed in. (This functionality is required in FSI problems; set is used to avoid double counting).

Reimplemented from oomph::FiniteElement.

176  {
177  // Loop over the node update data and add to the set
178  const unsigned n_geom_data = Geom_data_pt.size();
179  for (unsigned n = 0; n < n_geom_data; n++)
180  {
181  geometric_data_pt.insert(Geom_data_pt[n]);
182  }
183  }

References Geom_data_pt, and n.

◆ is_fill_in_jacobian_from_geometric_data_bypassed()

bool oomph::ElementWithMovingNodes::is_fill_in_jacobian_from_geometric_data_bypassed ( ) const
inline

Test whether the call to fill_in_jacobian_from_geometric_data is bypassed

288  {
290  }

References Bypass_fill_in_jacobian_from_geometric_data.

◆ method_for_shape_derivs()

int& oomph::ElementWithMovingNodes::method_for_shape_derivs ( )
inline

Access to method (enumerated flag) for determination of shape derivs.

269  {
271  }

References Method_for_shape_derivs.

Referenced by oomph::RefineableQElement< 3 >::build(), oomph::RefineableQElement< 1 >::build(), and oomph::RefineableQElement< 2 >::build().

◆ ngeom_data()

unsigned oomph::ElementWithMovingNodes::ngeom_data ( ) const
inlinevirtual

◆ ngeom_dof()

unsigned oomph::ElementWithMovingNodes::ngeom_dof ( ) const
inline

Number of geometric dofs.

120  {
121  return Ngeom_dof;
122  }

References Ngeom_dof.

◆ operator=()

void oomph::ElementWithMovingNodes::operator= ( const ElementWithMovingNodes )
delete

Broken assignment operator.

Member Data Documentation

◆ Bypass_fill_in_jacobian_from_geometric_data

bool oomph::ElementWithMovingNodes::Bypass_fill_in_jacobian_from_geometric_data
private

◆ Evaluate_dresidual_dnodal_coordinates_by_fd

bool oomph::ElementWithMovingNodes::Evaluate_dresidual_dnodal_coordinates_by_fd
private

Boolean to decide if shape derivatives are to be evaluated by fd (using FiniteElement::get_dresidual_dnodal_coordinates()) or analytically, using the overloaded version of this function in this class.

Referenced by are_dresidual_dnodal_coordinates_always_evaluated_by_fd(), disable_always_evaluate_dresidual_dnodal_coordinates_by_fd(), enable_always_evaluate_dresidual_dnodal_coordinates_by_fd(), and fill_in_jacobian_from_geometric_data().

◆ Geom_data_pt

Vector<Data*> oomph::ElementWithMovingNodes::Geom_data_pt
protected

Vector that stores pointers to all Data that affect the node update operations, i.e. the variables that can affect the position of the node.

Referenced by assign_all_generic_local_eqn_numbers(), complete_setup_of_dependencies(), describe_local_dofs(), fill_in_jacobian_from_geometric_data(), geometric_data_local_eqn(), get_dnodal_coordinates_dgeom_dofs(), identify_geometric_data(), and ngeom_data().

◆ Geometric_data_local_eqn

int** oomph::ElementWithMovingNodes::Geometric_data_local_eqn
private

Array to hold local eqn number information for the geometric Data variables

Referenced by assign_all_generic_local_eqn_numbers(), geometric_data_local_eqn(), and ~ElementWithMovingNodes().

◆ Method_for_shape_derivs

int oomph::ElementWithMovingNodes::Method_for_shape_derivs
private

Choose method for evaluation of shape derivatives (this takes one of the values in the enumeration)

Referenced by evaluate_shape_derivs_by_chain_rule(), evaluate_shape_derivs_by_direct_fd(), evaluate_shape_derivs_by_fastest_method(), fill_in_jacobian_from_geometric_data(), and method_for_shape_derivs().

◆ Ngeom_dof

unsigned oomph::ElementWithMovingNodes::Ngeom_dof
private

Number of geometric dofs (computed on the fly when equation numbers are set up)

Referenced by assign_all_generic_local_eqn_numbers(), fill_in_jacobian_from_geometric_data(), and ngeom_dof().


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