oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT > Class Template Reference

#include <steady_axisym_advection_diffusion_elements.h>

+ Inheritance diagram for oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >:

Public Types

typedef void(* SteadyAxisymAdvectionDiffusionPrescribedBetaFctPt) (const Vector< double > &x, double &beta)
 
typedef void(* SteadyAxisymAdvectionDiffusionPrescribedAlphaFctPt) (const Vector< double > &x, double &alpha)
 
- 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

 SteadyAxisymAdvectionDiffusionFluxElement (FiniteElement *const &bulk_el_pt, const int &face_index)
 
 SteadyAxisymAdvectionDiffusionFluxElement ()
 Broken empty constructor. More...
 
 SteadyAxisymAdvectionDiffusionFluxElement (const SteadyAxisymAdvectionDiffusionFluxElement &dummy)=delete
 Broken copy constructor. More...
 
SteadyAxisymAdvectionDiffusionPrescribedBetaFctPtbeta_fct_pt ()
 Broken assignment operator. More...
 
SteadyAxisymAdvectionDiffusionPrescribedAlphaFctPtalpha_fct_pt ()
 Access function for the prescribed-alpha function pointer. More...
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Add the element's contribution to its residual vector. More...
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 
void output (std::ostream &outfile)
 
void output (std::ostream &outfile, const unsigned &nplot)
 
- Public Member Functions inherited from oomph::FaceElement
 FaceElement ()
 Constructor: Initialise all appropriate member data. More...
 
virtual ~FaceElement ()
 Empty virtual destructor. More...
 
 FaceElement (const FaceElement &)=delete
 Broken copy constructor. More...
 
const unsignedboundary_number_in_bulk_mesh () const
 Broken assignment operator. More...
 
void set_boundary_number_in_bulk_mesh (const unsigned &b)
 Set function for the boundary number in bulk mesh. More...
 
double J_eulerian (const Vector< double > &s) const
 
double J_eulerian_at_knot (const unsigned &ipt) const
 
void check_J_eulerian_at_knots (bool &passed) const
 
double interpolated_x (const Vector< double > &s, const unsigned &i) const
 
double interpolated_x (const unsigned &t, const Vector< double > &s, const unsigned &i) const
 
void interpolated_x (const Vector< double > &s, Vector< double > &x) const
 
void interpolated_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) const
 
double interpolated_dxdt (const Vector< double > &s, const unsigned &i, const unsigned &t)
 
void interpolated_dxdt (const Vector< double > &s, const unsigned &t, Vector< double > &dxdt)
 
intnormal_sign ()
 
int normal_sign () const
 
intface_index ()
 
int face_index () const
 
const Vector< double > * tangent_direction_pt () const
 Public access function for the tangent direction pointer. More...
 
void set_tangent_direction (Vector< double > *tangent_direction_pt)
 Set the tangent direction vector. More...
 
void turn_on_warning_for_discontinuous_tangent ()
 
void turn_off_warning_for_discontinuous_tangent ()
 
void continuous_tangent_and_outer_unit_normal (const Vector< double > &s, Vector< Vector< double >> &tang_vec, Vector< double > &unit_normal) const
 
void continuous_tangent_and_outer_unit_normal (const unsigned &ipt, Vector< Vector< double >> &tang_vec, Vector< double > &unit_normal) const
 
void outer_unit_normal (const Vector< double > &s, Vector< double > &unit_normal) const
 Compute outer unit normal at the specified local coordinate. More...
 
void outer_unit_normal (const unsigned &ipt, Vector< double > &unit_normal) const
 Compute outer unit normal at ipt-th integration point. More...
 
FiniteElement *& bulk_element_pt ()
 Pointer to higher-dimensional "bulk" element. More...
 
FiniteElementbulk_element_pt () const
 Pointer to higher-dimensional "bulk" element (const version) More...
 
CoordinateMappingFctPtface_to_bulk_coordinate_fct_pt ()
 
CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt () const
 
BulkCoordinateDerivativesFctPtbulk_coordinate_derivatives_fct_pt ()
 
BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt () const
 
Vector< doublelocal_coordinate_in_bulk (const Vector< double > &s) const
 
void get_local_coordinate_in_bulk (const Vector< double > &s, Vector< double > &s_bulk) const
 
