28 #ifndef OOMPH_HELMHOLTZ_FLUX_ELEMENTS_HEADER
29 #define OOMPH_HELMHOLTZ_FLUX_ELEMENTS_HEADER
34 #include <oomph-lib-config.h>
38 #include "../generic/Qelements.h"
50 template<
class ELEMENT>
59 std::complex<double>& flux);
70 "Don't call empty constructor for HelmholtzFluxElement",
71 OOMPH_CURRENT_FUNCTION,
72 OOMPH_EXCEPTION_LOCATION);
111 residuals, jacobian, 1);
122 const unsigned&
i)
const
137 void output(std::ostream& outfile,
const unsigned& n_plot)
153 void output(FILE* file_pt,
const unsigned& n_plot)
177 unsigned n_node =
nnode();
183 for (
unsigned i = 0;
i < n_node;
i++)
201 unsigned n_node =
nnode();
207 for (
unsigned i = 0;
i < n_node;
i++)
224 flux = std::complex<double>(0.0, 0.0);
229 (*Flux_fct_pt)(x, flux);
245 const unsigned& flag);
267 template<
class ELEMENT>
280 ELEMENT* elem_pt =
dynamic_cast<ELEMENT*
>(bulk_el_pt);
282 if (elem_pt->dim() == 3)
291 throw OomphLibError(
"This flux element will not work correctly if "
292 "nodes are hanging\n",
293 OOMPH_CURRENT_FUNCTION,
294 OOMPH_EXCEPTION_LOCATION);
328 "Bulk element must inherit from HelmholtzEquations.";
330 "Nodes are one dimensional, but cannot cast the bulk element to\n";
331 error_string +=
"HelmholtzEquations<1>\n.";
332 error_string +=
"If you desire this functionality, you must "
333 "implement it yourself\n";
336 error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
356 "Bulk element must inherit from HelmholtzEquations.";
358 "Nodes are two dimensional, but cannot cast the bulk element to\n";
359 error_string +=
"HelmholtzEquations<2>\n.";
360 error_string +=
"If you desire this functionality, you must "
361 "implement it yourself\n";
364 error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
384 "Bulk element must inherit from HelmholtzEquations.";
385 error_string +=
"Nodes are three dimensional, but cannot cast the "
387 error_string +=
"HelmholtzEquations<3>\n.";
388 error_string +=
"If you desire this functionality, you must "
389 "implement it yourself\n";
392 error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
404 std::ostringstream error_stream;
405 error_stream <<
"Dimension of node is " <<
Dim
406 <<
". It should be 1,2, or 3!" << std::endl;
409 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
418 template<
class ELEMENT>
423 const unsigned& flag)
426 const unsigned n_node = nnode();
429 Shape psif(n_node), testf(n_node);
432 const unsigned n_intpt = integral_pt()->nweight();
438 int local_eqn_real = 0, local_eqn_imag = 0;
442 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
445 for (
unsigned i = 0;
i < (
Dim - 1);
i++)
447 s[
i] = integral_pt()->knot(ipt,
i);
451 double w = integral_pt()->weight(ipt);
455 double J = shape_and_test(
s, psif, testf);
464 for (
unsigned l = 0; l < n_node; l++)
466 for (
unsigned i = 0;
i <
Dim;
i++)
468 interpolated_x[
i] += nodal_position(l,
i) * psif[l];
473 std::complex<double> flux(0.0, 0.0);
474 get_flux(interpolated_x, flux);
478 for (
unsigned l = 0; l < n_node; l++)
480 local_eqn_real = nodal_local_eqn(l, U_index_helmholtz.real());
482 if (local_eqn_real >= 0)
485 residuals[local_eqn_real] -= flux.real() * testf[l] *
W;
491 local_eqn_imag = nodal_local_eqn(l, U_index_helmholtz.imag());
493 if (local_eqn_imag >= 0)
496 residuals[local_eqn_imag] -= flux.imag() * testf[l] *
W;
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
In a FaceElement, the "global" intrinsic coordinate of the element along the boundary,...
double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s....
double J_eulerian_at_knot(const unsigned &ipt) const
Return the Jacobian of the mapping from local to global coordinates at the ipt-th integration point O...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
A general Finite Element class.
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
unsigned nnode() const
Return the number of nodes.
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Function for building a lower dimensional FaceElement on the specified face of the FiniteElement....
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Return the geometric shape function at the ipt-th integration point.
bool has_hanging_nodes() const
Return boolean to indicate if any of the element's nodes are geometrically hanging.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
A class for all isoparametric elements that solve the Helmholtz equations.
virtual std::complex< unsigned > u_index_helmholtz() const
Broken assignment operator.
A class for elements that allow the imposition of an applied flux on the boundaries of Helmholtz elem...
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Specify the value of nodal zeta from the face geometry The "global" intrinsic coordinate of the eleme...
void output(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
void get_flux(const Vector< double > &x, std::complex< double > &flux)
Function to calculate the prescribed flux at a given spatial position.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
HelmholtzPrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
double shape_and_test_at_knot(const unsigned &ipt, Shape &psi, Shape &test) const
Function to compute the shape and test functions and to return the Jacobian of mapping between local ...
virtual void fill_in_generic_residual_contribution_helmholtz_flux(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Add the element's contribution to its residual vector. flag=1(or 0): do (or don't) compute the contri...
unsigned Dim
The spatial dimension of the problem.
HelmholtzFluxElement()
Broken empty constructor.
HelmholtzFluxElement(const HelmholtzFluxElement &dummy)=delete
Broken copy constructor.
double shape_and_test(const Vector< double > &s, Shape &psi, Shape &test) const
Function to compute the shape and test functions and to return the Jacobian of mapping between local ...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector.
void output(FILE *file_pt)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
HelmholtzPrescribedFluxFctPt & flux_fct_pt()
Broken assignment operator.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the element's contribution to its residual vector and its Jacobian matrix.
std::complex< unsigned > U_index_helmholtz
The index at which the real and imag part of the unknown is stored at the nodes.
void(* HelmholtzPrescribedFluxFctPt)(const Vector< double > &x, std::complex< double > &flux)
Function pointer to the prescribed-flux function fct(x,f(x)) – x is a Vector and the flux is a comple...
virtual std::complex< unsigned > u_index_helmholtz() const
Return the index at which the unknown value is stored.
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
An OomphLibError object which should be thrown when an run-time error is encountered....
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
unsigned Dim
Dimension of zeta tuples (set by get_dim_helper) – needed because we store the scalar coordinates in ...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...