oomph::SteadyAxisymAdvectionDiffusionEquations Class Referenceabstract

#include <steady_axisym_advection_diffusion_elements.h>

+ Inheritance diagram for oomph::SteadyAxisymAdvectionDiffusionEquations:

Public Types

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

 SteadyAxisymAdvectionDiffusionEquations ()
 
 SteadyAxisymAdvectionDiffusionEquations (const SteadyAxisymAdvectionDiffusionEquations &dummy)=delete
 Broken copy constructor. More...
 
virtual unsigned u_index_axisym_adv_diff () const
 Broken assignment operator. More...
 
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...
 
SteadyAxisymAdvectionDiffusionSourceFctPtsource_fct_pt ()
 Access function: Pointer to source function. More...
 
SteadyAxisymAdvectionDiffusionSourceFctPt source_fct_pt () const
 Access function: Pointer to source function. Const version. More...
 
SteadyAxisymAdvectionDiffusionWindFctPtwind_fct_pt ()
 Access function: Pointer to wind function. More...
 
SteadyAxisymAdvectionDiffusionWindFctPt 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...
 
virtual void get_source_axisym_adv_diff (const unsigned &ipt, const Vector< double > &x, double &source) const
 
virtual void get_wind_axisym_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] = \nabla u = \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)
 
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)
 
virtual void disable_ALE ()
 
virtual void enable_ALE ()
 
virtual unsigned required_nvalue (const unsigned &n) const
 
unsigned nnodal_position_type () const
 
bool has_hanging_nodes () const
 
unsigned nodal_dimension () const
 Return the required Eulerian dimension of the nodes in this element. More...
 
virtual 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 unsigned nscalar_paraview () const
 
virtual void scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
 
virtual std::string scalar_name_paraview (const unsigned &i) const
 
