oomph::AxisymAdvectionDiffusionEquations Class Referenceabstract

#include <axisym_advection_diffusion_elements.h>

+ Inheritance diagram for oomph::AxisymAdvectionDiffusionEquations:

Public Types

typedef void(* AxisymAdvectionDiffusionSourceFctPt) (const Vector< double > &x, double &f)
 
typedef void(* AxisymAdvectionDiffusionWindFctPt) (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 > &)
 

Public Member Functions

 AxisymAdvectionDiffusionEquations ()
 
 AxisymAdvectionDiffusionEquations (const AxisymAdvectionDiffusionEquations &dummy)=delete
 Broken copy constructor. More...
 
virtual unsigned u_index_axi_adv_diff () const
 Broken assignment operator. More...
 
double du_dt_axi_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 output (std::ostream &outfile)
 Output with default number of plot points. More...
 
void output (std::ostream &outfile, const unsigned &nplot)
 
void output (FILE *file_pt)
 C_style output with default number of plot points. More...
 
void output (FILE *file_pt, const unsigned &n_plot)
 
void output_fct (std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 Output exact soln: r,z,u_exact at nplot^2 plot points. More...
 
void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 Get error against and norm of exact solution. More...
 
AxisymAdvectionDiffusionSourceFctPtsource_fct_pt ()
 Access function: Pointer to source function. More...
 
AxisymAdvectionDiffusionSourceFctPt source_fct_pt () const
 Access function: Pointer to source function. Const version. More...
 
AxisymAdvectionDiffusionWindFctPtwind_fct_pt ()
 Access function: Pointer to wind function. More...
 
AxisymAdvectionDiffusionWindFctPt 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...
 
const doubled () const
 Peclet number multiplied by Strouhal number. More...
 
double *& d_pt ()
 Pointer to Peclet number multipled by Strouha number. More...
 
virtual void get_source_axi_adv_diff (const unsigned &ipt, const Vector< double > &x, double &source) const
 
virtual void get_wind_axi_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: [du/dr,du/dz]. 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)
 
double interpolated_u_axi_adv_diff (const Vector< double > &s) const
 Return FE representation of function value u(s) at local coordinate s. More...
 
virtual void dinterpolated_u_axi_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)
 
virtual unsigned required_nvalue (const unsigned &n) const
 
unsigned nnodal_position_type () const
 
bool has_hanging_nodes () const
 
unsigned nodal_dimension () const
 Return the required Eulerian dimension of the nodes in this element. More...
 
virtual 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, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
 Output a time-dependent exact solution over the element. More...
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, const SolutionFunctorBase &exact_soln) const
 Output a time-dependent exact solution over the element. More...
 
virtual void 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::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_abs_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error)
 
void integrate_fct (FiniteElement::SteadyExactSolutionFctPt integrand_fct_pt, Vector< double > &integral)
 Integrate Vector-valued function over element. More...
 
void integrate_fct (FiniteElement::UnsteadyExactSolutionFctPt integrand_fct_pt, const double &time, Vector< double > &integral)
 Integrate Vector-valued time-dep function over element. More...
 
virtual void build_face_element (const int &face_index, FaceElement *face_element_pt)
 
virtual unsigned 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

