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)
61 error_stream <<
"Bad arguments in routine plgndr1: x=" << x
62 <<
" but should be less than 1 in absolute value.\n";
82 for (
i = 2;
i <=
n;
i++)
97 double plgndr2(
const unsigned&
l,
const unsigned&
m,
const double& x)
105 if (std::fabs(x) > 1.0)
108 error_stream <<
"Bad arguments in routine plgndr2: x=" << x
109 <<
" but should be less than 1 in absolute value.\n";
127 for (
i = 1;
i <=
m;
i++)
133 if (
l ==
m)
return pmm;
169 template<
unsigned NNODE_1D>
186 const unsigned&
flag)
217 std::complex<double> interpolated_u(0.0, 0.0);
227 for (
unsigned j = 0;
j < 2;
j++)
233 const std::complex<double>
u_value(
241 for (
unsigned j = 0;
j < 2;
j++)
249 std::complex<double> source(0.0, 0.0);
281 for (
unsigned k = 0;
k < 2;
k++)
284 interpolated_dudx[
k].real() *
dtestdx(
l,
k) *
r * W;
301 for (
unsigned i = 0;
i < 2;
i++)
333 for (
unsigned k = 0;
k < 2;
k++)
336 interpolated_dudx[
k].imag() *
dtestdx(
l,
k) *
r * W;
353 for (
unsigned i = 0;
i < 2;
i++)
403 const unsigned&
nplot)
418 for (
unsigned i = 0;
i < 2;
i++)
422 outfile << u.real() <<
" " << u.imag() << std::endl;
443 const unsigned&
nplot)
458 for (
unsigned i = 0;
i < 2;
i++)
478 const unsigned&
nplot)
494 for (
unsigned i = 0;
i < 2;
i++)
499 for (
unsigned i = 0;
i < 2;
i++)
522 const unsigned&
nplot,
551 for (
unsigned i = 0;
i < 2;
i++)
577 const unsigned&
nplot,
606 for (
unsigned i = 0;
i < 2;
i++)
651 outfile <<
"ZONE" << std::endl;
660 for (
unsigned i = 0;
i < 2;
i++)
678 std::complex<double>
u_fe =
685 for (
unsigned i = 0;
i < 2;
i++)
729 for (
unsigned i = 0;
i < 2;
i++)
744 std::complex<double>
u_fe =
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.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
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 .
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.
virtual std::complex< unsigned > u_index_fourier_decomposed_helmholtz() const
Broken assignment operator.
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.
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.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
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...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
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.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).