oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations Class Reference

#include <refineable_gen_axisym_advection_diffusion_elements.h>

+ Inheritance diagram for oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations:

Public Member Functions

 RefineableGeneralisedAxisymAdvectionDiffusionEquations ()
 Empty Constructor. More...
 
 RefineableGeneralisedAxisymAdvectionDiffusionEquations (const RefineableGeneralisedAxisymAdvectionDiffusionEquations &dummy)=delete
 Broken copy constructor. More...
 
unsigned num_Z2_flux_terms ()
 Broken assignment operator. More...
 
void get_Z2_flux (const Vector< double > &s, Vector< double > &flux)
 
void get_interpolated_values (const Vector< double > &s, Vector< double > &values)
 
void get_interpolated_values (const unsigned &t, const Vector< double > &s, Vector< double > &values)
 
double geometric_jacobian (const Vector< double > &x)
 Fill in the geometric Jacobian, which in this case is r. More...
 
void further_build ()
 Further build: Copy source function pointer from father element. More...
 
- Public Member Functions inherited from oomph::RefineableElement
 RefineableElement ()
 
virtual ~RefineableElement ()
 Destructor, delete the allocated storage for the hanging equations. More...
 
 RefineableElement (const RefineableElement &)=delete
 Broken copy constructor. More...
 
void operator= (const RefineableElement &)=delete
 Broken assignment operator. More...
 
Treetree_pt ()
 Access function: Pointer to quadtree representation of this element. More...
 
void set_tree_pt (Tree *my_tree_pt)
 Set pointer to quadtree representation of this element. More...
 
virtual unsigned required_nsons () const
 
bool refinement_is_enabled ()
 Flag to indicate suppression of any refinement. More...
 
void disable_refinement ()
 Suppress of any refinement for this element. More...
 
void enable_refinement ()
 Emnable refinement for this element. More...
 
template<class ELEMENT >
void split (Vector< ELEMENT * > &son_pt) const
 
int local_hang_eqn (Node *const &node_pt, const unsigned &i)
 
virtual void build (Mesh *&mesh_pt, Vector< Node * > &new_node_pt, bool &was_already_built, std::ofstream &new_nodes_file)=0
 
void set_refinement_level (const int &refine_level)
 Set the refinement level. More...
 
unsigned refinement_level () const
 Return the Refinement level. More...
 
void select_for_refinement ()
 Select the element for refinement. More...
 
void deselect_for_refinement ()
 Deselect the element for refinement. More...
 
void select_sons_for_unrefinement ()
 Unrefinement will be performed by merging the four sons of this element. More...
 
void deselect_sons_for_unrefinement ()
 
bool to_be_refined ()
 Has the element been selected for refinement? More...
 
bool sons_to_be_unrefined ()
 Has the element been selected for unrefinement? More...
 
virtual void rebuild_from_sons (Mesh *&mesh_pt)=0
 
virtual void unbuild ()
 
virtual void deactivate_element ()
 
virtual bool nodes_built ()
 Return true if all the nodes have been built, false if not. More...
 
long number () const
 Element number (for debugging/plotting) More...
 
void set_number (const long &mynumber)
 Set element number (for debugging/plotting) More...
 
virtual unsigned ncont_interpolated_values () const =0
 
virtual Nodeinterpolating_node_pt (const unsigned &n, const int &value_id)
 
virtual double local_one_d_fraction_of_interpolating_node (const unsigned &n1d, const unsigned &i, const int &value_id)
 
virtual Nodeget_interpolating_node_at_local_coordinate (const Vector< double > &s, const int &value_id)
 
virtual unsigned ninterpolating_node (const int &value_id)
 
virtual unsigned ninterpolating_node_1d (const int &value_id)
 
virtual void interpolating_basis (const Vector< double > &s, Shape &psi, const int &value_id) const
 
virtual void check_integrity (double &max_error)=0
 
void identify_field_data_for_interactions (std::set< std::pair< Data *, unsigned >> &paired_field_data)
 
void assign_nodal_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual RefineableElementroot_element_pt ()
 
virtual RefineableElementfather_element_pt () const
 Return a pointer to the father element. More...
 
void get_father_at_refinement_level (unsigned &refinement_level, RefineableElement *&father_at_reflevel_pt)
 
