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

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

#include <Qspectral_elements.h>

+ Inheritance diagram for oomph::QSpectralElement< 3, 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<3,..> 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 &nplot)
 Output at nplot 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::BrickElementBase
 BrickElementBase ()
 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< 3, 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< 3, NNODE_1D >

General QSpectralElement class specialised to three spatial dimensions.

Constructor & Destructor Documentation

◆ QSpectralElement()

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

Constructor.

978  {
979  // There are NNODE_1D^3 nodes in this element
980  this->set_n_node(NNODE_1D * NNODE_1D * NNODE_1D);
981  Spectral_order.resize(3, NNODE_1D);
982  Nodal_spectral_order.resize(3, NNODE_1D);
983  // Set the elemental and nodal dimensions
984  this->set_dimension(3);
985  // Assign integral point
987  // Calculate the nodal positions for the shape functions
989  // OneDLegendreShapeParam::calculate_nodal_positions(NNODE_1D);
990  }
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< 3, NNODE_1D > integral
Assign the static integral.
Definition: Qspectral_elements.h:973
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< 3, NNODE_1D >::build_face_element ( const int face_index,
FaceElement face_element_pt 
)
virtual

Build the lower-dimensional FaceElement (an element of type QSpectralElement<2,NNODE_1D>). The face index takes six values corresponding to the six possible faces: -1 (Left) s[0] = -1.0 +1 (Right) s[0] = 1.0 -2 (Down) s[1] = -1.0 +2 (Up) s[1] = 1.0 -3 (Back) s[2] = -1.0 +3 (Front) s[2] = 1.0

Function to setup geometrical information for lower-dimensional FaceElements (which are of type QElement<3,NNODE_1D>).

Reimplemented from oomph::FiniteElement.

