FrontISTR  5.9.0
Large-scale structural analysis program with finit element method
hecmw_solver_las_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 
7  use hecmw_util
8  implicit none
9 
10  private
11 
12  public :: hecmw_matvec_44
15  public :: hecmw_matresid_44
16  public :: hecmw_rel_resid_l2_44
17  public :: hecmw_tvec_44
18  public :: hecmw_ttvec_44
19  public :: hecmw_ttmattvec_44
20  public :: hecmw_mat_diag_sr_44
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_44
32  !C***
33  !C
34  subroutine hecmw_matvec_44 (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_44(hecmesh, hecmat, x, y, wk, time_ax, tcomm)
53  deallocate(wk)
54  else
55  call hecmw_matvec_44_inner(hecmesh, hecmat, x, y, time_ax, tcomm)
56  endif
57 
58  if (present(commtime)) commtime = commtime + tcomm
59  end subroutine hecmw_matvec_44
60 
61  !C
62  !C***
63  !C*** hecmw_matvec_44_set_async
64  !C***
65  !C
66  subroutine hecmw_matvec_44_set_async (hecMAT)
67  use hecmw_util
68  implicit none
69  type (hecmwst_matrix), intent(in) :: hecmat
70 
71  end subroutine hecmw_matvec_44_set_async
72 
73  !C
74  !C***
75  !C*** hecmw_matvec_44_unset_async
76  !C***
77  !C
79  implicit none
80 
81  end subroutine hecmw_matvec_44_unset_async
82 
83  !C
84  !C***
85  !C*** hecmw_matvec_44_inner ( private subroutine )
86  !C***
87  !C
88  subroutine hecmw_matvec_44_inner (hecMESH, hecMAT, X, Y, time_Ax, COMMtime)
89  use hecmw_util
90  use m_hecmw_comm_f
92  use hecmw_jad_type
93  use hecmw_tuning_fx
94  !$ use omp_lib
95 
96  implicit none
97  type (hecmwst_local_mesh), intent(in) :: hecmesh
98  type (hecmwst_matrix), intent(in), target :: hecmat
99  real(kind=kreal), intent(in) :: x(:)
100  real(kind=kreal), intent(out) :: y(:)
101  real(kind=kreal), intent(inout) :: time_ax
102  real(kind=kreal), intent(inout), optional :: commtime
103 
104  real(kind=kreal) :: start_time, end_time, tcomm
105  integer(kind=kint) :: i, j, js, je, in
106  real(kind=kreal) :: yv1, yv2, yv3, yv4, x1, x2, x3, x4
107 
108  integer(kind=kint) :: n, np
109  integer(kind=kint), pointer :: indexl(:), iteml(:), indexu(:), itemu(:), indexa(:), itema(:)
110  real(kind=kreal), pointer :: al(:), au(:), d(:), a(:)
111 
112  ! added for turning >>>
113  integer, parameter :: numofblockperthread = 100
114  logical, save :: isfirst = .true.
115  integer, save :: numofthread = 1
116  integer, save, allocatable :: startpos(:), endpos(:)
117  integer(kind=kint), save :: sectorcachesize0, sectorcachesize1
118  integer(kind=kint) :: threadnum, blocknum, numofblock
119  integer(kind=kint) :: numofelement, elementcount, blockindex
120  real(kind=kreal) :: numofelementperblock
121  ! <<< added for turning
122 
123  if (hecmw_jad_is_initialized().ne.0) then
124  tcomm = 0.d0
125  start_time = hecmw_wtime()
126  call hecmw_jad_matvec(hecmesh, hecmat, x, y, tcomm)
127  end_time = hecmw_wtime()
128  time_ax = time_ax + end_time - start_time - tcomm
129  if (present(commtime)) commtime = commtime + tcomm
130  else
131 
132  n = hecmat%N
133  np = hecmat%NP
134  indexl => hecmat%indexL
135  indexu => hecmat%indexU
136  indexa => hecmat%indexA
137  iteml => hecmat%itemL
138  itemu => hecmat%itemU
139  itema => hecmat%itemA
140  al => hecmat%AL
141  au => hecmat%AU
142  d => hecmat%D
143  a => hecmat%A
144 
145  ! added for turning >>>
146 #ifndef _OPENACC
147  if (.not. isfirst) then
148  numofblock = numofthread * numofblockperthread
149  if (endpos(numofblock-1) .ne. n-1) then
150  deallocate(startpos, endpos)
151  isfirst = .true.
152  endif
153  endif
154  if (isfirst) then
155  !$ numOfThread = omp_get_max_threads()
156  numofblock = numofthread * numofblockperthread
157  allocate (startpos(0 : numofblock - 1), endpos(0 : numofblock - 1))
158  numofelement = n + indexl(n) + indexu(n)
159  numofelementperblock = dble(numofelement) / numofblock
160  blocknum = 0
161  elementcount = 0
162  startpos(blocknum) = 1
163  do i= 1, n
164  elementcount = elementcount + 1
165  elementcount = elementcount + (indexl(i) - indexl(i-1))
166  elementcount = elementcount + (indexu(i) - indexu(i-1))
167  if (elementcount > (blocknum + 1) * numofelementperblock) then
168  endpos(blocknum) = i
169  ! write(9000+hecMESH%my_rank,*) mod(blockNum, numOfThread), &
170  ! startPos(blockNum), endPos(blockNum)
171  blocknum = blocknum + 1
172  startpos(blocknum) = i + 1
173  if (blocknum == (numofblock - 1)) exit
174  endif
175  enddo
176  endpos(blocknum) = n
177  ! write(9000+hecMESH%my_rank,*) mod(blockNum, numOfThread), &
178  ! startPos(blockNum), endPos(blockNum)
179  ! for irregular data
180  do i= blocknum+1, numofblock-1
181  startpos(i) = n
182  endpos(i) = n-1
183  ! write(9000+hecMESH%my_rank,*) mod(i, numOfThread), &
184  ! startPos(i), endPos(i)
185  end do
186 
188  sectorcachesize0, sectorcachesize1)
189 
190  isfirst = .false.
191  endif
192 #endif
193  ! <<< added for turning
194 
195  start_time= hecmw_wtime()
196 
197  call hecmw_update_r (hecmesh, x, np, 4)
198 
199  ! endif
200  end_time= hecmw_wtime()
201  if (present(commtime)) commtime = commtime + end_time - start_time
202 
203  start_time = hecmw_wtime()
204 
205  !call fapp_start("loopInMatvec44", 1, 0)
206  !call start_collection("loopInMatvec44")
207 
208  !OCL CACHE_SECTOR_SIZE(sectorCacheSize0,sectorCacheSize1)
209  !OCL CACHE_SUBSECTOR_ASSIGN(X)
210 
211 #ifdef _OPENACC
212  !$acc kernels
213  !$acc loop independent
214  do i = 1, n
215  x1= x(4*i-3)
216  x2= x(4*i-2)
217  x3= x(4*i-1)
218  x4= x(4*i )
219  yv1= 0
220  yv2= 0
221  yv3= 0
222  yv4= 0
223  js= indexa(i) + 1
224  je= indexa(i+1)
225  do j= js, je
226  in = itema(j)
227  x1= x(4*in-3)
228  x2= x(4*in-2)
229  x3= x(4*in-1)
230  x4= x(4*in )
231  yv1= yv1 + a(16*j-15)*x1 + a(16*j-14)*x2 + a(16*j-13)*x3 + a(16*j-12)*x4
232  yv2= yv2 + a(16*j-11)*x1 + a(16*j-10)*x2 + a(16*j- 9)*x3 + a(16*j- 8)*x4
233  yv3= yv3 + a(16*j- 7)*x1 + a(16*j- 6)*x2 + a(16*j- 5)*x3 + a(16*j- 4)*x4
234  yv4= yv4 + a(16*j- 3)*x1 + a(16*j- 2)*x2 + a(16*j- 1)*x3 + a(16*j )*x4
235  enddo
236  y(4*i-3)= yv1
237  y(4*i-2)= yv2
238  y(4*i-1)= yv3
239  y(4*i )= yv4
240  enddo
241  !$acc end kernels
242 #else
243  !$OMP PARALLEL DEFAULT(NONE) &
244  !$OMP&PRIVATE(i,X1,X2,X3,X4,YV1,YV2,YV3,YV4,jS,jE,j,in,threadNum,blockNum,blockIndex) &
245  !$OMP&SHARED(D,AL,AU,indexL,itemL,indexU,itemU,X,Y,startPos,endPos,numOfThread,N,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  x1= x(4*i-3)
252  x2= x(4*i-2)
253  x3= x(4*i-1)
254  x4= x(4*i )
255  yv1= d(16*i-15)*x1 + d(16*i-14)*x2 + d(16*i-13)*x3 + d(16*i-12)*x4
256  yv2= d(16*i-11)*x1 + d(16*i-10)*x2 + d(16*i- 9)*x3 + d(16*i- 8)*x4
257  yv3= d(16*i- 7)*x1 + d(16*i- 6)*x2 + d(16*i- 5)*x3 + d(16*i- 4)*x4
258  yv4= d(16*i- 3)*x1 + d(16*i- 2)*x2 + d(16*i- 1)*x3 + d(16*i )*x4
259 
260  js= indexl(i-1) + 1
261  je= indexl(i )
262  do j= js, je
263  in = iteml(j)
264  x1= x(4*in-3)
265  x2= x(4*in-2)
266  x3= x(4*in-1)
267  x4= x(4*in )
268  yv1= yv1 + al(16*j-15)*x1 + al(16*j-14)*x2 + al(16*j-13)*x3 + al(16*j-12)*x4
269  yv2= yv2 + al(16*j-11)*x1 + al(16*j-10)*x2 + al(16*j- 9)*x3 + al(16*j- 8)*x4
270  yv3= yv3 + al(16*j- 7)*x1 + al(16*j- 6)*x2 + al(16*j- 5)*x3 + al(16*j- 4)*x4
271  yv4= yv4 + al(16*j- 3)*x1 + al(16*j- 2)*x2 + al(16*j- 1)*x3 + al(16*j )*x4
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  x1= x(4*in-3)
279  x2= x(4*in-2)
280  x3= x(4*in-1)
281  x4= x(4*in )
282  yv1= yv1 + au(16*j-15)*x1 + au(16*j-14)*x2 + au(16*j-13)*x3 + au(16*j-12)*x4
283  yv2= yv2 + au(16*j-11)*x1 + au(16*j-10)*x2 + au(16*j- 9)*x3 + au(16*j- 8)*x4
284  yv3= yv3 + au(16*j- 7)*x1 + au(16*j- 6)*x2 + au(16*j- 5)*x3 + au(16*j- 4)*x4
285  yv4= yv4 + au(16*j- 3)*x1 + au(16*j- 2)*x2 + au(16*j- 1)*x3 + au(16*j )*x4
286  enddo
287  y(4*i-3)= yv1
288  y(4*i-2)= yv2
289  y(4*i-1)= yv3
290  y(4*i )= yv4
291  enddo
292  enddo
293  !$OMP END PARALLEL
294 #endif
295 
296  !OCL END_CACHE_SUBSECTOR
297  !OCL END_CACHE_SECTOR_SIZE
298 
299  !call stop_collection("loopInMatvec44")
300  !call fapp_stop("loopInMatvec44", 1, 0)
301 
302  end_time = hecmw_wtime()
303  time_ax = time_ax + end_time - start_time
304 
305  endif
306  end subroutine hecmw_matvec_44_inner
307 
308  !C
309  !C***
310  !C*** hecmw_matresid_44
311  !C***
312  !C
313  subroutine hecmw_matresid_44 (hecMESH, hecMAT, X, B, R, time_Ax, COMMtime)
314  use hecmw_util
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_44 (hecmesh, hecmat, x, r, time_ax, tcomm)
328  if (present(commtime)) commtime = commtime + tcomm
329 #ifdef _OPENACC
330  !$acc kernels
331  !$acc loop independent
332 #else
333  !$omp parallel default(none),private(i),shared(hecMAT,R,B)
334  !$omp do
335 #endif
336  do i = 1, hecmat%N * 4
337  r(i) = b(i) - r(i)
338  enddo
339 #ifdef _OPENACC
340  !$acc end kernels
341 #else
342  !$omp end do
343  !$omp end parallel
344 #endif
345  end subroutine hecmw_matresid_44
346 
347  !C
348  !C***
349  !C*** hecmw_rel_resid_L2_44
350  !C***
351  !C
352  function hecmw_rel_resid_l2_44 (hecMESH, hecMAT, time_Ax, COMMtime)
353  use hecmw_util
355  implicit none
356  real(kind=kreal) :: hecmw_rel_resid_l2_44
357  type ( hecmwst_local_mesh ), intent(in) :: hecmesh
358  type ( hecmwst_matrix ), intent(in) :: hecmat
359  real(kind=kreal), intent(inout) :: time_ax
360  real(kind=kreal), intent(inout), optional :: commtime
361 
362  real(kind=kreal), allocatable :: r(:)
363  real(kind=kreal) :: bnorm2, rnorm2
364  real(kind=kreal) :: tcomm
365 
366  allocate(r(hecmat%NDOF*hecmat%NP))
367 
368  tcomm = 0.d0
369  call hecmw_innerproduct_r(hecmesh, hecmat%NDOF, &
370  hecmat%B, hecmat%B, bnorm2, tcomm)
371  if (bnorm2 == 0.d0) then
372  bnorm2 = 1.d0
373  endif
374  call hecmw_matresid_44(hecmesh, hecmat, hecmat%X, hecmat%B, r, time_ax, tcomm)
375  call hecmw_innerproduct_r(hecmesh, hecmat%NDOF, r, r, rnorm2, tcomm)
376  hecmw_rel_resid_l2_44 = sqrt(rnorm2 / bnorm2)
377 
378  if (present(commtime)) commtime = commtime + tcomm
379 
380  deallocate(r)
381  end function hecmw_rel_resid_l2_44
382 
383  !C
384  !C***
385  !C*** hecmw_Tvec_44
386  !C***
387  !C
388  subroutine hecmw_tvec_44 (hecMESH, X, Y, COMMtime)
389  use hecmw_util
390  use m_hecmw_comm_f
391  implicit none
392  type (hecmwst_local_mesh), intent(in) :: hecmesh
393  real(kind=kreal), intent(in) :: x(:)
394  real(kind=kreal), intent(out) :: y(:)
395  real(kind=kreal), intent(inout) :: commtime
396 
397  real(kind=kreal) :: start_time, end_time
398  integer(kind=kint) :: i, j, jj, k, kk
399 
400  start_time= hecmw_wtime()
401  call hecmw_update_r (hecmesh, x, hecmesh%n_node, 4)
402  end_time= hecmw_wtime()
403  commtime = commtime + end_time - start_time
404 
405 #ifdef _OPENACC
406  !$acc kernels
407  !$acc loop independent
408 #else
409  !$omp parallel default(none),private(i,k,kk,j,jj),shared(hecMESH,X,Y)
410  !$omp do
411 #endif
412  do i= 1, hecmesh%nn_internal * hecmesh%n_dof
413  y(i)= x(i)
414  enddo
415 #ifndef _OPENACC
416  !$omp end do
417 #endif
418 
419 #ifdef _OPENACC
420  !$acc loop independent
421 #else
422  !$omp do
423 #endif
424  outer: do i= 1, hecmesh%mpc%n_mpc
425  do j= hecmesh%mpc%mpc_index(i-1) + 1, hecmesh%mpc%mpc_index(i)
426  if (hecmesh%mpc%mpc_dof(j) > 4) cycle outer
427  enddo
428  k = hecmesh%mpc%mpc_index(i-1) + 1
429  kk = 4 * (hecmesh%mpc%mpc_item(k) - 1) + hecmesh%mpc%mpc_dof(k)
430  y(kk) = 0.d0
431  do j= hecmesh%mpc%mpc_index(i-1) + 2, hecmesh%mpc%mpc_index(i)
432  jj = 4 * (hecmesh%mpc%mpc_item(j) - 1) + hecmesh%mpc%mpc_dof(j)
433  y(kk) = y(kk) - hecmesh%mpc%mpc_val(j) * x(jj)
434  enddo
435  enddo outer
436 #ifdef _OPENACC
437  !$acc end kernels
438 #else
439  !$omp end do
440  !$omp end parallel
441 #endif
442 
443  end subroutine hecmw_tvec_44
444 
445  !C
446  !C***
447  !C*** hecmw_Ttvec_44
448  !C***
449  !C
450  subroutine hecmw_ttvec_44 (hecMESH, X, Y, COMMtime)
451  use hecmw_util
452  use m_hecmw_comm_f
453  implicit none
454  type (hecmwst_local_mesh), intent(in) :: hecmesh
455  real(kind=kreal), intent(in) :: x(:)
456  real(kind=kreal), intent(out) :: y(:)
457  real(kind=kreal), intent(inout) :: commtime
458 
459  real(kind=kreal) :: start_time, end_time
460  integer(kind=kint) :: i, j, jj, k, kk
461 
462  start_time= hecmw_wtime()
463  call hecmw_update_r (hecmesh, x, hecmesh%n_node, 4)
464  end_time= hecmw_wtime()
465  commtime = commtime + end_time - start_time
466 
467 #ifdef _OPENACC
468  !$acc kernels
469  !$acc loop independent
470 #else
471  !$omp parallel default(none),private(i,k,kk,j,jj),shared(hecMESH,X,Y)
472  !$omp do
473 #endif
474  do i= 1, hecmesh%nn_internal * hecmesh%n_dof
475  y(i)= x(i)
476  enddo
477 #ifndef _OPENACC
478  !$omp end do
479 #endif
480 
481 #ifdef _OPENACC
482  !$acc loop independent
483 #else
484  !$omp do
485 #endif
486  outer: do i= 1, hecmesh%mpc%n_mpc
487  do j= hecmesh%mpc%mpc_index(i-1) + 1, hecmesh%mpc%mpc_index(i)
488  if (hecmesh%mpc%mpc_dof(j) > 4) cycle outer
489  enddo
490  k = hecmesh%mpc%mpc_index(i-1) + 1
491  kk = 4 * (hecmesh%mpc%mpc_item(k) - 1) + hecmesh%mpc%mpc_dof(k)
492  y(kk) = 0.d0
493  do j= hecmesh%mpc%mpc_index(i-1) + 2, hecmesh%mpc%mpc_index(i)
494  jj = 4 * (hecmesh%mpc%mpc_item(j) - 1) + hecmesh%mpc%mpc_dof(j)
495 #ifdef _OPENACC
496  !$acc atomic update
497 #else
498  !omp atomic
499 #endif
500  y(jj) = y(jj) - hecmesh%mpc%mpc_val(j) * x(kk)
501  enddo
502  enddo outer
503 #ifdef _OPENACC
504  !$acc end kernels
505 #else
506  !$omp end do
507  !$omp end parallel
508 #endif
509 
510  end subroutine hecmw_ttvec_44
511 
512  !C
513  !C***
514  !C*** hecmw_TtmatTvec_44
515  !C***
516  !C
517  subroutine hecmw_ttmattvec_44 (hecMESH, hecMAT, X, Y, W, time_Ax, COMMtime)
518  use hecmw_util
519  implicit none
520  type (hecmwst_local_mesh), intent(in) :: hecmesh
521  type (hecmwst_matrix), intent(in) :: hecmat
522  real(kind=kreal), intent(in) :: x(:)
523  real(kind=kreal), intent(out) :: y(:), w(:)
524  real(kind=kreal), intent(inout) :: time_ax
525  real(kind=kreal), intent(inout) :: commtime
526 
527  call hecmw_tvec_44(hecmesh, x, y, commtime)
528  call hecmw_matvec_44_inner(hecmesh, hecmat, y, w, time_ax, commtime)
529  call hecmw_ttvec_44(hecmesh, w, y, commtime)
530 
531  end subroutine hecmw_ttmattvec_44
532 
533 
534  !C
535  !C***
536  !C*** hecmw_mat_diag_sr_44
537  !C***
538  !C
539  subroutine hecmw_mat_diag_sr_44(hecMESH, hecMAT, COMMtime)
540  use hecmw_util
541  use m_hecmw_comm_f
542  implicit none
543  type (hecmwst_local_mesh), intent(in) :: hecmesh
544  type (hecmwst_matrix), intent(inout), target :: hecmat
545  real(kind=kreal), intent(inout), optional :: commtime
546  real(kind=kreal), allocatable :: w(:,:)
547  real(kind=kreal), pointer :: d(:)
548  integer(kind=kint) :: ip
549  real(kind=kreal) :: start_time, end_time
550  allocate(w(4*hecmat%NP,4))
551  d => hecmat%D
552  do ip= 1, hecmat%N
553  w(4*ip-3,1)= d(16*ip-15); w(4*ip-3,2)= d(16*ip-14); w(4*ip-3,3)= d(16*ip-13); w(4*ip-3,4)= d(16*ip-12)
554  w(4*ip-2,1)= d(16*ip-11); w(4*ip-2,2)= d(16*ip-10); w(4*ip-2,3)= d(16*ip- 9); w(4*ip-2,4)= d(16*ip- 8)
555  w(4*ip-1,1)= d(16*ip- 7); w(4*ip-1,2)= d(16*ip- 6); w(4*ip-1,3)= d(16*ip- 5); w(4*ip-1,4)= d(16*ip- 4)
556  w(4*ip ,1)= d(16*ip- 3); w(4*ip ,2)= d(16*ip- 2); w(4*ip ,3)= d(16*ip- 1); w(4*ip ,4)= d(16*ip )
557  enddo
558  start_time= hecmw_wtime()
559  call hecmw_update_r (hecmesh, w(:,1), hecmat%NP, 4)
560  call hecmw_update_r (hecmesh, w(:,2), hecmat%NP, 4)
561  call hecmw_update_r (hecmesh, w(:,3), hecmat%NP, 4)
562  call hecmw_update_r (hecmesh, w(:,4), hecmat%NP, 4)
563  end_time= hecmw_wtime()
564  if (present(commtime)) commtime = commtime + end_time - start_time
565  do ip= hecmat%N+1, hecmat%NP
566  d(16*ip-15)= w(4*ip-3,1); d(16*ip-14)= w(4*ip-3,2); d(16*ip-13)= w(4*ip-3,3); d(16*ip-12)= w(4*ip-3,4)
567  d(16*ip-11)= w(4*ip-2,1); d(16*ip-10)= w(4*ip-2,2); d(16*ip- 9)= w(4*ip-2,3); d(16*ip- 8)= w(4*ip-2,4)
568  d(16*ip- 7)= w(4*ip-1,1); d(16*ip- 6)= w(4*ip-1,2); d(16*ip- 5)= w(4*ip-1,3); d(16*ip- 4)= w(4*ip-1,4)
569  d(16*ip- 3)= w(4*ip ,1); d(16*ip- 2)= w(4*ip ,2); d(16*ip- 1)= w(4*ip ,3); d(16*ip )= w(4*ip ,4)
570  enddo
571  deallocate(w)
572  end subroutine hecmw_mat_diag_sr_44
573 
574 end module hecmw_solver_las_44
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_44(hecMESH, hecMAT, X, B, R, time_Ax, COMMtime)
subroutine, public hecmw_mat_diag_sr_44(hecMESH, hecMAT, COMMtime)
subroutine, public hecmw_matvec_44(hecMESH, hecMAT, X, Y, time_Ax, COMMtime)
subroutine, public hecmw_ttvec_44(hecMESH, X, Y, COMMtime)
subroutine, public hecmw_matvec_44_unset_async
subroutine, public hecmw_matvec_44_set_async(hecMAT)
subroutine, public hecmw_tvec_44(hecMESH, X, Y, COMMtime)
real(kind=kreal) function, public hecmw_rel_resid_l2_44(hecMESH, hecMAT, time_Ax, COMMtime)
subroutine, public hecmw_ttmattvec_44(hecMESH, hecMAT, X, Y, W, time_Ax, COMMtime)
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)