28 #ifndef OOMPH_POISSON_FLUX_ELEMENTS_HEADER
29 #define OOMPH_POISSON_FLUX_ELEMENTS_HEADER
34 #include <oomph-lib-config.h>
38 #include "../generic/Qelements.h"
48 template<
class ELEMENT>
65 throw OomphLibError(
"Don't call empty constructor for PoissonFluxElement",
83 const unsigned&
i)
const
112 residuals, jacobian, 1);
118 const unsigned n_plot = 5;
123 void output(std::ostream& outfile,
const unsigned& nplot)
136 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
142 for (
unsigned i = 0;
i <
el_dim + 1;
i++)
145 outfile <<
x[
i] <<
" ";
149 outfile <<
flux << std::endl;
167 void output(FILE* file_pt,
const unsigned& n_plot)
182 unsigned n_node =
nnode();
188 for (
unsigned i = 0;
i < n_node;
i++)
206 unsigned n_node =
nnode();
212 for (
unsigned i = 0;
i < n_node;
i++)
234 (*Flux_fct_pt)(
x,
flux);
245 const unsigned& flag);
270 template<
class ELEMENT>
283 ELEMENT* elem_pt =
dynamic_cast<ELEMENT*
>(bulk_el_pt);
285 if (elem_pt->dim() == 3)
294 throw OomphLibError(
"This flux element will not work correctly if "
295 "nodes are hanging\n",
331 "Bulk element must inherit from PoissonEquations.";
333 "Nodes are one dimensional, but cannot cast the bulk element to\n";
334 error_string +=
"PoissonEquations<1>\n.";
335 error_string +=
"If you desire this functionality, you must "
336 "implement it yourself\n";
359 "Bulk element must inherit from PoissonEquations.";
361 "Nodes are two dimensional, but cannot cast the bulk element to\n";
362 error_string +=
"PoissonEquations<2>\n.";
363 error_string +=
"If you desire this functionality, you must "
364 "implement it yourself\n";
386 "Bulk element must inherit from PoissonEquations.";
387 error_string +=
"Nodes are three dimensional, but cannot cast the "
389 error_string +=
"PoissonEquations<3>\n.";
390 error_string +=
"If you desire this functionality, you must "
391 "implement it yourself\n";
406 std::ostringstream error_stream;
407 error_stream <<
"Dimension of node is " <<
Dim
408 <<
". It should be 1,2, or 3!" << std::endl;
420 template<
class ELEMENT>
425 const unsigned& flag)
428 const unsigned n_node = nnode();
431 Shape psif(n_node), testf(n_node);
434 const unsigned n_intpt = integral_pt()->nweight();
443 const unsigned u_index_poisson = U_index_poisson;
447 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
450 for (
unsigned i = 0;
i < (
Dim - 1);
i++)
452 s[
i] = integral_pt()->knot(ipt,
i);
456 double w = integral_pt()->weight(ipt);
460 double J = shape_and_test(
s, psif, testf);
469 for (
unsigned l = 0; l < n_node; l++)
472 for (
unsigned i = 0;
i <
Dim;
i++)
474 interpolated_x[
i] += nodal_position(l,
i) * psif[l];
485 for (
unsigned l = 0; l < n_node; l++)
487 local_eqn = nodal_local_eqn(l, u_index_poisson);
492 residuals[local_eqn] -=
flux * testf[l] *
W;
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
RowVector3d w
Definition: Matrix_resize_int.cpp:3
Definition: elements.h:4338
int & face_index()
Definition: elements.h:4626
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: elements.h:4497
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Definition: elements.h:4528
double J_eulerian(const Vector< double > &s) const
Definition: elements.cc:5242
double J_eulerian_at_knot(const unsigned &ipt) const
Definition: elements.cc:5328
Definition: elements.h:4998
Definition: elements.h:1313
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
virtual void output(std::ostream &outfile)
Definition: elements.h:3050
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Definition: elements.cc:5132
virtual void shape(const Vector< double > &s, Shape &psi) const =0
unsigned dim() const
Definition: elements.h:2611
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
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
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Definition: elements.h:3174
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Definition: elements.cc:3220
bool has_hanging_nodes() const
Definition: elements.h:2470
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
Definition: oomph_definitions.h:222
virtual unsigned u_index_poisson() const
Definition: poisson_elements.h:85
Definition: poisson_flux_elements.h:51
double shape_and_test(const Vector< double > &s, Shape &psi, Shape &test) const
Definition: poisson_flux_elements.h:177
void output(FILE *file_pt)
Definition: poisson_flux_elements.h:159
PoissonFluxElement()
Broken empty constructor.
Definition: poisson_flux_elements.h:63
PoissonFluxElement(const PoissonFluxElement &dummy)=delete
Broken copy constructor.
unsigned Dim
The spatial dimension of the problem.
Definition: poisson_flux_elements.h:252
void operator=(const PoissonFluxElement &)=delete
Broken assignment operator.
double shape_and_test_at_knot(const unsigned &ipt, Shape &psi, Shape &test) const
Definition: poisson_flux_elements.h:201
void get_flux(const Vector< double > &x, double &flux)
Definition: poisson_flux_elements.h:224
void(* PoissonPrescribedFluxFctPt)(const Vector< double > &x, double &flux)
Definition: poisson_flux_elements.h:55
void fill_in_generic_residual_contribution_poisson_flux(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Compute the element's residual vector and the (zero) Jacobian matrix.
Definition: poisson_flux_elements.h:422
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: poisson_flux_elements.h:107
void output(std::ostream &outfile, const unsigned &nplot)
Output function.
Definition: poisson_flux_elements.h:123
PoissonPrescribedFluxFctPt & flux_fct_pt()
Access function for the prescribed-flux function pointer.
Definition: poisson_flux_elements.h:89
unsigned U_index_poisson
The index at which the unknown is stored at the nodes.
Definition: poisson_flux_elements.h:255
void output(std::ostream &outfile)
Output function.
Definition: poisson_flux_elements.h:116
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector.
Definition: poisson_flux_elements.h:96
void output(FILE *file_pt, const unsigned &n_plot)
Definition: poisson_flux_elements.h:167
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: poisson_flux_elements.h:81
PoissonPrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
Definition: poisson_flux_elements.h:249
Definition: refineable_elements.h:97
RealScalar s
Definition: level1_cplx_impl.h:130
char char char int int * k
Definition: level2_impl.h:374
static const unsigned Dim
Problem dimension.
Definition: two_d_tilted_square.cc:62
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
@ W
Definition: quadtree.h:63
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
Definition: gen_axisym_advection_diffusion_elements.h:424
list x
Definition: plotDoE.py:28
Definition: indexed_view.cpp:20
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
unsigned el_dim
dimension
Definition: overloaded_cartesian_element_body.h:30