27#ifndef OOMPH_SPACE_TIME_UNSTEADY_HEAT_ELEMENTS_HEADER
28#define OOMPH_SPACE_TIME_UNSTEADY_HEAT_ELEMENTS_HEADER
32#include <oomph-lib-config.h>
74 template<
unsigned SPATIAL_DIM>
164 const unsigned&
nplot,
172 const unsigned&
nplot,
210 const unsigned&
nplot)
const
216 error_stream <<
"Space-time unsteady heat elements only store a single "
217 <<
"field so i must be 0 rather than " <<
i << std::endl;
246 const unsigned&
nplot,
253 error_stream <<
"Space-time unsteady heat elements only store a single "
254 <<
"field so i must be 0 rather than " <<
i << std::endl;
299 const unsigned&
nplot,
307 error_stream <<
"Space-time unsteady heat elements only store a single "
308 <<
"field so i must be 0 rather than " <<
i << std::endl;
368 error_stream <<
"These unsteady heat elements only store 1 field, \n"
369 <<
"but i is currently " <<
i << std::endl;
400 double& source)
const
412 (*Source_fct_pt)(
t, x, source);
522 double interpolated_u = 0.0;
532 return interpolated_u;
570 double interpolated_u = 0.0;
580 return interpolated_u;
664 const unsigned&
flag);
701 template<
unsigned SPATIAL_DIM,
unsigned NNODE_1D>
703 :
public virtual QElement<SPATIAL_DIM + 1, NNODE_1D>,
820 template<
unsigned SPATIAL_DIM,
unsigned NNODE_1D>
848 template<
unsigned SPATIAL_DIM,
unsigned NNODE_1D>
885 template<
unsigned SPATIAL_DIM,
unsigned NNODE_1D>
887 :
public virtual QElement<SPATIAL_DIM, NNODE_1D>
905 template<
unsigned NNODE_1D>
924 template<
class UNSTEADY_HEAT_ELEMENT>
947 error_stream <<
"SpaceTimeUnsteadyHeat elements only store a single "
948 <<
"field so fld must be 0 rather than " <<
fld
964 for (
unsigned j = 0;
j <
nnod;
j++)
995 error_stream <<
"SpaceTimeUnsteadyHeat elements only store a single "
996 <<
"field so fld must be 0 rather than " <<
fld
1033 error_stream <<
"SpaceTimeUnsteadyHeat elements only store a single "
1034 <<
"field so fld must be 0 rather than " <<
fld
1072 const unsigned&
fld,
1083 error_stream <<
"SpaceTimeUnsteadyHeat elements only store a single "
1084 <<
"field so fld must be 0 rather than " <<
fld
1106 double interpolated_u = 0.0;
1116 return interpolated_u;
1131 error_stream <<
"SpaceTimeUnsteadyHeat elements only store a single "
1132 <<
"field so fld must be 0 rather than " <<
fld
1142 return this->
nnode();
1157 error_stream <<
"SpaceTimeUnsteadyHeat elements only store a single "
1158 <<
"field so fld must be 0 rather than " <<
fld
1205 outfile << this->interpolated_u_ust_heat(
s) <<
" ";
1208 outfile << this->interpolated_du_dt_ust_heat(
s) <<
" ";
1215 for (
unsigned t = 1;
t < n_prev;
t++)
1229 for (
unsigned t = 1;
t < n_prev;
t++)
1232 outfile << this->interpolated_u_ust_heat(
t,
s) <<
" ";
1249 template<
class ELEMENT>
1262 template<
class ELEMENT>
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.
unsigned ntstorage() const
Return total number of doubles stored per value to record time history of each value (one for steady ...
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
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...
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual void local_coordinate_of_node(const unsigned &j, Vector< double > &s) const
Get local coordinates of node j in the element; vector sets its own size (broken virtual)
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 std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
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.
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 dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
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...
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
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...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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 DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
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.
Wrapper class for projectable elements. Adds "projectability" to the underlying ELEMENT.
SpaceTimeUnsteadyHeat upgraded to become projectable.
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Return interpolated field fld at local coordinate s, at time level t (t=0: present; t>0: history valu...
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. (Note: count includes current value!...
void output(std::ostream &outfile, const unsigned &nplot)
Output FE representation of soln: x,t,u or x,y,t,u at n_plot^(SPATIAL_DIM+1) plot points.
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Specify the values associated with field fld. The information is returned in a vector of pairs which ...
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (Note: count includes current value!)
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
unsigned nfields_for_projection()
Number of fields to be projected: Just one.
ProjectableUnsteadyHeatSpaceTimeElement()
Constructor [this was only required explicitly from gcc 4.5.2 onwards...].
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Return Jacobian of mapping and shape functions of field fld at local coordinate s.
QUnsteadyHeatSpaceTimeElement elements are quadrilateral/brick- shaped UnsteadyHeat elements with iso...
QUnsteadyHeatSpaceTimeElement(const QUnsteadyHeatSpaceTimeElement< SPATIAL_DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
double dshape_and_dtest_eulerian_at_knot_ust_heat(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
double dshape_and_dtest_eulerian_ust_heat(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, const unsigned &n_plot)
C-style output function: x,t,u or x,y,t,u at n_plot^(SPATIAL_DIM+1) plot points.
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,t,u_exact or x,y,t,u_exact at n_plot^(SPATIAL_...
QUnsteadyHeatSpaceTimeElement()
Constructor: Call constructors for QElement and SpaceTimeUnsteadyHeat equations.
void output(FILE *file_pt)
C-style output function: x,t,u or x,y,t,u.
void output(std::ostream &outfile)
Output function: x,t,u or x,y,t,u.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: x,t,u or x,y,t,u at n_plot^(SPATIAL_DIM+1) plot points.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: x,t,u_exact or x,y,t,u_exact at n_plot^(SPATIAL_DIM+1) plot po...
static const unsigned Initial_Nvalue
Static array of ints to hold number of variables at nodes: Initial_Nvalue[n].
unsigned required_nvalue(const unsigned &n) const
Required number of 'values' (pinned or dofs) at node n.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Base class so that we don't need to know the dimension just to set the source function!
void(* SpaceTimeUnsteadyHeatSourceFctPt)(const double &time, const Vector< double > &x, double &u)
Function pointer to source function fct(t,x,f(x,t)) – x is a Vector!
virtual SpaceTimeUnsteadyHeatSourceFctPt & source_fct_pt()=0
Access function: Pointer to source function.
A class for all isoparametric elements that solve the SpaceTimeUnsteadyHeat equations.
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^SPATIAL_DIM plot points.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error and norm against exact solution.
double *& alpha_pt()
Pointer to Alpha parameter (thermal inertia)
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute element residual Vector (wrapper)
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...
SpaceTimeUnsteadyHeatEquations(const SpaceTimeUnsteadyHeatEquations &dummy)=delete
Broken copy constructor.
void(* SpaceTimeUnsteadyHeatSourceFctPt)(const double &time, const Vector< double > &x, double &u)
Function pointer to source function fct(t,x,f(x,t)) – x is a Vector! DRAIG: Why is this here?...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute element residual Vector and element Jacobian matrix (wrapper)
double *& beta_pt()
Pointer to Beta parameter (thermal conductivity)
void scalar_value_fct_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specif...
double interpolated_u_ust_heat(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.
double * Beta_pt
Pointer to Beta parameter (thermal conductivity)
unsigned nscalar_paraview() const
Number of scalars/fields output by this element. Reimplements broken virtual function in base class.
SpaceTimeUnsteadyHeatSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
void compute_norm(double &norm)
Compute norm of FE solution.
double * Alpha_pt
Pointer to Alpha parameter (thermal inertia)
unsigned self_test()
Self-test: Return 0 for OK.
SpaceTimeUnsteadyHeatSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
static double Default_beta_parameter
Static default value for the Beta parameter (thermal conductivity): One for natural scaling.
std::string scalar_name_paraview(const unsigned &i) const
Name of the i-th scalar field. Default implementation returns V1 for the first one,...
virtual double dshape_and_dtest_eulerian_at_knot_ust_heat(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(std::ostream &outfile)
Output with default number of plot points.
double interpolated_du_dt_ust_heat(const Vector< double > &s) const
Return FE representation of function value du/dt(s) at local coordinate s.
virtual void fill_in_generic_residual_contribution_ust_heat(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Compute element residual Vector only (if flag=and/or element Jacobian matrix.
virtual double dshape_and_dtest_eulerian_ust_heat(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 coordinate s; return Jacobian of map...
void disable_ALE()
Disable ALE, i.e. assert the mesh is not moving – you do this at your own risk!
void output_element_paraview(std::ofstream &outfile, const unsigned &nplot)
C-style output FE representation of soln: x,y,u or x,y,z,u at nplot^SPATIAL_DIM plot points.
virtual unsigned u_index_ust_heat() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
SpaceTimeUnsteadyHeatSourceFctPt Source_fct_pt
Pointer to source function:
const double & beta() const
Beta parameter (thermal conductivity)
double interpolated_u_ust_heat(const unsigned &t, const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s at previous time t (t=0: presen...
void scalar_value_fct_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specif...
const double & alpha() const
Alpha parameter (thermal inertia)
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
double du_dt_ust_heat(const unsigned &n) const
Calculate du/dt at the n-th local node. Uses suitably interpolated value for hanging nodes.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i]=du/dx_i.
SpaceTimeUnsteadyHeatEquations()
Constructor: Initialises the Source_fct_pt to null and sets flag to use ALE formulation of the equati...
void enable_ALE()
(Re-)enable ALE, i.e. take possible mesh motion into account when evaluating the time-derivative....
static double Default_alpha_parameter
Static default value for the Alpha parameter (thermal inertia): One for natural scaling.
virtual void get_source_ust_heat(const double &t, const unsigned &ipt, const Vector< double > &x, double &source) const
Get source term at continous time t and (Eulerian) position x. Virtual so it can be overloaded in der...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
TAdvectionDiffusionReactionElement()
Constructor: Call constructors for TElement and AdvectionDiffusionReaction equations.
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).