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

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

#include <Qspectral_elements.h>

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

Public Member Functions

 QSpectralElement ()
 Constructor. More...
 
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 &n, 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 &n, Vector< double > &s_fraction)
 Get the local fractino of node j in the element. More...
 
double local_one_d_fraction_of_node (const unsigned &n1d, const unsigned &i)
 The local one-d fraction is the same. More...
 
void shape (const Vector< double > &s, Shape &psi) const
 Calculate the geometric shape functions at local coordinate s. More...
 
void dshape_local (const Vector< double > &s, Shape &psi, DShape &dpsids) const
 Derivatives of shape functions for specific QSpectralElement<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
 
unsigned nnode_1d () const
 Number of nodes along each element edge. 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 output (std::ostream &outfile)
 Output. More...
 
void output (std::ostream &outfile, const unsigned &n_plot)
 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
 
void build_face_element (const int &face_index, FaceElement *face_element_pt)
 
- Public Member Functions inherited from oomph::SpectralElement
 SpectralElement ()
 
virtual ~SpectralElement ()
 
Dataspectral_data_pt (const unsigned &i) const
 
virtual void describe_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)
 
unsigned nspectral () 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...
 
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_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...
 
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
 
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 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)
 
virtual unsigned nplot_points_paraview (const unsigned &nplot) const
 
virtual unsigned nsub_elements_paraview (const unsigned &nplot) const
 
