oomph::QElement< 2, NNODE_1D > Class Template Reference

General QElement class specialised to two spatial dimensions. More...

#include <Qelements.h>

+ Inheritance diagram for oomph::QElement< 2, NNODE_1D >:

Public Member Functions

 QElement ()
 Constructor. More...
 
 QElement (const QElement &)=delete
 Broken copy constructor. More...
 
void shape (const Vector< double > &s, Shape &psi) const
 Broken assignment operator. More...
 
void dshape_local (const Vector< double > &s, Shape &psi, DShape &dpsids) const
 Derivatives of shape functions for specific QElement<2,..> More...
 
void d2shape_local (const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
 
double invert_jacobian_mapping (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
double s_min () const
 Min. value of local coordinate. More...
 
double s_max () const
 Max. value of local coordinate. More...
 
unsigned nvertex_node () const
 Number of vertex nodes in the element. More...
 
Nodevertex_node_pt (const unsigned &j) const
 Pointer to the j-th vertex node in the element. More...
 
void local_coordinate_of_node (const unsigned &j, Vector< double > &s) const
 Get local coordinates of node j in the element; vector sets its own size. More...
 
void local_fraction_of_node (const unsigned &j, Vector< double > &s_fraction)
 Get the local fraction of node j in the element. More...
 
double local_one_d_fraction_of_node (const unsigned &n1d, const unsigned &i)
 
Nodeget_node_at_local_coordinate (const Vector< double > &s) const
 Get the node at the specified local coordinate. More...
 
unsigned nnode_1d () const
 Number of nodes along each element edge. More...
 
unsigned nplot_points_paraview (const unsigned &nplot) const
 
unsigned nsub_elements_paraview (const unsigned &nplot) const
 
void write_paraview_output_offset_information (std::ofstream &file_out, const unsigned &nplot, unsigned &counter) const
 
void write_paraview_type (std::ofstream &file_out, const unsigned &nplot) const
 
void write_paraview_offsets (std::ofstream &file_out, const unsigned &nplot, unsigned &offset_sum) const
 
void output (std::ostream &outfile)
 Output. More...
 
void output (std::ostream &outfile, const unsigned &n_plot)
 Output at n_plot points. More...
 
void output (FILE *file_pt)
 C-style output. More...
 
void output (FILE *file_pt, const unsigned &n_plot)
 C_style output at n_plot points. More...
 
void get_s_plot (const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &use_equally_spaced_interior_sample_points=false) const
 
std::string tecplot_zone_string (const unsigned &nplot) const
 
unsigned nplot_points (const unsigned &nplot) const
 
unsigned get_bulk_node_number (const int &face_index, const unsigned &i) const
 
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...
 
CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt (const int &face_index) const
 
BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt (const int &face_index) const
 
- Public Member Functions inherited from oomph::QuadElementBase
 QuadElementBase ()
 Constructor. Empty. More...
 
- Public Member Functions inherited from oomph::QElementBase
 QElementBase ()
 Constructor: Initialise pointers to macro element reference coords. More...
 
 QElementBase (const QElementBase &)=delete
 Broken copy constructor. More...
 
virtual ~QElementBase ()
 Broken assignment operator. More...
 
bool local_coord_is_valid (const Vector< double > &s)
 Check whether the local coordinate are valid or not. More...
 
void move_local_coord_back_into_element (Vector< double > &s) const
 
virtual void set_macro_elem_pt (MacroElement *macro_elem_pt)
 
doubles_macro_ll (const unsigned &i)
 
doubles_macro_ur (const unsigned &i)
 
double s_macro_ll (const unsigned &i) const
 
double s_macro_ur (const unsigned &i) const
 
void get_x_from_macro_element (const Vector< double > &s, Vector< double > &x) const
 
void get_x_from_macro_element (const unsigned &t, const Vector< double > &s, Vector< double > &x)
 
unsigned nnode_on_face () const
 
ElementGeometry::ElementGeometry element_geometry () const
 It's a Q element! More...
 
- Public Member Functions inherited from oomph::QElementGeometricBase
 QElementGeometricBase ()
 Empty default constructor. More...
 
 QElementGeometricBase (const QElementGeometricBase &)=delete
 Broken copy constructor. 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...
 
void get_centre_of_gravity_and_max_radius_in_terms_of_zeta (Vector< double > &cog, double &max_radius) const
 
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 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_at_knot (const unsigned &ipt, Shape &psi) const
 
virtual void dshape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids) 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_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...
 
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 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
 
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 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 void identify_geometric_data (std::set< Data * > &geometric_data_pt)
 
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)
 
