27 #ifndef OOMPH_GEN_AXISYM_ADV_DIFF_ELEMENTS_HEADER
28 #define OOMPH_GEN_AXISYM_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"
62 typedef void (*GeneralisedAxisymAdvectionDiffusionSourceFctPt)(
67 typedef void (*GeneralisedAxisymAdvectionDiffusionWindFctPt)(
110 virtual inline unsigned u_index_cons_axisym_adv_diff()
const
128 const unsigned u_nodal_index = u_index_cons_axisym_adv_diff();
133 for (
unsigned t = 0;
t < n_time;
t++)
169 void output(std::ostream& outfile,
const unsigned& nplot);
180 void output(FILE* file_pt,
const unsigned& n_plot);
185 const unsigned& nplot,
192 std::ostream& outfile,
193 const unsigned& nplot,
197 throw OomphLibError(
"There is no time-dependent output_fct() for "
198 "Advection Diffusion elements",
199 OOMPH_CURRENT_FUNCTION,
200 OOMPH_EXCEPTION_LOCATION);
219 "No time-dependent compute_error() for Advection Diffusion elements",
220 OOMPH_CURRENT_FUNCTION,
221 OOMPH_EXCEPTION_LOCATION);
284 const double&
pe()
const
313 double& source)
const
323 (*Source_fct_pt)(x, source);
342 for (
unsigned i = 0;
i < 3;
i++)
350 (*Wind_fct_pt)(x,
wind);
370 for (
unsigned i = 0;
i < 3;
i++)
378 (*Conserved_wind_fct_pt)(x,
wind);
400 for (
unsigned i = 0;
i < 3;
i++)
402 for (
unsigned j = 0; j < 3; j++)
418 (*Diff_fct_pt)(x,
D);
427 const unsigned n_node = this->
nnode();
430 const unsigned u_nodal_index = this->u_index_cons_axisym_adv_diff();
440 for (
unsigned j = 0; j < 2; j++)
446 for (
unsigned l = 0; l < n_node; l++)
448 const double u_value = this->
nodal_value(l, u_nodal_index);
450 flux[0] += u_value * dpsidx(l, 0);
451 flux[1] += u_value * dpsidx(l, 1);
460 const unsigned n_node =
nnode();
463 const unsigned u_nodal_index = u_index_cons_axisym_adv_diff();
475 double interpolated_u = 0.0;
480 for (
unsigned l = 0; l < n_node; l++)
483 const double u_value = this->
nodal_value(l, u_nodal_index);
484 interpolated_u += u_value * psi(l);
486 for (
unsigned j = 0; j < 2; j++)
489 interpolated_dudx[j] += u_value * dpsidx(l, j);
509 for (
unsigned i = 0;
i < 2;
i++)
512 for (
unsigned j = 0; j < 2; j++)
514 total_flux[
i] +=
D(
i, j) * interpolated_dudx[j];
516 total_flux[
i] -= conserved_wind[
i] * interpolated_u;
554 residuals, jacobian, mass_matrix, 2);
563 unsigned n_node =
nnode();
566 unsigned u_nodal_index = u_index_cons_axisym_adv_diff();
575 double interpolated_u = 0.0;
578 for (
unsigned l = 0; l < n_node; l++)
580 interpolated_u +=
nodal_value(l, u_nodal_index) * psi[l];
583 return (interpolated_u);
598 DShape& dtestdx)
const = 0;
607 DShape& dtestdx)
const = 0;
656 template<
unsigned NNODE_1D>
658 :
public virtual QElement<2, NNODE_1D>,
699 void output(std::ostream& outfile,
const unsigned& n_plot)
714 void output(FILE* file_pt,
const unsigned& n_plot)
722 const unsigned& n_plot,
726 outfile, n_plot, exact_soln_pt);
734 const unsigned& n_plot,
739 outfile, n_plot, time, exact_soln_pt);
772 template<
unsigned NNODE_1D>
781 double J = this->dshape_eulerian(
s, psi, dpsidx);
785 for (
unsigned i = 0;
i < NNODE_1D;
i++)
788 for (
unsigned j = 0; j < 2; j++)
790 dtestdx(
i, j) = dpsidx(
i, j);
805 template<
unsigned NNODE_1D>
815 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
837 template<
unsigned NNODE_1D>
839 :
public virtual QElement<1, 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.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
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 .
GeneralisedAxisymAdvectionDiffusionWindFctPt wind_fct_pt() const
Access function: Pointer to wind function. Const version.
GeneralisedAxisymAdvectionDiffusionDiffFctPt Diff_fct_pt
Pointer to diffusivity funciton.
double *& pe_st_pt()
Pointer to Peclet number multipled by Strouha number.
const double & pe_st() const
Peclet number multiplied by Strouhal number.
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
double * Pe_pt
Pointer to global Peclet number.
double interpolated_u_cons_axisym_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
void output(FILE *file_pt)
C_style output with default number of plot points.
GeneralisedAxisymAdvectionDiffusionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
double du_dt_cons_axisym_adv_diff(const unsigned &n) const
du/dt at local node n. Uses suitably interpolated value for hanging nodes.
virtual void get_conserved_wind_cons_axisym_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....
virtual double dshape_and_dtest_eulerian_at_knot_cons_axisym_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
GeneralisedAxisymAdvectionDiffusionDiffFctPt diff_fct_pt() const
Access function: Pointer to diffusion function. Const version.
virtual void get_diff_cons_axisym_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...
virtual void output_fct(std::ostream &outfile, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,z,,u_exact at nplot^2 plot points (dummy time-dependent version to keep intel co...
void get_total_flux(const Vector< double > &s, Vector< double > &total_flux) const
Get flux: .
GeneralisedAxisymAdvectionDiffusionWindFctPt & conserved_wind_fct_pt()
Access function: Pointer to additional (conservative) wind function.
static double Default_peclet_number
Static default value for the Peclet number.
GeneralisedAxisymAdvectionDiffusionWindFctPt conserved_wind_fct_pt() const
Access function: Pointer to additional (conservative) wind function. Const version.
virtual void get_wind_cons_axisym_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...
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,z,u_exact at nplot^2 plot points.
GeneralisedAxisymAdvectionDiffusionSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
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.
GeneralisedAxisymAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
double integrate_u()
Integrate the concentration over the element.
virtual void fill_in_generic_residual_contribution_cons_axisym_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.
void disable_ALE()
Disable ALE, i.e. assert the mesh is not moving – you do this at your own risk!
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
double *& pe_pt()
Pointer to Peclet number.
Function pointer to wind function Vector< double > & wind
GeneralisedAxisymAdvectionDiffusionDiffFctPt & diff_fct_pt()
Access function: Pointer to diffusion function.
unsigned self_test()
Self-test: Return 0 for OK.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
virtual double dshape_and_dtest_eulerian_cons_axisym_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...
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
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....
GeneralisedAxisymAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
virtual void get_source_cons_axisym_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...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
GeneralisedAxisymAdvectionDiffusionWindFctPt Conserved_wind_fct_pt
Pointer to additional (conservative) wind function:
Function pointer to a diffusivity function typedef void(* GeneralisedAxisymAdvectionDiffusionDiffFctPt)(const Vector< double > &x, DenseMatrix< double > &D)
Function pointer to source function double & f
const double & pe() const
Peclet number.
Broken copy constructor GeneralisedAxisymAdvectionDiffusionEquations(const GeneralisedAxisymAdvectionDiffusionEquations &dummy)=delete
GeneralisedAxisymAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
void output(std::ostream &outfile)
Output with default number of plot points.
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....
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: r,z,u_exact at n_plot^2 plot points.
void output(std::ostream &outfile)
Output function: r,z,u.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: r,z,u at n_plot^2 plot points.
double dshape_and_dtest_eulerian_at_knot_cons_axisym_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....
static const unsigned Initial_Nvalue
Static array of ints to hold number of variables at nodes: Initial_Nvalue[n].
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. r,z,u_exact at n_plot^2 plot points (Calls the s...
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
double dshape_and_dtest_eulerian_cons_axisym_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.
void output(FILE *file_pt)
C-style output function: r,z,u.
QGeneralisedAxisymAdvectionDiffusionElement(const QGeneralisedAxisymAdvectionDiffusionElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: r,z,u at n_plot^2 plot points.
QGeneralisedAxisymAdvectionDiffusionElement()
Constructor: Call constructors for QElement and Advection Diffusion equations.
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...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...