10#include "scaleFElib.h"
17 use scale_const,
only: &
56 private :: gen_jacobigaussquadraturepts
70 integer,
intent(in) :: nord
71 real(rp),
intent(in) :: x_lgl(nord+1)
72 real(rp),
intent(in) :: x(:)
73 real(rp) :: l(size(x), nord+1)
76 real(rp) :: p_lgl(nord+1,nord+1)
77 real(rp) :: p(size(x),nord+1)
78 real(rp) :: pr(size(x),nord+1)
90 if ( abs(x(i)-x_lgl(n)) < 1e-15_rp )
then
94 (x(i) - 1.0_rp)*(x(i) + 1.0_rp)*pr(i,nord+1) &
95 / (dble(nord*(nord+1))*p_lgl(n,nord+1)*(x(i) - x_lgl(n)))
129 integer,
intent(in) :: nord
130 real(rp),
intent(in) :: x_lgl(nord+1)
131 real(rp) :: lr(nord+1, nord+1)
134 real(rp) :: p(nord+1,nord+1)
143 if (k==1 .and. n==1)
then
144 lr(k,n) = - 0.25_rp*dble(nord*(nord+1))
145 else if (k==nord+1 .and. n==nord+1)
then
146 lr(k,n) = + 0.25_rp*dble(nord*(nord+1))
150 lr(k,n) = p(n,nord+1)/(p(k,nord+1)*(x_lgl(n) - x_lgl(k)))
154 lr_nn = lr_nn + lr(k,n)
172 integer,
intent(in) :: nord
173 real(rp),
intent(in) :: x(:)
174 real(rp),
intent(out) :: p(size(x), nord+1)
180 log_error(
"Polynominal_GenLegendrePoly",*)
"Nord must be larger than 0."
188 p(:,n+1) = ( dble(2*n-1)*x(:)*p(:,n) - dble(n-1)*p(:,n-1) )/dble(n)
203 integer,
intent(in) :: nord
204 real(rp),
intent(in) :: x(:)
205 real(rp) :: p(size(x), nord+1)
222 integer,
intent(in) :: nord
223 real(rp),
intent(in) :: x(:)
224 real(rp),
intent(in) :: p(:,:)
225 real(rp) :: gradp(size(x), nord+1)
231 log_error(
"Polynominal_GenDLegendrePoly",*)
"Nord must be larger than 0."
235 if (nord == 0)
return
239 gradp(:,n+1) = 2.0_rp*x(:)*gradp(:,n) - gradp(:,n-1) + p(:,n)
253 integer,
intent(in) :: nord
254 real(rp) :: pts(nord+1)
259 pts(1) = -1.0_rp; pts(nord+1) = 1.0_rp
262 call gen_jacobigaussquadraturepts( 1, 1, nord-2, pts(2:nord) )
274 integer,
intent(in) :: nord
275 real(rp) :: int_weight_lgl(nord+1)
277 real(rp) :: lglpts1d(nord+1)
278 real(rp) :: p1d_ori(nord+1, nord+1)
284 int_weight_lgl(:) = 2.0_rp/(dble(nord*(nord+1))*p1d_ori(:,nord+1)**2)
297 integer,
intent(in) :: nord
298 real(rp) :: pts(nord)
301 call gen_jacobigaussquadraturepts( 0, 0, nord-1, pts(:) )
313 integer,
intent(in) :: nord
314 real(rp) :: int_weight_gl(nord)
316 real(rp) :: glpts1d(nord)
317 real(rp) :: p1d_ori(nord, nord+1)
318 real(rp) :: dp1d_ori(nord, nord+1)
325 int_weight_gl(:) = 2.0_rp / ( (1.0_rp - glpts1d(:)**2) * dp1d_ori(:,nord+1)**2 )
334 subroutine gen_jacobigaussquadraturepts( alpha, beta, N, &
338 integer,
intent(in) :: alpha
339 integer,
intent(in) :: beta
340 integer,
intent(in) :: n
341 real(rp),
intent(out) :: x(n+1)
344 real(dp) :: d(n+1), e(n)
345 real(dp) :: work(2*(n+1)-2), z(n+1,n+1)
351 x(1) = - dble(alpha - beta) / dble(alpha + beta + 2)
356 h1(i+1) = dble( 2*i + alpha + beta )
360 d(i) = - dble(alpha**2 - beta**2) / (h1(i) * (h1(i) + 2d0))
363 e(i) = 2d0 / (h1(i) + 2d0) &
364 * sqrt( dble(i * (i + alpha + beta) * (i + alpha) * (i + beta)) &
365 / ((h1(i) + 1d0) * (h1(i) + 3d0)) )
367 if ( dble(alpha + beta) < 1d-16) d(1) = 0.0_rp
369 call dstev(
'Vectors', n+1, d, e, z, n+1, work, info )
373 end subroutine gen_jacobigaussquadraturepts
Module common / Polynominal.
real(rp) function, dimension(nord), public polynominal_gengausslegendreptintweight(nord)
A function to calculate the Gauss-Legendre (GL) weights.
real(rp) function, dimension(nord), public polynominal_gengausslegendrept(nord)
A function to calculate the Gauss-Legendre (GL) points.
real(rp) function, dimension(size(x), nord+1), public polynominal_genlagrangepoly(nord, x_lgl, x)
A function to obtain the Lagrange basis functions related to the Gauss-Legendre-Lobatto (GLL) points.
real(rp) function, dimension(nord+1), public polynominal_gengausslobattopt(nord)
A function to calculate the Legendre-Gauss-Lobatto (LGL) points.
real(rp) function, dimension(size(x), nord+1), public polynominal_genlegendrepoly(nord, x)
A function to obtain the values of Legendre polynomials which are evaluated at arbitrary points.
real(rp) function, dimension(size(x), nord+1), public polynominal_gendlegendrepoly(nord, x, p)
A function to obtain differential values of Legendre polynomials which are evaluated at arbitrary poi...
real(rp) function, dimension(nord+1), public polynominal_gengausslobattoptintweight(nord)
A function to calculate the Gauss-Lobbato weights.
real(rp) function, dimension(nord+1, nord+1), public polynominal_gendlagrangepoly_lglpt(nord, x_lgl)
A function to obtain the differential values of Lagrange basis functions at the GLL points.
subroutine, public polynominal_genlegendrepoly_sub(nord, x, p)
A function to obtain the values of Legendre polynomials which are evaluated at arbitrary points.