27 #ifndef OOMPH_POISSON_ELEMENTS_HEADER
28 #define OOMPH_POISSON_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
39 #include "../generic/projection.h"
40 #include "../generic/nodes.h"
41 #include "../generic/Qelements.h"
42 #include "../generic/oomph_utilities.h"
54 template<
unsigned DIM>
96 unsigned dim =
s.size();
102 for (
unsigned i = 0;
i <
dim;
i++)
121 const unsigned& nplot)
const
126 std::stringstream error_stream;
128 <<
"Poisson elements only store a single field so i must be 0 rather"
129 <<
" than " <<
i << std::endl;
136 for (
unsigned j = 0;
j < local_loop;
j++)
154 std::stringstream error_stream;
156 <<
"Poisson elements only store a single field so i must be 0 rather"
157 <<
" than " <<
i << std::endl;
163 return "Poisson solution";
170 std::ofstream& file_out,
172 const unsigned& nplot,
178 std::stringstream error_stream;
179 error_stream <<
"Poisson equation has only one field. Can't call "
180 <<
"this function for value " <<
i << std::endl;
197 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
216 const unsigned n_plot = 5;
222 void output(std::ostream& outfile,
const unsigned& n_plot);
227 const unsigned n_plot = 5;
233 void output(FILE* file_pt,
const unsigned& n_plot);
238 const unsigned& n_plot,
245 std::ostream& outfile,
246 const unsigned& n_plot,
251 "There is no time-dependent output_fct() for Poisson elements ",
275 "There is no time-dependent compute_error() for Poisson elements",
345 double source_pls = 0.0;
347 for (
unsigned i = 0;
i <
DIM;
i++)
351 gradient[
i] = (source_pls -
source) / eps_fd;
358 (*Source_fct_gradient_pt)(
x, gradient);
367 const unsigned n_node =
nnode();
380 for (
unsigned j = 0;
j <
DIM;
j++)
386 for (
unsigned l = 0; l < n_node; l++)
389 for (
unsigned j = 0;
j <
DIM;
j++)
403 const unsigned n_node =
nnode();
413 for (
unsigned i = 0;
i <
DIM;
i++)
415 for (
unsigned j = 0;
j < n_node;
j++)
417 dflux_dnodal_u[
i][
j] = dpsidx(
j,
i);
448 const unsigned n_node =
nnode();
460 double interpolated_u = 0.0;
463 for (
unsigned l = 0; l < n_node; l++)
465 interpolated_u += this->
nodal_value(l, u_nodal_index) * psi[l];
468 return (interpolated_u);
490 DShape& dtestdx)
const = 0;
500 DShape& dtestdx)
const = 0;
520 const unsigned& flag);
539 template<
unsigned DIM,
unsigned NNODE_1D>
576 void output(std::ostream& outfile,
const unsigned& n_plot)
592 void output(FILE* file_pt,
const unsigned& n_plot)
601 const unsigned& n_plot,
612 const unsigned& n_plot,
663 template<
unsigned DIM,
unsigned NNODE_1D>
672 const double J = this->dshape_eulerian(
s, psi, dpsidx);
689 template<
unsigned DIM,
unsigned NNODE_1D>
698 const double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
717 template<
unsigned DIM,
unsigned NNODE_1D>
730 const double J = this->dshape_eulerian_at_knot(
731 ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
736 d_dtestdx_dX = d_dpsidx_dX;
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>
803 std::stringstream error_stream;
804 error_stream <<
"Poisson elements only store a single field so fld "
806 <<
" than " << fld << std::endl;
813 unsigned nnod = this->
nnode();
817 for (
unsigned j = 0;
j < nnod;
j++)
820 data_values[
j] = std::make_pair(this->
node_pt(
j), fld);
840 std::stringstream error_stream;
841 error_stream <<
"Poisson elements only store a single field so fld "
843 <<
" than " << fld << std::endl;
867 std::stringstream error_stream;
868 error_stream <<
"Poisson elements only store a single field so fld "
870 <<
" than " << fld << std::endl;
875 unsigned n_dim = this->
dim();
876 unsigned n_node = this->
nnode();
878 DShape dpsidx(n_node, n_dim), dtestdx(n_node, n_dim);
880 this->dshape_and_dtest_eulerian_poisson(
s, psi, dpsidx,
test, dtestdx);
894 std::stringstream error_stream;
895 error_stream <<
"Poisson elements only store a single field so fld "
897 <<
" than " << fld << std::endl;
903 unsigned u_nodal_index = this->u_index_poisson();
906 unsigned n_node = this->
nnode();
913 double interpolated_u = 0.0;
916 for (
unsigned l = 0; l < n_node; l++)
918 interpolated_u += this->
nodal_value(l, u_nodal_index) * psi[l];
920 return interpolated_u;
930 std::stringstream error_stream;
931 error_stream <<
"Poisson elements only store a single field so fld "
933 <<
" than " << fld << std::endl;
938 return this->
nnode();
948 std::stringstream error_stream;
949 error_stream <<
"Poisson elements only store a single field so fld "
951 <<
" than " << fld << std::endl;
956 const unsigned u_nodal_index = this->u_index_poisson();
966 template<
class ELEMENT>
979 template<
class ELEMENT>
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
int data[]
Definition: Map_placement_new.cpp:1
unsigned ntstorage() const
Definition: nodes.cc:879
FaceGeometry()
Definition: poisson_elements.h:984
FaceGeometry()
Definition: poisson_elements.h:971
FaceGeometry()
Definition: poisson_elements.h:778
FaceGeometry()
Definition: poisson_elements.h:761
Definition: elements.h:4998
Definition: elements.h:1313
virtual unsigned nplot_points_paraview(const unsigned &nplot) const
Definition: elements.h:2862
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
virtual void shape(const Vector< double > &s, Shape &psi) const =0
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned dim() const
Definition: elements.h:2611
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Definition: elements.h:1765
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3298
static double Default_fd_jacobian_step
Definition: elements.h:1198
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
Definition: nodes.h:1022
Definition: oomph_definitions.h:222
Definition: elements.h:3439
Definition: poisson_elements.h:56
virtual void fill_in_generic_residual_contribution_poisson(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: poisson_elements.cc:49
void point_output_data(const Vector< double > &s, Vector< double > &data)
Definition: poisson_elements.h:93
unsigned self_test()
Self-test: Return 0 for OK.
Definition: poisson_elements.cc:298
virtual void get_source_gradient_poisson(const unsigned &ipt, const Vector< double > &x, Vector< double > &gradient) const
Definition: poisson_elements.h:331
virtual void get_source_poisson(const unsigned &ipt, const Vector< double > &x, double &source) const
Definition: poisson_elements.h:309
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
PoissonSourceFctGradientPt source_fct_gradient_pt() const
Access function: Pointer to gradient source function. Const version.
Definition: poisson_elements.h:299
void get_dflux_dnodal_u(const Vector< double > &s, Vector< Vector< double >> &dflux_dnodal_u) const
Definition: poisson_elements.h:399
void(* PoissonSourceFctPt)(const Vector< double > &x, double &f)
Definition: poisson_elements.h:60
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: poisson_elements.h:214
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
Definition: poisson_elements.h:424
void scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Definition: poisson_elements.h:119
PoissonEquations()
Constructor (must initialise the Source_fct_pt to null)
Definition: poisson_elements.h:70
PoissonSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
Definition: poisson_elements.h:287
void compute_norm(double &norm)
Compute norm of solution: square of the L2 norm.
Definition: poisson_elements.cc:447
PoissonSourceFctPt Source_fct_pt
Pointer to source function:
Definition: poisson_elements.h:523
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
Definition: poisson_elements.h:268
PoissonSourceFctGradientPt Source_fct_gradient_pt
Pointer to gradient of source function.
Definition: poisson_elements.h:526
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Definition: poisson_elements.cc:165
void(* PoissonSourceFctGradientPt)(const Vector< double > &x, Vector< double > &gradient)
Definition: poisson_elements.h:65
PoissonSourceFctGradientPt & source_fct_gradient_pt()
Access function: Pointer to gradient of source function.
Definition: poisson_elements.h:293
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: poisson_elements.h:244
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
unsigned nscalar_paraview() const
Definition: poisson_elements.h:112
virtual double interpolated_u_poisson(const Vector< double > &s) const
Definition: poisson_elements.h:445
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
Definition: poisson_elements.h:364
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: poisson_elements.cc:500
std::string scalar_name_paraview(const unsigned &i) const
Definition: poisson_elements.h:149
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: poisson_elements.h:435
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: poisson_elements.cc:400
virtual unsigned u_index_poisson() const
Definition: poisson_elements.h:85
PoissonSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: poisson_elements.h:281
virtual double dshape_and_dtest_eulerian_at_knot_poisson(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
void scalar_value_fct_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt) const
Definition: poisson_elements.h:169
void output(FILE *file_pt)
C_style output with default number of plot points.
Definition: poisson_elements.h:225
PoissonEquations(const PoissonEquations &dummy)=delete
Broken copy constructor.
Definition: projection.h:183
Poisson upgraded to become projectable.
Definition: poisson_elements.h:793
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
Definition: poisson_elements.h:925
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: poisson_elements.h:798
unsigned nfields_for_projection()
Number of fields to be projected: Just one.
Definition: poisson_elements.h:828
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: poisson_elements.h:860
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: poisson_elements.h:887
unsigned nhistory_values_for_coordinate_projection()
Definition: poisson_elements.h:853
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: poisson_elements.h:835
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: poisson_elements.h:943
Definition: Qelements.h:459
Definition: poisson_elements.h:542
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
Definition: poisson_elements.h:719
double dshape_and_dtest_eulerian_at_knot_poisson(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: poisson_elements.h:691
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: poisson_elements.h:611
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: poisson_elements.h:576
QPoissonElement(const QPoissonElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
unsigned required_nvalue(const unsigned &n) const
Definition: poisson_elements.h:561
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: poisson_elements.h:600
QPoissonElement()
Definition: poisson_elements.h:551
void output(std::ostream &outfile)
Definition: poisson_elements.h:568
double dshape_and_dtest_eulerian_poisson(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: poisson_elements.h:664
void output(FILE *file_pt)
Definition: poisson_elements.h:584
static const unsigned Initial_Nvalue
Definition: poisson_elements.h:546
void operator=(const QPoissonElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
void output(FILE *file_pt, const unsigned &n_plot)
Definition: poisson_elements.h:592
A Rank 4 Tensor class.
Definition: matrices.h:1701
A Rank 3 Tensor class.
Definition: matrices.h:1370
unsigned ntstorage() const
Definition: timesteppers.h:601
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
squared absolute sa ArrayBase::abs2 DOXCOMMA MatrixBase::cwiseAbs2 sa Eigen::abs2 DOXCOMMA Eigen::pow DOXCOMMA ArrayBase::square nearest sa Eigen::floor DOXCOMMA Eigen::ceil DOXCOMMA ArrayBase::round nearest integer not less than the given sa Eigen::floor DOXCOMMA ArrayBase::ceil not a number test
Definition: GlobalFunctions.h:109
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
void exact_soln(const double &time, const Vector< double > &x, Vector< double > &soln)
Definition: unstructured_two_d_curved.cc:301
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
int error
Definition: calibrate.py:297
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
list x
Definition: plotDoE.py:28
t
Definition: plotPSD.py:36
Definition: indexed_view.cpp:20
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2