virtual double dshape_and_dtest_eulerian_axi_adv_diff (const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
virtual double dshape_and_dtest_eulerian_at_knot_axi_adv_diff (const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
virtual void fill_in_generic_residual_contribution_axi_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_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void fill_in_contribution_to_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void fill_in_contribution_to_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void fill_in_contribution_to_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void fill_in_contribution_to_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void fill_in_contribution_to_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 

Protected Attributes

doublePe_pt
 Pointer to global Peclet number. More...
 
doublePeSt_pt
 Pointer to global Peclet number multiplied by Strouhal number. More...
 
doubleD_pt
 Pointer to global Diffusion parameter. More...
 
AxisymAdvectionDiffusionSourceFctPt Source_fct_pt
 Pointer to source function: More...
 
AxisymAdvectionDiffusionWindFctPt Wind_fct_pt
 Pointer to wind function: More...
 
bool ALE_is_disabled
 Boolean flag to indicate whether AlE formulation is disable. More...
 
- Protected Attributes inherited from oomph::FiniteElement
MacroElementMacro_elem_pt
 Pointer to the element's macro element (NULL by default) More...
 
- Protected Attributes inherited from oomph::GeomObject
unsigned NLagrangian
 Number of Lagrangian (intrinsic) coordinates. More...
 
unsigned Ndim
 Number of Eulerian coordinates. More...
 
TimeStepperGeom_object_time_stepper_pt
 

Static Private Attributes

static double Default_peclet_number = 0.0
 Static default value for the Peclet number. More...
 
static double Default_diffusion_parameter = 1.0
 Static default value for the Peclet number. More...
 

Additional Inherited Members

- 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

A class for all elements that solve the Advection Diffusion equations in a cylindrical polar coordinate system using isoparametric elements.

\[ Pe \mathbf{w}\cdot(\mathbf{x}) \nabla u = \nabla \cdot \left( \nabla u \right) + f(\mathbf{x}) \]

This contains the generic maths. Shape functions, geometric mapping etc. must get implemented in derived class.

Member Typedef Documentation

◆ AxisymAdvectionDiffusionSourceFctPt

typedef void(* oomph::AxisymAdvectionDiffusionEquations::AxisymAdvectionDiffusionSourceFctPt) (const Vector< double > &x, double &f)

Function pointer to source function fct(x,f(x)) – x is a Vector!

◆ AxisymAdvectionDiffusionWindFctPt

typedef void(* oomph::AxisymAdvectionDiffusionEquations::AxisymAdvectionDiffusionWindFctPt) (const Vector< double > &x, Vector< double > &wind)

Function pointer to wind function fct(x,w(x)) – x is a Vector!

Constructor & Destructor Documentation

◆ AxisymAdvectionDiffusionEquations() [1/2]

oomph::AxisymAdvectionDiffusionEquations::AxisymAdvectionDiffusionEquations ( )
inline

Constructor: Initialise the Source_fct_pt and Wind_fct_pt to null and set (pointer to) Peclet number to default

72  {
73  // Set pointer to Peclet number to the default value zero
77  }
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
Definition: axisym_advection_diffusion_elements.h:558
AxisymAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
Definition: axisym_advection_diffusion_elements.h:564
AxisymAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
Definition: axisym_advection_diffusion_elements.h:567
static double Default_peclet_number
Static default value for the Peclet number.
Definition: axisym_advection_diffusion_elements.h:574
static double Default_diffusion_parameter
Static default value for the Peclet number.
Definition: axisym_advection_diffusion_elements.h:577
bool ALE_is_disabled
Boolean flag to indicate whether AlE formulation is disable.
Definition: axisym_advection_diffusion_elements.h:570
double * D_pt
Pointer to global Diffusion parameter.
Definition: axisym_advection_diffusion_elements.h:561
double * Pe_pt
Pointer to global Peclet number.
Definition: axisym_advection_diffusion_elements.h:555

References D_pt, Default_diffusion_parameter, Default_peclet_number, Pe_pt, and PeSt_pt.

◆ AxisymAdvectionDiffusionEquations() [2/2]

oomph::AxisymAdvectionDiffusionEquations::AxisymAdvectionDiffusionEquations ( const AxisymAdvectionDiffusionEquations dummy)
delete

Broken copy constructor.

Member Function Documentation

◆ compute_error()

void oomph::AxisymAdvectionDiffusionEquations::compute_error ( std::ostream &  outfile,
FiniteElement::SteadyExactSolutionFctPt  exact_soln_pt,
double error,
double norm 
)
virtual

Get error against and norm of exact solution.

Validate against exact solution

Solution is provided via function pointer. Plot error at a given number of plot points.

Reimplemented from oomph::FiniteElement.

420  {
421  // Initialise
422  error = 0.0;
423  norm = 0.0;
424 
425  // Vector of local coordinates
426  Vector<double> s(2);
427 
428  // Vector for coordintes
429  Vector<double> x(2);
430 
431  // Find out how many nodes there are in the element
432  unsigned n_node = nnode();
433 
434  Shape psi(n_node);
435 
436  // Set the value of n_intpt
437  unsigned n_intpt = integral_pt()->nweight();
438 
439  // Tecplot header info
440  outfile << "ZONE" << std::endl;
441 
442  // Exact solution Vector (here a scalar)
443  Vector<double> exact_soln(1);
444 
445  // Loop over the integration points
446  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
447  {
448  // Assign values of s
449  for (unsigned i = 0; i < 2; i++)
450  {
451  s[i] = integral_pt()->knot(ipt, i);
452  }
453 
454  // Get the integral weight
455  double w = integral_pt()->weight(ipt);
456 
457  // Get jacobian of mapping
458  double J = J_eulerian(s);
459 
460  // Premultiply the weights and the Jacobian
461  double W = w * J;
462 
463  // Get x position as Vector
464  interpolated_x(s, x);
465 
466  // Get FE function value
467  double u_fe = interpolated_u_axi_adv_diff(s);
468 
469  // Get exact solution at this point
470  (*exact_soln_pt)(x, exact_soln);
471 
472  // Output x,y,...,error
473  for (unsigned i = 0; i < 2; i++)
474  {
475  outfile << x[i] << " ";
476  }
477  outfile << exact_soln[0] << " " << exact_soln[0] - u_fe << std::endl;
478 
479  // Add to error and norm
480  norm += exact_soln[0] * exact_soln[0] * x[0] * W;
481  error += (exact_soln[0] - u_fe) * (exact_soln[0] - u_fe) * x[0] * W;
482  }
483  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
RowVector3d w
Definition: Matrix_resize_int.cpp:3
double interpolated_u_axi_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: axisym_advection_diffusion_elements.h:437
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
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual double J_eulerian(const Vector< double > &s) const
Definition: elements.cc:4103
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
RealScalar s
Definition: level1_cplx_impl.h:130
void exact_soln(const double &time, const Vector< double > &x, Vector< double > &soln)
Definition: unstructured_two_d_curved.cc:301
int error
Definition: calibrate.py:297
@ W
Definition: quadtree.h:63
list x
Definition: plotDoE.py:28

References calibrate::error, ProblemParameters::exact_soln(), i, oomph::FiniteElement::integral_pt(), interpolated_u_axi_adv_diff(), oomph::FiniteElement::interpolated_x(), J, oomph::FiniteElement::J_eulerian(), oomph::Integral::knot(), oomph::FiniteElement::nnode(), oomph::Integral::nweight(), s, w, oomph::QuadTreeNames::W, oomph::Integral::weight(), and plotDoE::x.

◆ d()

const double& oomph::AxisymAdvectionDiffusionEquations::d ( ) const
inline

Peclet number multiplied by Strouhal number.

324  {
325  return *D_pt;
326  }

References D_pt.

Referenced by fill_in_generic_residual_contribution_axi_adv_diff().

◆ d_pt()

double*& oomph::AxisymAdvectionDiffusionEquations::d_pt ( )
inline

Pointer to Peclet number multipled by Strouha number.

330  {
331  return D_pt;
332  }

References D_pt.

◆ dinterpolated_u_axi_adv_diff_ddata()

virtual void oomph::AxisymAdvectionDiffusionEquations::dinterpolated_u_axi_adv_diff_ddata ( const Vector< double > &  s,
Vector< double > &  du_ddata,
Vector< unsigned > &  global_eqn_number 
)
inlinevirtual

Return derivative of u at point s with respect to all data that can affect its value. In addition, return the global equation numbers corresponding to the data. This is virtual so that it can be overloaded in the refineable version

473  {
474  // Find number of nodes
475  const unsigned n_node = nnode();
476 
477  // Get the nodal index at which the unknown is stored
478  const unsigned u_nodal_index = u_index_axi_adv_diff();
479 
480  // Local shape function
481  Shape psi(n_node);
482 
483  // Find values of shape function
484  shape(s, psi);
485 
486  // Find the number of dofs associated with interpolated u
487  unsigned n_u_dof = 0;
488  for (unsigned l = 0; l < n_node; l++)
489  {
490  int global_eqn = this->node_pt(l)->eqn_number(u_nodal_index);
491  // If it's positive add to the count
492  if (global_eqn >= 0)
493  {
494  ++n_u_dof;
495  }
496  }
497 
498  // Now resize the storage schemes
499  du_ddata.resize(n_u_dof, 0.0);
500  global_eqn_number.resize(n_u_dof, 0);
501 
502  // Loop over the nodes again and set the derivatives
503  unsigned count = 0;
504  for (unsigned l = 0; l < n_node; l++)
505  {
506  // Get the global equation number
507  int global_eqn = this->node_pt(l)->eqn_number(u_nodal_index);
508  // If it's positive
509  if (global_eqn >= 0)
510  {
511  // Set the global equation number
512  global_eqn_number[count] = global_eqn;
513  // Set the derivative with respect to the unknown
514  du_ddata[count] = psi[l];
515  // Increase the counter
516  ++count;
517  }
518  }
519  }
virtual unsigned u_index_axi_adv_diff() const
Broken assignment operator.
Definition: axisym_advection_diffusion_elements.h:98
long & eqn_number(const unsigned &i)
Return the equation number of the i-th stored variable.
Definition: nodes.h:367
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
virtual void shape(const Vector< double > &s, Shape &psi) const =0

References oomph::Data::eqn_number(), oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), s, oomph::FiniteElement::shape(), and u_index_axi_adv_diff().

Referenced by RefineableQAxisymCrouzeixRaviartBoussinesqElement::get_dbody_force_axi_nst_dexternal_element_data(), and QAxisymCrouzeixRaviartElementWithExternalElement::get_dbody_force_axi_nst_dexternal_element_data().

◆ disable_ALE()

void oomph::AxisymAdvectionDiffusionEquations::disable_ALE ( )
inlinevirtual

Disable ALE, i.e. assert the mesh is not moving – you do this at your own risk!

Reimplemented from oomph::FiniteElement.

134  {
135  ALE_is_disabled = true;
136  }

References ALE_is_disabled.

Referenced by oomph::RefineableBuoyantQAxisymCrouzeixRaviartElement::disable_ALE().

◆ dshape_and_dtest_eulerian_at_knot_axi_adv_diff()

virtual double oomph::AxisymAdvectionDiffusionEquations::dshape_and_dtest_eulerian_at_knot_axi_adv_diff ( const unsigned ipt,
Shape psi,
DShape dpsidx,
Shape test,
DShape dtestdx 
) const
protectedpure virtual

◆ dshape_and_dtest_eulerian_axi_adv_diff()

virtual double oomph::AxisymAdvectionDiffusionEquations::dshape_and_dtest_eulerian_axi_adv_diff ( const Vector< double > &  s,
Shape psi,
DShape dpsidx,
Shape test,
DShape dtestdx 
) const
protectedpure virtual

Shape/test functions and derivs w.r.t. to global coords at local coord. s; return Jacobian of mapping

Implemented in oomph::QAxisymAdvectionDiffusionElement< NNODE_1D >, and oomph::QAxisymAdvectionDiffusionElement< 3 >.

◆ du_dt_axi_adv_diff()

double oomph::AxisymAdvectionDiffusionEquations::du_dt_axi_adv_diff ( const unsigned n) const
inline

du/dt at local node n. Uses suitably interpolated value for hanging nodes.

107  {
108  // Get the data's timestepper
109  TimeStepper* time_stepper_pt = this->node_pt(n)->time_stepper_pt();
110 
111  // Initialise dudt
112  double dudt = 0.0;
113  // Loop over the timesteps, if there is a non Steady timestepper
114  if (!time_stepper_pt->is_steady())
115  {
116  // Find the index at which the variable is stored
117  const unsigned u_nodal_index = u_index_axi_adv_diff();
118 
119  // Number of timsteps (past & present)
120  const unsigned n_time = time_stepper_pt->ntstorage();
121 
122  for (unsigned t = 0; t < n_time; t++)
123  {
124  dudt +=
125  time_stepper_pt->weight(1, t) * nodal_value(t, n, u_nodal_index);
126  }
127  }
128  return dudt;
129  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
TimeStepper *& time_stepper_pt()
Definition: geom_objects.h:192
unsigned ntstorage() const
Definition: timesteppers.h:601
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
Definition: timesteppers.h:594
bool is_steady() const
Definition: timesteppers.h:389
t
Definition: plotPSD.py:36

References oomph::TimeStepper::is_steady(), n, oomph::FiniteElement::nodal_value(), oomph::FiniteElement::node_pt(), oomph::TimeStepper::ntstorage(), plotPSD::t, oomph::GeomObject::time_stepper_pt(), oomph::Data::time_stepper_pt(), u_index_axi_adv_diff(), and oomph::TimeStepper::weight().

Referenced by fill_in_generic_residual_contribution_axi_adv_diff(), and oomph::RefineableAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_axi_adv_diff().

◆ enable_ALE()

void oomph::AxisymAdvectionDiffusionEquations::enable_ALE ( )
inlinevirtual

(Re-)enable ALE, i.e. take possible mesh motion into account when evaluating the time-derivative. Note: By default, ALE is enabled, at the expense of possibly creating unnecessary work in problems where the mesh is, in fact, stationary.

Reimplemented from oomph::FiniteElement.

144  {
145  ALE_is_disabled = false;
146  }

References ALE_is_disabled.

Referenced by oomph::AxisymmetricQAdvectionCrouzeixRaviartElement::disable_ALE(), oomph::RefineableBuoyantQAxisymCrouzeixRaviartElement::enable_ALE(), and oomph::AxisymmetricQAdvectionCrouzeixRaviartElement::enable_ALE().

◆ fill_in_contribution_to_jacobian()

void oomph::AxisymAdvectionDiffusionEquations::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
inlinevirtual

Add the element's contribution to its residual vector and the element Jacobian matrix (wrapper)

Reimplemented from oomph::GeneralisedElement.

429  {
430  // Call the generic routine with the flag set to 1
432  residuals, jacobian, GeneralisedElement::Dummy_matrix, 1);
433  }
virtual void fill_in_generic_residual_contribution_axi_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: axisym_advection_diffusion_elements.cc:50
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227

References oomph::GeneralisedElement::Dummy_matrix, and fill_in_generic_residual_contribution_axi_adv_diff().

Referenced by oomph::RefineableBuoyantQAxisymCrouzeixRaviartElement::fill_in_contribution_to_jacobian(), RefineableQAxisymAdvectionDiffusionBoussinesqElement::fill_in_contribution_to_jacobian(), QAxisymAdvectionDiffusionElementWithExternalElement::fill_in_contribution_to_jacobian(), and oomph::AxisymmetricQAdvectionCrouzeixRaviartElement::fill_in_contribution_to_jacobian().

◆ fill_in_contribution_to_residuals()

void oomph::AxisymAdvectionDiffusionEquations::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inlinevirtual

◆ fill_in_generic_residual_contribution_axi_adv_diff()

void oomph::AxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_axi_adv_diff ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  mass_matrix,
unsigned  flag 
)
protectedvirtual

Add the element's contribution to its residual vector only (if flag=and/or element Jacobian matrix

Compute element residual Vector and/or element Jacobian matrix

flag=1: compute both flag=0: compute only residual Vector

Pure version without hanging nodes

Reimplemented in oomph::RefineableAxisymAdvectionDiffusionEquations.

55  {
56  // Find out how many nodes there are
57  const unsigned n_node = nnode();
58 
59  // Get the nodal index at which the unknown is stored
60  const unsigned u_nodal_index = u_index_axi_adv_diff();
61 
62  // Set up memory for the shape and test functions
63  Shape psi(n_node), test(n_node);
64  DShape dpsidx(n_node, 2), dtestdx(n_node, 2);
65 
66  // Set the value of n_intpt
67  const unsigned n_intpt = integral_pt()->nweight();
68 
69  // Set the Vector to hold local coordinates
70  Vector<double> s(2);
71 
72  // Get Physical Variables from Element
73  const double scaled_peclet = pe();
74 
75  // Get peclet number multiplied by Strouhal number
76  const double scaled_peclet_st = pe_st();
77 
78  // Integers used to store the local equation number and local unknown
79  // indices for the residuals and jacobians
80  int local_eqn = 0, local_unknown = 0;
81 
82  // Loop over the integration points
83  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
84  {
85  // Assign values of s
86  for (unsigned i = 0; i < 2; i++) s[i] = integral_pt()->knot(ipt, i);
87 
88  // Get the integral weight
89  double w = integral_pt()->weight(ipt);
90 
91  // Call the derivatives of the shape and test functions
93  ipt, psi, dpsidx, test, dtestdx);
94 
95  // Premultiply the weights and the Jacobian
96  double W = w * J;
97 
98  // Calculate local values of the solution and its derivatives
99  // Allocate
100  double interpolated_u = 0.0;
101  double dudt = 0.0;
102  Vector<double> interpolated_x(2, 0.0);
103  Vector<double> interpolated_dudx(2, 0.0);
104  Vector<double> mesh_velocity(2, 0.0);
105 
106  // Calculate function value and derivatives:
107  //-----------------------------------------
108  // Loop over nodes
109  for (unsigned l = 0; l < n_node; l++)
110  {
111  // Get the value at the node
112  double u_value = raw_nodal_value(l, u_nodal_index);
113  interpolated_u += u_value * psi(l);
114  dudt += du_dt_axi_adv_diff(l) * psi(l);
115  // Loop over the two coordinates directions
116  for (unsigned j = 0; j < 2; j++)
117  {
118  interpolated_x[j] += raw_nodal_position(l, j) * psi(l);
119  interpolated_dudx[j] += u_value * dpsidx(l, j);
120  }
121  }
122 
123  if (!ALE_is_disabled)
124  {
125  for (unsigned l = 0; l < n_node; l++)
126  {
127  for (unsigned j = 0; j < 2; j++)
128  {
129  mesh_velocity[j] += raw_dnodal_position_dt(l, j) * psi(l);
130  }
131  }
132  }
133 
134 
135  // Get source function
136  //-------------------
137  double source;
139 
140 
141  // Get wind three potential components
142  //---------------------------------------
143  Vector<double> wind(3);
145 
146  // Get the diffusion
147  double diff = this->d();
148 
149  // r is the first position component
150  double r = interpolated_x[0];
151 
152  // TEMPERATURE EQUATION (Neglected viscous dissipation)
153  // Assemble residuals and Jacobian
154  //--------------------------------
155 
156  // Loop over the test functions
157  for (unsigned l = 0; l < n_node; l++)
158  {
159  // Set the local equation number
160  local_eqn = nodal_local_eqn(l, u_nodal_index);
161 
162  /*IF it's not a boundary condition*/
163  if (local_eqn >= 0)
164  {
165  // Add body force/source term
166  residuals[local_eqn] -=
167  (scaled_peclet_st * dudt + source) * r * test(l) * W;
168 
169  // The Advection Diffusion bit itself
170  residuals[local_eqn] -=
171  // radial terms
172  (interpolated_dudx[0] *
173  (scaled_peclet * wind[0] * test(l) + diff * dtestdx(l, 0)) +
174  // azimuthal terms
175  (interpolated_dudx[1] *
176  (scaled_peclet * wind[1] * test(l) + diff * dtestdx(l, 1)))) *
177  r * W;
178 
179  if (!ALE_is_disabled)
180  {
181  residuals[local_eqn] += scaled_peclet_st *
182  (mesh_velocity[0] * interpolated_dudx[0] +
183  mesh_velocity[1] * interpolated_dudx[1]) *
184  test(l) * r * W;
185  }
186 
187  // Calculate the jacobian
188  //-----------------------
189  if (flag)
190  {
191  // Loop over the velocity shape functions again
192  for (unsigned l2 = 0; l2 < n_node; l2++)
193  {
194  // Set the number of the unknown
195  local_unknown = nodal_local_eqn(l2, u_nodal_index);
196 
197  // If at a non-zero degree of freedom add in the entry
198  if (local_unknown >= 0)
199  {
200  // Mass matrix term
201  jacobian(local_eqn, local_unknown) -=
202  scaled_peclet_st * test(l) * psi(l2) *
203  node_pt(l2)->time_stepper_pt()->weight(1, 0) * r * W;
204 
205  // add the mass matrix term
206  if (flag == 2)
207  {
208  mass_matrix(local_eqn, local_unknown) +=
209  scaled_peclet_st * test(l) * psi(l2) * r * W;
210  }
211 
212  // Assemble Jacobian term
213  jacobian(local_eqn, local_unknown) -=
214  // radial terms
215  (dpsidx(l2, 0) * (scaled_peclet * wind[0] * test(l) +
216  diff * dtestdx(l, 0)) +
217  // azimuthal terms
218  (dpsidx(l2, 1) * (scaled_peclet * wind[1] * test(l) +
219  diff * dtestdx(l, 1)))) *
220  r * W;
221 
222  if (!ALE_is_disabled)
223  {
224  jacobian(local_eqn, local_unknown) +=
225  scaled_peclet_st *
226  (mesh_velocity[0] * dpsidx(l2, 0) +
227  mesh_velocity[1] * dpsidx(l2, 1)) *
228  test(l) * r * W;
229  }
230  }
231  }
232  }
233  }
234  }
235 
236 
237  } // End of loop over integration points
238  }
virtual void get_wind_axi_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Definition: axisym_advection_diffusion_elements.h:359
double du_dt_axi_adv_diff(const unsigned &n) const
Definition: axisym_advection_diffusion_elements.h:106
const double & pe() const
Peclet number.
Definition: axisym_advection_diffusion_elements.h:299
virtual double dshape_and_dtest_eulerian_at_knot_axi_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
const double & d() const
Peclet number multiplied by Strouhal number.
Definition: axisym_advection_diffusion_elements.h:323
virtual void get_source_axi_adv_diff(const unsigned &ipt, const Vector< double > &x, double &source) const
Definition: axisym_advection_diffusion_elements.h:338
const double & pe_st() const
Peclet number multiplied by Strouhal number.
Definition: axisym_advection_diffusion_elements.h:311
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
double raw_dnodal_position_dt(const unsigned &n, const unsigned &i) const
Definition: elements.h:2256
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2576
double raw_nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.cc:1686
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
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
r
Definition: UniformPSDSelfTest.py:20
Definition: indexed_view.cpp:20
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References ALE_is_disabled, d(), dshape_and_dtest_eulerian_at_knot_axi_adv_diff(), du_dt_axi_adv_diff(), get_source_axi_adv_diff(), get_wind_axi_adv_diff(), i, oomph::FiniteElement::integral_pt(), oomph::FiniteElement::interpolated_x(), J, j, oomph::Integral::knot(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_local_eqn(), oomph::FiniteElement::node_pt(), oomph::Integral::nweight(), pe(), pe_st(), UniformPSDSelfTest::r, oomph::FiniteElement::raw_dnodal_position_dt(), oomph::FiniteElement::raw_nodal_position(), oomph::FiniteElement::raw_nodal_value(), s, TestProblem::source(), Eigen::test, oomph::Data::time_stepper_pt(), u_index_axi_adv_diff(), w, oomph::QuadTreeNames::W, oomph::Integral::weight(), and oomph::TimeStepper::weight().

Referenced by fill_in_contribution_to_jacobian(), and fill_in_contribution_to_residuals().

◆ get_flux()

void oomph::AxisymAdvectionDiffusionEquations::get_flux ( const Vector< double > &  s,
Vector< double > &  flux 
) const
inline

Get flux: [du/dr,du/dz].

381  {
382  // Find out how many nodes there are in the element
383  const unsigned n_node = nnode();
384 
385  // Get the nodal index at which the unknown is stored
386  const unsigned u_nodal_index = u_index_axi_adv_diff();
387 
388  // Set up memory for the shape and test functions
389  Shape psi(n_node);
390  DShape dpsidx(n_node, 2);
391 
392  // Call the derivatives of the shape and test functions
393  dshape_eulerian(s, psi, dpsidx);
394 
395  // Initialise to zero
396  for (unsigned j = 0; j < 2; j++)
397  {
398  flux[j] = 0.0;
399  }
400 
401  // Loop over nodes
402  for (unsigned l = 0; l < n_node; l++)
403  {
404  const double u_value = this->nodal_value(l, u_nodal_index);
405  // Add in the derivative directions
406  flux[0] += u_value * dpsidx(l, 0);
407  flux[1] += u_value * dpsidx(l, 1);
408  }
409  }
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3298
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59

References oomph::FiniteElement::dshape_eulerian(), ProblemParameters::flux(), j, oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), s, and u_index_axi_adv_diff().

