27 #ifndef OOMPH_HELMHOLTZ_ELEMENTS_HEADER
28 #define OOMPH_HELMHOLTZ_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
38 #include "../generic/projection.h"
39 #include "../generic/nodes.h"
40 #include "../generic/Qelements.h"
41 #include "../generic/oomph_utilities.h"
54 template<
unsigned DIM>
61 std::complex<double>&
f);
82 return std::complex<unsigned>(0, 1);
100 "Please set pointer to k_squared using access fct to pointer!",
120 const unsigned& nplot)
const
127 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
139 file_out << u.real() << std::endl;
144 file_out << u.imag() << std::endl;
149 std::stringstream error_stream;
151 <<
"Helmholtz elements only store 2 fields so i must be 0 or 1"
173 return "Imaginary part";
178 std::stringstream error_stream;
180 <<
"Helmholtz elements only store 2 fields so i must be 0 or 1"
195 const unsigned n_plot = 5;
201 void output(std::ostream& outfile,
const unsigned& n_plot);
210 const unsigned& n_plot);
215 const unsigned n_plot = 5;
221 void output(FILE* file_pt,
const unsigned& n_plot);
226 const unsigned& n_plot,
232 std::ostream& outfile,
233 const unsigned& n_plot,
238 "There is no time-dependent output_fct() for Helmholtz elements ",
251 const unsigned& n_plot,
270 "There is no time-dependent compute_error() for Helmholtz elements",
294 std::complex<double>&
source)
const
299 source = std::complex<double>(0.0, 0.0);
314 const unsigned n_node =
nnode();
325 const std::complex<double>
zero(0.0, 0.0);
326 for (
unsigned j = 0;
j <
DIM;
j++)
332 for (
unsigned l = 0; l < n_node; l++)
335 const std::complex<double> u_value(
339 for (
unsigned j = 0;
j <
DIM;
j++)
341 flux[
j] += u_value * dpsidx(l,
j);
373 const unsigned n_node =
nnode();
382 std::complex<double> interpolated_u(0.0, 0.0);
389 for (
unsigned l = 0; l < n_node; l++)
392 const std::complex<double> u_value(
396 interpolated_u += u_value * psi[l];
398 return interpolated_u;
414 DShape& dtestdx)
const = 0;
424 DShape& dtestdx)
const = 0;
431 const unsigned& flag);
450 template<
unsigned DIM,
unsigned NNODE_1D>
490 void output(std::ostream& outfile,
const unsigned& n_plot)
502 const unsigned& n_plot)
518 void output(FILE* file_pt,
const unsigned& n_plot)
527 const unsigned& n_plot,
541 const unsigned& n_plot,
545 outfile, phi, n_plot, exact_soln_pt);
553 const unsigned& n_plot,
591 template<
unsigned DIM,
unsigned NNODE_1D>
600 const double J = this->dshape_eulerian(
s, psi, dpsidx);
617 template<
unsigned DIM,
unsigned NNODE_1D>
626 const double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
647 template<
unsigned DIM,
unsigned NNODE_1D>
649 :
public virtual QElement<DIM - 1, NNODE_1D>
665 template<
unsigned NNODE_1D>
684 template<
class HELMHOLTZ_ELEMENT>
701 std::stringstream error_stream;
702 error_stream <<
"Helmholtz elements only store 2 fields so fld = "
703 << fld <<
" is illegal \n";
710 unsigned nnod = this->
nnode();
714 for (
unsigned j = 0;
j < nnod;
j++)
717 data_values[
j] = std::make_pair(this->
node_pt(
j), fld);
737 std::stringstream error_stream;
738 error_stream <<
"Helmholtz elements only store two fields so fld = "
739 << fld <<
" is illegal\n";
763 std::stringstream error_stream;
764 error_stream <<
"Helmholtz elements only store two fields so fld = "
765 << fld <<
" is illegal.\n";
770 unsigned n_dim = this->
dim();
771 unsigned n_node = this->
nnode();
773 DShape dpsidx(n_node, n_dim), dtestdx(n_node, n_dim);
774 double J = this->dshape_and_dtest_eulerian_helmholtz(
775 s, psi, dpsidx,
test, dtestdx);
789 std::stringstream error_stream;
790 error_stream <<
"Helmholtz elements only store two fields so fld = "
791 << fld <<
" is illegal\n";
797 std::complex<unsigned> complex_u_nodal_index = this->u_index_helmholtz();
798 unsigned u_nodal_index = 0;
801 u_nodal_index = complex_u_nodal_index.real();
805 u_nodal_index = complex_u_nodal_index.imag();
810 unsigned n_node = this->
nnode();
817 double interpolated_u = 0.0;
820 for (
unsigned l = 0; l < n_node; l++)
822 interpolated_u += this->
nodal_value(t, l, u_nodal_index) * psi[l];
824 return interpolated_u;
834 std::stringstream error_stream;
835 error_stream <<
"Helmholtz elements only store two fields so fld = "
836 << fld <<
" is illegal\n";
841 return this->
nnode();
851 std::stringstream error_stream;
852 error_stream <<
"Helmholtz elements only store two fields so fld = "
853 << fld <<
" is illegal\n";
858 std::complex<unsigned> complex_u_nodal_index = this->u_index_helmholtz();
859 unsigned u_nodal_index = 0;
862 u_nodal_index = complex_u_nodal_index.real();
866 u_nodal_index = complex_u_nodal_index.imag();
874 void output(std::ostream& outfile,
const unsigned& nplot)
885 template<
class ELEMENT>
898 template<
class ELEMENT>
AnnoyingScalar imag(const AnnoyingScalar &)
Definition: AnnoyingScalar.h:132
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
unsigned ntstorage() const
Definition: nodes.cc:879
FaceGeometry()
Definition: helmholtz_elements.h:903
FaceGeometry()
Definition: helmholtz_elements.h:890
FaceGeometry()
Definition: helmholtz_elements.h:672
FaceGeometry()
Definition: helmholtz_elements.h:654
Definition: elements.h:4998
Definition: elements.h:1313
virtual unsigned nplot_points_paraview(const unsigned &nplot) const
Definition: elements.h:2862
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual void shape(const Vector< double > &s, Shape &psi) const =0
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned dim() const
Definition: elements.h:2611
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Definition: elements.h:1765
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3298
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
Definition: helmholtz_elements.h:56
HelmholtzSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: helmholtz_elements.h:276
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: helmholtz_elements.cc:374
std::string scalar_name_paraview(const unsigned &i) const
Definition: helmholtz_elements.h:164
void get_flux(const Vector< double > &s, Vector< std::complex< double >> &flux) const
Get flux: flux[i] = du/dx_i for real and imag part.
Definition: helmholtz_elements.h:310
unsigned self_test()
Self-test: Return 0 for OK.
Definition: helmholtz_elements.cc:224
virtual void fill_in_generic_residual_contribution_helmholtz(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: helmholtz_elements.cc:49
virtual std::complex< unsigned > u_index_helmholtz() const
Broken assignment operator.
Definition: helmholtz_elements.h:80
unsigned nscalar_paraview() const
Definition: helmholtz_elements.h:111
HelmholtzEquations(const HelmholtzEquations &dummy)=delete
Broken copy constructor.
void output_real(std::ostream &outfile, const double &phi, const unsigned &n_plot)
Definition: helmholtz_elements.cc:294
virtual double dshape_and_dtest_eulerian_at_knot_helmholtz(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
HelmholtzEquations()
Constructor (must initialise the Source_fct_pt to null)
Definition: helmholtz_elements.h:65
double * K_squared_pt
Pointer to square of wavenumber.
Definition: helmholtz_elements.h:437
void output(FILE *file_pt)
C_style output with default number of plot points.
Definition: helmholtz_elements.h:213
void(* HelmholtzSourceFctPt)(const Vector< double > &x, std::complex< double > &f)
Definition: helmholtz_elements.h:60
virtual void get_source_helmholtz(const unsigned &ipt, const Vector< double > &x, std::complex< double > &source) const
Definition: helmholtz_elements.h:292
void scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Definition: helmholtz_elements.h:118
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
Definition: helmholtz_elements.h:348
HelmholtzSourceFctPt Source_fct_pt
Pointer to source function:
Definition: helmholtz_elements.h:434
std::complex< double > interpolated_u_helmholtz(const Vector< double > &s) const
Definition: helmholtz_elements.h:369
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: helmholtz_elements.cc:482
virtual double dshape_and_dtest_eulerian_helmholtz(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
void output_real_fct(std::ostream &outfile, const double &phi, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: helmholtz_elements.cc:429
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
Definition: helmholtz_elements.h:263
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: helmholtz_elements.h:231
double k_squared()
Get the square of wavenumber.
Definition: helmholtz_elements.h:94
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: helmholtz_elements.h:193
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: helmholtz_elements.h:359
double *& k_squared_pt()
Get pointer to square of wavenumber.
Definition: helmholtz_elements.h:87
HelmholtzSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
Definition: helmholtz_elements.h:282
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
Definition: nodes.h:1022
Definition: oomph_definitions.h:222
Definition: elements.h:3439
Definition: projection.h:183
Helmholtz upgraded to become projectable.
Definition: helmholtz_elements.h:687
ProjectableHelmholtzElement()
Definition: helmholtz_elements.h:691
unsigned nhistory_values_for_coordinate_projection()
Definition: helmholtz_elements.h:749
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: helmholtz_elements.h:756
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: helmholtz_elements.h:782
unsigned nfields_for_projection()
Number of fields to be projected: 2 (real and imag part)
Definition: helmholtz_elements.h:725
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: helmholtz_elements.h:732
void output(std::ostream &outfile, const unsigned &nplot)
Definition: helmholtz_elements.h:874
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
Definition: helmholtz_elements.h:829
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: helmholtz_elements.h:846
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: helmholtz_elements.h:696
Definition: Qelements.h:459
Definition: helmholtz_elements.h:453
void output(FILE *file_pt, const unsigned &n_plot)
Definition: helmholtz_elements.h:518
static const unsigned Initial_Nvalue
Definition: helmholtz_elements.h:457
void output_real(std::ostream &outfile, const double &phi, const unsigned &n_plot)
Definition: helmholtz_elements.h:500
void output_real_fct(std::ostream &outfile, const double &phi, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: helmholtz_elements.h:539
QHelmholtzElement(const QHelmholtzElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
QHelmholtzElement()
Definition: helmholtz_elements.h:462
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
Definition: helmholtz_elements.h:475
void output(std::ostream &outfile)
Definition: helmholtz_elements.h:482
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: helmholtz_elements.h:490
void output(FILE *file_pt)
Definition: helmholtz_elements.h:510
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: helmholtz_elements.h:552
double dshape_and_dtest_eulerian_at_knot_helmholtz(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: helmholtz_elements.h:619
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: helmholtz_elements.h:526
double dshape_and_dtest_eulerian_helmholtz(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: helmholtz_elements.h:592
unsigned ntstorage() const
Definition: timesteppers.h:601
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237
float real
Definition: datatypes.h:10
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
squared absolute sa ArrayBase::abs2 DOXCOMMA MatrixBase::cwiseAbs2 sa Eigen::abs2 DOXCOMMA Eigen::pow DOXCOMMA ArrayBase::square nearest sa Eigen::floor DOXCOMMA Eigen::ceil DOXCOMMA ArrayBase::round nearest integer not less than the given sa Eigen::floor DOXCOMMA ArrayBase::ceil not a number test
Definition: GlobalFunctions.h:109
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
int error
Definition: calibrate.py:297
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
list x
Definition: plotDoE.py:28
t
Definition: plotPSD.py:36
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
EIGEN_DONT_INLINE Scalar zero()
Definition: svd_common.h:232
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2