void output_paraview (std::ofstream &file_out, const unsigned &nplot) 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 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 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 ()
 
void face_node_number_error_check (const unsigned &i) const
 Range check for face node numbers. More...
 
- 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
 

Static Private Attributes

static Gauss< 2, NNODE_1D > Default_integration_scheme
 Assign the spatial integration scheme. More...
 

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
 
- 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 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 inherited from oomph::FiniteElement
MacroElementMacro_elem_pt
 Pointer to the element's macro element (NULL by default) More...
 
- Protected Attributes inherited from oomph::GeomObject
unsigned NLagrangian
 Number of Lagrangian (intrinsic) coordinates. More...
 
unsigned Ndim
 Number of Eulerian coordinates. More...
 
TimeStepperGeom_object_time_stepper_pt
 
- Static Protected Attributes inherited from oomph::FiniteElement
static const unsigned Default_Initial_Nvalue = 0
 Default value for the number of values at a node. More...
 
static const double Node_location_tolerance = 1.0e-14
 
static const unsigned N2deriv [] = {0, 1, 3, 6}
 
- Static Protected Attributes inherited from oomph::GeneralisedElement
static DenseMatrix< doubleDummy_matrix
 
static std::deque< double * > Dof_pt_deque
 

Detailed Description

template<unsigned NNODE_1D>
class oomph::QElement< 2, NNODE_1D >

General QElement class specialised to two spatial dimensions.

Constructor & Destructor Documentation

◆ QElement() [1/2]

template<unsigned NNODE_1D>
oomph::QElement< 2, NNODE_1D >::QElement ( )
inline

Constructor.

850  {
851  // There are NNODE_1D*NNODE_1D nodes in this element
852  this->set_n_node(NNODE_1D * NNODE_1D);
853  // Set the dimensions of the element and the nodes, by default, both 2D
854  set_dimension(2);
855  // Assign default (full) spatial integration scheme
857  }
void set_dimension(const unsigned &dim)
Definition: elements.h:1380
void set_n_node(const unsigned &n)
Definition: elements.h:1404
virtual void set_integration_scheme(Integral *const &integral_pt)
Set the spatial integration scheme.
Definition: elements.cc:3210
static Gauss< 2, NNODE_1D > Default_integration_scheme
Assign the spatial integration scheme.
Definition: Qelements.h:845

◆ QElement() [2/2]

template<unsigned NNODE_1D>
oomph::QElement< 2, NNODE_1D >::QElement ( const QElement< 2, NNODE_1D > &  )
delete

Broken copy constructor.

Member Function Documentation

◆ bulk_coordinate_derivatives_fct_pt()

template<unsigned NNODE_1D>
BulkCoordinateDerivativesFctPt oomph::QElement< 2, NNODE_1D >::bulk_coordinate_derivatives_fct_pt ( const int face_index) const
inlinevirtual

Get a pointer to the derivative of the mapping from face to bulk coordinates.

Reimplemented from oomph::FiniteElement.