Referenced by oomph::RefineableAxisymAdvectionDiffusionEquations::get_Z2_flux().

◆ get_source_axi_adv_diff()

virtual void oomph::AxisymAdvectionDiffusionEquations::get_source_axi_adv_diff ( const unsigned ipt,
const Vector< double > &  x,
double source 
) const
inlinevirtual

Get source term at (Eulerian) position x. This function is virtual to allow overloading in multi-physics problems where the strength of the source function might be determined by another system of equations

341  {
342  // If no source function has been set, return zero
343  if (Source_fct_pt == 0)
344  {
345  source = 0.0;
346  }
347  else
348  {
349  // Get source strength
350  (*Source_fct_pt)(x, source);
351  }
352  }

References TestProblem::source(), Source_fct_pt, and plotDoE::x.

Referenced by fill_in_generic_residual_contribution_axi_adv_diff(), and oomph::RefineableAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_axi_adv_diff().

◆ get_wind_axi_adv_diff()

virtual void oomph::AxisymAdvectionDiffusionEquations::get_wind_axi_adv_diff ( const unsigned ipt,
const Vector< double > &  s,
const Vector< double > &  x,
Vector< double > &  wind 
) const
inlinevirtual

Get wind at (Eulerian) position x and/or local coordinate s. This function is virtual to allow overloading in multi-physics problems where the wind function might be determined by another system of equations

