oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM > Class Template Referenceabstract

#include <advection_diffusion_reaction_elements.h>

+ Inheritance diagram for oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >:

Public Types

typedef void(* AdvectionDiffusionReactionSourceFctPt) (const Vector< double > &x, Vector< double > &f)
 
typedef void(* AdvectionDiffusionReactionReactionFctPt) (const Vector< double > &c, Vector< double > &R)
 Function pointer to reaction terms. More...
 
typedef void(* AdvectionDiffusionReactionReactionDerivFctPt) (const Vector< double > &c, DenseMatrix< double > &dRdC)
 Function pointer to derivative of reaction terms. More...
 
typedef void(* AdvectionDiffusionReactionWindFctPt) (const double &time, const Vector< double > &x, Vector< double > &wind)
 
- Public Types inherited from oomph::FiniteElement
typedef void(* SteadyExactSolutionFctPt) (const Vector< double > &, Vector< double > &)
 
typedef void(* UnsteadyExactSolutionFctPt) (const double &, const Vector< double > &, Vector< double > &)
 

Public Member Functions

 AdvectionDiffusionReactionEquations ()
 
 AdvectionDiffusionReactionEquations (const AdvectionDiffusionReactionEquations &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const AdvectionDiffusionReactionEquations &)=delete
 Broken assignment operator. More...
 
virtual unsigned c_index_adv_diff_react (const unsigned &i) const
 
double dc_dt_adv_diff_react (const unsigned &n, const unsigned &r) const
 
void dc_dt_adv_diff_react (const unsigned &n, Vector< double > &dc_dt) 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_norm (double &norm)
 Compute norm of the solution: sum of squares of L2 norms for reagents. More...
 
