oomph::GeneralisedAdvectionDiffusionEquations< DIM > Class Template Referenceabstract

#include <gen_advection_diffusion_elements.h>

+ Inheritance diagram for oomph::GeneralisedAdvectionDiffusionEquations< DIM >:

Public Types

typedef void(* GeneralisedAdvectionDiffusionSourceFctPt) (const Vector< double > &x, double &f)
 
typedef void(* GeneralisedAdvectionDiffusionWindFctPt) (const Vector< double > &x, Vector< double > &wind)
 
typedef void(* GeneralisedAdvectionDiffusionDiffFctPt) (const Vector< double > &x, DenseMatrix< double > &D)
 Funciton pointer to a diffusivity function. More...
 
- 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

 GeneralisedAdvectionDiffusionEquations ()
 
 GeneralisedAdvectionDiffusionEquations (const GeneralisedAdvectionDiffusionEquations &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const GeneralisedAdvectionDiffusionEquations &)=delete
 Broken assignment operator. More...
 
virtual unsigned u_index_cons_adv_diff () const
 
double du_dt_cons_adv_diff (const unsigned &n) const
 
void disable_ALE ()
 
void enable_ALE ()
 
void output (std::ostream &outfile)
 Output with default number of plot points. More...
 
void output (std::ostream &outfile, const unsigned &nplot)
 
void output (FILE *file_pt)
 C_style output with default number of plot points. More...
 
void output (FILE *file_pt, const unsigned &n_plot)
 
void output_fct (std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
 Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^DIM plot points. More...
 
virtual void output_fct (std::ostream &outfile, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
 
void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 Get error against and norm of exact solution. More...
 
void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Dummy, time dependent error checker. More...
 
double integrate_u ()
 Integrate the concentration over the element. More...
 
GeneralisedAdvectionDiffusionSourceFctPtsource_fct_pt ()
 Access function: Pointer to source function. More...
 
GeneralisedAdvectionDiffusionSourceFctPt source_fct_pt () const
 Access function: Pointer to source function. Const version. More...
 
GeneralisedAdvectionDiffusionWindFctPtwind_fct_pt ()
 Access function: Pointer to wind function. More...
 
GeneralisedAdvectionDiffusionWindFctPt wind_fct_pt () const
 Access function: Pointer to wind function. Const version. More...
 
GeneralisedAdvectionDiffusionWindFctPtconserved_wind_fct_pt ()
 Access function: Pointer to additional (conservative) wind function. More...
 
GeneralisedAdvectionDiffusionWindFctPt conserved_wind_fct_pt () const
 
GeneralisedAdvectionDiffusionDiffFctPtdiff_fct_pt ()
 Access function: Pointer to diffusion function. More...
 
GeneralisedAdvectionDiffusionDiffFctPt diff_fct_pt () const
 Access function: Pointer to diffusion function. Const version. More...
 
const doublepe () const
 Peclet number. More...
 
double *& pe_pt ()
 Pointer to Peclet number. More...
 
const doublepe_st () const
 Peclet number multiplied by Strouhal number. More...
 
double *& pe_st_pt ()
 Pointer to Peclet number multipled by Strouha number. More...
 
virtual void get_source_cons_adv_diff (const unsigned &ipt, const Vector< double > &x, double &source) const
 
virtual void get_wind_cons_adv_diff (const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
 
virtual void get_conserved_wind_cons_adv_diff (const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
 
virtual void get_diff_cons_adv_diff (const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, DenseMatrix< double > &D) const
 
void get_flux (const Vector< double > &s, Vector< double > &flux) const
 Get flux: \(\mbox{flux}[i] = \mbox{d}u / \mbox{d}x_i \). More...
 
void get_total_flux (const Vector< double > &s, Vector< double > &total_flux) const
 Get flux: \(\mbox{flux}[i] = \mbox{d}u / \mbox{d}x_i \). More...
 
void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Add the element's contribution to its residual vector (wrapper) More...
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
void fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
double interpolated_u_cons_adv_diff (const Vector< double > &s) const
 Return FE representation of function value u(s) at local coordinate s. More...
 
unsigned self_test ()
 Self-test: Return 0 for OK. More...
 
- Public Member Functions inherited from oomph::FiniteElement
void set_dimension (const unsigned &dim)
 
void set_nodal_dimension (const unsigned &nodal_dim)
 
void set_nnodal_position_type (const unsigned &nposition_type)
 Set the number of types required to interpolate the coordinate. More...
 
void set_n_node (const unsigned &n)
 
int nodal_local_eqn (const unsigned &n, const unsigned &i) const
 
double dJ_eulerian_at_knot (const unsigned &ipt, Shape &psi, DenseMatrix< double > &djacobian_dX) const
 
 FiniteElement ()
 Constructor. More...
 
virtual ~FiniteElement ()
 
 FiniteElement (const FiniteElement &)=delete
 Broken copy constructor. More...
 
virtual bool local_coord_is_valid (const Vector< double > &s)
 Broken assignment operator. More...
 
virtual void move_local_coord_back_into_element (Vector< double > &s) const
 
void get_centre_of_gravity_and_max_radius_in_terms_of_zeta (Vector< double > &cog, double &max_radius) const
 
virtual void local_coordinate_of_node (const unsigned &j, Vector< double > &s) const
 
virtual void local_fraction_of_node (const unsigned &j, Vector< double > &s_fraction)
 
virtual double local_one_d_fraction_of_node (const unsigned &n1d, const unsigned &i)
 
virtual void set_macro_elem_pt (MacroElement *macro_elem_pt)
 
MacroElementmacro_elem_pt ()
 Access function to pointer to macro element. More...
 
void get_x (const Vector< double > &s, Vector< double > &x) const
 
void get_x (const unsigned &t, const Vector< double > &s, Vector< double > &x)
 
virtual void get_x_from_macro_element (const Vector< double > &s, Vector< double > &x) const
 
virtual void get_x_from_macro_element (const unsigned &t, const Vector< double > &s, Vector< double > &x)
 
virtual void set_integration_scheme (Integral *const &integral_pt)
 Set the spatial integration scheme. More...
 
Integral *const & integral_pt () const
 Return the pointer to the integration scheme (const version) More...
 
virtual void shape (const Vector< double > &s, Shape &psi) const =0
 
virtual void shape_at_knot (const unsigned &ipt, Shape &psi) const
 
virtual void dshape_local (const Vector< double > &s, Shape &psi, DShape &dpsids) const
 
virtual void dshape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids) const
 
virtual void d2shape_local (const Vector< double > &s, Shape &psi, DShape &dpsids, DShape &d2psids) const
 
virtual void d2shape_local_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsids, DShape &d2psids) const
 
virtual double J_eulerian (const Vector< double > &s) const
 
virtual double J_eulerian_at_knot (const unsigned &ipt) const
 
void check_J_eulerian_at_knots (bool &passed) const
 
void check_jacobian (const double &jacobian) const
 
double dshape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx) const
 
virtual double dshape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsi, DenseMatrix< double > &djacobian_dX, RankFourTensor< double > &d_dpsidx_dX) const
 
double d2shape_eulerian (const Vector< double > &s, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual double d2shape_eulerian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
 
virtual void assign_nodal_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void describe_nodal_local_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void assign_all_generic_local_eqn_numbers (const bool &store_local_dof_pt)
 
Node *& node_pt (const unsigned &n)
 Return a pointer to the local node n. More...
 
Node *const & node_pt (const unsigned &n) const
 Return a pointer to the local node n (const version) More...
 
unsigned nnode () const
 Return the number of nodes. More...
 
virtual unsigned nnode_1d () const
 
double raw_nodal_position (const unsigned &n, const unsigned &i) const
 
double raw_nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &i) const
 
double raw_dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double raw_dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position (const unsigned &n, const unsigned &i) const
 
