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