FE-Project
Loading...
Searching...
No Matches
scale_polynominal.F90
Go to the documentation of this file.
1!> Module common / Polynominal
2!!
3!! @par Description
4!! A module to provide utilities for polynomials
5!!
6!! @par Reference
7!!
8!! @author Yuta Kawai, Team SCALE
9!!
10#include "scaleFElib.h"
12 !-----------------------------------------------------------------------------
13 !
14 !++ used modules
15 !
16 use scale_precision
17 use scale_const, only: &
18 pi => const_pi
19 use scale_io
20
21 !-----------------------------------------------------------------------------
22 implicit none
23 private
24 !-----------------------------------------------------------------------------
25 !
26 !++ Public procedure
27 !
28
32
35
38
41
42 !-----------------------------------------------------------------------------
43 !
44 !++ Public parameters & variables
45 !
46
47 !-----------------------------------------------------------------------------
48 !
49 !++ Private procedure
50 !
51
52 !-----------------------------------------------------------------------------
53 !
54 !++ Private parameters & variables
55 !
56 private :: gen_jacobigaussquadraturepts
57
58 !-----------------------------------------------------------------------------
59
60contains
61 !> A function to obtain the Lagrange basis functions related to the Gauss-Legendre-Lobatto (GLL) points
62 !!
63 !! @param Nord Order of Lagrange polynomial
64 !! @param x_lgl Positions of GLL points
65 !! @param x Positions where the Lagrange basis functions are evaluated
66!OCL SERIAL
67 function polynominal_genlagrangepoly(Nord, x_lgl, x) result(l)
68 implicit none
69
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)
74
75 integer :: n, i
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)
79
80 ! real(RP) :: w(Nord+1)
81 ! real(RP) :: int_w(Nord+1)
82 !---------------------------------------------------------------------------
83
84 p_lgl(:,:) = polynominal_genlegendrepoly(nord, x_lgl)
85 p(:,:) = polynominal_genlegendrepoly(nord, x)
86 pr(:,:) = polynominal_gendlegendrepoly(nord, x, p)
87
88 do n=1, nord+1
89 do i=1, size(x)
90 if ( abs(x(i)-x_lgl(n)) < 1e-15_rp ) then
91 l(i,n) = 1.0_rp
92 else
93 l(i,n) = &
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)))
96 end if
97 end do
98 end do
99
100 !- Calculate interpolation coefficient based on barycentric form
101 ! Eq. (3.46) in Wang, Huybrechs & Vandewalle (2012): Explicit barycentric weights for polynomial interpolation in the roots or extrema of classical orthogonal polynomials
102 ! int_w(:) = Polynominal_GenGaussLobattoPtIntWeight( Nord )
103 ! do n=1, Nord+1
104 ! w(n) = (-1)**mod(n-1,2) * sqrt(int_w(n))
105 ! end do
106 ! do n=1, Nord+1
107 ! do i=1, size(x)
108 ! if ( abs(x(i)-x_lgl(n)) < 1E-16_RP ) then
109 ! l(i,n) = 1.0_RP
110 ! else
111 ! l(i,n) = &
112 ! ( w(n) / ( x(i) - x_lgl(n) ) ) &
113 ! / sum( w(:) / ( x(i) - x_lgl(:) ) )
114 ! end if
115 ! end do
116 ! end do
117
118 return
119 end function polynominal_genlagrangepoly
120
121 !> A function to obtain the differential values of Lagrange basis functions at the GLL points
122 !!
123 !! @param Nord Order of Lagrange polynomial
124 !! @param x_lgl Positions of GLL points
125!OCL SERIAL
126 function polynominal_gendlagrangepoly_lglpt(Nord, x_lgl) result(lr)
127 implicit none
128
129 integer, intent(in) :: nord
130 real(rp), intent(in) :: x_lgl(nord+1)
131 real(rp) :: lr(nord+1, nord+1)
132
133 integer :: n, k
134 real(rp) :: p(nord+1,nord+1)
135 real(rp) :: lr_nn
136 !---------------------------------------------------------------------------
137
138 p(:,:) = polynominal_genlegendrepoly(nord, x_lgl)
139
140 do n=1, nord+1
141 lr_nn = 0.0_rp
142 do k=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))
147 else if (k==n) then
148 lr(k,n) = 0.0_rp
149 else
150 lr(k,n) = p(n,nord+1)/(p(k,nord+1)*(x_lgl(n) - x_lgl(k)))
151 end if
152
153 if ( k /= n ) then
154 lr_nn = lr_nn + lr(k,n)
155 end if
156 end do
157 lr(n,n) = - lr_nn
158 end do
159
160 return
162
163 !> A function to obtain the values of Legendre polynomials which are evaluated at arbitrary points.
164 !!
165 !! @param Nord Order of Lagrange polynomial
166 !! @param x Positions where the Legendre polynomials are evaluated
167 !! @param P Values of the Legendre polynomials at x
168!OCL SERIAL
169 subroutine polynominal_genlegendrepoly_sub(Nord, x, P)
170 implicit none
171
172 integer, intent(in) :: nord
173 real(rp), intent(in) :: x(:)
174 real(rp), intent(out) :: p(size(x), nord+1)
175
176 integer :: n
177 !---------------------------------------------------------------------------
178
179 if (nord < 0) then
180 log_error("Polynominal_GenLegendrePoly",*) "Nord must be larger than 0."
181 end if
182
183 p(:,1) = 1.0_rp
184 if (nord==0) return
185
186 p(:,2) = x(:)
187 do n=2, nord
188 p(:,n+1) = ( dble(2*n-1)*x(:)*p(:,n) - dble(n-1)*p(:,n-1) )/dble(n)
189 end do
190
191 return
193
194 !> A function to obtain the values of Legendre polynomials which are evaluated at arbitrary points.
195 !!
196 !! @param Nord Order of Lagrange polynomial
197 !! @param x Positions where the Legendre polynomials are evaluated
198 !! @param P Values of the Legendre polynomials at x
199!OCL SERIAL
200 function polynominal_genlegendrepoly(Nord, x) result(P)
201 implicit none
202
203 integer, intent(in) :: nord
204 real(rp), intent(in) :: x(:)
205 real(rp) :: p(size(x), nord+1)
206 !---------------------------------------------------------------------------
207
208 call polynominal_genlegendrepoly_sub( nord, x(:), & ! (in)
209 p(:,:) ) ! (out)
210 return
211 end function polynominal_genlegendrepoly
212
213 !> A function to obtain differential values of Legendre polynomials which are evaluated at arbitrary points.
214 !!
215 !! @param Nord Order of Lagrange polynomial
216 !! @param x Positions where the Legendre polynomials are evaluated
217 !! @param P Values of the Legendre polynomials at x
218!OCL SERIAL
219 function polynominal_gendlegendrepoly(Nord, x, P) result(GradP)
220 implicit none
221
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)
226
227 integer :: n
228 !---------------------------------------------------------------------------
229
230 if (nord < 0) then
231 log_error("Polynominal_GenDLegendrePoly",*) "Nord must be larger than 0."
232 end if
233
234 gradp(:,1) = 0.0_rp
235 if (nord == 0) return
236
237 gradp(:,2) = 1.0_rp
238 do n=2, nord
239 gradp(:,n+1) = 2.0_rp*x(:)*gradp(:,n) - gradp(:,n-1) + p(:,n)
240 end do
241
242 return
244
245 !> A function to calculate the Legendre-Gauss-Lobatto (LGL) points.
246 !!
247 !! @param Nord Order of Lagrange polynomial
248 !! @param pts Position of the LGL points
249!OCL SERIAL
250 function polynominal_gengausslobattopt(Nord) result(pts)
251 implicit none
252
253 integer, intent(in) :: nord
254 real(rp) :: pts(nord+1)
255
256 integer :: n1
257 !---------------------------------------------------------------------------
258
259 pts(1) = -1.0_rp; pts(nord+1) = 1.0_rp
260 if (nord==1) return
261
262 call gen_jacobigaussquadraturepts( 1, 1, nord-2, pts(2:nord) )
263 return
265
266 !> A function to calculate the Gauss-Lobbato weights.
267 !!
268 !! @param Nord Order of Lagrange polynomial
269 !! @param int_weight_lgl Gauss-Lobbato weights
270!OCL SERIAL
271 function polynominal_gengausslobattoptintweight(Nord) result(int_weight_lgl)
272 implicit none
273
274 integer, intent(in) :: nord
275 real(rp) :: int_weight_lgl(nord+1)
276
277 real(rp) :: lglpts1d(nord+1)
278 real(rp) :: p1d_ori(nord+1, nord+1)
279 !---------------------------------------------------------------------------
280
281 lglpts1d(:) = polynominal_gengausslobattopt( nord )
282 p1d_ori(:,:) = polynominal_genlegendrepoly( nord, lglpts1d )
283
284 int_weight_lgl(:) = 2.0_rp/(dble(nord*(nord+1))*p1d_ori(:,nord+1)**2)
285
286 return
288
289 !> A function to calculate the Gauss-Legendre (GL) points.
290 !!
291 !! @param Nord Order of the Legendre polynomial
292 !! @param pts Position of the GL points
293!OCL SERIAL
294 function polynominal_gengausslegendrept(Nord) result(pts)
295 implicit none
296
297 integer, intent(in) :: nord
298 real(rp) :: pts(nord)
299 !---------------------------------------------------------------------------
300
301 call gen_jacobigaussquadraturepts( 0, 0, nord-1, pts(:) )
302 return
304
305 !> A function to calculate the Gauss-Legendre (GL) weights.
306 !!
307 !! @param Nord Order of the Legendre polynomial
308 !! @param int_weight_gl Gauss-Legendre weights
309!OCL SERIAL
310 function polynominal_gengausslegendreptintweight(Nord) result(int_weight_gl)
311 implicit none
312
313 integer, intent(in) :: nord
314 real(rp) :: int_weight_gl(nord)
315
316 real(rp) :: glpts1d(nord)
317 real(rp) :: p1d_ori(nord, nord+1)
318 real(rp) :: dp1d_ori(nord, nord+1)
319 !---------------------------------------------------------------------------
320
321 glpts1d(:) = polynominal_gengausslegendrept( nord )
322 p1d_ori(:,:) = polynominal_genlegendrepoly( nord, glpts1d )
323 dp1d_ori(:,:) = polynominal_gendlegendrepoly( nord, glpts1d, p1d_ori )
324
325 int_weight_gl(:) = 2.0_rp / ( (1.0_rp - glpts1d(:)**2) * dp1d_ori(:,nord+1)**2 )
326
327 return
329
330 !- private -------------------------------
331
332 !> Calculate the N'th-order Gauss quadrature points and weights associated the Jacobi polynomial of type (alpha,beta).
333!OCL SERIAL
334 subroutine gen_jacobigaussquadraturepts( alpha, beta, N, &
335 x )
336
337 implicit none
338 integer, intent(in) :: alpha
339 integer, intent(in) :: beta
340 integer, intent(in) :: n
341 real(rp), intent(out) :: x(n+1)
342
343 integer :: i
344 real(dp) :: d(n+1), e(n)
345 real(dp) :: work(2*(n+1)-2), z(n+1,n+1)
346 real(dp) :: h1(n+1)
347 integer :: info
348 !--------------------------------------------------------------
349
350 if (n==0) then
351 x(1) = - dble(alpha - beta) / dble(alpha + beta + 2)
352 return
353 end if
354
355 do i=0, n
356 h1(i+1) = dble( 2*i + alpha + beta )
357 end do
358
359 do i=1, n+1
360 d(i) = - dble(alpha**2 - beta**2) / (h1(i) * (h1(i) + 2d0))
361 end do
362 do i=1, n
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)) )
366 end do
367 if ( dble(alpha + beta) < 1d-16) d(1) = 0.0_rp
368
369 call dstev( 'Vectors', n+1, d, e, z, n+1, work, info )
370 x(:) = d(:)
371
372 return
373 end subroutine gen_jacobigaussquadraturepts
374
375end module scale_polynominal
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.