void get_ds_bulk_ds_face (const Vector< double > &s, DenseMatrix< double > &dsbulk_dsface, unsigned &interior_direction) const
 
unsignedbulk_position_type (const unsigned &i)
 
const unsignedbulk_position_type (const unsigned &i) const
 
void bulk_node_number_resize (const unsigned &i)
 Resize the storage for the bulk node numbers. More...
 
unsignedbulk_node_number (const unsigned &n)
 
const unsignedbulk_node_number (const unsigned &n) const
 
void bulk_position_type_resize (const unsigned &i)
 Resize the storage for bulk_position_type to i entries. More...
 
unsignednbulk_value (const unsigned &n)
 
unsigned nbulk_value (const unsigned &n) const
 
void nbulk_value_resize (const unsigned &i)
 
void resize_nodes (Vector< unsigned > &nadditional_data_values)
 
void output_zeta (std::ostream &outfile, const unsigned &nplot)
 Output boundary coordinate zeta. 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
 
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...
 
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)
 
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 (FILE *file_pt)
 
virtual void output (FILE *file_pt, const unsigned &n_plot)
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 Output an exact solution over the element. More...
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
 Output a time-dependent exact solution over the element. More...
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, const SolutionFunctorBase &exact_soln) const
 Output a time-dependent exact solution over the element. More...
 
virtual void get_s_plot (const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
 
virtual std::string tecplot_zone_string (const unsigned &nplot) const
 
virtual void write_tecplot_zone_footer (std::ostream &outfile, const unsigned &nplot) const
 
virtual void write_tecplot_zone_footer (FILE *file_pt, const unsigned &nplot) const
 
virtual unsigned nplot_points (const unsigned &nplot) const
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 Calculate the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Calculate the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_abs_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error)
 
void integrate_fct (FiniteElement::SteadyExactSolutionFctPt integrand_fct_pt, Vector< double > &integral)
 Integrate Vector-valued function over element. More...
 
void integrate_fct (FiniteElement::UnsteadyExactSolutionFctPt integrand_fct_pt, const double &time, Vector< double > &integral)
 Integrate Vector-valued time-dep function over element. More...
 
virtual void build_face_element (const int &face_index, FaceElement *face_element_pt)
 
virtual unsigned self_test ()
 
virtual unsigned get_bulk_node_number (const int &face_index, const unsigned &i) const
 
virtual int face_outer_unit_normal_sign (const int &face_index) const
 Get the sign of the outer unit normal on the face given by face_index. More...
 
virtual unsigned nnode_on_face () const
 
void face_node_number_error_check (const unsigned &i) const
 Range check for face node numbers. More...
 
virtual CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt (const int &face_index) const
 
virtual BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt (const int &face_index) const
 
- Public Member Functions inherited from oomph::GeneralisedElement
 GeneralisedElement ()
 Constructor: Initialise all pointers and all values to zero. More...
 
virtual ~GeneralisedElement ()
 Virtual destructor to clean up any memory allocated by the object. More...
 
 GeneralisedElement (const GeneralisedElement &)=delete
 Broken copy constructor. More...
 
void operator= (const GeneralisedElement &)=delete
 Broken assignment operator. More...
 
Data *& internal_data_pt (const unsigned &i)
 Return a pointer to i-th internal data object. More...
 
Data *const & internal_data_pt (const unsigned &i) const
 Return a pointer to i-th internal data object (const version) More...
 
Data *& external_data_pt (const unsigned &i)
 Return a pointer to i-th external data object. More...
 
Data *const & external_data_pt (const unsigned &i) const
 Return a pointer to i-th external data object (const version) More...
 
unsigned long eqn_number (const unsigned &ieqn_local) const
 
int local_eqn_number (const unsigned long &ieqn_global) const
 
unsigned add_external_data (Data *const &data_pt, const bool &fd=true)
 
bool external_data_fd (const unsigned &i) const
 
void exclude_external_data_fd (const unsigned &i)
 
void include_external_data_fd (const unsigned &i)
 
void flush_external_data ()
 Flush all external data. More...
 
void flush_external_data (Data *const &data_pt)
 Flush the object addressed by data_pt from the external data array. More...
 
unsigned ninternal_data () const
 Return the number of internal data objects. More...
 
unsigned nexternal_data () const
 Return the number of external data objects. More...
 
