27 #ifndef OOMPH_GEN_ADV_DIFF_ELEMENTS_HEADER
28 #define OOMPH_GEN_ADV_DIFF_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
37 #include "../generic/nodes.h"
38 #include "../generic/Qelements.h"
39 #include "../generic/oomph_utilities.h"
50 template<
unsigned DIM>
121 for (
unsigned t = 0;
t < n_time;
t++)
157 void output(std::ostream& outfile,
const unsigned& nplot);
169 void output(FILE* file_pt,
const unsigned& n_plot);
174 const unsigned& nplot,
181 std::ostream& outfile,
182 const unsigned& nplot,
186 throw OomphLibError(
"There is no time-dependent output_fct() for "
187 "Advection Diffusion elements",
188 OOMPH_CURRENT_FUNCTION,
189 OOMPH_EXCEPTION_LOCATION);
208 "No time-dependent compute_error() for Advection Diffusion elements",
209 OOMPH_CURRENT_FUNCTION,
210 OOMPH_EXCEPTION_LOCATION);
273 const double&
pe()
const
302 double& source)
const
312 (*Source_fct_pt)(x, source);
329 for (
unsigned i = 0;
i < DIM;
i++)
337 (*Wind_fct_pt)(x, wind);
357 for (
unsigned i = 0;
i < DIM;
i++)
365 (*Conserved_wind_fct_pt)(x, wind);
384 for (
unsigned i = 0;
i < DIM;
i++)
386 for (
unsigned j = 0; j < DIM; j++)
402 (*Diff_fct_pt)(x,
D);
411 unsigned n_node =
nnode();
418 DShape dpsidx(n_node, DIM);
424 for (
unsigned j = 0; j < DIM; j++)
430 for (
unsigned l = 0; l < n_node; l++)
433 for (
unsigned j = 0; j < DIM; j++)
435 flux[j] +=
nodal_value(l, u_nodal_index) * dpsidx(l, j);
445 const unsigned n_node =
nnode();
452 DShape dpsidx(n_node, DIM);
460 double interpolated_u = 0.0;
465 for (
unsigned l = 0; l < n_node; l++)
468 const double u_value = this->
nodal_value(l, u_nodal_index);
469 interpolated_u += u_value * psi(l);
471 for (
unsigned j = 0; j < DIM; j++)
474 interpolated_dudx[j] += u_value * dpsidx(l, j);
491 for (
unsigned i = 0;
i < DIM;
i++)
494 for (
unsigned j = 0; j < DIM; j++)
496 total_flux[
i] +=
D(
i, j) * interpolated_dudx[j];
498 total_flux[
i] -= conserved_wind[
i] * interpolated_u;
536 residuals, jacobian, mass_matrix, 2);
544 unsigned n_node =
nnode();
556 double interpolated_u = 0.0;
559 for (
unsigned l = 0; l < n_node; l++)
561 interpolated_u +=
nodal_value(l, u_nodal_index) * psi[l];
564 return (interpolated_u);
579 DShape& dtestdx)
const = 0;
588 DShape& dtestdx)
const = 0;
637 template<
unsigned DIM,
unsigned NNODE_1D>
639 :
public virtual QElement<DIM, NNODE_1D>,
680 void output(std::ostream& outfile,
const unsigned& n_plot)
695 void output(FILE* file_pt,
const unsigned& n_plot)
703 const unsigned& n_plot,
707 outfile, n_plot, exact_soln_pt);
715 const unsigned& n_plot,
720 outfile, n_plot, time, exact_soln_pt);
753 template<
unsigned DIM,
unsigned NNODE_1D>
762 double J = this->dshape_eulerian(
s, psi, dpsidx);
766 for (
unsigned i = 0;
i < NNODE_1D;
i++)
769 for (
unsigned j = 0; j < DIM; j++)
771 dtestdx(
i, j) = dpsidx(
i, j);
786 template<
unsigned DIM,
unsigned NNODE_1D>
795 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
817 template<
unsigned DIM,
unsigned NNODE_1D>
819 :
public virtual QElement<DIM - 1, NNODE_1D>
837 template<
unsigned NNODE_1D>
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
A general Finite Element class.
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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 shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
unsigned nnode() const
Return the number of nodes.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
double nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. If the node is hanging, the appropriate interpolation is ...
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
A class for all elements that solve the Advection Diffusion equations in conservative form using isop...
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)
GeneralisedAdvectionDiffusionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
double *& pe_pt()
Pointer to Peclet number.
virtual double dshape_and_dtest_eulerian_at_knot_cons_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
virtual void get_wind_cons_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Get wind at (Eulerian) position x and/or local coordinate s. This function is virtual to allow overlo...
double * Pe_pt
Pointer to global Peclet number.
virtual void get_diff_cons_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, DenseMatrix< double > &D) const
Get diffusivity tensor at (Eulerian) position x and/or local coordinate s. This function is virtual t...
void disable_ALE()
Disable ALE, i.e. assert the mesh is not moving – you do this at your own risk!
void output(FILE *file_pt)
C_style output with default number of plot points.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
static double Default_peclet_number
Static default value for the Peclet number.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
virtual void get_source_cons_adv_diff(const unsigned &ipt, const Vector< double > &x, double &source) const
Get source term at (Eulerian) position x. This function is virtual to allow overloading in multi-phys...
GeneralisedAdvectionDiffusionEquations()
Constructor: Initialise the Source_fct_pt and Wind_fct_pt to null and set (pointer to) Peclet number ...
virtual void output_fct(std::ostream &outfile, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: x,y,u_exact or x,y,z,u_exact at nplot^DIM plot points (dummy time-dependent versio...
GeneralisedAdvectionDiffusionDiffFctPt Diff_fct_pt
Pointer to diffusivity funciton.
double *& pe_st_pt()
Pointer to Peclet number multipled by Strouha number.
GeneralisedAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
GeneralisedAdvectionDiffusionWindFctPt wind_fct_pt() const
Access function: Pointer to wind function. Const version.
GeneralisedAdvectionDiffusionWindFctPt Conserved_wind_fct_pt
Pointer to additional (conservative) wind function:
void(* GeneralisedAdvectionDiffusionDiffFctPt)(const Vector< double > &x, DenseMatrix< double > &D)
Funciton pointer to a diffusivity function.
const double & pe() const
Peclet number.
void output(std::ostream &outfile)
Output with default number of plot points.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
double du_dt_cons_adv_diff(const unsigned &n) const
du/dt at local node n. Uses suitably interpolated value for hanging nodes.
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.
unsigned self_test()
Self-test: Return 0 for OK.
void get_total_flux(const Vector< double > &s, Vector< double > &total_flux) const
Get flux: .
GeneralisedAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
virtual double dshape_and_dtest_eulerian_cons_adv_diff(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at local coord. s; return Jacobian of mapping...
GeneralisedAdvectionDiffusionDiffFctPt diff_fct_pt() const
Access function: Pointer to diffusion function. Const version.
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
GeneralisedAdvectionDiffusionWindFctPt conserved_wind_fct_pt() const
Access function: Pointer to additoinal (conservative) wind function. Const version.
virtual unsigned u_index_cons_adv_diff() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
const double & pe_st() const
Peclet number multiplied by Strouhal number.
GeneralisedAdvectionDiffusionDiffFctPt & diff_fct_pt()
Access function: Pointer to diffusion function.
GeneralisedAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
GeneralisedAdvectionDiffusionEquations(const GeneralisedAdvectionDiffusionEquations &dummy)=delete
Broken copy constructor.
void operator=(const GeneralisedAdvectionDiffusionEquations &)=delete
Broken assignment operator.
void(* GeneralisedAdvectionDiffusionWindFctPt)(const Vector< double > &x, Vector< double > &wind)
Function pointer to wind function fct(x,w(x)) – x is a Vector!
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.
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
GeneralisedAdvectionDiffusionWindFctPt & conserved_wind_fct_pt()
Access function: Pointer to additional (conservative) wind function.
double integrate_u()
Integrate the concentration over the element.
void(* GeneralisedAdvectionDiffusionSourceFctPt)(const Vector< double > &x, double &f)
Function pointer to source function fct(x,f(x)) – x is a Vector!
GeneralisedAdvectionDiffusionSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
void enable_ALE()
(Re-)enable ALE, i.e. take possible mesh motion into account when evaluating the time-derivative....
virtual void fill_in_generic_residual_contribution_cons_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Add the element's contribution to its residual vector only (if flag=and/or element Jacobian matrix.
double interpolated_u_cons_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
virtual void get_conserved_wind_cons_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Get additional (conservative) wind at (Eulerian) position x and/or local coordinate s....
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
An OomphLibError object which should be thrown when an run-time error is encountered....
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
QGeneralisedAdvectionDiffusionElement(const QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void operator=(const QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
void output(FILE *file_pt)
C-style output function: x,y,u or x,y,z,u.
static const unsigned Initial_Nvalue
Static array of ints to hold number of variables at nodes: Initial_Nvalue[n].
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: x,y,u or x,y,z,u at n_plot^DIM plot points.
double dshape_and_dtest_eulerian_at_knot_cons_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape, test functions & derivs. w.r.t. to global coords. at integration point ipt....
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output function for a time-dependent exact solution. x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot ...
double dshape_and_dtest_eulerian_cons_adv_diff(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
QGeneralisedAdvectionDiffusionElement()
Constructor: Call constructors for QElement and Advection Diffusion equations.
unsigned required_nvalue(const unsigned &n) const
Required # of ‘values’ (pinned or dofs) at node n.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points.
void output(std::ostream &outfile)
Output function: x,y,u or x,y,z,u.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: x,y,u or x,y,z,u at n_plot^DIM plot points.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
////////////////////////////////////////////////////////////////////// //////////////////////////////...
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
bool is_steady() const
Flag to indicate if a timestepper has been made steady (possibly temporarily to switch off time-depen...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...