double nodal_position (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double dnodal_position_dt (const unsigned &n, const unsigned &i) const
 Return the i-th component of nodal velocity: dx/dt at local node n. More...
 
double dnodal_position_dt (const unsigned &n, const unsigned &j, const unsigned &i) const
 
double nodal_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double nodal_position_gen (const unsigned &t, const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &n, const unsigned &k, const unsigned &i) const
 
double dnodal_position_gen_dt (const unsigned &j, const unsigned &n, const unsigned &k, const unsigned &i) const
 
virtual void get_dresidual_dnodal_coordinates (RankThreeTensor< double > &dresidual_dnodal_coordinates)
 
virtual unsigned required_nvalue (const unsigned &n) const
 
unsigned nnodal_position_type () const
 
bool has_hanging_nodes () const
 
unsigned nodal_dimension () const
 Return the required Eulerian dimension of the nodes in this element. More...
 
virtual unsigned nvertex_node () const
 
virtual Nodevertex_node_pt (const unsigned &j) const
 
virtual Nodeconstruct_node (const unsigned &n)
 
virtual Nodeconstruct_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
virtual Nodeconstruct_boundary_node (const unsigned &n)
 
virtual Nodeconstruct_boundary_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 
int get_node_number (Node *const &node_pt) const
 
virtual Nodeget_node_at_local_coordinate (const Vector< double > &s) const
 
double raw_nodal_value (const unsigned &n, const unsigned &i) const
 
double raw_nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &n, const unsigned &i) const
 
double nodal_value (const unsigned &t, const unsigned &n, const unsigned &i) const
 
unsigned dim () const
 
virtual ElementGeometry::ElementGeometry element_geometry () const
 Return the geometry type of the element (either Q or T usually). More...
 
virtual double interpolated_x (const Vector< double > &s, const unsigned &i) const
 Return FE interpolated coordinate x[i] at local coordinate s. More...
 
virtual double interpolated_x (const unsigned &t, const Vector< double > &s, const unsigned &i) const
 
virtual void interpolated_x (const Vector< double > &s, Vector< double > &x) const
 Return FE interpolated position x[] at local coordinate s as Vector. More...
 
virtual void interpolated_x (const unsigned &t, const Vector< double > &s, Vector< double > &x) const
 
virtual double interpolated_dxdt (const Vector< double > &s, const unsigned &i, const unsigned &t)
 
virtual void interpolated_dxdt (const Vector< double > &s, const unsigned &t, Vector< double > &dxdt)
 
unsigned ngeom_data () const
 
Datageom_data_pt (const unsigned &j)
 
void position (const Vector< double > &zeta, Vector< double > &r) const
 
