28#ifndef OOMPH_PML_HELMHOLTZ_FLUX_ELEMENTS_HEADER
29#define OOMPH_PML_HELMHOLTZ_FLUX_ELEMENTS_HEADER
34#include <oomph-lib-config.h>
50 template<
class ELEMENT>
64 "Don't call empty constructor for PMLHelmholtzPowerElement",
88 const unsigned&
i)
const
183 std::complex<double>
dphi_dn(0.0, 0.0);
185 bulk_dim, std::complex<double>(0.0, 0.0));
286 std::complex<double> flux(0.0, 0.0);
325 std::complex<double>
dphi_dn(0.0, 0.0);
327 bulk_dim, std::complex<double>(0.0, 0.0));
393 template<
class ELEMENT>
412 throw OomphLibError(
"This flux element will not work correctly if "
413 "nodes are hanging\n",
451 "Bulk element must inherit from PMLHelmholtzEquations.";
453 "Nodes are one dimensional, but cannot cast the bulk element to\n";
455 error_string +=
"If you desire this functionality, you must "
456 "implement it yourself\n";
479 "Bulk element must inherit from PMLHelmholtzEquations.";
481 "Nodes are two dimensional, but cannot cast the bulk element to\n";
483 error_string +=
"If you desire this functionality, you must "
484 "implement it yourself\n";
507 "Bulk element must inherit from PMLHelmholtzEquations.";
508 error_string +=
"Nodes are three dimensional, but cannot cast the "
511 error_string +=
"If you desire this functionality, you must "
512 "implement it yourself\n";
529 <<
". It should be 1,2, or 3!" << std::endl;
549 template<
class ELEMENT>
558 std::complex<double>& flux);
569 "Don't call empty constructor for PMLHelmholtzFluxElement",
616 const unsigned&
i)
const
692 if (local_eqn_number >= 0)
708 if (local_eqn_number >= 0)
777 flux = std::complex<double>(0.0, 0.0);
782 (*Flux_fct_pt)(x, flux);
798 const unsigned&
flag);
820 template<
class ELEMENT>
839 throw OomphLibError(
"This flux element will not work correctly if "
840 "nodes are hanging\n",
881 "Bulk element must inherit from PMLHelmholtzEquations.";
883 "Nodes are one dimensional, but cannot cast the bulk element to\n";
885 error_string +=
"If you desire this functionality, you must "
886 "implement it yourself\n";
909 "Bulk element must inherit from PMLHelmholtzEquations.";
911 "Nodes are two dimensional, but cannot cast the bulk element to\n";
913 error_string +=
"If you desire this functionality, you must "
914 "implement it yourself\n";
937 "Bulk element must inherit from PMLHelmholtzEquations.";
938 error_string +=
"Nodes are three dimensional, but cannot cast the "
941 error_string +=
"If you desire this functionality, you must "
942 "implement it yourself\n";
959 <<
". It should be 1,2, or 3!" << std::endl;
971 template<
class ELEMENT>
976 const unsigned&
flag)
998 for (
unsigned i = 0;
i < (Dim - 1);
i++)
1019 for (
unsigned i = 0;
i < Dim;
i++)
1026 std::complex<double> flux(0.0, 0.0);
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
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)
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
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,...
Vector< double > local_coordinate_in_bulk(const Vector< double > &s) const
Return vector of local coordinates in bulk element, given the local coordinates in this FaceElement.
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s. Overloaded to get information from bulk...
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
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...
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
A general Finite Element class.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
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...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
unsigned nnode() const
Return the number of nodes.
double nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. If the node is hanging, the appropriate interpolation is ...
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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.
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number.
int local_eqn_number(const unsigned long &ieqn_global) const
Return the local equation number corresponding to the ieqn_global-th global equation number....
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
An OomphLibError object which should be thrown when an run-time error is encountered....
A class for elements that allow the imposition of an applied flux on the boundaries of PMLHelmholtz e...
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.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
PMLHelmholtzPrescribedFluxFctPt & flux_fct_pt()
Broken assignment operator.
std::complex< unsigned > U_index_helmholtz
The index at which the real and imag part of the unknown is stored at the nodes.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
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 ...
virtual std::complex< unsigned > u_index_helmholtz() const
Return the index at which the unknown value is stored.
PMLHelmholtzFluxElement(const PMLHelmholtzFluxElement &dummy)=delete
Broken copy constructor.
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...
unsigned Dim
The spatial dimension of the problem.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: real and imag...
PMLHelmholtzPrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned > > &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
void get_flux(const Vector< double > &x, std::complex< double > &flux)
Function to calculate the prescribed flux at a given spatial position.
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...
void output(FILE *file_pt)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
PMLHelmholtzFluxElement()
Broken empty constructor.
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 ...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector.
void(* PMLHelmholtzPrescribedFluxFctPt)(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...
void output(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
A class for elements that allow the post-processing of radiated power and flux on the boundaries of P...
std::complex< double > global_flux_contribution(std::ofstream &outfile)
Compute the element's contribution to the integral of the flux over the artificial boundary....
double global_power_contribution(std::ofstream &outfile)
Compute the element's contribution to the time-averaged radiated power over the artificial boundary....
std::complex< unsigned > U_index_helmholtz
The index at which the real and imag part of the unknown is stored at the nodes.
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Broken assignment operator.
double global_power_contribution()
Compute the element's contribution to the time-averaged radiated power over the artificial boundary.
unsigned Dim
The spatial dimension of the problem.
PMLHelmholtzPowerElement()
Broken empty constructor.
virtual std::complex< unsigned > u_index_helmholtz() const
Return the index at which the unknown value is stored.
PMLHelmholtzPowerElement(const PMLHelmholtzPowerElement &dummy)=delete
Broken copy constructor.
std::complex< double > global_flux_contribution()
Compute the element's contribution to the time-averaged flux.
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...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
TAdvectionDiffusionReactionElement()
Constructor: Call constructors for TElement and AdvectionDiffusionReaction equations.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).