virtual void output (const unsigned &t, std::ostream &outfile, const unsigned &n_plot) const
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, 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_adv_diff (const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
virtual double dshape_and_dtest_eulerian_at_knot_adv_diff (const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
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)
 
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...
 
SteadyAxisymAdvectionDiffusionSourceFctPt Source_fct_pt
 Pointer to source function: More...
 
SteadyAxisymAdvectionDiffusionWindFctPt Wind_fct_pt
 Pointer to wind function: 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...
 

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 Steady Axisymmetric Advection Diffusion equations 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

◆ SteadyAxisymAdvectionDiffusionSourceFctPt

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

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

◆ SteadyAxisymAdvectionDiffusionWindFctPt

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

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

Constructor & Destructor Documentation

◆ SteadyAxisymAdvectionDiffusionEquations() [1/2]

oomph::SteadyAxisymAdvectionDiffusionEquations::SteadyAxisymAdvectionDiffusionEquations ( )
inline

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

68  : Source_fct_pt(0), Wind_fct_pt(0)
69  {
70  // Set pointer to Peclet number to the default value zero
72  }
SteadyAxisymAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
Definition: steady_axisym_advection_diffusion_elements.h:400
static double Default_peclet_number
Static default value for the Peclet number.
Definition: steady_axisym_advection_diffusion_elements.h:407
SteadyAxisymAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
Definition: steady_axisym_advection_diffusion_elements.h:403
double * Pe_pt
Pointer to global Peclet number.
Definition: steady_axisym_advection_diffusion_elements.h:397

References Default_peclet_number, and Pe_pt.

◆ SteadyAxisymAdvectionDiffusionEquations() [2/2]

oomph::SteadyAxisymAdvectionDiffusionEquations::SteadyAxisymAdvectionDiffusionEquations ( const SteadyAxisymAdvectionDiffusionEquations dummy)
delete

Broken copy constructor.

Member Function Documentation

◆ compute_error()

void oomph::SteadyAxisymAdvectionDiffusionEquations::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.

360  {
361  // Initialise
362  error = 0.0;
363  norm = 0.0;
364 
365  // Vector of local coordinates
366  Vector<double> s(2);
367 
368  // Vector for coordintes
369  Vector<double> x(2);
370 
371  // Find out how many nodes there are in the element
372  unsigned n_node = nnode();
373 
374  Shape psi(n_node);
375 
376  // Set the value of n_intpt
377  unsigned n_intpt = integral_pt()->nweight();
378 
379  // Tecplot header info
380  outfile << "ZONE" << std::endl;
381 
382  // Exact solution Vector (here a scalar)
383  Vector<double> exact_soln(1);
384 
385  // Loop over the integration points
386  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
387  {
388  // Assign values of s
389  for (unsigned i = 0; i < 2; i++)
390  {
391  s[i] = integral_pt()->knot(ipt, i);
392  }
393 
394  // Get the integral weight
395  double w = integral_pt()->weight(ipt);
396 
397  // Get jacobian of mapping
398  double J = J_eulerian(s);
399 
400  // Premultiply the weights and the Jacobian
401  double W = w * J;
402 
403  // Get x position as Vector
404  interpolated_x(s, x);
405 
406  // Get FE function value
407  double u_fe = interpolated_u_adv_diff(s);
408 
409  // Get exact solution at this point
410  (*exact_soln_pt)(x, exact_soln);
411 
412  // Output x,y,...,error
413  for (unsigned i = 0; i < 2; i++)
414  {
415  outfile << x[i] << " ";
416  }
417  outfile << exact_soln[0] << " " << exact_soln[0] - u_fe << std::endl;
418 
419  // Add to error and norm
420  norm += exact_soln[0] * exact_soln[0] * W;
421  error += (exact_soln[0] - u_fe) * (exact_soln[0] - u_fe) * W;
422  }
423  }
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
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.
double interpolated_u_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: steady_axisym_advection_diffusion_elements.h:279
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_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.

◆ dinterpolated_u_adv_diff_ddata()

virtual void oomph::SteadyAxisymAdvectionDiffusionEquations::dinterpolated_u_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

315  {
316  // Find number of nodes
317  const unsigned n_node = nnode();
318 
319  // Get the nodal index at which the unknown is stored
320  const unsigned u_nodal_index = u_index_axisym_adv_diff();
321 
322  // Local shape function
323  Shape psi(n_node);
324 
325  // Find values of shape function
326  shape(s, psi);
327 
328  // Find the number of dofs associated with interpolated u
329  unsigned n_u_dof = 0;
330  for (unsigned l = 0; l < n_node; l++)
331  {
332  int global_eqn = this->node_pt(l)->eqn_number(u_nodal_index);
333  // If it's positive add to the count
334  if (global_eqn >= 0)
335  {
336  ++n_u_dof;
337  }
338  }
339 
340  // Now resize the storage schemes
341  du_ddata.resize(n_u_dof, 0.0);
342  global_eqn_number.resize(n_u_dof, 0);
343 
344  // Loop over the nodes again and set the derivatives
345  unsigned count = 0;
346  for (unsigned l = 0; l < n_node; l++)
347  {
348  // Get the global equation number
349  int global_eqn = this->node_pt(l)->eqn_number(u_nodal_index);
350  // If it's positive
351  if (global_eqn >= 0)
352  {
353  // Set the global equation number
354  global_eqn_number[count] = global_eqn;
355  // Set the derivative with respect to the unknown
356  du_ddata[count] = psi[l];
357  // Increase the counter
358  ++count;
359  }
360  }
361  }
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
virtual unsigned u_index_axisym_adv_diff() const
Broken assignment operator.
Definition: steady_axisym_advection_diffusion_elements.h:94

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

◆ dshape_and_dtest_eulerian_adv_diff()

virtual double oomph::SteadyAxisymAdvectionDiffusionEquations::dshape_and_dtest_eulerian_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::QSteadyAxisymAdvectionDiffusionElement< NNODE_1D >.

◆ dshape_and_dtest_eulerian_at_knot_adv_diff()

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

Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of mapping

Implemented in oomph::QSteadyAxisymAdvectionDiffusionElement< NNODE_1D >.

Referenced by fill_in_generic_residual_contribution_adv_diff().

◆ fill_in_contribution_to_jacobian()

void oomph::SteadyAxisymAdvectionDiffusionEquations::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::FiniteElement.

271  {
272  // Call the generic routine with the flag set to 1
274  residuals, jacobian, GeneralisedElement::Dummy_matrix, 1);
275  }
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
virtual void fill_in_generic_residual_contribution_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: steady_axisym_advection_diffusion_elements.cc:46

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

◆ fill_in_contribution_to_residuals()

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

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

Reimplemented from oomph::GeneralisedElement.