void output_paraview (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_output_offset_information (std::ofstream &file_out, const unsigned &nplot, unsigned &counter) const
 
virtual void write_paraview_type (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_offsets (std::ofstream &file_out, const unsigned &nplot, unsigned &offset_sum) const
 
virtual unsigned nscalar_paraview () const
 
virtual void scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
 
virtual std::string scalar_name_paraview (const unsigned &i) const
 
virtual void output (const unsigned &t, std::ostream &outfile, const unsigned &n_plot) const
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 Output an exact solution over the element. More...
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
 Output a time-dependent exact solution over the element. More...
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, const SolutionFunctorBase &exact_soln) const
 Output a time-dependent exact solution over the element. More...
 
virtual void 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 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...
 
void face_node_number_error_check (const unsigned &i) const
 Range check for face node numbers. More...
 
virtual CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt (const int &face_index) const
 
virtual BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt (const int &face_index) const
 
- Public Member Functions inherited from oomph::GeneralisedElement
 GeneralisedElement ()
 Constructor: Initialise all pointers and all values to zero. More...
 
virtual ~GeneralisedElement ()
 Virtual destructor to clean up any memory allocated by the object. More...
 
 GeneralisedElement (const GeneralisedElement &)=delete
 Broken copy constructor. More...
 
void operator= (const GeneralisedElement &)=delete
 Broken assignment operator. More...
 
Data *& internal_data_pt (const unsigned &i)
 Return a pointer to i-th internal data object. More...
 
Data *const & internal_data_pt (const unsigned &i) const
 Return a pointer to i-th internal data object (const version) More...
 
Data *& external_data_pt (const unsigned &i)
 Return a pointer to i-th external data object. More...
 
Data *const & external_data_pt (const unsigned &i) const
 Return a pointer to i-th external data object (const version) More...
 
unsigned long eqn_number (const unsigned &ieqn_local) const
 
int local_eqn_number (const unsigned long &ieqn_global) const
 
unsigned add_external_data (Data *const &data_pt, const bool &fd=true)
 
bool external_data_fd (const unsigned &i) const
 
void exclude_external_data_fd (const unsigned &i)
 
void include_external_data_fd (const unsigned &i)
 
void flush_external_data ()
 Flush all external data. More...
 
void flush_external_data (Data *const &data_pt)
 Flush the object addressed by data_pt from the external data array. More...
 
unsigned ninternal_data () const
 Return the number of internal data objects. More...
 
unsigned nexternal_data () const
 Return the number of external data objects. More...
 
unsigned ndof () const
 Return the number of equations/dofs in the element. More...
 
void dof_vector (const unsigned &t, Vector< double > &dof)
 Return the vector of dof values at time level t. More...
 
void dof_pt_vector (Vector< double * > &dof_pt)
 Return the vector of pointers to dof values. More...
 
void set_internal_data_time_stepper (const unsigned &i, TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void assign_internal_eqn_numbers (unsigned long &global_number, Vector< double * > &Dof_pt)
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
void add_internal_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt)
 
virtual void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void complete_setup_of_dependencies ()
 
virtual void get_residuals (Vector< double > &residuals)
 
virtual void get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
virtual void get_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void get_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void get_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void get_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void get_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void get_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void get_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void get_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_norm (double &norm)
 
virtual unsigned ndof_types () const
 
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 
- Public Member Functions inherited from oomph::GeomObject
 GeomObject ()
 Default constructor. More...
 
 GeomObject (const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim, TimeStepper *time_stepper_pt)
 
 GeomObject (const GeomObject &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const GeomObject &)=delete
 Broken assignment operator. More...
 
virtual ~GeomObject ()
 (Empty) destructor More...
 
unsigned nlagrangian () const
 Access function to # of Lagrangian coordinates. More...
 
unsigned ndim () const
 Access function to # of Eulerian coordinates. More...
 
void set_nlagrangian_and_ndim (const unsigned &n_lagrangian, const unsigned &n_dim)
 Set # of Lagrangian and Eulerian coordinates. More...
 
TimeStepper *& time_stepper_pt ()
 
TimeSteppertime_stepper_pt () const
 
virtual void position (const double &t, const Vector< double > &zeta, Vector< double > &r) const
 
virtual void dposition (const Vector< double > &zeta, DenseMatrix< double > &drdzeta) const
 
virtual void d2position (const Vector< double > &zeta, RankThreeTensor< double > &ddrdzeta) const
 
virtual void d2position (const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const
 
- Public 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...
 

Static Private Attributes

static GaussLobattoLegendre< 2, NNODE_1D > integral
 Assign the static integral. 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::SpectralElement
Vector< Data * > * Spectral_data_pt
 Additional Storage for shared spectral data. More...
 
Vector< unsignedSpectral_order
 Vector that represents the spectral order in each dimension. More...
 
Vector< unsignedNodal_spectral_order
 Vector that represents the nodal spectral order. More...
 
DenseMatrix< intSpectral_local_eqn
 Local equation numbers for the spectral degrees of freedom. More...
 
- 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::QSpectralElement< 2, NNODE_1D >

General QSpectralElement class specialised to two spatial dimensions.

Constructor & Destructor Documentation

◆ QSpectralElement()

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

Constructor.

663  {
664  // There are NNODE_1D^2 nodes in this element
665  this->set_n_node(NNODE_1D * NNODE_1D);
666  Spectral_order.resize(2, NNODE_1D);
667  Nodal_spectral_order.resize(2, NNODE_1D);
668  // Set the elemental and nodal dimensions
669  this->set_dimension(2);
670  // Assign integral pointer
672  // Calculate the nodal positions for the shape functions
674  // OneDLegendreShapeParam::calculate_nodal_positions(NNODE_1D);
675  }
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 void calculate_nodal_positions()
Static function used to populate the stored positions.
Definition: shape.h:1241
static GaussLobattoLegendre< 2, NNODE_1D > integral
Assign the static integral.
Definition: Qspectral_elements.h:658
Vector< unsigned > Nodal_spectral_order
Vector that represents the nodal spectral order.
Definition: Qspectral_elements.h:162
Vector< unsigned > Spectral_order
Vector that represents the spectral order in each dimension.
Definition: Qspectral_elements.h:159

References oomph::OneDimensionalLegendreShape< NNODE_1D >::calculate_nodal_positions().

Member Function Documentation

◆ build_face_element()

template<unsigned NNODE_1D>
void oomph::QSpectralElement< 2, NNODE_1D >::build_face_element ( const int face_index,
FaceElement face_element_pt 
)
virtual

Build the lower-dimensional FaceElement (an element of type QSpectralElement<1,NNODE_1D>). The face index takes four values corresponding to the four possible faces: -1 (West) s[0] = -1.0 +1 (East) s[0] = 1.0 -2 (South) s[1] = -1.0 +2 (North) s[1] = 1.0

Function to setup geometrical information for lower-dimensional FaceElements (which are of type QSpectralElement<0,1>).

Reimplemented from oomph::FiniteElement.

307  {
308  // Set the nodal dimension from the "bulk"
309  face_element_pt->set_nodal_dimension(this->node_pt(0)->ndim());
310 
311  // Set the pointer to the "bulk" element
312  face_element_pt->bulk_element_pt() = this;
313 
314 #ifdef OOMPH_HAS_MPI
315  // If the bulk element is halo then the face element must be too
316  // if (this->is_halo())
317  {
318  face_element_pt->set_halo(Non_halo_proc_ID);
319  }
320 #endif
321 
322  // Resize the storage for the original number of values at
323  // NNODE_1D nodes of the FaceElement
324  face_element_pt->nbulk_value_resize(NNODE_1D);
325 
326  // Resize the storage for the bulk node numbers corresponding
327  // to the NNODE_1D nodes of the FaceElement
328  face_element_pt->bulk_node_number_resize(NNODE_1D);
329 
330  // Set the face index in the face element
331  // The faces are
332  // +1 East
333  // -1 West
334  // +2 North
335  // -2 South
336 
337  // Set the face index in the face element
338  face_element_pt->face_index() = face_index;
339 
340  // Now set up the node pointers
341  // =================================
342  // The convention here is that interior_tangent X tangent X tangent
343  // is the OUTWARD normal
344  switch (face_index)
345  {
346  unsigned bulk_number;
347  // West face, normal sign is positive
348  case (-1):
349  // Set the pointer to the bulk coordinate translation scheme
350  face_element_pt->face_to_bulk_coordinate_fct_pt() =
352 
353  // Set the pointer to the derivative mappings
354  face_element_pt->bulk_coordinate_derivatives_fct_pt() =
356 
357  for (unsigned i = 0; i < NNODE_1D; i++)
358  {
359  bulk_number = i * NNODE_1D;
360  face_element_pt->node_pt(i) = this->node_pt(bulk_number);
361  face_element_pt->bulk_node_number(i) = bulk_number;
362  face_element_pt->normal_sign() = 1;
363  // Set the number of values originally stored at this node
364  face_element_pt->nbulk_value(i) = this->required_nvalue(bulk_number);
365  }
366  break;
367  // South face, normal sign is positive
368  case (-2):
369  // Set the pointer to the bulk coordinate translation scheme
370  face_element_pt->face_to_bulk_coordinate_fct_pt() =
372 
373  // Set the pointer to the derivative mappings
374  face_element_pt->bulk_coordinate_derivatives_fct_pt() =
376 
377  for (unsigned i = 0; i < NNODE_1D; i++)
378  {
379  bulk_number = i;
380  face_element_pt->node_pt(i) = this->node_pt(bulk_number);
381  face_element_pt->bulk_node_number(i) = bulk_number;
382  face_element_pt->normal_sign() = 1;
383  // Set the number of values originally stored at this node
384  face_element_pt->nbulk_value(i) = this->required_nvalue(bulk_number);
385  }
386  break;
387  // East face, normal sign is negative
388  case (1):
389  // Set the pointer to the bulk coordinate translation scheme
390  face_element_pt->face_to_bulk_coordinate_fct_pt() =
392 
393  // Set the pointer to the derivative mappings
394  face_element_pt->bulk_coordinate_derivatives_fct_pt() =
396 
397  for (unsigned i = 0; i < NNODE_1D; i++)
398  {
399  bulk_number = NNODE_1D * i + NNODE_1D - 1;
400  face_element_pt->node_pt(i) = this->node_pt(bulk_number);
401  face_element_pt->bulk_node_number(i) = bulk_number;
402  face_element_pt->normal_sign() = -1;
403  // Set the number of values originally stored at this node
404  face_element_pt->nbulk_value(i) = this->required_nvalue(bulk_number);
405  }
406  break;
407  // North face, normal sign is negative
408  case (2):
409  // Set the pointer to the bulk coordinate translation scheme
410  face_element_pt->face_to_bulk_coordinate_fct_pt() =
412 
413  // Set the pointer to the derivative mappings
414  face_element_pt->bulk_coordinate_derivatives_fct_pt() =
416 
417  for (unsigned i = 0; i < NNODE_1D; i++)
418  {
419  bulk_number = NNODE_1D * (NNODE_1D - 1) + i;
420  face_element_pt->node_pt(i) = this->node_pt(bulk_number);
421  face_element_pt->bulk_node_number(i) = bulk_number;
422  face_element_pt->normal_sign() = -1;
423  // Set the number of values originally stored at this node
424  face_element_pt->nbulk_value(i) = this->required_nvalue(bulk_number);
425  }
426  break;
427 
428  // Now cover the other cases
429  default:
430  std::ostringstream error_message;
431  error_message
432  << "Face index should only take the values +/- 1 or +/- 2,"
433  << " not " << face_index << std::endl;
434  throw OomphLibError(error_message.str(),
437  }
438  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
virtual unsigned required_nvalue(const unsigned &n) const
Definition: elements.h:2455
unsigned ndim() const
Access function to # of Eulerian coordinates.
Definition: geom_objects.h:177
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
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
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References oomph::FaceElement::bulk_coordinate_derivatives_fct_pt(), oomph::FaceElement::bulk_element_pt(), oomph::FaceElement::bulk_node_number(), oomph::FaceElement::bulk_node_number_resize(), oomph::QElement2FaceToBulkCoordinates::face0(), oomph::QElement2FaceToBulkCoordinates::face1(), oomph::QElement2FaceToBulkCoordinates::face2(), oomph::QElement2FaceToBulkCoordinates::face3(), oomph::FaceElement::face_index(), oomph::FaceElement::face_to_bulk_coordinate_fct_pt(), oomph::QElement2BulkCoordinateDerivatives::faces0(), oomph::QElement2BulkCoordinateDerivatives::faces1(), i, oomph::FaceElement::nbulk_value(), oomph::FaceElement::nbulk_value_resize(), oomph::FiniteElement::node_pt(), oomph::FaceElement::normal_sign(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::FiniteElement::set_nodal_dimension().

◆ d2shape_local()

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

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

938  {
939  std::ostringstream error_message;
940  error_message
941  << "\nd2shpe_local currently not implemented for this element\n";
942  throw OomphLibError(
943  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
944 
945  /* //Call the shape functions and derivatives */
946  /* OneDimensionalShape<NNODE_1D> psi1(s[0]); */
947  /* OneDimensionalDShape<NNODE_1D> dpsi1ds(s[0]); */
948  /* OneDimensionalD2Shape<NNODE_1D> d2psi1ds(s[0]); */
949 
950  /* //Loop over shape functions in element and assign to psi */
951  /* for(unsigned l=0;l<NNODE_1D;l++) */
952  /* { */
953  /* psi[l] = psi1[l]; */
954  /* dpsids[l][0] = dpsi1ds[l]; */
955  /* d2psids[l][0] = d2psi1ds[l]; */
956  /* } */
957  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ dshape_local()

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

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

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

Reimplemented from oomph::FiniteElement.

901  {
902  // Call the shape functions and derivatives
903  OneDimensionalLegendreShape<NNODE_1D> psi1(s[0]), psi2(s[1]);
904  OneDimensionalLegendreDShape<NNODE_1D> dpsi1ds(s[0]), dpsi2ds(s[1]);
905 
906  // OneDLegendreShapeParam psi1(NNODE_1D,s[0]), psi2(NNODE_1D,s[1]);
907  // OneDLegendreDShapeParam dpsi1ds(NNODE_1D,s[0]), dpsi2ds(NNODE_1D,s[1]);
908 
909  // Index for the shape functions
910  unsigned index = 0;
911  // Loop over shape functions in element
912  for (unsigned i = 0; i < NNODE_1D; i++)
913  {
914  for (unsigned j = 0; j < NNODE_1D; j++)
915  {
916  // Assign the values
917  dpsids(index, 0) = psi2[i] * dpsi1ds[j];
918  dpsids(index, 1) = dpsi2ds[i] * psi1[j];
919  psi[index] = psi2[i] * psi1[j];
920  // Increment the index
921  ++index;
922  }
923  }
924  }
RealScalar s
Definition: level1_cplx_impl.h:130
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References i, j, and s.

◆ get_s_plot()

template<unsigned NNODE_1D>
void oomph::QSpectralElement< 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.

813  {
814  if (nplot > 1)
815  {
816  unsigned i0 = i % nplot;
817  unsigned i1 = (i - i0) / nplot;
818 
819  s[0] = -1.0 + 2.0 * double(i0) / double(nplot - 1);
820  s[1] = -1.0 + 2.0 * double(i1) / double(nplot - 1);
821  if (use_equally_spaced_interior_sample_points)
822  {
823  double range = 2.0;
824  double dx_new = range / double(nplot);
825  double range_new = double(nplot - 1) * dx_new;
826  s[0] = -1.0 + 0.5 * dx_new + range_new * (1.0 + s[0]) / range;
827  s[1] = -1.0 + 0.5 * dx_new + range_new * (1.0 + s[1]) / range;
828  }
829  }
830  else
831  {
832  s[0] = 0.0;
833  s[1] = 0.0;
834  }
835  }

References i, and s.

◆ invert_jacobian_mapping()

template<unsigned NNODE_1D>
double oomph::QSpectralElement< 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 one-dimensional element, so use the 1D version.

Reimplemented from oomph::FiniteElement.

777  {
778  return FiniteElement::invert_jacobian<2>(jacobian, inverse_jacobian);
779  }

◆ local_coordinate_of_node()

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

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

Reimplemented from oomph::FiniteElement.

729  {
730  s.resize(2);
731  unsigned n0 = n % NNODE_1D;
732  unsigned n1 = unsigned(double(n) / double(NNODE_1D));
735  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
static double nodal_position(const unsigned &n)
Definition: shape.h:1250

References n, oomph::OneDimensionalLegendreShape< NNODE_1D >::nodal_position(), and s.

◆ local_fraction_of_node()

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

Get the local fractino of node j in the element.

Reimplemented from oomph::FiniteElement.

739  {
740  this->local_coordinate_of_node(n, s_fraction);
741  // Resize
742  s_fraction[0] = 0.5 * (s_fraction[0] + 1.0);
743  s_fraction[1] = 0.5 * (s_fraction[1] + 1.0);
744  }
void local_coordinate_of_node(const unsigned &n, Vector< double > &s) const
Get local coordinates of node j in the element; vector sets its own size.
Definition: Qspectral_elements.h:728

References n.

◆ local_one_d_fraction_of_node()

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

The local one-d fraction is the same.

Reimplemented from oomph::FiniteElement.

748  {
749  return 0.5 *
751  }

References oomph::OneDimensionalLegendreShape< NNODE_1D >::nodal_position().

◆ nnode_1d()

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

Number of nodes along each element edge.

Reimplemented from oomph::FiniteElement.

784  {
785  return NNODE_1D;
786  }

◆ nplot_points()

template<unsigned NNODE_1D>
unsigned oomph::QSpectralElement< 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.

850  {
851  unsigned DIM = 2;
852  unsigned np = 1;
853  for (unsigned i = 0; i < DIM; i++)
854  {
855  np *= nplot;
856  }
857  return np;
858  }
#define DIM
Definition: linearised_navier_stokes_elements.h:44

References DIM, and i.

◆ nvertex_node()

template<unsigned NNODE_1D>
unsigned oomph::QSpectralElement< 2, NNODE_1D >::nvertex_node ( ) const
inlinevirtual

Number of vertex nodes in the element.

Implements oomph::QuadElementBase.

691  {
692  return 4;
693  }

◆ output() [1/4]

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

C-style output.

Reimplemented from oomph::FiniteElement.

790  {
791  FiniteElement::output(file_pt);
792  }
virtual void output(std::ostream &outfile)
Definition: elements.h:3050

References oomph::FiniteElement::output().

◆ output() [2/4]

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

C_style output at n_plot points.

Reimplemented from oomph::FiniteElement.

796  {
797  FiniteElement::output(file_pt, n_plot);
798  }

References oomph::FiniteElement::output().

◆ output() [3/4]

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

Output.

The output function for general 1D QSpectralElements.

Reimplemented from oomph::FiniteElement.

231  {
232  // Tecplot header info
233  outfile << "ZONE I=" << NNODE_1D << ", J=" << NNODE_1D << std::endl;
234 
235  // Find the dimension of the nodes
236  unsigned n_dim = this->nodal_dimension();
237 
238  // Loop over element nodes
239  for (unsigned l2 = 0; l2 < NNODE_1D; l2++)
240  {
241  for (unsigned l1 = 0; l1 < NNODE_1D; l1++)
242  {
243  unsigned l = l2 * NNODE_1D + l1;
244 
245  // Loop over the dimensions and output the position
246  for (unsigned i = 0; i < n_dim; i++)
247  {
248  outfile << this->node_pt(l)->x(i) << " ";
249  }
250  // Find out how many data values at the node
251  unsigned initial_nvalue = this->node_pt(l)->nvalue();
252  // Loop over the data and output whether pinned or not
253  for (unsigned i = 0; i < initial_nvalue; i++)
254  {
255  outfile << this->node_pt(l)->is_pinned(i) << " ";
256  }
257  outfile << std::endl;
258  }
259  }
260  outfile << std::endl;
261  }
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
Definition: nodes.h:483
bool is_pinned(const unsigned &i) const
Test whether the i-th variable is pinned (1: true; 0: false).
Definition: nodes.h:417
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
Definition: elements.h:2484
double & x(const unsigned &i)
Return the i-th nodal coordinate.
Definition: nodes.h:1060

References i.

◆ output() [4/4]

template<unsigned NNODE_1D>
void oomph::QSpectralElement< 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.

270  {
271  // Local variables
272  Vector<double> s(2);
273 
274  // Tecplot header info
275  outfile << "ZONE I=" << n_plot << ", J=" << n_plot << std::endl;
276 
277  // Find the dimension of the first node
278  unsigned n_dim = this->nodal_dimension();
279 
280  // Loop over plot points
281  for (unsigned l2 = 0; l2 < n_plot; l2++)
282  {
283  s[1] = -1.0 + l2 * 2.0 / (n_plot - 1);
284  for (unsigned l1 = 0; l1 < n_plot; l1++)
285  {
286  s[0] = -1.0 + l1 * 2.0 / (n_plot - 1);
287 
288  // Output the coordinates
289  for (unsigned i = 0; i < n_dim; i++)
290  {
291  outfile << this->interpolated_x(s, i) << " ";
292  }
293  outfile << std::endl;
294  }
295  }
296  outfile << std::endl;
297  }
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

References i, and s.

◆ s_max()

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

Max. value of local coordinate.

Reimplemented from oomph::FiniteElement.

685  {
686  return 1.0;
687  }

◆ s_min()

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

Min. value of local coordinate.

Reimplemented from oomph::FiniteElement.

679  {
680  return -1.0;
681  }

◆ shape()

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

Calculate the geometric shape functions at local coordinate s.

Shape function for specific QSpectralElement<2,..>

Implements oomph::FiniteElement.

878  {
879  // Call the OneDimensional Shape functions
880  OneDimensionalLegendreShape<NNODE_1D> psi1(s[0]), psi2(s[1]);
881  // OneDLegendreShapeParam psi1(NNODE_1D,s[0]), psi2(NNODE_1D,s[1]);
882 
883  // Now let's loop over the nodal points in the element
884  // and copy the values back in
885  for (unsigned i = 0; i < NNODE_1D; i++)
886  {
887  for (unsigned j = 0; j < NNODE_1D; j++)
888  {
889  psi(NNODE_1D * i + j) = psi2[i] * psi1[j];
890  }
891  }
892  }

References i, j, and s.

◆ tecplot_zone_string()

template<unsigned NNODE_1D>
std::string oomph::QSpectralElement< 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.

841  {
842  std::ostringstream header;
843  header << "ZONE I=" << nplot << ", J=" << nplot << "\n";
844  return header.str();
845  }

◆ vertex_node_pt()

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

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

Implements oomph::QuadElementBase.

697  {
698  unsigned n_node_1d = nnode_1d();
699  Node* nod_pt;
700  switch (j)
701  {
702  case 0:
703  nod_pt = this->node_pt(0);
704  break;
705  case 1:
706  nod_pt = this->node_pt(n_node_1d - 1);
707  break;
708  case 2:
709  nod_pt = this->node_pt(n_node_1d * (n_node_1d - 1));
710  break;
711  case 3:
712  nod_pt = this->node_pt(n_node_1d * n_node_1d - 1);
713  break;
714  default:
715  std::ostringstream error_message;
716  error_message << "Vertex node number is " << j
717  << " but must be from 0 to 3\n";
718 
719  throw OomphLibError(error_message.str(),
722  }
723  return nod_pt;
724  }
unsigned nnode_1d() const
Number of nodes along each element edge.
Definition: Qspectral_elements.h:783

References j, OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Member Data Documentation

◆ integral

template<unsigned NNODE_1D>
GaussLobattoLegendre< 2, NNODE_1D > oomph::QSpectralElement< 2, NNODE_1D >::integral
staticprivate

Assign the static integral.

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

/////////////////////////////////////////////////////////////////////// 2D QLegendreElements ///////////////////////////////////////////////////////////////////////


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