FE-Project
Loading...
Searching...
No Matches
Data Types | Functions/Subroutines
scale_linalgebra Module Reference

module common / Linear algebra More...

Data Types

interface  linalgebra_solvelineq
 

Functions/Subroutines

real(rp) function, dimension(size(a, 1), size(a, 2)), public linalgebra_inv (a)
 
subroutine, public linalgebra_lu (a_lu, ipiv)
 
subroutine, public linalgebra_solvelineq_gmres (a, b, x, m, restart_num, conv_crit)
 

Detailed Description

module common / Linear algebra

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

Function/Subroutine Documentation

◆ linalgebra_inv()

real(rp) function, dimension(size(a,1),size(a,2)), public scale_linalgebra::linalgebra_inv ( real(rp), dimension(:,:), intent(in) a)

Definition at line 66 of file scale_linalgebra.F90.

67
68 real(RP), intent(in) :: A(:,:)
69 real(RP) :: Ainv(size(A,1),size(A,2))
70
71 real(RP) :: work(size(A,1))
72 integer :: ipiv(size(A,1))
73 integer :: n, info
74
75 !---------------------------------------------------------------------------
76
77 ainv(:,:) = a
78 n = size(a,1)
79
80 call dgetrf(n, n, ainv, n, ipiv, info)
81 if (info /=0 ) then
82 log_error("linalgebra_inv",*) "Matrix is singular"
83 call prc_abort
84 end if
85
86 call dgetri(n, ainv, n, ipiv, work, n, info)
87 if (info /=0 ) then
88 log_error("linalgebra_inv",*) "Matrix inversion is failed"
89 call prc_abort
90 end if
91
92 return

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

◆ linalgebra_lu()

subroutine, public scale_linalgebra::linalgebra_lu ( real(rp), dimension(:,:), intent(inout) a_lu,
integer, dimension(size(a_lu,1)), intent(out) ipiv )

Definition at line 160 of file scale_linalgebra.F90.

161
162 real(RP), intent(inout) :: A_lu(:,:)
163 integer, intent(out) :: ipiv(size(A_lu,1))
164
165 integer :: n
166 integer :: info
167
168 !---------------------------------------------------------------------------
169 n = size(a_lu,1)
170
171 call dgetrf(n, n, a_lu, n, ipiv, info)
172 if (info /=0 ) then
173 log_error("linalgebra_LU",*) "Matrix is singular"
174 call prc_abort
175 end if
176
177 return

Referenced by scale_atm_dyn_dgm_nonhydro3d_hevi_gmres::atm_dyn_dgm_nonhydro3d_hevi_cal_vi(), get_pmatd_lu(), and scale_atm_dyn_dgm_hydrostatic::hydrostaic_build_rho_xyz::hydrostaic_build_rho_xyz_moist().

◆ linalgebra_solvelineq_gmres()

subroutine, public scale_linalgebra::linalgebra_solvelineq_gmres ( type(sparsemat), intent(in) a,
real(rp), dimension(:), intent(in) b,
real(rp), dimension(:), intent(inout) x,
integer, intent(in) m,
integer, intent(in) restart_num,
real(rp), intent(in) conv_crit )

Definition at line 180 of file scale_linalgebra.F90.

181
183 implicit none
184
185 type(sparsemat), intent(in) :: A
186 real(RP), intent(in) :: b(:)
187 real(RP), intent(inout) :: x(:)
188 integer, intent(in) :: m
189 integer, intent(in) ::restart_num
190 real(RP), intent(in) :: CONV_CRIT
191
192 real(RP) :: x0(size(x))
193 real(RP) :: w(size(x))
194 real(RP) :: v(size(x),m+1)
195 real(RP) :: z(size(x),m+1)
196 real(RP) :: Av(size(x))
197 real(RP) :: g(size(x)+1)
198 real(RP) :: r0_l2, r0_l2_new
199 real(RP) :: h(m+1,m)
200 real(RP) :: r(m,m)
201 real(RP) :: c(m), s(m)
202 real(RP) :: y(m)
203 real(RP) :: tmp, tmp1, tmp2
204 integer :: i, j
205 integer :: restart_i
206
207 type(sparsemat) :: M_PC
208
209 !---------------------------------------------------------------------------
210
211 call precondstep_ilu0_constructmat(a, m_pc)
212 x0(:) = x
213
214 do restart_i=1, restart_num
215
216 if (restart_i == 1) then
217 call sparsemat_matmul(a, x, av)
218 v(:,1) = b - av
219 r0_l2 = sqrt(sum(v(:,1)**2))
220 end if
221 v(:,1) = v(:,1)/r0_l2
222 g(1) = r0_l2
223
224 do j=1, m
225
226 call precondstep_ilu0_solve(m_pc, v(:,j), z(:,j))
227
228 call sparsemat_matmul(a, z(:,j), w)
229
230 do i=1, j
231 h(i,j) = sum(w*v(:,i))
232 w(:) = w - h(i,j)*v(:,i)
233 end do
234 h(j+1,j) = sqrt(sum(w**2))
235 v(:,j+1) = w(:)/h(j+1,j)
236
237 r(1,j) = h(1,j)
238 do i=1, j-1
239 tmp1 = c(i)*r(i,j) + s(i)*h(i+1,j)
240 tmp2 = -s(i)*r(i,j) + c(i)*h(i+1,j)
241 r(i,j) = tmp1
242 r(i+1,j) = tmp2
243 end do
244
245 tmp = sqrt(r(j,j)**2 + h(j+1,j)**2)
246 c(j) = r(j,j)/tmp
247 s(j) = h(j+1,j)/tmp
248
249 g(j+1) = -s(j)*g(j)
250 g(j) = c(j)*g(j)
251 r(j,j) = c(j)*r(j,j) + s(j)*h(j+1,j)
252 end do
253
254 y(m) = g(m)/r(m,m)
255 do i=m-1,1,-1
256 y(i) = (g(i) - sum(r(i,i+1:m)*y(i+1:m)))/r(i,i)
257 end do
258
259 x(:) = x + matmul(z(:,1:m),y)
260
261 !write(*,*) "x:", x
262 !write(*,*) "x0:", x0
263
264 call sparsemat_matmul(a, x, av)
265 v(:,1) = b - av
266 r0_l2_new = sqrt(sum(v(:,1)**2))
267
268 !write(*,*) "r:", r0_l2, "->", r0_l2_new
269 if (r0_l2_new < conv_crit) exit
270 r0_l2 = r0_l2_new
271 end do
272
273 return
module common / sparsemat

References scale_sparsemat::sparsemat_storage_typeid_csr.