30#ifndef OOMPH_BOUSSINESQ_ELEMENTS_HEADER
31#define OOMPH_BOUSSINESQ_ELEMENTS_HEADER
36#include "advection_diffusion.h"
37#include "navier_stokes.h"
64 template<
unsigned DIM>
102 const double&
ra()
const
135 "This function hasn't been implemented for this element",
148 const unsigned&
nplot)
const
151 "This function hasn't been implemented for this element",
191 for (
unsigned i = 0;
i <
DIM;
i++)
197 for (
unsigned i = 0;
i <
DIM;
i++)
229 const unsigned&
Nplot,
239 const unsigned&
Nplot,
312 for (
unsigned i = 0;
i <
DIM;
i++)
336#ifdef USE_FD_JACOBIAN_FOR_BUOYANT_Q_ELEMENT
365#ifdef USE_OFF_DIAGONAL_FD_JACOBIAN_FOR_BUOYANT_Q_ELEMENT
375 for (
unsigned i = 0;
i <
DIM;
i++)
405 for (
unsigned i = 0;
i <
DIM;
i++)
426 for (
unsigned m = 0;
m <
n_dof;
m++)
479 for (
unsigned m = 0;
m <
n_dof;
m++)
492 for (
unsigned j = 0;
j <
DIM;
j++)
575 for (
unsigned i = 0;
i <
DIM;
i++)
594 double Ra = this->
ra();
595 double Pe = this->
pe();
624 for (
unsigned j = 0;
j <
DIM;
j++)
640 for (
unsigned i = 0;
i <
DIM;
i++)
754 template<
unsigned int DIM>
756 :
public virtual QElement<DIM - 1, 3>
799 template<
unsigned DIM>
832 const double&
ra()
const
867 "This function hasn't been implemented for this element",
880 const unsigned&
nplot)
const
883 "This function hasn't been implemented for this element",
921 for (
unsigned i = 0;
i <
DIM;
i++)
927 for (
unsigned i = 0;
i <
DIM;
i++)
958 const unsigned&
Nplot,
968 const unsigned&
Nplot,
1025 for (
unsigned i = 0;
i <
DIM;
i++)
1058 for (
unsigned i = 0;
i <
DIM;
i++)
1232 const unsigned&
ipt,
1245 for (
unsigned i = 0;
i <
DIM;
i++)
1269#ifdef USE_FD_JACOBIAN_FOR_REFINEABLE_BUOYANT_Q_ELEMENT
1314 for (
unsigned i = 0;
i <
DIM;
i++)
1333 double Ra = this->
ra();
1334 double Pe = this->
pe();
1363 for (
unsigned j = 0;
j <
DIM;
j++)
1416 for (
unsigned i = 0;
i <
DIM;
i++)
1588 double RefineableBuoyantQCrouzeixRaviartElement<
1589 2>::Default_Physical_Constant_Value = 0.0;
1593 3>::Default_Physical_Constant_Value = 0.0;
A class for all elements that solve the Advection Diffusion equations using isoparametric elements.
double interpolated_u_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the element's contribution to its residual vector and the element Jacobian matrix (wrapper)
void enable_ALE()
(Re-)enable ALE, i.e. take possible mesh motion into account when evaluating the time-derivative....
const double & pe() const
Peclet number.
void disable_ALE()
Disable ALE, i.e. assert the mesh is not moving – you do this at your own risk!
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
A class that solves the Boussinesq approximation of the Navier–Stokes and energy equations by couplin...
double * Ra_pt
Pointer to a private data member, the Rayleigh number.
double *& ra_pt()
Access function for the pointer to the Rayleigh number.
std::string scalar_name_paraview(const unsigned &i) const
Name of the i-th scalar field. Default implementation returns V1 for the first one,...
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Validate against exact solution. Solution is provided via function pointer. Plot at a given number of...
static double Default_Physical_Constant_Value
The static default value of the Rayleigh number.
void scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Write values of the i-th scalar field at the plot points. Broken virtual. Needs to be implemented for...
void output(std::ostream &outfile)
Overload the standard output function with the broken default.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: Broken default.
unsigned nscalar_paraview() const
Number of scalars/fields output by this element. Broken virtual. Needs to be implemented for each new...
void disable_ALE()
Final override for disable ALE.
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Validate against exact solution at given time Solution is provided via function pointer....
void fill_in_off_diagonal_jacobian_blocks_by_fd(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Helper function to get the off-diagonal blocks of the Jacobian matrix by finite differences.
unsigned required_nvalue(const unsigned &n) const
The required number of values stored at the nodes is the sum of the required values of the two single...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual vector and the Jacobian matrix. Jacobian is computed by finite-differe...
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Add the element's contribution to its residuals vector, jacobian matrix and mass matrix.
void output(std::ostream &outfile, const unsigned &nplot)
Output function: Output x, y, u, v, p, theta at Nplot^DIM plot points.
void unfix_pressure(const unsigned &p_dof)
Unpin p_dof-th pressure dof.
BuoyantQCrouzeixRaviartElement()
Constructor: call the underlying constructors and initialise the pointer to the Rayleigh number to po...
void output(FILE *file_pt)
C-style output function: Broken default.
void fill_in_off_diagonal_jacobian_blocks_analytic(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Helper function to get the off-diagonal blocks of the Jacobian matrix analytically.
const double & ra() const
Access function for the Rayleigh number (const version)
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Calculate the element's contribution to the residual vector. Recall that fill_in_* functions MUST NOT...
unsigned u_index_adv_diff() const
Overload the index at which the temperature variable is stored. We choose to store it after the fluid...
void output_fct(std::ostream &outfile, const unsigned &Nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: Broken default.
void get_wind_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Overload the wind function in the advection-diffusion equations. This provides the coupling from the ...
void output_fct(std::ostream &outfile, const unsigned &Nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output function for a time-dependent exact solution: Broken default.
void get_body_force_nst(const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &result)
Overload the body force in the Navier-Stokes equations This provides the coupling from the advection-...
void enable_ALE()
Final override for enable ALE.
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
void unpin(const unsigned &i)
Unpin the i-th stored variable.
double * value_pt(const unsigned &i) const
Return the pointer to the i-the stored value. Typically this is required when direct access to the st...
virtual void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)=0
Z2 'flux' terms for Z2 error estimation.
virtual unsigned num_Z2_flux_terms()=0
Number of 'flux' terms for Z2 error estimation.
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the elemental contribution to the jacobian matrix. and the residuals vector. Note that this funct...
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output an exact solution over the element.
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned nnode() const
Return the number of nodes.
virtual void disable_ALE()
This is an empty function that establishes a uniform interface for all (derived) elements that involv...
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.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
virtual void enable_ALE()
(Re-)enable ALE, i.e. take possible mesh motion into account when evaluating the time-derivative....
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n but do NOT take hanging nodes into account.
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Add tecplot zone "footer" to output stream (when plotting nplot points in each "coordinate direction"...
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
static double Default_fd_jacobian_step
Double used for the default finite difference step in elemental jacobian calculations.
virtual void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Add the elemental contribution to the jacobian matrix, mass matrix and the residuals vector....
virtual void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the elemental contribution to the residuals vector. Note that this function will NOT initialise t...
unsigned ndof() const
Return the number of equations/dofs in the element.
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
virtual void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the elemental contribution to the jacobian matrix. and the residuals vector. Note that this funct...
Class that contains data for hanging nodes.
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
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.
A class for elements that solve the cartesian Navier–Stokes equations, templated by the dimension DIM...
void disable_ALE()
Disable ALE, i.e. assert the mesh is not moving – you do this at your own risk!
virtual unsigned u_index_nst(const unsigned &i) const
Return the index at which the i-th unknown velocity component is stored. The default value,...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute the element's residual Vector.
void interpolated_u_nst(const Vector< double > &s, Vector< double > &veloc) const
Compute vector of FE interpolated velocity u at local coordinate s.
void enable_ALE()
(Re-)enable ALE, i.e. take possible mesh motion into account when evaluating the time-derivative....
virtual double interpolated_p_nst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual Vector and the jacobian matrix Virtual function can be overloaded by h...
const Vector< double > & g() const
Vector of gravitational components.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
bool is_hanging() const
Test whether the node is geometrically hanging.
HangInfo *const & hanging_pt() const
Return pointer to hanging node data (this refers to the geometric hanging node status) (const version...
An OomphLibError object which should be thrown when an run-time error is encountered....
Point element has just a single node and a single shape function which is identically equal to one.
QAdvectionDiffusionElement elements are linear/quadrilateral/brick-shaped Advection Diffusion element...
Crouzeix_Raviart elements are Navier–Stokes elements with quadratic interpolation for velocities and ...
double dshape_and_dtest_eulerian_at_knot_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Velocity shape and test functions and their derivs w.r.t. to global coords at ipt-th integation point...
unsigned P_nst_internal_index
Internal index that indicates at which internal data the pressure is stored.
A version of the Advection Diffusion equations that can be used with non-uniform mesh refinement....
void further_build()
Further build: Copy source function pointer from father element.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from AdvectionDiffusion equations.
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Get the function value u in Vector. Note: Given the generality of the interface (this function is usu...
A RefineableElement class that solves the Boussinesq approximation of the Navier–Stokes and energy eq...
void rebuild_from_sons(Mesh *&mesh_pt)
Call the rebuild_from_sons functions for each of the constituent multi-physics elements.
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Add the element's contribution to its residuals vector, jacobian matrix and mass matrix.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Validate against exact solution. Solution is provided via function pointer. Plot at a given number of...
void get_wind_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Overload the wind function in the advection-diffusion equations. This provides the coupling from the ...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual Vector and the jacobian matrix using full finite differences,...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the constituent elements' contribution to the residual vector.
RefineableBuoyantQCrouzeixRaviartElement()
Constructor: call the underlying constructors and initialise the pointer to the Rayleigh number to ad...
unsigned nrecovery_order()
The recovery order is that of the NavierStokes elements.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get the Z2 flux by concatenating the fluxes from the fluid and the advection diffusion elements.
unsigned nscalar_paraview() const
Number of scalars/fields output by this element. Broken virtual. Needs to be implemented for each new...
static double Default_Physical_Constant_Value
The static default value of the Rayleigh number.
void output(std::ostream &outfile, const unsigned &nplot)
Output function: x,y,u or x,y,z,u at Nplot^DIM plot points.
unsigned ncompound_fluxes()
The number of compound fluxes is two (one for the fluid and one for the temperature)
const double & ra() const
Access function for the Rayleigh number (const version)
void further_build()
Call the underlying single-physics element's further_build() functions and make sure that the pointer...
unsigned num_Z2_flux_terms()
The number of Z2 flux terms is the same as that in the fluid element plus that in the advection-diffu...
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: Broken default.
double * Ra_pt
Pointer to a new physical variable, the Rayleigh number.
unsigned required_nvalue(const unsigned &n) const
The required number of values stored at the nodes is the sum of the required values of the two single...
void get_body_force_nst(const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &result)
Overload the body force in the navier-stokes equations This provides the coupling from the advection-...
void fill_in_off_diagonal_jacobian_blocks_analytic(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the contribution of the off-diagonal blocks analytically.
void enable_ALE()
Final override for enable ALE.
void disable_ALE()
Final override for disable ALE.
unsigned nvertex_node() const
Number of vertex nodes in the element is obtained from the geometric element.
void output(FILE *file_pt)
C-style output function: Broken default.
unsigned ncont_interpolated_values() const
The total number of continously interpolated values is DIM+1 (DIM fluid velocities and one temperatur...
void scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Write values of the i-th scalar field at the plot points. Broken virtual. Needs to be implemented for...
void output_fct(std::ostream &outfile, const unsigned &Nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: Broken default.
void output_fct(std::ostream &outfile, const unsigned &Nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output function for a time-dependent exact solution. Broken default.
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Get the continuously interpolated values at the local coordinate s. We choose to put the fluid veloci...
void further_setup_hanging_nodes()
The additional hanging node information must be set up for both single-physics elements.
double *& ra_pt()
Access function for the pointer to the Rayleigh number.
void output(std::ostream &outfile)
Overload the standard output function with the broken default.
unsigned u_index_adv_diff() const
Overload the index at which the temperature variable is stored. We choose to store is after the fluid...
void get_Z2_compound_flux_indices(Vector< unsigned > &flux_index)
Fill in which flux components are associated with the fluid measure and which are associated with the...
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Validate against exact solution at given time Solution is provided via function pointer....
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element, Call the geometric element's function.
std::string scalar_name_paraview(const unsigned &i) const
Name of the i-th scalar field. Default implementation returns V1 for the first one,...
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Get all continuously interpolated values at the local coordinate s at time level t (t=0: present; t>0...
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
int local_hang_eqn(Node *const &node_pt, const unsigned &i)
Access function that returns the local equation number for the hanging node variables (values stored ...
Refineable version of QAdvectionDiffusionElement. Inherit from the standard QAdvectionDiffusionElemen...
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Refineable version of Crouzeix Raviart elements. Generic class definitions.
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Get the function value u in Vector. Note: Given the generality of the interface (this function is usu...
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
void further_build()
Further build for Crouzeix_Raviart interpolates the internal pressure dofs from father element: Make ...
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: Reconstruct pressure from the (merged) sons This must be specialised for each dime...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
std::string to_string(T object, unsigned float_precision=8)
Conversion function that should work for anything with operator<< defined (at least all basic types).
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).