virtual void initial_setup (Tree *const &adopted_father_pt=0, const unsigned &initial_p_order=0)
 
virtual void pre_build (Mesh *&mesh_pt, Vector< Node * > &new_node_pt)
 Pre-build the element. More...
 
virtual void setup_hanging_nodes (Vector< std::ofstream * > &output_stream)
 
virtual void further_setup_hanging_nodes ()
 
void get_dresidual_dnodal_coordinates (RankThreeTensor< double > &dresidual_dnodal_coordinates)
 
unsigned nshape_controlling_nodes ()
 
std::map< Node *, unsignedshape_controlling_node_lookup ()
 
- 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 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 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 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_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 (std::ostream &outfile)
 
virtual void output (std::ostream &outfile, const unsigned &n_plot)
 
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
 
- Public Member Functions inherited from oomph::ElementWithZ2ErrorEstimator
 ElementWithZ2ErrorEstimator ()
 Default empty constructor. More...
 
 ElementWithZ2ErrorEstimator (const ElementWithZ2ErrorEstimator &)=delete
 Broken copy constructor. More...
 
void operator= (const ElementWithZ2ErrorEstimator &)=delete
 Broken assignment operator. More...
 
virtual unsigned ncompound_fluxes ()
 
virtual void compute_exact_Z2_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_flux_pt, double &error, double &norm)
 
virtual void get_Z2_compound_flux_indices (Vector< unsigned > &flux_index)
 
virtual unsigned nvertex_node () const =0
 Number of vertex nodes in the element. More...
 
virtual Nodevertex_node_pt (const unsigned &j) const =0
 
virtual unsigned nrecovery_order ()=0
 Order of recovery shape functions. More...
 

Protected Member Functions

void fill_in_generic_residual_contribution_cons_axisym_adv_diff (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
 
- Protected Member Functions inherited from oomph::RefineableElement
void assemble_local_to_eulerian_jacobian (const DShape &dpsids, DenseMatrix< double > &jacobian) const
 
void assemble_local_to_eulerian_jacobian2 (const DShape &d2psids, DenseMatrix< double > &jacobian2) const
 
void assemble_eulerian_base_vectors (const DShape &dpsids, DenseMatrix< double > &interpolated_G) const
 
double local_to_eulerian_mapping_diagonal (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
void assign_hanging_local_eqn_numbers (const bool &store_local_dof_pt)
 Assign the local equation numbers for hanging node variables. More...
 
virtual void fill_in_jacobian_from_nodal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
- Protected Member Functions inherited from oomph::FiniteElement
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 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
 
void fill_in_jacobian_from_nodal_by_fd (DenseMatrix< double > &jacobian)
 
virtual void update_before_nodal_fd ()
 
virtual void reset_after_nodal_fd ()
 
virtual void update_in_nodal_fd (const unsigned &i)
 
virtual void reset_in_nodal_fd (const unsigned &i)
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Zero-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 One-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Two-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
- Protected Member Functions inherited from oomph::GeneralisedElement
unsigned add_internal_data (Data *const &data_pt, const bool &fd=true)
 
bool internal_data_fd (const unsigned &i) const
 
void exclude_internal_data_fd (const unsigned &i)
 
void include_internal_data_fd (const unsigned &i)
 
void clear_global_eqn_numbers ()
 
void add_global_eqn_numbers (std::deque< unsigned long > const &global_eqn_numbers, std::deque< double * > const &global_dof_pt)
 
virtual void assign_internal_and_external_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void assign_additional_local_eqn_numbers ()
 
int internal_local_eqn (const unsigned &i, const unsigned &j) const
 
int external_local_eqn (const unsigned &i, const unsigned &j)
 
virtual void fill_in_contribution_to_residuals (Vector< double > &residuals)
 
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)
 

Additional Inherited Members

- Public Types inherited from oomph::FiniteElement
typedef void(* SteadyExactSolutionFctPt) (const Vector< double > &, Vector< double > &)
 
typedef void(* UnsteadyExactSolutionFctPt) (const double &, const Vector< double > &, Vector< double > &)
 
- Static Public Member Functions inherited from oomph::RefineableElement
static doublemax_integrity_tolerance ()
 Max. allowed discrepancy in element integrity check. More...
 
- Static Public Attributes inherited from oomph::FiniteElement
static double Tolerance_for_singular_jacobian = 1.0e-16
 Tolerance below which the jacobian is considered singular. More...
 
static bool Accept_negative_jacobian = false
 
static bool Suppress_output_while_checking_for_inverted_elements
 
- Static Public Attributes inherited from oomph::GeneralisedElement
static bool Suppress_warning_about_repeated_internal_data
 
static bool Suppress_warning_about_repeated_external_data = true
 
static double Default_fd_jacobian_step = 1.0e-8
 
- Static Protected Member Functions inherited from oomph::RefineableElement
static void check_value_id (const int &n_continuously_interpolated_values, const int &value_id)
 
- Protected Attributes inherited from oomph::RefineableElement
TreeTree_pt
 A pointer to a general tree object. More...
 
unsigned Refine_level
 Refinement level. More...
 
bool To_be_refined
 Flag for refinement. More...
 
bool Refinement_is_enabled
 Flag to indicate suppression of any refinement. More...
 
bool Sons_to_be_unrefined
 Flag for unrefinement. More...
 
long Number
 Global element number – for plotting/validation purposes. 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 Protected Attributes inherited from oomph::RefineableElement
static double Max_integrity_tolerance = 1.0e-8
 Max. allowed discrepancy in element integrity check. More...
 
- 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 version of the GeneralisedAxisymAdvectionDiffusion equations that can be used with non-uniform mesh refinement. In essence, the class overloads the fill_in_generic_residual_contribution_cons_axisym_adv_diff() function so that contributions from hanging nodes (or alternatively in-compatible function values) are taken into account.

Constructor & Destructor Documentation

◆ RefineableGeneralisedAxisymAdvectionDiffusionEquations() [1/2]

oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations::RefineableGeneralisedAxisymAdvectionDiffusionEquations ( )
inline

Empty Constructor.

65  {
66  }
ElementWithZ2ErrorEstimator()
Default empty constructor.
Definition: error_estimator.h:82
RefineableElement()
Definition: refineable_elements.h:188

◆ RefineableGeneralisedAxisymAdvectionDiffusionEquations() [2/2]

oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations::RefineableGeneralisedAxisymAdvectionDiffusionEquations ( const RefineableGeneralisedAxisymAdvectionDiffusionEquations dummy)
delete

Broken copy constructor.

Member Function Documentation

◆ fill_in_generic_residual_contribution_cons_axisym_adv_diff()

void oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_cons_axisym_adv_diff ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  mass_matrix,
unsigned  flag 
)
protected

