FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
hecmw_solver_las_33.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
8  implicit none
9 
10  private
11 
12  public :: hecmw_matvec_33
15  public :: hecmw_matresid_33
16  public :: hecmw_rel_resid_l2_33
17  public :: hecmw_tvec_33
18  public :: hecmw_ttvec_33
19  public :: hecmw_ttmattvec_33
20  public :: hecmw_mat_diag_sr_33
21 
22  ! ! for communication hiding in matvec
23  ! integer(kind=kint), save, allocatable :: index_o(:), item_o(:)
24  ! real(kind=kreal), save, allocatable :: A_o(:)
25  logical, save :: async_matvec_flg = .false.
26 
27 contains
28 
29  !C
30  !C***
31  !C*** hecmw_matvec_33
32  !C***
33  !C
34  subroutine hecmw_matvec_33 (hecMESH, hecMAT, X, Y, time_Ax, COMMtime)
37  implicit none
38  type (hecmwst_local_mesh), intent(in) :: hecmesh
39  type (hecmwst_matrix), intent(in), target :: hecmat
40  real(kind=kreal), intent(in) :: x(:)
41  real(kind=kreal), intent(out) :: y(:)
42  real(kind=kreal), intent(inout) :: time_ax
43  real(kind=kreal), intent(inout), optional :: commtime
44 
45  real(kind=kreal) :: tcomm
46  real(kind=kreal), allocatable :: wk(:)
47 
48  tcomm = 0.d0
49 
50  if (hecmw_mat_get_flag_mpcmatvec(hecmat) /= 0) then
51  allocate(wk(hecmat%NP * hecmat%NDOF))
52  call hecmw_ttmattvec_33(hecmesh, hecmat, x, y, wk, time_ax, tcomm)
53  deallocate(wk)
54  else
55  call hecmw_matvec_33_inner(hecmesh, hecmat, x, y, time_ax, tcomm)
56  endif
57 
58  if (present(commtime)) commtime = commtime + tcomm
59  end subroutine hecmw_matvec_33
60 
61  !C
62  !C***
63  !C*** hecmw_matvec_33_set_async
64  !C***
65  !C
66  subroutine hecmw_matvec_33_set_async (hecMAT)
68  implicit none
69  type (hecmwst_matrix), intent(in) :: hecmat
70  ! integer(kind=kint) :: i, j, jS, jE, idx, in
71 
72  ! allocate(index_o(0:hecMAT%N))
73  ! index_o(0) = 0
74  ! do i = 1, hecMAT%N
75  ! jS= hecMAT%indexU(i-1) + 1
76  ! jE= hecMAT%indexU(i )
77  ! idx = index_o(i-1)
78  ! do j= jS, jE
79  ! in = hecMAT%itemU(j)
80  ! if (in <= hecMAT%N) cycle
81  ! idx = idx + 1
82  ! enddo
83  ! index_o(i) = idx
84  ! enddo
85  ! allocate(item_o(idx))
86  ! allocate(A_o(idx*9))
87  ! do i = 1, hecMAT%N
88  ! jS= hecMAT%indexU(i-1) + 1
89  ! jE= hecMAT%indexU(i )
90  ! idx = index_o(i-1)
91  ! do j= jS, jE
92  ! in = hecMAT%itemU(j)
93  ! if (in <= hecMAT%N) cycle
94  ! idx = idx + 1
95  ! item_o(idx) = hecMAT%itemU(j) - hecMAT%N
96  ! A_o(9*idx-8:9*idx) = hecMAT%AU(9*j-8:9*j)
97  ! enddo
98  ! enddo
99  ! async_matvec_flg = .true.
100  end subroutine hecmw_matvec_33_set_async
101 
102  !C
103  !C***
104  !C*** hecmw_matvec_33_unset_async
105  !C***
106  !C
107  subroutine hecmw_matvec_33_unset_async
108  implicit none
109  ! if (allocated(index_o)) deallocate(index_o)
110  ! if (allocated(item_o)) deallocate(item_o)
111  ! if (allocated(A_o)) deallocate(A_o)
112  ! async_matvec_flg = .false.
113  end subroutine hecmw_matvec_33_unset_async
114 
115  !C
116  !C***
117  !C*** hecmw_matvec_33_inner ( private subroutine )
118  !C***
119  !C
120  subroutine hecmw_matvec_33_inner (hecMESH, hecMAT, X, Y, time_Ax, COMMtime)
121  use hecmw_util
122  use m_hecmw_comm_f
124  use hecmw_jad_type
125  use hecmw_tuning_fx
126  !$ use omp_lib
127 
128  implicit none
129  type (hecmwst_local_mesh), intent(in) :: hecmesh
130  type (hecmwst_matrix), intent(in), target :: hecmat
131  real(kind=kreal), intent(in) :: x(:)
132  real(kind=kreal), intent(out) :: y(:)
133  real(kind=kreal), intent(inout) :: time_ax
134  real(kind=kreal), intent(inout), optional :: commtime
135 
136  real(kind=kreal) :: start_time, end_time, tcomm
137  integer(kind=kint) :: i, j, js, je, in
138  real(kind=kreal) :: yv1, yv2, yv3, x1, x2, x3
139 
140  integer(kind=kint) :: n, np
141  integer(kind=kint), pointer :: indexl(:), iteml(:), indexu(:), itemu(:)
142  real(kind=kreal), pointer :: al(:), au(:), d(:)
143 
144  ! added for turning >>>
145  integer, parameter :: numofblockperthread = 100
146  logical, save :: isfirst = .true.
147  integer, save :: numofthread = 1
148  integer, save, allocatable :: startpos(:), endpos(:)
149  integer(kind=kint), save :: sectorcachesize0, sectorcachesize1
150  integer(kind=kint) :: threadnum, blocknum, numofblock
151  integer(kind=kint) :: numofelement, elementcount, blockindex
152  real(kind=kreal) :: numofelementperblock
153  ! <<< added for turning
154 
155  if (hecmw_jad_is_initialized().ne.0) then
156  tcomm = 0.d0
157  start_time = hecmw_wtime()
158  call hecmw_jad_matvec(hecmesh, hecmat, x, y, tcomm)
159  end_time = hecmw_wtime()
160  time_ax = time_ax + end_time - start_time - tcomm
161  if (present(commtime)) commtime = commtime + tcomm
162  else
163 
164  n = hecmat%N
165  np = hecmat%NP
166  indexl => hecmat%indexL
167  indexu => hecmat%indexU
168  iteml => hecmat%itemL
169  itemu => hecmat%itemU
170  al => hecmat%AL
171  au => hecmat%AU
172  d => hecmat%D
173 
174  ! added for turning >>>
175  if (.not. isfirst) then
176  numofblock = numofthread * numofblockperthread
177  if (endpos(numofblock-1) .ne. n-1) then
178  deallocate(startpos, endpos)
179  isfirst = .true.
180  endif
181  endif
182  if (isfirst) then
183  !$ numOfThread = omp_get_max_threads()
184  numofblock = numofthread * numofblockperthread
185  allocate (startpos(0 : numofblock - 1), endpos(0 : numofblock - 1))
186  numofelement = n + indexl(n) + indexu(n)
187  numofelementperblock = dble(numofelement) / numofblock
188  blocknum = 0
189  elementcount = 0
190  startpos(blocknum) = 1
191  do i= 1, n
192  elementcount = elementcount + 1
193  elementcount = elementcount + (indexl(i) - indexl(i-1))
194  elementcount = elementcount + (indexu(i) - indexu(i-1))
195  if (elementcount > (blocknum + 1) * numofelementperblock) then
196  endpos(blocknum) = i
197  ! write(9000+hecMESH%my_rank,*) mod(blockNum, numOfThread), &
198  ! startPos(blockNum), endPos(blockNum)
199  blocknum = blocknum + 1
200  startpos(blocknum) = i + 1
201  if (blocknum == (numofblock - 1)) exit
202  endif
203  enddo
204  endpos(blocknum) = n
205  ! write(9000+hecMESH%my_rank,*) mod(blockNum, numOfThread), &
206  ! startPos(blockNum), endPos(blockNum)
207  ! for irregular data
208  do i= blocknum+1, numofblock-1
209  startpos(i) = n
210  endpos(i) = n-1
211  ! write(9000+hecMESH%my_rank,*) mod(i, numOfThread), &
212  ! startPos(i), endPos(i)
213  end do
214 
216  sectorcachesize0, sectorcachesize1)
217 
218  isfirst = .false.
219  endif
220  ! <<< added for turning
221 
222  start_time= hecmw_wtime()
223  ! if (async_matvec_flg) then
224  ! call hecmw_update_3_R_async (hecMESH, X, NP, ireq)
225  ! else
226  call hecmw_update_r (hecmesh, x, np, 3)
227  ! endif
228  end_time= hecmw_wtime()
229  if (present(commtime)) commtime = commtime + end_time - start_time
230 
231  start_time = hecmw_wtime()
232 
233  !call fapp_start("loopInMatvec33", 1, 0)
234  !call start_collection("loopInMatvec33")
235 
236  !OCL CACHE_SECTOR_SIZE(sectorCacheSize0,sectorCacheSize1)
237  !OCL CACHE_SUBSECTOR_ASSIGN(X)
238 
239  !$OMP PARALLEL DEFAULT(NONE) &
240  !$OMP&PRIVATE(i,X1,X2,X3,YV1,YV2,YV3,jS,jE,j,in,threadNum,blockNum,blockIndex) &
241  !$OMP&SHARED(D,AL,AU,indexL,itemL,indexU,itemU,X,Y,startPos,endPos,numOfThread,N,async_matvec_flg)
242  threadnum = 0
243  !$ threadNum = omp_get_thread_num()
244  do blocknum = 0 , numofblockperthread - 1
245  blockindex = blocknum * numofthread + threadnum
246  do i = startpos(blockindex), endpos(blockindex)
247  x1= x(3*i-2)
248  x2= x(3*i-1)
249  x3= x(3*i )
250  yv1= d(9*i-8)*x1 + d(9*i-7)*x2 + d(9*i-6)*x3
251  yv2= d(9*i-5)*x1 + d(9*i-4)*x2 + d(9*i-3)*x3
252  yv3= d(9*i-2)*x1 + d(9*i-1)*x2 + d(9*i )*x3
253 
254  js= indexl(i-1) + 1
255  je= indexl(i )
256  do j= js, je
257  in = iteml(j)
258  x1= x(3*in-2)
259  x2= x(3*in-1)
260  x3= x(3*in )
261  yv1= yv1 + al(9*j-8)*x1 + al(9*j-7)*x2 + al(9*j-6)*x3
262  yv2= yv2 + al(9*j-5)*x1 + al(9*j-4)*x2 + al(9*j-3)*x3
263  yv3= yv3 + al(9*j-2)*x1 + al(9*j-1)*x2 + al(9*j )*x3
264  enddo
265  js= indexu(i-1) + 1
266  je= indexu(i )
267  do j= js, je
268  in = itemu(j)
269  ! if (async_matvec_flg .and. in > N) cycle
270  x1= x(3*in-2)
271  x2= x(3*in-1)
272  x3= x(3*in )
273  yv1= yv1 + au(9*j-8)*x1 + au(9*j-7)*x2 + au(9*j-6)*x3
274  yv2= yv2 + au(9*j-5)*x1 + au(9*j-4)*x2 + au(9*j-3)*x3
275  yv3= yv3 + au(9*j-2)*x1 + au(9*j-1)*x2 + au(9*j )*x3
276  enddo
277  y(3*i-2)= yv1
278  y(3*i-1)= yv2
279  y(3*i )= yv3
280  enddo
281  enddo
282  !$OMP END PARALLEL
283 
284  !OCL END_CACHE_SUBSECTOR
285  !OCL END_CACHE_SECTOR_SIZE
286 
287  !call stop_collection("loopInMatvec33")
288  !call fapp_stop("loopInMatvec33", 1, 0)
289 
290  end_time = hecmw_wtime()
291  time_ax = time_ax + end_time - start_time
292 
293  ! if (async_matvec_flg) then
294  ! START_TIME= HECMW_WTIME()
295  ! call hecmw_update_3_R_wait (hecMESH, ireq)
296  ! END_TIME= HECMW_WTIME()
297  ! if (present(COMMtime)) COMMtime = COMMtime + END_TIME - START_TIME
298 
299  ! START_TIME = hecmw_Wtime()
300 
301  ! do i = 1, N
302  ! jS= index_o(i-1) + 1
303  ! jE= index_o(i )
304  ! if (jS > jE) cycle
305  ! YV1= 0.d0
306  ! YV2= 0.d0
307  ! YV3= 0.d0
308  ! do j=jS, jE
309  ! in = item_o(j)
310  ! X1= X(3*(N+in)-2)
311  ! X2= X(3*(N+in)-1)
312  ! X3= X(3*(N+in) )
313  ! YV1= YV1 + A_o(9*j-8)*X1 + A_o(9*j-7)*X2 + A_o(9*j-6)*X3
314  ! YV2= YV2 + A_o(9*j-5)*X1 + A_o(9*j-4)*X2 + A_o(9*j-3)*X3
315  ! YV3= YV3 + A_o(9*j-2)*X1 + A_o(9*j-1)*X2 + A_o(9*j )*X3
316  ! enddo
317  ! Y(3*i-2)= Y(3*i-2)+YV1
318  ! Y(3*i-1)= Y(3*i-1)+YV2
319  ! Y(3*i )= Y(3*i )+YV3
320  ! enddo
321 
322  ! END_TIME = hecmw_Wtime()
323  ! time_Ax = time_Ax + END_TIME - START_TIME
324  ! endif
325 
326  endif
327  end subroutine hecmw_matvec_33_inner
328 
329  !C
330  !C***
331  !C*** hecmw_matresid_33
332  !C***
333  !C
334  subroutine hecmw_matresid_33 (hecMESH, hecMAT, X, B, R, time_Ax, COMMtime)
336  implicit none
337  type (hecmwst_local_mesh), intent(in) :: hecmesh
338  type (hecmwst_matrix), intent(in) :: hecmat
339  real(kind=kreal), intent(in) :: x(:), b(:)
340  real(kind=kreal), intent(out) :: r(:)
341  real(kind=kreal), intent(inout) :: time_ax
342  real(kind=kreal), intent(inout), optional :: commtime
343 
344  integer(kind=kint) :: i
345  real(kind=kreal) :: tcomm
346 
347  tcomm = 0.d0
348  call hecmw_matvec_33 (hecmesh, hecmat, x, r, time_ax, tcomm)
349  if (present(commtime)) commtime = commtime + tcomm
350  !$omp parallel default(none),private(i),shared(hecMAT,R,B)
351  !$omp do
352  do i = 1, hecmat%N * 3
353  r(i) = b(i) - r(i)
354  enddo
355  !$omp end do
356  !$omp end parallel
357  end subroutine hecmw_matresid_33
358 
359  !C
360  !C***
361  !C*** hecmw_rel_resid_L2_33
362  !C***
363  !C
364  function hecmw_rel_resid_l2_33 (hecMESH, hecMAT, time_Ax, COMMtime)
367  implicit none
368  real(kind=kreal) :: hecmw_rel_resid_l2_33
369  type ( hecmwst_local_mesh ), intent(in) :: hecmesh
370  type ( hecmwst_matrix ), intent(in) :: hecmat
371  real(kind=kreal), intent(inout) :: time_ax
372  real(kind=kreal), intent(inout), optional :: commtime
373 
374  real(kind=kreal), allocatable :: r(:)
375  real(kind=kreal) :: bnorm2, rnorm2
376  real(kind=kreal) :: tcomm
377 
378  allocate(r(hecmat%NDOF*hecmat%NP))
379 
380  tcomm = 0.d0
381  call hecmw_innerproduct_r(hecmesh, hecmat%NDOF, &
382  hecmat%B, hecmat%B, bnorm2, tcomm)
383  if (bnorm2 == 0.d0) then
384  bnorm2 = 1.d0
385  endif
386  call hecmw_matresid_33(hecmesh, hecmat, hecmat%X, hecmat%B, r, time_ax, tcomm)
387  call hecmw_innerproduct_r(hecmesh, hecmat%NDOF, r, r, rnorm2, tcomm)
388  hecmw_rel_resid_l2_33 = sqrt(rnorm2 / bnorm2)
389 
390  if (present(commtime)) commtime = commtime + tcomm
391 
392  deallocate(r)
393  end function hecmw_rel_resid_l2_33
394 
395  !C
396  !C***
397  !C*** hecmw_Tvec_33
398  !C***
399  !C
400  subroutine hecmw_tvec_33 (hecMESH, X, Y, COMMtime)
402  use m_hecmw_comm_f
403  implicit none
404  type (hecmwst_local_mesh), intent(in) :: hecmesh
405  real(kind=kreal), intent(in) :: x(:)
406  real(kind=kreal), intent(out) :: y(:)
407  real(kind=kreal), intent(inout) :: commtime
408 
409  real(kind=kreal) :: start_time, end_time
410  integer(kind=kint) :: i, j, jj, k, kk
411 
412  start_time= hecmw_wtime()
413  call hecmw_update_r (hecmesh, x, hecmesh%n_node, 3)
414  end_time= hecmw_wtime()
415  commtime = commtime + end_time - start_time
416 
417  !$omp parallel default(none),private(i,k,kk,j,jj),shared(hecMESH,X,Y)
418  !$omp do
419  do i= 1, hecmesh%nn_internal * hecmesh%n_dof
420  y(i)= x(i)
421  enddo
422  !$omp end do
423 
424  !$omp do
425  outer: do i= 1, hecmesh%mpc%n_mpc
426  do j= hecmesh%mpc%mpc_index(i-1) + 1, hecmesh%mpc%mpc_index(i)
427  if (hecmesh%mpc%mpc_dof(j) > 3) cycle outer
428  enddo
429  k = hecmesh%mpc%mpc_index(i-1) + 1
430  kk = 3 * (hecmesh%mpc%mpc_item(k) - 1) + hecmesh%mpc%mpc_dof(k)
431  y(kk) = 0.d0
432  do j= hecmesh%mpc%mpc_index(i-1) + 2, hecmesh%mpc%mpc_index(i)
433  jj = 3 * (hecmesh%mpc%mpc_item(j) - 1) + hecmesh%mpc%mpc_dof(j)
434  y(kk) = y(kk) - hecmesh%mpc%mpc_val(j) * x(jj)
435  enddo
436  enddo outer
437  !$omp end do
438  !$omp end parallel
439 
440  end subroutine hecmw_tvec_33
441 
442  !C
443  !C***
444  !C*** hecmw_Ttvec_33
445  !C***
446  !C
447  subroutine hecmw_ttvec_33 (hecMESH, X, Y, COMMtime)
449  use m_hecmw_comm_f
450  implicit none
451  type (hecmwst_local_mesh), intent(in) :: hecmesh
452  real(kind=kreal), intent(in) :: x(:)
453  real(kind=kreal), intent(out) :: y(:)
454  real(kind=kreal), intent(inout) :: commtime
455 
456  real(kind=kreal) :: start_time, end_time
457  integer(kind=kint) :: i, j, jj, k, kk
458 
459  start_time= hecmw_wtime()
460  call hecmw_update_r (hecmesh, x, hecmesh%n_node, 3)
461  end_time= hecmw_wtime()
462  commtime = commtime + end_time - start_time
463 
464  !$omp parallel default(none),private(i,k,kk,j,jj),shared(hecMESH,X,Y)
465  !$omp do
466  do i= 1, hecmesh%nn_internal * hecmesh%n_dof
467  y(i)= x(i)
468  enddo
469  !$omp end do
470 
471  !$omp do
472  outer: do i= 1, hecmesh%mpc%n_mpc
473  do j= hecmesh%mpc%mpc_index(i-1) + 1, hecmesh%mpc%mpc_index(i)
474  if (hecmesh%mpc%mpc_dof(j) > 3) cycle outer
475  enddo
476  k = hecmesh%mpc%mpc_index(i-1) + 1
477  kk = 3 * (hecmesh%mpc%mpc_item(k) - 1) + hecmesh%mpc%mpc_dof(k)
478  y(kk) = 0.d0
479  do j= hecmesh%mpc%mpc_index(i-1) + 2, hecmesh%mpc%mpc_index(i)
480  jj = 3 * (hecmesh%mpc%mpc_item(j) - 1) + hecmesh%mpc%mpc_dof(j)
481  !$omp atomic
482  y(jj) = y(jj) - hecmesh%mpc%mpc_val(j) * x(kk)
483  enddo
484  enddo outer
485  !$omp end do
486  !$omp end parallel
487 
488  end subroutine hecmw_ttvec_33
489 
490  !C
491  !C***
492  !C*** hecmw_TtmatTvec_33
493  !C***
494  !C
495  subroutine hecmw_ttmattvec_33 (hecMESH, hecMAT, X, Y, W, time_Ax, COMMtime)
497  implicit none
498  type (hecmwst_local_mesh), intent(in) :: hecmesh
499  type (hecmwst_matrix), intent(in) :: hecmat
500  real(kind=kreal), intent(in) :: x(:)
501  real(kind=kreal), intent(out) :: y(:), w(:)
502  real(kind=kreal), intent(inout) :: time_ax
503  real(kind=kreal), intent(inout) :: commtime
504 
505  call hecmw_tvec_33(hecmesh, x, y, commtime)
506  call hecmw_matvec_33_inner(hecmesh, hecmat, y, w, time_ax, commtime)
507  call hecmw_ttvec_33(hecmesh, w, y, commtime)
508 
509  end subroutine hecmw_ttmattvec_33
510 
511 
512  !C
513  !C***
514  !C*** hecmw_mat_diag_sr_33
515  !C***
516  !C
517  subroutine hecmw_mat_diag_sr_33(hecMESH, hecMAT, COMMtime)
519  use m_hecmw_comm_f
520  implicit none
521  type (hecmwst_local_mesh), intent(in) :: hecmesh
522  type (hecmwst_matrix), intent(inout), target :: hecmat
523  real(kind=kreal), intent(inout), optional :: commtime
524  real(kind=kreal), allocatable :: w(:,:)
525  real(kind=kreal), pointer :: d(:)
526  integer(kind=kint) :: ip
527  real(kind=kreal) :: start_time, end_time
528  allocate(w(3*hecmat%NP,3))
529  d => hecmat%D
530  do ip= 1, hecmat%N
531  w(3*ip-2,1)= d(9*ip-8); w(3*ip-2,2)= d(9*ip-7); w(3*ip-2,3)= d(9*ip-6)
532  w(3*ip-1,1)= d(9*ip-5); w(3*ip-1,2)= d(9*ip-4); w(3*ip-1,3)= d(9*ip-3)
533  w(3*ip ,1)= d(9*ip-2); w(3*ip ,2)= d(9*ip-1); w(3*ip ,3)= d(9*ip )
534  enddo
535  start_time= hecmw_wtime()
536  call hecmw_update_r (hecmesh, w(:,1), hecmat%NP, 3)
537  call hecmw_update_r (hecmesh, w(:,2), hecmat%NP, 3)
538  call hecmw_update_r (hecmesh, w(:,3), hecmat%NP, 3)
539  end_time= hecmw_wtime()
540  if (present(commtime)) commtime = commtime + end_time - start_time
541  do ip= hecmat%N+1, hecmat%NP
542  d(9*ip-8)= w(3*ip-2,1); d(9*ip-7)= w(3*ip-2,2); d(9*ip-6)= w(3*ip-2,3)
543  d(9*ip-5)= w(3*ip-1,1); d(9*ip-4)= w(3*ip-1,2); d(9*ip-3)= w(3*ip-1,3)
544  d(9*ip-2)= w(3*ip ,1); d(9*ip-1)= w(3*ip ,2); d(9*ip )= w(3*ip ,3)
545  enddo
546  deallocate(w)
547  end subroutine hecmw_mat_diag_sr_33
548 
549 end module hecmw_solver_las_33
hecmw_solver_las_33::hecmw_matvec_33
subroutine, public hecmw_matvec_33(hecMESH, hecMAT, X, Y, time_Ax, COMMtime)
Definition: hecmw_solver_las_33.f90:35
hecmw_solver_las_33::hecmw_ttvec_33
subroutine, public hecmw_ttvec_33(hecMESH, X, Y, COMMtime)
Definition: hecmw_solver_las_33.f90:448
hecmw_matrix_misc::hecmw_mat_get_flag_mpcmatvec
integer(kind=kint) function, public hecmw_mat_get_flag_mpcmatvec(hecMAT)
Definition: hecmw_matrix_misc.f90:643
hecmw_util::hecmw_wtime
real(kind=kreal) function hecmw_wtime()
Definition: hecmw_util_f.F90:549
hecmw_solver_las_33::hecmw_matresid_33
subroutine, public hecmw_matresid_33(hecMESH, hecMAT, X, B, R, time_Ax, COMMtime)
Definition: hecmw_solver_las_33.f90:335
hecmw_solver_las_33::hecmw_matvec_33_unset_async
subroutine, public hecmw_matvec_33_unset_async
Definition: hecmw_solver_las_33.f90:108
hecmw_jad_type::hecmw_jad_is_initialized
integer(kind=kint) function, public hecmw_jad_is_initialized()
Definition: hecmw_jadm.f90:56
hecmw_solver_las_33::hecmw_matvec_33_set_async
subroutine, public hecmw_matvec_33_set_async(hecMAT)
Definition: hecmw_solver_las_33.f90:67
hecmw_jad_type
Jagged Diagonal Matrix storage for vector processors. Original code was provided by JAMSTEC.
Definition: hecmw_jadm.f90:8
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_solver_misc::hecmw_innerproduct_r
subroutine hecmw_innerproduct_r(hecMESH, ndof, X, Y, sum, COMMtime)
Definition: hecmw_solver_misc.f90:47
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
hecmw_util::hecmwst_local_mesh
Definition: hecmw_util_f.F90:234
m_hecmw_comm_f
Definition: hecmw_comm_f.F90:6
hecmw_solver_las_33::hecmw_tvec_33
subroutine, public hecmw_tvec_33(hecMESH, X, Y, COMMtime)
Definition: hecmw_solver_las_33.f90:401
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
hecmw_solver_las_33::hecmw_rel_resid_l2_33
real(kind=kreal) function, public hecmw_rel_resid_l2_33(hecMESH, hecMAT, time_Ax, COMMtime)
Definition: hecmw_solver_las_33.f90:365
hecmw_solver_las_33::hecmw_ttmattvec_33
subroutine, public hecmw_ttmattvec_33(hecMESH, hecMAT, X, Y, W, time_Ax, COMMtime)
Definition: hecmw_solver_las_33.f90:496
hecmw_matrix_misc
Definition: hecmw_matrix_misc.f90:6
hecmw_solver_misc
Definition: hecmw_solver_misc.f90:6
hecmw_solver_las_33
Definition: hecmw_solver_las_33.f90:6
hecmw_solver_las_33::hecmw_mat_diag_sr_33
subroutine, public hecmw_mat_diag_sr_33(hecMESH, hecMAT, COMMtime)
Definition: hecmw_solver_las_33.f90:518
m_hecmw_comm_f::hecmw_update_r
subroutine hecmw_update_r(hecMESH, val, n, m)
Definition: hecmw_comm_f.F90:683
hecmw_jad_type::hecmw_jad_matvec
subroutine, public hecmw_jad_matvec(hecMESH, hecMAT, X, Y, COMMtime)
Definition: hecmw_jadm.f90:61
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444