28 #ifndef OOMPH_PML_FOURIER_DECOMPOSED_HELMHOLTZ_FLUX_ELEMENTS_HEADER
29 #define OOMPH_PML_FOURIER_DECOMPOSED_HELMHOLTZ_FLUX_ELEMENTS_HEADER
34 #include <oomph-lib-config.h>
41 #include "../generic/Qelements.h"
52 template<
class ELEMENT>
73 "PMLFourierDecomposedHelmholtzFluxElement",
74 OOMPH_CURRENT_FUNCTION,
75 OOMPH_EXCEPTION_LOCATION);
116 residuals, jacobian, 1);
128 void output(std::ostream& outfile,
const unsigned& n_plot)
144 void output(FILE* file_pt,
const unsigned& n_plot)
155 return std::complex<unsigned>(
170 unsigned n_node =
nnode();
176 for (
unsigned i = 0;
i < n_node;
i++)
194 unsigned n_node =
nnode();
200 for (
unsigned i = 0;
i < n_node;
i++)
217 flux = std::complex<double>(0.0, 0.0);
222 (*Flux_fct_pt)(x, flux);
238 const unsigned& flag);
257 template<
class ELEMENT>
260 const int& face_index)
280 std::string error_string =
"Bulk element must inherit from "
281 "PMLFourierDecomposedHelmholtzEquations.";
283 error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
297 template<
class ELEMENT>
302 const unsigned& flag)
305 const unsigned n_node = nnode();
308 Shape psif(n_node), testf(n_node);
311 const unsigned n_intpt = integral_pt()->nweight();
317 int local_eqn_real = 0, local_eqn_imag = 0;
321 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
324 for (
unsigned i = 0;
i < 1;
i++)
326 s[
i] = integral_pt()->knot(ipt,
i);
330 double w = integral_pt()->weight(ipt);
334 double J = shape_and_test(
s, psif, testf);
343 for (
unsigned l = 0; l < n_node; l++)
346 for (
unsigned i = 0;
i < 2;
i++)
348 interpolated_x[
i] += nodal_position(l,
i) * psif[l];
353 double r = interpolated_x[0];
356 std::complex<double> flux(0.0, 0.0);
357 get_flux(interpolated_x, flux);
361 for (
unsigned l = 0; l < n_node; l++)
364 nodal_local_eqn(l, U_index_pml_fourier_decomposed_helmholtz.real());
367 if (local_eqn_real >= 0)
370 residuals[local_eqn_real] -= flux.real() * testf[l] * r *
W;
377 nodal_local_eqn(l, U_index_pml_fourier_decomposed_helmholtz.imag());
380 if (local_eqn_imag >= 0)
383 residuals[local_eqn_imag] -= flux.imag() * testf[l] * r *
W;
406 template<
class ELEMENT>
421 "PMLFourierDecomposedHelmholtzPowerMonitorElement",
422 OOMPH_CURRENT_FUNCTION,
423 OOMPH_EXCEPTION_LOCATION);
442 const unsigned&
i)
const
457 void output(std::ostream& outfile,
const unsigned& n_plot)
472 void output(FILE* file_pt,
const unsigned& n_plot)
482 return std::complex<unsigned>(
494 std::ofstream outfile;
506 ELEMENT* bulk_elem_pt =
dynamic_cast<ELEMENT*
>(this->
bulk_element_pt());
509 unsigned nnode_bulk = bulk_elem_pt->nnode();
510 const unsigned n_node_local = this->
nnode();
513 const unsigned bulk_dim = bulk_elem_pt->dim();
514 const unsigned local_dim = this->
node_pt(0)->
ndim();
517 Shape psi(n_node_local);
520 Shape psi_bulk(nnode_bulk);
521 DShape dpsi_bulk_dx(nnode_bulk, bulk_dim);
534 if (outfile.is_open())
541 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
544 for (
unsigned i = 0;
i < (local_dim - 1);
i++)
566 (void)bulk_elem_pt->dshape_eulerian(s_bulk, psi_bulk, dpsi_bulk_dx);
570 std::complex<double> dphi_dn(0.0, 0.0);
572 std::complex<double> interpolated_phi(0.0, 0.0);
578 for (
unsigned l = 0; l < nnode_bulk; l++)
581 const std::complex<double> phi_value(
582 bulk_elem_pt->nodal_value(
584 bulk_elem_pt->u_index_pml_fourier_decomposed_helmholtz().real()),
585 bulk_elem_pt->nodal_value(
587 bulk_elem_pt->u_index_pml_fourier_decomposed_helmholtz().imag()));
590 for (
unsigned i = 0;
i < bulk_dim;
i++)
592 interpolated_dphidx[
i] += phi_value * dpsi_bulk_dx(l,
i);
596 for (
unsigned l = 0; l < n_node_local; l++)
599 const std::complex<double> phi_value(
603 interpolated_phi += phi_value * psi(l);
607 for (
unsigned i = 0;
i < bulk_dim;
i++)
609 dphi_dn += interpolated_dphidx[
i] * unit_normal[
i];
613 double integrand = (interpolated_phi.real() * dphi_dn.imag() -
614 interpolated_phi.imag() * dphi_dn.real());
617 double theta = atan2(x[0], x[1]);
619 if (outfile.is_open())
621 outfile << x[0] <<
" " << x[1] <<
" " << theta <<
" " << integrand
643 unsigned nnod =
nnode();
644 for (
unsigned i = 0;
i < nnod;
i++)
664 unsigned n_node =
nnode();
670 for (
unsigned i = 0;
i < n_node;
i++)
672 for (
unsigned j = 0; j < 1; j++)
675 dtest_ds(
i, j) = dpsi_ds(
i, j);
696 template<
class ELEMENT>
719 std::string error_string =
"Bulk element must inherit from "
720 "PMLFourierDecomposedHelmholtzEquations.";
722 error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
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.
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
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...
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.
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...
unsigned nnode() const
Return the number of nodes.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Function to compute the geometric shape functions and derivatives w.r.t. local coordinates at local c...
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.
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....
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
virtual std::complex< unsigned > u_index_pml_fourier_decomposed_helmholtz() const
Broken assignment operator.
A class for elements that allow the imposition of an applied flux on the boundaries of Fourier decomp...
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
virtual std::complex< unsigned > u_index_pml_fourier_decomposed_helmholtz() const
Return the index at which the unknown value is stored. (Real/imag part gives real index of real/imag ...
void(* PMLFourierDecomposedHelmholtzPrescribedFluxFctPt)(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 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...
std::complex< unsigned > U_index_pml_fourier_decomposed_helmholtz
The index at which the real and imag part of the unknown is stored at the nodes.
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 ...
PMLFourierDecomposedHelmholtzFluxElement(const PMLFourierDecomposedHelmholtzFluxElement &dummy)=delete
Broken copy constructor.
PMLFourierDecomposedHelmholtzFluxElement()
Broken empty constructor.
virtual void fill_in_generic_residual_contribution_pml_fourier_decomposed_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(std::ostream &outfile)
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 ...
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 get_flux(const Vector< double > &x, std::complex< double > &flux)
Function to calculate the prescribed flux at a given spatial position.
PMLFourierDecomposedHelmholtzPrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
PMLFourierDecomposedHelmholtzPrescribedFluxFctPt & flux_fct_pt()
Broken assignment operator.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
////////////////////////////////////////////////////////////////// //////////////////////////////////...
PMLFourierDecomposedHelmholtzPowerMonitorElement()
Broken empty constructor.
std::complex< unsigned > U_index_pml_fourier_decomposed_helmholtz
The index at which the real and imag part of the unknown is stored at the nodes.
void output(FILE *file_pt)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
double d_shape_and_test_local(const Vector< double > &s, Shape &psi, Shape &test, DShape &dpsi_ds, DShape &dtest_ds) const
Function to compute the shape, test functions and derivates and to return the Jacobian of mapping bet...
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Broken assignment operator.
void output(FILE *file_pt, const unsigned &n_plot)
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...
double global_power_contribution(std::ofstream &outfile)
Compute the element's contribution to the time-averaged radiated power over the artificial boundary....
double shape_and_test(const Vector< double > &s, Shape &psi, Shape &test) const
Function to compute the test functions and to return the Jacobian of mapping between local and global...
virtual std::complex< unsigned > u_index_pml_fourier_decomposed_helmholtz() const
Return the index at which the real/imag unknown value is stored.
void output(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
double global_power_contribution()
Compute the element's contribution to the time-averaged radiated power over the artificial boundary....
PMLFourierDecomposedHelmholtzPowerMonitorElement(const PMLFourierDecomposedHelmholtzPowerMonitorElement &dummy)=delete
Broken copy constructor.
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.
const double Pi
50 digits from maple
//////////////////////////////////////////////////////////////////// ////////////////////////////////...