FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_precond_SSOR_44.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_precond_SSOR_44
9 !C***
10 !C
12  use hecmw_util
17  !$ use omp_lib
18 
19  private
20 
24 
25  integer(kind=kint) :: N
26  real(kind=kreal), pointer :: d(:) => null()
27  real(kind=kreal), pointer :: al(:) => null()
28  real(kind=kreal), pointer :: au(:) => null()
29  integer(kind=kint), pointer :: indexL(:) => null()
30  integer(kind=kint), pointer :: indexU(:) => null()
31  integer(kind=kint), pointer :: itemL(:) => null()
32  integer(kind=kint), pointer :: itemU(:) => null()
33  real(kind=kreal), pointer :: alu(:) => null()
34 
35  integer(kind=kint) :: NColor
36  integer(kind=kint), pointer :: COLORindex(:) => null()
37  integer(kind=kint), pointer :: perm(:) => null()
38  integer(kind=kint), pointer :: iperm(:) => null()
39 
40  logical, save :: isFirst = .true.
41 
42  logical, save :: INITIALIZED = .false.
43 
44 contains
45 
46  subroutine hecmw_precond_ssor_44_setup(hecMAT)
47  implicit none
48  type(hecmwst_matrix), intent(inout) :: hecmat
49  integer(kind=kint ) :: npl, npu, npcl, npcu
50  real (kind=kreal), allocatable :: cd(:)
51  integer(kind=kint ) :: ncolor_in
52  real (kind=kreal) :: sigma_diag
53  real (kind=kreal) :: alutmp(4,4), pw(4)
54  integer(kind=kint ) :: ii, i, j, k
55  integer(kind=kint ) :: nthreads = 1
56  integer(kind=kint ), allocatable :: perm_tmp(:)
57  !real (kind=kreal) :: t0
58 
59  !t0 = hecmw_Wtime()
60  !write(*,*) 'DEBUG: SSOR setup start', hecmw_Wtime()-t0
61 
62  if (initialized) then
63  if (hecmat%Iarray(98) == 1) then ! need symbolic and numerical setup
64  call hecmw_precond_ssor_44_clear(hecmat)
65  else if (hecmat%Iarray(97) == 1) then ! need numerical setup only
66  call hecmw_precond_ssor_44_clear(hecmat) ! TEMPORARY
67  else
68  return
69  endif
70  endif
71 
72  !$ nthreads = omp_get_max_threads()
73 
74  n = hecmat%N
75  ! N = hecMAT%NP
76  ncolor_in = hecmw_mat_get_ncolor_in(hecmat)
77  sigma_diag = hecmw_mat_get_sigma_diag(hecmat)
78 
79  if (nthreads == 1) then
80  ncolor = 1
81  allocate(colorindex(0:1), perm(n), iperm(n))
82  colorindex(0) = 0
83  colorindex(1) = n
84  do i=1,n
85  perm(i) = i
86  iperm(i) = i
87  end do
88  else
89  allocate(colorindex(0:n), perm_tmp(n), perm(n), iperm(n))
90  call hecmw_matrix_ordering_rcm(n, hecmat%indexL, hecmat%itemL, &
91  hecmat%indexU, hecmat%itemU, perm_tmp, iperm)
92  !write(*,*) 'DEBUG: RCM ordering done', hecmw_Wtime()-t0
93  call hecmw_matrix_ordering_mc(n, hecmat%indexL, hecmat%itemL, &
94  hecmat%indexU, hecmat%itemU, perm_tmp, &
95  ncolor_in, ncolor, colorindex, perm, iperm)
96  !write(*,*) 'DEBUG: MC ordering done', hecmw_Wtime()-t0
97  deallocate(perm_tmp)
98 
99  !call write_debug_info
100  endif
101 
102  npl = hecmat%indexL(n)
103  npu = hecmat%indexU(n)
104  allocate(indexl(0:n), indexu(0:n), iteml(npl), itemu(npu))
105  call hecmw_matrix_reorder_profile(n, perm, iperm, &
106  hecmat%indexL, hecmat%indexU, hecmat%itemL, hecmat%itemU, &
107  indexl, indexu, iteml, itemu)
108  !write(*,*) 'DEBUG: reordering profile done', hecmw_Wtime()-t0
109 
110  !call check_ordering
111 
112  allocate(d(16*n), al(16*npl), au(16*npu))
113  call hecmw_matrix_reorder_values(n, 4, perm, iperm, &
114  hecmat%indexL, hecmat%indexU, hecmat%itemL, hecmat%itemU, &
115  hecmat%AL, hecmat%AU, hecmat%D, &
116  indexl, indexu, iteml, itemu, al, au, d)
117  !write(*,*) 'DEBUG: reordering values done', hecmw_Wtime()-t0
118 
119  call hecmw_matrix_reorder_renum_item(n, perm, indexl, iteml)
120  call hecmw_matrix_reorder_renum_item(n, perm, indexu, itemu)
121 
122  allocate(alu(16*n))
123  alu = 0.d0
124 
125  do ii= 1, 16*n
126  alu(ii) = d(ii)
127  enddo
128 
129  do ii= 1, n
130  alutmp(1,1)= alu(16*ii-15) * sigma_diag
131  alutmp(1,2)= alu(16*ii-14)
132  alutmp(1,3)= alu(16*ii-13)
133  alutmp(1,4)= alu(16*ii-12)
134  alutmp(2,1)= alu(16*ii-11)
135  alutmp(2,2)= alu(16*ii-10) * sigma_diag
136  alutmp(2,3)= alu(16*ii- 9)
137  alutmp(2,4)= alu(16*ii- 8)
138  alutmp(3,1)= alu(16*ii- 7)
139  alutmp(3,2)= alu(16*ii- 6)
140  alutmp(3,3)= alu(16*ii- 5) * sigma_diag
141  alutmp(3,4)= alu(16*ii- 4)
142  alutmp(4,1)= alu(16*ii- 3)
143  alutmp(4,2)= alu(16*ii- 2)
144  alutmp(4,3)= alu(16*ii- 1)
145  alutmp(4,4)= alu(16*ii ) * sigma_diag
146 
147  do k= 1, 4
148  alutmp(k,k)= 1.d0/alutmp(k,k)
149  do i= k+1, 4
150  alutmp(i,k)= alutmp(i,k) * alutmp(k,k)
151  do j= k+1, 4
152  pw(j)= alutmp(i,j) - alutmp(i,k)*alutmp(k,j)
153  enddo
154  do j= k+1, 4
155  alutmp(i,j)= pw(j)
156  enddo
157  enddo
158  enddo
159  alu(16*ii-15)= alutmp(1,1)
160  alu(16*ii-14)= alutmp(1,2)
161  alu(16*ii-13)= alutmp(1,3)
162  alu(16*ii-12)= alutmp(1,4)
163  alu(16*ii-11)= alutmp(2,1)
164  alu(16*ii-10)= alutmp(2,2)
165  alu(16*ii- 9)= alutmp(2,3)
166  alu(16*ii- 8)= alutmp(2,4)
167  alu(16*ii- 7)= alutmp(3,1)
168  alu(16*ii- 6)= alutmp(3,2)
169  alu(16*ii- 5)= alutmp(3,3)
170  alu(16*ii- 4)= alutmp(3,4)
171  alu(16*ii- 3)= alutmp(4,1)
172  alu(16*ii- 2)= alutmp(4,2)
173  alu(16*ii- 1)= alutmp(4,3)
174  alu(16*ii )= alutmp(4,4)
175 
176  enddo
177 
178  isfirst = .true.
179 
180  initialized = .true.
181  hecmat%Iarray(98) = 0 ! symbolic setup done
182  hecmat%Iarray(97) = 0 ! numerical setup done
183 
184  !write(*,*) 'DEBUG: SSOR setup done', hecmw_Wtime()-t0
185 
186  end subroutine hecmw_precond_ssor_44_setup
187 
188  subroutine hecmw_precond_ssor_44_apply(ZP)
190  implicit none
191  real(kind=kreal), intent(inout) :: zp(:)
192  integer(kind=kint) :: ic, i, iold, j, isl, iel, isu, ieu, k
193  real(kind=kreal) :: sw1, sw2, sw3, sw4, x1, x2, x3, x4
194 
195  ! added for turning >>>
196  integer(kind=kint), parameter :: numofblockperthread = 100
197  integer(kind=kint), save :: numofthread = 1, numofblock
198  integer(kind=kint), save, allocatable :: ictoblockindex(:)
199  integer(kind=kint), save, allocatable :: blockindextocolorindex(:)
200  integer(kind=kint), save :: sectorcachesize0, sectorcachesize1
201  integer(kind=kint) :: blockindex, elementcount, numofelement, ii
202  real(kind=kreal) :: numofelementperblock
203  integer(kind=kint) :: my_rank
204 
205  if (isfirst) then
206  !$ numOfThread = omp_get_max_threads()
207  numofblock = numofthread * numofblockperthread
208  if (allocated(ictoblockindex)) deallocate(ictoblockindex)
209  if (allocated(blockindextocolorindex)) deallocate(blockindextocolorindex)
210  allocate (ictoblockindex(0:ncolor), &
211  blockindextocolorindex(0:numofblock + ncolor))
212  numofelement = n + indexl(n) + indexu(n)
213  numofelementperblock = dble(numofelement) / numofblock
214  blockindex = 0
215  ictoblockindex = -1
216  ictoblockindex(0) = 0
217  blockindextocolorindex = -1
218  blockindextocolorindex(0) = 0
219  my_rank = hecmw_comm_get_rank()
220  ! write(9000+my_rank,*) &
221  ! '# numOfElementPerBlock =', numOfElementPerBlock
222  ! write(9000+my_rank,*) &
223  ! '# ic, blockIndex, colorIndex, elementCount'
224  do ic = 1, ncolor
225  elementcount = 0
226  ii = 1
227  do i = colorindex(ic-1)+1, colorindex(ic)
228  elementcount = elementcount + 1
229  elementcount = elementcount + (indexl(i) - indexl(i-1))
230  elementcount = elementcount + (indexu(i) - indexu(i-1))
231  if (elementcount > ii * numofelementperblock &
232  .or. i == colorindex(ic)) then
233  ii = ii + 1
234  blockindex = blockindex + 1
235  blockindextocolorindex(blockindex) = i
236  ! write(9000+my_rank,*) ic, blockIndex, &
237  ! blockIndexToColorIndex(blockIndex), elementCount
238  endif
239  enddo
240  ictoblockindex(ic) = blockindex
241  enddo
242  numofblock = blockindex
243 
245  sectorcachesize0, sectorcachesize1 )
246 
247  isfirst = .false.
248  endif
249  ! <<< added for turning
250 
251  !call start_collection("loopInPrecond44")
252 
253  !OCL CACHE_SECTOR_SIZE(sectorCacheSize0,sectorCacheSize1)
254  !OCL CACHE_SUBSECTOR_ASSIGN(ZP)
255 
256  !$omp parallel default(none) &
257  !$omp&shared(NColor,indexL,itemL,indexU,itemU,AL,AU,D,ALU,perm,&
258  !$omp& ZP,icToBlockIndex,blockIndexToColorIndex) &
259  !$omp&private(SW1,SW2,SW3,SW4,X1,X2,X3,X4,ic,i,iold,isL,ieL,isU,ieU,j,k,blockIndex)
260 
261  !C-- FORWARD
262  do ic=1,ncolor
263  !$omp do schedule (static, 1)
264  do blockindex = ictoblockindex(ic-1)+1, ictoblockindex(ic)
265  do i = blockindextocolorindex(blockindex-1)+1, &
266  blockindextocolorindex(blockindex)
267  ! do i = startPos(threadNum, ic), endPos(threadNum, ic)
268  iold = perm(i)
269  sw1= zp(4*iold-3)
270  sw2= zp(4*iold-2)
271  sw3= zp(4*iold-1)
272  sw4= zp(4*iold )
273  isl= indexl(i-1)+1
274  iel= indexl(i)
275  do j= isl, iel
276  !k= perm(itemL(j))
277  k= iteml(j)
278  x1= zp(4*k-3)
279  x2= zp(4*k-2)
280  x3= zp(4*k-1)
281  x4= zp(4*k )
282  sw1= sw1 - al(16*j-15)*x1 - al(16*j-14)*x2 - al(16*j-13)*x3 - al(16*j-12)*x4
283  sw2= sw2 - al(16*j-11)*x1 - al(16*j-10)*x2 - al(16*j- 9)*x3 - al(16*j- 8)*x4
284  sw3= sw3 - al(16*j- 7)*x1 - al(16*j- 6)*x2 - al(16*j- 5)*x3 - al(16*j- 4)*x4
285  sw4= sw4 - al(16*j- 3)*x1 - al(16*j- 2)*x2 - al(16*j- 1)*x3 - al(16*j- 0)*x4
286  enddo ! j
287 
288  x1= sw1
289  x2= sw2
290  x3= sw3
291  x4= sw4
292  x2= x2 - alu(16*i-11)*x1
293  x3= x3 - alu(16*i- 7)*x1 - alu(16*i- 6)*x2
294  x4= x4 - alu(16*i- 3)*x1 - alu(16*i- 2)*x2 - alu(16*i- 1)*x3
295 
296  x4= alu(16*i )* x4
297  x3= alu(16*i- 5)*( x3 - alu(16*i- 4)*x4 )
298  x2= alu(16*i-10)*( x2 - alu(16*i- 8)*x4 - alu(16*i- 9)*x3 )
299  x1= alu(16*i-15)*( x1 - alu(16*i-12)*x4 - alu(16*i-13)*x3 - alu(16*i-14)*x2 )
300 
301  zp(4*iold-3)= x1
302  zp(4*iold-2)= x2
303  zp(4*iold-1)= x3
304  zp(4*iold )= x4
305  enddo ! i
306  enddo ! blockIndex
307  !$omp end do
308  enddo ! ic
309 
310  !C-- BACKWARD
311  do ic=ncolor, 1, -1
312  !$omp do schedule (static, 1)
313  do blockindex = ictoblockindex(ic), ictoblockindex(ic-1)+1, -1
314  do i = blockindextocolorindex(blockindex), &
315  blockindextocolorindex(blockindex-1)+1, -1
316  ! do blockIndex = icToBlockIndex(ic-1)+1, icToBlockIndex(ic)
317  ! do i = blockIndexToColorIndex(blockIndex-1)+1, &
318  ! blockIndexToColorIndex(blockIndex)
319  ! do i = endPos(threadNum, ic), startPos(threadNum, ic), -1
320  sw1= 0.d0
321  sw2= 0.d0
322  sw3= 0.d0
323  sw4= 0.d0
324  isu= indexu(i-1) + 1
325  ieu= indexu(i)
326  do j= ieu, isu, -1
327  !k= perm(itemU(j))
328  k= itemu(j)
329  x1= zp(4*k-3)
330  x2= zp(4*k-2)
331  x3= zp(4*k-1)
332  x4= zp(4*k )
333  sw1= sw1 + au(16*j-15)*x1 + au(16*j-14)*x2 + au(16*j-13)*x3 + au(16*j-12)*x4
334  sw2= sw2 + au(16*j-11)*x1 + au(16*j-10)*x2 + au(16*j- 9)*x3 + au(16*j- 8)*x4
335  sw3= sw3 + au(16*j- 7)*x1 + au(16*j- 6)*x2 + au(16*j- 5)*x3 + au(16*j- 4)*x4
336  sw4= sw4 + au(16*j- 3)*x1 + au(16*j- 2)*x2 + au(16*j- 1)*x3 + au(16*j- 0)*x4
337 
338  enddo ! j
339 
340  x1= sw1
341  x2= sw2
342  x3= sw3
343  x4= sw4
344  x2= x2 - alu(16*i-11)*x1
345  x3= x3 - alu(16*i- 7)*x1 - alu(16*i- 6)*x2
346  x4= x4 - alu(16*i- 3)*x1 - alu(16*i- 2)*x2 - alu(16*i- 1)*x3
347  x4= alu(16*i )* x4
348  x3= alu(16*i- 5)*( x3 - alu(16*i- 4)*x4 )
349  x2= alu(16*i-10)*( x2 - alu(16*i- 8)*x4 - alu(16*i- 9)*x3 )
350  x1= alu(16*i-15)*( x1 - alu(16*i-12)*x4 - alu(16*i-13)*x3 - alu(16*i-14)*x2 )
351  iold = perm(i)
352  zp(4*iold-3)= zp(4*iold-3) - x1
353  zp(4*iold-2)= zp(4*iold-2) - x2
354  zp(4*iold-1)= zp(4*iold-1) - x3
355  zp(4*iold )= zp(4*iold ) - x4
356  enddo ! i
357  enddo ! blockIndex
358  !$omp end do
359  enddo ! ic
360  !$omp end parallel
361 
362  !OCL END_CACHE_SUBSECTOR
363  !OCL END_CACHE_SECTOR_SIZE
364 
365  !call stop_collection("loopInPrecond44")
366 
367  end subroutine hecmw_precond_ssor_44_apply
368 
369  subroutine hecmw_precond_ssor_44_clear(hecMAT)
370  implicit none
371  type(hecmwst_matrix), intent(inout) :: hecmat
372  integer(kind=kint ) :: nthreads = 1
373  !$ nthreads = omp_get_max_threads()
374  if (associated(colorindex)) deallocate(colorindex)
375  if (associated(perm)) deallocate(perm)
376  if (associated(iperm)) deallocate(iperm)
377  if (associated(alu)) deallocate(alu)
378  if (nthreads >= 1) then
379  if (associated(d)) deallocate(d)
380  if (associated(al)) deallocate(al)
381  if (associated(au)) deallocate(au)
382  if (associated(indexl)) deallocate(indexl)
383  if (associated(indexu)) deallocate(indexu)
384  if (associated(iteml)) deallocate(iteml)
385  if (associated(itemu)) deallocate(itemu)
386  end if
387  nullify(colorindex)
388  nullify(perm)
389  nullify(iperm)
390  nullify(alu)
391  nullify(d)
392  nullify(al)
393  nullify(au)
394  nullify(indexl)
395  nullify(indexu)
396  nullify(iteml)
397  nullify(itemu)
398  end subroutine hecmw_precond_ssor_44_clear
399 
400  subroutine write_debug_info
401  implicit none
402  integer(kind=kint) :: my_rank, ic, in
403  my_rank = hecmw_comm_get_rank()
404  !--------------------> debug: shizawa
405  if (my_rank.eq.0) then
406  write(*,*) 'DEBUG: Output fort.19000+myrank and fort.29000+myrank for coloring information'
407  endif
408  write(19000+my_rank,'(a)') '#NCOLORTot'
409  write(19000+my_rank,*) ncolor
410  write(19000+my_rank,'(a)') '#ic COLORindex(ic-1)+1 COLORindex(ic)'
411  do ic=1,ncolor
412  write(19000+my_rank,*) ic, colorindex(ic-1)+1,colorindex(ic)
413  enddo ! ic
414  write(29000+my_rank,'(a)') '#n_node'
415  write(29000+my_rank,*) n
416  write(29000+my_rank,'(a)') '#in OLDtoNEW(in) NEWtoOLD(in)'
417  do in=1,n
418  write(29000+my_rank,*) in, iperm(in), perm(in)
419  if (perm(iperm(in)) .ne. in) then
420  write(29000+my_rank,*) '** WARNING **: NEWtoOLD and OLDtoNEW: ',in
421  endif
422  enddo
423  end subroutine write_debug_info
424 
425  subroutine check_ordering
426  implicit none
427  integer(kind=kint) :: ic, i, j, k
428  integer(kind=kint), allocatable :: iicolor(:)
429  ! check color dependence of neighbouring nodes
430  if (ncolor.gt.1) then
431  allocate(iicolor(n))
432  do ic=1,ncolor
433  do i= colorindex(ic-1)+1, colorindex(ic)
434  iicolor(i) = ic
435  enddo ! i
436  enddo ! ic
437  ! FORWARD: L-part
438  do ic=1,ncolor
439  do i= colorindex(ic-1)+1, colorindex(ic)
440  do j= indexl(i-1)+1, indexl(i)
441  k= iteml(j)
442  if (iicolor(i).eq.iicolor(k)) then
443  write(*,*) .eq.'** ERROR **: L-part: iicolor(i)iicolor(k)',i,k,iicolor(i)
444  endif
445  enddo ! j
446  enddo ! i
447  enddo ! ic
448  ! BACKWARD: U-part
449  do ic=ncolor, 1, -1
450  do i= colorindex(ic), colorindex(ic-1)+1, -1
451  do j= indexu(i-1)+1, indexu(i)
452  k= itemu(j)
453  if (iicolor(i).eq.iicolor(k)) then
454  write(*,*) .eq.'** ERROR **: U-part: iicolor(i)iicolor(k)',i,k,iicolor(i)
455  endif
456  enddo ! j
457  enddo ! i
458  enddo ! ic
459  deallocate(iicolor)
460  endif ! if (NColor.gt.1)
461  !--------------------< debug: shizawa
462  end subroutine check_ordering
463 
464 end module hecmw_precond_ssor_44
hecmw_precond_ssor_44
Definition: hecmw_precond_SSOR_44.f90:11
m_hecmw_matrix_ordering_cm
Definition: hecmw_matrix_ordering_CM.f90:6
hecmw_matrix_reorder::hecmw_matrix_reorder_values
subroutine, public hecmw_matrix_reorder_values(N, NDOF, perm, iperm, indexL, indexU, itemL, itemU, AL, AU, D, indexLp, indexUp, itemLp, itemUp, ALp, AUp, Dp)
Definition: hecmw_matrix_reorder.f90:71
hecmw_precond_ssor_44::hecmw_precond_ssor_44_apply
subroutine, public hecmw_precond_ssor_44_apply(ZP)
Definition: hecmw_precond_SSOR_44.f90:189
hecmw_matrix_reorder::hecmw_matrix_reorder_profile
subroutine, public hecmw_matrix_reorder_profile(N, perm, iperm, indexL, indexU, itemL, itemU, indexLp, indexUp, itemLp, itemUp)
Definition: hecmw_matrix_reorder.f90:22
hecmw_matrix_reorder::hecmw_matrix_reorder_renum_item
subroutine, public hecmw_matrix_reorder_renum_item(N, perm, indexXp, itemXp)
Definition: hecmw_matrix_reorder.f90:144
hecmw_matrix_misc::hecmw_mat_get_ncolor_in
integer(kind=kint) function, public hecmw_mat_get_ncolor_in(hecMAT)
Definition: hecmw_matrix_misc.f90:531
hecmw_tuning_fx
Definition: hecmw_tuning_fx.f90:6
hecmw_tuning_fx::hecmw_tuning_fx_calc_sector_cache
subroutine, public hecmw_tuning_fx_calc_sector_cache(N, NDOF, sectorCacheSize0, sectorCacheSize1)
Definition: hecmw_tuning_fx.f90:25
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
m_hecmw_matrix_ordering_mc::hecmw_matrix_ordering_mc
subroutine, public hecmw_matrix_ordering_mc(N, indexL, itemL, indexU, itemU, perm_cur, ncolor_in, ncolor_out, COLORindex, perm, iperm)
Definition: hecmw_matrix_ordering_MC.f90:17
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
hecmw_matrix_misc::hecmw_mat_get_sigma_diag
real(kind=kreal) function, public hecmw_mat_get_sigma_diag(hecMAT)
Definition: hecmw_matrix_misc.f90:694
m_hecmw_matrix_ordering_cm::hecmw_matrix_ordering_rcm
subroutine, public hecmw_matrix_ordering_rcm(N, indexL, itemL, indexU, itemU, perm, iperm)
Definition: hecmw_matrix_ordering_CM.f90:67
hecmw_precond_ssor_44::hecmw_precond_ssor_44_clear
subroutine, public hecmw_precond_ssor_44_clear(hecMAT)
Definition: hecmw_precond_SSOR_44.f90:370
hecmw_precond_ssor_44::hecmw_precond_ssor_44_setup
subroutine, public hecmw_precond_ssor_44_setup(hecMAT)
Definition: hecmw_precond_SSOR_44.f90:47
m_hecmw_matrix_ordering_mc
Definition: hecmw_matrix_ordering_MC.f90:6
hecmw_matrix_misc
Definition: hecmw_matrix_misc.f90:6
hecmw_matrix_reorder
Definition: hecmw_matrix_reorder.f90:6
hecmw_util::hecmw_comm_get_rank
integer(kind=kint) function hecmw_comm_get_rank()
Definition: hecmw_util_f.F90:582
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444