27#ifndef OOMPH_POISSON_ELEMENTS_HEADER
28#define OOMPH_POISSON_ELEMENTS_HEADER
33#include <oomph-lib-config.h>
54 template<
unsigned DIM>
102 for (
unsigned i = 0;
i <
dim;
i++)
121 const unsigned&
nplot)
const
128 <<
"Poisson elements only store a single field so i must be 0 rather"
129 <<
" than " <<
i << std::endl;
156 <<
"Poisson elements only store a single field so i must be 0 rather"
157 <<
" than " <<
i << std::endl;
163 return "Poisson solution";
172 const unsigned&
nplot,
179 error_stream <<
"Poisson equation has only one field. Can't call "
180 <<
"this function for value " <<
i << std::endl;
216 const unsigned n_plot = 5;
227 const unsigned n_plot = 5;
251 "There is no time-dependent output_fct() for Poisson elements ",
275 "There is no time-dependent compute_error() for Poisson elements",
311 double& source)
const
321 (*Source_fct_pt)(x, source);
347 for (
unsigned i = 0;
i <
DIM;
i++)
358 (*Source_fct_gradient_pt)(x,
gradient);
380 for (
unsigned j = 0;
j <
DIM;
j++)
389 for (
unsigned j = 0;
j <
DIM;
j++)
413 for (
unsigned i = 0;
i <
DIM;
i++)
460 double interpolated_u = 0.0;
468 return (interpolated_u);
520 const unsigned&
flag);
539 template<
unsigned DIM,
unsigned NNODE_1D>
663 template<
unsigned DIM,
unsigned NNODE_1D>
689 template<
unsigned DIM,
unsigned NNODE_1D>
717 template<
unsigned DIM,
unsigned NNODE_1D>
754 template<
unsigned DIM,
unsigned NNODE_1D>
756 :
public virtual QElement<DIM - 1, NNODE_1D>
772 template<
unsigned NNODE_1D>
790 template<
class POISSON_ELEMENT>
804 error_stream <<
"Poisson elements only store a single field so fld "
806 <<
" than " <<
fld << std::endl;
817 for (
unsigned j = 0;
j <
nnod;
j++)
841 error_stream <<
"Poisson elements only store a single field so fld "
843 <<
" than " <<
fld << std::endl;
868 error_stream <<
"Poisson elements only store a single field so fld "
870 <<
" than " <<
fld << std::endl;
895 error_stream <<
"Poisson elements only store a single field so fld "
897 <<
" than " <<
fld << std::endl;
913 double interpolated_u = 0.0;
920 return interpolated_u;
931 error_stream <<
"Poisson elements only store a single field so fld "
933 <<
" than " <<
fld << std::endl;
938 return this->
nnode();
949 error_stream <<
"Poisson elements only store a single field so fld "
951 <<
" than " <<
fld << std::endl;
966 template<
class ELEMENT>
979 template<
class ELEMENT>
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
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...
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 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-...
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
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...
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.
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
static double Default_fd_jacobian_step
Double used for the default finite difference step in elemental jacobian calculations.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
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.
A class for all isoparametric elements that solve the Poisson equations.
virtual void fill_in_generic_residual_contribution_poisson(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Compute element residual Vector only (if flag=and/or element Jacobian matrix.
void point_output_data(const Vector< double > &s, Vector< double > &data)
Output solution in data vector at local cordinates s: x,y [,z], u.
unsigned self_test()
Self-test: Return 0 for OK.
virtual void get_source_gradient_poisson(const unsigned &ipt, const Vector< double > &x, Vector< double > &gradient) const
Get gradient of source term at (Eulerian) position x. This function is virtual to allow overloading i...
virtual void get_source_poisson(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...
virtual double dshape_and_dtest_eulerian_at_knot_poisson(const unsigned &ipt, Shape &psi, DShape &dpsidx, RankFourTensor< double > &d_dpsidx_dX, Shape &test, DShape &dtestdx, RankFourTensor< double > &d_dtestdx_dX, DenseMatrix< double > &djacobian_dX) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
void get_dflux_dnodal_u(const Vector< double > &s, Vector< Vector< double > > &dflux_dnodal_u) const
Get derivative of flux w.r.t. to nodal values: dflux_dnodal_u[i][j] = d ( du/dx_i ) / dU_j.
PoissonSourceFctGradientPt source_fct_gradient_pt() const
Access function: Pointer to gradient source function. Const version.
void(* PoissonSourceFctPt)(const Vector< double > &x, double &f)
Function pointer to source function fct(x,f(x)) – x is a Vector!
void output(std::ostream &outfile)
Output with default number of plot points.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its 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...
PoissonEquations()
Constructor (must initialise the Source_fct_pt to null)
PoissonSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
PoissonSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
void compute_norm(double &norm)
Compute norm of solution: square of the L2 norm.
PoissonSourceFctPt Source_fct_pt
Pointer to source function:
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
PoissonSourceFctGradientPt Source_fct_gradient_pt
Pointer to gradient of source function.
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Compute derivatives of elemental residual vector with respect to nodal coordinates....
void(* PoissonSourceFctGradientPt)(const Vector< double > &x, Vector< double > &gradient)
Function pointer to gradient of source function fct(x,g(x)) – x is a Vector!
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points (dummy time-dependent versi...
void operator=(const PoissonEquations &)=delete
Broken assignment operator.
virtual double dshape_and_dtest_eulerian_poisson(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...
PoissonSourceFctGradientPt & source_fct_gradient_pt()
Access function: Pointer to gradient of source function.
unsigned nscalar_paraview() const
Number of scalars/fields output by this element. Reimplements broken virtual function in base class.
virtual double interpolated_u_poisson(const Vector< double > &s) const
Return FE representation of function value u_poisson(s) at local coordinate s.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
std::string scalar_name_paraview(const unsigned &i) const
Name of the i-th scalar field. Default implementation returns V1 for the first one,...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the element's contribution to its residual vector and element Jacobian matrix (wrapper)
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points.
virtual unsigned u_index_poisson() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
virtual double dshape_and_dtest_eulerian_at_knot_poisson(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 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...
void output(FILE *file_pt)
C_style output with default number of plot points.
PoissonEquations(const PoissonEquations &dummy)=delete
Broken copy constructor.
Wrapper class for projectable elements. Adds "projectability" to the underlying ELEMENT.
Poisson upgraded to become projectable.
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.
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.
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 ...
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_coordinate_projection()
Number of positional history values (Note: count includes current value!)
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field (includes current value!)
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
QPoissonElement elements are linear/quadrilateral/brick-shaped Poisson elements with isoparametric in...
double dshape_and_dtest_eulerian_at_knot_poisson(const unsigned &ipt, Shape &psi, DShape &dpsidx, RankFourTensor< double > &d_dpsidx_dX, Shape &test, DShape &dtestdx, RankFourTensor< double > &d_dtestdx_dX, DenseMatrix< double > &djacobian_dX) const
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
double dshape_and_dtest_eulerian_at_knot_poisson(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 ...
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.
QPoissonElement(const QPoissonElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
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.
QPoissonElement()
Constructor: Call constructors for QElement and Poisson equations.
void output(std::ostream &outfile)
Output function: x,y,u or x,y,z,u.
double dshape_and_dtest_eulerian_poisson(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: x,y,u or x,y,z,u.
static const unsigned Initial_Nvalue
Static int that holds the number of variables at nodes: always the same.
void operator=(const QPoissonElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
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.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
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)
A slight extension to the standard template vector class so that we can include "graceful" array rang...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).