Reimplemented in oomph::AxisymmetricQAdvectionCrouzeixRaviartElement, QAxisymAdvectionDiffusionElementWithExternalElement, RefineableQAxisymAdvectionDiffusionBoussinesqElement, oomph::RefineableBuoyantQAxisymCrouzeixRaviartElement, QAxisymAdvectionDiffusionElementWithExternalElement, RefineableQAxisymAdvectionDiffusionBoussinesqElement, and oomph::RefineableBuoyantQAxisymCrouzeixRaviartElement.

363  {
364  // If no wind function has been set, return zero
365  if (Wind_fct_pt == 0)
366  {
367  for (unsigned i = 0; i < 3; i++)
368  {
369  wind[i] = 0.0;
370  }
371  }
372  else
373  {
374  // Get wind
375  (*Wind_fct_pt)(x, wind);
376  }
377  }

References i, Wind_fct_pt, and plotDoE::x.

Referenced by fill_in_generic_residual_contribution_axi_adv_diff(), oomph::RefineableAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_axi_adv_diff(), output(), and scalar_value_paraview().

◆ interpolated_u_axi_adv_diff()

double oomph::AxisymAdvectionDiffusionEquations::interpolated_u_axi_adv_diff ( const Vector< double > &  s) const
inline

Return FE representation of function value u(s) at local coordinate s.

