FE-Project
Loading...
Searching...
No Matches
scale_polynominal Module Reference

Module common / Polynominal. More...

Functions/Subroutines

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, 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.
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 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(nord+1), public polynominal_gengausslobattoptintweight (nord)
 A function to calculate the Gauss-Lobbato weights.
real(rp) function, dimension(nord), public polynominal_gengausslegendrept (nord)
 A function to calculate the Gauss-Legendre (GL) points.
real(rp) function, dimension(nord), public polynominal_gengausslegendreptintweight (nord)
 A function to calculate the Gauss-Legendre (GL) weights.

Detailed Description

Module common / Polynominal.

Description
A module to provide utilities for polynomials
Reference
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ polynominal_genlagrangepoly()

real(rp) function, dimension(size(x), nord+1), public scale_polynominal::polynominal_genlagrangepoly ( integer, intent(in) nord,
real(rp), dimension(nord+1), intent(in) x_lgl,
real(rp), dimension(:), intent(in) x )

A function to obtain the Lagrange basis functions related to the Gauss-Legendre-Lobatto (GLL) points.

Parameters
NordOrder of Lagrange polynomial
x_lglPositions of GLL points
xPositions where the Lagrange basis functions are evaluated

Definition at line 67 of file scale_polynominal.F90.

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

References polynominal_gendlegendrepoly(), and polynominal_genlegendrepoly().

Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), scale_element_quadrilateral::quadrilateralelement_init(), and scale_element_siacfilter::siac_filter_init().

◆ polynominal_gendlagrangepoly_lglpt()

real(rp) function, dimension(nord+1, nord+1), public scale_polynominal::polynominal_gendlagrangepoly_lglpt ( integer, intent(in) nord,
real(rp), dimension(nord+1), intent(in) x_lgl )

A function to obtain the differential values of Lagrange basis functions at the GLL points.

Parameters
NordOrder of Lagrange polynomial
x_lglPositions of GLL points

Definition at line 126 of file scale_polynominal.F90.

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

References polynominal_genlegendrepoly().

Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), and scale_element_quadrilateral::quadrilateralelement_init().

◆ polynominal_genlegendrepoly_sub()

subroutine, public scale_polynominal::polynominal_genlegendrepoly_sub ( integer, intent(in) nord,
real(rp), dimension(:), intent(in) x,
real(rp), dimension(size(x), nord+1), intent(out) p )

A function to obtain the values of Legendre polynomials which are evaluated at arbitrary points.

Parameters
NordOrder of Lagrange polynomial
xPositions where the Legendre polynomials are evaluated
PValues of the Legendre polynomials at x

Definition at line 169 of file scale_polynominal.F90.

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

Referenced by scale_file_common_meshfield::file_common_meshfield_put_field1d_cartesbuf(), scale_file_common_meshfield::file_common_meshfield_put_field2d_cartesbuf(), scale_file_common_meshfield::file_common_meshfield_put_field3d_cartesbuf(), and polynominal_genlegendrepoly().

◆ polynominal_genlegendrepoly()

real(rp) function, dimension(size(x), nord+1), public scale_polynominal::polynominal_genlegendrepoly ( integer, intent(in) nord,
real(rp), dimension(:), intent(in) x )

A function to obtain the values of Legendre polynomials which are evaluated at arbitrary points.

Parameters
NordOrder of Lagrange polynomial
xPositions where the Legendre polynomials are evaluated
PValues of the Legendre polynomials at x

Definition at line 200 of file scale_polynominal.F90.

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

References polynominal_genlegendrepoly_sub().

Referenced by scale_file_common_meshfield::file_common_meshfield_put_field2d_cubedsphere_cartesbuf(), scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), polynominal_gendlagrangepoly_lglpt(), polynominal_gengausslegendreptintweight(), polynominal_gengausslobattoptintweight(), polynominal_genlagrangepoly(), scale_element_quadrilateral::quadrilateralelement_init(), and scale_element_siacfilter::siac_filter_init().

◆ polynominal_gendlegendrepoly()

real(rp) function, dimension(size(x), nord+1), public scale_polynominal::polynominal_gendlegendrepoly ( integer, intent(in) nord,
real(rp), dimension(:), intent(in) x,
real(rp), dimension(:,:), intent(in) p )

A function to obtain differential values of Legendre polynomials which are evaluated at arbitrary points.

Parameters
NordOrder of Lagrange polynomial
xPositions where the Legendre polynomials are evaluated
PValues of the Legendre polynomials at x

Definition at line 219 of file scale_polynominal.F90.

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

Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), polynominal_gengausslegendreptintweight(), polynominal_genlagrangepoly(), and scale_element_quadrilateral::quadrilateralelement_init().

◆ polynominal_gengausslobattopt()

real(rp) function, dimension(nord+1), public scale_polynominal::polynominal_gengausslobattopt ( integer, intent(in) nord)

A function to calculate the Legendre-Gauss-Lobatto (LGL) points.

Parameters
NordOrder of Lagrange polynomial
ptsPosition of the LGL points

Definition at line 250 of file scale_polynominal.F90.

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

Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), polynominal_gengausslobattoptintweight(), and scale_element_quadrilateral::quadrilateralelement_init().

◆ polynominal_gengausslobattoptintweight()

real(rp) function, dimension(nord+1), public scale_polynominal::polynominal_gengausslobattoptintweight ( integer, intent(in) nord)

A function to calculate the Gauss-Lobbato weights.

Parameters
NordOrder of Lagrange polynomial
int_weight_lglGauss-Lobbato weights

Definition at line 271 of file scale_polynominal.F90.

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

References polynominal_gengausslobattopt(), and polynominal_genlegendrepoly().

Referenced by scale_atm_dyn_dgm_trcadvect3d_heve::atm_dyn_dgm_trcadvect3d_heve_init(), scale_atm_phy_mp_dgm_common::atm_phy_mp_dgm_common_gen_intweight(), scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), and scale_element_quadrilateral::quadrilateralelement_init().

◆ polynominal_gengausslegendrept()

real(rp) function, dimension(nord), public scale_polynominal::polynominal_gengausslegendrept ( integer, intent(in) nord)

A function to calculate the Gauss-Legendre (GL) points.

Parameters
NordOrder of the Legendre polynomial
ptsPosition of the GL points

Definition at line 294 of file scale_polynominal.F90.

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

Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), polynominal_gengausslegendreptintweight(), and scale_element_quadrilateral::quadrilateralelement_init().

◆ polynominal_gengausslegendreptintweight()

real(rp) function, dimension(nord), public scale_polynominal::polynominal_gengausslegendreptintweight ( integer, intent(in) nord)

A function to calculate the Gauss-Legendre (GL) weights.

Parameters
NordOrder of the Legendre polynomial
int_weight_glGauss-Legendre weights

Definition at line 310 of file scale_polynominal.F90.

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

References polynominal_gendlegendrepoly(), polynominal_gengausslegendrept(), and polynominal_genlegendrepoly().

Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), and scale_element_quadrilateral::quadrilateralelement_init().