10#include "scaleFElib.h"
17 use scale_const,
only: &
56 private :: gen_jacobigaussquadraturepts
67 integer,
intent(in) :: nord
68 real(rp),
intent(in) :: x_lgl(nord+1)
69 real(rp),
intent(in) :: x(:)
70 real(rp) :: l(size(x), nord+1)
73 real(rp) :: p_lgl(nord+1,nord+1)
74 real(rp) :: p(size(x),nord+1)
75 real(rp) :: pr(size(x),nord+1)
84 if ( abs(x(i)-x_lgl(n)) < 1e-15_rp )
then
88 (x(i) - 1.0_rp)*(x(i) + 1.0_rp)*pr(i,nord+1) &
89 / (dble(nord*(nord+1))*p_lgl(n,nord+1)*(x(i) - x_lgl(n)))
103 integer,
intent(in) :: nord
104 real(rp),
intent(in) :: x_lgl(nord+1)
105 real(rp) :: lr(nord+1, nord+1)
108 real(rp) :: p(nord+1,nord+1)
117 if (k==1 .and. n==1)
then
118 lr(k,n) = - 0.25_rp*dble(nord*(nord+1))
119 else if (k==nord+1 .and. n==nord+1)
then
120 lr(k,n) = + 0.25_rp*dble(nord*(nord+1))
124 lr(k,n) = p(n,nord+1)/(p(k,nord+1)*(x_lgl(n) - x_lgl(k)))
128 lr_nn = lr_nn + lr(k,n)
143 integer,
intent(in) :: nord
144 real(rp),
intent(in) :: x(:)
145 real(rp),
intent(out) :: p(size(x), nord+1)
151 log_error(
"Polynominal_GenLegendrePoly",*)
"Nord must be larger than 0."
159 p(:,n+1) = ( dble(2*n-1)*x(:)*p(:,n) - dble(n-1)*p(:,n-1) )/dble(n)
171 integer,
intent(in) :: nord
172 real(rp),
intent(in) :: x(:)
173 real(rp) :: p(size(x), nord+1)
187 integer,
intent(in) :: nord
188 real(rp),
intent(in) :: x(:)
189 real(rp),
intent(in) :: p(:,:)
190 real(rp) :: gradp(size(x), nord+1)
196 log_error(
"Polynominal_GenDLegendrePoly",*)
"Nord must be larger than 0."
200 if (nord == 0)
return
204 gradp(:,n+1) = 2.0_rp*x(:)*gradp(:,n) - gradp(:,n-1) + p(:,n)
216 integer,
intent(in) :: nord
217 real(rp) :: pts(nord+1)
222 pts(1) = -1.0_rp; pts(nord+1) = 1.0_rp
225 call gen_jacobigaussquadraturepts( 1, 1, nord-2, pts(2:nord) )
235 integer,
intent(in) :: nord
236 real(rp) :: int_weight_lgl(nord+1)
238 real(rp) :: lglpts1d(nord+1)
239 real(rp) :: p1d_ori(nord+1, nord+1)
245 int_weight_lgl(:) = 2.0_rp/(dble(nord*(nord+1))*p1d_ori(:,nord+1)**2)
256 integer,
intent(in) :: nord
257 real(rp) :: pts(nord)
260 call gen_jacobigaussquadraturepts( 0, 0, nord-1, pts(:) )
270 integer,
intent(in) :: nord
271 real(rp) :: int_weight_gl(nord)
273 real(rp) :: glpts1d(nord)
274 real(rp) :: p1d_ori(nord, nord+1)
275 real(rp) :: dp1d_ori(nord, nord+1)
282 int_weight_gl(:) = 2.0_rp / ( (1.0_rp - glpts1d(:)**2) * dp1d_ori(:,nord+1)**2 )
291 subroutine gen_jacobigaussquadraturepts( alpha, beta, N, &
295 integer,
intent(in) :: alpha
296 integer,
intent(in) :: beta
297 integer,
intent(in) :: n
298 real(rp),
intent(out) :: x(n+1)
301 real(dp) :: d(n+1), e(n)
302 real(dp) :: work(2*(n+1)-2), z(n+1,n+1)
308 x(1) = - dble(alpha - beta) / dble(alpha + beta + 2)
313 h1(i+1) = dble( 2*i + alpha + beta )
317 d(i) = - dble(alpha**2 - beta**2) / (h1(i) * (h1(i) + 2d0))
320 e(i) = 2d0 / (h1(i) + 2d0) &
321 * sqrt( dble(i * (i + alpha + beta) * (i + alpha) * (i + beta)) &
322 / ((h1(i) + 1d0) * (h1(i) + 3d0)) )
324 if ( dble(alpha + beta) < 1d-16) d(1) = 0.0_rp
326 call dstev(
'Vectors', n+1, d, e, z, n+1, work, info )
330 end subroutine gen_jacobigaussquadraturepts
module common / Polynominal
real(rp) function, dimension(nord), public polynominal_gengausslegendreptintweight(nord)
A function to calcuate the Gauss-Legendre weights.
real(rp) function, dimension(nord), public polynominal_gengausslegendrept(nord)
A function to calcuate the Gauss-Legendre points.
real(rp) function, dimension(size(x), nord+1), public polynominal_genlagrangepoly(nord, x_lgl, x)
A function to obtain the values of Lagrange basis functions which are evaluated over aribitary points...
real(rp) function, dimension(nord+1), public polynominal_gengausslobattopt(nord)
A function to calcuate the Legendre-Gauss-Lobtatto (LGL) points.
real(rp) function, dimension(size(x), nord+1), public polynominal_genlegendrepoly(nord, x)
A function to obtain the values of Legendre polynominals which are evaluated at aribitary points.
real(rp) function, dimension(size(x), nord+1), public polynominal_gendlegendrepoly(nord, x, p)
A function to obtain differential values of Legendre polynominals which are evaluated at aribitary po...
real(rp) function, dimension(nord+1), public polynominal_gengausslobattoptintweight(nord)
A function to calcuate 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 which are evaluated over ari...
subroutine, public polynominal_genlegendrepoly_sub(nord, x, p)
A function to obtain the values of Legendre polynominals which are evaluated at aribitary points.