void position (const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const
 
void dposition_dt (const Vector< double > &zeta, const unsigned &t, Vector< double > &drdt)
 
virtual double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 
void interpolated_zeta (const Vector< double > &s, Vector< double > &zeta) const
 
void locate_zeta (const Vector< double > &zeta, GeomObject *&geom_object_pt, Vector< double > &s, const bool &use_coordinate_as_initial_guess=false)
 
virtual void node_update ()
 
virtual void identify_field_data_for_interactions (std::set< std::pair< Data *, unsigned >> &paired_field_data)
 
virtual void identify_geometric_data (std::set< Data * > &geometric_data_pt)
 
virtual double s_min () const
 Min value of local coordinate. More...
 
virtual double s_max () const
 Max. value of local coordinate. More...
 
double size () const
 
virtual double compute_physical_size () const
 
virtual void point_output_data (const Vector< double > &s, Vector< double > &data)
 
void point_output (std::ostream &outfile, const Vector< double > &s)
 
virtual unsigned nplot_points_paraview (const unsigned &nplot) const
 
virtual unsigned nsub_elements_paraview (const unsigned &nplot) const
 
void output_paraview (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_output_offset_information (std::ofstream &file_out, const unsigned &nplot, unsigned &counter) const
 
virtual void write_paraview_type (std::ofstream &file_out, const unsigned &nplot) const
 
virtual void write_paraview_offsets (std::ofstream &file_out, const unsigned &nplot, unsigned &offset_sum) const
 
virtual unsigned nscalar_paraview () const
 
virtual void scalar_value_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
 
virtual void scalar_value_fct_paraview (std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
 
virtual std::string scalar_name_paraview (const unsigned &i) const
 
virtual void output (const unsigned &t, std::ostream &outfile, const unsigned &n_plot) const
 
virtual void output_fct (std::ostream &outfile, const unsigned &n_plot, const double &time, const SolutionFunctorBase &exact_soln) const
 Output a time-dependent exact solution over the element. More...
 
virtual void get_s_plot (const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
 
virtual std::string tecplot_zone_string (const unsigned &nplot) const
 
virtual void write_tecplot_zone_footer (std::ostream &outfile, const unsigned &nplot) const
 
virtual void write_tecplot_zone_footer (FILE *file_pt, const unsigned &nplot) const
 
virtual unsigned nplot_points (const unsigned &nplot) const
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 Calculate the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Calculate the norm of the error and that of the exact solution. More...
 
virtual void compute_error (FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, Vector< double > &error, Vector< double > &norm)
 
virtual void compute_abs_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error)
 
void integrate_fct (FiniteElement::SteadyExactSolutionFctPt integrand_fct_pt, Vector< double > &integral)
 Integrate Vector-valued function over element. More...
 
void integrate_fct (FiniteElement::UnsteadyExactSolutionFctPt integrand_fct_pt, const double &time, Vector< double > &integral)
 Integrate Vector-valued time-dep function over element. More...
 
virtual void build_face_element (const int &face_index, FaceElement *face_element_pt)
 
virtual unsigned get_bulk_node_number (const int &face_index, const unsigned &i) const
 
virtual int face_outer_unit_normal_sign (const int &face_index) const
 Get the sign of the outer unit normal on the face given by face_index. More...
 
virtual unsigned nnode_on_face () const
 
void face_node_number_error_check (const unsigned &i) const
 Range check for face node numbers. More...
 
virtual CoordinateMappingFctPt face_to_bulk_coordinate_fct_pt (const int &face_index) const
 
virtual BulkCoordinateDerivativesFctPt bulk_coordinate_derivatives_fct_pt (const int &face_index) const
 
- Public Member Functions inherited from oomph::GeneralisedElement
 GeneralisedElement ()
 Constructor: Initialise all pointers and all values to zero. More...
 
virtual ~GeneralisedElement ()
 Virtual destructor to clean up any memory allocated by the object. More...
 
 GeneralisedElement (const GeneralisedElement &)=delete
 Broken copy constructor. More...
 
void operator= (const GeneralisedElement &)=delete
 Broken assignment operator. More...
 
Data *& internal_data_pt (const unsigned &i)
 Return a pointer to i-th internal data object. More...
 
Data *const & internal_data_pt (const unsigned &i) const
 Return a pointer to i-th internal data object (const version) More...
 
Data *& external_data_pt (const unsigned &i)
 Return a pointer to i-th external data object. More...
 
Data *const & external_data_pt (const unsigned &i) const
 Return a pointer to i-th external data object (const version) More...
 
unsigned long eqn_number (const unsigned &ieqn_local) const
 
int local_eqn_number (const unsigned long &ieqn_global) const
 
unsigned add_external_data (Data *const &data_pt, const bool &fd=true)
 
bool external_data_fd (const unsigned &i) const
 
void exclude_external_data_fd (const unsigned &i)
 
void include_external_data_fd (const unsigned &i)
 
void flush_external_data ()
 Flush all external data. More...
 
void flush_external_data (Data *const &data_pt)
 Flush the object addressed by data_pt from the external data array. More...
 
unsigned ninternal_data () const
 Return the number of internal data objects. More...
 
unsigned nexternal_data () const
 Return the number of external data objects. More...
 
unsigned ndof () const
 Return the number of equations/dofs in the element. More...
 
void dof_vector (const unsigned &t, Vector< double > &dof)
 Return the vector of dof values at time level t. More...
 
void dof_pt_vector (Vector< double * > &dof_pt)
 Return the vector of pointers to dof values. More...
 
void set_internal_data_time_stepper (const unsigned &i, TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void assign_internal_eqn_numbers (unsigned long &global_number, Vector< double * > &Dof_pt)
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
void add_internal_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt)
 
virtual void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void complete_setup_of_dependencies ()
 
virtual void get_residuals (Vector< double > &residuals)
 
virtual void get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
virtual void get_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void get_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void get_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void get_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void get_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void get_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void get_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void get_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_norm (double &norm)
 
virtual unsigned ndof_types () const
 
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 
- Public Member Functions inherited from oomph::GeomObject
 GeomObject ()
 Default constructor. More...
 
 GeomObject (const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim, TimeStepper *time_stepper_pt)
 
 GeomObject (const GeomObject &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const GeomObject &)=delete
 Broken assignment operator. More...
 
virtual ~GeomObject ()
 (Empty) destructor More...
 
unsigned nlagrangian () const
 Access function to # of Lagrangian coordinates. More...
 
unsigned ndim () const
 Access function to # of Eulerian coordinates. More...
 
void set_nlagrangian_and_ndim (const unsigned &n_lagrangian, const unsigned &n_dim)
 Set # of Lagrangian and Eulerian coordinates. More...
 
TimeStepper *& time_stepper_pt ()
 
TimeSteppertime_stepper_pt () const
 
virtual void position (const double &t, const Vector< double > &zeta, Vector< double > &r) const
 
virtual void dposition (const Vector< double > &zeta, DenseMatrix< double > &drdzeta) const
 
virtual void d2position (const Vector< double > &zeta, RankThreeTensor< double > &ddrdzeta) const
 
virtual void d2position (const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const
 

Protected Member Functions

virtual double dshape_and_dtest_eulerian_cons_adv_diff (const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
virtual double dshape_and_dtest_eulerian_at_knot_cons_adv_diff (const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
virtual void fill_in_generic_residual_contribution_cons_adv_diff (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
 
- Protected Member Functions inherited from oomph::FiniteElement
virtual void assemble_local_to_eulerian_jacobian (const DShape &dpsids, DenseMatrix< double > &jacobian) const
 
virtual void assemble_local_to_eulerian_jacobian2 (const DShape &d2psids, DenseMatrix< double > &jacobian2) const
 
virtual void assemble_eulerian_base_vectors (const DShape &dpsids, DenseMatrix< double > &interpolated_G) const
 
template<unsigned DIM>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual double invert_jacobian_mapping (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual double local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
double local_to_eulerian_mapping (const DShape &dpsids, DenseMatrix< double > &inverse_jacobian) const
 
virtual double local_to_eulerian_mapping_diagonal (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
virtual void dJ_eulerian_dnodal_coordinates (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<unsigned DIM>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
virtual void d_dshape_eulerian_dnodal_coordinates (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<unsigned DIM>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
virtual void transform_derivatives (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const
 
void transform_derivatives_diagonal (const DenseMatrix< double > &inverse_jacobian, DShape &dbasis) const
 
virtual void transform_second_derivatives (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<unsigned DIM>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<unsigned DIM>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
virtual void fill_in_jacobian_from_nodal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
void fill_in_jacobian_from_nodal_by_fd (DenseMatrix< double > &jacobian)
 
virtual void update_before_nodal_fd ()
 
virtual void reset_after_nodal_fd ()
 
virtual void update_in_nodal_fd (const unsigned &i)
 
virtual void reset_in_nodal_fd (const unsigned &i)
 
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Zero-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 One-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Two-d specialisation of function to calculate inverse of jacobian mapping. More...
 
template<>
double invert_jacobian (const DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void dJ_eulerian_dnodal_coordinates_templated_helper (const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void d_dshape_eulerian_dnodal_coordinates_templated_helper (const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
 
template<>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_template (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
template<>
void transform_second_derivatives_diagonal (const DenseMatrix< double > &jacobian, const DenseMatrix< double > &inverse_jacobian, const DenseMatrix< double > &jacobian2, DShape &dbasis, DShape &d2basis) const
 
- Protected Member Functions inherited from oomph::GeneralisedElement
unsigned add_internal_data (Data *const &data_pt, const bool &fd=true)
 
bool internal_data_fd (const unsigned &i) const
 
void exclude_internal_data_fd (const unsigned &i)
 
void include_internal_data_fd (const unsigned &i)
 
void clear_global_eqn_numbers ()
 
void add_global_eqn_numbers (std::deque< unsigned long > const &global_eqn_numbers, std::deque< double * > const &global_dof_pt)
 
virtual void assign_internal_and_external_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void assign_additional_local_eqn_numbers ()
 
int internal_local_eqn (const unsigned &i, const unsigned &j) const
 
int external_local_eqn (const unsigned &i, const unsigned &j)
 
void fill_in_jacobian_from_internal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_internal_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
virtual void update_before_internal_fd ()
 
virtual void reset_after_internal_fd ()
 
virtual void update_in_internal_fd (const unsigned &i)
 
virtual void reset_in_internal_fd (const unsigned &i)
 
virtual void update_before_external_fd ()
 
virtual void reset_after_external_fd ()
 
virtual void update_in_external_fd (const unsigned &i)
 
virtual void reset_in_external_fd (const unsigned &i)
 
virtual void fill_in_contribution_to_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void fill_in_contribution_to_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void fill_in_contribution_to_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void fill_in_contribution_to_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void fill_in_contribution_to_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void fill_in_contribution_to_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 

Protected Attributes

doublePe_pt
 Pointer to global Peclet number. More...
 
doublePeSt_pt
 Pointer to global Peclet number multiplied by Strouhal number. More...
 
GeneralisedAdvectionDiffusionSourceFctPt Source_fct_pt
 Pointer to source function: More...
 
GeneralisedAdvectionDiffusionWindFctPt Wind_fct_pt
 Pointer to wind function: More...
 
GeneralisedAdvectionDiffusionWindFctPt Conserved_wind_fct_pt
 Pointer to additional (conservative) wind function: More...
 
GeneralisedAdvectionDiffusionDiffFctPt Diff_fct_pt
 Pointer to diffusivity funciton. More...
 
bool ALE_is_disabled
 
- Protected Attributes inherited from oomph::FiniteElement
MacroElementMacro_elem_pt
 Pointer to the element's macro element (NULL by default) More...
 
- Protected Attributes inherited from oomph::GeomObject
unsigned NLagrangian
 Number of Lagrangian (intrinsic) coordinates. More...
 
unsigned Ndim
 Number of Eulerian coordinates. More...
 
TimeStepperGeom_object_time_stepper_pt
 

Static Private Attributes

static double Default_peclet_number
 Static default value for the Peclet number. More...
 

Additional Inherited Members

- Static Public Attributes inherited from oomph::FiniteElement
static double Tolerance_for_singular_jacobian = 1.0e-16
 Tolerance below which the jacobian is considered singular. More...
 
static bool Accept_negative_jacobian = false
 
static bool Suppress_output_while_checking_for_inverted_elements
 
- Static Public Attributes inherited from oomph::GeneralisedElement
static bool Suppress_warning_about_repeated_internal_data
 
static bool Suppress_warning_about_repeated_external_data = true
 
static double Default_fd_jacobian_step = 1.0e-8
 
- Static Protected Attributes inherited from oomph::FiniteElement
static const unsigned Default_Initial_Nvalue = 0
 Default value for the number of values at a node. More...
 
static const double Node_location_tolerance = 1.0e-14
 
static const unsigned N2deriv [] = {0, 1, 3, 6}
 
- Static Protected Attributes inherited from oomph::GeneralisedElement
static DenseMatrix< doubleDummy_matrix
 
static std::deque< double * > Dof_pt_deque
 

Detailed Description

template<unsigned DIM>
class oomph::GeneralisedAdvectionDiffusionEquations< DIM >

A class for all elements that solve the Advection Diffusion equations in conservative form using isoparametric elements.

\[ \frac{\partial}{\partial x_{i}}\left( Pe w_{i}(x_{k}) u - D_{ij}(x_{k})\frac{\partial u}{\partial x_{j}}\right) = f(x_{j}) \]

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

Member Typedef Documentation

◆ GeneralisedAdvectionDiffusionDiffFctPt

template<unsigned DIM>
typedef void(* oomph::GeneralisedAdvectionDiffusionEquations< DIM >::GeneralisedAdvectionDiffusionDiffFctPt) (const Vector< double > &x, DenseMatrix< double > &D)

Funciton pointer to a diffusivity function.

◆ GeneralisedAdvectionDiffusionSourceFctPt

template<unsigned DIM>
typedef void(* oomph::GeneralisedAdvectionDiffusionEquations< DIM >::GeneralisedAdvectionDiffusionSourceFctPt) (const Vector< double > &x, double &f)

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

◆ GeneralisedAdvectionDiffusionWindFctPt

template<unsigned DIM>
typedef void(* oomph::GeneralisedAdvectionDiffusionEquations< DIM >::GeneralisedAdvectionDiffusionWindFctPt) (const Vector< double > &x, Vector< double > &wind)

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

Constructor & Destructor Documentation

◆ GeneralisedAdvectionDiffusionEquations() [1/2]

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

72  : Source_fct_pt(0),
73  Wind_fct_pt(0),
75  Diff_fct_pt(0),
76  ALE_is_disabled(false)
77  {
78  // Set Peclet number to default
80  // Set Peclet Strouhal number to default
82  }
double * Pe_pt
Pointer to global Peclet number.
Definition: gen_advection_diffusion_elements.h:599
static double Default_peclet_number
Static default value for the Peclet number.
Definition: gen_advection_diffusion_elements.h:623
GeneralisedAdvectionDiffusionDiffFctPt Diff_fct_pt
Pointer to diffusivity funciton.
Definition: gen_advection_diffusion_elements.h:614
GeneralisedAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
Definition: gen_advection_diffusion_elements.h:608
GeneralisedAdvectionDiffusionWindFctPt Conserved_wind_fct_pt
Pointer to additional (conservative) wind function:
Definition: gen_advection_diffusion_elements.h:611
GeneralisedAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
Definition: gen_advection_diffusion_elements.h:605
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
Definition: gen_advection_diffusion_elements.h:602
bool ALE_is_disabled
Definition: gen_advection_diffusion_elements.h:619

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Default_peclet_number, oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Pe_pt, and oomph::GeneralisedAdvectionDiffusionEquations< DIM >::PeSt_pt.

◆ GeneralisedAdvectionDiffusionEquations() [2/2]

Broken copy constructor.

Member Function Documentation

◆ compute_error() [1/2]

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::compute_error ( std::ostream &  outfile,
FiniteElement::SteadyExactSolutionFctPt  exact_soln_pt,
double error,
double norm 
)
virtual

Get error against and norm of exact solution.

Validate against exact solution

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

Reimplemented from oomph::FiniteElement.

454  {
455  // Initialise
456  error = 0.0;
457  norm = 0.0;
458 
459  // Vector of local coordinates
460  Vector<double> s(DIM);
461 
462  // Vector for coordintes
463  Vector<double> x(DIM);
464 
465  // Find out how many nodes there are in the element
466  unsigned n_node = nnode();
467 
468  Shape psi(n_node);
469 
470  // Set the value of n_intpt
471  unsigned n_intpt = integral_pt()->nweight();
472 
473  // Tecplot header info
474  outfile << "ZONE" << std::endl;
475 
476  // Exact solution Vector (here a scalar)
477  Vector<double> exact_soln(1);
478 
479  // Loop over the integration points
480  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
481  {
482  // Assign values of s
483  for (unsigned i = 0; i < DIM; i++)
484  {
485  s[i] = integral_pt()->knot(ipt, i);
486  }
487 
488  // Get the integral weight
489  double w = integral_pt()->weight(ipt);
490 
491  // Get jacobian of mapping
492  double J = J_eulerian(s);
493 
494  // Premultiply the weights and the Jacobian
495  double W = w * J;
496 
497  // Get x position as Vector
498  interpolated_x(s, x);
499 
500  // Get FE function value
501  double u_fe = interpolated_u_cons_adv_diff(s);
502 
503  // Get exact solution at this point
504  (*exact_soln_pt)(x, exact_soln);
505 
506  // Output x,y,...,error
507  for (unsigned i = 0; i < DIM; i++)
508  {
509  outfile << x[i] << " ";
510  }
511  outfile << exact_soln[0] << " " << exact_soln[0] - u_fe << std::endl;
512 
513  // Add to error and norm
514  norm += exact_soln[0] * exact_soln[0] * W;
515  error += (exact_soln[0] - u_fe) * (exact_soln[0] - u_fe) * W;
516  }
517  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
RowVector3d w
Definition: Matrix_resize_int.cpp:3
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual double J_eulerian(const Vector< double > &s) const
Definition: elements.cc:4103
double interpolated_u_cons_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: gen_advection_diffusion_elements.h:541
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
void exact_soln(const double &time, const Vector< double > &x, Vector< double > &soln)
Definition: unstructured_two_d_curved.cc:301
int error
Definition: calibrate.py:297
@ W
Definition: quadtree.h:63
list x
Definition: plotDoE.py:28

References DIM, calibrate::error, ProblemParameters::exact_soln(), i, J, s, w, oomph::QuadTreeNames::W, and plotDoE::x.

◆ compute_error() [2/2]

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::compute_error ( std::ostream &  outfile,
FiniteElement::UnsteadyExactSolutionFctPt  exact_soln_pt,
const double time,
double error,
double norm 
)
inlinevirtual

Dummy, time dependent error checker.

Reimplemented from oomph::FiniteElement.

206  {
207  throw OomphLibError(
208  "No time-dependent compute_error() for Advection Diffusion elements",
211  }
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ conserved_wind_fct_pt() [1/2]

template<unsigned DIM>
GeneralisedAdvectionDiffusionWindFctPt& oomph::GeneralisedAdvectionDiffusionEquations< DIM >::conserved_wind_fct_pt ( )
inline

Access function: Pointer to additional (conservative) wind function.

247  {
248  return Conserved_wind_fct_pt;
249  }

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Conserved_wind_fct_pt.

Referenced by oomph::RefineableGeneralisedAdvectionDiffusionEquations< DIM >::further_build().

◆ conserved_wind_fct_pt() [2/2]

template<unsigned DIM>
GeneralisedAdvectionDiffusionWindFctPt oomph::GeneralisedAdvectionDiffusionEquations< DIM >::conserved_wind_fct_pt ( ) const
inline

Access function: Pointer to additoinal (conservative) wind function. Const version

256  {
257  return Conserved_wind_fct_pt;
258  }

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Conserved_wind_fct_pt.

◆ diff_fct_pt() [1/2]

Access function: Pointer to diffusion function.

262  {
263  return Diff_fct_pt;
264  }

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Diff_fct_pt.

Referenced by oomph::RefineableGeneralisedAdvectionDiffusionEquations< DIM >::further_build().

◆ diff_fct_pt() [2/2]

template<unsigned DIM>
GeneralisedAdvectionDiffusionDiffFctPt oomph::GeneralisedAdvectionDiffusionEquations< DIM >::diff_fct_pt ( ) const
inline

Access function: Pointer to diffusion function. Const version.

268  {
269  return Diff_fct_pt;
270  }

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Diff_fct_pt.

◆ disable_ALE()

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::disable_ALE ( )
inlinevirtual

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

Reimplemented from oomph::FiniteElement.

133  {
134  ALE_is_disabled = true;
135  }

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::ALE_is_disabled.

◆ dshape_and_dtest_eulerian_at_knot_cons_adv_diff()

template<unsigned DIM>
virtual double oomph::GeneralisedAdvectionDiffusionEquations< DIM >::dshape_and_dtest_eulerian_at_knot_cons_adv_diff ( const unsigned ipt,
Shape psi,
DShape dpsidx,
Shape test,
DShape dtestdx 
) const
protectedpure virtual

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

Implemented in oomph::QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D >.

◆ dshape_and_dtest_eulerian_cons_adv_diff()

template<unsigned DIM>
virtual double oomph::GeneralisedAdvectionDiffusionEquations< DIM >::dshape_and_dtest_eulerian_cons_adv_diff ( const Vector< double > &  s,
Shape psi,
DShape dpsidx,
Shape test,
DShape dtestdx 
) const
protectedpure virtual

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

Implemented in oomph::QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D >.

◆ du_dt_cons_adv_diff()

template<unsigned DIM>
double oomph::GeneralisedAdvectionDiffusionEquations< DIM >::du_dt_cons_adv_diff ( const unsigned n) const
inline

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

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

References oomph::TimeStepper::is_steady(), n, oomph::FiniteElement::nodal_value(), oomph::FiniteElement::node_pt(), oomph::TimeStepper::ntstorage(), plotPSD::t, oomph::GeomObject::time_stepper_pt(), oomph::Data::time_stepper_pt(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::u_index_cons_adv_diff(), and oomph::TimeStepper::weight().

◆ enable_ALE()

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::enable_ALE ( )
inlinevirtual

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

Reimplemented from oomph::FiniteElement.

143  {
144  ALE_is_disabled = false;
145  }

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::ALE_is_disabled.

◆ fill_in_contribution_to_jacobian()

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
inlinevirtual

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

Reimplemented from oomph::GeneralisedElement.

520  {
521  // Call the generic routine with the flag set to 1
523  residuals, jacobian, GeneralisedElement::Dummy_matrix, 1);
524  }
virtual void fill_in_generic_residual_contribution_cons_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: gen_advection_diffusion_elements.cc:49
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227

References oomph::GeneralisedElement::Dummy_matrix, and oomph::GeneralisedAdvectionDiffusionEquations< DIM >::fill_in_generic_residual_contribution_cons_adv_diff().

◆ fill_in_contribution_to_jacobian_and_mass_matrix()

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::fill_in_contribution_to_jacobian_and_mass_matrix ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  mass_matrix 
)
inlinevirtual

Add the element's contribution to its residuals vector, jacobian matrix and mass matrix

Reimplemented from oomph::GeneralisedElement.

533  {
534  // Call the generic routine with the flag set to 2
536  residuals, jacobian, mass_matrix, 2);
537  }

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::fill_in_generic_residual_contribution_cons_adv_diff().

◆ fill_in_contribution_to_residuals()

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::fill_in_contribution_to_residuals ( Vector< double > &  residuals)
inlinevirtual

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

Reimplemented from oomph::GeneralisedElement.

505  {
506  // Call the generic residuals function with flag set to 0 and using
507  // a dummy matrix
509  residuals,
512  0);
513  }

References oomph::GeneralisedElement::Dummy_matrix, and oomph::GeneralisedAdvectionDiffusionEquations< DIM >::fill_in_generic_residual_contribution_cons_adv_diff().

◆ fill_in_generic_residual_contribution_cons_adv_diff()

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::fill_in_generic_residual_contribution_cons_adv_diff ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  mass_matrix,
unsigned  flag 
)
protectedvirtual

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

Compute element residual Vector and/or element Jacobian matrix

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

Pure version without hanging nodes

Reimplemented in oomph::RefineableGeneralisedAdvectionDiffusionEquations< DIM >.

54  {
55  // Find out how many nodes there are
56  const unsigned n_node = nnode();
57 
58  // Get the nodal index at which the unknown is stored
59  const unsigned u_nodal_index = u_index_cons_adv_diff();
60 
61  // Set up memory for the shape and test functions
62  Shape psi(n_node), test(n_node);
63  DShape dpsidx(n_node, DIM), dtestdx(n_node, DIM);
64 
65  // Set the value of n_intpt
66  const unsigned n_intpt = integral_pt()->nweight();
67 
68  // Set the Vector to hold local coordinates
69  Vector<double> s(DIM);
70 
71  // Get Peclet number
72  const double peclet = pe();
73 
74  // Get the Peclet*Strouhal number
75  const double peclet_st = pe_st();
76 
77  // Integers used to store the local equation number and local unknown
78  // indices for the residuals and jacobians
79  int local_eqn = 0, local_unknown = 0;
80 
81  // Loop over the integration points
82  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
83  {
84  // Assign values of s
85  for (unsigned i = 0; i < DIM; i++) s[i] = integral_pt()->knot(ipt, i);
86 
87  // Get the integral weight
88  double w = integral_pt()->weight(ipt);
89 
90  // Call the derivatives of the shape and test functions
92  ipt, psi, dpsidx, test, dtestdx);
93 
94  // Premultiply the weights and the Jacobian
95  double W = w * J;
96 
97  // Calculate local values of the solution and its derivatives
98  // Allocate
99  double interpolated_u = 0.0;
100  double dudt = 0.0;
101  Vector<double> interpolated_x(DIM, 0.0);
102  Vector<double> interpolated_dudx(DIM, 0.0);
103  Vector<double> mesh_velocity(DIM, 0.0);
104 
105 
106  // Calculate function value and derivatives:
107  //-----------------------------------------
108  // Loop over nodes
109  for (unsigned l = 0; l < n_node; l++)
110  {
111  // Get the value at the node
112  double u_value = raw_nodal_value(l, u_nodal_index);
113  interpolated_u += u_value * psi(l);
114  dudt += du_dt_cons_adv_diff(l) * psi(l);
115  // Loop over directions
116  for (unsigned j = 0; j < DIM; j++)
117  {
118  interpolated_x[j] += raw_nodal_position(l, j) * psi(l);
119  interpolated_dudx[j] += u_value * dpsidx(l, j);
120  }
121  }
122 
123  // Mesh velocity?
124  if (!ALE_is_disabled)
125  {
126  for (unsigned l = 0; l < n_node; l++)
127  {
128  for (unsigned j = 0; j < DIM; j++)
129  {
130  mesh_velocity[j] += raw_dnodal_position_dt(l, j) * psi(l);
131  }
132  }
133  }
134 
135 
136  // Get source function
137  //-------------------
138  double source;
140 
141 
142  // Get wind
143  //--------
144  Vector<double> wind(DIM);
146 
147  // Get the conserved wind (non-divergence free)
148  Vector<double> conserved_wind(DIM);
149  get_conserved_wind_cons_adv_diff(ipt, s, interpolated_x, conserved_wind);
150 
151 
152  // Get diffusivity tensor
155 
156  // Assemble residuals and Jacobian
157  //--------------------------------
158 
159  // Loop over the test functions
160  for (unsigned l = 0; l < n_node; l++)
161  {
162  // Set the local equation number
163  local_eqn = nodal_local_eqn(l, u_nodal_index);
164 
165  /*IF it's not a boundary condition*/
166  if (local_eqn >= 0)
167  {
168  // Add body force/source term and time derivative
169  residuals[local_eqn] -= (peclet_st * dudt + source) * test(l) * W;
170 
171  // The Generalised Advection Diffusion bit itself
172  for (unsigned k = 0; k < DIM; k++)
173  {
174  // Terms that multiply the test function
175  // divergence-free wind
176  double tmp = peclet * wind[k];
177  // If the mesh is moving need to subtract the mesh velocity
178  if (!ALE_is_disabled)
179  {
180  tmp -= peclet_st * mesh_velocity[k];
181  }
182  tmp *= interpolated_dudx[k];
183 
184  // Terms that multiply the derivative of the test function
185  // Advective term
186  double tmp2 = -conserved_wind[k] * interpolated_u;
187  // Now the diuffusive term
188  for (unsigned j = 0; j < DIM; j++)
189  {
190  tmp2 += D(k, j) * interpolated_dudx[j];
191  }
192  // Now construct the contribution to the residuals
193  residuals[local_eqn] -= (tmp * test(l) + tmp2 * dtestdx(l, k)) * W;
194  }
195 
196  // Calculate the jacobian
197  //-----------------------
198  if (flag)
199  {
200  // Loop over the velocity shape functions again
201  for (unsigned l2 = 0; l2 < n_node; l2++)
202  {
203  // Set the number of the unknown
204  local_unknown = nodal_local_eqn(l2, u_nodal_index);
205 
206  // If at a non-zero degree of freedom add in the entry
207  if (local_unknown >= 0)
208  {
209  // Mass matrix term
210  jacobian(local_eqn, local_unknown) -=
211  peclet_st * test(l) * psi(l2) *
212  node_pt(l2)->time_stepper_pt()->weight(1, 0) * W;
213 
214  // Add the mass matrix term
215  if (flag == 2)
216  {
217  mass_matrix(local_eqn, local_unknown) +=
218  peclet_st * test(l) * psi(l2) * W;
219  }
220 
221  // Add contribution to Elemental Matrix
222  for (unsigned k = 0; k < DIM; k++)
223  {
224  // Temporary term used in assembly
225  double tmp = -peclet * wind[k];
226  if (!ALE_is_disabled)
227  {
228  tmp -= peclet_st * mesh_velocity[k];
229  }
230  tmp *= dpsidx(l2, k);
231 
232  double tmp2 = -conserved_wind[k] * psi(l2);
233  // Now the diffusive term
234  for (unsigned j = 0; j < DIM; j++)
235  {
236  tmp2 += D(k, j) * dpsidx(l2, j);
237  }
238 
239  // Now assemble Jacobian term
240  jacobian(local_eqn, local_unknown) -=
241  (tmp * test(l) + tmp2 * dtestdx(l, k)) * W;
242  }
243  }
244  }
245  }
246  }
247  }
248 
249 
250  } // End of loop over integration points
251  }
dominoes D
Definition: Domino.cpp:55
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
double raw_dnodal_position_dt(const unsigned &n, const unsigned &i) const
Definition: elements.h:2256
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2576
double raw_nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.cc:1686
virtual double dshape_and_dtest_eulerian_at_knot_cons_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
virtual void get_wind_cons_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Definition: gen_advection_diffusion_elements.h:321
virtual void get_diff_cons_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, DenseMatrix< double > &D) const
Definition: gen_advection_diffusion_elements.h:376
virtual void get_source_cons_adv_diff(const unsigned &ipt, const Vector< double > &x, double &source) const
Definition: gen_advection_diffusion_elements.h:300
const double & pe() const
Peclet number.
Definition: gen_advection_diffusion_elements.h:273
double du_dt_cons_adv_diff(const unsigned &n) const
Definition: gen_advection_diffusion_elements.h:105
const double & pe_st() const
Peclet number multiplied by Strouhal number.
Definition: gen_advection_diffusion_elements.h:285
virtual void get_conserved_wind_cons_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Definition: gen_advection_diffusion_elements.h:348
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
Definition: indexed_view.cpp:20
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::ALE_is_disabled, D, DIM, i, J, j, oomph::pe(), oomph::pe_st(), s, TestProblem::source(), Eigen::test, w, and oomph::QuadTreeNames::W.

Referenced by oomph::GeneralisedAdvectionDiffusionEquations< DIM >::fill_in_contribution_to_jacobian(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::fill_in_contribution_to_jacobian_and_mass_matrix(), and oomph::GeneralisedAdvectionDiffusionEquations< DIM >::fill_in_contribution_to_residuals().

◆ get_conserved_wind_cons_adv_diff()

template<unsigned DIM>
virtual void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::get_conserved_wind_cons_adv_diff ( const unsigned ipt,
const Vector< double > &  s,
const Vector< double > &  x,
Vector< double > &  wind 
) const
inlinevirtual

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

353  {
354  // If no wind function has been set, return zero
355  if (Conserved_wind_fct_pt == 0)
356  {
357  for (unsigned i = 0; i < DIM; i++)
358  {
359  wind[i] = 0.0;
360  }
361  }
362  else
363  {
364  // Get wind
365  (*Conserved_wind_fct_pt)(x, wind);
366  }
367  }

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Conserved_wind_fct_pt, DIM, i, and plotDoE::x.

Referenced by oomph::GeneralisedAdvectionDiffusionEquations< DIM >::get_total_flux().

◆ get_diff_cons_adv_diff()

template<unsigned DIM>
virtual void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::get_diff_cons_adv_diff ( const unsigned ipt,
const Vector< double > &  s,
const Vector< double > &  x,
DenseMatrix< double > &  D 
) const
inlinevirtual

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

380  {
381  // If no wind function has been set, return identity
382  if (Diff_fct_pt == 0)
383  {
384  for (unsigned i = 0; i < DIM; i++)
385  {
386  for (unsigned j = 0; j < DIM; j++)
387  {
388  if (i == j)
389  {
390  D(i, j) = 1.0;
391  }
392  else
393  {
394  D(i, j) = 0.0;
395  }
396  }
397  }
398  }
399  else
400  {
401  // Get diffusivity tensor
402  (*Diff_fct_pt)(x, D);
403  }
404  }

References D, oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Diff_fct_pt, DIM, i, j, and plotDoE::x.

Referenced by oomph::GeneralisedAdvectionDiffusionEquations< DIM >::get_total_flux().

◆ get_flux()

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::get_flux ( const Vector< double > &  s,
Vector< double > &  flux 
) const
inline

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

409  {
410  // Find out how many nodes there are in the element
411  unsigned n_node = nnode();
412 
413  // Get the nodal index at which the unknown is stored
414  unsigned u_nodal_index = u_index_cons_adv_diff();
415 
416  // Set up memory for the shape and test functions
417  Shape psi(n_node);
418  DShape dpsidx(n_node, DIM);
419 
420  // Call the derivatives of the shape and test functions
421  dshape_eulerian(s, psi, dpsidx);
422 
423  // Initialise to zero
424  for (unsigned j = 0; j < DIM; j++)
425  {
426  flux[j] = 0.0;
427  }
428 
429  // Loop over nodes
430  for (unsigned l = 0; l < n_node; l++)
431  {
432  // Loop over derivative directions
433  for (unsigned j = 0; j < DIM; j++)
434  {
435  flux[j] += nodal_value(l, u_nodal_index) * dpsidx(l, j);
436  }
437  }
438  }
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3298
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59

References DIM, oomph::FiniteElement::dshape_eulerian(), ProblemParameters::flux(), j, oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), s, and oomph::GeneralisedAdvectionDiffusionEquations< DIM >::u_index_cons_adv_diff().

Referenced by oomph::RefineableGeneralisedAdvectionDiffusionEquations< DIM >::get_Z2_flux().

◆ get_source_cons_adv_diff()

template<unsigned DIM>
virtual void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::get_source_cons_adv_diff ( const unsigned ipt,
const Vector< double > &  x,
double source 
) const
inlinevirtual

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

303  {
304  // If no source function has been set, return zero
305  if (Source_fct_pt == 0)
306  {
307  source = 0.0;
308  }
309  else
310  {
311  // Get source strength
312  (*Source_fct_pt)(x, source);
313  }
314  }

References TestProblem::source(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Source_fct_pt, and plotDoE::x.

◆ get_total_flux()

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::get_total_flux ( const Vector< double > &  s,
Vector< double > &  total_flux 
) const
inline

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

443  {
444  // Find out how many nodes there are in the element
445  const unsigned n_node = nnode();
446 
447  // Get the nodal index at which the unknown is stored
448  const unsigned u_nodal_index = u_index_cons_adv_diff();
449 
450  // Set up memory for the shape and test functions
451  Shape psi(n_node);
452  DShape dpsidx(n_node, DIM);
453 
454  // Call the derivatives of the shape and test functions
455  dshape_eulerian(s, psi, dpsidx);
456 
457  // Storage for the Eulerian position
458  Vector<double> interpolated_x(DIM, 0.0);
459  // Storage for the concentration
460  double interpolated_u = 0.0;
461  // Storage for the derivatives of the concentration
462  Vector<double> interpolated_dudx(DIM, 0.0);
463 
464  // Loop over nodes
465  for (unsigned l = 0; l < n_node; l++)
466  {
467  // Get the value at the node
468  const double u_value = this->nodal_value(l, u_nodal_index);
469  interpolated_u += u_value * psi(l);
470  // Loop over directions
471  for (unsigned j = 0; j < DIM; j++)
472  {
473  interpolated_x[j] += this->nodal_position(l, j) * psi(l);
474  interpolated_dudx[j] += u_value * dpsidx(l, j);
475  }
476  }
477 
478  // Dummy integration point
479  unsigned ipt = 0;
480 
481  // Get the conserved wind (non-divergence free)
482  Vector<double> conserved_wind(DIM);
483  get_conserved_wind_cons_adv_diff(ipt, s, interpolated_x, conserved_wind);
484 
485  // Get diffusivity tensor
488 
489  // Calculate the total flux made up of the diffusive flux
490  // and the conserved wind
491  for (unsigned i = 0; i < DIM; i++)
492  {
493  total_flux[i] = 0.0;
494  for (unsigned j = 0; j < DIM; j++)
495  {
496  total_flux[i] += D(i, j) * interpolated_dudx[j];
497  }
498  total_flux[i] -= conserved_wind[i] * interpolated_u;
499  }
500  }
double nodal_position(const unsigned &n, const unsigned &i) const
Definition: elements.h:2317

References D, DIM, oomph::FiniteElement::dshape_eulerian(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::get_conserved_wind_cons_adv_diff(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::get_diff_cons_adv_diff(), i, oomph::FiniteElement::interpolated_x(), j, oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_position(), oomph::FiniteElement::nodal_value(), s, and oomph::GeneralisedAdvectionDiffusionEquations< DIM >::u_index_cons_adv_diff().

◆ get_wind_cons_adv_diff()

template<unsigned DIM>
virtual void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::get_wind_cons_adv_diff ( const unsigned ipt,
const Vector< double > &  s,
const Vector< double > &  x,
Vector< double > &  wind 
) const
inlinevirtual

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

325  {
326  // If no wind function has been set, return zero
327  if (Wind_fct_pt == 0)
328  {
329  for (unsigned i = 0; i < DIM; i++)
330  {
331  wind[i] = 0.0;
332  }
333  }
334  else
335  {
336  // Get wind
337  (*Wind_fct_pt)(x, wind);
338  }
339  }

References DIM, i, oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Wind_fct_pt, and plotDoE::x.

◆ integrate_u()

template<unsigned DIM>
double oomph::GeneralisedAdvectionDiffusionEquations< DIM >::integrate_u

Integrate the concentration over the element.

Calculate the integrated value of the unknown over the element

525  {
526  // Initialise
527  double sum = 0.0;
528 
529  // Vector of local coordinates
530  Vector<double> s(DIM);
531 
532  // Find out how many nodes there are in the element
533  const unsigned n_node = nnode();
534 
535  // Find the index at which the concentration is stored
536  const unsigned u_nodal_index = this->u_index_cons_adv_diff();
537 
538  // Allocate memory for the shape functions
539  Shape psi(n_node);
540 
541  // Set the value of n_intpt
542  const unsigned n_intpt = integral_pt()->nweight();
543 
544  // Loop over the integration points
545  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
546  {
547  // Get the integral weight
548  const double w = integral_pt()->weight(ipt);
549 
550  // Get the shape functions
551  this->shape_at_knot(ipt, psi);
552 
553  // Calculate the concentration
554  double interpolated_u = 0.0;
555  for (unsigned l = 0; l < n_node; l++)
556  {
557  interpolated_u += this->nodal_value(l, u_nodal_index) * psi(l);
558  }
559 
560  // Get jacobian of mapping
561  const double J = J_eulerian_at_knot(ipt);
562 
563  // Add the values to the sum
564  sum += interpolated_u * w * J;
565  }
566 
567  // return the sum
568  return sum;
569  }
virtual double J_eulerian_at_knot(const unsigned &ipt) const
Definition: elements.cc:4168
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Definition: elements.cc:3220

References DIM, J, s, and w.

◆ interpolated_u_cons_adv_diff()

template<unsigned DIM>
double oomph::GeneralisedAdvectionDiffusionEquations< DIM >::interpolated_u_cons_adv_diff ( const Vector< double > &  s) const
inline

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

542  {
543  // Find number of nodes
544  unsigned n_node = nnode();
545 
546  // Get the nodal index at which the unknown is stored
547  unsigned u_nodal_index = u_index_cons_adv_diff();
548 
549  // Local shape function
550  Shape psi(n_node);
551 
552  // Find values of shape function
553  shape(s, psi);
554 
555  // Initialise value of u
556  double interpolated_u = 0.0;
557 
558  // Loop over the local nodes and sum
559  for (unsigned l = 0; l < n_node; l++)
560  {
561  interpolated_u += nodal_value(l, u_nodal_index) * psi[l];
562  }
563 
564  return (interpolated_u);
565  }
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::GeneralisedAdvectionDiffusionEquations< DIM >::u_index_cons_adv_diff().

◆ operator=()

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::operator= ( const GeneralisedAdvectionDiffusionEquations< DIM > &  )
delete

Broken assignment operator.

◆ output() [1/4]

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::output ( FILE *  file_pt)
inlinevirtual

C_style output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D >.

162  {
163  unsigned n_plot = 5;
164  output(file_pt, n_plot);
165  }
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: gen_advection_diffusion_elements.h:149

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::output().

◆ output() [2/4]

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::output ( FILE *  file_pt,
const unsigned nplot 
)
virtual

C-style output FE representation of soln: x,y,u or x,y,z,u at n_plot^DIM plot points

C-style output function:

x,y,u or x,y,z,u

nplot points in each coordinate direction

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D >.

363  {
364  // Vector of local coordinates
365  Vector<double> s(DIM);
366 
367  // Tecplot header info
368  fprintf(file_pt, "%s", tecplot_zone_string(nplot).c_str());
369 
370  // Loop over plot points
371  unsigned num_plot_points = nplot_points(nplot);
372  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
373  {
374  // Get local coordinates of plot point
375  get_s_plot(iplot, nplot, s);
376 
377  for (unsigned i = 0; i < DIM; i++)
378  {
379  fprintf(file_pt, "%g ", interpolated_x(s, i));
380  }
381  fprintf(file_pt, "%g \n", interpolated_u_cons_adv_diff(s));
382  }
383 
384  // Write tecplot footer (e.g. FE connectivity lists)
385  write_tecplot_zone_footer(file_pt, nplot);
386  }
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Definition: elements.h:3174

References DIM, i, and s.

◆ output() [3/4]

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

Output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D >.

150  {
151  unsigned nplot = 5;
152  output(outfile, nplot);
153  }

Referenced by oomph::GeneralisedAdvectionDiffusionEquations< DIM >::output(), and oomph::QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D >::output().

◆ output() [4/4]

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::output ( std::ostream &  outfile,
const unsigned nplot 
)
virtual

Output FE representation of soln: x,y,u or x,y,z,u at nplot^DIM plot points

Output function:

x,y,u,w_x,w_y or x,y,z,u,w_x,w_y,w_z

nplot points in each coordinate direction

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D >.

290  {
291  // Vector of local coordinates
292  Vector<double> s(DIM);
293 
294 
295  // Tecplot header info
296  outfile << tecplot_zone_string(nplot);
297 
298  const unsigned n_node = this->nnode();
299  Shape psi(n_node);
300  DShape dpsidx(n_node, DIM);
301 
302  // Loop over plot points
303  unsigned num_plot_points = nplot_points(nplot);
304  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
305  {
306  // Get local coordinates of plot point
307  get_s_plot(iplot, nplot, s);
308 
309  // Get Eulerian coordinate of plot point
310  Vector<double> x(DIM);
311  interpolated_x(s, x);
312 
313  for (unsigned i = 0; i < DIM; i++)
314  {
315  outfile << x[i] << " ";
316  }
317  outfile << interpolated_u_cons_adv_diff(s) << " ";
318 
319  // Get the gradients
320  (void)this->dshape_eulerian(s, psi, dpsidx);
321  Vector<double> interpolated_dudx(DIM, 0.0);
322  for (unsigned n = 0; n < n_node; n++)
323  {
324  const double u_ = this->nodal_value(n, 0);
325  for (unsigned i = 0; i < DIM; i++)
326  {
327  interpolated_dudx[i] += u_ * dpsidx(n, i);
328  }
329  }
330 
331  for (unsigned i = 0; i < DIM; i++)
332  {
333  outfile << interpolated_dudx[i] << " ";
334  }
335 
336  // Get the wind
337  Vector<double> wind(DIM);
338  // Dummy integration point variable needed
339  unsigned ipt = 0;
340  get_wind_cons_adv_diff(ipt, s, x, wind);
341  for (unsigned i = 0; i < DIM; i++)
342  {
343  outfile << wind[i] << " ";
344  }
345  outfile << std::endl;
346  }
347 
348  // Write tecplot footer (e.g. FE connectivity lists)
349  write_tecplot_zone_footer(outfile, nplot);
350  }

References DIM, i, n, s, and plotDoE::x.

◆ output_fct() [1/2]

template<unsigned DIM>
virtual void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::output_fct ( std::ostream &  outfile,
const unsigned nplot,
const double time,
FiniteElement::UnsteadyExactSolutionFctPt  exact_soln_pt 
)
inlinevirtual

Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^DIM plot points (dummy time-dependent version to keep intel compiler happy)

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D >.

185  {
186  throw OomphLibError("There is no time-dependent output_fct() for "
187  "Advection Diffusion elements",
190  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ output_fct() [2/2]

template<unsigned DIM>
void oomph::GeneralisedAdvectionDiffusionEquations< DIM >::output_fct ( std::ostream &  outfile,
const unsigned nplot,
FiniteElement::SteadyExactSolutionFctPt  exact_soln_pt 
)
virtual

Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^DIM plot points.

Output exact solution

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

x,y,u_exact or x,y,z,u_exact

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D >.

402  {
403  // Vector of local coordinates
404  Vector<double> s(DIM);
405 
406  // Vector for coordintes
407  Vector<double> x(DIM);
408 
409  // Tecplot header info
410  outfile << tecplot_zone_string(nplot);
411 
412  // Exact solution Vector (here a scalar)
413  Vector<double> exact_soln(1);
414 
415  // Loop over plot points
416  unsigned num_plot_points = nplot_points(nplot);
417  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
418  {
419  // Get local coordinates of plot point
420  get_s_plot(iplot, nplot, s);
421 
422  // Get x position as Vector
423  interpolated_x(s, x);
424 
425  // Get exact solution at this point
426  (*exact_soln_pt)(x, exact_soln);
427 
428  // Output x,y,...,u_exact
429  for (unsigned i = 0; i < DIM; i++)
430  {
431  outfile << x[i] << " ";
432  }
433  outfile << exact_soln[0] << std::endl;
434  }
435 
436  // Write tecplot footer (e.g. FE connectivity lists)
437  write_tecplot_zone_footer(outfile, nplot);
438  }

References DIM, ProblemParameters::exact_soln(), i, s, and plotDoE::x.

Referenced by oomph::QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D >::output_fct().

◆ pe()

template<unsigned DIM>
const double& oomph::GeneralisedAdvectionDiffusionEquations< DIM >::pe ( ) const
inline

Peclet number.

274  {
275  return *Pe_pt;
276  }

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Pe_pt.

◆ pe_pt()

◆ pe_st()

template<unsigned DIM>
const double& oomph::GeneralisedAdvectionDiffusionEquations< DIM >::pe_st ( ) const
inline

Peclet number multiplied by Strouhal number.

286  {
287  return *PeSt_pt;
288  }

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::PeSt_pt.

◆ pe_st_pt()

template<unsigned DIM>
double*& oomph::GeneralisedAdvectionDiffusionEquations< DIM >::pe_st_pt ( )
inline

Pointer to Peclet number multipled by Strouha number.

292  {
293  return PeSt_pt;
294  }

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::PeSt_pt.

Referenced by oomph::RefineableGeneralisedAdvectionDiffusionEquations< DIM >::further_build().

◆ self_test()

template<unsigned DIM>
unsigned oomph::GeneralisedAdvectionDiffusionEquations< DIM >::self_test
virtual

Self-test: Return 0 for OK.

Reimplemented from oomph::FiniteElement.

259  {
260  bool passed = true;
261 
262  // Check lower-level stuff
263  if (FiniteElement::self_test() != 0)
264  {
265  passed = false;
266  }
267 
268  // Return verdict
269  if (passed)
270  {
271  return 0;
272  }
273  else
274  {
275  return 1;
276  }
277  }
virtual unsigned self_test()
Definition: elements.cc:4440

References oomph::FiniteElement::self_test().

◆ source_fct_pt() [1/2]

◆ source_fct_pt() [2/2]

template<unsigned DIM>
GeneralisedAdvectionDiffusionSourceFctPt oomph::GeneralisedAdvectionDiffusionEquations< DIM >::source_fct_pt ( ) const
inline

Access function: Pointer to source function. Const version.

226  {
227  return Source_fct_pt;
228  }

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Source_fct_pt.

◆ u_index_cons_adv_diff()

template<unsigned DIM>
virtual unsigned oomph::GeneralisedAdvectionDiffusionEquations< DIM >::u_index_cons_adv_diff ( ) const
inlinevirtual

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

99  {
100  return 0;
101  }

Referenced by oomph::GeneralisedAdvectionDiffusionEquations< DIM >::du_dt_cons_adv_diff(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::get_flux(), oomph::RefineableGeneralisedAdvectionDiffusionEquations< DIM >::get_interpolated_values(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::get_total_flux(), and oomph::GeneralisedAdvectionDiffusionEquations< DIM >::interpolated_u_cons_adv_diff().

◆ wind_fct_pt() [1/2]

◆ wind_fct_pt() [2/2]

template<unsigned DIM>
GeneralisedAdvectionDiffusionWindFctPt oomph::GeneralisedAdvectionDiffusionEquations< DIM >::wind_fct_pt ( ) const
inline

Access function: Pointer to wind function. Const version.

240  {
241  return Wind_fct_pt;
242  }

References oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Wind_fct_pt.

Member Data Documentation

◆ ALE_is_disabled

template<unsigned DIM>
bool oomph::GeneralisedAdvectionDiffusionEquations< DIM >::ALE_is_disabled
protected

Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed. Only set to false if you're sure that the mesh is stationary.

Referenced by oomph::GeneralisedAdvectionDiffusionEquations< DIM >::disable_ALE(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::enable_ALE(), and oomph::RefineableGeneralisedAdvectionDiffusionEquations< DIM >::further_build().

◆ Conserved_wind_fct_pt

◆ Default_peclet_number

template<unsigned DIM>
double oomph::GeneralisedAdvectionDiffusionEquations< DIM >::Default_peclet_number
staticprivate
Initial value:
=
0.0

Static default value for the Peclet number.

2D GeneralisedAdvection Diffusion elements

Default value for Peclet number

Referenced by oomph::GeneralisedAdvectionDiffusionEquations< DIM >::GeneralisedAdvectionDiffusionEquations().

◆ Diff_fct_pt

◆ Pe_pt

◆ PeSt_pt

◆ Source_fct_pt

◆ Wind_fct_pt


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