oomph::TTaylorHoodElement< DIM > Class Template Reference

#include <Tnavier_stokes_elements.h>

+ Inheritance diagram for oomph::TTaylorHoodElement< DIM >:

Public Member Functions

 TTaylorHoodElement ()
 Constructor, no internal data points. More...
 
 TTaylorHoodElement (const TTaylorHoodElement< DIM > &dummy)=delete
 Broken copy constructor. More...
 
virtual unsigned required_nvalue (const unsigned &n) const
 Broken assignment operator. More...
 
void pshape_nst (const Vector< double > &s, Shape &psi) const
 Test whether the pressure dof p_dof hanging or not? More...
 
void pshape_nst (const Vector< double > &s, Shape &psi, Shape &test) const
 Pressure shape and test functions at local coordinte s. More...
 
unsigned p_index_nst ()
 Which nodal value represents the pressure? More...
 
int p_local_eqn (const unsigned &n) const
 Pointer to n_p-th pressure node. More...
 
double p_nst (const unsigned &n_p) const
 
double p_nst (const unsigned &t, const unsigned &n_p) const
 
int p_nodal_index_nst () const
 Set the value at which the pressure is stored in the nodes. More...
 
unsigned npres_nst () const
 Return number of pressure values. More...
 
void fix_pressure (const unsigned &p_dof, const double &p_value)
 Pin p_dof-th pressure dof and set it to value specified by p_value. More...
 
void build_fp_press_adv_diff_robin_bc_element (const unsigned &face_index)
 
void identify_load_data (std::set< std::pair< Data *, unsigned >> &paired_load_data)
 
void identify_pressure_data (std::set< std::pair< Data *, unsigned >> &paired_pressure_data)
 
void output (std::ostream &outfile)
 Redirect output to NavierStokesEquations output. More...
 
void output (std::ostream &outfile, const unsigned &nplot)
 Redirect output to NavierStokesEquations output. More...
 
void output (FILE *file_pt)
 Redirect output to NavierStokesEquations output. More...
 
void output (FILE *file_pt, const unsigned &n_plot)
 Redirect output to NavierStokesEquations output. More...
 
unsigned nrecovery_order ()
 
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...
 
unsigned num_Z2_flux_terms ()
 Number of 'flux' terms for Z2 error estimation. More...
 
void get_Z2_flux (const Vector< double > &s, Vector< double > &flux)
 
unsigned ndof_types () const
 
void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 
unsigned npres_nst () const
 
unsigned npres_nst () const
 
void pshape_nst (const Vector< double > &s, Shape &psi) const
 
void pshape_nst (const Vector< double > &s, Shape &psi) const
 
- Public Member Functions inherited from oomph::TElement< DIM, 3 >
const unsigned Node_on_face [3][2]
 Assign the nodal translation schemes. More...
 
const unsigned Node_on_face [3][3]
 
const unsigned Node_on_face [3][4]
 
const unsigned Node_on_face [4][3]
 Assign the nodal translation scheme for linear elements. More...
 
const unsigned Node_on_face [4][6]
 Assign the nodal translation scheme for quadratic elements. More...
 
- Public Member Functions inherited from oomph::NavierStokesEquations< DIM >
 NavierStokesEquations ()
 
const doublere () const
 Reynolds number. More...
 
const doublere_st () const
 Product of Reynolds and Strouhal number (=Womersley number) More...
 
double *& re_pt ()
 Pointer to Reynolds number. More...
 
double *& re_st_pt ()
 Pointer to product of Reynolds and Strouhal number (=Womersley number) More...
 
const doubleviscosity_ratio () const
 
double *& viscosity_ratio_pt ()
 Pointer to Viscosity Ratio. More...
 
const doubledensity_ratio () const
 
double *& density_ratio_pt ()
 Pointer to Density ratio. More...
 
const doublere_invfr () const
 Global inverse Froude number. More...
 
double *& re_invfr_pt ()
 Pointer to global inverse Froude number. More...
 
const Vector< double > & g () const
 Vector of gravitational components. More...
 
Vector< double > *& g_pt ()
 Pointer to Vector of gravitational components. More...
 
NavierStokesBodyForceFctPtbody_force_fct_pt ()
 Access function for the body-force pointer. More...
 
NavierStokesBodyForceFctPt body_force_fct_pt () const
 Access function for the body-force pointer. Const version. More...
 
NavierStokesSourceFctPtsource_fct_pt ()
 Access function for the source-function pointer. More...
 
NavierStokesSourceFctPt source_fct_pt () const
 Access function for the source-function pointer. Const version. More...
 
NavierStokesPressureAdvDiffSourceFctPtsource_fct_for_pressure_adv_diff ()
 
NavierStokesPressureAdvDiffSourceFctPt source_fct_for_pressure_adv_diff () const
 
intpinned_fp_pressure_eqn ()
 
double u_nst (const unsigned &n, const unsigned &i) const
 
double u_nst (const unsigned &t, const unsigned &n, const unsigned &i) const
 
virtual unsigned u_index_nst (const unsigned &i) const
 
unsigned n_u_nst () const
 
double du_dt_nst (const unsigned &n, const unsigned &i) const
 
void disable_ALE ()
 
void enable_ALE ()
 
double pressure_integral () const
 Integral of pressure over element. More...
 
double dissipation () const
 Return integral of dissipation over element. More...
 
double dissipation (const Vector< double > &s) const
 Return dissipation at local coordinate s. More...
 
void get_vorticity (const Vector< double > &s, Vector< double > &vorticity) const
 Compute the vorticity vector at local coordinate s. More...
 
void get_vorticity (const Vector< double > &s, double &vorticity) const
 Compute the scalar vorticity at local coordinate s (2D) More...
 
double kin_energy () const
 Get integral of kinetic energy over element. More...
 
double d_kin_energy_dt () const
 Get integral of time derivative of kinetic energy over element. More...
 
void strain_rate (const Vector< double > &s, DenseMatrix< double > &strain_rate) const
 Strain-rate tensor: 1/2 (du_i/dx_j + du_j/dx_i) More...
 
void get_traction (const Vector< double > &s, const Vector< double > &N, Vector< double > &traction)
 
void get_traction (const Vector< double > &s, const Vector< double > &N, Vector< double > &traction_p, Vector< double > &traction_visc_n, Vector< double > &traction_visc_t)
 
void get_load (const Vector< double > &s, const Vector< double > &N, Vector< double > &load)
 
void get_pressure_and_velocity_mass_matrix_diagonal (Vector< double > &press_mass_diag, Vector< double > &veloc_mass_diag, const unsigned &which_one=0)
 
unsigned nscalar_paraview () const
 
void scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
 
void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
 
std::string scalar_name_paraview (const unsigned &i) const
 
void full_output (std::ostream &outfile)
 
void full_output (std::ostream &outfile, const unsigned &nplot)
 
void output_veloc (std::ostream &outfile, const unsigned &nplot, const unsigned &t)
 
void output_vorticity (std::ostream &outfile, const unsigned &nplot)
 