438  {
439  // Find number of nodes
440  const unsigned n_node = nnode();
441 
442  // Get the nodal index at which the unknown is stored
443  const unsigned u_nodal_index = u_index_axi_adv_diff();
444 
445  // Local shape function
446  Shape psi(n_node);
447 
448  // Find values of shape function
449  shape(s, psi);
450 
451  // Initialise value of u
452  double interpolated_u = 0.0;
453 
454  // Loop over the local nodes and sum
455  for (unsigned l = 0; l < n_node; l++)
456  {
457  interpolated_u += nodal_value(l, u_nodal_index) * psi[l];
458  }
459 
460  return (interpolated_u);
461  }

References oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), s, oomph::FiniteElement::shape(), and u_index_axi_adv_diff().

Referenced by compute_error(), oomph::RefineableBuoyantQAxisymCrouzeixRaviartElement::get_body_force_axi_nst(), oomph::RefineableBuoyantQAxisymCrouzeixRaviartElement::output(), oomph::AxisymmetricQAdvectionCrouzeixRaviartElement::output(), output(), and scalar_value_paraview().

◆ nscalar_paraview()

unsigned oomph::AxisymAdvectionDiffusionEquations::nscalar_paraview ( ) const
inlinevirtual

Number of scalars/fields output by this element. Reimplements broken virtual function in base class.

