31 #ifndef OOMPH_GEN_ADV_DIFF_ELEMENTS_HEADER 32 #define OOMPH_GEN_ADV_DIFF_ELEMENTS_HEADER 37 #include <oomph-lib-config.h> 41 #include "../generic/nodes.h" 42 #include "../generic/Qelements.h" 43 #include "../generic/oomph_utilities.h" 59 template <
unsigned DIM>
129 const unsigned n_time = time_stepper_pt->
ntstorage();
131 for(
unsigned t=0;
t<n_time;
t++)
166 void output(std::ostream &outfile,
const unsigned &nplot);
178 void output(FILE* file_pt,
const unsigned &n_plot);
182 void output_fct(std::ostream &outfile,
const unsigned &nplot,
189 virtual void output_fct(std::ostream &outfile,
const unsigned &nplot,
194 "There is no time-dependent output_fct() for Advection Diffusion elements",
195 OOMPH_CURRENT_FUNCTION,
196 OOMPH_EXCEPTION_LOCATION);
203 exact_soln_pt,
double& error,
double& norm);
210 const double& time,
double& error,
double& norm)
213 "No time-dependent compute_error() for Advection Diffusion elements",
214 OOMPH_CURRENT_FUNCTION,
215 OOMPH_EXCEPTION_LOCATION);
278 double& source)
const 285 (*Source_fct_pt)(x,source);
302 for(
unsigned i=0;
i<DIM;
i++) {wind[
i]= 0.0;}
307 (*Wind_fct_pt)(x,wind);
327 for(
unsigned i=0;
i<DIM;
i++) {wind[
i]= 0.0;}
332 (*Conserved_wind_fct_pt)(x,wind);
351 for(
unsigned i=0;
i<DIM;
i++)
353 for(
unsigned j=0;j<DIM;j++)
355 if(
i==j) {
D(
i,j) = 1.0;}
else {
D(
i,j) = 0.0;}
371 unsigned n_node =
nnode();
378 DShape dpsidx(n_node,DIM);
384 for(
unsigned j=0;j<DIM;j++) {flux[j] = 0.0;}
387 for(
unsigned l=0;l<n_node;l++)
390 for(
unsigned j=0;j<DIM;j++)
392 flux[j] +=
nodal_value(l,u_nodal_index)*dpsidx(l,j);
402 const unsigned n_node =
nnode();
409 DShape dpsidx(n_node,DIM);
417 double interpolated_u = 0.0;
422 for(
unsigned l=0;l<n_node;l++)
425 const double u_value = this->
nodal_value(l,u_nodal_index);
426 interpolated_u += u_value*psi(l);
428 for(
unsigned j=0;j<DIM;j++)
431 interpolated_dudx[j] += u_value*dpsidx(l,j);
448 for(
unsigned i=0;
i<DIM;
i++)
451 for(
unsigned j=0;j<DIM;j++)
453 total_flux[
i] +=
D(
i,j)*interpolated_dudx[j];
455 total_flux[
i] -= conserved_wind[
i]*interpolated_u;
490 jacobian,mass_matrix,2);
498 unsigned n_node =
nnode();
510 double interpolated_u = 0.0;
513 for(
unsigned l=0;l<n_node;l++)
515 interpolated_u +=
nodal_value(l,u_nodal_index)*psi[l];
518 return(interpolated_u);
595 template <
unsigned DIM,
unsigned NNODE_1D>
597 public virtual QElement<DIM,NNODE_1D>,
632 {
return Initial_Nvalue;}
641 void output(std::ostream &outfile,
const unsigned &n_plot)
654 void output(FILE* file_pt,
const unsigned &n_plot)
661 void output_fct(std::ostream &outfile,
const unsigned &n_plot,
670 void output_fct(std::ostream &outfile,
const unsigned &n_plot,
676 output_fct(outfile,n_plot,time,exact_soln_pt);
711 template<
unsigned DIM,
unsigned NNODE_1D>
724 for(
unsigned i=0;
i<NNODE_1D;
i++)
727 for(
unsigned j=0;j<DIM;j++)
729 dtestdx(
i,j) = dpsidx(
i,j);
745 template<
unsigned DIM,
unsigned NNODE_1D>
778 template<
unsigned DIM,
unsigned NNODE_1D>
780 public virtual QElement<DIM-1,NNODE_1D>
801 template<
unsigned NNODE_1D>
GeneralisedAdvectionDiffusionEquations()
Constructor: Initialise the Source_fct_pt and Wind_fct_pt to null and set (pointer to) Peclet number ...
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
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.
void output(FILE *file_pt)
C_style output with default number of plot points.
void(* GeneralisedAdvectionDiffusionWindFctPt)(const Vector< double > &x, Vector< double > &wind)
Function pointer to wind function fct(x,w(x)) – x is a Vector!
double du_dt_cons_adv_diff(const unsigned &n) const
du/dt at local node n.
virtual double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Return the geometric shape functions and also first derivatives w.r.t. global coordinates at the ipt-...
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.
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...
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_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.
void enable_ALE()
(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.
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...
unsigned required_nvalue(const unsigned &n) const
Required # of `values' (pinned or dofs) at node n.
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
double interpolated_u_cons_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
double integrate_u()
Integrate the concentration over the element.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement. ...
A general Finite Element class.
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(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as ...
void get_total_flux(const Vector< double > &s, Vector< double > &total_flux) const
Get flux: .
GeneralisedAdvectionDiffusionDiffFctPt & diff_fct_pt()
Access function: Pointer to diffusion function.
QGeneralisedAdvectionDiffusionElement(const QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D > &dummy)
Broken copy constructor.
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 operator=(const QGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D > &)
Broken assignment operator.
A class for all elements that solve the Advection Diffusion equations in conservative form using isop...
GeneralisedAdvectionDiffusionDiffFctPt diff_fct_pt() const
Access function: Pointer to diffusion function. Const version.
void operator=(const GeneralisedAdvectionDiffusionEquations &)
Broken assignment operator.
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
GeneralisedAdvectionDiffusionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
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...
const double & pe_st() const
Peclet number multiplied by Strouhal number.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as ...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
void output(FILE *file_pt)
C-style output function: x,y,u or x,y,z,u.
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...
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed...
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...
double * Pe_pt
Pointer to global Peclet number.
void(* GeneralisedAdvectionDiffusionSourceFctPt)(const Vector< double > &x, double &f)
Function pointer to source function fct(x,f(x)) – x is a Vector!
const double & pe() const
Peclet number.
GeneralisedAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
GeneralisedAdvectionDiffusionWindFctPt wind_fct_pt() const
Access function: Pointer to wind function. Const version.
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 ...
void(* GeneralisedAdvectionDiffusionDiffFctPt)(const Vector< double > &x, DenseMatrix< double > &D)
Funciton pointer to a diffusivity function.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
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...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
GeneralisedAdvectionDiffusionSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
void output(std::ostream &outfile)
Output function: x,y,u or x,y,z,u.
GeneralisedAdvectionDiffusionWindFctPt Conserved_wind_fct_pt
Pointer to additional (conservative) wind function:
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...
static double Default_peclet_number
Static default value for the Peclet number.
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
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. Return Jacobian.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement. ...
unsigned self_test()
Self-test: Return 0 for OK.
QGeneralisedAdvectionDiffusionElement elements are linear/quadrilateral/brick-shaped Advection Diffus...
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 ...
void broken_assign(const std::string &class_name)
Issue error message and terminate execution.
void output(std::ostream &outfile)
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 const unsigned Initial_Nvalue
Static array of ints to hold number of variables at nodes: Initial_Nvalue[n].
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...
GeneralisedAdvectionDiffusionDiffFctPt Diff_fct_pt
Pointer to diffusivity funciton.
GeneralisedAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
GeneralisedAdvectionDiffusionEquations(const GeneralisedAdvectionDiffusionEquations &dummy)
Broken copy constructor.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
double *& pe_pt()
Pointer to Peclet number.
double *& pe_st_pt()
Pointer to Peclet number multipled by Strouha number.
GeneralisedAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
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...
bool is_steady() const
Flag to indicate if a timestepper has been made steady (possibly temporarily to switch off time-depen...
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) ...
QGeneralisedAdvectionDiffusionElement()
Constructor: Call constructors for QElement and Advection Diffusion equations.
unsigned nnode() const
Return the number of nodes.
Base class for time-stepping schemes. Timestepper provides an approximation of the temporal derivativ...
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
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...
GeneralisedAdvectionDiffusionWindFctPt conserved_wind_fct_pt() const
Access function: Pointer to additoinal (conservative) wind function. Const version.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
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.
GeneralisedAdvectionDiffusionWindFctPt & conserved_wind_fct_pt()
Access function: Pointer to additional (conservative) wind function.