FE-Project
Loading...
Searching...
No Matches
Functions/Subroutines
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 values of Lagrange basis functions which are evaluated over aribitary 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 which are evaluated over aribitary points.
 
subroutine, public polynominal_genlegendrepoly_sub (nord, x, p)
 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_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 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(nord+1), public polynominal_gengausslobattoptintweight (nord)
 A function to calcuate the Gauss-Lobbato weights.
 
real(rp) function, dimension(nord), public polynominal_gengausslegendrept (nord)
 A function to calcuate the Gauss-Legendre points.
 
real(rp) function, dimension(nord), public polynominal_gengausslegendreptintweight (nord)
 A function to calcuate the Gauss-Legendre weights.
 

Detailed Description

module common / Polynominal

Description
A module to provide utilities for polynominal
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 values of Lagrange basis functions which are evaluated over aribitary points.

Definition at line 64 of file scale_polynominal.F90.

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

References polynominal_gendlegendrepoly(), and polynominal_genlegendrepoly().

Referenced by scale_element_hexahedral::hexhedralelement_init(), scale_element_line::lineelement_init(), and scale_element_quadrilateral::quadrilateralelement_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 which are evaluated over aribitary points.

Definition at line 100 of file scale_polynominal.F90.

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

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 polynominals which are evaluated at aribitary points.

Definition at line 140 of file scale_polynominal.F90.

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

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 polynominals which are evaluated at aribitary points.

Definition at line 168 of file scale_polynominal.F90.

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

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(), and scale_element_quadrilateral::quadrilateralelement_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 polynominals which are evaluated at aribitary points.

Definition at line 184 of file scale_polynominal.F90.

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

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 calcuate the Legendre-Gauss-Lobtatto (LGL) points.

Definition at line 213 of file scale_polynominal.F90.

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

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 calcuate the Gauss-Lobbato weights.

Definition at line 232 of file scale_polynominal.F90.

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

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 calcuate the Gauss-Legendre points.

Definition at line 253 of file scale_polynominal.F90.

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

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 calcuate the Gauss-Legendre weights.

Definition at line 267 of file scale_polynominal.F90.

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

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().