void output_fct (std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 
void output_fct (std::ostream &outfile, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
 
void compute_norm (double &norm)
 Compute norm of solution: square of the L2 norm of the velocities. More...
 
void compute_norm (Vector< double > &norm)
 Compute the vector norm of the FEM solution. More...
 
void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 
void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 
void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Compute the element's residual Vector. More...
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
void fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
void fill_in_contribution_to_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 Compute the element's residual Vector. More...
 
void fill_in_contribution_to_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
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)
 
void fill_in_pressure_advection_diffusion_residuals (Vector< double > &residuals)
 
void fill_in_pressure_advection_diffusion_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
void pin_all_non_pressure_dofs (std::map< Data *, std::vector< int >> &eqn_number_backup)
 Pin all non-pressure dofs and backup eqn numbers. More...
 
void output_pressure_advection_diffusion_robin_elements (std::ostream &outfile)
 
void delete_pressure_advection_diffusion_robin_elements ()
 
virtual void get_dresidual_dnodal_coordinates (RankThreeTensor< double > &dresidual_dnodal_coordinates)
 
void interpolated_u_nst (const Vector< double > &s, Vector< double > &veloc) const
 Compute vector of FE interpolated velocity u at local coordinate s. More...
 
double interpolated_u_nst (const Vector< double > &s, const unsigned &i) const
 Return FE interpolated velocity u[i] at local coordinate s. More...
 
double interpolated_u_nst (const unsigned &t, const Vector< double > &s, const unsigned &i) const
 
virtual void dinterpolated_u_nst_ddata (const Vector< double > &s, const unsigned &i, Vector< double > &du_ddata, Vector< unsigned > &global_eqn_number)
 
virtual double interpolated_p_nst (const Vector< double > &s) const
 Return FE interpolated pressure at local coordinate s. More...
 
double interpolated_p_nst (const unsigned &t, const Vector< double > &s) const
 Return FE interpolated pressure at local coordinate s at time level t. More...
 
double interpolated_dudx_nst (const Vector< double > &s, const unsigned &i, const unsigned &j) const
 
void point_output_data (const Vector< double > &s, Vector< double > &data)
 
void get_vorticity (const Vector< double > &s, Vector< double > &vorticity) const
 
void get_vorticity (const Vector< double > &s, double &vorticity) const
 
void get_vorticity (const Vector< double > &s, Vector< double > &vorticity) const
 Compute 3D vorticity vector at local coordinate s. More...
 
- Public Member Functions inherited from oomph::FSIFluidElement
 FSIFluidElement ()
 Constructor. More...
 
 FSIFluidElement (const FSIFluidElement &)=delete
 Broken copy constructor. More...
 
void operator= (const FSIFluidElement &)=delete
 Broken assignment operator. More...
 
- Public Member Functions inherited from oomph::FiniteElement
void set_dimension (const unsigned &dim)
 
void set_nodal_dimension (const unsigned &nodal_dim)
 
void set_nnodal_position_type (const unsigned &nposition_type)
 Set the number of types required to interpolate the coordinate. More...
 
void set_n_node (const unsigned &n)
 
int nodal_local_eqn (const unsigned &n, const unsigned &i) const
 
double dJ_eulerian_at_knot (const unsigned &ipt, Shape &psi, DenseMatrix< double > &djacobian_dX) const
 
 FiniteElement ()
 Constructor. More...
 
virtual ~FiniteElement ()
 
 FiniteElement (const FiniteElement &)=delete
 Broken copy constructor. More...
 
virtual bool local_coord_is_valid (const Vector< double > &s)
 Broken assignment operator. More...
 
virtual void move_local_coord_back_into_element (Vector< double > &s) const
 
void get_centre_of_gravity_and_max_radius_in_terms_of_zeta (Vector< double > &cog, double &max_radius) const
 
virtual void local_coordinate_of_node (const unsigned &j, Vector< double > &s) const
 
virtual void local_fraction_of_node (const unsigned &j, Vector< double > &s_fraction)
 
virtual double local_one_d_fraction_of_node (const unsigned &n1d, const unsigned &i)
 
virtual void set_macro_elem_pt (MacroElement *macro_elem_pt)
 
MacroElementmacro_elem_pt ()
 Access function to pointer to macro element. More...
 
void get_x (const Vector< double > &s, Vector< double > &x) const
 
void get_x (const unsigned &t, const Vector< double > &s, Vector< double > &x)
 
virtual void get_x_from_macro_element (const Vector< double > &s, Vector< double > &x) const
 
virtual void get_x_from_macro_element (const unsigned &t, const Vector< double > &s, Vector< double > &x)
 
virtual void set_integration_scheme (Integral *const &integral_pt)
 Set the spatial integration scheme. More...
 
Integral *const & integral_pt () const
 Return the pointer to the integration scheme (const version) More...
 
virtual void shape (const Vector< double > &s, Shape &psi) const =0
 
virtual void shape_at_knot (const unsigned &ipt, Shape &psi) const
 
virtual void dshape_local (const Vector< double > &s, Shape &psi, DShape &dpsids) const
 
virtual void dshape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids) const
 
virtual void d2shape_local (const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
 
virtual void d2shape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids, DShape &d2psids) const
 
virtual double J_eulerian (const Vector< double > &s) const
 
virtual double J_eulerian_at_knot (const unsigned &ipt) const
 
void check_J_eulerian_at_knots (bool &passed) const
 
void check_jacobian (const double &jacobian) const
 
double dshape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsi, DenseMatrix< double > &djacobian_dX, RankFourTensor< double > &d_dpsidx_dX) const
 
double d2shape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual double d2shape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual void assign_nodal_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void describe_nodal_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void assign_all_generic_local_eqn_numbers (const bool &store_local_dof_pt)
 
Node *& node_pt (const unsigned &n)
 Return a pointer to the local node n. More...
 
Node *const & node_pt (const unsigned &n) const
 Return a pointer to the local node n (const version) More...
 
unsigned nnode () const
 Return the number of nodes. More...
 
virtual unsigned nnode_1d () const
 
double raw_nodal_position (const unsigned &n, const unsigned &i) const
 
double raw_nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position (const unsigned &n, const unsigned &i) const
 
double nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double dnodal_position_dt (const unsigned &n, const unsigned &i) const
 Return the i-th component of nodal velocity: dx/dt at local node n. More...
 
double dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
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 ElementGeometry::ElementGeometry element_geometry () const
 Return the geometry type of the element (either Q or T usually). More...
 
virtual double interpolated_x (const Vector< double > &s, const unsigned &i) const
 Return FE interpolated coordinate x[i] at local coordinate s. More...
 
virtual double interpolated_x (const unsigned &t, const Vector< double > &s, const unsigned &i) const
 
virtual void interpolated_x (const Vector< double > &s, Vector< double > &x) const
 Return FE interpolated position x[] at local coordinate s as Vector. More...
 
virtual void interpolated_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) const
 
virtual double interpolated_dxdt (const Vector< double > &s, const unsigned &i, const unsigned &t)
 
virtual void interpolated_dxdt (const Vector< double > &s, const unsigned &t, Vector< double > &dxdt)
 
unsigned ngeom_data () const
 
Datageom_data_pt (const unsigned &j)
 
void position (const Vector< double > &zeta, Vector< double > &r) const
 
void position (const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const
 
void dposition_dt (const Vector< double > &zeta, const unsigned &t, Vector< double > &drdt)
 
virtual double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 
void interpolated_zeta (const Vector< double > &s, Vector< double > &zeta) const
 
void locate_zeta (const Vector< double > &zeta, GeomObject *&geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
 
virtual void node_update ()
 
virtual void identify_field_data_for_interactions (std::set< std::pair< Data *, unsigned >> &paired_field_data)
 
virtual void identify_geometric_data (std::set< Data * > &geometric_data_pt)
 
virtual double s_min () const
 Min value of local coordinate. More...
 
virtual double s_max () const
 Max. value of local coordinate. More...
 
double size () const
 
virtual double compute_physical_size () const
 
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 void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) 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, const double &time, const SolutionFunctorBase &exact_soln) const
 Output a time-dependent exact solution over the element. More...
 