Add the element's contribution to the elemental residual vector and/or Jacobian matrix flag=1: compute both flag=0: compute only residual vector

Add the element's contribution to the elemental residual vector and/or elemental jacobian matrix. This function overloads the standard version so that the possible presence of hanging nodes is taken into account.

42  {
43  // Find out how many nodes there are in the element
44  unsigned n_node = nnode();
45 
46  // Get the nodal index at which the unknown is stored
47  unsigned u_nodal_index = this->u_index_cons_axisym_adv_diff();
48 
49  // Set up memory for the shape and test functions
50  Shape psi(n_node), test(n_node);
51  DShape dpsidx(n_node, 2), dtestdx(n_node, 2);
52 
53  // Set the value of n_intpt
54  unsigned n_intpt = integral_pt()->nweight();
55 
56  // Set the Vector to hold local coordinates
57  Vector<double> s(2);
58 
59  // Get Peclet number
60  double peclet = this->pe();
61 
62  // Get the Peclet multiplied by the Strouhal number
63  double peclet_st = this->pe_st();
64 
65  // Integers used to store the local equation number and local unknown
66  // indices for the residuals and jacobians
67  int local_eqn = 0, local_unknown = 0;
68 
69  // Local storage for pointers to hang_info objects
70  HangInfo *hang_info_pt = 0, *hang_info2_pt = 0;
71 
72  // Local variable to determine the ALE stuff
73  bool ALE_is_disabled_flag = this->ALE_is_disabled;
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 function, initialise to zero
92  double dudt = 0.0;
93  double interpolated_u = 0.0;
94 
95  // These need to be a Vector to be ANSI C++, initialise to zero
96  Vector<double> interpolated_x(2, 0.0);
97  Vector<double> interpolated_dudx(2, 0.0);
98  Vector<double> mesh_velocity(2, 0.0);
99 
100  // Calculate function value and derivatives:
101  //-----------------------------------------
102 
103  // Loop over nodes
104  for (unsigned l = 0; l < n_node; l++)
105  {
106  // Get the value at the node
107  double u_value = this->nodal_value(l, u_nodal_index);
108  interpolated_u += u_value * psi(l);
109  dudt += this->du_dt_cons_axisym_adv_diff(l) * psi(l);
110  // Loop over directions
111  for (unsigned j = 0; j < 2; j++)
112  {
113  interpolated_x[j] += nodal_position(l, j) * psi(l);
114  interpolated_dudx[j] += u_value * dpsidx(l, j);
115  }
116  }
117 
118  // Get the mesh velocity, if required
119  if (!ALE_is_disabled_flag)
120  {
121  for (unsigned l = 0; l < n_node; l++)
122  {
123  // Loop over directions
124  for (unsigned j = 0; j < 2; j++)
125  {
126  mesh_velocity[j] += dnodal_position_dt(l, j) * psi(l);
127  }
128  }
129  }
130 
131 
132  // Get body force
133  double source;
135 
136 
137  // Get wind
138  //--------
139  Vector<double> wind(3);
141 
142  // Get the conserved wind (non-divergence free)
143  Vector<double> conserved_wind(3);
145  ipt, s, interpolated_x, conserved_wind);
146 
147  // Get diffusivity tensor
148  DenseMatrix<double> D(3, 3);
150 
151  // r is the first component of position
152  double r = interpolated_x[0];
153 
154  // Assemble residuals and Jacobian
155  //================================
156 
157  // Loop over the nodes for the test functions
158  for (unsigned l = 0; l < n_node; l++)
159  {
160  // Local variables to store the number of master nodes and
161  // the weight associated with the shape function if the node is hanging
162  unsigned n_master = 1;
163  double hang_weight = 1.0;
164  // Local bool (is the node hanging)
165  bool is_node_hanging = this->node_pt(l)->is_hanging();
166 
167  // If the node is hanging, get the number of master nodes
168  if (is_node_hanging)
169  {
170  hang_info_pt = this->node_pt(l)->hanging_pt();
171  n_master = hang_info_pt->nmaster();
172  }
173  // Otherwise there is just one master node, the node itself
174  else
175  {
176  n_master = 1;
177  }
178 
179  // Loop over the number of master nodes
180  for (unsigned m = 0; m < n_master; m++)
181  {
182  // Get the local equation number and hang_weight
183  // If the node is hanging
184  if (is_node_hanging)
185  {
186  // Read out the local equation from the master node
187  local_eqn = this->local_hang_eqn(hang_info_pt->master_node_pt(m),
188  u_nodal_index);
189  // Read out the weight from the master node
190  hang_weight = hang_info_pt->master_weight(m);
191  }
192  // If the node is not hanging
193  else
194  {
195  // The local equation number comes from the node itself
196  local_eqn = this->nodal_local_eqn(l, u_nodal_index);
197  // The hang weight is one
198  hang_weight = 1.0;
199  }
200 
201  // If the nodal equation is not a boundary conditino
202  if (local_eqn >= 0)
203  {
204  // Add du/dt and body force/source term here
205  residuals[local_eqn] -=
206  (peclet_st * dudt + source) * test(l) * r * W * hang_weight;
207 
208  // The Mesh velocity and GeneralisedAxisymAdvection--Diffusion bit
209  for (unsigned k = 0; k < 2; k++)
210  {
211  // Terms that multiply the test function
212  double tmp = peclet * wind[k];
213  // If the mesh is moving need to subtract the mesh velocity
214  if (!ALE_is_disabled_flag)
215  {
216  tmp -= peclet_st * mesh_velocity[k];
217  }
218  tmp *= interpolated_dudx[k];
219 
220  // Terms that multiply the derivative of the test function
221  // Advective term
222  double tmp2 = -conserved_wind[k] * interpolated_u;
223  // Now the diffusive term
224  for (unsigned j = 0; j < 2; j++)
225  {
226  tmp2 += D(k, j) * interpolated_dudx[j];
227  }
228  // Now construct the contribution to the residuals
229  residuals[local_eqn] -=
230  (tmp * test(l) + tmp2 * dtestdx(l, k)) * r * W * hang_weight;
231  }
232 
233  // Calculate the Jacobian
234  if (flag)
235  {
236  // Local variables to store the number of master nodes
237  // and the weights associated with each hanging node
238  unsigned n_master2 = 1;
239  double hang_weight2 = 1.0;
240  // Loop over the nodes for the variables
241  for (unsigned l2 = 0; l2 < n_node; l2++)
242  {
243  // Local bool (is the node hanging)
244  bool is_node2_hanging = this->node_pt(l2)->is_hanging();
245  // If the node is hanging, get the number of master nodes
246  if (is_node2_hanging)
247  {
248  hang_info2_pt = this->node_pt(l2)->hanging_pt();
249  n_master2 = hang_info2_pt->nmaster();
250  }
251  // Otherwise there is one master node, the node itself
252  else
253  {
254  n_master2 = 1;
255  }
256 
257  // Loop over the master nodes
258  for (unsigned m2 = 0; m2 < n_master2; m2++)
259  {
260  // Get the local unknown and weight
261  // If the node is hanging
262  if (is_node2_hanging)
263  {
264  // Read out the local unknown from the master node
265  local_unknown = this->local_hang_eqn(
266  hang_info2_pt->master_node_pt(m2), u_nodal_index);
267  // Read out the hanging weight from the master node
268  hang_weight2 = hang_info2_pt->master_weight(m2);
269  }
270  // If the node is not hanging
271  else
272  {
273  // The local unknown number comes from the node
274  local_unknown = this->nodal_local_eqn(l2, u_nodal_index);
275  // The hang weight is one
276  hang_weight2 = 1.0;
277  }
278 
279  // If the unknown is not pinned
280  if (local_unknown >= 0)
281  {
282  // Add contribution to Elemental Matrix
283  // Mass matrix du/dt term
284  jacobian(local_eqn, local_unknown) -=
285  peclet_st * test(l) * psi(l2) *
286  this->node_pt(l2)->time_stepper_pt()->weight(1, 0) * r *
287  W * hang_weight * hang_weight2;
288 
289  // Add contribution to mass matrix
290  if (flag == 2)
291  {
292  mass_matrix(local_eqn, local_unknown) +=
293  peclet_st * test(l) * psi(l2) * r * W * hang_weight *
294  hang_weight2;
295  }
296 
297  // Add contribution to Elemental Matrix
298  for (unsigned k = 0; k < 2; k++)
299  {
300  // Temporary term used in assembly
301  double tmp = peclet * wind[k];
302  if (!ALE_is_disabled_flag)
303  {
304  tmp -= peclet_st * mesh_velocity[k];
305  }
306  tmp *= dpsidx(l2, k);
307 
308  double tmp2 = -conserved_wind[k] * psi(l2);
309  // Now the diffusive term
310  for (unsigned j = 0; j < 2; j++)
311  {
312  tmp2 += D(k, j) * dpsidx(l2, j);
313  }
314 
315  // Now assemble Jacobian term
316  jacobian(local_eqn, local_unknown) -=
317  (tmp * test(l) + tmp2 * dtestdx(l, k)) * W * r *
318  hang_weight * hang_weight2;
319  }
320  }
321  } // End of loop over master nodes
322  } // End of loop over nodes
323  } // End of Jacobian calculation
324 
325  } // End of non-zero equation
326 
327  } // End of loop over the master nodes for residual
328  } // End of loop over nodes
329 
330  } // End of loop over integration points
331  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
dominoes D
Definition: Domino.cpp:55
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
RowVector3d w
Definition: Matrix_resize_int.cpp:3
MatrixType m2(n_dims)
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
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
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 dnodal_position_dt(const unsigned &n, const unsigned &i) const
Return the i-th component of nodal velocity: dx/dt at local node n.
Definition: elements.h:2333
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.
HangInfo *const & hanging_pt() const
Definition: nodes.h:1228
bool is_hanging() const
Test whether the node is geometrically hanging.
Definition: nodes.h:1285
int local_hang_eqn(Node *const &node_pt, const unsigned &i)
Definition: refineable_elements.h:278
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
Definition: timesteppers.h:594
RealScalar s
Definition: level1_cplx_impl.h:130
int * m
Definition: level2_cplx_impl.h:294
char char char int int * k
Definition: level2_impl.h:374
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > tmp
Definition: level3_impl.h:365
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
@ W
Definition: quadtree.h:63
virtual unsigned u_index_cons_axisym_adv_diff() const
Definition: gen_axisym_advection_diffusion_elements.h:110
virtual double dshape_and_dtest_eulerian_at_knot_cons_axisym_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
virtual void get_diff_cons_axisym_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, DenseMatrix< double > &D) const
Definition: gen_axisym_advection_diffusion_elements.h:389
const double & pe_st() const
Peclet number multiplied by Strouhal number.
Definition: gen_axisym_advection_diffusion_elements.h:296
virtual void get_conserved_wind_cons_axisym_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Definition: gen_axisym_advection_diffusion_elements.h:361
virtual void get_source_cons_axisym_adv_diff(const unsigned &ipt, const Vector< double > &x, double &source) const
Definition: gen_axisym_advection_diffusion_elements.h:311
virtual void get_wind_cons_axisym_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Definition: gen_axisym_advection_diffusion_elements.h:332
double du_dt_cons_axisym_adv_diff(const unsigned &n) const
Definition: gen_axisym_advection_diffusion_elements.h:117
bool ALE_is_disabled
Definition: gen_axisym_advection_diffusion_elements.h:638
const double & pe() const
Peclet number.
Definition: gen_axisym_advection_diffusion_elements.h:284
Definition: indexed_view.cpp:20
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::ALE_is_disabled, D, oomph::FiniteElement::dnodal_position_dt(), oomph::dshape_and_dtest_eulerian_at_knot_cons_axisym_adv_diff(), oomph::du_dt_cons_axisym_adv_diff(), oomph::get_conserved_wind_cons_axisym_adv_diff(), oomph::get_diff_cons_axisym_adv_diff(), oomph::get_source_cons_axisym_adv_diff(), oomph::get_wind_cons_axisym_adv_diff(), oomph::Node::hanging_pt(), i, oomph::FiniteElement::integral_pt(), oomph::FiniteElement::interpolated_x(), oomph::Node::is_hanging(), J, j, k, oomph::Integral::knot(), oomph::RefineableElement::local_hang_eqn(), m, m2(), oomph::HangInfo::master_node_pt(), oomph::HangInfo::master_weight(), oomph::HangInfo::nmaster(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_local_eqn(), oomph::FiniteElement::nodal_position(), oomph::FiniteElement::nodal_value(), oomph::FiniteElement::node_pt(), oomph::Integral::nweight(), oomph::pe(), oomph::pe_st(), UniformPSDSelfTest::r, s, TestProblem::source(), Eigen::test, oomph::Data::time_stepper_pt(), tmp, oomph::u_index_cons_axisym_adv_diff(), w, oomph::QuadTreeNames::W, oomph::Integral::weight(), and oomph::TimeStepper::weight().

◆ further_build()

void oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations::further_build ( )
inlinevirtual

Further build: Copy source function pointer from father element.

Reimplemented from oomph::RefineableElement.

171  {
173  cast_father_element_pt =
175  this->father_element_pt());
176 
177  // Set the values of the pointers from the father
178  this->Source_fct_pt = cast_father_element_pt->source_fct_pt();
179  this->Wind_fct_pt = cast_father_element_pt->wind_fct_pt();
180  this->Conserved_wind_fct_pt =
181  cast_father_element_pt->conserved_wind_fct_pt();
182  this->Diff_fct_pt = cast_father_element_pt->diff_fct_pt();
183  this->Pe_pt = cast_father_element_pt->pe_pt();
184  this->PeSt_pt = cast_father_element_pt->pe_st_pt();
185 
186  // Set the ALE status
187  this->ALE_is_disabled = cast_father_element_pt->ALE_is_disabled;
188  }
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
Definition: refineable_elements.h:539
RefineableGeneralisedAxisymAdvectionDiffusionEquations()
Empty Constructor.
Definition: refineable_gen_axisym_advection_diffusion_elements.h:61
GeneralisedAxisymAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
Definition: gen_axisym_advection_diffusion_elements.h:624
GeneralisedAxisymAdvectionDiffusionDiffFctPt Diff_fct_pt
Pointer to diffusivity funciton.
Definition: gen_axisym_advection_diffusion_elements.h:633
GeneralisedAxisymAdvectionDiffusionWindFctPt Conserved_wind_fct_pt
Pointer to additional (conservative) wind function:
Definition: gen_axisym_advection_diffusion_elements.h:630
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
Definition: gen_axisym_advection_diffusion_elements.h:621
GeneralisedAxisymAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
Definition: gen_axisym_advection_diffusion_elements.h:627
double * Pe_pt
Pointer to global Peclet number.
Definition: gen_axisym_advection_diffusion_elements.h:618