Reimplemented from oomph::FiniteElement.

152  {
153  return 4;
154  }

Referenced by oomph::AxisymmetricQAdvectionCrouzeixRaviartElement::nscalar_paraview().

◆ output() [1/4]

void oomph::AxisymAdvectionDiffusionEquations::output ( FILE *  file_pt)
inlinevirtual

C_style output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QAxisymAdvectionDiffusionElement< NNODE_1D >, and oomph::QAxisymAdvectionDiffusionElement< 3 >.

248  {
249  unsigned n_plot = 5;
250  output(file_pt, n_plot);
251  }
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: axisym_advection_diffusion_elements.h:235

References output().

◆ output() [2/4]

void oomph::AxisymAdvectionDiffusionEquations::output ( FILE *  file_pt,
const unsigned nplot 
)
virtual

C-style output FE representation of soln: r,z,u at n_plot^2 plot points

C-style output function:

r,z,u

nplot points in each coordinate direction

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QAxisymAdvectionDiffusionElement< NNODE_1D >, and oomph::QAxisymAdvectionDiffusionElement< 3 >.

329  {
330  // Vector of local coordinates
331  Vector<double> s(2);
332 
333  // Tecplot header info
334  fprintf(file_pt, "%s", tecplot_zone_string(nplot).c_str());
335 
336  // Loop over plot points
337  unsigned num_plot_points = nplot_points(nplot);
338  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
339  {
340  // Get local coordinates of plot point
341  get_s_plot(iplot, nplot, s);
342 
343  for (unsigned i = 0; i < 2; i++)
344  {
345  fprintf(file_pt, "%g ", interpolated_x(s, i));
346  }
347 
348  fprintf(file_pt, "%g \n", interpolated_u_axi_adv_diff(s));
349  }
350 
351  // Write tecplot footer (e.g. FE connectivity lists)
352  write_tecplot_zone_footer(file_pt, nplot);
353  }
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 oomph::FiniteElement::get_s_plot(), i, interpolated_u_axi_adv_diff(), oomph::FiniteElement::interpolated_x(), oomph::FiniteElement::nplot_points(), s, oomph::FiniteElement::tecplot_zone_string(), and oomph::FiniteElement::write_tecplot_zone_footer().

◆ output() [3/4]

void oomph::AxisymAdvectionDiffusionEquations::output ( std::ostream &  outfile)
inlinevirtual

Output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QAxisymAdvectionDiffusionElement< NNODE_1D >, and oomph::QAxisymAdvectionDiffusionElement< 3 >.

236  {
237  unsigned nplot = 5;
238  output(outfile, nplot);
239  }

Referenced by output(), and oomph::QAxisymAdvectionDiffusionElement< NNODE_1D >::output().

◆ output() [4/4]

void oomph::AxisymAdvectionDiffusionEquations::output ( std::ostream &  outfile,
const unsigned nplot 
)
virtual

Output FE representation of soln: r,z,u at nplot^2 plot points

Output function:

r,z,u,w_r,w_z

nplot points in each coordinate direction

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QAxisymAdvectionDiffusionElement< NNODE_1D >, and oomph::QAxisymAdvectionDiffusionElement< 3 >.