1203  {
1204  if (face_index == 1)
1205  {
1207  }
1208  else if (face_index == -1)
1209  {
1211  }
1212  else if (face_index == -2)
1213  {
1215  }
1216  else if (face_index == 2)
1217  {
1219  }
1220  else
1221  {
1222  std::string err = "Face index should be in {-1, +1}.";
1223  throw OomphLibError(
1225  }
1226  }
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
void faces1(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for the north and south faces, along which s1 is fixed.
Definition: Qelement_face_coordinate_translation_schemes.cc:131
void faces0(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for the east and west faces, along which s0 is fixed.
Definition: Qelement_face_coordinate_translation_schemes.cc:116
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References oomph::QElement2BulkCoordinateDerivatives::faces0(), oomph::QElement2BulkCoordinateDerivatives::faces1(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::Global_string_for_annotation::string().

◆ d2shape_local()

template<unsigned NNODE_1D>
void oomph::QElement< 2, NNODE_1D >::d2shape_local ( const Vector< double > &  s,
Shape psi,
DShape dpsids,
DShape d2psids 
) const
virtual

Compute the geometric shape functions, derivatives and second derivatives w.r.t. local coordinates at local coordinate s d2psids(i,0) = \( \partial^2 \psi_j / \partial s_0^2 \) d2psids(i,1) = \( \partial^2 \psi_j / \partial s_1^2 \) d2psids(i,2) = \( \partial^2 \psi_j / \partial s_0 \partial s_1 \)

Second derivatives of shape functions for specific QElement<2,..>: d2psids(i,0) = \( \partial^2 \psi_j / \partial s_0^2 \) d2psids(i,1) = \( \partial^2 \psi_j / \partial s_1^2 \) d2psids(i,2) = \( \partial^2 \psi_j / \partial s_0 \partial s_1 \)

Reimplemented from oomph::FiniteElement.

390  {
391  // Local storage
392  double Psi[2][NNODE_1D];
393  double DPsi[2][NNODE_1D];
394  double D2Psi[2][NNODE_1D];
395  // Index for the assembly process
396  unsigned index = 0;
397 
398  // Call the shape functions and derivatives
399  OneDimLagrange::shape<NNODE_1D>(s[0], Psi[0]);
400  OneDimLagrange::shape<NNODE_1D>(s[1], Psi[1]);
401  OneDimLagrange::dshape<NNODE_1D>(s[0], DPsi[0]);
402  OneDimLagrange::dshape<NNODE_1D>(s[1], DPsi[1]);
403  OneDimLagrange::d2shape<NNODE_1D>(s[0], D2Psi[0]);
404  OneDimLagrange::d2shape<NNODE_1D>(s[1], D2Psi[1]);
405 
406  // Loop over shape functions in element
407  for (unsigned i = 0; i < NNODE_1D; i++)
408  {
409  for (unsigned j = 0; j < NNODE_1D; j++)
410  {
411  // Assign the values
412  psi[index] = Psi[1][i] * Psi[0][j];
413  // First derivatives
414  dpsids(index, 0) = Psi[1][i] * DPsi[0][j];
415  dpsids(index, 1) = DPsi[1][i] * Psi[0][j];
416  // Second derivatives
417  // N.B. index 2 is the mixed derivative
418  d2psids(index, 0) = Psi[1][i] * D2Psi[0][j];
419  d2psids(index, 1) = D2Psi[1][i] * Psi[0][j];
420  d2psids(index, 2) = DPsi[1][i] * DPsi[0][j];
421  // Increment the index
422  ++index;
423  }
424  }
425  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
RealScalar s
Definition: level1_cplx_impl.h:130
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References i, j, and s.

◆ dshape_local()

template<unsigned NNODE_1D>
void oomph::QElement< 2, NNODE_1D >::dshape_local ( const Vector< double > &  s,
Shape psi,
DShape dpsids 
) const
virtual

Derivatives of shape functions for specific QElement<2,..>

Compute the geometric shape functions and derivatives w.r.t. local coordinates at local coordinate s

Reimplemented from oomph::FiniteElement.

351  {
352  // Local storage
353  double Psi[2][NNODE_1D];
354  double DPsi[2][NNODE_1D];
355  unsigned index = 0;
356 
357  // Call the shape functions and derivatives
358  OneDimLagrange::shape<NNODE_1D>(s[0], Psi[0]);
359  OneDimLagrange::shape<NNODE_1D>(s[1], Psi[1]);
360  OneDimLagrange::dshape<NNODE_1D>(s[0], DPsi[0]);
361  OneDimLagrange::dshape<NNODE_1D>(s[1], DPsi[1]);
362 
363  // Loop over shape functions in element
364  for (unsigned i = 0; i < NNODE_1D; i++)
365  {
366  for (unsigned j = 0; j < NNODE_1D; j++)
367  {
368  // Assign the values
369  dpsids(index, 0) = Psi[1][i] * DPsi[0][j];
370  dpsids(index, 1) = DPsi[1][i] * Psi[0][j];
371  psi[index] = Psi[1][i] * Psi[0][j];
372  // Increment the index
373  ++index;
374  }
375  }
376  }

References i, j, and s.

◆ face_outer_unit_normal_sign()

template<unsigned NNODE_1D>
int oomph::QElement< 2, NNODE_1D >::face_outer_unit_normal_sign ( const int face_index) const
inlinevirtual

Get the sign of the outer unit normal on the face given by face_index.

Reimplemented from oomph::FiniteElement.

1153  {
1154  if (face_index < 0)
1155  {
1156  return 1;
1157  }
1158  else if (face_index > 0)
1159  {
1160  return -1;
1161  }
1162  else
1163  {
1164  std::string err = "Face index should be one of {-1, +1, -2, +2}.";
1165  throw OomphLibError(
1167  }
1168  }

References OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::Global_string_for_annotation::string().

◆ face_to_bulk_coordinate_fct_pt()

template<unsigned NNODE_1D>
CoordinateMappingFctPt oomph::QElement< 2, NNODE_1D >::face_to_bulk_coordinate_fct_pt ( const int face_index) const
inlinevirtual

Get a pointer to the function mapping face coordinates to bulk coordinates

Reimplemented from oomph::FiniteElement.

1174  {
1175  if (face_index == 1)
1176  {
1178  }
1179  else if (face_index == -1)
1180  {
1182  }
1183  else if (face_index == -2)
1184  {
1186  }
1187  else if (face_index == 2)
1188  {
1190  }
1191  else
1192  {
1193  std::string err = "Face index should be in {-1, +1}.";
1194  throw OomphLibError(
1196  }
1197  }
void face2(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the east face (s0 = 1.0)
Definition: Qelement_face_coordinate_translation_schemes.cc:93
void face1(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the south face (s1 = -1.0)
Definition: Qelement_face_coordinate_translation_schemes.cc:86
void face0(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the west face (s0 = -1.0)
Definition: Qelement_face_coordinate_translation_schemes.cc:79
void face3(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the north face (s1 = 1.0)
Definition: Qelement_face_coordinate_translation_schemes.cc:100

References oomph::QElement2FaceToBulkCoordinates::face0(), oomph::QElement2FaceToBulkCoordinates::face1(), oomph::QElement2FaceToBulkCoordinates::face2(), oomph::QElement2FaceToBulkCoordinates::face3(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::Global_string_for_annotation::string().

◆ get_bulk_node_number()

template<unsigned NNODE_1D>
unsigned oomph::QElement< 2, NNODE_1D >::get_bulk_node_number ( const int face_index,
const unsigned i 
) const
inlinevirtual

Get the number of the ith node on face face_index (in the bulk node vector).

Reimplemented from oomph::FiniteElement.

1122  {
1124 
1125  const unsigned nn1d = nnode_1d();
1126 
1127  if (face_index == -1)
1128  {
1129  return i * nn1d;
1130  }
1131  else if (face_index == +1)
1132  {
1133  return nn1d * i + nn1d - 1;
1134  }
1135  else if (face_index == -2)
1136  {
1137  return i;
1138  }
1139  else if (face_index == +2)
1140  {
1141  return nn1d * (nn1d - 1) + i;
1142  }
1143  else
1144  {
1145  std::string err = "Face index should be in {-1, +1, -2, +2}.";
1146  throw OomphLibError(
1148  }
1149  }
void face_node_number_error_check(const unsigned &i) const
Range check for face node numbers.
Definition: elements.h:3395
unsigned nnode_1d() const
Number of nodes along each element edge.
Definition: Qelements.h:979

References i, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::Global_string_for_annotation::string().

◆ get_node_at_local_coordinate()

template<unsigned NNODE_1D>
Node * oomph::QElement< 2, NNODE_1D >::get_node_at_local_coordinate ( const Vector< double > &  s) const
virtual

Get the node at the specified local coordinate.

Return the node at the specified local coordinate.

Reimplemented from oomph::FiniteElement.

267  {
268  // Load the tolerance into a local variable
270  // There are two possible indices.
271  Vector<int> index(2);
272  // Loop over the coordinates and determine the indices
273  for (unsigned i = 0; i < 2; i++)
274  {
275  // If we are at the lower limit, the index is zero
276  if (std::fabs(s[i] + 1.0) < tol)
277  {
278  index[i] = 0;
279  }
280  // If we are at the upper limit, the index is the number of nodes minus 1
281  else if (std::fabs(s[i] - 1.0) < tol)
282  {
283  index[i] = NNODE_1D - 1;
284  }
285  // Otherwise, we have to calculate the index in general
286  else
287  {
288  // For uniformly spaced nodes the node number would be
289  double float_index = 0.5 * (1.0 + s[i]) * (NNODE_1D - 1);
290  // Convert to an integer by taking the floor (rounding down)
291  index[i] = static_cast<int>(std::floor(float_index));
292  // What is the excess. This should be safe because the
293  // we have rounded down
294  double excess = float_index - index[i];
295  // If the excess is bigger than our tolerance there is no node,
296  // return null
297  // Note that we test at both lower and upper ends.
298  if ((excess > tol) && ((1.0 - excess) > tol))
299  {
300  return 0;
301  }
302  // If we are at the upper end (i.e. the system has actually rounded up)
303  // we need to add one to the index
304  if ((1.0 - excess) <= tol)
305  {
306  index[i] += 1;
307  }
308  }
309  }
310  // If we've got here we have a node, so let's return a pointer to it
311  return node_pt(index[0] + NNODE_1D * index[1]);
312  }
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
static const double Node_location_tolerance
Definition: elements.h:1374
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 floor(const bfloat16 &a)
Definition: BFloat16.h:643
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117

References boost::multiprecision::fabs(), Eigen::bfloat16_impl::floor(), i, oomph::FiniteElement::Node_location_tolerance, and s.

◆ get_s_plot()

template<unsigned NNODE_1D>
void oomph::QElement< 2, NNODE_1D >::get_s_plot ( const unsigned i,
const unsigned nplot,
Vector< double > &  s,
const bool use_equally_spaced_interior_sample_points = false 
) const
inlinevirtual

Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direction).

Reimplemented from oomph::FiniteElement.

1072  {
1073  if (nplot > 1)
1074  {
1075  unsigned i0 = i % nplot;
1076  unsigned i1 = (i - i0) / nplot;
1077 
1078  s[0] = -1.0 + 2.0 * double(i0) / double(nplot - 1);
1079  s[1] = -1.0 + 2.0 * double(i1) / double(nplot - 1);
1080  if (use_equally_spaced_interior_sample_points)
1081  {
1082  double range = 2.0;
1083  double dx_new = range / double(nplot);
1084  double range_new = double(nplot - 1) * dx_new;
1085  s[0] = -1.0 + 0.5 * dx_new + range_new * (1.0 + s[0]) / range;
1086  s[1] = -1.0 + 0.5 * dx_new + range_new * (1.0 + s[1]) / range;
1087  }
1088  }
1089  else
1090  {
1091  s[0] = 0.0;
1092  s[1] = 0.0;
1093  }
1094  }

References i, and s.

Referenced by oomph::RefineableBuoyantQAxisymCrouzeixRaviartElement::output(), oomph::AxisymmetricQAdvectionCrouzeixRaviartElement::output(), and oomph::RefineableBuoyantQSphericalCrouzeixRaviartElement::output().

◆ invert_jacobian_mapping()

template<unsigned NNODE_1D>
double oomph::QElement< 2, NNODE_1D >::invert_jacobian_mapping ( const DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  inverse_jacobian 
) const
inlinevirtual

Overload the template-free interface for the calculation of inverse jacobian matrix. This is a two-dimensional element, so use the two-d version.

Reimplemented from oomph::FiniteElement.

889  {
890  return FiniteElement::invert_jacobian<2>(jacobian, inverse_jacobian);
891  }

◆ local_coordinate_of_node()

template<unsigned NNODE_1D>
void oomph::QElement< 2, NNODE_1D >::local_coordinate_of_node ( const unsigned j,
Vector< double > &  s 
) const
inlinevirtual

Get local coordinates of node j in the element; vector sets its own size.

Reimplemented from oomph::FiniteElement.

946  {
947  s.resize(2);
948  unsigned j0 = j % NNODE_1D;
949  unsigned j1 = j / NNODE_1D;
950  const double S_min = this->s_min();
951  const double S_range = this->s_max() - S_min;
952  s[0] = S_min + double(j0) / double(NNODE_1D - 1) * S_range;
953  s[1] = S_min + double(j1) / double(NNODE_1D - 1) * S_range;
954  }
double s_min() const
Min. value of local coordinate.
Definition: Qelements.h:894
double s_max() const
Max. value of local coordinate.
Definition: Qelements.h:900

References j, and s.

◆ local_fraction_of_node()

template<unsigned NNODE_1D>
void oomph::QElement< 2, NNODE_1D >::local_fraction_of_node ( const unsigned j,
Vector< double > &  s_fraction 
)
inlinevirtual

Get the local fraction of node j in the element.

Reimplemented from oomph::FiniteElement.

958  {
959  s_fraction.resize(2);
960  unsigned j0 = j % NNODE_1D;
961  unsigned j1 = j / NNODE_1D;
962  s_fraction[0] = double(j0) / double(NNODE_1D - 1);
963  s_fraction[1] = double(j1) / double(NNODE_1D - 1);
964  }

References j.

◆ local_one_d_fraction_of_node()

template<unsigned NNODE_1D>
double oomph::QElement< 2, NNODE_1D >::local_one_d_fraction_of_node ( const unsigned n1d,
const unsigned i 
)
inlinevirtual

This function returns the local fraction of ant nodes in the n-th positoin in a one dimensional expansion along the i-th local coordinate

Reimplemented from oomph::FiniteElement.

970  {
971  // It's just the value of the node divided by the number of 1-D nodes
972  return double(n1d) / double(NNODE_1D - 1);
973  }

◆ nnode_1d()

template<unsigned NNODE_1D>
unsigned oomph::QElement< 2, NNODE_1D >::nnode_1d ( ) const
inlinevirtual

Number of nodes along each element edge.

Reimplemented from oomph::FiniteElement.

980  {
981  return NNODE_1D;
982  }

◆ nplot_points()

template<unsigned NNODE_1D>
unsigned oomph::QElement< 2, NNODE_1D >::nplot_points ( const unsigned nplot) const
inlinevirtual

Return total number of plot points (when plotting nplot points in each "coordinate direction)

Reimplemented from oomph::FiniteElement.

1108  {
1109  unsigned DIM = 2;
1110  unsigned np = 1;
1111  for (unsigned i = 0; i < DIM; i++)
1112  {
1113  np *= nplot;
1114  }
1115  return np;
1116  };
#define DIM
Definition: linearised_navier_stokes_elements.h:44

References DIM, and i.

Referenced by oomph::RefineableBuoyantQAxisymCrouzeixRaviartElement::output(), oomph::AxisymmetricQAdvectionCrouzeixRaviartElement::output(), and oomph::RefineableBuoyantQSphericalCrouzeixRaviartElement::output().

◆ nplot_points_paraview()

template<unsigned NNODE_1D>
unsigned oomph::QElement< 2, NNODE_1D >::nplot_points_paraview ( const unsigned nplot) const
inlinevirtual

Return the number of actual plot points for paraview plot with parameter nplot.

Reimplemented from oomph::FiniteElement.

987  {
988  return nplot * nplot;
989  }

◆ nsub_elements_paraview()

template<unsigned NNODE_1D>
unsigned oomph::QElement< 2, NNODE_1D >::nsub_elements_paraview ( const unsigned nplot) const
inlinevirtual

Return the number of local sub-elements for paraview plot with parameter nplot.

Reimplemented from oomph::FiniteElement.

994  {
995  return (nplot - 1) * (nplot - 1);
996  }

◆ nvertex_node()

◆ output() [1/4]

template<unsigned NNODE_1D>
void oomph::QElement< 2, NNODE_1D >::output ( FILE *  file_pt)
virtual

◆ output() [2/4]

template<unsigned NNODE_1D>
void oomph::QElement< 2, NNODE_1D >::output ( FILE *  file_pt,
const unsigned n_plot 
)
virtual

C_style output at n_plot points.

C-style output function for n_plot points in each coordinate direction.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QTimeHarmonicFourierDecomposedLinearElasticityElement< NNODE_1D >, oomph::QSteadyAxisymAdvectionDiffusionElement< NNODE_1D >, oomph::QSphericalAdvectionDiffusionElement< NNODE_1D >, oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >, oomph::SolidQElement< 2, NNODE_1D >, and oomph::QAxisymAdvectionDiffusionElement< 3 >.

487  {
488  // Local variables
489  Vector<double> s(2);
490 
491  // Find the dimension of the nodes
492  unsigned n_dim = this->nodal_dimension();
493 
494  // Tecplot header info
495  // outfile << "ZONE I=" << n_plot << ", J=" << n_plot << std::endl;
496  fprintf(file_pt, "ZONE I=%i, J=%i\n", n_plot, n_plot);
497 
498  // Loop over element nodes
499  for (unsigned l2 = 0; l2 < n_plot; l2++)
500  {
501  s[1] = -1.0 + l2 * 2.0 / (n_plot - 1);
502  for (unsigned l1 = 0; l1 < n_plot; l1++)
503  {
504  s[0] = -1.0 + l1 * 2.0 / (n_plot - 1);
505 
506  // Output the coordinates
507  for (unsigned i = 0; i < n_dim; i++)
508  {
509  // outfile << interpolated_x(s,i) << " " ;
510  fprintf(file_pt, "%g ", interpolated_x(s, i));
511  }
512  // outfile << std::endl;
513  fprintf(file_pt, "\n");
514  }
515  }
516  // outfile << std::endl;
517  fprintf(file_pt, "\n");
518  }
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
Definition: elements.h:2484

References i, and s.

◆ output() [3/4]

◆ output() [4/4]

template<unsigned NNODE_1D>
void oomph::QElement< 2, NNODE_1D >::output ( std::ostream &  outfile,
const unsigned n_plot 
)
virtual

Output at n_plot points.

The output function for n_plot points in each coordinate direction.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QYoungLaplaceElement< NNODE_1D >, oomph::QTimeHarmonicFourierDecomposedLinearElasticityElement< NNODE_1D >, oomph::QSteadyAxisymAdvectionDiffusionElement< NNODE_1D >, oomph::QSphericalAdvectionDiffusionElement< NNODE_1D >, oomph::QPMLFourierDecomposedHelmholtzElement< NNODE_1D >, oomph::SolidQElement< 2, NNODE_1D >, and oomph::QAxisymAdvectionDiffusionElement< 3 >.

443  {
444  // Local variables
445  Vector<double> s(2);
446 
447  // Tecplot header info
448  outfile << "ZONE I=" << n_plot << ", J=" << n_plot << std::endl;
449 
450  // Find the dimension of the nodes
451  unsigned n_dim = this->nodal_dimension();
452 
453  // Loop over plot points
454  for (unsigned l2 = 0; l2 < n_plot; l2++)
455  {
456  s[1] = -1.0 + l2 * 2.0 / (n_plot - 1);
457  for (unsigned l1 = 0; l1 < n_plot; l1++)
458  {
459  s[0] = -1.0 + l1 * 2.0 / (n_plot - 1);
460 
461  // Output the coordinates
462  for (unsigned i = 0; i < n_dim; i++)
463  {
464  outfile << interpolated_x(s, i) << " ";
465  }
466  outfile << std::endl;
467  }
468  }
469  outfile << std::endl;
470  }

References i, and s.

◆ s_max()

template<unsigned NNODE_1D>
double oomph::QElement< 2, NNODE_1D >::s_max ( ) const
inlinevirtual

Max. value of local coordinate.

Reimplemented from oomph::FiniteElement.

901  {
902  return 1.0;
903  }

◆ s_min()

template<unsigned NNODE_1D>
double oomph::QElement< 2, NNODE_1D >::s_min ( ) const
inlinevirtual

Min. value of local coordinate.

Reimplemented from oomph::FiniteElement.

895  {
896  return -1.0;
897  }

◆ shape()

template<unsigned NNODE_1D>
void oomph::QElement< 2, NNODE_1D >::shape ( const Vector< double > &  s,
Shape psi 
) const
virtual

Broken assignment operator.

Calculate the geometric shape functions at local coordinate s

Shape function for specific QElement<2,..>

Implements oomph::FiniteElement.

321  {
322  // Local storage
323  double Psi[2][NNODE_1D];
324  // Call the OneDimensional Shape functions
325  OneDimLagrange::shape<NNODE_1D>(s[0], Psi[0]);
326  OneDimLagrange::shape<NNODE_1D>(s[1], Psi[1]);
327  // Index for the shape functions
328  unsigned index = 0;
329 
330  // Now let's loop over the nodal points in the element
331  // s1 is the "x" coordinate, s2 the "y"
332  for (unsigned i = 0; i < NNODE_1D; i++)
333  {
334  for (unsigned j = 0; j < NNODE_1D; j++)
335  {
336  // Multiply the two 1D functions together to get the 2D function
337  psi[index] = Psi[1][i] * Psi[0][j];
338  // Incremenet the index
339  ++index;
340  }
341  }
342  }

References i, j, and s.

◆ tecplot_zone_string()

template<unsigned NNODE_1D>
std::string oomph::QElement< 2, NNODE_1D >::tecplot_zone_string ( const unsigned nplot) const
inlinevirtual

Return string for tecplot zone header (when plotting nplot points in each "coordinate direction)

Reimplemented from oomph::FiniteElement.

1099  {
1100  std::ostringstream header;
1101  header << "ZONE I=" << nplot << ", J=" << nplot << "\n";
1102  return header.str();
1103  }

Referenced by oomph::RefineableBuoyantQAxisymCrouzeixRaviartElement::output(), oomph::AxisymmetricQAdvectionCrouzeixRaviartElement::output(), and oomph::RefineableBuoyantQSphericalCrouzeixRaviartElement::output().

◆ vertex_node_pt()

template<unsigned NNODE_1D>
Node* oomph::QElement< 2, NNODE_1D >::vertex_node_pt ( const unsigned j) const
inlinevirtual

Pointer to the j-th vertex node in the element.

Implements oomph::QuadElementBase.

Reimplemented in oomph::RefineableQYoungLaplaceElement< NNODE_1D >, oomph::RefineableQSphericalAdvectionDiffusionElement< NNODE_1D >, oomph::RefineableQSphericalAdvectionDiffusionElement< 3 >, and oomph::RefineableQAxisymAdvectionDiffusionElement< 3 >.

914  {
915  unsigned n_node_1d = nnode_1d();
916  Node* nod_pt;
917  switch (j)
918  {
919  case 0:
920  nod_pt = node_pt(0);
921  break;
922  case 1:
923  nod_pt = node_pt(n_node_1d - 1);
924  break;
925  case 2:
926  nod_pt = node_pt(n_node_1d * (n_node_1d - 1));
927  break;
928  case 3:
929  nod_pt = node_pt(n_node_1d * n_node_1d - 1);
930  break;
931  default:
932  std::ostringstream error_message;
933  error_message << "Vertex node number is " << j
934  << " but must be from 0 to 3\n";
935 
936  throw OomphLibError(error_message.str(),
939  }
940  return nod_pt;
941  }

References j, OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::RefineableQAxisymAdvectionDiffusionElement< NNODE_1D >::vertex_node_pt(), oomph::RefineableQSphericalAdvectionDiffusionElement< NNODE_1D >::vertex_node_pt(), and oomph::RefineableQYoungLaplaceElement< NNODE_1D >::vertex_node_pt().

◆ write_paraview_offsets()

template<unsigned NNODE_1D>
void oomph::QElement< 2, NNODE_1D >::write_paraview_offsets ( std::ofstream &  file_out,
const unsigned nplot,
unsigned offset_sum 
) const
inlinevirtual

Return the offsets for the paraview sub-elements. Needs to be implemented for each new geometric element type; see http://www.vtk.org/VTK/img/file-formats.pdf

Reimplemented from oomph::FiniteElement.

1041  {
1042  // Loop over all local elements and add its offset to the overall
1043  // offset_sum
1044  unsigned local_loop = nsub_elements_paraview(nplot);
1045  for (unsigned i = 0; i < local_loop; i++)
1046  {
1047  offset_sum += 4;
1048  file_out << offset_sum << std::endl;
1049  }
1050  }
unsigned nsub_elements_paraview(const unsigned &nplot) const
Definition: Qelements.h:993

References i.

◆ write_paraview_output_offset_information()

template<unsigned NNODE_1D>
void oomph::QElement< 2, NNODE_1D >::write_paraview_output_offset_information ( std::ofstream &  file_out,
const unsigned nplot,
unsigned counter 
) const
inlinevirtual

Fill in the offset information for paraview plot. Needs to be implemented for each new geometric element type; see http://www.vtk.org/VTK/img/file-formats.pdf

Reimplemented from oomph::FiniteElement.

1004  {
1005  // Number of local elements we want to plot over
1006  unsigned plot = nsub_elements_paraview(nplot);
1007 
1008  // loops over the i-th local element in parent element
1009  for (unsigned i = 0; i < plot; i++)
1010  {
1011  unsigned d = (i - (i % (nplot - 1))) / (nplot - 1);
1012 
1013  file_out << i % (nplot - 1) + d * nplot + counter << " "
1014  << i % (nplot - 1) + 1 + d * nplot + counter << " "
1015  << i % (nplot - 1) + 1 + (d + 1) * nplot + counter << " "
1016  << i % (nplot - 1) + (d + 1) * nplot + counter << std::endl;
1017  }
1018  counter += nplot_points_paraview(nplot);
1019  }
unsigned nplot_points_paraview(const unsigned &nplot) const
Definition: Qelements.h:986
void plot()
Plot.
Definition: sphere_scattering.cc:180

References i, and PlanarWave::plot().

◆ write_paraview_type()

template<unsigned NNODE_1D>
void oomph::QElement< 2, NNODE_1D >::write_paraview_type ( std::ofstream &  file_out,
const unsigned nplot 
) const
inlinevirtual

Return the paraview element type. Needs to be implemented for each new geometric element type; see http://www.vtk.org/VTK/img/file-formats.pdf Use type "VTK_QUAD" (== 9) for 2D quad elements

Reimplemented from oomph::FiniteElement.

1027  {
1028  unsigned local_loop = nsub_elements_paraview(nplot);
1029  for (unsigned i = 0; i < local_loop; i++)
1030  {
1031  file_out << "9" << std::endl;
1032  }
1033  }

References i.

Member Data Documentation

◆ Default_integration_scheme

template<unsigned NNODE_1D>
Gauss< 2, NNODE_1D > oomph::QElement< 2, NNODE_1D >::Default_integration_scheme
staticprivate

Assign the spatial integration scheme.

Default integration rule: Gaussian integration of same 'order' as the element


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