unsigned ndof () const
 Return the number of equations/dofs in the element. More...
 
void dof_vector (const unsigned &t, Vector< double > &dof)
 Return the vector of dof values at time level t. More...
 
void dof_pt_vector (Vector< double * > &dof_pt)
 Return the vector of pointers to dof values. More...
 
void set_internal_data_time_stepper (const unsigned &i, TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void assign_internal_eqn_numbers (unsigned long &global_number, Vector< double * > &Dof_pt)
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
void add_internal_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt)
 
virtual void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void complete_setup_of_dependencies ()
 
virtual void get_residuals (Vector< double > &residuals)
 
virtual void get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
virtual void get_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void get_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void get_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void get_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void get_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void get_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void get_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void get_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_norm (double &norm)
 
virtual unsigned ndof_types () const
 
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 
- Public Member Functions inherited from oomph::GeomObject
 GeomObject ()
 Default constructor. More...
 
 GeomObject (const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim, TimeStepper *time_stepper_pt)
 
 GeomObject (const GeomObject &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const GeomObject &)=delete
 Broken assignment operator. More...
 
virtual ~GeomObject ()
 (Empty) destructor More...
 
unsigned nlagrangian () const
 Access function to # of Lagrangian coordinates. More...
 
unsigned ndim () const
 Access function to # of Eulerian coordinates. More...
 
void set_nlagrangian_and_ndim (const unsigned &n_lagrangian, const unsigned &n_dim)
 Set # of Lagrangian and Eulerian coordinates. More...
 
TimeStepper *& time_stepper_pt ()
 
TimeSteppertime_stepper_pt () const
 
virtual void position (const double &t, const Vector< double > &zeta, Vector< double > &r) const
 
virtual void dposition (const Vector< double > &zeta, DenseMatrix< double > &drdzeta) const
 
virtual void d2position (const Vector< double > &zeta, RankThreeTensor< double > &ddrdzeta) const
 
virtual void d2position (const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const
 

Protected Member Functions

double shape_and_test (const Vector< double > &s, Shape &psi, Shape &test) const
 
double shape_and_test_at_knot (const unsigned &ipt, Shape &psi, Shape &test) const
 
void get_beta (const Vector< double > &x, double &beta)
 
void get_alpha (const Vector< double > &x, double &alpha)
 
- Protected Member Functions inherited from oomph::FaceElement
void add_additional_values (const Vector< unsigned > &nadditional_values, const unsigned &id)
 
- 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)
 

Private Member Functions

