28 #ifndef OOMPH_AXISYM_ADV_DIFF_ELEMENTS_HEADER
29 #define OOMPH_AXISYM_ADV_DIFF_ELEMENTS_HEADER
34 #include <oomph-lib-config.h>
38 #include "../generic/nodes.h"
39 #include "../generic/Qelements.h"
40 #include "../generic/refineable_elements.h"
41 #include "../generic/oomph_utilities.h"
122 for (
unsigned t = 0;
t < n_time;
t++)
160 const unsigned& nplot)
const
167 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
186 file_out << wind[
i] << std::endl;
196 std::stringstream error_stream;
197 error_stream <<
"Advection Diffusion Elements only store 4 fields "
200 OOMPH_CURRENT_FUNCTION,
201 OOMPH_EXCEPTION_LOCATION);
219 return "Advection Diffusion";
224 std::stringstream error_stream;
225 error_stream <<
"Advection Diffusion Elements only store 4 fields "
228 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
243 void output(std::ostream& outfile,
const unsigned& nplot);
255 void output(FILE* file_pt,
const unsigned& n_plot);
260 const unsigned& nplot,
299 inline const double&
pe()
const
323 inline const double&
d()
const
340 double& source)
const
350 (*Source_fct_pt)(x, source);
367 for (
unsigned i = 0;
i < 3;
i++)
375 (*Wind_fct_pt)(x, wind);
383 const unsigned n_node =
nnode();
396 for (
unsigned j = 0; j < 2; j++)
402 for (
unsigned l = 0; l < n_node; l++)
404 const double u_value = this->
nodal_value(l, u_nodal_index);
406 flux[0] += u_value * dpsidx(l, 0);
407 flux[1] += u_value * dpsidx(l, 1);
440 const unsigned n_node =
nnode();
452 double interpolated_u = 0.0;
455 for (
unsigned l = 0; l < n_node; l++)
457 interpolated_u +=
nodal_value(l, u_nodal_index) * psi[l];
460 return (interpolated_u);
475 const unsigned n_node =
nnode();
487 unsigned n_u_dof = 0;
488 for (
unsigned l = 0; l < n_node; l++)
499 du_ddata.resize(n_u_dof, 0.0);
500 global_eqn_number.resize(n_u_dof, 0);
504 for (
unsigned l = 0; l < n_node; l++)
512 global_eqn_number[count] = global_eqn;
514 du_ddata[count] = psi[l];
533 DShape& dtestdx)
const = 0;
542 DShape& dtestdx)
const = 0;
594 template<
unsigned NNODE_1D>
596 :
public virtual QElement<2, NNODE_1D>,
636 void output(std::ostream& outfile,
const unsigned& n_plot)
651 void output(FILE* file_pt,
const unsigned& n_plot)
659 const unsigned& n_plot,
663 outfile, n_plot, exact_soln_pt);
697 template<
unsigned NNODE_1D>
698 double QAxisymAdvectionDiffusionElement<
706 double J = this->dshape_eulerian(
s, psi, dpsidx);
710 for (
unsigned i = 0;
i < NNODE_1D;
i++)
713 for (
unsigned j = 0; j < 2; j++)
715 dtestdx(
i, j) = dpsidx(
i, j);
731 template<
unsigned NNODE_1D>
740 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
754 template<
unsigned NNODE_1D>
756 :
public virtual QElement<1, NNODE_1D>
A class for all elements that solve the Advection Diffusion equations in a cylindrical polar coordina...
unsigned nscalar_paraview() const
Number of scalars/fields output by this element. Reimplements broken virtual function in base class.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
AxisymAdvectionDiffusionSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
virtual void fill_in_generic_residual_contribution_axi_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 * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
double *& d_pt()
Pointer to Peclet number multipled by Strouha number.
void disable_ALE()
Disable ALE, i.e. assert the mesh is not moving – you do this at your own risk!
unsigned self_test()
Self-test: Return 0 for OK.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: [du/dr,du/dz].
AxisymAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
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(* AxisymAdvectionDiffusionSourceFctPt)(const Vector< double > &x, double &f)
Function pointer to source function fct(x,f(x)) – x is a Vector!
void(* AxisymAdvectionDiffusionWindFctPt)(const Vector< double > &x, Vector< double > &wind)
Function pointer to wind function fct(x,w(x)) – x is a Vector!
virtual void get_wind_axi_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...
virtual unsigned u_index_axi_adv_diff() const
Broken assignment operator.
double du_dt_axi_adv_diff(const unsigned &n) const
du/dt at local node n. Uses suitably interpolated value for hanging nodes.
const double & pe() const
Peclet number.
virtual double dshape_and_dtest_eulerian_at_knot_axi_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 ...
void output(FILE *file_pt)
C_style output with default number of plot points.
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.
AxisymAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
virtual void dinterpolated_u_axi_adv_diff_ddata(const Vector< double > &s, Vector< double > &du_ddata, Vector< unsigned > &global_eqn_number)
Return derivative of u at point s with respect to all data that can affect its value....
static double Default_peclet_number
Static default value for the Peclet number.
double *& pe_st_pt()
Pointer to Peclet number multipled by Strouha number.
double interpolated_u_axi_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
virtual double dshape_and_dtest_eulerian_axi_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.
AxisymAdvectionDiffusionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
const double & d() const
Peclet number multiplied by Strouhal 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. Needs to be implemented for each new specif...
double *& pe_pt()
Pointer to Peclet number.
AxisymAdvectionDiffusionEquations()
Constructor: Initialise the Source_fct_pt and Wind_fct_pt to null and set (pointer to) Peclet number ...
AxisymAdvectionDiffusionWindFctPt wind_fct_pt() const
Access function: Pointer to wind function. Const version.
static double Default_diffusion_parameter
Static default value for the Peclet number.
bool ALE_is_disabled
Boolean flag to indicate whether AlE formulation is disable.
AxisymAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
virtual void get_source_axi_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...
AxisymAdvectionDiffusionEquations(const AxisymAdvectionDiffusionEquations &dummy)=delete
Broken copy constructor.
double * D_pt
Pointer to global Diffusion parameter.
void output(std::ostream &outfile)
Output with default number of plot points.
void enable_ALE()
(Re-)enable ALE, i.e. take possible mesh motion into account when evaluating the time-derivative....
std::string scalar_name_paraview(const unsigned &i) const
Name of the i-th scalar field. Default implementation returns V1 for the first one,...
const double & pe_st() const
Peclet number multiplied by Strouhal number.
double * Pe_pt
Pointer to global Peclet number.
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
long & eqn_number(const unsigned &i)
Return the equation number of the i-th stored variable.
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.
virtual unsigned nplot_points_paraview(const unsigned &nplot) const
Return the number of actual plot points for paraview plot with parameter nplot. Broken virtual; can b...
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 .
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...
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...
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.
QAxisymAdvectionDiffusionElement()
Constructor: Call constructors for QElement and Advection Diffusion equations.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: r,z,u at n_plot^2 plot points.
static const unsigned Initial_Nvalue
Static array of ints to hold number of variables at nodes: Initial_Nvalue[n].
void output(std::ostream &outfile)
Output function: r,z,u.
double dshape_and_dtest_eulerian_at_knot_axi_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(FILE *file_pt, const unsigned &n_plot)
C-style output function: r,z,u at n_plot^2 plot points.
QAxisymAdvectionDiffusionElement(const QAxisymAdvectionDiffusionElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
double dshape_and_dtest_eulerian_axi_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.
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
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...
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
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).
//////////////////////////////////////////////////////////////////// ////////////////////////////////...