30 #ifndef OOMPH_UNSTEADY_HEAT_FLUX_ELEMENTS_HEADER
31 #define OOMPH_UNSTEADY_HEAT_FLUX_ELEMENTS_HEADER
35 #include <oomph-lib-config.h>
42 #include "../generic/Qelements.h"
53 template<
class ELEMENT>
103 residuals, jacobian, 1);
114 const unsigned&
i)
const
128 void output(std::ostream& outfile,
const unsigned& n_plot)
143 void output(FILE* file_pt,
const unsigned& n_plot)
157 unsigned n_node =
nnode();
163 for (
unsigned i = 0;
i < n_node;
i++)
184 (*Flux_fct_pt)(time,
x,
flux);
214 template<
class ELEMENT>
227 ELEMENT* elem_pt =
dynamic_cast<ELEMENT*
>(bulk_el_pt);
230 if (elem_pt->dim() == 3)
239 throw OomphLibError(
"This flux element will not work correctly if "
240 "nodes are hanging\n",
276 "Bulk element must inherit from UnsteadyHeatEquations.";
278 "Nodes are one dimensional, but cannot cast the bulk element to\n";
279 error_string +=
"UnsteadyHeatEquations<1>\n.";
280 error_string +=
"If you desire this functionality, you must "
281 "implement it yourself\n";
304 "Bulk element must inherit from UnsteadyHeatEquations.";
306 "Nodes are two dimensional, but cannot cast the bulk element to\n";
307 error_string +=
"UnsteadyHeatEquations<2>\n.";
308 error_string +=
"If you desire this functionality, you must "
309 "implement it yourself\n";
331 "Bulk element must inherit from UnsteadyHeatEquations.";
332 error_string +=
"Nodes are three dimensional, but cannot cast the "
334 error_string +=
"UnsteadyHeatEquations<3>\n.";
335 error_string +=
"If you desire this functionality, you must "
336 "implement it yourself\n";
351 std::ostringstream error_stream;
352 error_stream <<
"Dimension of node is " <<
Dim
353 <<
". It should be 1,2, or 3!" << std::endl;
365 template<
class ELEMENT>
371 const unsigned n_node = nnode();
374 double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
377 Shape psif(n_node), testf(n_node);
380 const unsigned n_intpt = integral_pt()->nweight();
389 const unsigned u_index_ust_heat = U_index_ust_heat;
393 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
396 for (
unsigned i = 0;
i < (
Dim - 1);
i++)
398 s[
i] = integral_pt()->knot(ipt,
i);
402 double w = integral_pt()->weight(ipt);
406 double J = shape_and_test(
s, psif, testf);
415 for (
unsigned i = 0;
i <
Dim;
i++)
417 interpolated_x[
i] = 0.0;
421 for (
unsigned l = 0; l < n_node; l++)
424 for (
unsigned i = 0;
i <
Dim;
i++)
426 interpolated_x[
i] += nodal_position(l,
i) * psif[l];
437 for (
unsigned l = 0; l < n_node; l++)
439 local_eqn = nodal_local_eqn(l, u_index_ust_heat);
444 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 J_eulerian(const Vector< double > &s) const
Definition: elements.cc:5242
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 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
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
Definition: unsteady_heat_elements.h:72
virtual unsigned u_index_ust_heat() const
Broken assignment operator.
Definition: unsteady_heat_elements.h:112
Definition: unsteady_heat_flux_elements.h:56
double shape_and_test(const Vector< double > &s, Shape &psi, Shape &test) const
Definition: unsteady_heat_flux_elements.h:152
unsigned Dim
The spatial dimension of the problem.
Definition: unsteady_heat_flux_elements.h:200
void get_flux(const double &time, const Vector< double > &x, double &flux)
Definition: unsteady_heat_flux_elements.h:174
UnsteadyHeatPrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
Definition: unsteady_heat_flux_elements.h:197
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Definition: unsteady_heat_flux_elements.h:112
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute the element residual vector.
Definition: unsteady_heat_flux_elements.h:88
UnsteadyHeatPrescribedFluxFctPt & flux_fct_pt()
Broken assignment operator.
Definition: unsteady_heat_flux_elements.h:82
void output(FILE *file_pt)
Definition: unsteady_heat_flux_elements.h:135
UnsteadyHeatFluxElement(FiniteElement *const &bulk_el_pt, const int &face_index)
Definition: unsteady_heat_flux_elements.h:215
unsigned U_index_ust_heat
The index at which the unknown is stored at the nodes.
Definition: unsteady_heat_flux_elements.h:203
void output(FILE *file_pt, const unsigned &n_plot)
Definition: unsteady_heat_flux_elements.h:143
void(* UnsteadyHeatPrescribedFluxFctPt)(const double &time, const Vector< double > &x, double &flux)
Definition: unsteady_heat_flux_elements.h:60
UnsteadyHeatFluxElement(const UnsteadyHeatFluxElement &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile)
Definition: unsteady_heat_flux_elements.h:121
void fill_in_generic_residual_contribution_ust_heat_flux(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Compute the element's residual vector and the (zero) Jacobian matrix.
Definition: unsteady_heat_flux_elements.h:367
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: unsteady_heat_flux_elements.h:128
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual vector and its Jacobian matrix.
Definition: unsteady_heat_flux_elements.h:98
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
void output(std::ostream &outfile, const unsigned &nplot)
Overload output function.
Definition: overloaded_element_body.h:490