276  {
277  // Vector of local coordinates
278  Vector<double> s(2);
279 
280  // Tecplot header info
281  outfile << tecplot_zone_string(nplot);
282 
283  // Loop over plot points
284  unsigned num_plot_points = nplot_points(nplot);
285  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
286  {
287  // Get local coordinates of plot point
288  get_s_plot(iplot, nplot, s);
289 
290  // Get Eulerian coordinate of plot point
291  Vector<double> x(2);
292  interpolated_x(s, x);
293 
294  for (unsigned i = 0; i < 2; i++)
295  {
296  outfile << x[i] << " ";
297  }
298 
299  // Output concentration
300  outfile << this->interpolated_u_axi_adv_diff(s) << " ";
301 
302  // Get the wind
303  Vector<double> wind(3);
304  // Dummy ipt argument needed... ?
305  unsigned ipt = 0;
306  get_wind_axi_adv_diff(ipt, s, x, wind);
307  for (unsigned i = 0; i < 3; i++)
308  {
309  outfile << wind[i] << " ";
310  }
311  outfile << std::endl;
312  }
313 
314  // Write tecplot footer (e.g. FE connectivity lists)
315  write_tecplot_zone_footer(outfile, nplot);
316  }

References oomph::FiniteElement::get_s_plot(), get_wind_axi_adv_diff(), i, interpolated_u_axi_adv_diff(), oomph::FiniteElement::interpolated_x(), oomph::FiniteElement::nplot_points(), s, oomph::FiniteElement::tecplot_zone_string(), oomph::FiniteElement::write_tecplot_zone_footer(), and plotDoE::x.

◆ output_fct()

void oomph::AxisymAdvectionDiffusionEquations::output_fct ( std::ostream &  outfile,
const unsigned nplot,
FiniteElement::SteadyExactSolutionFctPt  exact_soln_pt 
)
virtual

Output exact soln: r,z,u_exact at nplot^2 plot points.

Output exact solution

Solution is provided via function pointer. Plot at a given number of plot points.

r,z,u_exact

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QAxisymAdvectionDiffusionElement< NNODE_1D >, and oomph::QAxisymAdvectionDiffusionElement< 3 >.

368  {
369  // Vector of local coordinates
370  Vector<double> s(2);
371 
372  // Vector for coordintes
373  Vector<double> x(2);
374 
375  // Tecplot header info
376  outfile << tecplot_zone_string(nplot);
377 
378  // Exact solution Vector (here a scalar)
379  Vector<double> exact_soln(1);
380 
381  // Loop over plot points
382  unsigned num_plot_points = nplot_points(nplot);
383  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
384  {
385  // Get local coordinates of plot point
386  get_s_plot(iplot, nplot, s);
387 
388  // Get x position as Vector
389  interpolated_x(s, x);
390 
391  // Get exact solution at this point
392  (*exact_soln_pt)(x, exact_soln);
393 
394  // Output x,y,...,u_exact
395  for (unsigned i = 0; i < 2; i++)
396  {
397  outfile << x[i] << " ";
398  }
399  outfile << exact_soln[0] << std::endl;
400  }
401 
402  // Write tecplot footer (e.g. FE connectivity lists)
403  write_tecplot_zone_footer(outfile, nplot);
404  }

References ProblemParameters::exact_soln(), oomph::FiniteElement::get_s_plot(), i, oomph::FiniteElement::interpolated_x(), oomph::FiniteElement::nplot_points(), s, oomph::FiniteElement::tecplot_zone_string(), oomph::FiniteElement::write_tecplot_zone_footer(), and plotDoE::x.

Referenced by oomph::QAxisymAdvectionDiffusionElement< NNODE_1D >::output_fct().

◆ pe()

◆ pe_pt()

double*& oomph::AxisymAdvectionDiffusionEquations::pe_pt ( )
inline

Pointer to Peclet number.

306  {
307  return Pe_pt;
308  }

References Pe_pt.

Referenced by oomph::RefineableAxisymAdvectionDiffusionEquations::further_build().

◆ pe_st()

const double& oomph::AxisymAdvectionDiffusionEquations::pe_st ( ) const
inline

◆ pe_st_pt()

double*& oomph::AxisymAdvectionDiffusionEquations::pe_st_pt ( )
inline

Pointer to Peclet number multipled by Strouha number.

318  {
319  return PeSt_pt;
320  }

References PeSt_pt.

Referenced by oomph::RefineableAxisymAdvectionDiffusionEquations::further_build().

◆ scalar_name_paraview()

std::string oomph::AxisymAdvectionDiffusionEquations::scalar_name_paraview ( const unsigned i) const
inlinevirtual

Name of the i-th scalar field. Default implementation returns V1 for the first one, V2 for the second etc. Can (should!) be overloaded with more meaningful names in specific elements.

Reimplemented from oomph::FiniteElement.