References oomph::ALE_is_disabled, oomph::Conserved_wind_fct_pt, oomph::Diff_fct_pt, oomph::RefineableElement::father_element_pt(), oomph::Pe_pt, oomph::PeSt_pt, oomph::Source_fct_pt, and oomph::Wind_fct_pt.

◆ geometric_jacobian()

double oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations::geometric_jacobian ( const Vector< double > &  x)
inlinevirtual

Fill in the geometric Jacobian, which in this case is r.

Reimplemented from oomph::ElementWithZ2ErrorEstimator.

164  {
165  return x[0];
166  }
list x
Definition: plotDoE.py:28

References plotDoE::x.

◆ get_interpolated_values() [1/2]

void oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations::get_interpolated_values ( const unsigned t,
const Vector< double > &  s,
Vector< double > &  values 
)
inlinevirtual

Get the function value u in Vector. Note: Given the generality of the interface (this function is usually called from black-box documentation or interpolation routines), the values Vector sets its own size in here.

Implements oomph::RefineableElement.

136  {
137  // Set size of Vector:
138  values.resize(1);
139 
140  // Find out how many nodes there are
141  const unsigned n_node = nnode();
142 
143  // Find the nodal index at which the unknown is stored
144  const unsigned u_nodal_index = this->u_index_cons_axisym_adv_diff();
145 
146  // Shape functions
147  Shape psi(n_node);
148 
149  // Find values of shape function
150  shape(s, psi);
151 
152  // Initialise the value of u
153  values[0] = 0.0;
154 
155  // Calculate value
156  for (unsigned l = 0; l < n_node; l++)
157  {
158  values[0] += this->nodal_value(t, l, u_nodal_index) * psi[l];
159  }
160  }
virtual void shape(const Vector< double > &s, Shape &psi) const =0

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

