35 namespace Legendre_functions_helper
42 if (l == 0)
return 1.0;
50 double plgndr1(
const unsigned& n,
const double& x)
58 if (std::fabs(x) > 1.0)
60 std::ostringstream error_stream;
61 error_stream <<
"Bad arguments in routine plgndr1: x=" << x
62 <<
" but should be less than 1 in absolute value.\n";
64 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
82 for (
i = 2;
i <= n;
i++)
84 pmm2 = (x * (2 *
i - 1) * pmm1 - (
i - 1) * pmm) /
i;
97 double plgndr2(
const unsigned& l,
const unsigned& m,
const double& x)
100 double fact, pmm, pmmp1, somx2;
105 if (std::fabs(x) > 1.0)
107 std::ostringstream error_stream;
108 error_stream <<
"Bad arguments in routine plgndr2: x=" << x
109 <<
" but should be less than 1 in absolute value.\n";
111 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
125 somx2 = sqrt((1.0 - x) * (1.0 + x));
127 for (
i = 1;
i <= m;
i++)
129 pmm *= -fact * somx2;
133 if (l == m)
return pmm;
138 pmmp1 = x * (2 * m + 1) * pmm;
146 for (ll = m + 2; ll <= l; ll++)
148 pll = (x * (2 * ll - 1) * pmmp1 - (ll + m - 1) * pmm) / (ll - m);
169 template<
unsigned NNODE_1D>
186 const unsigned& flag)
189 const unsigned n_node =
nnode();
192 Shape psi(n_node), test(n_node);
193 DShape dpsidx(n_node, 2), dtestdx(n_node, 2);
199 int local_eqn_real = 0, local_unknown_real = 0;
200 int local_eqn_imag = 0, local_unknown_imag = 0;
203 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
210 ipt, psi, dpsidx, test, dtestdx);
217 std::complex<double> interpolated_u(0.0, 0.0);
224 for (
unsigned l = 0; l < n_node; l++)
227 for (
unsigned j = 0; j < 2; j++)
233 const std::complex<double> u_value(
238 interpolated_u += u_value * psi(l);
241 for (
unsigned j = 0; j < 2; j++)
243 interpolated_dudx[j] += u_value * dpsidx(l, j);
249 std::complex<double> source(0.0, 0.0);
255 double n_squared = n * n;
261 for (
unsigned l = 0; l < n_node; l++)
271 if (local_eqn_real >= 0)
275 residuals[local_eqn_real] +=
277 ((-n_squared / rr) +
k_squared()) * interpolated_u.real()) *
281 for (
unsigned k = 0; k < 2; k++)
283 residuals[local_eqn_real] +=
284 interpolated_dudx[k].real() * dtestdx(l, k) * r *
W;
292 for (
unsigned l2 = 0; l2 < n_node; l2++)
298 if (local_unknown_real >= 0)
301 for (
unsigned i = 0;
i < 2;
i++)
303 jacobian(local_eqn_real, local_unknown_real) +=
304 dpsidx(l2,
i) * dtestdx(l,
i) * r *
W;
308 jacobian(local_eqn_real, local_unknown_real) -=
309 ((-n_squared / rr) +
k_squared()) * psi(l2) * test(l) * r *
W;
324 if (local_eqn_imag >= 0)
327 residuals[local_eqn_imag] +=
329 ((-n_squared / rr) +
k_squared()) * interpolated_u.imag()) *
333 for (
unsigned k = 0; k < 2; k++)
335 residuals[local_eqn_imag] +=
336 interpolated_dudx[k].imag() * dtestdx(l, k) * r *
W;
344 for (
unsigned l2 = 0; l2 < n_node; l2++)
350 if (local_unknown_imag >= 0)
353 for (
unsigned i = 0;
i < 2;
i++)
355 jacobian(local_eqn_imag, local_unknown_imag) +=
356 dpsidx(l2,
i) * dtestdx(l,
i) * r *
W;
359 jacobian(local_eqn_imag, local_unknown_imag) -=
360 ((-n_squared / rr) +
k_squared()) * psi(l2) * test(l) * r *
W;
403 const unsigned& nplot)
413 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
418 for (
unsigned i = 0;
i < 2;
i++)
422 outfile << u.real() <<
" " << u.imag() << std::endl;
443 const unsigned& nplot)
453 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
458 for (
unsigned i = 0;
i < 2;
i++)
462 outfile << u.real() * cos(phi) + u.imag() * sin(phi) << std::endl;
478 const unsigned& nplot)
488 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
494 for (
unsigned i = 0;
i < 2;
i++)
499 for (
unsigned i = 0;
i < 2;
i++)
503 fprintf(file_pt,
"%g ", u.real());
504 fprintf(file_pt,
"%g \n", u.imag());
521 std::ostream& outfile,
522 const unsigned& nplot,
539 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
548 (*exact_soln_pt)(x, exact_soln);
551 for (
unsigned i = 0;
i < 2;
i++)
553 outfile << x[
i] <<
" ";
555 outfile << exact_soln[0] <<
" " << exact_soln[1] << std::endl;
575 std::ostream& outfile,
577 const unsigned& nplot,
594 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
603 (*exact_soln_pt)(x, exact_soln);
606 for (
unsigned i = 0;
i < 2;
i++)
608 outfile << x[
i] <<
" ";
610 outfile << exact_soln[0] * cos(phi) + exact_soln[1] * sin(phi)
627 std::ostream& outfile,
643 unsigned n_node =
nnode();
651 outfile <<
"ZONE" << std::endl;
657 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
660 for (
unsigned i = 0;
i < 2;
i++)
678 std::complex<double> u_fe =
682 (*exact_soln_pt)(x, exact_soln);
685 for (
unsigned i = 0;
i < 2;
i++)
687 outfile << x[
i] <<
" ";
689 outfile << exact_soln[0] <<
" " << exact_soln[1] <<
" "
690 << exact_soln[0] - u_fe.real() <<
" "
691 << exact_soln[1] - u_fe.imag() << std::endl;
695 (exact_soln[0] * exact_soln[0] + exact_soln[1] * exact_soln[1]) *
W;
696 error += ((exact_soln[0] - u_fe.real()) * (exact_soln[0] - u_fe.real()) +
697 (exact_soln[1] - u_fe.imag()) * (exact_soln[1] - u_fe.imag())) *
718 unsigned n_node =
nnode();
726 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
729 for (
unsigned i = 0;
i < 2;
i++)
744 std::complex<double> u_fe =
748 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.
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 void fill_in_generic_residual_contribution_fourier_decomposed_helmholtz(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Compute element residual Vector only (if flag=and/or element Jacobian matrix.
unsigned self_test()
Self-test: Return 0 for OK.
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 compute_norm(double &norm)
Compute norm of fe solution.
double k_squared()
Get the square of wavenumber.
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.
virtual double dshape_and_dtest_eulerian_at_knot_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 ...
void output(std::ostream &outfile)
Output with default number of plot points.
std::complex< double > interpolated_u_fourier_decomposed_helmholtz(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
virtual std::complex< unsigned > u_index_fourier_decomposed_helmholtz() const
Broken assignment operator.
virtual void get_source_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...
int fourier_wavenumber()
Get the Fourier wavenumber.
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.
An OomphLibError object which should be thrown when an run-time error is encountered....
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
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.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...