oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D > Class Template Reference

#include <supg_advection_diffusion_elements.h>

+ Inheritance diagram for oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D >:

Public Member Functions

 QSUPGAdvectionDiffusionElement ()
 
double get_Tau_SUPG ()
 Get stabilisation parameter for the element. More...
 
void switch_off_stabilisation ()
 Set stabilisation parameter for the element to zero. More...
 
void compute_stabilisation_parameter ()
 Compute stabilisation parameter for the element. More...
 
void output (std::ostream &outfile, const unsigned &nplot)
 
void output (std::ostream &outfile)
 Output at default number of plot points. More...
 
void output (FILE *file_pt)
 C-style output. More...
 
void output (FILE *file_pt, const unsigned &n_plot)
 C_style output at n_plot points. More...
 
- Public Member Functions inherited from oomph::QAdvectionDiffusionElement< DIM, NNODE_1D >
 QAdvectionDiffusionElement ()
 
 QAdvectionDiffusionElement (const QAdvectionDiffusionElement< DIM, NNODE_1D > &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const QAdvectionDiffusionElement< DIM, NNODE_1D > &)=delete
 Broken assignment operator. More...
 
unsigned required_nvalue (const unsigned &n) const
 
void output_fct (std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 
void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
 
- Public Member Functions inherited from oomph::AdvectionDiffusionEquations< DIM >
 AdvectionDiffusionEquations ()
 
 AdvectionDiffusionEquations (const AdvectionDiffusionEquations &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const AdvectionDiffusionEquations &)=delete
 Broken assignment operator. More...
 
virtual unsigned u_index_adv_diff () const
 
double du_dt_adv_diff (const unsigned &n) const
 
void disable_ALE ()
 
void enable_ALE ()
 
unsigned nscalar_paraview () const
 
void scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
 
std::string scalar_name_paraview (const unsigned &i) const
 
void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 Get error against and norm of exact solution. More...
 
void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Dummy, time dependent error checker. More...
 
AdvectionDiffusionSourceFctPtsource_fct_pt ()
 Access function: Pointer to source function. More...
 
AdvectionDiffusionSourceFctPt source_fct_pt () const
 Access function: Pointer to source function. Const version. More...
 
AdvectionDiffusionWindFctPtwind_fct_pt ()
 Access function: Pointer to wind function. More...
 
AdvectionDiffusionWindFctPt wind_fct_pt () const
 Access function: Pointer to wind function. Const version. More...
 
const doublepe () const
 Peclet number. More...
 
double *& pe_pt ()
 Pointer to Peclet number. More...
 
const doublepe_st () const
 Peclet number multiplied by Strouhal number. More...
 
double *& pe_st_pt ()
 Pointer to Peclet number multipled by Strouha number. More...
 
virtual void get_source_adv_diff (const unsigned &ipt, const Vector< double > &x, double &source) const
 
virtual void get_wind_adv_diff (const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
 
void get_flux (const Vector< double > &s, Vector< double > &flux) const
 Get flux: \(\mbox{flux}[i] = \mbox{d}u / \mbox{d}x_i \). More...
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Add the element's contribution to its residual vector (wrapper) 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)
 
double interpolated_u_adv_diff (const Vector< double > &s) const
 Return FE representation of function value u(s) at local coordinate s. More...
 
virtual void dinterpolated_u_adv_diff_ddata (const Vector< double > &s, Vector< double > &du_ddata, Vector< unsigned > &global_eqn_number)
 
unsigned self_test ()
 Self-test: Return 0 for OK. 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
 
virtual void get_dresidual_dnodal_coordinates (RankThreeTensor< double > &dresidual_dnodal_coordinates)
 
unsigned nnodal_position_type () const
 
bool has_hanging_nodes () const
 
unsigned nodal_dimension () const
 Return the required Eulerian dimension of the nodes in this element. More...
 
virtual unsigned nvertex_node () const
 
virtual Nodevertex_node_pt (const unsigned &j) const
 
virtual Nodeconstruct_node (const unsigned &n)
 
virtual Nodeconstruct_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
virtual Nodeconstruct_boundary_node (const unsigned &n)
 
virtual Nodeconstruct_boundary_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
int get_node_number (Node *const &node_pt) const
 
virtual Nodeget_node_at_local_coordinate (const Vector< double > &s) const
 
double raw_nodal_value (const unsigned &n, const unsigned &i) const
 
double raw_nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
unsigned dim () const
 
virtual ElementGeometry::ElementGeometry element_geometry () const
 Return the geometry type of the element (either Q or T usually). More...
 
virtual double interpolated_x (const Vector< double > &s, const unsigned &i) const
 Return FE interpolated coordinate x[i] at local coordinate s. More...
 
virtual double interpolated_x (const unsigned &t, const Vector< double > &s, const unsigned &i) const
 
virtual void interpolated_x (const Vector< double > &s, Vector< double > &x) const
 Return FE interpolated position x[] at local coordinate s as Vector. More...
 
virtual void interpolated_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) const
 
virtual double interpolated_dxdt (const Vector< double > &s, const unsigned &i, const unsigned &t)
 
virtual void interpolated_dxdt (const Vector< double > &s, const unsigned &t, Vector< double > &dxdt)
 
unsigned ngeom_data () const
 
Datageom_data_pt (const unsigned &j)
 
void position (const Vector< double > &zeta, Vector< double > &r) const
 
void position (const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const
 
void dposition_dt (const Vector< double > &zeta, const unsigned &t, Vector< double > &drdt)
 
virtual double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 
void interpolated_zeta (const Vector< double > &s, Vector< double > &zeta) const
 
void locate_zeta (const Vector< double > &zeta, GeomObject *&geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
 
virtual void node_update ()
 
virtual void identify_field_data_for_interactions (std::set< std::pair< Data *, unsigned >> &paired_field_data)
 
virtual 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
 
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 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 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, 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, 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 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)
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_norm (double &norm)
 
virtual unsigned ndof_types () const
 
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 
- Public Member Functions inherited from oomph::GeomObject
 GeomObject ()
 Default constructor. More...
 
 GeomObject (const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim, TimeStepper *time_stepper_pt)
 
 GeomObject (const GeomObject &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const GeomObject &)=delete
 Broken assignment operator. More...
 
virtual ~GeomObject ()
 (Empty) destructor More...
 
unsigned nlagrangian () const
 Access function to # of Lagrangian coordinates. More...
 
unsigned ndim () const
 Access function to # of Eulerian coordinates. More...
 
void set_nlagrangian_and_ndim (const unsigned &n_lagrangian, const unsigned &n_dim)
 Set # of Lagrangian and Eulerian coordinates. More...
 
TimeStepper *& time_stepper_pt ()
 
TimeSteppertime_stepper_pt () const
 
virtual void position (const double &t, const Vector< double > &zeta, Vector< double > &r) const
 
virtual void dposition (const Vector< double > &zeta, DenseMatrix< double > &drdzeta) const
 
virtual void d2position (const Vector< double > &zeta, RankThreeTensor< double > &ddrdzeta) const
 
virtual void d2position (const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const
 

Protected Member Functions

double dshape_and_dtest_eulerian_adv_diff (const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
 
double dshape_and_dtest_eulerian_at_knot_adv_diff (const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
 
- Protected Member Functions inherited from oomph::AdvectionDiffusionEquations< DIM >
virtual void fill_in_generic_residual_contribution_adv_diff (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
 
- Protected Member Functions inherited from oomph::FiniteElement
virtual void assemble_local_to_eulerian_jacobian (const DShape &dpsids, DenseMatrix< double > &jacobian) const
 
virtual void assemble_local_to_eulerian_jacobian2 (const DShape &d2psids, DenseMatrix< double > &jacobian2) const
 
virtual void assemble_eulerian_base_vectors (const DShape &dpsids, DenseMatrix< double > &interpolated_G) const
 
template<unsigned DIM>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual double invert_jacobian_mapping (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual double local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
double local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &inverse_jacobian) const
 
virtual double local_to_eulerian_mapping_diagonal (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual void dJ_eulerian_dnodal_coordinates (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<unsigned DIM>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
virtual void d_dshape_eulerian_dnodal_coordinates (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<unsigned DIM>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
virtual void transform_derivatives (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const
 
void transform_derivatives_diagonal (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const
 
virtual void transform_second_derivatives (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<unsigned DIM>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<unsigned DIM>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
virtual void fill_in_jacobian_from_nodal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
void fill_in_jacobian_from_nodal_by_fd (DenseMatrix< double > &jacobian)
 
virtual void update_before_nodal_fd ()
 
virtual void reset_after_nodal_fd ()
 
virtual void update_in_nodal_fd (const unsigned &i)
 
virtual void reset_in_nodal_fd (const unsigned &i)
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Zero-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 One-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Two-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
- Protected Member Functions inherited from oomph::GeneralisedElement
unsigned add_internal_data (Data *const &data_pt, const bool &fd=true)
 
bool internal_data_fd (const unsigned &i) const
 
void exclude_internal_data_fd (const unsigned &i)
 
void include_internal_data_fd (const unsigned &i)
 
void clear_global_eqn_numbers ()
 
void add_global_eqn_numbers (std::deque< unsigned long > const &global_eqn_numbers, std::deque< double * > const &global_dof_pt)
 
virtual void assign_internal_and_external_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void assign_additional_local_eqn_numbers ()
 
int internal_local_eqn (const unsigned &i, const unsigned &j) const
 
int external_local_eqn (const unsigned &i, const unsigned &j)
 
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_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

double Tau_SUPG
 SUPG stabilisation parameter. More...
 
- Protected Attributes inherited from oomph::AdvectionDiffusionEquations< DIM >
doublePe_pt
 Pointer to global Peclet number. More...
 
doublePeSt_pt
 Pointer to global Peclet number multiplied by Strouhal number. More...
 
AdvectionDiffusionSourceFctPt Source_fct_pt
 Pointer to source function: More...
 
AdvectionDiffusionWindFctPt Wind_fct_pt
 Pointer to wind function: More...
 
bool ALE_is_disabled
 
- 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
 

Additional Inherited Members

- Public Types inherited from oomph::AdvectionDiffusionEquations< DIM >
typedef void(* AdvectionDiffusionSourceFctPt) (const Vector< double > &x, double &f)
 
typedef void(* AdvectionDiffusionWindFctPt) (const Vector< double > &x, Vector< double > &wind)
 
- 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
 
- 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 DIM, unsigned NNODE_1D>
class oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D >

QSUPGAdvectionDiffusionElement<DIM,NNODE_1D> elements are SUPG-stabilised Advection Diffusion elements with NNODE_1D nodal points in each coordinate direction. Inherits from QAdvectionDiffusionElement and overwrites their test functions

Constructor & Destructor Documentation

◆ QSUPGAdvectionDiffusionElement()

template<unsigned DIM, unsigned NNODE_1D>
oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D >::QSUPGAdvectionDiffusionElement ( )
inline

Constructor: Call constructors for underlying QAdvectionDiffusion element. Initialise stabilisation parameter to zero

50  : QAdvectionDiffusionElement<DIM, NNODE_1D>()
51  {
52  Tau_SUPG = 0.0;
53  }
double Tau_SUPG
SUPG stabilisation parameter.
Definition: supg_advection_diffusion_elements.h:250

References oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D >::Tau_SUPG.

Member Function Documentation

◆ compute_stabilisation_parameter()

template<unsigned DIM, unsigned NNODE_1D>
void oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D >::compute_stabilisation_parameter ( )
inline

Compute stabilisation parameter for the element.

71  {
72  // Find out how many nodes there are
73  unsigned n_node = this->nnode();
74 
75  // Set up memory for the shape functions and their derivatives
76  Shape psi(n_node), test(n_node);
77  DShape dpsidx(n_node, DIM);
78 
79  // Evaluate everything at the element centroid
80  Vector<double> s(DIM, 0.0);
81 
82  // Call the geometrical shape functions and derivatives
83  (void)QElement<DIM, NNODE_1D>::dshape_eulerian(s, psi, dpsidx);
84 
85  // Calculate Eulerian coordinates
86  Vector<double> interpolated_x(DIM, 0.0);
87 
88  // Loop over nodes
89  for (unsigned l = 0; l < n_node; l++)
90  {
91  // Loop over directions (we're in 2D)
92  for (unsigned j = 0; j < DIM; j++)
93  {
94  interpolated_x[j] += this->nodal_position(l, j) * psi[l];
95  }
96  }
97 
98  // Element size: Choose the max. diagonal
99  double h = 0;
100  if (DIM == 1)
101  {
102  h = std::fabs(this->vertex_node_pt(1)->x(0) -
103  this->vertex_node_pt(0)->x(0));
104  }
105  else if (DIM == 2)
106  {
107  h =
108  pow(this->vertex_node_pt(3)->x(0) - this->vertex_node_pt(0)->x(0),
109  2) +
110  pow(this->vertex_node_pt(3)->x(1) - this->vertex_node_pt(0)->x(1), 2);
111  double h1 =
112  pow(this->vertex_node_pt(2)->x(0) - this->vertex_node_pt(1)->x(0),
113  2) +
114  pow(this->vertex_node_pt(2)->x(1) - this->vertex_node_pt(1)->x(1), 2);
115  if (h1 > h) h = h1;
116  h = sqrt(h);
117  }
118  else if (DIM == 3)
119  {
120  // diagonals are from nodes 0-7, 1-6, 2-5, 3-4
121  unsigned n1 = 0;
122  unsigned n2 = 7;
123  for (unsigned i = 0; i < 4; i++)
124  {
125  double h1 =
126  pow(this->vertex_node_pt(n1)->x(0) - this->vertex_node_pt(n2)->x(0),
127  2) +
128  pow(this->vertex_node_pt(n1)->x(1) - this->vertex_node_pt(n2)->x(1),
129  2) +
130  pow(this->vertex_node_pt(n1)->x(2) - this->vertex_node_pt(n2)->x(2),
131  2);
132  if (h1 > h) h = h1;
133  n1++;
134  n2--;
135  }
136  h = sqrt(h);
137  }
138 
139  // Get wind
140  Vector<double> wind(DIM);
141  // Dummy ipt argument?
142  unsigned ipt = 0;
143  this->get_wind_adv_diff(ipt, s, interpolated_x, wind);
144  double abs_wind = 0;
145  for (unsigned j = 0; j < DIM; j++)
146  {
147  abs_wind += wind[j] * wind[j];
148  }
149  abs_wind = sqrt(abs_wind);
150 
151  // Mesh Peclet number
152  double Pe_mesh = 0.5 * this->pe() * h * abs_wind;
153 
154  // Stabilisation parameter
155  if (Pe_mesh > 1.0)
156  {
157  Tau_SUPG = h / (2.0 * abs_wind) * (1.0 - 1.0 / Pe_mesh);
158  }
159  else
160  {
161  Tau_SUPG = 0.0;
162  }
163  }
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
virtual void get_wind_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Definition: advection_diffusion_elements.h:366
const double & pe() const
Peclet number.
Definition: advection_diffusion_elements.h:318
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
virtual Node * vertex_node_pt(const unsigned &j) const
Definition: elements.h:2500
double nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.h:2317
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625
squared absolute sa ArrayBase::abs2 DOXCOMMA MatrixBase::cwiseAbs2 sa Eigen::abs2 DOXCOMMA Eigen::pow DOXCOMMA ArrayBase::square nearest sa Eigen::floor DOXCOMMA Eigen::ceil DOXCOMMA ArrayBase::round nearest integer not less than the given sa Eigen::floor DOXCOMMA ArrayBase::ceil not a number test
Definition: GlobalFunctions.h:109
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117
list x
Definition: plotDoE.py:28
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References DIM, boost::multiprecision::fabs(), oomph::AdvectionDiffusionEquations< DIM >::get_wind_adv_diff(), i, oomph::FiniteElement::interpolated_x(), j, oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_position(), oomph::AdvectionDiffusionEquations< DIM >::pe(), Eigen::bfloat16_impl::pow(), s, sqrt(), oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D >::Tau_SUPG, Eigen::test, oomph::FiniteElement::vertex_node_pt(), and plotDoE::x.

◆ dshape_and_dtest_eulerian_adv_diff()

template<unsigned DIM, unsigned NNODE_1D>
double oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D >::dshape_and_dtest_eulerian_adv_diff ( const Vector< double > &  s,
Shape psi,
DShape dpsidx,
Shape test,
DShape dtestdx 
) const
protectedvirtual

Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.

QSUPGAdvectionDiffusionElement<DIM,NNODE_1D> elements are SUPG-stabilised Advection Diffusion elements with NNODE_1D nodal points in each coordinate direction. Inherits from QAdvectionDiffusionElement and overwrites their test functions Define the shape functions and test functions and derivatives w.r.t. global coordinates and return Jacobian of mapping.

SUPG stabilisation: Petrov-Galerkin, i.e. test functions \( \ne \) shape functions

Reimplemented from oomph::QAdvectionDiffusionElement< DIM, NNODE_1D >.

54  {
55  // Call the geometrical shape functions and derivatives
56  double J = QElement<DIM, NNODE_1D>::dshape_eulerian(s, psi, dpsidx);
57 
58  // Find out how many nodes there are
59  unsigned n_node = this->nnode();
60 
61  // Calculate Eulerian coordinates
62  Vector<double> interpolated_x(DIM, 0.0);
63 
64  // Loop over nodes
65  for (unsigned l = 0; l < n_node; l++)
66  {
67  // Loop over directions
68  for (unsigned j = 0; j < DIM; j++)
69  {
70  interpolated_x[j] += this->nodal_position(l, j) * psi[l];
71  }
72  }
73 
74  // Get wind
75  Vector<double> wind(DIM);
76  // Dummy ipt argument
77  unsigned ipt = 0;
78  this->get_wind_adv_diff(ipt, s, interpolated_x, wind);
79 
80  // Loop over the test functions and derivatives and set them equal to the
81  // shape functions + add stabilisation
82  for (unsigned j = 0; j < n_node; j++)
83  {
84  test[j] = psi[j];
85 
86  for (unsigned i = 0; i < DIM; i++)
87  {
88  dtestdx(j, i) = dpsidx(j, i);
89  test[j] += Tau_SUPG * wind[i] * dpsidx(j, i);
90  }
91  }
92 
93  // Return the jacobian
94  return J;
95  }
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
Definition: indexed_view.cpp:20

References DIM, i, J, j, and s.

◆ dshape_and_dtest_eulerian_at_knot_adv_diff()

template<unsigned DIM, unsigned NNODE_1D>
double oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D >::dshape_and_dtest_eulerian_at_knot_adv_diff ( const unsigned ipt,
Shape psi,
DShape dpsidx,
Shape test,
DShape dtestdx 
) const
protectedvirtual

Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.

Define the shape functions and test functions and derivatives w.r.t. global coordinates and return Jacobian of mapping.

SUPG stabilisation: Petrov-Galerkin, i.e. test functions \( \ne \) shape functions

Reimplemented from oomph::QAdvectionDiffusionElement< DIM, NNODE_1D >.

112  {
113  // Call the geometrical shape functions and derivatives
114  double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
115 
116  // Find out how many nodes there are
117  unsigned n_node = this->nnode();
118 
119  // Calculate Eulerian coordinates
120  Vector<double> interpolated_x(DIM, 0.0);
121 
122  // Loop over nodes
123  for (unsigned l = 0; l < n_node; l++)
124  {
125  // Loop over directions
126  for (unsigned j = 0; j < DIM; j++)
127  {
128  interpolated_x[j] += this->nodal_position(l, j) * psi(l);
129  }
130  }
131 
132  // Find the dimension of the element
133  unsigned Dim = this->dim();
134  // Storage for the local coordinates of the integration point
135  Vector<double> s(Dim);
136  // Set the local coordinate
137  for (unsigned i = 0; i < Dim; i++)
138  {
139  s[i] = this->integral_pt()->knot(ipt, i);
140  }
141 
142  // Get wind
143  Vector<double> wind(DIM);
144  this->get_wind_adv_diff(ipt, s, interpolated_x, wind);
145 
146  // Loop over the test functions and derivatives and set them equal to the
147  // shape functions + add stabilisation
148  for (unsigned j = 0; j < n_node; j++)
149  {
150  test(j) = psi(j);
151  for (unsigned i = 0; i < DIM; i++)
152  {
153  dtestdx(j, i) = dpsidx(j, i);
154  test(j) += Tau_SUPG * wind[i] * dpsidx(j, i);
155  }
156  }
157 
158 
159  // Return the jacobian
160  return J;
161  }
unsigned dim() const
Definition: elements.h:2611
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3325
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
static const unsigned Dim
Problem dimension.
Definition: two_d_tilted_square.cc:62

References Global_Variables::Dim, DIM, i, J, j, s, and Eigen::test.

◆ get_Tau_SUPG()

template<unsigned DIM, unsigned NNODE_1D>
double oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D >::get_Tau_SUPG ( )
inline

Get stabilisation parameter for the element.

57  {
58  return Tau_SUPG;
59  }

References oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D >::Tau_SUPG.

◆ output() [1/4]

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

C-style output.

Reimplemented from oomph::QAdvectionDiffusionElement< DIM, NNODE_1D >.

220  {
221  FiniteElement::output(file_pt);
222  }
virtual void output(std::ostream &outfile)
Definition: elements.h:3050

References oomph::FiniteElement::output().

◆ output() [2/4]

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

C_style output at n_plot points.

Reimplemented from oomph::QAdvectionDiffusionElement< DIM, NNODE_1D >.

226  {
227  FiniteElement::output(file_pt, n_plot);
228  }

References oomph::FiniteElement::output().

◆ output() [3/4]

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

Output at default number of plot points.

Reimplemented from oomph::QAdvectionDiffusionElement< DIM, NNODE_1D >.

214  {
215  FiniteElement::output(outfile);
216  }

References oomph::FiniteElement::output().

◆ output() [4/4]

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

Output function: x,y,u,w_x,w_y,tau_supg or x,y,z,u,w_x,w_y,w_z,tau_supg nplot points in each coordinate direction

Reimplemented from oomph::QAdvectionDiffusionElement< DIM, NNODE_1D >.

170  {
171  // Vector of local coordinates
172  Vector<double> s(DIM);
173 
174  // Tecplot header info
175  outfile << this->tecplot_zone_string(nplot);
176 
177  // Loop over plot points
178  unsigned num_plot_points = this->nplot_points(nplot);
179  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
180  {
181  // Get local coordinates of plot point
182  this->get_s_plot(iplot, nplot, s);
183 
184  // Get Eulerian coordinate of plot point
185  Vector<double> x(DIM);
186  this->interpolated_x(s, x);
187 
188  for (unsigned i = 0; i < DIM; i++)
189  {
190  outfile << x[i] << " ";
191  }
192  outfile << this->interpolated_u_adv_diff(s) << " ";
193 
194  // Get the wind
195  Vector<double> wind(DIM);
196  // Dummy ipt argument
197  unsigned ipt = 0;
198  this->get_wind_adv_diff(ipt, s, x, wind);
199  for (unsigned i = 0; i < DIM; i++)
200  {
201  outfile << wind[i] << " ";
202  }
203 
204  // Output stabilisation parameter
205  outfile << Tau_SUPG << std::endl;
206  }
207 
208  // Write tecplot footer (e.g. FE connectivity lists)
209  this->write_tecplot_zone_footer(outfile, nplot);
210  }
double interpolated_u_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: advection_diffusion_elements.h:458
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Definition: elements.h:3174

References DIM, oomph::FiniteElement::get_s_plot(), oomph::AdvectionDiffusionEquations< DIM >::get_wind_adv_diff(), i, oomph::AdvectionDiffusionEquations< DIM >::interpolated_u_adv_diff(), oomph::FiniteElement::interpolated_x(), oomph::FiniteElement::nplot_points(), s, oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D >::Tau_SUPG, oomph::FiniteElement::tecplot_zone_string(), oomph::FiniteElement::write_tecplot_zone_footer(), and plotDoE::x.

◆ switch_off_stabilisation()

template<unsigned DIM, unsigned NNODE_1D>
void oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D >::switch_off_stabilisation ( )
inline

Set stabilisation parameter for the element to zero.

64  {
65  Tau_SUPG = 0.0;
66  }

References oomph::QSUPGAdvectionDiffusionElement< DIM, NNODE_1D >::Tau_SUPG.

Member Data Documentation

◆ Tau_SUPG


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