virtual void get_s_plot (const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
 
virtual std::string tecplot_zone_string (const unsigned &nplot) const
 
virtual void write_tecplot_zone_footer (std::ostream &outfile, const unsigned &nplot) const
 
virtual void write_tecplot_zone_footer (FILE *file_pt, const unsigned &nplot) const
 
virtual unsigned nplot_points (const unsigned &nplot) const
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, 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, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_abs_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error)
 
void integrate_fct (FiniteElement::SteadyExactSolutionFctPt integrand_fct_pt, Vector< double > &integral)
 Integrate Vector-valued function over element. More...
 
void integrate_fct (FiniteElement::UnsteadyExactSolutionFctPt integrand_fct_pt, const double &time, Vector< double > &integral)
 Integrate Vector-valued time-dep function over element. More...
 
virtual void build_face_element (const int &face_index, FaceElement *face_element_pt)
 
virtual unsigned self_test ()
 
virtual unsigned get_bulk_node_number (const int &face_index, const unsigned &i) const
 
virtual int face_outer_unit_normal_sign (const int &face_index) const
 Get the sign of the outer unit normal on the face given by face_index. More...
 
virtual unsigned nnode_on_face () const
 
void face_node_number_error_check (const unsigned &i) const
 Range check for face node numbers. More...
 
virtual CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt (const int &face_index) const
 
virtual BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt (const int &face_index) const
 
- Public Member Functions inherited from oomph::GeneralisedElement
 GeneralisedElement ()
 Constructor: Initialise all pointers and all values to zero. More...
 
virtual ~GeneralisedElement ()
 Virtual destructor to clean up any memory allocated by the object. More...
 
 GeneralisedElement (const GeneralisedElement &)=delete
 Broken copy constructor. More...
 
void operator= (const GeneralisedElement &)=delete
 Broken assignment operator. More...
 
Data *& internal_data_pt (const unsigned &i)
 Return a pointer to i-th internal data object. More...
 
Data *const & internal_data_pt (const unsigned &i) const
 Return a pointer to i-th internal data object (const version) More...
 
Data *& external_data_pt (const unsigned &i)
 Return a pointer to i-th external data object. More...
 
Data *const & external_data_pt (const unsigned &i) const
 Return a pointer to i-th external data object (const version) More...
 
unsigned long eqn_number (const unsigned &ieqn_local) const
 
int local_eqn_number (const unsigned long &ieqn_global) const
 
unsigned add_external_data (Data *const &data_pt, const bool &fd=true)
 
bool external_data_fd (const unsigned &i) const
 
void exclude_external_data_fd (const unsigned &i)
 
void include_external_data_fd (const unsigned &i)
 
void flush_external_data ()
 Flush all external data. More...
 
void flush_external_data (Data *const &data_pt)
 Flush the object addressed by data_pt from the external data array. More...
 
unsigned ninternal_data () const
 Return the number of internal data objects. More...
 
unsigned nexternal_data () const
 Return the number of external data objects. More...
 
unsigned ndof () const
 Return the number of equations/dofs in the element. More...
 
void dof_vector (const unsigned &t, Vector< double > &dof)
 Return the vector of dof values at time level t. More...
 
void dof_pt_vector (Vector< double * > &dof_pt)
 Return the vector of pointers to dof values. More...
 
void set_internal_data_time_stepper (const unsigned &i, TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void assign_internal_eqn_numbers (unsigned long &global_number, Vector< double * > &Dof_pt)
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
void add_internal_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt)
 
virtual void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void complete_setup_of_dependencies ()
 
virtual void get_residuals (Vector< double > &residuals)
 
virtual void get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
virtual void get_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void get_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void get_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void get_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void get_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void get_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void get_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void get_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
- 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::TemplateFreeNavierStokesEquationsBase
 TemplateFreeNavierStokesEquationsBase ()
 Constructor (empty) More...
 
virtual ~TemplateFreeNavierStokesEquationsBase ()
 Virtual destructor (empty) More...
 
- Public Member Functions inherited from oomph::NavierStokesElementWithDiagonalMassMatrices
 NavierStokesElementWithDiagonalMassMatrices ()
 Empty constructor. More...
 
virtual ~NavierStokesElementWithDiagonalMassMatrices ()
 Virtual destructor. More...
 
 NavierStokesElementWithDiagonalMassMatrices (const NavierStokesElementWithDiagonalMassMatrices &)=delete
 Broken copy constructor. More...
 
void operator= (const NavierStokesElementWithDiagonalMassMatrices &)=delete
 Broken assignment operator. More...
 
- Public Member Functions inherited from oomph::ElementWithZ2ErrorEstimator
 ElementWithZ2ErrorEstimator ()
 Default empty constructor. More...
 
 ElementWithZ2ErrorEstimator (const ElementWithZ2ErrorEstimator &)=delete
 Broken copy constructor. More...
 
void operator= (const ElementWithZ2ErrorEstimator &)=delete
 Broken assignment operator. More...
 
virtual unsigned ncompound_fluxes ()
 
virtual void compute_exact_Z2_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_flux_pt, double &error, double &norm)
 
virtual void get_Z2_compound_flux_indices (Vector< unsigned > &flux_index)
 
virtual double geometric_jacobian (const Vector< double > &x)
 

Static Public Attributes

static const unsigned TEMPLATE_PARAMETER_DIM = DIM
 Publicly exposed template parameter. More...
 
static const unsigned TEMPLATE_PARAMETER_NNODE_1D = 3
 Publicly exposed template parameter. More...
 
- Static Public Attributes inherited from oomph::NavierStokesEquations< DIM >
static Vector< doubleGamma
 Vector to decide whether the stress-divergence form is used or not. More...
 
- 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

