28 #ifndef OOMPH_ADV_DIFF_FLUX_ELEMENTS_HEADER
29 #define OOMPH_ADV_DIFF_FLUX_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
37 #include "../generic/Qelements.h"
52 template<
class ELEMENT>
73 "Don't call empty constructor for AdvectionDiffusionFluxElement",
109 residuals, jacobian, 1);
119 const unsigned&
i)
const
133 void output(std::ostream& outfile,
const unsigned& nplot)
148 unsigned n_node =
nnode();
154 for (
unsigned i = 0;
i < n_node;
i++)
172 unsigned n_node =
nnode();
178 for (
unsigned i = 0;
i < n_node;
i++)
200 (*Flux_fct_pt)(
x,
flux);
231 template<
class ELEMENT>
244 ELEMENT* elem_pt =
dynamic_cast<ELEMENT*
>(bulk_el_pt);
246 if (elem_pt->dim() == 3)
255 throw OomphLibError(
"This flux element will not work correctly if "
256 "nodes are hanging\n",
293 "Bulk element must inherit from AdvectionDiffusionEquations.";
295 "Nodes are one dimensional, but cannot cast the bulk element to\n";
296 error_string +=
"AdvectionDiffusionEquations<1>\n.";
297 error_string +=
"If you desire this functionality, you must "
298 "implement it yourself\n";
321 "Bulk element must inherit from AdvectionDiffusionEquations.";
323 "Nodes are two dimensional, but cannot cast the bulk element to\n";
324 error_string +=
"AdvectionDiffusionEquations<2>\n.";
325 error_string +=
"If you desire this functionality, you must "
326 "implement it yourself\n";
348 "Bulk element must inherit from AdvectionDiffusionEquations.";
349 error_string +=
"Nodes are three dimensional, but cannot cast the "
351 error_string +=
"AdvectionDiffusionEquations<3>\n.";
352 error_string +=
"If you desire this functionality, you must "
353 "implement it yourself\n";
368 std::ostringstream error_stream;
369 error_stream <<
"Dimension of node is " <<
Dim
370 <<
". It should be 1,2, or 3!" << std::endl;
382 template<
class ELEMENT>
388 const unsigned n_node = nnode();
391 Shape psif(n_node), testf(n_node);
394 const unsigned n_intpt = integral_pt()->nweight();
404 const unsigned u_index_adv_diff = U_index_adv_diff;
409 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
412 for (
unsigned i = 0;
i < (
Dim - 1);
i++)
414 s[
i] = integral_pt()->knot(ipt,
i);
418 double w = integral_pt()->weight(ipt);
422 double J = shape_and_test(
s, psif, testf);
431 for (
unsigned l = 0; l < n_node; l++)
434 for (
unsigned i = 0;
i <
Dim;
i++)
436 interpolated_x[
i] += nodal_position(l,
i) * psif[l];
447 for (
unsigned l = 0; l < n_node; l++)
450 local_eqn = nodal_local_eqn(l, u_index_adv_diff);
455 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: advection_diffusion_elements.h:52
virtual unsigned u_index_adv_diff() const
Definition: advection_diffusion_elements.h:91
Definition: advection_diffusion_flux_elements.h:55
void output(std::ostream &outfile)
Definition: advection_diffusion_flux_elements.h:126
double shape_and_test_at_knot(const unsigned &ipt, Shape &psi, Shape &test) const
Definition: advection_diffusion_flux_elements.h:167
AdvectionDiffusionFluxElement(const AdvectionDiffusionFluxElement &dummy)=delete
Broken copy constructor.
void operator=(const AdvectionDiffusionFluxElement &)=delete
Broken assignment operator.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: advection_diffusion_flux_elements.h:104
void output(std::ostream &outfile, const unsigned &nplot)
Definition: advection_diffusion_flux_elements.h:133
unsigned U_index_adv_diff
The index at which the unknown is stored at the nodes.
Definition: advection_diffusion_flux_elements.h:218
double shape_and_test(const Vector< double > &s, Shape &psi, Shape &test) const
Definition: advection_diffusion_flux_elements.h:143
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector.
Definition: advection_diffusion_flux_elements.h:93
AdvectionDiffusionPrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
Definition: advection_diffusion_flux_elements.h:212
AdvectionDiffusionPrescribedFluxFctPt & flux_fct_pt()
Access function for the prescribed-flux function pointer.
Definition: advection_diffusion_flux_elements.h:86
unsigned Dim
The spatial dimension of the problem.
Definition: advection_diffusion_flux_elements.h:215
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: advection_diffusion_flux_elements.h:117
void get_flux(const Vector< double > &x, double &flux)
Definition: advection_diffusion_flux_elements.h:190
void fill_in_generic_residual_contribution_adv_diff_flux(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Compute the element's residual vector and the (zero) Jacobian matrix.
Definition: advection_diffusion_flux_elements.h:384
void(* AdvectionDiffusionPrescribedFluxFctPt)(const Vector< double > &x, double &flux)
Definition: advection_diffusion_flux_elements.h:59
AdvectionDiffusionFluxElement()
Broken empty constructor.
Definition: advection_diffusion_flux_elements.h:70
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 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 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 nnode() const
Return the number of nodes.
Definition: elements.h:2210
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
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