FrontISTR  5.9.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)
35  use hecmw_util
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)
67  use hecmw_util
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
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(:), indexa(:), itema(:)
142  real(kind=kreal), pointer :: al(:), au(:), d(:), a(:)
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  indexa => hecmat%indexA
169  iteml => hecmat%itemL
170  itemu => hecmat%itemU
171  itema => hecmat%itemA
172  al => hecmat%AL
173  au => hecmat%AU
174  d => hecmat%D
175  a => hecmat%A
176 
177  ! added for turning >>>
178 #ifndef _OPENACC
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 
220  sectorcachesize0, sectorcachesize1)
221 
222  isfirst = .false.
223  endif
224 #endif
225  ! <<< added for turning
226 
227  start_time= hecmw_wtime()
228  ! if (async_matvec_flg) then
229  ! call hecmw_update_3_R_async (hecMESH, X, NP, ireq)
230  ! else
231  call hecmw_update_r (hecmesh, x, np, 3)
232  ! endif
233  end_time= hecmw_wtime()
234  if (present(commtime)) commtime = commtime + end_time - start_time
235 
236  start_time = hecmw_wtime()
237 
238  !call fapp_start("loopInMatvec33", 1, 0)
239  !call start_collection("loopInMatvec33")
240 
241  !OCL CACHE_SECTOR_SIZE(sectorCacheSize0,sectorCacheSize1)
242  !OCL CACHE_SUBSECTOR_ASSIGN(X)
243 
244 #ifdef _OPENACC
245  !$acc kernels
246  !$acc loop independent
247  do i = 1, n
248  x1= x(3*i-2)
249  x2= x(3*i-1)
250  x3= x(3*i )
251  yv1= 0
252  yv2= 0
253  yv3= 0
254  js= indexa(i) + 1
255  je= indexa(i+1)
256  do j = js, je
257  in = itema(j)
258  x1= x(3*in-2)
259  x2= x(3*in-1)
260  x3= x(3*in )
261  yv1= yv1 + a(9*j-8)*x1 + a(9*j-7)*x2 + a(9*j-6)*x3
262  yv2= yv2 + a(9*j-5)*x1 + a(9*j-4)*x2 + a(9*j-3)*x3
263  yv3= yv3 + a(9*j-2)*x1 + a(9*j-1)*x2 + a(9*j )*x3
264  enddo
265  y(3*i-2)= yv1
266  y(3*i-1)= yv2
267  y(3*i )= yv3
268  enddo
269  !$acc end kernels
270 #else
271  !$OMP PARALLEL DEFAULT(NONE) &
272  !$OMP&PRIVATE(i,X1,X2,X3,YV1,YV2,YV3,jS,jE,j,in,threadNum,blockNum,blockIndex) &
273  !$OMP&SHARED(D,AL,AU,indexL,itemL,indexU,itemU,X,Y,startPos,endPos,numOfThread,N,async_matvec_flg)
274  threadnum = 0
275  !$ threadNum = omp_get_thread_num()
276  do blocknum = 0 , numofblockperthread - 1
277  blockindex = blocknum * numofthread + threadnum
278  do i = startpos(blockindex), endpos(blockindex)
279  x1= x(3*i-2)
280  x2= x(3*i-1)
281  x3= x(3*i )
282  yv1= d(9*i-8)*x1 + d(9*i-7)*x2 + d(9*i-6)*x3
283  yv2= d(9*i-5)*x1 + d(9*i-4)*x2 + d(9*i-3)*x3
284  yv3= d(9*i-2)*x1 + d(9*i-1)*x2 + d(9*i )*x3
285 
286  js= indexl(i-1) + 1
287  je= indexl(i )
288  do j= js, je
289  in = iteml(j)
290  x1= x(3*in-2)
291  x2= x(3*in-1)
292  x3= x(3*in )
293  yv1= yv1 + al(9*j-8)*x1 + al(9*j-7)*x2 + al(9*j-6)*x3
294  yv2= yv2 + al(9*j-5)*x1 + al(9*j-4)*x2 + al(9*j-3)*x3
295  yv3= yv3 + al(9*j-2)*x1 + al(9*j-1)*x2 + al(9*j )*x3
296  enddo
297  js= indexu(i-1) + 1
298  je= indexu(i )
299  do j= js, je
300  in = itemu(j)
301  ! if (async_matvec_flg .and. in > N) cycle
302  x1= x(3*in-2)
303  x2= x(3*in-1)
304  x3= x(3*in )
305  yv1= yv1 + au(9*j-8)*x1 + au(9*j-7)*x2 + au(9*j-6)*x3
306  yv2= yv2 + au(9*j-5)*x1 + au(9*j-4)*x2 + au(9*j-3)*x3
307  yv3= yv3 + au(9*j-2)*x1 + au(9*j-1)*x2 + au(9*j )*x3
308  enddo
309  y(3*i-2)= yv1
310  y(3*i-1)= yv2
311  y(3*i )= yv3
312  enddo
313  enddo
314  !$OMP END PARALLEL
315 #endif
316 
317  !OCL END_CACHE_SUBSECTOR
318  !OCL END_CACHE_SECTOR_SIZE
319 
320  !call stop_collection("loopInMatvec33")
321  !call fapp_stop("loopInMatvec33", 1, 0)
322 
323  end_time = hecmw_wtime()
324  time_ax = time_ax + end_time - start_time
325 
326  ! if (async_matvec_flg) then
327  ! START_TIME= HECMW_WTIME()
328  ! call hecmw_update_3_R_wait (hecMESH, ireq)
329  ! END_TIME= HECMW_WTIME()
330  ! if (present(COMMtime)) COMMtime = COMMtime + END_TIME - START_TIME
331 
332  ! START_TIME = hecmw_Wtime()
333 
334  ! do i = 1, N
335  ! jS= index_o(i-1) + 1
336  ! jE= index_o(i )
337  ! if (jS > jE) cycle
338  ! YV1= 0.d0
339  ! YV2= 0.d0
340  ! YV3= 0.d0
341  ! do j=jS, jE
342  ! in = item_o(j)
343  ! X1= X(3*(N+in)-2)
344  ! X2= X(3*(N+in)-1)
345  ! X3= X(3*(N+in) )
346  ! YV1= YV1 + A_o(9*j-8)*X1 + A_o(9*j-7)*X2 + A_o(9*j-6)*X3
347  ! YV2= YV2 + A_o(9*j-5)*X1 + A_o(9*j-4)*X2 + A_o(9*j-3)*X3
348  ! YV3= YV3 + A_o(9*j-2)*X1 + A_o(9*j-1)*X2 + A_o(9*j )*X3
349  ! enddo
350  ! Y(3*i-2)= Y(3*i-2)+YV1
351  ! Y(3*i-1)= Y(3*i-1)+YV2
352  ! Y(3*i )= Y(3*i )+YV3
353  ! enddo
354 
355  ! END_TIME = hecmw_Wtime()
356  ! time_Ax = time_Ax + END_TIME - START_TIME
357  ! endif
358 
359  endif
360  end subroutine hecmw_matvec_33_inner
361 
362  !C
363  !C***
364  !C*** hecmw_matresid_33
365  !C***
366  !C
367  subroutine hecmw_matresid_33 (hecMESH, hecMAT, X, B, R, time_Ax, COMMtime)
368  use hecmw_util
369  implicit none
370  type (hecmwst_local_mesh), intent(in) :: hecmesh
371  type (hecmwst_matrix), intent(in) :: hecmat
372  real(kind=kreal), intent(in) :: x(:), b(:)
373  real(kind=kreal), intent(out) :: r(:)
374  real(kind=kreal), intent(inout) :: time_ax
375  real(kind=kreal), intent(inout), optional :: commtime
376 
377  integer(kind=kint) :: i
378  real(kind=kreal) :: tcomm
379 
380  tcomm = 0.d0
381  call hecmw_matvec_33 (hecmesh, hecmat, x, r, time_ax, tcomm)
382  if (present(commtime)) commtime = commtime + tcomm
383 #ifdef _OPENACC
384  !$acc kernels
385  !$acc loop independent
386 #else
387  !$omp parallel default(none),private(i),shared(hecMAT,R,B)
388  !$omp do
389 #endif
390  do i = 1, hecmat%N * 3
391  r(i) = b(i) - r(i)
392  enddo
393 #ifdef _OPENACC
394  !$acc end kernels
395 #else
396  !$omp end do
397  !$omp end parallel
398 #endif
399  end subroutine hecmw_matresid_33
400 
401  !C
402  !C***
403  !C*** hecmw_rel_resid_L2_33
404  !C***
405  !C
406  function hecmw_rel_resid_l2_33 (hecMESH, hecMAT, time_Ax, COMMtime)
407  use hecmw_util
409  implicit none
410  real(kind=kreal) :: hecmw_rel_resid_l2_33
411  type ( hecmwst_local_mesh ), intent(in) :: hecmesh
412  type ( hecmwst_matrix ), intent(in) :: hecmat
413  real(kind=kreal), intent(inout) :: time_ax
414  real(kind=kreal), intent(inout), optional :: commtime
415 
416  real(kind=kreal), allocatable :: r(:)
417  real(kind=kreal) :: bnorm2, rnorm2
418  real(kind=kreal) :: tcomm
419 
420  allocate(r(hecmat%NDOF*hecmat%NP))
421 
422  tcomm = 0.d0
423  call hecmw_innerproduct_r(hecmesh, hecmat%NDOF, &
424  hecmat%B, hecmat%B, bnorm2, tcomm)
425  if (bnorm2 == 0.d0) then
426  bnorm2 = 1.d0
427  endif
428  call hecmw_matresid_33(hecmesh, hecmat, hecmat%X, hecmat%B, r, time_ax, tcomm)
429  call hecmw_innerproduct_r(hecmesh, hecmat%NDOF, r, r, rnorm2, tcomm)
430  hecmw_rel_resid_l2_33 = sqrt(rnorm2 / bnorm2)
431 
432  if (present(commtime)) commtime = commtime + tcomm
433 
434  deallocate(r)
435  end function hecmw_rel_resid_l2_33
436 
437  !C
438  !C***
439  !C*** hecmw_Tvec_33
440  !C***
441  !C
442  subroutine hecmw_tvec_33 (hecMESH, X, Y, COMMtime)
443  use hecmw_util
444  use m_hecmw_comm_f
445  implicit none
446  type (hecmwst_local_mesh), intent(in) :: hecmesh
447  real(kind=kreal), intent(in) :: x(:)
448  real(kind=kreal), intent(out) :: y(:)
449  real(kind=kreal), intent(inout) :: commtime
450 
451  real(kind=kreal) :: start_time, end_time
452  integer(kind=kint) :: i, j, jj, k, kk
453 
454  start_time= hecmw_wtime()
455  call hecmw_update_r (hecmesh, x, hecmesh%n_node, 3)
456  end_time= hecmw_wtime()
457  commtime = commtime + end_time - start_time
458 
459 #ifdef _OPENACC
460  !$acc kernels
461  !$acc loop independent private(i)
462 #else
463  !$omp parallel default(none),private(i,k,kk,j,jj),shared(hecMESH,X,Y)
464  !$omp do
465 #endif
466  do i= 1, hecmesh%nn_internal * hecmesh%n_dof
467  y(i)= x(i)
468  enddo
469 #ifndef _OPENACC
470  !$omp end do
471 #endif
472 
473 #ifdef _OPENACC
474  !$acc loop independent private(i,k,kk,j,jj)
475 #else
476  !$omp do
477 #endif
478  outer: do i= 1, hecmesh%mpc%n_mpc
479  do j= hecmesh%mpc%mpc_index(i-1) + 1, hecmesh%mpc%mpc_index(i)
480  if (hecmesh%mpc%mpc_dof(j) > 3) cycle outer
481  enddo
482  k = hecmesh%mpc%mpc_index(i-1) + 1
483  kk = 3 * (hecmesh%mpc%mpc_item(k) - 1) + hecmesh%mpc%mpc_dof(k)
484  y(kk) = 0.d0
485  do j= hecmesh%mpc%mpc_index(i-1) + 2, hecmesh%mpc%mpc_index(i)
486  jj = 3 * (hecmesh%mpc%mpc_item(j) - 1) + hecmesh%mpc%mpc_dof(j)
487  y(kk) = y(kk) - hecmesh%mpc%mpc_val(j) * x(jj)
488  enddo
489  enddo outer
490 #ifdef _OPENACC
491  !$acc end kernels
492 #else
493  !$omp end do
494  !$omp end parallel
495 #endif
496 
497  end subroutine hecmw_tvec_33
498 
499  !C
500  !C***
501  !C*** hecmw_Ttvec_33
502  !C***
503  !C
504  subroutine hecmw_ttvec_33 (hecMESH, X, Y, COMMtime)
505  use hecmw_util
506  use m_hecmw_comm_f
507  implicit none
508  type (hecmwst_local_mesh), intent(in) :: hecmesh
509  real(kind=kreal), intent(in) :: x(:)
510  real(kind=kreal), intent(out) :: y(:)
511  real(kind=kreal), intent(inout) :: commtime
512 
513  real(kind=kreal) :: start_time, end_time
514  integer(kind=kint) :: i, j, jj, k, kk
515 
516  start_time= hecmw_wtime()
517  call hecmw_update_r (hecmesh, x, hecmesh%n_node, 3)
518  end_time= hecmw_wtime()
519  commtime = commtime + end_time - start_time
520 
521 #ifdef _OPENACC
522  !$acc kernels
523  !$acc loop independent private(i)
524 #else
525  !$omp parallel default(none),private(i,k,kk,j,jj),shared(hecMESH,X,Y)
526  !$omp do
527 #endif
528  do i= 1, hecmesh%nn_internal * hecmesh%n_dof
529  y(i)= x(i)
530  enddo
531 #ifndef _OPENACC
532  !$omp end do
533 #endif
534 
535 #ifdef _OPENACC
536  !$acc loop independent private(i,k,kk,j,jj)
537 #else
538  !$omp do
539 #endif
540  outer: do i= 1, hecmesh%mpc%n_mpc
541  do j= hecmesh%mpc%mpc_index(i-1) + 1, hecmesh%mpc%mpc_index(i)
542  if (hecmesh%mpc%mpc_dof(j) > 3) cycle outer
543  enddo
544  k = hecmesh%mpc%mpc_index(i-1) + 1
545  kk = 3 * (hecmesh%mpc%mpc_item(k) - 1) + hecmesh%mpc%mpc_dof(k)
546  y(kk) = 0.d0
547  do j= hecmesh%mpc%mpc_index(i-1) + 2, hecmesh%mpc%mpc_index(i)
548  jj = 3 * (hecmesh%mpc%mpc_item(j) - 1) + hecmesh%mpc%mpc_dof(j)
549 #ifdef _OPENACC
550  !$acc atomic update
551 #else
552  !$omp atomic
553 #endif
554  y(jj) = y(jj) - hecmesh%mpc%mpc_val(j) * x(kk)
555  !$acc end atomic
556  enddo
557  enddo outer
558 #ifdef _OPENACC
559  !$acc end kernels
560 #else
561  !$omp end do
562  !$omp end parallel
563 #endif
564 
565  end subroutine hecmw_ttvec_33
566 
567  !C
568  !C***
569  !C*** hecmw_TtmatTvec_33
570  !C***
571  !C
572  subroutine hecmw_ttmattvec_33 (hecMESH, hecMAT, X, Y, W, time_Ax, COMMtime)
573  use hecmw_util
574  implicit none
575  type (hecmwst_local_mesh), intent(in) :: hecmesh
576  type (hecmwst_matrix), intent(in) :: hecmat
577  real(kind=kreal), intent(in) :: x(:)
578  real(kind=kreal), intent(out) :: y(:), w(:)
579  real(kind=kreal), intent(inout) :: time_ax
580  real(kind=kreal), intent(inout) :: commtime
581 
582  call hecmw_tvec_33(hecmesh, x, y, commtime)
583  call hecmw_matvec_33_inner(hecmesh, hecmat, y, w, time_ax, commtime)
584  call hecmw_ttvec_33(hecmesh, w, y, commtime)
585 
586  end subroutine hecmw_ttmattvec_33
587 
588 
589  !C
590  !C***
591  !C*** hecmw_mat_diag_sr_33
592  !C***
593  !C
594  subroutine hecmw_mat_diag_sr_33(hecMESH, hecMAT, COMMtime)
595  use hecmw_util
596  use m_hecmw_comm_f
597  implicit none
598  type (hecmwst_local_mesh), intent(in) :: hecmesh
599  type (hecmwst_matrix), intent(inout), target :: hecmat
600  real(kind=kreal), intent(inout), optional :: commtime
601  real(kind=kreal), allocatable :: w(:,:)
602  real(kind=kreal), pointer :: d(:)
603  integer(kind=kint) :: ip
604  real(kind=kreal) :: start_time, end_time
605  allocate(w(3*hecmat%NP,3))
606  d => hecmat%D
607  do ip= 1, hecmat%N
608  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)
609  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)
610  w(3*ip ,1)= d(9*ip-2); w(3*ip ,2)= d(9*ip-1); w(3*ip ,3)= d(9*ip )
611  enddo
612  start_time= hecmw_wtime()
613  call hecmw_update_r (hecmesh, w(:,1), hecmat%NP, 3)
614  call hecmw_update_r (hecmesh, w(:,2), hecmat%NP, 3)
615  call hecmw_update_r (hecmesh, w(:,3), hecmat%NP, 3)
616  end_time= hecmw_wtime()
617  if (present(commtime)) commtime = commtime + end_time - start_time
618  do ip= hecmat%N+1, hecmat%NP
619  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)
620  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)
621  d(9*ip-2)= w(3*ip ,1); d(9*ip-1)= w(3*ip ,2); d(9*ip )= w(3*ip ,3)
622  enddo
623  deallocate(w)
624  end subroutine hecmw_mat_diag_sr_33
625 
626 end module hecmw_solver_las_33
Jagged Diagonal Matrix storage for vector processors. Original code was provided by JAMSTEC.
Definition: hecmw_jadm.f90:8
subroutine, public hecmw_jad_matvec(hecMESH, hecMAT, X, Y, COMMtime)
Definition: hecmw_jadm.f90:61
integer(kind=kint) function, public hecmw_jad_is_initialized()
Definition: hecmw_jadm.f90:56
integer(kind=kint) function, public hecmw_mat_get_flag_mpcmatvec(hecMAT)
subroutine, public hecmw_matresid_33(hecMESH, hecMAT, X, B, R, time_Ax, COMMtime)
subroutine, public hecmw_matvec_33(hecMESH, hecMAT, X, Y, time_Ax, COMMtime)
subroutine, public hecmw_mat_diag_sr_33(hecMESH, hecMAT, COMMtime)
subroutine, public hecmw_matvec_33_unset_async
subroutine, public hecmw_ttmattvec_33(hecMESH, hecMAT, X, Y, W, time_Ax, COMMtime)
subroutine, public hecmw_ttvec_33(hecMESH, X, Y, COMMtime)
real(kind=kreal) function, public hecmw_rel_resid_l2_33(hecMESH, hecMAT, time_Ax, COMMtime)
subroutine, public hecmw_tvec_33(hecMESH, X, Y, COMMtime)
subroutine, public hecmw_matvec_33_set_async(hecMAT)
subroutine hecmw_innerproduct_r(hecMESH, ndof, X, Y, sum, COMMtime)
subroutine, public hecmw_tuning_fx_calc_sector_cache(N, NDOF, sectorCacheSize0, sectorCacheSize1)
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)