double dshape_and_dtest_eulerian_nst (const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
 
double dshape_and_dtest_eulerian_at_knot_nst (const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
 
double dshape_and_dtest_eulerian_at_knot_nst (const unsigned &ipt, Shape &psi, DShape &dpsidx, RankFourTensor< double > &d_dpsidx_dX, Shape &test, DShape &dtestdx, RankFourTensor< double > &d_dtestdx_dX, DenseMatrix< double > &djacobian_dX) const
 
virtual double dpshape_and_dptest_eulerian_nst (const Vector< double > &s, Shape &ppsi, DShape &dppsidx, Shape &ptest, DShape &dptestdx) const
 
void unpin_all_nodal_pressure_dofs ()
 Unpin all pressure dofs. More...
 
void pin_all_nodal_pressure_dofs ()
 Pin all nodal pressure dofs. More...
 
void unpin_proper_nodal_pressure_dofs ()
 Unpin the proper nodal pressure dofs. More...
 
const unsigned Pconv [3]
 
const unsigned Pconv [4]
 
double dpshape_and_dptest_eulerian_nst (const Vector< double > &s, Shape &ppsi, DShape &dppsidx, Shape &ptest, DShape &dptestdx) const
 
double dpshape_and_dptest_eulerian_nst (const Vector< double > &s, Shape &ppsi, DShape &dppsidx, Shape &ptest, DShape &dptestdx) const
 
double dshape_and_dtest_eulerian_at_knot_nst (const unsigned &ipt, Shape &psi, DShape &dpsidx, RankFourTensor< double > &d_dpsidx_dX, Shape &test, DShape &dtestdx, RankFourTensor< double > &d_dtestdx_dX, DenseMatrix< double > &djacobian_dX) const
 
double dshape_and_dtest_eulerian_at_knot_nst (const unsigned &ipt, Shape &psi, DShape &dpsidx, RankFourTensor< double > &d_dpsidx_dX, Shape &test, DShape &dtestdx, RankFourTensor< double > &d_dtestdx_dX, DenseMatrix< double > &djacobian_dX) const
 
- Protected Member Functions inherited from oomph::NavierStokesEquations< DIM >
virtual void get_body_force_nst (const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &result)
 
virtual void get_body_force_gradient_nst (const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, DenseMatrix< double > &d_body_force_dx)
 
virtual double get_source_nst (const double &time, const unsigned &ipt, const Vector< double > &x)
 
virtual void get_source_gradient_nst (const double &time, const unsigned &ipt, const Vector< double > &x, Vector< double > &gradient)
 
virtual void fill_in_generic_residual_contribution_nst (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
 
virtual void fill_in_generic_pressure_advection_diffusion_contribution_nst (Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
 
virtual void fill_in_generic_dresidual_contribution_nst (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam, unsigned flag)
 
void fill_in_contribution_to_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
- Protected Member Functions inherited from oomph::FiniteElement
virtual void assemble_local_to_eulerian_jacobian (const DShape &dpsids, DenseMatrix< double > &jacobian) const
 
virtual void assemble_local_to_eulerian_jacobian2 (const DShape &d2psids, DenseMatrix< double > &jacobian2) const
 
virtual void assemble_eulerian_base_vectors (const DShape &dpsids, DenseMatrix< double > &interpolated_G) const
 
template<unsigned DIM>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual double invert_jacobian_mapping (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual double local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
double local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &inverse_jacobian) const
 
virtual double local_to_eulerian_mapping_diagonal (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual void dJ_eulerian_dnodal_coordinates (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<unsigned DIM>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
virtual void d_dshape_eulerian_dnodal_coordinates (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<unsigned DIM>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
virtual void transform_derivatives (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const
 
void transform_derivatives_diagonal (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const
 
virtual void transform_second_derivatives (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<unsigned DIM>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<unsigned DIM>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
virtual void fill_in_jacobian_from_nodal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
void fill_in_jacobian_from_nodal_by_fd (DenseMatrix< double > &jacobian)
 
virtual void update_before_nodal_fd ()
 
virtual void reset_after_nodal_fd ()
 
virtual void update_in_nodal_fd (const unsigned &i)
 
virtual void reset_in_nodal_fd (const unsigned &i)
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Zero-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 One-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Two-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
- Protected Member Functions inherited from oomph::GeneralisedElement
unsigned add_internal_data (Data *const &data_pt, const bool &fd=true)
 
bool internal_data_fd (const unsigned &i) const
 
void exclude_internal_data_fd (const unsigned &i)
 
void include_internal_data_fd (const unsigned &i)
 
void clear_global_eqn_numbers ()
 
void add_global_eqn_numbers (std::deque< unsigned long > const &global_eqn_numbers, std::deque< double * > const &global_dof_pt)
 
virtual void assign_internal_and_external_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void assign_additional_local_eqn_numbers ()
 
int internal_local_eqn (const unsigned &i, const unsigned &j) const
 
int external_local_eqn (const unsigned &i, const unsigned &j)
 
void fill_in_jacobian_from_internal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_internal_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
virtual void update_before_internal_fd ()
 
virtual void reset_after_internal_fd ()
 
virtual void update_in_internal_fd (const unsigned &i)
 
virtual void reset_in_internal_fd (const unsigned &i)
 
virtual void update_before_external_fd ()
 
virtual void reset_after_external_fd ()
 
virtual void update_in_external_fd (const unsigned &i)
 
virtual void reset_in_external_fd (const unsigned &i)
 
virtual void fill_in_contribution_to_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_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)
 

Static Protected Attributes

static const unsigned Pconv []
 
- 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
 

Private Member Functions

const unsigned Initial_Nvalue [6]
 
const unsigned Initial_Nvalue [10]
 

Static Private Attributes

static const unsigned Initial_Nvalue []
 Static array of ints to hold number of variables at node. More...
 

Additional Inherited Members

- Public Types inherited from oomph::NavierStokesEquations< DIM >
typedef void(* NavierStokesBodyForceFctPt) (const double &time, const Vector< double > &x, Vector< double > &body_force)
 
typedef double(* NavierStokesSourceFctPt) (const double &time, const Vector< double > &x)
 
typedef double(* NavierStokesPressureAdvDiffSourceFctPt) (const Vector< double > &x)
 
- Public Types inherited from oomph::FiniteElement
typedef void(* SteadyExactSolutionFctPt) (const Vector< double > &, Vector< double > &)
 
typedef void(* UnsteadyExactSolutionFctPt) (const double &, const Vector< double > &, Vector< double > &)
 
- Protected Attributes inherited from oomph::NavierStokesEquations< DIM >
doubleViscosity_Ratio_pt
 
doubleDensity_Ratio_pt
 
doubleRe_pt
 Pointer to global Reynolds number. More...
 
doubleReSt_pt
 Pointer to global Reynolds number x Strouhal number (=Womersley) More...
 
doubleReInvFr_pt
 
Vector< double > * G_pt
 Pointer to global gravity Vector. More...
 
NavierStokesBodyForceFctPt Body_force_fct_pt
 Pointer to body force function. More...
 
NavierStokesSourceFctPt Source_fct_pt
 Pointer to volumetric source function. More...
 
NavierStokesPressureAdvDiffSourceFctPt Press_adv_diff_source_fct_pt
 
bool ALE_is_disabled
 
Vector< FpPressureAdvDiffRobinBCElementBase * > Pressure_advection_diffusion_robin_element_pt
 
int Pinned_fp_pressure_eqn
 
- 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
 

Detailed Description

template<unsigned DIM>
class oomph::TTaylorHoodElement< DIM >

///////////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////////// Taylor–Hood elements are Navier–Stokes elements with quadratic interpolation for velocities and positions and continous linear pressure interpolation

Constructor & Destructor Documentation

◆ TTaylorHoodElement() [1/2]

template<unsigned DIM>
oomph::TTaylorHoodElement< DIM >::TTaylorHoodElement ( )
inline

Constructor, no internal data points.

810 : TElement<DIM, 3>(), NavierStokesEquations<DIM>() {}

◆ TTaylorHoodElement() [2/2]

template<unsigned DIM>
oomph::TTaylorHoodElement< DIM >::TTaylorHoodElement ( const TTaylorHoodElement< DIM > &  dummy)
delete

Broken copy constructor.

Member Function Documentation

◆ build_fp_press_adv_diff_robin_bc_element()

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::build_fp_press_adv_diff_robin_bc_element ( const unsigned face_index)
inlinevirtual

Build FaceElements that apply the Robin boundary condition to the pressure advection diffusion problem required by Fp preconditioner

Implements oomph::NavierStokesEquations< DIM >.

890  {
892  new FpPressureAdvDiffRobinBCElement<TTaylorHoodElement<DIM>>(
893  this, face_index));
894  }
Vector< FpPressureAdvDiffRobinBCElementBase * > Pressure_advection_diffusion_robin_element_pt
Definition: navier_stokes_elements.h:475

References oomph::NavierStokesEquations< DIM >::Pressure_advection_diffusion_robin_element_pt.

◆ dpshape_and_dptest_eulerian_nst() [1/3]

template<unsigned DIM>
virtual double oomph::TTaylorHoodElement< DIM >::dpshape_and_dptest_eulerian_nst ( const Vector< double > &  s,
Shape ppsi,
DShape dppsidx,
Shape ptest,
DShape dptestdx 
) const
protectedvirtual

Compute the pressure shape and test functions and derivatives w.r.t. global coords at local coordinate s. Return Jacobian of mapping between local and global coordinates.

Implements oomph::NavierStokesEquations< DIM >.

◆ dpshape_and_dptest_eulerian_nst() [2/3]

double oomph::TTaylorHoodElement< 2 >::dpshape_and_dptest_eulerian_nst ( const Vector< double > &  s,
Shape ppsi,
DShape dppsidx,
Shape ptest,
DShape dptestdx 
) const
inlineprotectedvirtual

2D : Pressure shape and test functions and derivs w.r.t. to Eulerian coords. Return Jacobian of mapping between local and global coordinates.

Implements oomph::NavierStokesEquations< DIM >.

1147  {
1148  ppsi[0] = s[0];
1149  ppsi[1] = s[1];
1150  ppsi[2] = 1.0 - s[0] - s[1];
1151 
1152  dppsidx(0, 0) = 1.0;
1153  dppsidx(0, 1) = 0.0;
1154 
1155  dppsidx(1, 0) = 0.0;
1156  dppsidx(1, 1) = 1.0;
1157 
1158  dppsidx(2, 0) = -1.0;
1159  dppsidx(2, 1) = -1.0;
1160 
1161  // Allocate memory for the inverse 2x2 jacobian
1162  DenseMatrix<double> inverse_jacobian(2);
1163 
1164 
1165  // Get the values of the shape functions and their local derivatives
1166  Shape psi(6);
1167  DShape dpsi(6, 2);
1168  dshape_local(s, psi, dpsi);
1169 
1170  // Now calculate the inverse jacobian
1171  const double det = local_to_eulerian_mapping(dpsi, inverse_jacobian);
1172 
1173  // Now set the values of the derivatives to be derivs w.r.t. to the
1174  // Eulerian coordinates
1175  transform_derivatives(inverse_jacobian, dppsidx);
1176 
1177  // Test functions are shape functions
1178  ptest = ppsi;
1179  dptestdx = dppsidx;
1180 
1181  // Return the determinant of the jacobian
1182  return det;
1183  }
virtual double local_to_eulerian_mapping(const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
Definition: elements.h:1508
virtual void transform_derivatives(const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const
Definition: elements.cc:2833
virtual void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Definition: elements.h:1981
RealScalar s
Definition: level1_cplx_impl.h:130

References s.

◆ dpshape_and_dptest_eulerian_nst() [3/3]

double oomph::TTaylorHoodElement< 3 >::dpshape_and_dptest_eulerian_nst ( const Vector< double > &  s,
Shape ppsi,
DShape dppsidx,
Shape ptest,
DShape dptestdx 
) const
inlineprotectedvirtual

3D : Pressure shape and test functions and derivs w.r.t. to Eulerian coords. Return Jacobian of mapping between local and global coordinates.

Implements oomph::NavierStokesEquations< DIM >.

1198  {
1199  ppsi[0] = s[0];
1200  ppsi[1] = s[1];
1201  ppsi[2] = s[2];
1202  ppsi[3] = 1.0 - s[0] - s[1] - s[2];
1203 
1204  dppsidx(0, 0) = 1.0;
1205  dppsidx(0, 1) = 0.0;
1206  dppsidx(0, 2) = 0.0;
1207 
1208  dppsidx(1, 0) = 0.0;
1209  dppsidx(1, 1) = 1.0;
1210  dppsidx(1, 2) = 0.0;
1211 
1212  dppsidx(2, 0) = 0.0;
1213  dppsidx(2, 1) = 0.0;
1214  dppsidx(2, 2) = 1.0;
1215 
1216  dppsidx(3, 0) = -1.0;
1217  dppsidx(3, 1) = -1.0;
1218  dppsidx(3, 2) = -1.0;
1219 
1220 
1221  // Get the values of the shape functions and their local derivatives
1222  Shape psi(10);
1223  DShape dpsi(10, 3);
1224  dshape_local(s, psi, dpsi);
1225 
1226  // Allocate memory for the inverse 3x3 jacobian
1227  DenseMatrix<double> inverse_jacobian(3);
1228 
1229  // Now calculate the inverse jacobian
1230  const double det = local_to_eulerian_mapping(dpsi, inverse_jacobian);
1231 
1232  // Now set the values of the derivatives to be derivs w.r.t. to the
1233  // Eulerian coordinates
1234  transform_derivatives(inverse_jacobian, dppsidx);
1235 
1236  // Test functions are shape functions
1237  ptest = ppsi;
1238  dptestdx = dppsidx;
1239 
1240  // Return the determinant of the jacobian
1241  return det;
1242  }

References s.

◆ dshape_and_dtest_eulerian_at_knot_nst() [1/4]

template<unsigned DIM>
double oomph::TTaylorHoodElement< DIM >::dshape_and_dtest_eulerian_at_knot_nst ( const unsigned ipt,
Shape psi,
DShape dpsidx,
RankFourTensor< double > &  d_dpsidx_dX,
Shape test,
DShape dtestdx,
RankFourTensor< double > &  d_dtestdx_dX,
DenseMatrix< double > &  djacobian_dX 
) const
inlineprotectedvirtual

Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of mapping (J). Also compute derivatives of dpsidx, dtestdx and J w.r.t. nodal coordinates.

Implements oomph::NavierStokesEquations< DIM >.

◆ dshape_and_dtest_eulerian_at_knot_nst() [2/4]

double oomph::TTaylorHoodElement< 2 >::dshape_and_dtest_eulerian_at_knot_nst ( const unsigned ipt,
Shape psi,
DShape dpsidx,
RankFourTensor< double > &  d_dpsidx_dX,
Shape test,
DShape dtestdx,
RankFourTensor< double > &  d_dtestdx_dX,
DenseMatrix< double > &  djacobian_dX 
) const
inlineprotectedvirtual

2D : Define the shape functions (psi) and test functions (test) and their derivatives w.r.t. global coordinates (dpsidx and dtestdx) and return Jacobian of mapping (J). Additionally compute the derivatives of dpsidx, dtestdx and J w.r.t. nodal coordinates.

Galerkin: Test functions = shape functions

Implements oomph::NavierStokesEquations< DIM >.

1264  {
1265  // Call the geometrical shape functions and derivatives
1266  const double J = this->dshape_eulerian_at_knot(
1267  ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
1268 
1269  // Loop over the test functions and derivatives and set them equal to the
1270  // shape functions
1271  for (unsigned i = 0; i < 6; i++)
1272  {
1273  test[i] = psi[i];
1274 
1275  for (unsigned k = 0; k < 2; k++)
1276  {
1277  dtestdx(i, k) = dpsidx(i, k);
1278 
1279  for (unsigned p = 0; p < 2; p++)
1280  {
1281  for (unsigned q = 0; q < 6; q++)
1282  {
1283  d_dtestdx_dX(p, q, i, k) = d_dpsidx_dX(p, q, i, k);
1284  }
1285  }
1286  }
1287  }
1288 
1289  // Return the jacobian
1290  return J;
1291  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
float * p
Definition: Tutorial_Map_using.cpp:9
virtual double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3325
char char char int int * k
Definition: level2_impl.h:374
EIGEN_DEVICE_FUNC const Scalar & q
Definition: SpecialFunctionsImpl.h:2019
Definition: indexed_view.cpp:20

References i, J, k, p, and Eigen::numext::q.

◆ dshape_and_dtest_eulerian_at_knot_nst() [3/4]

double oomph::TTaylorHoodElement< 3 >::dshape_and_dtest_eulerian_at_knot_nst ( const unsigned ipt,
Shape psi,
DShape dpsidx,
RankFourTensor< double > &  d_dpsidx_dX,
Shape test,
DShape dtestdx,
RankFourTensor< double > &  d_dtestdx_dX,
DenseMatrix< double > &  djacobian_dX 
) const
inlineprotectedvirtual

3D : Define the shape functions (psi) and test functions (test) and their derivatives w.r.t. global coordinates (dpsidx and dtestdx) and return Jacobian of mapping (J). Additionally compute the derivatives of dpsidx, dtestdx and J w.r.t. nodal coordinates.

Galerkin: Test functions = shape functions

Implements oomph::NavierStokesEquations< DIM >.

1313  {
1314  // Call the geometrical shape functions and derivatives
1315  const double J = this->dshape_eulerian_at_knot(
1316  ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
1317 
1318  // Loop over the test functions and derivatives and set them equal to the
1319  // shape functions
1320  for (unsigned i = 0; i < 10; i++)
1321  {
1322  test[i] = psi[i];
1323 
1324  for (unsigned k = 0; k < 3; k++)
1325  {
1326  dtestdx(i, k) = dpsidx(i, k);
1327 
1328  for (unsigned p = 0; p < 3; p++)
1329  {
1330  for (unsigned q = 0; q < 10; q++)
1331  {
1332  d_dtestdx_dX(p, q, i, k) = d_dpsidx_dX(p, q, i, k);
1333  }
1334  }
1335  }
1336  }
1337 
1338  // Return the jacobian
1339  return J;
1340  }

References i, J, k, p, and Eigen::numext::q.

◆ dshape_and_dtest_eulerian_at_knot_nst() [4/4]

template<unsigned DIM>
double oomph::TTaylorHoodElement< DIM >::dshape_and_dtest_eulerian_at_knot_nst ( const unsigned ipt,
Shape psi,
DShape dpsidx,
Shape test,
DShape dtestdx 
) const
inlineprotectedvirtual

Velocity shape and test functions and their derivs w.r.t. to global coords at local coordinate s (taken from geometry) Return Jacobian of mapping between local and global coordinates.

Derivatives of the shape functions and test functions w.r.t to global (Eulerian) coordinates. Return Jacobian of mapping between local and global coordinates.

Implements oomph::NavierStokesEquations< DIM >.

1125  {
1126  // Call the geometrical shape functions and derivatives
1127  double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
1128  // Test functions are the shape functions
1129  test = psi;
1130  dtestdx = dpsidx;
1131  // Return the jacobian
1132  return J;
1133  }

References J.

◆ dshape_and_dtest_eulerian_nst()

template<unsigned DIM>
double oomph::TTaylorHoodElement< DIM >::dshape_and_dtest_eulerian_nst ( const Vector< double > &  s,
Shape psi,
DShape dpsidx,
Shape test,
DShape dtestdx 
) const
inlineprotectedvirtual

Velocity shape and test functions and their derivs w.r.t. to global coords at local coordinate s (taken from geometry) Return Jacobian of mapping between local and global coordinates.

2D : Derivatives of the shape functions and test functions w.r.t to global (Eulerian) coordinates. Return Jacobian of mapping between local and global coordinates.

Implements oomph::NavierStokesEquations< DIM >.

1102  {
1103  // Call the geometrical shape functions and derivatives
1104  double J = this->dshape_eulerian(s, psi, dpsidx);
1105  // Test functions are the shape functions
1106  test = psi;
1107  dtestdx = dpsidx;
1108  // Return the jacobian
1109  return J;
1110  }
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3298

References J, and s.

◆ fix_pressure()

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::fix_pressure ( const unsigned p_dof,
const double p_value 
)
inlinevirtual

Pin p_dof-th pressure dof and set it to value specified by p_value.

Implements oomph::NavierStokesEquations< DIM >.

880  {
881  this->node_pt(Pconv[p_dof])->pin(DIM);
882  this->node_pt(Pconv[p_dof])->set_value(DIM, p_value);
883  }
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
static const unsigned Pconv[]
Definition: Tnavier_stokes_elements.h:756
#define DIM
Definition: linearised_navier_stokes_elements.h:44

References DIM, oomph::FiniteElement::node_pt(), oomph::Data::pin(), and oomph::Data::set_value().

◆ get_dof_numbers_for_unknowns()

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::get_dof_numbers_for_unknowns ( std::list< std::pair< unsigned long, unsigned >> &  dof_lookup_list) const
inlinevirtual

Create a list of pairs for all unknowns in this element, so that the first entry in each pair contains the global equation number of the unknown, while the second one contains the number of the "DOF type" that this unknown is associated with. (Function can obviously only be called if the equation numbering scheme has been set up.) Velocity=0; Pressure=1

Reimplemented from oomph::GeneralisedElement.

Reimplemented in oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.

1026  {
1027  // number of nodes
1028  unsigned n_node = this->nnode();
1029 
1030  // temporary pair (used to store dof lookup prior to being added to list)
1031  std::pair<unsigned, unsigned> dof_lookup;
1032 
1033  // loop over the nodes
1034  for (unsigned n = 0; n < n_node; n++)
1035  {
1036  // find the number of Navier Stokes values at this node
1037  unsigned nv = this->required_nvalue(n);
1038 
1039  // loop over these values
1040  for (unsigned v = 0; v < nv; v++)
1041  {
1042  // determine local eqn number
1043  int local_eqn_number = this->nodal_local_eqn(n, v);
1044 
1045  // ignore pinned values - far away degrees of freedom resulting
1046  // from hanging nodes can be ignored since these are be dealt
1047  // with by the element containing their master nodes
1048  if (local_eqn_number >= 0)
1049  {
1050  // store dof lookup in temporary pair: Global equation number
1051  // is the first entry in pair
1052  dof_lookup.first = this->eqn_number(local_eqn_number);
1053 
1054  // set dof numbers: Dof number is the second entry in pair
1055  dof_lookup.second = v;
1056 
1057  // add to list
1058  dof_lookup_list.push_front(dof_lookup);
1059  }
1060  }
1061  }
1062  }
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
unsigned long eqn_number(const unsigned &ieqn_local) const
Definition: elements.h:704
int local_eqn_number(const unsigned long &ieqn_global) const
Definition: elements.h:726
virtual unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
Definition: Tnavier_stokes_elements.h:821

References oomph::GeneralisedElement::eqn_number(), oomph::GeneralisedElement::local_eqn_number(), n, oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_local_eqn(), oomph::TTaylorHoodElement< DIM >::required_nvalue(), and v.

◆ get_Z2_flux()

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::get_Z2_flux ( const Vector< double > &  s,
Vector< double > &  flux 
)
inlinevirtual

Get 'flux' for Z2 error recovery: Upper triangular entries in strain rate tensor.

Implements oomph::ElementWithZ2ErrorEstimator.

Reimplemented in oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.

971  {
972 #ifdef PARANOID
973  unsigned num_entries = DIM + (DIM * (DIM - 1)) / 2;
974  if (flux.size() < num_entries)
975  {
976  std::ostringstream error_message;
977  error_message << "The flux vector has the wrong number of entries, "
978  << flux.size() << ", whereas it should be at least "
979  << num_entries << std::endl;
980  throw OomphLibError(error_message.str(),
983  }
984 #endif
985 
986  // Get strain rate matrix
987  DenseMatrix<double> strainrate(DIM);
988  this->strain_rate(s, strainrate);
989 
990  // Pack into flux Vector
991  unsigned icount = 0;
992 
993  // Start with diagonal terms
994  for (unsigned i = 0; i < DIM; i++)
995  {
996  flux[icount] = strainrate(i, i);
997  icount++;
998  }
999 
1000  // Off diagonals row by row
1001  for (unsigned i = 0; i < DIM; i++)
1002  {
1003  for (unsigned j = i + 1; j < DIM; j++)
1004  {
1005  flux[icount] = strainrate(i, j);
1006  icount++;
1007  }
1008  }
1009  }
void strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate) const
Strain-rate tensor: 1/2 (du_i/dx_j + du_j/dx_i)
Definition: navier_stokes_elements.cc:1185
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
#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 DIM, ProblemParameters::flux(), i, j, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::NavierStokesEquations< DIM >::strain_rate().

◆ identify_load_data()

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::identify_load_data ( std::set< std::pair< Data *, unsigned >> &  paired_load_data)
virtual

Add to the set paired_load_data pairs containing

  • the pointer to a Data object and
  • the index of the value in that Data object

for all values (pressures, velocities) that affect the load computed in the get_load(...) function.

Implements oomph::FSIFluidElement.

214  {
215  // Loop over the nodes
216  unsigned n_node = this->nnode();
217  for (unsigned n = 0; n < n_node; n++)
218  {
219  // Loop over the velocity components and add pointer to their data
220  // and indices to the vectors
221  for (unsigned i = 0; i < DIM; i++)
222  {
223  paired_load_data.insert(std::make_pair(this->node_pt(n), i));
224  }
225  }
226 
227  // Add the pressure data
228  identify_pressure_data(paired_load_data);
229  }
void identify_pressure_data(std::set< std::pair< Data *, unsigned >> &paired_pressure_data)
Definition: Tnavier_stokes_elements.cc:241

References DIM, i, and n.

◆ identify_pressure_data()

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::identify_pressure_data ( std::set< std::pair< Data *, unsigned >> &  paired_load_data)
virtual

Add to the set paired_pressure_data pairs containing

  • the pointer to a Data object and
  • the index of the value in that Data object

for all pressure values that affect the load computed in the get_load(...) function.

Add to the set paired_load_data pairs containing

  • the pointer to a Data object and
  • the index of the value in that Data object

for all values (pressures, velocities) that affect the load computed in the get_load(...) function.

Implements oomph::FSIFluidElement.

243  {
244  // Loop over the pressure data
245  unsigned n_pres = npres_nst();
246  for (unsigned l = 0; l < n_pres; l++)
247  {
248  // The DIMth entry in each nodal data is the pressure, which
249  // affects the traction
250  paired_load_data.insert(std::make_pair(this->node_pt(Pconv[l]), DIM));
251  }
252  }
unsigned npres_nst() const
Return number of pressure values.

References DIM.

◆ Initial_Nvalue() [1/2]

const unsigned oomph::TTaylorHoodElement< 3 >::Initial_Nvalue
private

◆ Initial_Nvalue() [2/2]

const unsigned oomph::TTaylorHoodElement< 2 >::Initial_Nvalue
private

//////////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////////////

◆ ndof_types()

template<unsigned DIM>
unsigned oomph::TTaylorHoodElement< DIM >::ndof_types ( ) const
inlinevirtual

The number of "DOF types" that degrees of freedom in this element are sub-divided into: Velocity and pressure.

Reimplemented from oomph::GeneralisedElement.

Reimplemented in oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.

1014  {
1015  return DIM + 1;
1016  }

References DIM.

◆ npres_nst() [1/3]

template<unsigned DIM>
unsigned oomph::TTaylorHoodElement< DIM >::npres_nst ( ) const
virtual

Return number of pressure values.

Implements oomph::NavierStokesEquations< DIM >.

◆ npres_nst() [2/3]

unsigned oomph::TTaylorHoodElement< 2 >::npres_nst ( ) const
inlinevirtual

2D : Number of pressure values

Implements oomph::NavierStokesEquations< DIM >.

1074  {
1075  return 3;
1076  }

◆ npres_nst() [3/3]

unsigned oomph::TTaylorHoodElement< 3 >::npres_nst ( ) const
inlinevirtual

3D : Number of pressure values

Implements oomph::NavierStokesEquations< DIM >.

1084  {
1085  return 4;
1086  }

◆ nrecovery_order()

template<unsigned DIM>
unsigned oomph::TTaylorHoodElement< DIM >::nrecovery_order ( )
inlinevirtual

Order of recovery shape functions for Z2 error estimation: Same order as shape functions.

Implements oomph::ElementWithZ2ErrorEstimator.

Reimplemented in oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.

944  {
945  return 2;
946  }

◆ num_Z2_flux_terms()

template<unsigned DIM>
unsigned oomph::TTaylorHoodElement< DIM >::num_Z2_flux_terms ( )
inlinevirtual

Number of 'flux' terms for Z2 error estimation.

Implements oomph::ElementWithZ2ErrorEstimator.

Reimplemented in oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.

963  {
964  // DIM diagonal strain rates, DIM(DIM -1) /2 off diagonal rates
965  return DIM + (DIM * (DIM - 1)) / 2;
966  }

References DIM.

◆ nvertex_node()

template<unsigned DIM>
unsigned oomph::TTaylorHoodElement< DIM >::nvertex_node ( ) const
inlinevirtual

Number of vertex nodes in the element.

Implements oomph::ElementWithZ2ErrorEstimator.

Reimplemented in oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.

950  {
951  return DIM + 1;
952  }

References DIM.

◆ output() [1/4]

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::output ( FILE *  file_pt)
inlinevirtual

Redirect output to NavierStokesEquations output.

Reimplemented from oomph::NavierStokesEquations< DIM >.

Reimplemented in oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.

931  {
933  }
void output(std::ostream &outfile)
Definition: navier_stokes_elements.h:1155

References oomph::NavierStokesEquations< DIM >::output().

◆ output() [2/4]

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::output ( FILE *  file_pt,
const unsigned n_plot 
)
inlinevirtual

Redirect output to NavierStokesEquations output.

Reimplemented from oomph::NavierStokesEquations< DIM >.

Reimplemented in oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.

937  {
938  NavierStokesEquations<DIM>::output(file_pt, n_plot);
939  }

References oomph::NavierStokesEquations< DIM >::output().

◆ output() [3/4]

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::output ( std::ostream &  outfile)
inlinevirtual

Redirect output to NavierStokesEquations output.

Reimplemented from oomph::NavierStokesEquations< DIM >.

Reimplemented in oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.

919  {
921  }

References oomph::NavierStokesEquations< DIM >::output().

◆ output() [4/4]

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::output ( std::ostream &  outfile,
const unsigned nplot 
)
inlinevirtual

Redirect output to NavierStokesEquations output.

Reimplemented from oomph::NavierStokesEquations< DIM >.

Reimplemented in oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.

925  {
926  NavierStokesEquations<DIM>::output(outfile, nplot);
927  }

References oomph::NavierStokesEquations< DIM >::output().

◆ p_index_nst()

template<unsigned DIM>
unsigned oomph::TTaylorHoodElement< DIM >::p_index_nst ( )
inline

Which nodal value represents the pressure?

841  {
842  return DIM;
843  }

References DIM.

◆ p_local_eqn()

template<unsigned DIM>
int oomph::TTaylorHoodElement< DIM >::p_local_eqn ( const unsigned n) const
inlinevirtual

Pointer to n_p-th pressure node.

Return the local equation numbers for the pressure values.

Implements oomph::TemplateFreeNavierStokesEquationsBase.

851  {
852  return this->nodal_local_eqn(Pconv[n], DIM);
853  }

References DIM, n, and oomph::FiniteElement::nodal_local_eqn().

◆ p_nodal_index_nst()

template<unsigned DIM>
int oomph::TTaylorHoodElement< DIM >::p_nodal_index_nst ( ) const
inlinevirtual

Set the value at which the pressure is stored in the nodes.

Reimplemented from oomph::NavierStokesEquations< DIM >.

871  {
872  return static_cast<int>(DIM);
873  }

References DIM.

◆ p_nst() [1/2]

template<unsigned DIM>
double oomph::TTaylorHoodElement< DIM >::p_nst ( const unsigned n_p) const
inlinevirtual

Access function for the pressure values at local pressure node n_p (const version)

Implements oomph::NavierStokesEquations< DIM >.

858  {
859  return this->nodal_value(Pconv[n_p], DIM);
860  }
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593

References DIM, and oomph::FiniteElement::nodal_value().

◆ p_nst() [2/2]

template<unsigned DIM>
double oomph::TTaylorHoodElement< DIM >::p_nst ( const unsigned t,
const unsigned n_p 
) const
inlinevirtual

Access function for the pressure values at local pressure node n_p (const version)

Implements oomph::NavierStokesEquations< DIM >.

865  {
866  return this->nodal_value(t, Pconv[n_p], DIM);
867  }

References DIM, oomph::FiniteElement::nodal_value(), and oomph::TTaylorHoodElement< DIM >::Pconv.

◆ Pconv() [1/2]

const unsigned oomph::TTaylorHoodElement< 2 >::Pconv
protected

◆ Pconv() [2/2]

const unsigned oomph::TTaylorHoodElement< 3 >::Pconv
protected

◆ pin_all_nodal_pressure_dofs()

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::pin_all_nodal_pressure_dofs
protected

Pin all nodal pressure dofs.

Pin all nodal pressure dofs, incl the mid-face/side ones where they have been allocated (e.g. in the refineable version of this element).

171  {
172  // Loop over all nodes and pin pressure
173  unsigned n_node = this->nnode();
174  for (unsigned n = 0; n < n_node; n++)
175  {
176  if (this->node_pt(n)->nvalue() == DIM + 1)
177  {
178  this->node_pt(n)->pin(DIM);
179  }
180  }
181  }

References DIM, and n.

◆ pshape_nst() [1/4]

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::pshape_nst ( const Vector< double > &  s,
Shape psi 
) const
inlinevirtual

Test whether the pressure dof p_dof hanging or not?

Pressure shape functions at local coordinate s

Implements oomph::NavierStokesEquations< DIM >.

◆ pshape_nst() [2/4]

void oomph::TTaylorHoodElement< 2 >::pshape_nst ( const Vector< double > &  s,
Shape psi 
) const
inlinevirtual

2D : Pressure shape functions

Implements oomph::NavierStokesEquations< DIM >.

1350  {
1351  psi[0] = s[0];
1352  psi[1] = s[1];
1353  psi[2] = 1.0 - s[0] - s[1];
1354  }

References s.

◆ pshape_nst() [3/4]

void oomph::TTaylorHoodElement< 3 >::pshape_nst ( const Vector< double > &  s,
Shape psi 
) const
inlinevirtual

3D : Pressure shape functions

Implements oomph::NavierStokesEquations< DIM >.

1363  {
1364  psi[0] = s[0];
1365  psi[1] = s[1];
1366  psi[2] = s[2];
1367  psi[3] = 1.0 - s[0] - s[1] - s[2];
1368  }

References s.

◆ pshape_nst() [4/4]

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::pshape_nst ( const Vector< double > &  s,
Shape psi,
Shape test 
) const
inlinevirtual

Pressure shape and test functions at local coordinte s.

Pressure shape and test functions.

Implements oomph::NavierStokesEquations< DIM >.

1378  {
1379  // Call the pressure shape functions
1380  this->pshape_nst(s, psi);
1381  // Test functions are shape functions
1382  test = psi;
1383  }
void pshape_nst(const Vector< double > &s, Shape &psi) const
Test whether the pressure dof p_dof hanging or not?

References s.

◆ required_nvalue()

template<unsigned DIM>
virtual unsigned oomph::TTaylorHoodElement< DIM >::required_nvalue ( const unsigned n) const
inlinevirtual

Broken assignment operator.

Number of values (pinned or dofs) required at node n. Can be overwritten for hanging node version

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.

822  {
823  return Initial_Nvalue[n];
824  }
static const unsigned Initial_Nvalue[]
Static array of ints to hold number of variables at node.
Definition: Tnavier_stokes_elements.h:751

References oomph::TTaylorHoodElement< DIM >::Initial_Nvalue, and n.

Referenced by oomph::TTaylorHoodElement< DIM >::get_dof_numbers_for_unknowns().

◆ unpin_all_nodal_pressure_dofs()

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::unpin_all_nodal_pressure_dofs
protected

Unpin all pressure dofs.

Unpin all pressure dofs, incl the mid-face/side ones where they have been allocated (e.g. in the refineable version of this element).

151  {
152  unsigned n_node = this->nnode();
153  // loop over nodes
154  for (unsigned l = 0; l < n_node; l++)
155  {
156  if (this->node_pt(l)->nvalue() == DIM + 1)
157  {
158  // unpin pressure dof
159  this->node_pt(l)->unpin(DIM);
160  }
161  }
162  }
void unpin(const unsigned &i)
Unpin the i-th stored variable.
Definition: nodes.h:391

References DIM.

◆ unpin_proper_nodal_pressure_dofs()

template<unsigned DIM>
void oomph::TTaylorHoodElement< DIM >::unpin_proper_nodal_pressure_dofs
protected

Unpin the proper nodal pressure dofs.

Unpin the proper nodal pressure dofs which are not hanging.

188  {
189  // Loop over all pressure nodes and unpin if they're not hanging
190  unsigned n_pres = npres_nst();
191  for (unsigned l = 0; l < n_pres; l++)
192  {
193  Node* nod_pt = this->node_pt(Pconv[l]);
194  if (!nod_pt->is_hanging(DIM))
195  {
196  nod_pt->unpin(DIM);
197  }
198  }
199  }

References DIM, oomph::Node::is_hanging(), and oomph::Data::unpin().

◆ vertex_node_pt()

template<unsigned DIM>
Node* oomph::TTaylorHoodElement< DIM >::vertex_node_pt ( const unsigned j) const
inlinevirtual

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

Implements oomph::ElementWithZ2ErrorEstimator.

Reimplemented in oomph::PseudoSolidNodeUpdateElement< TTaylorHoodElement< 2 >, TPVDElement< 2, 3 > >.

956  {
957  return node_pt(j);
958  }

References j, and oomph::FiniteElement::node_pt().

Member Data Documentation

◆ Initial_Nvalue

template<unsigned DIM>
const unsigned oomph::TTaylorHoodElement< DIM >::Initial_Nvalue[]
staticprivate

Static array of ints to hold number of variables at node.

Referenced by oomph::TTaylorHoodElement< DIM >::required_nvalue().

◆ Pconv

template<unsigned DIM>
const unsigned oomph::TTaylorHoodElement< DIM >::Pconv[]
staticprotected

Static array of ints to hold conversion from pressure node numbers to actual node numbers

Referenced by oomph::TTaylorHoodElement< DIM >::p_nst().

◆ TEMPLATE_PARAMETER_DIM

template<unsigned DIM>
const unsigned oomph::TTaylorHoodElement< DIM >::TEMPLATE_PARAMETER_DIM = DIM
static

Publicly exposed template parameter.

◆ TEMPLATE_PARAMETER_NNODE_1D

template<unsigned DIM>
const unsigned oomph::TTaylorHoodElement< DIM >::TEMPLATE_PARAMETER_NNODE_1D = 3
static

Publicly exposed template parameter.


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