FE-Project
Loading...
Searching...
No Matches
scale_polynominal.F90
Go to the documentation of this file.
1
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
63!OCL SERIAL
64 function polynominal_genlagrangepoly(Nord, x_lgl, x) result(l)
65 implicit none
66
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)
71
72 integer :: n, i
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)
76 !---------------------------------------------------------------------------
77
78 p_lgl(:,:) = polynominal_genlegendrepoly(nord, x_lgl)
79 p(:,:) = polynominal_genlegendrepoly(nord, x)
80 pr(:,:) = polynominal_gendlegendrepoly(nord, x, p)
81
82 do n=1, nord+1
83 do i=1, size(x)
84 if ( abs(x(i)-x_lgl(n)) < 1e-15_rp ) then
85 l(i,n) = 1.0_rp
86 else
87 l(i,n) = &
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)))
90 end if
91 end do
92 end do
93
94 return
96
99!OCL SERIAL
100 function polynominal_gendlagrangepoly_lglpt(Nord, x_lgl) result(lr)
101 implicit none
102
103 integer, intent(in) :: nord
104 real(rp), intent(in) :: x_lgl(nord+1)
105 real(rp) :: lr(nord+1, nord+1)
106
107 integer :: n, k
108 real(rp) :: p(nord+1,nord+1)
109 real(rp) :: lr_nn
110 !---------------------------------------------------------------------------
111
112 p(:,:) = polynominal_genlegendrepoly(nord, x_lgl)
113
114 do n=1, nord+1
115 lr_nn = 0.0_rp
116 do k=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))
121 else if (k==n) then
122 lr(k,n) = 0.0_rp
123 else
124 lr(k,n) = p(n,nord+1)/(p(k,nord+1)*(x_lgl(n) - x_lgl(k)))
125 end if
126
127 if ( k /= n ) then
128 lr_nn = lr_nn + lr(k,n)
129 end if
130 end do
131 lr(n,n) = - lr_nn
132 end do
133
134 return
136
139!OCL SERIAL
140 subroutine polynominal_genlegendrepoly_sub(Nord, x, P)
141 implicit none
142
143 integer, intent(in) :: nord
144 real(rp), intent(in) :: x(:)
145 real(rp), intent(out) :: p(size(x), nord+1)
146
147 integer :: n
148 !---------------------------------------------------------------------------
149
150 if (nord < 0) then
151 log_error("Polynominal_GenLegendrePoly",*) "Nord must be larger than 0."
152 end if
153
154 p(:,1) = 1.0_rp
155 if (nord==0) return
156
157 p(:,2) = x(:)
158 do n=2, nord
159 p(:,n+1) = ( dble(2*n-1)*x(:)*p(:,n) - dble(n-1)*p(:,n-1) )/dble(n)
160 end do
161
162 return
164
167!OCL SERIAL
168 function polynominal_genlegendrepoly(Nord, x) result(P)
169 implicit none
170
171 integer, intent(in) :: nord
172 real(rp), intent(in) :: x(:)
173 real(rp) :: p(size(x), nord+1)
174 !---------------------------------------------------------------------------
175
176 call polynominal_genlegendrepoly_sub( nord, x(:), & ! (in)
177 p(:,:) ) ! (out)
178 return
179 end function polynominal_genlegendrepoly
180
183!OCL SERIAL
184 function polynominal_gendlegendrepoly(Nord, x, P) result(GradP)
185 implicit none
186
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)
191
192 integer :: n
193 !---------------------------------------------------------------------------
194
195 if (nord < 0) then
196 log_error("Polynominal_GenDLegendrePoly",*) "Nord must be larger than 0."
197 end if
198
199 gradp(:,1) = 0.0_rp
200 if (nord == 0) return
201
202 gradp(:,2) = 1.0_rp
203 do n=2, nord
204 gradp(:,n+1) = 2.0_rp*x(:)*gradp(:,n) - gradp(:,n-1) + p(:,n)
205 end do
206
207 return
209
212!OCL SERIAL
213 function polynominal_gengausslobattopt(Nord) result(pts)
214 implicit none
215
216 integer, intent(in) :: nord
217 real(rp) :: pts(nord+1)
218
219 integer :: n1
220 !---------------------------------------------------------------------------
221
222 pts(1) = -1.0_rp; pts(nord+1) = 1.0_rp
223 if (nord==1) return
224
225 call gen_jacobigaussquadraturepts( 1, 1, nord-2, pts(2:nord) )
226 return
228
231!OCL SERIAL
232 function polynominal_gengausslobattoptintweight(Nord) result(int_weight_lgl)
233 implicit none
234
235 integer, intent(in) :: nord
236 real(rp) :: int_weight_lgl(nord+1)
237
238 real(rp) :: lglpts1d(nord+1)
239 real(rp) :: p1d_ori(nord+1, nord+1)
240 !---------------------------------------------------------------------------
241
242 lglpts1d(:) = polynominal_gengausslobattopt( nord )
243 p1d_ori(:,:) = polynominal_genlegendrepoly( nord, lglpts1d )
244
245 int_weight_lgl(:) = 2.0_rp/(dble(nord*(nord+1))*p1d_ori(:,nord+1)**2)
246
247 return
249
252!OCL SERIAL
253 function polynominal_gengausslegendrept(Nord) result(pts)
254 implicit none
255
256 integer, intent(in) :: nord
257 real(rp) :: pts(nord)
258 !---------------------------------------------------------------------------
259
260 call gen_jacobigaussquadraturepts( 0, 0, nord-1, pts(:) )
261 return
263
266!OCL SERIAL
267 function polynominal_gengausslegendreptintweight(Nord) result(int_weight_gl)
268 implicit none
269
270 integer, intent(in) :: nord
271 real(rp) :: int_weight_gl(nord)
272
273 real(rp) :: glpts1d(nord)
274 real(rp) :: p1d_ori(nord, nord+1)
275 real(rp) :: dp1d_ori(nord, nord+1)
276 !---------------------------------------------------------------------------
277
278 glpts1d(:) = polynominal_gengausslegendrept( nord )
279 p1d_ori(:,:) = polynominal_genlegendrepoly( nord, glpts1d )
280 dp1d_ori(:,:) = polynominal_gendlegendrepoly( nord, glpts1d, p1d_ori )
281
282 int_weight_gl(:) = 2.0_rp / ( (1.0_rp - glpts1d(:)**2) * dp1d_ori(:,nord+1)**2 )
283
284 return
286
287 !- private -------------------------------
288
290!OCL SERIAL
291 subroutine gen_jacobigaussquadraturepts( alpha, beta, N, &
292 x )
293
294 implicit none
295 integer, intent(in) :: alpha
296 integer, intent(in) :: beta
297 integer, intent(in) :: n
298 real(rp), intent(out) :: x(n+1)
299
300 integer :: i
301 real(dp) :: d(n+1), e(n)
302 real(dp) :: work(2*(n+1)-2), z(n+1,n+1)
303 real(dp) :: h1(n+1)
304 integer :: info
305 !--------------------------------------------------------------
306
307 if (n==0) then
308 x(1) = - dble(alpha - beta) / dble(alpha + beta + 2)
309 return
310 end if
311
312 do i=0, n
313 h1(i+1) = dble( 2*i + alpha + beta )
314 end do
315
316 do i=1, n+1
317 d(i) = - dble(alpha**2 - beta**2) / (h1(i) * (h1(i) + 2d0))
318 end do
319 do i=1, n
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)) )
323 end do
324 if ( dble(alpha + beta) < 1d-16) d(1) = 0.0_rp
325
326 call dstev( 'Vectors', n+1, d, e, z, n+1, work, info )
327 x(:) = d(:)
328
329 return
330 end subroutine gen_jacobigaussquadraturepts
331
332end module scale_polynominal
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.