FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
hecmw_solver_CG.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
5 
6 !C
7 !C***
8 !C*** module hecmw_solver_CG
9 !C***
10 !C
12 
13  public :: hecmw_solve_cg
14 
15 contains
16  !C
17  !C*** CG_nn
18  !C
19  subroutine hecmw_solve_cg( hecMESH, hecMAT, ITER, RESID, error, &
20  & Tset, Tsol, Tcomm )
21 
22  use hecmw_util
24  use m_hecmw_comm_f
29  use hecmw_precond
30  use hecmw_jad_type
32 
33  implicit none
34 
35  type(hecmwst_local_mesh) :: hecmesh
36  type(hecmwst_matrix) :: hecmat
37  integer(kind=kint ), intent(inout):: iter, error
38  real (kind=kreal), intent(inout):: resid, tset, tsol, tcomm
39 
40  integer(kind=kint ) :: n, np, ndof, nndof
41  integer(kind=kint ) :: my_rank
42  integer(kind=kint ) :: iterlog, timelog
43  real(kind=kreal), pointer :: b(:), x(:)
44 
45  real(kind=kreal), dimension(:,:), allocatable :: ww
46 
47  integer(kind=kint), parameter :: r= 1
48  integer(kind=kint), parameter :: z= 2
49  integer(kind=kint), parameter :: q= 2
50  integer(kind=kint), parameter :: p= 3
51  integer(kind=kint), parameter :: wk= 4
52 
53  integer(kind=kint ) :: maxit
54 
55  ! local variables
56  real (kind=kreal) :: tol
57  integer(kind=kint )::i
58  real (kind=kreal)::s_time,s1_time,e_time,e1_time, start_time, end_time
59  real (kind=kreal)::bnrm2
60  real (kind=kreal)::rho,rho1,beta,c1,alpha,dnrm2
61  real (kind=kreal)::t_max,t_min,t_avg,t_sd
62  integer(kind=kint) ::estcond
63  real (kind=kreal), allocatable::d(:),e(:)
64  real (kind=kreal)::alpha1
65  integer(kind=kint) :: n_indef_precond
66 
67  integer(kind=kint), parameter :: n_iter_recompute_r= 50
68 
69  call hecmw_barrier(hecmesh)
70  s_time= hecmw_wtime()
71 
72  !C===
73  !C +-------+
74  !C | INIT. |
75  !C +-------+
76  !C===
77  n = hecmat%N
78  np = hecmat%NP
79  ndof = hecmat%NDOF
80  nndof = n * ndof
81  my_rank = hecmesh%my_rank
82  x => hecmat%X
83  b => hecmat%B
84 
85  iterlog = hecmw_mat_get_iterlog( hecmat )
86  timelog = hecmw_mat_get_timelog( hecmat )
87  maxit = hecmw_mat_get_iter( hecmat )
88  tol = hecmw_mat_get_resid( hecmat )
89  estcond = hecmw_mat_get_estcond( hecmat )
90 
91  error = 0
92  n_indef_precond = 0
93  rho1 = 0.0d0
94  alpha1 = 0.0d0
95  beta = 0.0d0
96 
97  allocate (ww(ndof*np, 4))
98  ww = 0.d0
99 
100  !C
101  !C-- SCALING
102  call hecmw_solver_scaling_fw(hecmesh, hecmat, tcomm)
103 
104  !C
105  !C-- matrix integration for OpenACC
106  !C
107  !C @note:
108  !C Combine hecMAT%AL, D, and AU into a single matrix for GPU execution.
109  !C This is a no-op for CPU builds.
110  call hecmw_mat_integrate(hecmat)
111 
112  if (hecmw_mat_get_usejad(hecmat).ne.0) then
113  call hecmw_jad_init(hecmat)
114  endif
115 
116  if (estcond /= 0 .and. hecmesh%my_rank == 0) then
117  allocate(d(maxit),e(maxit-1))
118  endif
119  !C===
120  !C +----------------------+
121  !C | SETUP PRECONDITIONER |
122  !C +----------------------+
123  !C===
124  call hecmw_precond_setup(hecmat, hecmesh, 1)
125 
126  !C===
127  !C +---------------------+
128  !C | {r0}= {b} - [A]{x0} |
129  !C +---------------------+
130  !C===
131  call hecmw_matresid(hecmesh, hecmat, x, b, ww(:,r), tcomm)
132 
133  !C-- compute ||{b}||
134  call hecmw_innerproduct_r(hecmesh, ndof, b, b, bnrm2, tcomm)
135  if (bnrm2.eq.0.d0) then
136  iter = 0
137  maxit = 0
138  resid = 0.d0
139  x = 0.d0
140  endif
141 
142  e_time = hecmw_wtime()
143  if (timelog.eq.2) then
144  call hecmw_time_statistics(hecmesh, e_time - s_time, &
145  t_max, t_min, t_avg, t_sd)
146  if (hecmesh%my_rank.eq.0) then
147  write(*,*) 'Time solver setup'
148  write(*,*) ' Max :',t_max
149  write(*,*) ' Min :',t_min
150  write(*,*) ' Avg :',t_avg
151  write(*,*) ' Std Dev :',t_sd
152  endif
153  tset = t_max
154  else
155  tset = e_time - s_time
156  endif
157 
158  tcomm = 0.d0
159  call hecmw_barrier(hecmesh)
160  s1_time = hecmw_wtime()
161  !C
162  !C************************************************* Conjugate Gradient Iteration start
163  !C
164  do iter = 1, maxit
165 
166  !C===
167  !C +----------------+
168  !C | {z}= [Minv]{r} |
169  !C +----------------+
170  !C===
171  call hecmw_precond_apply(hecmesh, hecmat, ww(:,r), ww(:,z), ww(:,wk), tcomm)
172 
173 
174  !C===
175  !C +---------------+
176  !C | {RHO}= {r}{z} |
177  !C +---------------+
178  !C===
179  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,r), ww(:,z), rho, tcomm)
180  ! if RHO is NaN or Inf then no converge
181  if (rho == 0.d0) then
182  ! converged due to RHO==0
183  exit
184  elseif (iter > 1 .and. rho*rho1 <= 0) then
185  n_indef_precond = n_indef_precond + 1
186  if (n_indef_precond >= 3) then
187  ! diverged due to indefinite preconditioner
189  exit
190  endif
191  elseif (rho /= rho) then ! RHO is NaN
193  exit
194  endif
195 
196  !C===
197  !C +-----------------------------+
198  !C | {p} = {z} if ITER=1 |
199  !C | BETA= RHO / RHO1 otherwise |
200  !C +-----------------------------+
201  !C===
202  if ( iter.eq.1 ) then
203  call hecmw_copy_r(nndof, ww(:,z), ww(:,p))
204  else
205  beta = rho / rho1
206  call hecmw_xpay_r(nndof, beta, ww(:,z), ww(:,p))
207  endif
208 
209  !C===
210  !C +--------------+
211  !C | {q}= [A] {p} |
212  !C +--------------+
213  !C===
214  call hecmw_matvec(hecmesh, hecmat, ww(:,p), ww(:,q), tcomm)
215 
216  !C===
217  !C +---------------------+
218  !C | ALPHA= RHO / {p}{q} |
219  !C +---------------------+
220  !C===
221  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,p), ww(:,q), c1, tcomm)
222  ! if C1 is NaN or Inf, no converge
223  if (c1 <= 0) then
224  ! diverged due to indefinite or negative definite matrix
226  exit
227  elseif (c1 /= c1) then ! C1 is NaN
229  exit
230  endif
231 
232  alpha= rho / c1
233 
234  !C===
235  !C +----------------------+
236  !C | {x}= {x} + ALPHA*{p} |
237  !C | {r}= {r} - ALPHA*{q} |
238  !C +----------------------+
239  !C===
240  call hecmw_axpy_r(nndof, alpha, ww(:,p), x)
241 
242  if ( mod(iter,n_iter_recompute_r)==0 ) then
243  call hecmw_matresid(hecmesh, hecmat, x, b, ww(:,r), tcomm)
244  else
245  call hecmw_axpy_r(nndof, -alpha, ww(:,q), ww(:,r))
246  endif
247 
248  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,r), ww(:,r), dnrm2, tcomm)
249 
250  resid= dsqrt(dnrm2/bnrm2)
251 
252  !C##### ITERATION HISTORY
253  if (my_rank.eq.0.and.iterlog.eq.1) write (*,'(i7, 1pe16.6)') iter, resid
254  !C#####
255 
256  if (estcond /= 0 .and. hecmesh%my_rank == 0) then
257  if (iter == 1) then
258  d(1) = 1.d0 / alpha
259  else
260  d(iter) = 1.d0 / alpha + beta / alpha1
261  e(iter-1) = sqrt(beta) / alpha1
262  endif
263  alpha1 = alpha
264  if (mod(iter,estcond) == 0) call hecmw_estimate_condition_cg(iter, d, e)
265  endif
266 
267  if ( resid.le.tol ) then
268  if ( mod(iter,n_iter_recompute_r)==0 ) exit
269  !C----- recompute R to make sure it is really converged
270  call hecmw_matresid(hecmesh, hecmat, x, b, ww(:,r), tcomm)
271  call hecmw_innerproduct_r(hecmesh, ndof, ww(:,r), ww(:,r), dnrm2, tcomm)
272  resid= dsqrt(dnrm2/bnrm2)
273  if ( resid.le.tol ) exit
274  endif
275  if ( iter .eq.maxit ) error = hecmw_solver_error_noconv_maxit
276 
277  rho1 = rho
278 
279  enddo
280  !C
281  !C************************************************* Conjugate Gradient Iteration end
282  !C
283  call hecmw_solver_scaling_bk(hecmat)
284  !C
285  !C-- INTERFACE data EXCHANGE
286  !C
287  start_time= hecmw_wtime()
288  call hecmw_update_r (hecmesh, x, hecmat%NP, hecmat%NDOF)
289  end_time = hecmw_wtime()
290  tcomm = tcomm + end_time - start_time
291 
292  deallocate (ww)
293  !call hecmw_precond_clear(hecMAT)
294 
295  if (hecmw_mat_get_usejad(hecmat).ne.0) then
296  call hecmw_jad_finalize(hecmat)
297  endif
298 
299  if (estcond /= 0 .and. error == 0 .and. hecmesh%my_rank == 0) then
300  call hecmw_estimate_condition_cg(iter, d, e)
301  deallocate(d, e)
302  endif
303 
304  e1_time = hecmw_wtime()
305  if (timelog.eq.2) then
306  call hecmw_time_statistics(hecmesh, e1_time - s1_time, &
307  t_max, t_min, t_avg, t_sd)
308  if (hecmesh%my_rank.eq.0) then
309  write(*,*) 'Time solver iterations'
310  write(*,*) ' Max :',t_max
311  write(*,*) ' Min :',t_min
312  write(*,*) ' Avg :',t_avg
313  write(*,*) ' Std Dev :',t_sd
314  endif
315  tsol = t_max
316  else
317  tsol = e1_time - s1_time
318  endif
319 
320  end subroutine hecmw_solve_cg
321 
322 end module hecmw_solver_cg
subroutine, public hecmw_estimate_condition_cg(ITER, D, E)
Jagged Diagonal Matrix storage for vector processors. Original code was provided by JAMSTEC.
Definition: hecmw_jadm.f90:8
subroutine, public hecmw_jad_init(hecMAT)
Definition: hecmw_jadm.f90:29
subroutine, public hecmw_jad_finalize(hecMAT)
Definition: hecmw_jadm.f90:42
subroutine, public hecmw_mat_integrate(hecMAT)
Integrate matrix components into a single array for efficient access.
real(kind=kreal) function, public hecmw_mat_get_resid(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_iterlog(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_timelog(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_estcond(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_usejad(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_iter(hecMAT)
subroutine, public hecmw_precond_setup(hecMAT, hecMESH, sym)
subroutine, public hecmw_precond_apply(hecMESH, hecMAT, R, Z, ZP, COMMtime)
subroutine, public hecmw_solve_cg(hecMESH, hecMAT, ITER, RESID, error, Tset, Tsol, Tcomm)
subroutine, public hecmw_matresid(hecMESH, hecMAT, X, B, R, COMMtime)
subroutine, public hecmw_matvec(hecMESH, hecMAT, X, Y, COMMtime)
subroutine hecmw_xpay_r(n, alpha, X, Y)
subroutine hecmw_innerproduct_r(hecMESH, ndof, X, Y, sum, COMMtime)
subroutine hecmw_axpy_r(n, alpha, X, Y)
subroutine hecmw_copy_r(n, X, Y)
subroutine hecmw_time_statistics(hecMESH, time, t_max, t_min, t_avg, t_sd)
subroutine, public hecmw_solver_scaling_fw(hecMESH, hecMAT, COMMtime)
subroutine, public hecmw_solver_scaling_bk(hecMAT)
I/O and Utility.
Definition: hecmw_util_f.F90:7
integer(kind=4), parameter kreal
real(kind=kreal) function hecmw_wtime()
subroutine hecmw_update_r(hecMESH, val, n, m)
subroutine hecmw_barrier(hecMESH)
integer(kind=kint), parameter hecmw_solver_error_diverge_pc
integer(kind=kint), parameter hecmw_solver_error_diverge_nan
integer(kind=kint), parameter hecmw_solver_error_noconv_maxit
integer(kind=kint), parameter hecmw_solver_error_diverge_mat