31 #ifndef OOMPH_SPECTRAL_POISSON_ELEMENTS_HEADER 32 #define OOMPH_SPECTRAL_POISSON_ELEMENTS_HEADER 36 #include <oomph-lib-config.h> 41 #include "../generic/Qspectral_elements.h" 53 template <
unsigned DIM,
unsigned NNODE_1D>
101 void output(std::ostream &outfile,
const unsigned &n_plot)
113 void output(FILE* file_pt,
const unsigned &n_plot)
119 void output_fct(std::ostream &outfile,
const unsigned &n_plot,
128 void output_fct(std::ostream &outfile,
const unsigned &n_plot,
174 template<
unsigned DIM,
unsigned NNODE_1D>
188 unsigned nnod=this->
nnode();
189 for(
unsigned i=0;
i<nnod;
i++)
192 for(
unsigned j=0;j<DIM;j++)
194 dtestdx(
i,j) = dpsidx(
i,j);
208 template<
unsigned DIM,
unsigned NNODE_1D>
236 template<
unsigned DIM,
unsigned NNODE_1D>
250 djacobian_dX,d_dpsidx_dX);
255 d_dtestdx_dX = d_dpsidx_dX;
273 template<
unsigned DIM,
unsigned NNODE_1D>
290 template<
unsigned NNODE_1D>
void output(FILE *file_pt)
C-style output function: x,y,u or x,y,z,u.
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
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(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as ...
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement. ...
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(std::ostream &outfile)
Output function: x,y,u or x,y,z,u.
void output(std::ostream &outfile)
Output with default number of plot points.
QSpectralPoissonElement()
Constructor: Call constructors for QSpectralElement and Poisson equations.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as ...
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.
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.
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 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 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 ...
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement. ...
QSpectralPoissonElement(const QSpectralPoissonElement< DIM, NNODE_1D > &dummy)
Broken copy constructor.
unsigned nnode() const
Return the number of nodes.
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, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: x,y,u_exact or x,y,z,u_exact at n_plot^DIM plot points.
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. Return Jacobian.