210  {
211  // Winds
212  if (i < 3)
213  {
214  return "Wind " + StringConversion::to_string(i);
215  }
216  // Advection Diffusion field
217  else if (i == 3)
218  {
219  return "Advection Diffusion";
220  }
221  // Never get here
222  else
223  {
224  std::stringstream error_stream;
225  error_stream << "Advection Diffusion Elements only store 4 fields "
226  << std::endl;
227  throw OomphLibError(
228  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
229  // Dummy return
230  return " ";
231  }
232  }
std::string to_string(T object, unsigned float_precision=8)
Definition: oomph_utilities.h:189
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References i, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::StringConversion::to_string().

◆ scalar_value_paraview()

void oomph::AxisymAdvectionDiffusionEquations::scalar_value_paraview ( std::ofstream &  file_out,
const unsigned i,
const unsigned nplot 
) const
inlinevirtual

Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specific element type.

Reimplemented from oomph::FiniteElement.

161  {
162  // Vector of local coordinates
163  Vector<double> s(2);
164 
165  // Loop over plot points
166  unsigned num_plot_points = nplot_points_paraview(nplot);
167  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
168  {
169  // Get local coordinates of plot point
170  get_s_plot(iplot, nplot, s);
171 
172  // Get Eulerian coordinate of plot point
173  Vector<double> x(2);
174  interpolated_x(s, x);
175 
176  // Winds
177  if (i < 3)
178  {
179  // Get the wind
180  Vector<double> wind(3);
181 
182  // Dummy ipt argument needed... ?
183  unsigned ipt = 0;
184  get_wind_axi_adv_diff(ipt, s, x, wind);
185 
186  file_out << wind[i] << std::endl;
187  }
188  // Advection Diffusion
189  else if (i == 3)
190  {
191  file_out << this->interpolated_u_axi_adv_diff(s) << std::endl;
192  }
193  // Never get here
194  else
195  {
196  std::stringstream error_stream;
197  error_stream << "Advection Diffusion Elements only store 4 fields "
198  << std::endl;
199  throw OomphLibError(error_stream.str(),
202  }
203  }
204  }
virtual unsigned nplot_points_paraview(const unsigned &nplot) const
Definition: elements.h:2862

References oomph::FiniteElement::get_s_plot(), get_wind_axi_adv_diff(), i, interpolated_u_axi_adv_diff(), oomph::FiniteElement::interpolated_x(), oomph::FiniteElement::nplot_points_paraview(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, s, and plotDoE::x.

Referenced by oomph::AxisymmetricQAdvectionCrouzeixRaviartElement::scalar_value_paraview().

◆ self_test()

unsigned oomph::AxisymAdvectionDiffusionEquations::self_test ( )
virtual

Self-test: Return 0 for OK.

Reimplemented from oomph::FiniteElement.

246  {
247  bool passed = true;
248 
249  // Check lower-level stuff
250  if (FiniteElement::self_test() != 0)
251  {
252  passed = false;
253  }
254 
255  // Return verdict
256  if (passed)
257  {
258  return 0;
259  }
260  else
261  {
262  return 1;
263  }
264  }
virtual unsigned self_test()
Definition: elements.cc:4440

References oomph::FiniteElement::self_test().

◆ source_fct_pt() [1/2]

AxisymAdvectionDiffusionSourceFctPt& oomph::AxisymAdvectionDiffusionEquations::source_fct_pt ( )
inline

Access function: Pointer to source function.

271  {
272  return Source_fct_pt;
273  }

References Source_fct_pt.

Referenced by oomph::RefineableAxisymAdvectionDiffusionEquations::further_build().

◆ source_fct_pt() [2/2]

AxisymAdvectionDiffusionSourceFctPt oomph::AxisymAdvectionDiffusionEquations::source_fct_pt ( ) const
inline

Access function: Pointer to source function. Const version.

278  {
279  return Source_fct_pt;
280  }

References Source_fct_pt.

◆ u_index_axi_adv_diff()

virtual unsigned oomph::AxisymAdvectionDiffusionEquations::u_index_axi_adv_diff ( ) const
inlinevirtual

Broken assignment operator.

Return the index at which the unknown value is stored. The default value, 0, is appropriate for single-physics problems, when there is only one variable, the value that satisfies the spherical advection-diffusion equation. In derived multi-physics elements, this function should be overloaded to reflect the chosen storage scheme. Note that these equations require that the unknown is always stored at the same index at each node.

Reimplemented in oomph::AxisymmetricQAdvectionCrouzeixRaviartElement, oomph::RefineableBuoyantQAxisymCrouzeixRaviartElement, and oomph::RefineableBuoyantQAxisymCrouzeixRaviartElement.

99  {
100  return 0;
101  }

Referenced by oomph::RefineableAxisymAdvectionDiffusionEquations::dinterpolated_u_adv_diff_ddata(), dinterpolated_u_axi_adv_diff_ddata(), du_dt_axi_adv_diff(), fill_in_generic_residual_contribution_axi_adv_diff(), oomph::RefineableAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_axi_adv_diff(), get_flux(), oomph::RefineableAxisymAdvectionDiffusionEquations::get_interpolated_values(), and interpolated_u_axi_adv_diff().

◆ wind_fct_pt() [1/2]

AxisymAdvectionDiffusionWindFctPt& oomph::AxisymAdvectionDiffusionEquations::wind_fct_pt ( )
inline

Access function: Pointer to wind function.

285  {
286  return Wind_fct_pt;
287  }

References Wind_fct_pt.

Referenced by oomph::RefineableAxisymAdvectionDiffusionEquations::further_build().

◆ wind_fct_pt() [2/2]

AxisymAdvectionDiffusionWindFctPt oomph::AxisymAdvectionDiffusionEquations::wind_fct_pt ( ) const
inline

Access function: Pointer to wind function. Const version.

292  {
293  return Wind_fct_pt;
294  }

References Wind_fct_pt.

Member Data Documentation

◆ ALE_is_disabled

bool oomph::AxisymAdvectionDiffusionEquations::ALE_is_disabled
protected

◆ D_pt

double* oomph::AxisymAdvectionDiffusionEquations::D_pt
protected

Pointer to global Diffusion parameter.

Referenced by AxisymAdvectionDiffusionEquations(), d(), and d_pt().

◆ Default_diffusion_parameter

double oomph::AxisymAdvectionDiffusionEquations::Default_diffusion_parameter = 1.0
staticprivate

Static default value for the Peclet number.

Default Diffusion coefficient.

Referenced by AxisymAdvectionDiffusionEquations().

◆ Default_peclet_number

double oomph::AxisymAdvectionDiffusionEquations::Default_peclet_number = 0.0
staticprivate

Static default value for the Peclet number.

Default value for Peclet number.

Referenced by AxisymAdvectionDiffusionEquations().

◆ Pe_pt

double* oomph::AxisymAdvectionDiffusionEquations::Pe_pt
protected

◆ PeSt_pt

double* oomph::AxisymAdvectionDiffusionEquations::PeSt_pt
protected

Pointer to global Peclet number multiplied by Strouhal number.

Referenced by AxisymAdvectionDiffusionEquations(), oomph::RefineableAxisymAdvectionDiffusionEquations::further_build(), pe_st(), and pe_st_pt().

◆ Source_fct_pt

AxisymAdvectionDiffusionSourceFctPt oomph::AxisymAdvectionDiffusionEquations::Source_fct_pt
protected

◆ Wind_fct_pt

AxisymAdvectionDiffusionWindFctPt oomph::AxisymAdvectionDiffusionEquations::Wind_fct_pt
protected

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