◆ get_interpolated_values() [2/2]

void oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations::get_interpolated_values ( const Vector< double > &  s,
Vector< double > &  values 
)
inlinevirtual

Get the function value u in Vector. Note: Given the generality of the interface (this function is usually called from black-box documentation or interpolation routines), the values Vector sets its own size in here.

Reimplemented from oomph::RefineableElement.

103  {
104  // Set size of Vector: u
105  values.resize(1);
106 
107  // Find number of nodes
108  const unsigned n_node = nnode();
109 
110  // Find the index at which the unknown is stored
111  const unsigned u_nodal_index = this->u_index_cons_axisym_adv_diff();
112 
113  // Local shape function
114  Shape psi(n_node);
115 
116  // Find values of shape function
117  shape(s, psi);
118 
119  // Initialise value of u
120  values[0] = 0.0;
121 
122  // Loop over the local nodes and sum
123  for (unsigned l = 0; l < n_node; l++)
124  {
125  values[0] += this->nodal_value(l, u_nodal_index) * psi[l];
126  }
127  }

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

◆ get_Z2_flux()

void oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations::get_Z2_flux ( const Vector< double > &  s,
Vector< double > &  flux 
)
inlinevirtual

Get 'flux' for Z2 error recovery: Standard flux.from GeneralisedAxisymAdvectionDiffusion equations

Implements oomph::ElementWithZ2ErrorEstimator.

92  {
93  this->get_flux(s, flux);
94  }
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
Definition: gen_axisym_advection_diffusion_elements.h:424

References ProblemParameters::flux(), and oomph::get_flux().

◆ num_Z2_flux_terms()

unsigned oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations::num_Z2_flux_terms ( )
inlinevirtual

Broken assignment operator.

Number of 'flux' terms for Z2 error estimation

Implements oomph::ElementWithZ2ErrorEstimator.

85  {
86  return 2;
87  }

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