540  {
541  oomph_info << " WARNING UNTESTED CODE" << std::endl;
542 
543  // Set the nodal dimension from the "bulk"
544  face_element_pt->set_nodal_dimension(this->node_pt(0)->ndim());
545 
546  // Set the pointer to the orginal "bulk" element
547  face_element_pt->bulk_element_pt() = this;
548 
549 #ifdef OOMPH_HAS_MPI
550  // If the bulk element is halo then the face element must be too
551  // if (this->is_halo())
552  {
553  face_element_pt->set_halo(Non_halo_proc_ID);
554  }
555 #endif
556 
557  // Resize storage for the number of values originally stored
558  // at the face element's NNODE_1D*NNODE_1D nodes.
559  face_element_pt->nbulk_value_resize(NNODE_1D * NNODE_1D);
560 
561  // Set the face index in the element
562  // The faces are
563  // -3 : BACK (OLD: Bottom
564  // -2 : DOWN (OLD: Front
565  // -1 : LEFT (OLD: Left Side
566  // 1 : RIGHT (OLD: Right Side
567  // 2 : UP (OLD: Back
568  // 3 : FRONT (OLD: Top
569 
570  face_element_pt->face_index() = face_index;
571 
572  // Now set up the node pointers and the normal vectors
573  switch (face_index)
574  {
575  // BACK
576  //-----
577  case -3:
578  // Set the pointer to the bulk coordinate translation scheme
579  face_element_pt->face_to_bulk_coordinate_fct_pt() =
581 
582  // Set the pointer to the derivative mappings
583  face_element_pt->bulk_coordinate_derivatives_fct_pt() =
585 
586  // Copy nodes
587  for (unsigned i = 0; i < (NNODE_1D * NNODE_1D); i++)
588  {
589  face_element_pt->node_pt(i) = this->node_pt(i);
590  }
591  // Outer unit normal is negative of cross product of two in plane
592  // tangent vectors
593  face_element_pt->normal_sign() = -1;
594 
595  break;
596 
597  // FRONT
598  //------
599  case 3:
600 
601  // Set the pointer to the bulk coordinate translation scheme
602  face_element_pt->face_to_bulk_coordinate_fct_pt() =
604 
605  // Set the pointer to the derivative mappings
606  face_element_pt->bulk_coordinate_derivatives_fct_pt() =
608 
609  // Copy nodes
610  for (unsigned i = 0; i < (NNODE_1D * NNODE_1D); i++)
611  {
612  face_element_pt->node_pt(i) =
613  this->node_pt(i + (NNODE_1D * NNODE_1D) * (NNODE_1D - 1));
614  }
615  // Outer unit normal is cross product of two in plane
616  // tangent vectors
617  face_element_pt->normal_sign() = 1;
618 
619 
620  break;
621 
622  // DOWN:
623  //------
624  case -2:
625 
626  {
627  // Set the pointer to the bulk coordinate translation scheme
628  face_element_pt->face_to_bulk_coordinate_fct_pt() =
630 
631  // Set the pointer to the derivative mappings
632  face_element_pt->bulk_coordinate_derivatives_fct_pt() =
634 
635  // Copy nodes
636  unsigned count = 0;
637  for (unsigned i = 0; i < NNODE_1D; i++)
638  {
639  for (unsigned j = 0; j < NNODE_1D; j++)
640  {
641  face_element_pt->node_pt(count) =
642  this->node_pt(j + i * (NNODE_1D * NNODE_1D));
643  count++;
644  }
645  }
646 
647  // Outer unit normal is cross product of two in plane
648  // tangent vectors
649  face_element_pt->normal_sign() = 1;
650  }
651  break;
652 
653 
654  // UP:
655  //----
656  case 2:
657 
658  {
659  // Set the pointer to the bulk coordinate translation scheme
660  face_element_pt->face_to_bulk_coordinate_fct_pt() =
662 
663  // Set the pointer to the derivative mappings
664  face_element_pt->bulk_coordinate_derivatives_fct_pt() =
666 
667  // Copy nodes
668  unsigned count = 0;
669  for (unsigned i = 0; i < NNODE_1D; i++)
670  {
671  for (unsigned j = 0; j < NNODE_1D; j++)
672  {
673  face_element_pt->node_pt(count) = this->node_pt(
674  j + i * (NNODE_1D * NNODE_1D) + (NNODE_1D * (NNODE_1D - 1)));
675  count++;
676  }
677  }
678 
679  // Outer unit normal is negative of cross product of two in plane
680  // tangent vectors
681  face_element_pt->normal_sign() = -1;
682  }
683  break;
684 
685  // LEFT:
686  //------
687  case -1:
688 
689  {
690  // Set the pointer to the bulk coordinate translation scheme
691  face_element_pt->face_to_bulk_coordinate_fct_pt() =
693 
694  // Set the pointer to the derivative mappings
695  face_element_pt->bulk_coordinate_derivatives_fct_pt() =
697 
698  // Copy nodes
699  unsigned count = 0;
700  for (unsigned i = 0; i < NNODE_1D; i++)
701  {
702  for (unsigned j = 0; j < NNODE_1D; j++)
703  {
704  unsigned jj = j * NNODE_1D + i * (NNODE_1D * NNODE_1D);
705  face_element_pt->node_pt(count) = this->node_pt(jj);
706  count++;
707  }
708  }
709 
710  // Outer unit normal is negative of cross product of two in plane
711  // tangent vectors
712  face_element_pt->normal_sign() = -1;
713  }
714  break;
715 
716 
717  // RIGHT:
718  //-------
719  case 1:
720 
721  {
722  // Set the pointer to the bulk coordinate translation scheme
723  face_element_pt->face_to_bulk_coordinate_fct_pt() =
725 
726  // Set the pointer to the derivative mappings
727  face_element_pt->bulk_coordinate_derivatives_fct_pt() =
729 
730  // Copy nodes
731  unsigned count = 0;
732  for (unsigned i = 0; i < NNODE_1D; i++)
733  {
734  for (unsigned j = 0; j < NNODE_1D; j++)
735  {
736  unsigned jj =
737  j * NNODE_1D + i * (NNODE_1D * NNODE_1D) + (NNODE_1D - 1);
738  face_element_pt->node_pt(count) = this->node_pt(jj);
739  count++;
740  }
741  }
742 
743  // Outer unit normal is cross product of two in plane
744  // tangent vectors
745  face_element_pt->normal_sign() = 1;
746  }
747  break;
748 
749 
750  // Cover all other cases
751  default:
752  std::ostringstream error_message;
753  error_message
754  << "Face index should only take the values +/- 1, +/- 2 or +/- 3,"
755  << " not " << face_index << std::endl;
756  throw OomphLibError(error_message.str(),
759  } // end switch
760  }
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
unsigned ndim() const
Access function to # of Eulerian coordinates.
Definition: geom_objects.h:177
void faces2(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for the left and right faces, along which s2 is fixed.
Definition: Qelement_face_coordinate_translation_schemes.cc:249
void faces0(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for the back and front faces, along which s0 is fixed.
Definition: Qelement_face_coordinate_translation_schemes.cc:210
void faces1(const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction)
Function for the up and down faces, along which s1 is fixed.
Definition: Qelement_face_coordinate_translation_schemes.cc:230
void face4(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the up face (s1 = 1.0)
Definition: Qelement_face_coordinate_translation_schemes.cc:185
void face1(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the down face (s1 = -1.0)
Definition: Qelement_face_coordinate_translation_schemes.cc:161
void face2(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the back face (s2 = -1.0)
Definition: Qelement_face_coordinate_translation_schemes.cc:169
void face0(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the left face (s0 = -1.0)
Definition: Qelement_face_coordinate_translation_schemes.cc:153
void face5(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the front face (s2 = 1.0)
Definition: Qelement_face_coordinate_translation_schemes.cc:193
void face3(const Vector< double > &s, Vector< double > &s_bulk)
The translation scheme for the right face (s0 = 1.0)
Definition: Qelement_face_coordinate_translation_schemes.cc:177
OomphInfo oomph_info
Definition: oomph_definitions.cc:319
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::FaceElement::bulk_coordinate_derivatives_fct_pt(), oomph::FaceElement::bulk_element_pt(), oomph::QElement3FaceToBulkCoordinates::face0(), oomph::QElement3FaceToBulkCoordinates::face1(), oomph::QElement3FaceToBulkCoordinates::face2(), oomph::QElement3FaceToBulkCoordinates::face3(), oomph::QElement3FaceToBulkCoordinates::face4(), oomph::QElement3FaceToBulkCoordinates::face5(), oomph::FaceElement::face_index(), oomph::FaceElement::face_to_bulk_coordinate_fct_pt(), oomph::QElement3BulkCoordinateDerivatives::faces0(), oomph::QElement3BulkCoordinateDerivatives::faces1(), oomph::QElement3BulkCoordinateDerivatives::faces2(), i, j, oomph::FaceElement::nbulk_value_resize(), oomph::FiniteElement::node_pt(), oomph::FaceElement::normal_sign(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::oomph_info, and oomph::FiniteElement::set_nodal_dimension().

◆ d2shape_local()

template<unsigned NNODE_1D>
void oomph::QSpectralElement< 3, 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_2^2 \) d2psids(i,3) = \( \partial ^2 \psi_j / \partial s_0 \partial s_1 \) d2psids(i,4) = \( \partial ^2 \psi_j / \partial s_0 \partial s_2 \) d2psids(i,5) = \( \partial ^2 \psi_j / \partial s_1 \partial s_2 \)

Second derivatives of shape functions for specific QSpectralElement<3,..> 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_2^2 \) d2psids(i,3) = \( \partial ^2 \psi_j / \partial s_0 \partial s_1 \) d2psids(i,4) = \( \partial ^2 \psi_j / \partial s_0 \partial s_2 \) d2psids(i,5) = \( \partial ^2 \psi_j / \partial s_1 \partial s_2 \)

Reimplemented from oomph::FiniteElement.

1296  {
1297  std::ostringstream error_message;
1298  error_message
1299  << "\nd2shpe_local currently not implemented for this element\n";
1300  throw OomphLibError(
1301  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1302 
1303  /* //Call the shape functions and derivatives */
1304  /* OneDimensionalShape<NNODE_1D> psi1(s[0]); */
1305  /* OneDimensionalDShape<NNODE_1D> dpsi1ds(s[0]); */
1306  /* OneDimensionalD2Shape<NNODE_1D> d2psi1ds(s[0]); */
1307 
1308  /* //Loop over shape functions in element and assign to psi */
1309  /* for(unsigned l=0;l<NNODE_1D;l++) */
1310  /* { */
1311  /* psi[l] = psi1[l]; */
1312  /* dpsids[l][0] = dpsi1ds[l]; */
1313  /* d2psids[l][0] = d2psi1ds[l]; */
1314  /* } */
1315  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ dshape_local()

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

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

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

Reimplemented from oomph::FiniteElement.

1248  {
1249  // Call the shape functions and derivatives
1250  // OneDLegendreShapeParam psi1(NNODE_1D,s[0]), psi2(NNODE_1D,s[1]),
1251  // psi3(NNODE_1D,s[2]);
1252  // OneDLegendreDShapeParam dpsi1ds(NNODE_1D,s[0]), dpsi2ds(NNODE_1D,s[1]),
1253  // dpsi3ds(NNODE_1D,s[2]);
1254  OneDimensionalLegendreShape<NNODE_1D> psi1(s[0]), psi2(s[1]), psi3(s[2]);
1255  OneDimensionalLegendreDShape<NNODE_1D> dpsi1ds(s[0]), dpsi2ds(s[1]),
1256  dpsi3ds(s[2]);
1257 
1258 
1259  // Index for the shape functions
1260  unsigned index = 0;
1261 
1262  // Loop over shape functions in element
1263  for (unsigned i = 0; i < NNODE_1D; i++)
1264  {
1265  for (unsigned j = 0; j < NNODE_1D; j++)
1266  {
1267  for (unsigned k = 0; k < NNODE_1D; k++)
1268  {
1269  // Assign the values
1270  dpsids(index, 0) = psi3[i] * psi2[j] * dpsi1ds[k];
1271  dpsids(index, 1) = psi3[i] * dpsi2ds[j] * psi1[k];
1272  dpsids(index, 2) = dpsi3ds[i] * psi2[j] * psi1[k];
1273  psi[index] = psi3[i] * psi2[j] * psi1[k];
1274  // Increment the index
1275  ++index;
1276  }
1277  }
1278  }
1279  }
RealScalar s
Definition: level1_cplx_impl.h:130
char char char int int * k
Definition: level2_impl.h:374

References i, j, k, and s.

◆ get_s_plot()

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

1147  {
1148  if (nplot > 1)
1149  {
1150  unsigned i01 = i % (nplot * nplot);
1151  unsigned i0 = i01 % nplot;
1152  unsigned i1 = (i01 - i0) / nplot;
1153  unsigned i2 = (i - i01) / (nplot * nplot);
1154 
1155  s[0] = -1.0 + 2.0 * double(i0) / double(nplot - 1);
1156  s[1] = -1.0 + 2.0 * double(i1) / double(nplot - 1);
1157  s[2] = -1.0 + 2.0 * double(i2) / double(nplot - 1);
1158  if (use_equally_spaced_interior_sample_points)
1159  {
1160  double range = 2.0;
1161  double dx_new = range / double(nplot);
1162  double range_new = double(nplot - 1) * dx_new;
1163  s[0] = -1.0 + 0.5 * dx_new + range_new * (1.0 + s[0]) / range;
1164  s[1] = -1.0 + 0.5 * dx_new + range_new * (1.0 + s[1]) / range;
1165  s[2] = -1.0 + 0.5 * dx_new + range_new * (1.0 + s[2]) / range;
1166  }
1167  }
1168  else
1169  {
1170  s[0] = 0.0;
1171  s[1] = 0.0;
1172  s[2] = 0.0;
1173  }
1174  }

References i, and s.

◆ invert_jacobian_mapping()

template<unsigned NNODE_1D>
double oomph::QSpectralElement< 3, 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 3D version.

Reimplemented from oomph::FiniteElement.

1112  {
1113  return FiniteElement::invert_jacobian<3>(jacobian, inverse_jacobian);
1114  }

◆ local_coordinate_of_node()

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

1057  {
1058  s.resize(3);
1059  unsigned n0 = n % NNODE_1D;
1060  unsigned n1 = unsigned(double(n) / double(NNODE_1D)) % NNODE_1D;
1061  unsigned n2 = unsigned(double(n) / double(NNODE_1D * NNODE_1D));
1065  }
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< 3, 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.

1069  {
1070  this->local_coordinate_of_node(n, s_fraction);
1071  // Resize
1072  s_fraction[0] = 0.5 * (s_fraction[0] + 1.0);
1073  s_fraction[1] = 0.5 * (s_fraction[1] + 1.0);
1074  s_fraction[2] = 0.5 * (s_fraction[2] + 1.0);
1075  }
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:1056

References n.

◆ local_one_d_fraction_of_node()

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

1079  {
1080  return 0.5 *
1082  }

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

◆ nnode_1d()

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

Number of nodes along each element edge.

Reimplemented from oomph::FiniteElement.

1118  {
1119  return NNODE_1D;
1120  }

◆ nplot_points()

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

1190  {
1191  unsigned DIM = 3;
1192  unsigned np = 1;
1193  for (unsigned i = 0; i < DIM; i++)
1194  {
1195  np *= nplot;
1196  }
1197  return np;
1198  }
#define DIM
Definition: linearised_navier_stokes_elements.h:44

References DIM, and i.

◆ nvertex_node()

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

Number of vertex nodes in the element.

Implements oomph::BrickElementBase.

1006  {
1007  return 8;
1008  }

◆ output() [1/4]

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

C-style output.

Reimplemented from oomph::FiniteElement.

1124  {
1125  FiniteElement::output(file_pt);
1126  }
virtual void output(std::ostream &outfile)
Definition: elements.h:3050

References oomph::FiniteElement::output().

◆ output() [2/4]

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

C_style output at n_plot points.

Reimplemented from oomph::FiniteElement.

1130  {
1131  FiniteElement::output(file_pt, n_plot);
1132  }

References oomph::FiniteElement::output().

◆ output() [3/4]

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

Output.

The output function for general 1D QSpectralElements.

Reimplemented from oomph::FiniteElement.

456  {
457  // Tecplot header info
458  outfile << "ZONE I=" << NNODE_1D << ", J=" << NNODE_1D << ", K=" << NNODE_1D
459  << std::endl;
460 
461  // Find the dimension of the nodes
462  unsigned n_dim = this->nodal_dimension();
463 
464  // Loop over element nodes
465  for (unsigned l3 = 0; l3 < NNODE_1D; l3++)
466  {
467  for (unsigned l2 = 0; l2 < NNODE_1D; l2++)
468  {
469  for (unsigned l1 = 0; l1 < NNODE_1D; l1++)
470  {
471  unsigned l = l3 * NNODE_1D * NNODE_1D + l2 * NNODE_1D + l1;
472 
473  // Loop over the dimensions and output the position
474  for (unsigned i = 0; i < n_dim; i++)
475  {
476  outfile << this->node_pt(l)->x(i) << " ";
477  }
478  // Find out how many data values at the node
479  unsigned initial_nvalue = this->node_pt(l)->nvalue();
480  // Loop over the data and output whether pinned or not
481  for (unsigned i = 0; i < initial_nvalue; i++)
482  {
483  outfile << this->node_pt(l)->is_pinned(i) << " ";
484  }
485  outfile << std::endl;
486  }
487  }
488  }
489  outfile << std::endl;
490  }
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< 3, NNODE_1D >::output ( std::ostream &  outfile,
const unsigned nplot 
)
virtual

Output at nplot points.

The output function for n_plot points in each coordinate direction.

Reimplemented from oomph::FiniteElement.

498  {
499  // Local variables
500  Vector<double> s(3);
501 
502  // Tecplot header info
503  outfile << "ZONE I=" << n_plot << ", J=" << n_plot << ", K=" << n_plot
504  << std::endl;
505 
506  // Find the dimension of the first node
507  unsigned n_dim = this->nodal_dimension();
508 
509  // Loop over element nodes
510  for (unsigned l3 = 0; l3 < n_plot; l3++)
511  {
512  s[2] = -1.0 + l3 * 2.0 / (n_plot - 1);
513  for (unsigned l2 = 0; l2 < n_plot; l2++)
514  {
515  s[1] = -1.0 + l2 * 2.0 / (n_plot - 1);
516  for (unsigned l1 = 0; l1 < n_plot; l1++)
517  {
518  s[0] = -1.0 + l1 * 2.0 / (n_plot - 1);
519 
520  // Output the coordinates
521  for (unsigned i = 0; i < n_dim; i++)
522  {
523  outfile << this->interpolated_x(s, i) << " ";
524  }
525  outfile << std::endl;
526  }
527  }
528  }
529  outfile << std::endl;
530  }
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< 3, NNODE_1D >::s_max ( ) const
inlinevirtual

Max. value of local coordinate.

Reimplemented from oomph::FiniteElement.

1000  {
1001  return 1.0;
1002  }

◆ s_min()

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

Min. value of local coordinate.

Reimplemented from oomph::FiniteElement.

994  {
995  return -1.0;
996  }

◆ shape()

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

Implements oomph::FiniteElement.

1220  {
1221  // Call the OneDimensional Shape functions
1222  OneDimensionalLegendreShape<NNODE_1D> psi1(s[0]), psi2(s[1]), psi3(s[2]);
1223  // OneDLegendreShapeParam psi1(NNODE_1D,s[0]), psi2(NNODE_1D,s[1]),
1224  // psi3(NNODE_1D,s[2]);
1225 
1226  // Now let's loop over the nodal points in the element
1227  // and copy the values back in
1228  for (unsigned i = 0; i < NNODE_1D; i++)
1229  {
1230  for (unsigned j = 0; j < NNODE_1D; j++)
1231  {
1232  for (unsigned k = 0; k < NNODE_1D; k++)
1233  {
1234  psi(NNODE_1D * NNODE_1D * i + NNODE_1D * j + k) =
1235  psi3[i] * psi2[j] * psi1[k];
1236  }
1237  }
1238  }
1239  }

References i, j, k, and s.

◆ tecplot_zone_string()

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

1180  {
1181  std::ostringstream header;
1182  header << "ZONE I=" << nplot << ", J=" << nplot << ", K=" << nplot
1183  << "\n";
1184  return header.str();
1185  }

◆ vertex_node_pt()

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

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

Implements oomph::BrickElementBase.

1012  {
1013  unsigned n_node_1d = nnode_1d();
1014  Node* nod_pt;
1015  switch (j)
1016  {
1017  case 0:
1018  nod_pt = this->node_pt(0);
1019  break;
1020  case 1:
1021  nod_pt = this->node_pt(n_node_1d - 1);
1022  break;
1023  case 2:
1024  nod_pt = this->node_pt(n_node_1d * (n_node_1d - 1));
1025  break;
1026  case 3:
1027  nod_pt = this->node_pt(n_node_1d * n_node_1d - 1);
1028  break;
1029  case 4:
1030  nod_pt = this->node_pt(n_node_1d * n_node_1d * (n_node_1d - 1));
1031  break;
1032  case 5:
1033  nod_pt = this->node_pt(n_node_1d * n_node_1d * (n_node_1d - 1) +
1034  (n_node_1d - 1));
1035  break;
1036  case 6:
1037  nod_pt = this->node_pt(n_node_1d * n_node_1d * n_node_1d - n_node_1d);
1038  break;
1039  case 7:
1040  nod_pt = this->node_pt(n_node_1d * n_node_1d * n_node_1d - 1);
1041  break;
1042  default:
1043  std::ostringstream error_message;
1044  error_message << "Vertex node number is " << j
1045  << " but must be from 0 to 7\n";
1046 
1047  throw OomphLibError(error_message.str(),
1050  }
1051  return nod_pt;
1052  }
unsigned nnode_1d() const
Number of nodes along each element edge.
Definition: Qspectral_elements.h:1117

References j, OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Member Data Documentation

◆ integral

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

Assign the static integral.

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

/////////////////////////////////////////////////////////////////////// 3D QLegendreElements ///////////////////////////////////////////////////////////////////////


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