39 namespace Legendre_functions_helper
46 if (l == 0)
return 1.0;
54 double plgndr1(
const unsigned& n,
const double& x)
62 if (std::fabs(x) > 1.0)
64 std::ostringstream error_stream;
65 error_stream <<
"Bad arguments in routine plgndr1: x=" << x
66 <<
" but should be less than 1 in absolute value.\n";
68 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
86 for (
i = 2;
i <= n;
i++)
88 pmm2 = (x * (2 *
i - 1) * pmm1 - (
i - 1) * pmm) /
i;
101 double plgndr2(
const unsigned& l,
const unsigned& m,
const double& x)
104 double fact, pmm, pmmp1, somx2;
109 if (std::fabs(x) > 1.0)
111 std::ostringstream error_stream;
112 error_stream <<
"Bad arguments in routine plgndr2: x=" << x
113 <<
" but should be less than 1 in absolute value.\n";
115 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
129 somx2 = sqrt((1.0 - x) * (1.0 + x));
131 for (
i = 1;
i <= m;
i++)
133 pmm *= -fact * somx2;
137 if (l == m)
return pmm;
142 pmmp1 = x * (2 * m + 1) * pmm;
150 for (ll = m + 2; ll <= l; ll++)
152 pll = (x * (2 * ll - 1) * pmmp1 - (ll + m - 1) * pmm) / (ll - m);
173 template<
unsigned NNODE_1D>
195 const unsigned& flag)
198 const unsigned n_node =
nnode();
201 Shape psi(n_node), test(n_node);
202 DShape dpsidx(n_node, 2), dtestdx(n_node, 2);
208 int local_eqn_real = 0, local_unknown_real = 0;
209 int local_eqn_imag = 0, local_unknown_imag = 0;
212 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
220 ipt, psi, dpsidx, test, dtestdx);
227 std::complex<double> interpolated_u(0.0, 0.0);
234 for (
unsigned l = 0; l < n_node; l++)
237 for (
unsigned j = 0; j < 2; j++)
243 const std::complex<double> u_value(
249 interpolated_u += u_value * psi(l);
252 for (
unsigned j = 0; j < 2; j++)
254 interpolated_dudx[j] += u_value * dpsidx(l, j);
260 std::complex<double> source(0.0, 0.0);
264 double n_squared = n * n;
270 std::complex<double> pml_k_squared_factor =
271 std::complex<double>(1.0, 0.0);
286 std::complex<double> complex_r = std::complex<double>(1.0, 0.0);
291 std::complex<double> pml_k_squared_jacobian =
292 -pml_k_squared_factor *
293 (
k_squared() - n_squared / (complex_r * complex_r)) * complex_r *
W;
297 for (
unsigned k = 0; k < 2; k++)
299 pml_laplace_jacobian[k] = pml_laplace_factor[k] * complex_r *
W;
305 std::complex<double> pml_k_squared_residual =
306 pml_k_squared_jacobian * interpolated_u;
310 for (
unsigned k = 0; k < 2; k++)
312 pml_laplace_residual[k] =
313 pml_laplace_jacobian[k] * interpolated_dudx[k];
319 for (
unsigned l = 0; l < n_node; l++)
331 if (local_eqn_real >= 0)
334 residuals[local_eqn_real] += pml_k_squared_residual.real() * test(l);
337 for (
unsigned k = 0; k < 2; k++)
339 residuals[local_eqn_real] +=
340 pml_laplace_residual[k].real() * dtestdx(l, k);
348 for (
unsigned l2 = 0; l2 < n_node; l2++)
357 if (local_unknown_real >= 0)
360 jacobian(local_eqn_real, local_unknown_real) +=
361 pml_k_squared_jacobian.real() * psi(l2) * test(l);
364 for (
unsigned k = 0; k < 2; k++)
366 jacobian(local_eqn_real, local_unknown_real) +=
367 pml_laplace_jacobian[k].real() * dpsidx(l2, k) *
372 if (local_unknown_imag >= 0)
375 jacobian(local_eqn_real, local_unknown_imag) +=
376 -pml_k_squared_jacobian.imag() * psi(l2) * test(l);
379 for (
unsigned k = 0; k < 2; k++)
381 jacobian(local_eqn_real, local_unknown_imag) +=
382 -pml_laplace_jacobian[k].imag() * dpsidx(l2, k) *
395 if (local_eqn_imag >= 0)
398 residuals[local_eqn_imag] += pml_k_squared_residual.imag() * test(l);
401 for (
unsigned k = 0; k < 2; k++)
403 residuals[local_eqn_imag] +=
404 pml_laplace_residual[k].imag() * dtestdx(l, k);
412 for (
unsigned l2 = 0; l2 < n_node; l2++)
420 if (local_unknown_imag >= 0)
423 jacobian(local_eqn_imag, local_unknown_imag) +=
424 pml_k_squared_jacobian.real() * psi(l2) * test(l);
427 for (
unsigned k = 0; k < 2; k++)
429 jacobian(local_eqn_imag, local_unknown_imag) +=
430 pml_laplace_jacobian[k].real() * dpsidx(l2, k) *
435 if (local_unknown_real >= 0)
438 jacobian(local_eqn_imag, local_unknown_real) +=
439 pml_k_squared_jacobian.imag() * psi(l2) * test(l);
442 for (
unsigned k = 0; k < 2; k++)
444 jacobian(local_eqn_imag, local_unknown_real) +=
445 pml_laplace_jacobian[k].imag() * dpsidx(l2, k) *
490 const unsigned& nplot)
500 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
504 std::complex<double> u(
506 for (
unsigned i = 0;
i < 2;
i++)
510 outfile << u.real() <<
" " << u.imag() << std::endl;
530 std::ostream& outfile,
const double& phi,
const unsigned& nplot)
540 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
544 std::complex<double> u(
546 for (
unsigned i = 0;
i < 2;
i++)
550 outfile << u.real() * cos(phi) + u.imag() * sin(phi) << std::endl;
566 const unsigned& nplot)
576 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
580 std::complex<double> u(
583 for (
unsigned i = 0;
i < 2;
i++)
588 for (
unsigned i = 0;
i < 2;
i++)
592 fprintf(file_pt,
"%g ", u.real());
593 fprintf(file_pt,
"%g \n", u.imag());
610 std::ostream& outfile,
611 const unsigned& nplot,
628 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
637 (*exact_soln_pt)(x, exact_soln);
640 for (
unsigned i = 0;
i < 2;
i++)
642 outfile << x[
i] <<
" ";
644 outfile << exact_soln[0] <<
" " << exact_soln[1] << std::endl;
664 std::ostream& outfile,
666 const unsigned& nplot,
683 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
692 (*exact_soln_pt)(x, exact_soln);
695 for (
unsigned i = 0;
i < 2;
i++)
697 outfile << x[
i] <<
" ";
699 outfile << exact_soln[0] * cos(phi) + exact_soln[1] * sin(phi)
716 std::ostream& outfile,
732 unsigned n_node =
nnode();
740 outfile <<
"ZONE" << std::endl;
746 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
749 for (
unsigned i = 0;
i < 2;
i++)
767 std::complex<double> u_fe =
771 (*exact_soln_pt)(x, exact_soln);
774 for (
unsigned i = 0;
i < 2;
i++)
776 outfile << x[
i] <<
" ";
778 outfile << exact_soln[0] <<
" " << exact_soln[1] <<
" "
779 << exact_soln[0] - u_fe.real() <<
" "
780 << exact_soln[1] - u_fe.imag() << std::endl;
784 (exact_soln[0] * exact_soln[0] + exact_soln[1] * exact_soln[1]) *
W;
785 error += ((exact_soln[0] - u_fe.real()) * (exact_soln[0] - u_fe.real()) +
786 (exact_soln[1] - u_fe.imag()) * (exact_soln[1] - u_fe.imag())) *
807 unsigned n_node =
nnode();
815 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
818 for (
unsigned i = 0;
i < 2;
i++)
833 std::complex<double> u_fe =
837 norm += (u_fe.real() * u_fe.real() + u_fe.imag() * u_fe.imag()) *
W;
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
virtual double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s.
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
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 nnode() const
Return the number of nodes.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n but do NOT take hanging nodes into account.
double raw_nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. Do not use the hanging node representation....
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Add tecplot zone "footer" to output stream (when plotting nplot points in each "coordinate direction"...
virtual unsigned self_test()
Self-test: Check inversion of element & do self-test for GeneralisedElement. Return 0 if OK.
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.
std::complex< double > interpolated_u_pml_fourier_decomposed_helmholtz(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
int pml_fourier_wavenumber()
Get the Fourier wavenumber.
virtual void fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Compute element residual Vector only (if flag=and/or element Jacobian matrix.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
void output_real(std::ostream &outfile, const double &phi, const unsigned &n_plot)
Output function for real part of full time-dependent solution u = Re( (u_r +i u_i) exp(-i omega t) at...
void output(std::ostream &outfile)
Output with default number of plot points.
void compute_complex_r(const unsigned &ipt, const Vector< double > &x, std::complex< double > &complex_r)
Compute complex variable r at position x[0] and integration point ipt.
static double Default_Physical_Constant_Value
Static default value for the physical constants (initialised to zero)
unsigned self_test()
Self-test: Return 0 for OK.
double k_squared()
Get k squared.
virtual void get_source_pml_fourier_decomposed_helmholtz(const unsigned &ipt, const Vector< double > &x, std::complex< double > &source) const
Get source term at (Eulerian) position x. This function is virtual to allow overloading in multi-phys...
void compute_norm(double &norm)
Compute norm of fe solution.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,z,u_re_exact,u_im_exact at n_plot^2 plot points.
void compute_pml_coefficients(const unsigned &ipt, const Vector< double > &x, Vector< std::complex< double >> &pml_laplace_factor, std::complex< double > &pml_k_squared_factor)
Compute pml coefficients at position x and integration point ipt. pml_laplace_factor is used in the r...
void output_real_fct(std::ostream &outfile, const double &phi, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for real part of full time-dependent fct u = Re( (u_r +i u_i) exp(-i omega t) at phas...
virtual double dshape_and_dtest_eulerian_at_knot_pml_fourier_decomposed_helmholtz(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
static BermudezPMLMappingAndTransformedCoordinate Default_pml_mapping_and_transformed_coordinate
Static so that the class doesn't need to instantiate a new default everytime it uses it.
virtual std::complex< unsigned > u_index_pml_fourier_decomposed_helmholtz() const
Broken assignment operator.
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
static const unsigned Initial_Nvalue
Static int that holds the number of variables at nodes: always the same.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
double factorial(const unsigned &l)
Factorial.
double plgndr1(const unsigned &n, const double &x)
Legendre polynomials depending on one parameter.
double plgndr2(const unsigned &l, const unsigned &m, const double &x)
Legendre polynomials depending on two parameters.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...