FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
hecmw_precond.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 
7  use hecmw_util
16  implicit none
17 
18  private
19  public :: hecmw_precond_setup
20  public :: hecmw_precond_clear
21  public :: hecmw_precond_apply
23  public :: hecmw_precond_get_timer
24 
25  real(kind=kreal) :: time_precond = 0.d0
26 
27 contains
28 
29  subroutine hecmw_precond_setup(hecMAT, hecMESH, sym)
30  implicit none
31  type (hecmwst_matrix), intent(inout) :: hecmat
32  type (hecmwst_local_mesh), intent(inout) :: hecmesh
33  integer(kind=kint) :: sym
34 
35  if (hecmw_mat_get_iterpremax( hecmat ).le.0) return
36 
37  select case(hecmw_mat_get_precond( hecmat ))
38  case(1,2)
39  call hecmw_precond_ssor_setup(hecmat)
40  case(3)
41  call hecmw_precond_diag_setup(hecmat)
42  case(5)
43  call hecmw_precond_ml_setup(hecmat, hecmesh, sym)
44  case(10,11,12)
45  call hecmw_precond_bilu_setup(hecmat)
46  case(20)
47  call hecmw_precond_sainv_setup(hecmat)
48  case(21)
49  call hecmw_precond_rif_setup(hecmat)
50  case default
51  write (*,'(/a )')'#### HEC-MW-SOLVER-E-1001'
52  write (*,'( a/)')' inconsistent solver/preconditioning'
54  end select
55  end subroutine hecmw_precond_setup
56 
57  subroutine hecmw_precond_clear(hecMAT)
58  implicit none
59  type (hecmwst_matrix), intent(inout) :: hecmat
60 
61  if (hecmw_mat_get_iterpremax( hecmat ).le.0) return
62 
63  select case(hecmw_mat_get_precond( hecmat ))
64  case(1,2)
65  call hecmw_precond_ssor_clear(hecmat)
66  case(3)
67  call hecmw_precond_diag_clear(hecmat%NDOF)
68  case(5)
69  call hecmw_precond_ml_clear(hecmat%NDOF)
70  case(10:12)
71  call hecmw_precond_bilu_clear(hecmat%NDOF)
72  case(20)
73  call hecmw_precond_sainv_clear(hecmat%NDOF)
74  case(21)
75  call hecmw_precond_rif_clear(hecmat%NDOF)
76  case default
77  end select
78 
79  end subroutine hecmw_precond_clear
80 
81  subroutine hecmw_precond_apply(hecMESH, hecMAT, R, Z, ZP, COMMtime)
82  implicit none
83  type (hecmwst_local_mesh), intent(inout) :: hecmesh
84  type (hecmwst_matrix), intent(inout) :: hecmat
85  real(kind=kreal), intent(inout) :: r(:)
86  real(kind=kreal), intent(inout) :: z(:), zp(:)
87  real(kind=kreal), intent(inout) :: commtime
88  integer(kind=kint ) :: i, n, np, nndof, npndof
89  integer(kind=kint) :: iterpremax, iterpre
90  real(kind=kreal) :: start_time, end_time
91 
92  start_time = hecmw_wtime()
93 
94  n = hecmat%N
95  np = hecmat%NP
96  nndof = n * hecmat%NDOF
97  npndof = np * hecmat%NDOF
98 
99  if (hecmw_mat_get_iterpremax( hecmat ).le.0) then
100  !$acc kernels
101  !$acc loop independent
102  do i= 1, nndof
103  z(i)= r(i)
104  enddo
105  !$acc end kernels
106  return
107  endif
108 
109  !C {z}= [Minv]{r}
110  !$acc kernels
111  !$acc loop independent
112  do i= 1, nndof
113  zp(i)= r(i)
114  enddo
115  !$acc end kernels
116 
117  !$acc kernels
118  !$acc loop independent
119  do i= nndof+1, npndof
120  zp(i) = 0.d0
121  enddo
122  !$acc end kernels
123 
124  !$acc kernels
125  !$acc loop independent
126  do i= 1, npndof
127  z(i)= 0.d0
128  enddo
129  !$acc end kernels
130 
131  iterpremax = hecmw_mat_get_iterpremax( hecmat )
132  do iterpre= 1, iterpremax
133 
134  select case(hecmw_mat_get_precond( hecmat ))
135  case(1,2)
136  call hecmw_precond_ssor_apply(zp,hecmat%NDOF)
137  case(3)
138  call hecmw_precond_diag_apply(zp,hecmat%NDOF)
139  case(5)
140  call hecmw_precond_ml_apply(zp,hecmat%NDOF)
141  case(10:12)
142  call hecmw_precond_bilu_apply(zp,hecmat%NDOF)
143  case(20)
144  call hecmw_precond_sainv_apply(r,zp,hecmat%NDOF)
145  case(21)
146  call hecmw_precond_rif_apply(zp,hecmat%NDOF)
147  case default
148  end select
149 
150  !C-- additive Schwartz
151  !$acc kernels
152  !$acc loop independent
153  do i= 1, hecmat%N * hecmat%NDOF
154  z(i)= z(i) + zp(i)
155  enddo
156  !$acc end kernels
157  if (iterpre.eq.iterpremax) exit
158 
159  !C-- {ZP} = {R} - [A] {Z}
160  call hecmw_matresid (hecmesh, hecmat, z, r, zp, commtime)
161  enddo
162 
163  end_time = hecmw_wtime()
164  time_precond = time_precond + end_time - start_time
165  end subroutine hecmw_precond_apply
166 
168  implicit none
169  time_precond = 0.d0
170  end subroutine hecmw_precond_clear_timer
171 
173  implicit none
174  real(kind=kreal) :: hecmw_precond_get_timer
175  hecmw_precond_get_timer = time_precond
176  end function hecmw_precond_get_timer
177 
178 end module hecmw_precond
integer(kind=kint) function, public hecmw_mat_get_iterpremax(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_precond(hecMAT)
subroutine, public hecmw_precond_bilu_setup(hecMAT)
subroutine, public hecmw_precond_bilu_apply(ZP, NDOF)
subroutine, public hecmw_precond_bilu_clear(NDOF)
subroutine, public hecmw_precond_diag_setup(hecMAT)
subroutine, public hecmw_precond_diag_clear(NDOF)
subroutine, public hecmw_precond_diag_apply(ZP, NDOF)
subroutine, public hecmw_precond_ml_clear(NDOF)
subroutine, public hecmw_precond_ml_setup(hecMAT, hecMESH, sym)
subroutine, public hecmw_precond_ml_apply(ZP, NDOF)
subroutine, public hecmw_precond_rif_clear(NDOF)
subroutine, public hecmw_precond_rif_apply(ZP, NDOF)
subroutine, public hecmw_precond_rif_setup(hecMAT)
subroutine, public hecmw_precond_sainv_setup(hecMAT)
subroutine, public hecmw_precond_sainv_clear(NDOF)
subroutine, public hecmw_precond_sainv_apply(R, ZP, NDOF)
subroutine, public hecmw_precond_ssor_apply(ZP, NDOF)
subroutine, public hecmw_precond_ssor_setup(hecMAT)
subroutine, public hecmw_precond_ssor_clear(hecMAT)
subroutine, public hecmw_precond_clear_timer
real(kind=kreal) function, public hecmw_precond_get_timer()
subroutine, public hecmw_precond_clear(hecMAT)
subroutine, public hecmw_precond_setup(hecMAT, hecMESH, sym)
subroutine, public hecmw_precond_apply(hecMESH, hecMAT, R, Z, ZP, COMMtime)
subroutine, public hecmw_matresid(hecMESH, hecMAT, X, B, R, COMMtime)
I/O and Utility.
Definition: hecmw_util_f.F90:7
integer(kind=kint) function hecmw_comm_get_comm()
integer(kind=4), parameter kreal
subroutine hecmw_abort(comm)
real(kind=kreal) function hecmw_wtime()