void compute_error (std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
 Get error against and norm of exact solution. More...
 
void compute_error (std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
 Dummy, time dependent error checker. More...
 
AdvectionDiffusionReactionSourceFctPtsource_fct_pt ()
 Access function: Pointer to source function. More...
 
AdvectionDiffusionReactionSourceFctPt source_fct_pt () const
 Access function: Pointer to source function. Const version. More...
 
AdvectionDiffusionReactionWindFctPtwind_fct_pt ()
 Access function: Pointer to wind function. More...
 
AdvectionDiffusionReactionWindFctPt wind_fct_pt () const
 Access function: Pointer to reaction function. Const version. More...
 
AdvectionDiffusionReactionReactionFctPtreaction_fct_pt ()
 Access function: Pointer to reaction function. More...
 
AdvectionDiffusionReactionReactionFctPt reaction_fct_pt () const
 Access function: Pointer to reaction function. Const version. More...
 
AdvectionDiffusionReactionReactionDerivFctPtreaction_deriv_fct_pt ()
 Access function: Pointer to reaction derivatives function. More...
 
AdvectionDiffusionReactionReactionDerivFctPt reaction_deriv_fct_pt () const
 Access function: Pointer to reaction function. Const version. More...
 
const Vector< double > & diff () const
 Vector of diffusion coefficients. More...
 
Vector< double > *& diff_pt ()
 Pointer to vector of diffusion coefficients. More...
 
const Vector< double > & tau () const
 Vector of dimensionless timescales. More...
 
Vector< double > *& tau_pt ()
 Pointer to vector of dimensionless timescales. More...
 
virtual void get_source_adv_diff_react (const unsigned &ipt, const Vector< double > &x, Vector< double > &source) const
 
virtual void get_wind_adv_diff_react (const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
 
virtual void get_reaction_adv_diff_react (const unsigned &ipt, const Vector< double > &C, Vector< double > &R) const
 
virtual void get_reaction_deriv_adv_diff_react (const unsigned &ipt, const Vector< double > &C, DenseMatrix< double > &dRdC) const
 
void get_flux (const Vector< double > &s, Vector< double > &flux) const
 Get flux: \(\mbox{flux}[DIM r + i] = \mbox{d}C_{r} / \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_c_adv_diff_react (const Vector< double > &s, const unsigned &i) const
 Return FE representation of function value c_i(s) at local coordinate s. More...
 
unsigned self_test ()
 Self-test: Return 0 for OK. More...
 
void integrate_reagents (Vector< double > &C) const
 Return the integrated reagent concentrations. 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 unsigned ndof_types () const
 
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 
- Public Member Functions inherited from oomph::GeomObject
 GeomObject ()
 Default constructor. More...
 
 GeomObject (const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim)
 
 GeomObject (const unsigned &nlagrangian, const unsigned &ndim, TimeStepper *time_stepper_pt)
 
 GeomObject (const GeomObject &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const GeomObject &)=delete
 Broken assignment operator. More...
 
virtual ~GeomObject ()
 (Empty) destructor More...
 
unsigned nlagrangian () const
 Access function to # of Lagrangian coordinates. More...
 
unsigned ndim () const
 Access function to # of Eulerian coordinates. More...
 
void set_nlagrangian_and_ndim (const unsigned &n_lagrangian, const unsigned &n_dim)
 Set # of Lagrangian and Eulerian coordinates. More...
 
TimeStepper *& time_stepper_pt ()
 
TimeSteppertime_stepper_pt () const
 
virtual void position (const double &t, const Vector< double > &zeta, Vector< double > &r) const
 
virtual void dposition (const Vector< double > &zeta, DenseMatrix< double > &drdzeta) const
 
virtual void d2position (const Vector< double > &zeta, RankThreeTensor< double > &ddrdzeta) const
 
virtual void d2position (const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const
 

Protected Member Functions

virtual double dshape_and_dtest_eulerian_adv_diff_react (const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
virtual double dshape_and_dtest_eulerian_at_knot_adv_diff_react (const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
 
virtual void fill_in_generic_residual_contribution_adv_diff_react (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

Vector< double > * Diff_pt
 Pointer to global diffusion coefficients. More...
 
Vector< double > * Tau_pt
 Pointer to global timescales. More...
 
AdvectionDiffusionReactionSourceFctPt Source_fct_pt
 Pointer to source function: More...
 
AdvectionDiffusionReactionWindFctPt Wind_fct_pt
 Pointer to wind function: More...
 
AdvectionDiffusionReactionReactionFctPt Reaction_fct_pt
 Pointer to reaction function. More...
 
AdvectionDiffusionReactionReactionDerivFctPt Reaction_deriv_fct_pt
 Pointer to reaction derivatives. 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 Protected Attributes

static const unsigned N_reagent
 
- 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
 

Static Private Attributes

static Vector< doubleDefault_dimensionless_number
 Static default value for the dimensionless numbers. 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
 

Detailed Description

template<unsigned NREAGENT, unsigned DIM>
class oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >

A class for all elements that solve the Advection Diffusion Reaction equations using isoparametric elements.

\[ \tau_{i} \frac{\partial C_{i}}{\partial t} + w_{j} \frac{\partial C_{i}}{\partial x_{j}} = D_{i}\frac{\partial^2 C_{i}}{\partial x_j^2} = - R_{i}(C_{i}) + f_{i} \]

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

Member Typedef Documentation

◆ AdvectionDiffusionReactionReactionDerivFctPt

template<unsigned NREAGENT, unsigned DIM>
typedef void(* oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::AdvectionDiffusionReactionReactionDerivFctPt) (const Vector< double > &c, DenseMatrix< double > &dRdC)

Function pointer to derivative of reaction terms.

◆ AdvectionDiffusionReactionReactionFctPt

template<unsigned NREAGENT, unsigned DIM>
typedef void(* oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::AdvectionDiffusionReactionReactionFctPt) (const Vector< double > &c, Vector< double > &R)

Function pointer to reaction terms.

◆ AdvectionDiffusionReactionSourceFctPt

template<unsigned NREAGENT, unsigned DIM>
typedef void(* oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::AdvectionDiffusionReactionSourceFctPt) (const Vector< double > &x, Vector< double > &f)

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

◆ AdvectionDiffusionReactionWindFctPt

template<unsigned NREAGENT, unsigned DIM>
typedef void(* oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::AdvectionDiffusionReactionWindFctPt) (const double &time, const Vector< double > &x, Vector< double > &wind)

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

Constructor & Destructor Documentation

◆ AdvectionDiffusionReactionEquations() [1/2]

template<unsigned NREAGENT, unsigned DIM>
oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::AdvectionDiffusionReactionEquations ( )
inline

Constructor: Initialise the Source_fct_pt, Wind_fct_pt, Reaction_fct_pt to null and initialise the dimensionless timescale and diffusion ratios

80  : Source_fct_pt(0),
81  Wind_fct_pt(0),
82  Reaction_fct_pt(0),
84  ALE_is_disabled(false)
85  {
86  // Set diffusion coefficients to default
88  // Set timescales to default
90  }
AdvectionDiffusionReactionSourceFctPt Source_fct_pt
Pointer to source function:
Definition: advection_diffusion_reaction_elements.h:628
bool ALE_is_disabled
Definition: advection_diffusion_reaction_elements.h:642
AdvectionDiffusionReactionReactionDerivFctPt Reaction_deriv_fct_pt
Pointer to reaction derivatives.
Definition: advection_diffusion_reaction_elements.h:637
Vector< double > * Tau_pt
Pointer to global timescales.
Definition: advection_diffusion_reaction_elements.h:625
AdvectionDiffusionReactionWindFctPt Wind_fct_pt
Pointer to wind function:
Definition: advection_diffusion_reaction_elements.h:631
static Vector< double > Default_dimensionless_number
Static default value for the dimensionless numbers.
Definition: advection_diffusion_reaction_elements.h:646
AdvectionDiffusionReactionReactionFctPt Reaction_fct_pt
Pointer to reaction function.
Definition: advection_diffusion_reaction_elements.h:634
Vector< double > * Diff_pt
Pointer to global diffusion coefficients.
Definition: advection_diffusion_reaction_elements.h:622

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Default_dimensionless_number, oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Diff_pt, and oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Tau_pt.

◆ AdvectionDiffusionReactionEquations() [2/2]

template<unsigned NREAGENT, unsigned DIM>
oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::AdvectionDiffusionReactionEquations ( const AdvectionDiffusionReactionEquations< NREAGENT, DIM > &  dummy)
delete

Broken copy constructor.

Member Function Documentation

◆ c_index_adv_diff_react()

template<unsigned NREAGENT, unsigned DIM>
virtual unsigned oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::c_index_adv_diff_react ( const unsigned i) const
inlinevirtual

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

Reimplemented in oomph::RefineableDoubleBuoyantQCrouzeixRaviartElement< DIM >, oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >, oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >, and oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >.

107  {
108  return i;
109  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9

References i.

Referenced by oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::dc_dt_adv_diff_react(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::get_flux(), oomph::RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM >::get_interpolated_values(), and oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::interpolated_c_adv_diff_react().

◆ compute_error() [1/2]

template<unsigned NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, 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.

575  {
576  // Initialise
577  error = 0.0;
578  norm = 0.0;
579 
580  // Vector of local coordinates
581  Vector<double> s(DIM);
582 
583  // Vector for coordintes
584  Vector<double> x(DIM);
585 
586  // Find out how many nodes there are in the element
587  unsigned n_node = nnode();
588 
589  Shape psi(n_node);
590 
591  // Set the value of n_intpt
592  unsigned n_intpt = integral_pt()->nweight();
593 
594  // Tecplot header info
595  outfile << "ZONE" << std::endl;
596 
597  // Exact solution Vector (here a scalar)
598  Vector<double> exact_soln(1);
599 
600  // Loop over the integration points
601  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
602  {
603  // Assign values of s
604  for (unsigned i = 0; i < DIM; i++)
605  {
606  s[i] = integral_pt()->knot(ipt, i);
607  }
608 
609  // Get the integral weight
610  double w = integral_pt()->weight(ipt);
611 
612  // Get jacobian of mapping
613  double J = J_eulerian(s);
614 
615  // Premultiply the weights and the Jacobian
616  double W = w * J;
617 
618  // Get x position as Vector
619  interpolated_x(s, x);
620 
621  // Get FE function value
622  double u_fe = interpolated_c_adv_diff_react(s, 0);
623 
624  // Get exact solution at this point
625  (*exact_soln_pt)(x, exact_soln);
626 
627  // Output x,y,...,error
628  for (unsigned i = 0; i < DIM; i++)
629  {
630  outfile << x[i] << " ";
631  }
632  outfile << exact_soln[0] << " " << exact_soln[0] - u_fe << std::endl;
633 
634  // Add to error and norm
635  norm += exact_soln[0] * exact_soln[0] * W;
636  error += (exact_soln[0] - u_fe) * (exact_soln[0] - u_fe) * W;
637  }
638  }
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
RowVector3d w
Definition: Matrix_resize_int.cpp:3
double interpolated_c_adv_diff_react(const Vector< double > &s, const unsigned &i) const
Return FE representation of function value c_i(s) at local coordinate s.
Definition: advection_diffusion_reaction_elements.h:556
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual double J_eulerian(const Vector< double > &s) const
Definition: elements.cc:4103
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
RealScalar s
Definition: level1_cplx_impl.h:130
#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 NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, 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.

256  {
257  throw OomphLibError(
258  "No time-dependent compute_error() for Advection Diffusion elements",
261  }
#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.

◆ compute_norm()

template<unsigned NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::compute_norm ( double norm)
virtual

Compute norm of the solution: sum of squares of L2 norms for reagents.

Compute norm of the solution: sum of squares of the L2 norm for each reagent

Reimplemented from oomph::GeneralisedElement.

279  {
280  // Initialise
281  norm = 0.0;
282 
283  // Vector of local coordinates
284  Vector<double> s(DIM);
285 
286  // Solution
287  double c = 0.0;
288 
289  // Find out how many nodes there are in the element
290  unsigned n_node = this->nnode();
291 
292  Shape psi(n_node);
293 
294  // Set the value of n_intpt
295  unsigned n_intpt = this->integral_pt()->nweight();
296 
297  // Loop over the integration points
298  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
299  {
300  // Assign values of s
301  for (unsigned i = 0; i < DIM; i++)
302  {
303  s[i] = this->integral_pt()->knot(ipt, i);
304  }
305 
306  // Get the integral weight
307  double w = this->integral_pt()->weight(ipt);
308 
309  // Get jacobian of mapping
310  double J = this->J_eulerian(s);
311 
312  // Premultiply the weights and the Jacobian
313  double W = w * J;
314 
315  // Loop over all reagents
316  for (unsigned r = 0; r < NREAGENT; ++r)
317  {
318  // Get FE function value
320  // Add to norm
321  norm += c * c * W;
322  }
323  }
324  }
r
Definition: UniformPSDSelfTest.py:20
int c
Definition: calibrate.py:100

References calibrate::c, DIM, i, J, UniformPSDSelfTest::r, s, w, and oomph::QuadTreeNames::W.

◆ dc_dt_adv_diff_react() [1/2]

template<unsigned NREAGENT, unsigned DIM>
double oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::dc_dt_adv_diff_react ( const unsigned n,
const unsigned r 
) const
inline

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

114  {
115  // Get the data's timestepper
116  TimeStepper* time_stepper_pt = this->node_pt(n)->time_stepper_pt();
117 
118  // Initialise dudt
119  double dudt = 0.0;
120  // Loop over the timesteps, if there is a non Steady timestepper
121  if (!time_stepper_pt->is_steady())
122  {
123  // Find the index at which the variable is stored
124  const unsigned c_nodal_index = c_index_adv_diff_react(r);
125 
126  // Number of timsteps (past & present)
127  const unsigned n_time = time_stepper_pt->ntstorage();
128 
129  for (unsigned t = 0; t < n_time; t++)
130  {
131  dudt +=
132  time_stepper_pt->weight(1, t) * nodal_value(t, n, c_nodal_index);
133  }
134  }
135  return dudt;
136  }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
virtual unsigned c_index_adv_diff_react(const unsigned &i) const
Definition: advection_diffusion_reaction_elements.h:106
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
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::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::c_index_adv_diff_react(), oomph::TimeStepper::is_steady(), n, oomph::FiniteElement::nodal_value(), oomph::FiniteElement::node_pt(), oomph::TimeStepper::ntstorage(), UniformPSDSelfTest::r, plotPSD::t, oomph::GeomObject::time_stepper_pt(), oomph::Data::time_stepper_pt(), and oomph::TimeStepper::weight().

◆ dc_dt_adv_diff_react() [2/2]

template<unsigned NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::dc_dt_adv_diff_react ( const unsigned n,
Vector< double > &  dc_dt 
) const
inline

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

142  {
143  // Get the data's timestepper
144  TimeStepper* time_stepper_pt = this->node_pt(n)->time_stepper_pt();
145 
146  // Initialise to zero
147  for (unsigned r = 0; r < NREAGENT; r++)
148  {
149  dc_dt[r] = 0.0;
150  }
151 
152  // Loop over the timesteps, if there is a non Steady timestepper
153  if (!time_stepper_pt->is_steady())
154  {
155  // Number of timsteps (past & present)
156  const unsigned n_time = time_stepper_pt->ntstorage();
157  // Local storage (cache) for the weights
158  double weight[n_time];
159  for (unsigned t = 0; t < n_time; t++)
160  {
161  weight[t] = time_stepper_pt->weight(1, t);
162  }
163 
164  // Loop over the reagents
165  for (unsigned r = 0; r < NREAGENT; r++)
166  {
167  // Find the index at which the variable is stored
168  const unsigned c_nodal_index = c_index_adv_diff_react(r);
169 
170  for (unsigned t = 0; t < n_time; t++)
171  {
172  dc_dt[r] += weight[t] * nodal_value(t, n, c_nodal_index);
173  }
174  }
175  }
176  }

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::c_index_adv_diff_react(), oomph::TimeStepper::is_steady(), n, oomph::FiniteElement::nodal_value(), oomph::FiniteElement::node_pt(), oomph::TimeStepper::ntstorage(), UniformPSDSelfTest::r, plotPSD::t, oomph::GeomObject::time_stepper_pt(), oomph::Data::time_stepper_pt(), and oomph::TimeStepper::weight().

◆ diff()

template<unsigned NREAGENT, unsigned DIM>
const Vector<double>& oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::diff ( ) const
inline

Vector of diffusion coefficients.

316  {
317  return *Diff_pt;
318  }

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Diff_pt.

◆ diff_pt()

template<unsigned NREAGENT, unsigned DIM>
Vector<double>*& oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::diff_pt ( )
inline

◆ disable_ALE()

template<unsigned NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, 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.

181  {
182  ALE_is_disabled = true;
183  }

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::ALE_is_disabled.

Referenced by oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >::disable_ALE(), and oomph::RefineableDoubleBuoyantQCrouzeixRaviartElement< DIM >::disable_ALE().

◆ dshape_and_dtest_eulerian_adv_diff_react()

template<unsigned NREAGENT, unsigned DIM>
virtual double oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::dshape_and_dtest_eulerian_adv_diff_react ( 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::TAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, oomph::QAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, and oomph::QAdvectionDiffusionReactionElement< 2, DIM, 3 >.

◆ dshape_and_dtest_eulerian_at_knot_adv_diff_react()

template<unsigned NREAGENT, unsigned DIM>
virtual double oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::dshape_and_dtest_eulerian_at_knot_adv_diff_react ( 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::TAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, oomph::QAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, and oomph::QAdvectionDiffusionReactionElement< 2, DIM, 3 >.

◆ enable_ALE()

template<unsigned NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, 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.

190  {
191  ALE_is_disabled = false;
192  }

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::ALE_is_disabled.

Referenced by oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >::enable_ALE(), and oomph::RefineableDoubleBuoyantQCrouzeixRaviartElement< DIM >::enable_ALE().

◆ fill_in_contribution_to_jacobian()

template<unsigned NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, 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.

535  {
536  // Call the generic routine with the flag set to 1
538  residuals, jacobian, GeneralisedElement::Dummy_matrix, 1);
539  }
virtual void fill_in_generic_residual_contribution_adv_diff_react(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: advection_diffusion_reaction_elements.cc:56
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227

References oomph::GeneralisedElement::Dummy_matrix, and oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::fill_in_generic_residual_contribution_adv_diff_react().

Referenced by oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >::fill_in_contribution_to_jacobian(), and oomph::RefineableDoubleBuoyantQCrouzeixRaviartElement< DIM >::fill_in_contribution_to_jacobian().

◆ fill_in_contribution_to_jacobian_and_mass_matrix()

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

◆ fill_in_contribution_to_residuals()

◆ fill_in_generic_residual_contribution_adv_diff_react()

template<unsigned NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::fill_in_generic_residual_contribution_adv_diff_react ( 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 and/or the mass matrix

flag=2: compute Jacobian, residuals and mass matrix flag=1: compute Jacobian and residuals flag=0: compute only residual Vector

Pure version without hanging nodes

Reimplemented in oomph::RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM >.

61  {
62  // Find out how many nodes there are
63  const unsigned n_node = nnode();
64 
65  // Get the nodal index at which the unknown is stored
66  unsigned c_nodal_index[NREAGENT];
67  for (unsigned r = 0; r < NREAGENT; r++)
68  {
69  c_nodal_index[r] = c_index_adv_diff_react(r);
70  }
71 
72  // Set up memory for the shape and test functions
73  Shape psi(n_node), test(n_node);
74  DShape dpsidx(n_node, DIM), dtestdx(n_node, DIM);
75 
76  // Set the value of n_intpt
77  unsigned n_intpt = integral_pt()->nweight();
78 
79  // Set the Vector to hold local coordinates
80  Vector<double> s(DIM);
81 
82  // Get diffusion coefficients
83  Vector<double> D = diff();
84 
85  // Get the timescales
86  Vector<double> T = tau();
87 
88  // Integers used to store the local equation number and local unknown
89  // indices for the residuals and jacobians
90  int local_eqn = 0, local_unknown = 0;
91 
92  // Loop over the integration points
93  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
94  {
95  // Assign values of s
96  for (unsigned i = 0; i < DIM; i++) s[i] = integral_pt()->knot(ipt, i);
97 
98  // Get the integral weight
99  double w = integral_pt()->weight(ipt);
100 
101  // Call the derivatives of the shape and test functions
103  ipt, psi, dpsidx, test, dtestdx);
104 
105  // Premultiply the weights and the Jacobian
106  double W = w * J;
107 
108  // Calculate local values of the solution and its derivatives
109  // Allocate
110  Vector<double> interpolated_c(NREAGENT, 0.0);
111  Vector<double> dcdt(NREAGENT, 0.0);
112  Vector<double> interpolated_x(DIM, 0.0);
113  DenseMatrix<double> interpolated_dcdx(NREAGENT, DIM, 0.0);
114  Vector<double> mesh_velocity(DIM, 0.0);
115 
116 
117  // Calculate function value and derivatives:
118  // Loop over nodes
119  for (unsigned l = 0; l < n_node; l++)
120  {
121  // Loop over directions to calculate the position
122  for (unsigned j = 0; j < DIM; j++)
123  {
124  interpolated_x[j] += raw_nodal_position(l, j) * psi(l);
125  }
126 
127  // Loop over the unknown reagents
128  for (unsigned r = 0; r < NREAGENT; r++)
129  {
130  // Get the value at the node
131  const double c_value = raw_nodal_value(l, c_nodal_index[r]);
132 
133  // Calculate the interpolated value
134  interpolated_c[r] += c_value * psi(l);
135  dcdt[r] += dc_dt_adv_diff_react(l, r) * psi(l);
136 
137  // Loop over directions to calculate the derivatie
138  for (unsigned j = 0; j < DIM; j++)
139  {
140  interpolated_dcdx(r, j) += c_value * dpsidx(l, j);
141  }
142  }
143  }
144 
145  // Mesh velocity?
146  if (!ALE_is_disabled)
147  {
148  for (unsigned l = 0; l < n_node; l++)
149  {
150  for (unsigned j = 0; j < DIM; j++)
151  {
152  mesh_velocity[j] += raw_dnodal_position_dt(l, j) * psi(l);
153  }
154  }
155  }
156 
157 
158  // Get source function
159  Vector<double> source(NREAGENT);
161 
162 
163  // Get wind
164  Vector<double> wind(DIM);
166 
167  // Get reaction terms
168  Vector<double> R(NREAGENT);
169  get_reaction_adv_diff_react(ipt, interpolated_c, R);
170 
171  // If we are getting the jacobian the get the derivative terms
172  DenseMatrix<double> dRdC(NREAGENT);
173  if (flag)
174  {
175  get_reaction_deriv_adv_diff_react(ipt, interpolated_c, dRdC);
176  }
177 
178 
179  // Assemble residuals and Jacobian
180  //--------------------------------
181 
182  // Loop over the test functions
183  for (unsigned l = 0; l < n_node; l++)
184  {
185  // Loop over the reagents
186  for (unsigned r = 0; r < NREAGENT; r++)
187  {
188  // Set the local equation number
189  local_eqn = nodal_local_eqn(l, c_nodal_index[r]);
190 
191  /*IF it's not a boundary condition*/
192  if (local_eqn >= 0)
193  {
194  // Add body force/source/reaction term and time derivative
195  residuals[local_eqn] -=
196  (T[r] * dcdt[r] + source[r] + R[r]) * test(l) * W;
197 
198  // The Advection Diffusion bit itself
199  for (unsigned k = 0; k < DIM; k++)
200  {
201  // Terms that multiply the test function
202  double tmp = wind[k];
203  // If the mesh is moving need to subtract the mesh velocity
204  if (!ALE_is_disabled)
205  {
206  tmp -= T[r] * mesh_velocity[k];
207  }
208  // Now construct the contribution to the residuals
209  residuals[local_eqn] -= interpolated_dcdx(r, k) *
210  (tmp * test(l) + D[r] * dtestdx(l, k)) *
211  W;
212  }
213 
214  // Calculate the jacobian
215  //-----------------------
216  if (flag)
217  {
218  // Loop over the velocity shape functions again
219  for (unsigned l2 = 0; l2 < n_node; l2++)
220  {
221  // Loop over the reagents again
222  for (unsigned r2 = 0; r2 < NREAGENT; r2++)
223  {
224  // Set the number of the unknown
225  local_unknown = nodal_local_eqn(l2, c_nodal_index[r2]);
226 
227  // If at a non-zero degree of freedom add in the entry
228  if (local_unknown >= 0)
229  {
230  // Diagonal terms (i.e. the basic equations)
231  if (r2 == r)
232  {
233  // Mass matrix term
234  jacobian(local_eqn, local_unknown) -=
235  T[r] * test(l) * psi(l2) *
236  node_pt(l2)->time_stepper_pt()->weight(1, 0) * W;
237 
238  // Add the mass matrix term
239  if (flag == 2)
240  {
241  mass_matrix(local_eqn, local_unknown) +=
242  T[r] * test(l) * psi(l2) * W;
243  }
244 
245  // Add contribution to Elemental Matrix
246  for (unsigned i = 0; i < DIM; i++)
247  {
248  // Temporary term used in assembly
249  double tmp = wind[i];
250  if (!ALE_is_disabled) tmp -= T[r] * mesh_velocity[i];
251  // Now assemble Jacobian term
252  jacobian(local_eqn, local_unknown) -=
253  dpsidx(l2, i) *
254  (tmp * test(l) + D[r] * dtestdx(l, i)) * W;
255  }
256 
257  } // End of diagonal terms
258 
259  // Now add the cross-reaction terms
260  jacobian(local_eqn, local_unknown) -=
261  dRdC(r, r2) * psi(l2) * test(l) * W;
262  }
263  }
264  }
265  } // End of jacobian
266  }
267  } // End of loop over reagents
268  } // End of loop over nodes
269  } // End of loop over integration points
270  }
dominoes D
Definition: Domino.cpp:55
@ R
Definition: StatisticsVector.h:21
virtual void get_reaction_deriv_adv_diff_react(const unsigned &ipt, const Vector< double > &C, DenseMatrix< double > &dRdC) const
Definition: advection_diffusion_reaction_elements.h:417
virtual void get_source_adv_diff_react(const unsigned &ipt, const Vector< double > &x, Vector< double > &source) const
Definition: advection_diffusion_reaction_elements.h:342
virtual void get_wind_adv_diff_react(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Definition: advection_diffusion_reaction_elements.h:366
const Vector< double > & tau() const
Vector of dimensionless timescales.
Definition: advection_diffusion_reaction_elements.h:327
virtual double dshape_and_dtest_eulerian_at_knot_adv_diff_react(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
const Vector< double > & diff() const
Vector of diffusion coefficients.
Definition: advection_diffusion_reaction_elements.h:315
double dc_dt_adv_diff_react(const unsigned &n, const unsigned &r) const
Definition: advection_diffusion_reaction_elements.h:113
virtual void get_reaction_adv_diff_react(const unsigned &ipt, const Vector< double > &C, Vector< double > &R) const
Definition: advection_diffusion_reaction_elements.h:395
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
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, k, UniformPSDSelfTest::r, R, s, TestProblem::source(), Eigen::test, tmp, w, and oomph::QuadTreeNames::W.

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

◆ get_flux()

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

Get flux: \(\mbox{flux}[DIM r + i] = \mbox{d}C_{r} / \mbox{d}x_i \).

482  {
483  // Find out how many nodes there are in the element
484  const unsigned n_node = nnode();
485 
486  // Set up memory for the shape and test functions
487  Shape psi(n_node);
488  DShape dpsidx(n_node, DIM);
489 
490  // Call the derivatives of the shape and test functions
491  dshape_eulerian(s, psi, dpsidx);
492 
493  // Initialise to zero
494  for (unsigned j = 0; j < DIM * NREAGENT; j++)
495  {
496  flux[j] = 0.0;
497  }
498 
499  // Loop over the reagent terms
500  for (unsigned r = 0; r < NREAGENT; r++)
501  {
502  unsigned c_nodal_index = c_index_adv_diff_react(r);
503 
504  // Loop over derivative directions
505  for (unsigned j = 0; j < DIM; j++)
506  {
507  unsigned index = r * DIM + j;
508  // Loop over the nodes
509  for (unsigned l = 0; l < n_node; l++)
510  {
511  flux[index] += nodal_value(l, c_nodal_index) * dpsidx(l, j);
512  }
513  }
514  }
515  }
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3298
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::c_index_adv_diff_react(), DIM, oomph::FiniteElement::dshape_eulerian(), ProblemParameters::flux(), j, oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), UniformPSDSelfTest::r, and s.

Referenced by oomph::RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM >::get_Z2_flux(), and oomph::TAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >::get_Z2_flux().

◆ get_reaction_adv_diff_react()

template<unsigned NREAGENT, unsigned DIM>
virtual void oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::get_reaction_adv_diff_react ( const unsigned ipt,
const Vector< double > &  C,
Vector< double > &  R 
) const
inlinevirtual

Get reaction as a function of the given reagent concentrations This function is virtual to allow overloading in multi-physics problems where the reaction function might be determined by another system of equations

Reimplemented in oomph::RefineableDoubleBuoyantQCrouzeixRaviartElement< DIM >, and oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >.

398  {
399  // If no wind function has been set, return zero
400  if (Reaction_fct_pt == 0)
401  {
402  for (unsigned r = 0; r < NREAGENT; r++)
403  {
404  R[r] = 0.0;
405  }
406  }
407  else
408  {
409  // Get reaction terms
410  (*Reaction_fct_pt)(C, R);
411  }
412  }
Matrix< Scalar, Dynamic, Dynamic > C
Definition: bench_gemm.cpp:49

References UniformPSDSelfTest::r, R, and oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Reaction_fct_pt.

◆ get_reaction_deriv_adv_diff_react()

template<unsigned NREAGENT, unsigned DIM>
virtual void oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::get_reaction_deriv_adv_diff_react ( const unsigned ipt,
const Vector< double > &  C,
DenseMatrix< double > &  dRdC 
) const
inlinevirtual

Get the derivatives of the reaction terms with respect to the concentration variables. If no explicit function pointer is set, these will be calculated by finite differences

Reimplemented in oomph::RefineableDoubleBuoyantQCrouzeixRaviartElement< DIM >, and oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >.

421  {
422  // If no reaction pointer set, return zero
423  if (Reaction_fct_pt == 0)
424  {
425  for (unsigned r = 0; r < NREAGENT; r++)
426  {
427  for (unsigned p = 0; p < NREAGENT; p++)
428  {
429  dRdC(r, p) = 0.0;
430  }
431  }
432  }
433  else
434  {
435  // If no function pointer get finite differences
436  if (Reaction_deriv_fct_pt == 0)
437  {
438  // Local copy of the unknowns
439  Vector<double> C_local = C;
440  // Finite differences
441  Vector<double> R(NREAGENT), R_plus(NREAGENT), R_minus(NREAGENT);
442  // Get the initial reaction terms
443  //(*Reaction_fct_pt)(C,R);
444  const double fd_step = GeneralisedElement::Default_fd_jacobian_step;
445  // Now loop over all the reagents
446  for (unsigned p = 0; p < NREAGENT; p++)
447  {
448  // Store the old value
449  double old_var = C_local[p];
450  // Increment the value
451  C_local[p] += fd_step;
452  // Get the new values
453  (*Reaction_fct_pt)(C_local, R_plus);
454 
455  // Reset the value
456  C_local[p] = old_var;
457  // Decrement the value
458  C_local[p] -= fd_step;
459  // Get the new values
460  (*Reaction_fct_pt)(C_local, R_minus);
461 
462  // Assemble the column of the jacobian
463  for (unsigned r = 0; r < NREAGENT; r++)
464  {
465  dRdC(r, p) = (R_plus[r] - R_minus[r]) / (2.0 * fd_step);
466  }
467 
468  // Reset the value
469  C_local[p] = old_var;
470  }
471  }
472  // Otherwise get the terms from the function
473  else
474  {
475  (*Reaction_deriv_fct_pt)(C, dRdC);
476  }
477  }
478  }
float * p
Definition: Tutorial_Map_using.cpp:9
static double Default_fd_jacobian_step
Definition: elements.h:1198

References oomph::GeneralisedElement::Default_fd_jacobian_step, p, UniformPSDSelfTest::r, R, oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Reaction_deriv_fct_pt, and oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Reaction_fct_pt.

◆ get_source_adv_diff_react()

template<unsigned NREAGENT, unsigned DIM>
virtual void oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::get_source_adv_diff_react ( const unsigned ipt,
const Vector< double > &  x,
Vector< 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

345  {
346  // If no source function has been set, return zero
347  if (Source_fct_pt == 0)
348  {
349  for (unsigned r = 0; r < NREAGENT; r++)
350  {
351  source[r] = 0.0;
352  }
353  }
354  else
355  {
356  // Get source strength
357  (*Source_fct_pt)(x, source);
358  }
359  }

References UniformPSDSelfTest::r, TestProblem::source(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Source_fct_pt, and plotDoE::x.

◆ get_wind_adv_diff_react()

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

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

Reimplemented in oomph::RefineableDoubleBuoyantQCrouzeixRaviartElement< DIM >, oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >, oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >, and oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >.

370  {
371  // If no wind function has been set, return zero
372  if (Wind_fct_pt == 0)
373  {
374  for (unsigned i = 0; i < DIM; i++)
375  {
376  wind[i] = 0.0;
377  }
378  }
379  else
380  {
381  // Get continuous time from timestepper of first node
382  double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
383 
384  // Get wind
385  (*Wind_fct_pt)(time, x, wind);
386  }
387  }
Time *const & time_pt() const
Access function for the pointer to time (const version)
Definition: timesteppers.h:572
double & time()
Return the current value of the continuous time.
Definition: timesteppers.h:123

References DIM, i, oomph::FiniteElement::node_pt(), oomph::Time::time(), oomph::TimeStepper::time_pt(), oomph::Data::time_stepper_pt(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Wind_fct_pt, and plotDoE::x.

◆ integrate_reagents()

template<unsigned NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::integrate_reagents ( Vector< double > &  C) const

Return the integrated reagent concentrations.

Integrate the reagent concentrations over the element.

359  {
360  // Find out how many nodes there are
361  const unsigned n_node = nnode();
362 
363  // Get the nodal index at which the unknown is stored
364  unsigned c_nodal_index[NREAGENT];
365  for (unsigned r = 0; r < NREAGENT; r++)
366  {
367  c_nodal_index[r] = c_index_adv_diff_react(r);
368  }
369 
370  // Set up memory for the shape and test functions
371  Shape psi(n_node);
372  DShape dpsidx(n_node, DIM);
373 
374  // Set the value of n_intpt
375  unsigned n_intpt = integral_pt()->nweight();
376 
377  // Loop over the integration points
378  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
379  {
380  // Get the integral weight
381  double w = integral_pt()->weight(ipt);
382 
383  // Call the derivatives of the shape and test functions
384  double J = dshape_eulerian_at_knot(ipt, psi, dpsidx);
385 
386  // Premultiply the weights and the Jacobian
387  double W = w * J;
388 
389  // Calculate local values of the solution and its derivatives
390  // Allocate
391  Vector<double> interpolated_c(NREAGENT, 0.0);
392 
393  // Calculate function value and derivatives:
394  // Loop over nodes
395  for (unsigned l = 0; l < n_node; l++)
396  {
397  // Loop over the unknown reagents
398  for (unsigned r = 0; r < NREAGENT; r++)
399  {
400  // Get the value at the node
401  const double c_value = raw_nodal_value(l, c_nodal_index[r]);
402 
403  // Calculate the interpolated value
404  interpolated_c[r] += c_value * psi(l);
405  }
406  }
407 
408  for (unsigned r = 0; r < NREAGENT; r++)
409  {
410  C[r] += interpolated_c[r] * W;
411  }
412  } // End of loop over integration points
413  }
virtual double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3325
Definition: matrices.h:74

References DIM, J, UniformPSDSelfTest::r, w, and oomph::QuadTreeNames::W.

◆ interpolated_c_adv_diff_react()

template<unsigned NREAGENT, unsigned DIM>
double oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::interpolated_c_adv_diff_react ( const Vector< double > &  s,
const unsigned i 
) const
inline

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

558  {
559  // Find number of nodes
560  unsigned n_node = nnode();
561 
562  // Get the nodal index at which the unknown is stored
563  unsigned c_nodal_index = c_index_adv_diff_react(i);
564 
565  // Local shape function
566  Shape psi(n_node);
567 
568  // Find values of shape function
569  shape(s, psi);
570 
571  // Initialise value of u
572  double interpolated_c = 0.0;
573 
574  // Loop over the local nodes and sum
575  for (unsigned l = 0; l < n_node; l++)
576  {
577  interpolated_c += nodal_value(l, c_nodal_index) * psi[l];
578  }
579 
580  return (interpolated_c);
581  }
virtual void shape(const Vector< double > &s, Shape &psi) const =0

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::c_index_adv_diff_react(), i, oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_value(), s, and oomph::FiniteElement::shape().

Referenced by oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >::output(), and oomph::RefineableDoubleBuoyantQCrouzeixRaviartElement< DIM >::output().

◆ operator=()

template<unsigned NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::operator= ( const AdvectionDiffusionReactionEquations< NREAGENT, DIM > &  )
delete

Broken assignment operator.

◆ output() [1/4]

template<unsigned NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::output ( FILE *  file_pt)
inlinevirtual

C_style output with default number of plot points.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::TAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, oomph::QAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, and oomph::QAdvectionDiffusionReactionElement< 2, DIM, 3 >.

209  {
210  unsigned n_plot = 5;
211  output(file_pt, n_plot);
212  }
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: advection_diffusion_reaction_elements.h:196

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::output().

◆ output() [2/4]

template<unsigned NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, 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::TAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, oomph::QAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, and oomph::QAdvectionDiffusionReactionElement< 2, DIM, 3 >.

481  {
482  // Vector of local coordinates
483  Vector<double> s(DIM);
484 
485  // Tecplot header info
486  fprintf(file_pt, "%s", tecplot_zone_string(nplot).c_str());
487 
488  // Loop over plot points
489  unsigned num_plot_points = nplot_points(nplot);
490  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
491  {
492  // Get local coordinates of plot point
493  get_s_plot(iplot, nplot, s);
494 
495  for (unsigned i = 0; i < DIM; i++)
496  {
497  fprintf(file_pt, "%g ", interpolated_x(s, i));
498  }
499  for (unsigned i = 0; i < NREAGENT; i++)
500  {
501  fprintf(file_pt, "%g \n", interpolated_c_adv_diff_react(s, i));
502  }
503  }
504 
505  // Write tecplot footer (e.g. FE connectivity lists)
506  write_tecplot_zone_footer(file_pt, nplot);
507  }
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]

◆ output() [4/4]

template<unsigned NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, 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::TAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, oomph::QAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, and oomph::QAdvectionDiffusionReactionElement< 2, DIM, 3 >.

426  {
427  // Vector of local coordinates
428  Vector<double> s(DIM);
429 
430 
431  // Tecplot header info
432  outfile << tecplot_zone_string(nplot);
433 
434  // Loop over plot points
435  unsigned num_plot_points = nplot_points(nplot);
436  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
437  {
438  // Get local coordinates of plot point
439  get_s_plot(iplot, nplot, s);
440 
441  // Get Eulerian coordinate of plot point
442  Vector<double> x(DIM);
443  interpolated_x(s, x);
444 
445  for (unsigned i = 0; i < DIM; i++)
446  {
447  outfile << x[i] << " ";
448  }
449  for (unsigned i = 0; i < NREAGENT; i++)
450  {
451  outfile << interpolated_c_adv_diff_react(s, i) << " ";
452  }
453 
454  // Get the wind
455  Vector<double> wind(DIM);
456  // Dummy integration point needed
457  unsigned ipt = 0;
458  get_wind_adv_diff_react(ipt, s, x, wind);
459  for (unsigned i = 0; i < DIM; i++)
460  {
461  outfile << wind[i] << " ";
462  }
463  outfile << std::endl;
464  }
465 
466  // Write tecplot footer (e.g. FE connectivity lists)
467  write_tecplot_zone_footer(outfile, nplot);
468  }

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

◆ output_fct() [1/2]

template<unsigned NREAGENT, unsigned DIM>
virtual void oomph::AdvectionDiffusionReactionEquations< NREAGENT, 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::RefineableDoubleBuoyantQCrouzeixRaviartElement< DIM >, oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >, oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >, oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >, oomph::TAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, oomph::QAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, and oomph::QAdvectionDiffusionReactionElement< 2, DIM, 3 >.

232  {
233  throw OomphLibError("There is no time-dependent output_fct() for "
234  "Advection Diffusion elements",
237  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ output_fct() [2/2]

template<unsigned NREAGENT, unsigned DIM>
void oomph::AdvectionDiffusionReactionEquations< NREAGENT, 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::TAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, oomph::QAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >, and oomph::QAdvectionDiffusionReactionElement< 2, DIM, 3 >.

523  {
524  // Vector of local coordinates
525  Vector<double> s(DIM);
526 
527  // Vector for coordintes
528  Vector<double> x(DIM);
529 
530  // Tecplot header info
531  outfile << tecplot_zone_string(nplot);
532 
533  // Exact solution Vector (here a scalar)
534  Vector<double> exact_soln(1);
535 
536  // Loop over plot points
537  unsigned num_plot_points = nplot_points(nplot);
538  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
539  {
540  // Get local coordinates of plot point
541  get_s_plot(iplot, nplot, s);
542 
543  // Get x position as Vector
544  interpolated_x(s, x);
545 
546  // Get exact solution at this point
547  (*exact_soln_pt)(x, exact_soln);
548 
549  // Output x,y,...,u_exact
550  for (unsigned i = 0; i < DIM; i++)
551  {
552  outfile << x[i] << " ";
553  }
554  outfile << exact_soln[0] << std::endl;
555  }
556 
557  // Write tecplot footer (e.g. FE connectivity lists)
558  write_tecplot_zone_footer(outfile, nplot);
559  }

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

Referenced by oomph::QAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >::output_fct(), and oomph::TAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D >::output_fct().

◆ reaction_deriv_fct_pt() [1/2]

template<unsigned NREAGENT, unsigned DIM>
AdvectionDiffusionReactionReactionDerivFctPt& oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::reaction_deriv_fct_pt ( )
inline

◆ reaction_deriv_fct_pt() [2/2]

template<unsigned NREAGENT, unsigned DIM>
AdvectionDiffusionReactionReactionDerivFctPt oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::reaction_deriv_fct_pt ( ) const
inline

Access function: Pointer to reaction function. Const version.

310  {
311  return Reaction_deriv_fct_pt;
312  }

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Reaction_deriv_fct_pt.

◆ reaction_fct_pt() [1/2]

template<unsigned NREAGENT, unsigned DIM>
AdvectionDiffusionReactionReactionFctPt& oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::reaction_fct_pt ( )
inline

◆ reaction_fct_pt() [2/2]

template<unsigned NREAGENT, unsigned DIM>
AdvectionDiffusionReactionReactionFctPt oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::reaction_fct_pt ( ) const
inline

Access function: Pointer to reaction function. Const version.

298  {
299  return Reaction_fct_pt;
300  }

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Reaction_fct_pt.

◆ self_test()

template<unsigned NREAGENT, unsigned DIM>
unsigned oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::self_test
virtual

Self-test: Return 0 for OK.

Reimplemented from oomph::FiniteElement.

332  {
333  bool passed = true;
334 
335  // Check lower-level stuff
336  if (FiniteElement::self_test() != 0)
337  {
338  passed = false;
339  }
340 
341  // Return verdict
342  if (passed)
343  {
344  return 0;
345  }
346  else
347  {
348  return 1;
349  }
350  }
virtual unsigned self_test()
Definition: elements.cc:4440

References oomph::FiniteElement::self_test().

◆ source_fct_pt() [1/2]

template<unsigned NREAGENT, unsigned DIM>
AdvectionDiffusionReactionSourceFctPt& oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::source_fct_pt ( )
inline

◆ source_fct_pt() [2/2]

template<unsigned NREAGENT, unsigned DIM>
AdvectionDiffusionReactionSourceFctPt oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::source_fct_pt ( ) const
inline

Access function: Pointer to source function. Const version.

273  {
274  return Source_fct_pt;
275  }

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Source_fct_pt.

◆ tau()

template<unsigned NREAGENT, unsigned DIM>
const Vector<double>& oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::tau ( ) const
inline

Vector of dimensionless timescales.

328  {
329  return *Tau_pt;
330  }

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Tau_pt.

◆ tau_pt()

template<unsigned NREAGENT, unsigned DIM>
Vector<double>*& oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::tau_pt ( )
inline

Pointer to vector of dimensionless timescales.

334  {
335  return Tau_pt;
336  }

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Tau_pt.

Referenced by oomph::RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM >::further_build().

◆ wind_fct_pt() [1/2]

template<unsigned NREAGENT, unsigned DIM>
AdvectionDiffusionReactionWindFctPt& oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::wind_fct_pt ( )
inline

◆ wind_fct_pt() [2/2]

template<unsigned NREAGENT, unsigned DIM>
AdvectionDiffusionReactionWindFctPt oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::wind_fct_pt ( ) const
inline

Access function: Pointer to reaction function. Const version.

286  {
287  return Wind_fct_pt;
288  }

References oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Wind_fct_pt.

Member Data Documentation

◆ ALE_is_disabled

template<unsigned NREAGENT, unsigned DIM>
bool oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::ALE_is_disabled
protected

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

Referenced by oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::disable_ALE(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::enable_ALE(), and oomph::RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM >::further_build().

◆ Default_dimensionless_number

template<unsigned NREAGENT, unsigned DIM>
Vector< double > oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::Default_dimensionless_number
staticprivate

Static default value for the dimensionless numbers.

2D Advection Diffusion elements

Default value for Peclet number

Referenced by oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::AdvectionDiffusionReactionEquations().

◆ Diff_pt

◆ N_reagent

template<unsigned NREAGENT, unsigned DIM>
const unsigned oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::N_reagent
staticprotected
Initial value:
=
NREAGENT

◆ Reaction_deriv_fct_pt

◆ Reaction_fct_pt

◆ Source_fct_pt

◆ Tau_pt

◆ Wind_fct_pt


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