256  {
257  // Call the generic residuals function with flag set to 0 and using
258  // a dummy matrix
260  residuals,
263  0);
264  }

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

◆ fill_in_generic_residual_contribution_adv_diff()

void oomph::SteadyAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_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

51  {
52  // Find out how many nodes there are
53  unsigned n_node = nnode();
54 
55  // Get the nodal index at which the unknown is stored
56  unsigned u_nodal_index = u_index_axisym_adv_diff();
57 
58  // Set up memory for the shape and test functions
59  Shape psi(n_node), test(n_node);
60  DShape dpsidx(n_node, 2), dtestdx(n_node, 2);
61 
62  // Set the value of n_intpt
63  unsigned n_intpt = integral_pt()->nweight();
64 
65  // Set the Vector to hold local coordinates
66  Vector<double> s(2);
67 
68  // Get Physical Variables from Element
69  double scaled_peclet = pe();
70 
71  // Integers used to store the local equation number and local unknown
72  // indices for the residuals and jacobians
73  int local_eqn = 0, local_unknown = 0;
74 
75  // Loop over the integration points
76  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
77  {
78  // Assign values of s
79  for (unsigned i = 0; i < 2; i++) s[i] = integral_pt()->knot(ipt, i);
80 
81  // Get the integral weight
82  double w = integral_pt()->weight(ipt);
83 
84  // Call the derivatives of the shape and test functions
86  ipt, psi, dpsidx, test, dtestdx);
87 
88  // Premultiply the weights and the Jacobian
89  double W = w * J;
90 
91  // Calculate local values of the solution and its derivatives
92  // Allocate
93  double interpolated_u = 0.0;
94  Vector<double> interpolated_x(2, 0.0);
95  Vector<double> interpolated_dudx(2, 0.0);
96 
97  // Calculate function value and derivatives:
98  //-----------------------------------------
99  // Loop over nodes
100  for (unsigned l = 0; l < n_node; l++)
101  {
102  // Get the value at the node
103  double u_value = raw_nodal_value(l, u_nodal_index);
104  interpolated_u += u_value * psi(l);
105  // Loop over the two coordinates directions
106  for (unsigned j = 0; j < 2; j++)
107  {
108  interpolated_x[j] += raw_nodal_position(l, j) * psi(l);
109  interpolated_dudx[j] += u_value * dpsidx(l, j);
110  }
111  }
112 
113  // Get source function
114  //-------------------
115  double source;
117 
118 
119  // Get wind
120  //--------
121  Vector<double> wind(3);
123 
124  // r is the first position component
125  double r = interpolated_x[0];
126 
127  // TEMPERATURE EQUATION (Neglected viscous dissipation)
128  // Assemble residuals and Jacobian
129  //--------------------------------
130 
131  // Loop over the test functions
132  for (unsigned l = 0; l < n_node; l++)
133  {
134  // Set the local equation number
135  local_eqn = nodal_local_eqn(l, u_nodal_index);
136 
137  /*IF it's not a boundary condition*/
138  if (local_eqn >= 0)
139  {
140  // Add body force/source term
141  residuals[local_eqn] += r * source * test(l) * W;
142 
143  // The Advection Diffusion bit itself
144  for (unsigned k = 0; k < 2; k++)
145  {
146  // Construct the contribution to the residuals
147  residuals[local_eqn] -=
148  r * interpolated_dudx[k] *
149  (scaled_peclet * wind[k] * test(l) + dtestdx(l, k)) * W;
150  }
151 
152  // Calculate the jacobian
153  //-----------------------
154  if (flag)
155  {
156  // Loop over the velocity shape functions again
157  for (unsigned l2 = 0; l2 < n_node; l2++)
158  {
159  // Set the number of the unknown
160  local_unknown = nodal_local_eqn(l2, u_nodal_index);
161 
162  // If at a non-zero degree of freedom add in the entry
163  if (local_unknown >= 0)
164  {
165  // Add contribution to Elemental Matrix
166  for (unsigned i = 0; i < 2; i++)
167  {
168  // Assemble Jacobian term
169  jacobian(local_eqn, local_unknown) -=
170  r * dpsidx(l2, i) *
171  (scaled_peclet * wind[i] * test(l) + dtestdx(l, i)) * W;
172  }
173  }
174  }
175  }
176  }
177  }
178 
179 
180  } // End of loop over integration points
181  }
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
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
virtual double dshape_and_dtest_eulerian_at_knot_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
virtual void get_wind_axisym_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Definition: steady_axisym_advection_diffusion_elements.h:200
const double & pe() const
Peclet number.
Definition: steady_axisym_advection_diffusion_elements.h:164
virtual void get_source_axisym_adv_diff(const unsigned &ipt, const Vector< double > &x, double &source) const
Definition: steady_axisym_advection_diffusion_elements.h:179
char char char int int * k
Definition: level2_impl.h:374
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 dshape_and_dtest_eulerian_at_knot_adv_diff(), get_source_axisym_adv_diff(), get_wind_axisym_adv_diff(), i, oomph::FiniteElement::integral_pt(), oomph::FiniteElement::interpolated_x(), J, j, k, oomph::Integral::knot(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_local_eqn(), oomph::Integral::nweight(), pe(), UniformPSDSelfTest::r, oomph::FiniteElement::raw_nodal_position(), oomph::FiniteElement::raw_nodal_value(), s, TestProblem::source(), Eigen::test, u_index_axisym_adv_diff(), w, oomph::QuadTreeNames::W, and oomph::Integral::weight().

Referenced by fill_in_contribution_to_jacobian(), and fill_in_contribution_to_residuals().

◆ get_flux()

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

Get flux: \( \mbox{flux}[i] = \nabla u = \mbox{d}u / \mbox{d}x_i \).

222  {
223  // Find out how many nodes there are in the element
224  unsigned n_node = nnode();
225 
226  // Get the nodal index at which the unknown is stored
227  unsigned u_nodal_index = u_index_axisym_adv_diff();
228 
229  // Set up memory for the shape and test functions
230  Shape psi(n_node);
231  DShape dpsidx(n_node, 2);
232 
233  // Call the derivatives of the shape and test functions
234  dshape_eulerian(s, psi, dpsidx);
235 
236  // Initialise to zero
237  for (unsigned j = 0; j < 2; j++)
238  {
239  flux[j] = 0.0;
240  }
241 
242  // Loop over nodes
243  for (unsigned l = 0; l < n_node; l++)
244  {
245  // Loop over derivative directions
246  for (unsigned j = 0; j < 2; j++)
247  {
248  flux[j] += nodal_value(l, u_nodal_index) * dpsidx(l, j);
249  }
250  }
251  }
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
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_axisym_adv_diff().

◆ get_source_axisym_adv_diff()

virtual void oomph::SteadyAxisymAdvectionDiffusionEquations::get_source_axisym_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

182  {
183  // If no source function has been set, return zero
184  if (Source_fct_pt == 0)
185  {
186  source = 0.0;
187  }
188  else
189  {
190  // Get source strength
191  (*Source_fct_pt)(x, source);
192  }
193  }

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

Referenced by fill_in_generic_residual_contribution_adv_diff().

◆ get_wind_axisym_adv_diff()

virtual void oomph::SteadyAxisymAdvectionDiffusionEquations::get_wind_axisym_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

204  {
205  // If no wind function has been set, return zero
206  if (Wind_fct_pt == 0)
207  {
208  for (unsigned i = 0; i < 2; i++)
209  {
210  wind[i] = 0.0;
211  }
212  }
213  else
214  {
215  // Get wind
216  (*Wind_fct_pt)(x, wind);
217  }
218  }

References i, Wind_fct_pt, and plotDoE::x.

Referenced by fill_in_generic_residual_contribution_adv_diff(), and output().

◆ interpolated_u_adv_diff()

double oomph::SteadyAxisymAdvectionDiffusionEquations::interpolated_u_adv_diff ( const Vector< double > &  s) const
inline

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

280  {
281  // Find number of nodes
282  unsigned n_node = nnode();
283 
284  // Get the nodal index at which the unknown is stored
285  unsigned u_nodal_index = u_index_axisym_adv_diff();
286 
287  // Local shape function
288  Shape psi(n_node);
289 
290  // Find values of shape function
291  shape(s, psi);
292 
293  // Initialise value of u
294  double interpolated_u = 0.0;
295 
296  // Loop over the local nodes and sum
297  for (unsigned l = 0; l < n_node; l++)
298  {
299  interpolated_u += nodal_value(l, u_nodal_index) * psi[l];
300  }
301 
302  return (interpolated_u);
303  }

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

Referenced by compute_error(), and output().

◆ output() [1/4]

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

C_style output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QSteadyAxisymAdvectionDiffusionElement< NNODE_1D >.

113  {
114  unsigned n_plot = 5;
115  output(file_pt, n_plot);
116  }
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: steady_axisym_advection_diffusion_elements.h:100

References output().

◆ output() [2/4]

void oomph::SteadyAxisymAdvectionDiffusionEquations::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::QSteadyAxisymAdvectionDiffusionElement< NNODE_1D >.

270  {
271  // Vector of local coordinates
272  Vector<double> s(2);
273 
274  // Tecplot header info
275  fprintf(file_pt, "%s", tecplot_zone_string(nplot).c_str());
276 
277  // Loop over plot points
278  unsigned num_plot_points = nplot_points(nplot);
279  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
280  {
281  // Get local coordinates of plot point
282  get_s_plot(iplot, nplot, s);
283 
284  for (unsigned i = 0; i < 2; i++)
285  {
286  fprintf(file_pt, "%g ", interpolated_x(s, i));
287  }
288  fprintf(file_pt, "%g \n", interpolated_u_adv_diff(s));
289  }
290 
291  // Write tecplot footer (e.g. FE connectivity lists)
292  write_tecplot_zone_footer(file_pt, nplot);
293  }
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_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::SteadyAxisymAdvectionDiffusionEquations::output ( std::ostream &  outfile)
inlinevirtual

Output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QSteadyAxisymAdvectionDiffusionElement< NNODE_1D >.

101  {
102  unsigned nplot = 5;
103  output(outfile, nplot);
104  }

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

◆ output() [4/4]

void oomph::SteadyAxisymAdvectionDiffusionEquations::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::QSteadyAxisymAdvectionDiffusionElement< NNODE_1D >.

219  {
220  // Vector of local coordinates
221  Vector<double> s(2);
222 
223  // Tecplot header info
224  outfile << tecplot_zone_string(nplot);
225 
226  // Loop over plot points
227  unsigned num_plot_points = nplot_points(nplot);
228  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
229  {
230  // Get local coordinates of plot point
231  get_s_plot(iplot, nplot, s);
232 
233  // Get Eulerian coordinate of plot point
234  Vector<double> x(2);
235  interpolated_x(s, x);
236 
237  for (unsigned i = 0; i < 2; i++)
238  {
239  outfile << x[i] << " ";
240  }
241  outfile << interpolated_u_adv_diff(s) << " ";
242 
243  // Get the wind
244  Vector<double> wind(2);
245  // Dummy ipt argument needed... ?
246  unsigned ipt = 0;
247  get_wind_axisym_adv_diff(ipt, s, x, wind);
248  for (unsigned i = 0; i < 2; i++)
249  {
250  outfile << wind[i] << " ";
251  }
252  outfile << std::endl;
253  }
254 
255  // Write tecplot footer (e.g. FE connectivity lists)
256  write_tecplot_zone_footer(outfile, nplot);
257  }

References oomph::FiniteElement::get_s_plot(), get_wind_axisym_adv_diff(), i, interpolated_u_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::SteadyAxisymAdvectionDiffusionEquations::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::QSteadyAxisymAdvectionDiffusionElement< NNODE_1D >.

308  {
309  // Vector of local coordinates
310  Vector<double> s(2);
311 
312  // Vector for coordintes
313  Vector<double> x(2);
314 
315  // Tecplot header info
316  outfile << tecplot_zone_string(nplot);
317 
318  // Exact solution Vector (here a scalar)
319  Vector<double> exact_soln(1);
320 
321  // Loop over plot points
322  unsigned num_plot_points = nplot_points(nplot);
323  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
324  {
325  // Get local coordinates of plot point
326  get_s_plot(iplot, nplot, s);
327 
328  // Get x position as Vector
329  interpolated_x(s, x);
330 
331  // Get exact solution at this point
332  (*exact_soln_pt)(x, exact_soln);
333 
334  // Output x,y,...,u_exact
335  for (unsigned i = 0; i < 2; i++)
336  {
337  outfile << x[i] << " ";
338  }
339  outfile << exact_soln[0] << std::endl;
340  }
341 
342  // Write tecplot footer (e.g. FE connectivity lists)
343  write_tecplot_zone_footer(outfile, nplot);
344  }

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::QSteadyAxisymAdvectionDiffusionElement< NNODE_1D >::output_fct().

◆ pe()

const double& oomph::SteadyAxisymAdvectionDiffusionEquations::pe ( ) const
inline

Peclet number.

165  {
166  return *Pe_pt;
167  }

References Pe_pt.

Referenced by fill_in_generic_residual_contribution_adv_diff().

◆ pe_pt()

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

Pointer to Peclet number.

171  {
172  return Pe_pt;
173  }

References Pe_pt.

◆ self_test()

unsigned oomph::SteadyAxisymAdvectionDiffusionEquations::self_test ( )
virtual

Self-test: Return 0 for OK.

Reimplemented from oomph::FiniteElement.

189  {
190  bool passed = true;
191 
192  // Check lower-level stuff
193  if (FiniteElement::self_test() != 0)
194  {
195  passed = false;
196  }
197 
198  // Return verdict
199  if (passed)
200  {
201  return 0;
202  }
203  else
204  {
205  return 1;
206  }
207  }
virtual unsigned self_test()
Definition: elements.cc:4440

References oomph::FiniteElement::self_test().

◆ source_fct_pt() [1/2]

SteadyAxisymAdvectionDiffusionSourceFctPt& oomph::SteadyAxisymAdvectionDiffusionEquations::source_fct_pt ( )
inline

Access function: Pointer to source function.

136  {
137  return Source_fct_pt;
138  }

References Source_fct_pt.

◆ source_fct_pt() [2/2]

SteadyAxisymAdvectionDiffusionSourceFctPt oomph::SteadyAxisymAdvectionDiffusionEquations::source_fct_pt ( ) const
inline

Access function: Pointer to source function. Const version.

143  {
144  return Source_fct_pt;
145  }

References Source_fct_pt.

◆ u_index_axisym_adv_diff()

virtual unsigned oomph::SteadyAxisymAdvectionDiffusionEquations::u_index_axisym_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 steady axisymmetric 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.

95  {
96  return 0;
97  }

Referenced by dinterpolated_u_adv_diff_ddata(), fill_in_generic_residual_contribution_adv_diff(), get_flux(), interpolated_u_adv_diff(), and oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::SteadyAxisymAdvectionDiffusionFluxElement().

◆ wind_fct_pt() [1/2]

SteadyAxisymAdvectionDiffusionWindFctPt& oomph::SteadyAxisymAdvectionDiffusionEquations::wind_fct_pt ( )
inline

Access function: Pointer to wind function.

150  {
151  return Wind_fct_pt;
152  }

References Wind_fct_pt.

◆ wind_fct_pt() [2/2]

SteadyAxisymAdvectionDiffusionWindFctPt oomph::SteadyAxisymAdvectionDiffusionEquations::wind_fct_pt ( ) const
inline

Access function: Pointer to wind function. Const version.

157  {
158  return Wind_fct_pt;
159  }

References Wind_fct_pt.

Member Data Documentation

◆ Default_peclet_number

double oomph::SteadyAxisymAdvectionDiffusionEquations::Default_peclet_number = 0.0
staticprivate

Static default value for the Peclet number.

2D Steady Axisymmetric Advection Diffusion elements

Default value for Peclet number

Referenced by SteadyAxisymAdvectionDiffusionEquations().

◆ Pe_pt

double* oomph::SteadyAxisymAdvectionDiffusionEquations::Pe_pt
protected

Pointer to global Peclet number.

Referenced by pe(), pe_pt(), and SteadyAxisymAdvectionDiffusionEquations().

◆ Source_fct_pt

SteadyAxisymAdvectionDiffusionSourceFctPt oomph::SteadyAxisymAdvectionDiffusionEquations::Source_fct_pt
protected

Pointer to source function:

Referenced by get_source_axisym_adv_diff(), and source_fct_pt().

◆ Wind_fct_pt

SteadyAxisymAdvectionDiffusionWindFctPt oomph::SteadyAxisymAdvectionDiffusionEquations::Wind_fct_pt
protected

Pointer to wind function:

Referenced by get_wind_axisym_adv_diff(), and wind_fct_pt().


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