void fill_in_generic_residual_contribution_adv_diff_flux (Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
 

Private Attributes

SteadyAxisymAdvectionDiffusionPrescribedBetaFctPt Beta_fct_pt
 Function pointer to the (global) prescribed-beta function. More...
 
SteadyAxisymAdvectionDiffusionPrescribedAlphaFctPt Alpha_fct_pt
 Function pointer to the (global) prescribed-alpha function. More...
 
unsigned U_index_adv_diff
 The index at which the unknown is stored at the nodes. 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
 
- Protected Attributes inherited from oomph::FaceElement
unsigned Boundary_number_in_bulk_mesh
 The boundary number in the bulk mesh to which this element is attached. More...
 
FiniteElementBulk_element_pt
 Pointer to the associated higher-dimensional "bulk" element. More...
 
Vector< unsignedBulk_node_number
 
Vector< unsignedNbulk_value
 
Vector< double > * Tangent_direction_pt
 
- Protected Attributes inherited from oomph::FiniteElement
MacroElementMacro_elem_pt
 Pointer to the element's macro element (NULL by default) More...
 
- Protected Attributes inherited from oomph::GeomObject
unsigned NLagrangian
 Number of Lagrangian (intrinsic) coordinates. More...
 
unsigned Ndim
 Number of Eulerian coordinates. More...
 
TimeStepperGeom_object_time_stepper_pt
 
- Static Protected Attributes inherited from oomph::FiniteElement
static const unsigned Default_Initial_Nvalue = 0
 Default value for the number of values at a node. More...
 
static const double Node_location_tolerance = 1.0e-14
 
static const unsigned N2deriv [] = {0, 1, 3, 6}
 
- Static Protected Attributes inherited from oomph::GeneralisedElement
static DenseMatrix< doubleDummy_matrix
 
static std::deque< double * > Dof_pt_deque
 

Detailed Description

template<class ELEMENT>
class oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >

///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// A class for elements that allow the imposition of an applied Robin boundary condition on the boundaries of Steady Axisymmnetric Advection Diffusion Flux elements.

\[ -\Delta u \cdot \mathbf{n} + \alpha(r,z) u = \beta(r,z) \]

The element geometry is obtained from the FaceGeometry<ELEMENT> policy class.

Member Typedef Documentation

◆ SteadyAxisymAdvectionDiffusionPrescribedAlphaFctPt

template<class ELEMENT >
typedef void(* oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::SteadyAxisymAdvectionDiffusionPrescribedAlphaFctPt) (const Vector< double > &x, double &alpha)

Function pointer to the prescribed-alpha function fct(x,alpha(x)) – x is a Vector!

◆ SteadyAxisymAdvectionDiffusionPrescribedBetaFctPt

template<class ELEMENT >
typedef void(* oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::SteadyAxisymAdvectionDiffusionPrescribedBetaFctPt) (const Vector< double > &x, double &beta)

Function pointer to the prescribed-beta function fct(x,beta(x)) – x is a Vector!

Constructor & Destructor Documentation

◆ SteadyAxisymAdvectionDiffusionFluxElement() [1/3]

template<class ELEMENT >
oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::SteadyAxisymAdvectionDiffusionFluxElement ( FiniteElement *const &  bulk_el_pt,
const int face_index 
)

Constructor, takes the pointer to the "bulk" element and the index of the face to be created

//////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// Constructor, takes the pointer to the "bulk" element and the index of the face to be created

820  : FaceGeometry<ELEMENT>(), FaceElement()
821 
822  {
823  // Let the bulk element build the FaceElement, i.e. setup the pointers
824  // to its nodes (by referring to the appropriate nodes in the bulk
825  // element), etc.
826  bulk_el_pt->build_face_element(face_index, this);
827 
828 
829 #ifdef PARANOID
830  {
831  // Check that the element is not a refineable 3d element
832  ELEMENT* elem_pt = dynamic_cast<ELEMENT*>(bulk_el_pt);
833  // If it's three-d
834  if (elem_pt->dim() == 3)
835  {
836  // Is it refineable
837  RefineableElement* ref_el_pt =
838  dynamic_cast<RefineableElement*>(elem_pt);
839  if (ref_el_pt != 0)
840  {
841  if (this->has_hanging_nodes())
842  {
843  throw OomphLibError("This flux element will not work correctly if "
844  "nodes are hanging\n",
847  }
848  }
849  }
850  }
851 #endif
852 
853  // Initialise the prescribed-beta function pointer to zero
854  Beta_fct_pt = 0;
855 
856  // Set up U_index_adv_diff. Initialise to zero, which probably won't change
857  // in most cases, oh well, the price we pay for generality
858  U_index_adv_diff = 0;
859 
860  // Cast to the appropriate AdvectionDiffusionEquation so that we can
861  // find the index at which the variable is stored
862  // We assume that the dimension of the full problem is the same
863  // as the dimension of the node, if this is not the case you will have
864  // to write custom elements, sorry
865  SteadyAxisymAdvectionDiffusionEquations* eqn_pt =
866  dynamic_cast<SteadyAxisymAdvectionDiffusionEquations*>(bulk_el_pt);
867 
868  // If the cast has failed die
869  if (eqn_pt == 0)
870  {
871  std::string error_string = "Bulk element must inherit from "
872  "SteadyAxisymAdvectionDiffusionEquations.";
873  error_string +=
874  "Nodes are two dimensional, but cannot cast the bulk element to\n";
875  error_string += "SteadyAxisymAdvectionDiffusionEquations<2>\n.";
876  error_string +=
877  "If you desire this functionality, you must implement it yourself\n";
878 
879  throw OomphLibError(
881  }
882  else
883  {
884  // Read the index from the (cast) bulk element.
885  U_index_adv_diff = eqn_pt->u_index_axisym_adv_diff();
886  }
887  }
int & face_index()
Definition: elements.h:4626
FaceElement()
Constructor: Initialise all appropriate member data.
Definition: elements.h:4446
bool has_hanging_nodes() const
Definition: elements.h:2470
unsigned U_index_adv_diff
The index at which the unknown is stored at the nodes.
Definition: steady_axisym_advection_diffusion_elements.h:801
SteadyAxisymAdvectionDiffusionPrescribedBetaFctPt Beta_fct_pt
Function pointer to the (global) prescribed-beta function.
Definition: steady_axisym_advection_diffusion_elements.h:795
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::Beta_fct_pt, oomph::FiniteElement::build_face_element(), oomph::FaceElement::face_index(), oomph::FiniteElement::has_hanging_nodes(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::Global_string_for_annotation::string(), oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::U_index_adv_diff, and oomph::SteadyAxisymAdvectionDiffusionEquations::u_index_axisym_adv_diff().

◆ SteadyAxisymAdvectionDiffusionFluxElement() [2/3]

Broken empty constructor.

630  {
631  throw OomphLibError("Don't call empty constructor for "
632  "SteadyAxisymAdvectionDiffusionFluxElement",
635  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ SteadyAxisymAdvectionDiffusionFluxElement() [3/3]

Broken copy constructor.

Member Function Documentation

◆ alpha_fct_pt()

Access function for the prescribed-alpha function pointer.

653  {
654  return Alpha_fct_pt;
655  }
SteadyAxisymAdvectionDiffusionPrescribedAlphaFctPt Alpha_fct_pt
Function pointer to the (global) prescribed-alpha function.
Definition: steady_axisym_advection_diffusion_elements.h:798

References oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::Alpha_fct_pt.

Referenced by AxisymFreeSurfaceNozzleAdvDiffRobinProblem< ELEMENT >::AxisymFreeSurfaceNozzleAdvDiffRobinProblem(), and Wind_fct_pt().

◆ beta_fct_pt()

Broken assignment operator.

Access function for the prescribed-beta function pointer

647  {
648  return Beta_fct_pt;
649  }

References oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::Beta_fct_pt.

Referenced by AxisymFreeSurfaceNozzleAdvDiffRobinProblem< ELEMENT >::AxisymFreeSurfaceNozzleAdvDiffRobinProblem(), and Wind_fct_pt().

◆ fill_in_contribution_to_jacobian()

template<class ELEMENT >
void oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
inlinevirtual

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

Reimplemented from oomph::FiniteElement.

672  {
673  // Call the generic routine with the flag set to 1
675  residuals, jacobian, 1);
676  }
void fill_in_generic_residual_contribution_adv_diff_flux(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: steady_axisym_advection_diffusion_elements.h:897

References oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::fill_in_generic_residual_contribution_adv_diff_flux().

◆ fill_in_contribution_to_residuals()

template<class ELEMENT >
void oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inlinevirtual

Add the element's contribution to its residual vector.

Reimplemented from oomph::GeneralisedElement.

660  {
661  // Call the generic residuals function with flag set to 0
662  // using a dummy matrix
664  residuals, GeneralisedElement::Dummy_matrix, 0);
665  }
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227

References oomph::GeneralisedElement::Dummy_matrix, and oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::fill_in_generic_residual_contribution_adv_diff_flux().

◆ fill_in_generic_residual_contribution_adv_diff_flux()

template<class ELEMENT >
void oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::fill_in_generic_residual_contribution_adv_diff_flux ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
unsigned  flag 
)
private

Add the element's contribution to its residual vector. flag=1(or 0): do (or don't) compute the Jacobian as well.

Compute the element's residual vector and the (zero) Jacobian matrix for the Robin boundary condition:

\[ \Delta u \cdot \mathbf{n} + \alpha (\mathbf{x}) = \beta (\mathbf{x}) \]

899  {
900  // Find out how many nodes there are
901  const unsigned n_node = nnode();
902 
903  // Locally cache the index at which the variable is stored
904  const unsigned u_index_axisym_adv_diff = U_index_adv_diff;
905 
906  // Set up memory for the shape and test functions
907  Shape psif(n_node), testf(n_node);
908 
909  // Set the value of n_intpt
910  const unsigned n_intpt = integral_pt()->nweight();
911 
912  // Set the Vector to hold local coordinates
913  Vector<double> s(1);
914 
915  // Integers used to store the local equation number and local unknown
916  // indices for the residuals and jacobians
917  int local_eqn = 0, local_unknown = 0;
918 
919  // Loop over the integration points
920  //--------------------------------
921  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
922  {
923  // Assign values of s
924  for (unsigned i = 0; i < 1; i++)
925  {
926  s[i] = integral_pt()->knot(ipt, i);
927  }
928 
929  // Get the integral weight
930  double w = integral_pt()->weight(ipt);
931 
932  // Find the shape and test functions and return the Jacobian
933  // of the mapping
934  double J = shape_and_test(s, psif, testf);
935 
936  // Premultiply the weights and the Jacobian
937  double W = w * J;
938 
939  // Calculate local values of the solution and its derivatives
940  // Allocate
941  double interpolated_u = 0.0;
942  Vector<double> interpolated_x(2, 0.0);
943 
944  // Calculate position
945  for (unsigned l = 0; l < n_node; l++)
946  {
947  // Get the value at the node
948  double u_value = raw_nodal_value(l, u_index_axisym_adv_diff);
949  interpolated_u += u_value * psif(l);
950  // Loop over coordinate direction
951  for (unsigned i = 0; i < 2; i++)
952  {
953  interpolated_x[i] += nodal_position(l, i) * psif(l);
954  }
955  }
956 
957  // Get the imposed beta (beta=flux when alpha=0.0)
958  double beta;
960 
961  // Get the imposed alpha
962  double alpha;
964 
965  // r is the first position component
966  double r = interpolated_x[0];
967 
968  // Now add to the appropriate equations
969 
970  // Loop over the test functions
971  for (unsigned l = 0; l < n_node; l++)
972  {
973  // Set the local equation number
974  local_eqn = nodal_local_eqn(l, u_index_axisym_adv_diff);
975  /*IF it's not a boundary condition*/
976  if (local_eqn >= 0)
977  {
978  // Add the prescribed beta terms
979  residuals[local_eqn] -=
980  r * (beta - alpha * interpolated_u) * testf(l) * W;
981 
982  // Calculate the Jacobian
983  //----------------------
984  if ((flag) && (alpha != 0.0))
985  {
986  // Loop over the velocity shape functions again
987  for (unsigned l2 = 0; l2 < n_node; l2++)
988  {
989  // Set the number of the unknown
990  local_unknown = nodal_local_eqn(l2, u_index_axisym_adv_diff);
991 
992  // If at a non-zero degree of freedom add in the entry
993  if (local_unknown >= 0)
994  {
995  jacobian(local_eqn, local_unknown) +=
996  r * alpha * psif[l2] * testf[l] * W;
997  }
998  }
999  }
1000  }
1001  } // end loop over test functions
1002 
1003  } // end loop over integration points
1004 
1005  } // end fill_in_generic_residual_contribution_adv_diff_flux
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_x(const Vector< double > &s, const unsigned &i) const
Definition: elements.h:4528
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
double nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.h:2317
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2576
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.
void get_alpha(const Vector< double > &x, double &alpha)
Definition: steady_axisym_advection_diffusion_elements.h:773
void get_beta(const Vector< double > &x, double &beta)
Definition: steady_axisym_advection_diffusion_elements.h:757
double shape_and_test(const Vector< double > &s, Shape &psi, Shape &test) const
Definition: steady_axisym_advection_diffusion_elements.h:711
RealScalar s
Definition: level1_cplx_impl.h:130
RealScalar alpha
Definition: level1_cplx_impl.h:151
Scalar beta
Definition: level2_cplx_impl.h:36
r
Definition: UniformPSDSelfTest.py:20
@ W
Definition: quadtree.h:63

References alpha, beta, i, J, UniformPSDSelfTest::r, s, w, and oomph::QuadTreeNames::W.

Referenced by oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::fill_in_contribution_to_jacobian(), and oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::fill_in_contribution_to_residuals().

◆ get_alpha()

template<class ELEMENT >
void oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::get_alpha ( const Vector< double > &  x,
double alpha 
)
inlineprotected

Function to calculate the prescribed alpha at a given spatial position

774  {
775  // If the function pointer is zero return zero
776  if (Alpha_fct_pt == 0)
777  {
778  alpha = 0.0;
779  }
780  // Otherwise call the function
781  else
782  {
783  (*Alpha_fct_pt)(x, alpha);
784  }
785  }
list x
Definition: plotDoE.py:28

References alpha, oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::Alpha_fct_pt, and plotDoE::x.

◆ get_beta()

template<class ELEMENT >
void oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::get_beta ( const Vector< double > &  x,
double beta 
)
inlineprotected

Function to calculate the prescribed beta at a given spatial position

758  {
759  // If the function pointer is zero return zero
760  if (Beta_fct_pt == 0)
761  {
762  beta = 0.0;
763  }
764  // Otherwise call the function
765  else
766  {
767  (*Beta_fct_pt)(x, beta);
768  }
769  }

References beta, oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::Beta_fct_pt, and plotDoE::x.

◆ output() [1/2]

template<class ELEMENT >
void oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::output ( std::ostream &  outfile)
inlinevirtual

Output function – forward to broken version in FiniteElement until somebody decides what exactly they want to plot here...

Reimplemented from oomph::FiniteElement.

695  {
696  FiniteElement::output(outfile);
697  }
virtual void output(std::ostream &outfile)
Definition: elements.h:3050

References oomph::FiniteElement::output().

◆ output() [2/2]

template<class ELEMENT >
void oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::output ( std::ostream &  outfile,
const unsigned nplot 
)
inlinevirtual

Output function – forward to broken version in FiniteElement until somebody decides what exactly they want to plot here...

Reimplemented from oomph::FiniteElement.

702  {
703  FiniteElement::output(outfile, nplot);
704  }

References oomph::FiniteElement::output().

◆ shape_and_test()

template<class ELEMENT >
double oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::shape_and_test ( const Vector< double > &  s,
Shape psi,
Shape test 
) const
inlineprotected

Function to compute the shape and test functions and to return the Jacobian of mapping between local and global (Eulerian) coordinates

714  {
715  // Find number of nodes
716  unsigned n_node = nnode();
717 
718  // Get the shape functions
719  shape(s, psi);
720 
721  // Set the test functions to be the same as the shape functions
722  for (unsigned i = 0; i < n_node; i++)
723  {
724  test[i] = psi[i];
725  }
726 
727  // Return the value of the jacobian
728  return J_eulerian(s);
729  }
double J_eulerian(const Vector< double > &s) const
Definition: elements.cc:5242
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Definition: indexed_view.cpp:20

References i, oomph::FaceElement::J_eulerian(), oomph::FiniteElement::nnode(), s, and oomph::FiniteElement::shape().

◆ shape_and_test_at_knot()

template<class ELEMENT >
double oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::shape_and_test_at_knot ( const unsigned ipt,
Shape psi,
Shape test 
) const
inlineprotected

Function to compute the shape and test functions and to return the Jacobian of mapping between local and global (Eulerian) coordinates

738  {
739  // Find number of nodes
740  unsigned n_node = nnode();
741 
742  // Get the shape functions
743  shape_at_knot(ipt, psi);
744 
745  // Set the test functions to be the same as the shape functions
746  for (unsigned i = 0; i < n_node; i++)
747  {
748  test[i] = psi[i];
749  }
750 
751  // Return the value of the jacobian
752  return J_eulerian_at_knot(ipt);
753  }
double J_eulerian_at_knot(const unsigned &ipt) const
Definition: elements.cc:5328
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Definition: elements.cc:3220

References i, oomph::FaceElement::J_eulerian_at_knot(), oomph::FiniteElement::nnode(), and oomph::FiniteElement::shape_at_knot().

◆ zeta_nodal()

template<class ELEMENT >
double oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::zeta_nodal ( const unsigned n,
const unsigned k,
const unsigned i 
) const
inlinevirtual

Specify the value of nodal zeta from the face geometry The "global" intrinsic coordinate of the element when viewed as part of a geometric object should be given by the FaceElement representation, by default (needed to break indeterminacy if bulk element is SolidElement)

Reimplemented from oomph::FaceElement.

687  {
688  return FaceElement::zeta_nodal(n, k, i);
689  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: elements.h:4497
char char char int int * k
Definition: level2_impl.h:374

References i, k, n, and oomph::FaceElement::zeta_nodal().

Member Data Documentation

◆ Alpha_fct_pt

◆ Beta_fct_pt

◆ U_index_adv_diff

template<class ELEMENT >
unsigned oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::U_index_adv_diff
private

The index at which the unknown is stored at the nodes.

Referenced by oomph::SteadyAxisymAdvectionDiffusionFluxElement< ELEMENT >::SteadyAxisymAdvectionDiffusionFluxElement().


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