FrontISTR  5.7.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)
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)
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(:)
110  real(kind=kreal), pointer :: al(:), au(:), d(:)
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  iteml => hecmat%itemL
137  itemu => hecmat%itemU
138  al => hecmat%AL
139  au => hecmat%AU
140  d => hecmat%D
141 
142  ! added for turning >>>
143  if (.not. isfirst) then
144  numofblock = numofthread * numofblockperthread
145  if (endpos(numofblock-1) .ne. n-1) then
146  deallocate(startpos, endpos)
147  isfirst = .true.
148  endif
149  endif
150  if (isfirst) then
151  !$ numOfThread = omp_get_max_threads()
152  numofblock = numofthread * numofblockperthread
153  allocate (startpos(0 : numofblock - 1), endpos(0 : numofblock - 1))
154  numofelement = n + indexl(n) + indexu(n)
155  numofelementperblock = dble(numofelement) / numofblock
156  blocknum = 0
157  elementcount = 0
158  startpos(blocknum) = 1
159  do i= 1, n
160  elementcount = elementcount + 1
161  elementcount = elementcount + (indexl(i) - indexl(i-1))
162  elementcount = elementcount + (indexu(i) - indexu(i-1))
163  if (elementcount > (blocknum + 1) * numofelementperblock) then
164  endpos(blocknum) = i
165  ! write(9000+hecMESH%my_rank,*) mod(blockNum, numOfThread), &
166  ! startPos(blockNum), endPos(blockNum)
167  blocknum = blocknum + 1
168  startpos(blocknum) = i + 1
169  if (blocknum == (numofblock - 1)) exit
170  endif
171  enddo
172  endpos(blocknum) = n
173  ! write(9000+hecMESH%my_rank,*) mod(blockNum, numOfThread), &
174  ! startPos(blockNum), endPos(blockNum)
175  ! for irregular data
176  do i= blocknum+1, numofblock-1
177  startpos(i) = n
178  endpos(i) = n-1
179  ! write(9000+hecMESH%my_rank,*) mod(i, numOfThread), &
180  ! startPos(i), endPos(i)
181  end do
182 
184  sectorcachesize0, sectorcachesize1)
185 
186  isfirst = .false.
187  endif
188  ! <<< added for turning
189 
190  start_time= hecmw_wtime()
191 
192  call hecmw_update_r (hecmesh, x, np, 4)
193 
194  ! endif
195  end_time= hecmw_wtime()
196  if (present(commtime)) commtime = commtime + end_time - start_time
197 
198  start_time = hecmw_wtime()
199 
200  !call fapp_start("loopInMatvec44", 1, 0)
201  !call start_collection("loopInMatvec44")
202 
203  !OCL CACHE_SECTOR_SIZE(sectorCacheSize0,sectorCacheSize1)
204  !OCL CACHE_SUBSECTOR_ASSIGN(X)
205 
206  !$OMP PARALLEL DEFAULT(NONE) &
207  !$OMP&PRIVATE(i,X1,X2,X3,X4,YV1,YV2,YV3,YV4,jS,jE,j,in,threadNum,blockNum,blockIndex) &
208  !$OMP&SHARED(D,AL,AU,indexL,itemL,indexU,itemU,X,Y,startPos,endPos,numOfThread,N,async_matvec_flg)
209  threadnum = 0
210  !$ threadNum = omp_get_thread_num()
211  do blocknum = 0 , numofblockperthread - 1
212  blockindex = blocknum * numofthread + threadnum
213  do i = startpos(blockindex), endpos(blockindex)
214  x1= x(4*i-3)
215  x2= x(4*i-2)
216  x3= x(4*i-1)
217  x4= x(4*i )
218  yv1= d(16*i-15)*x1 + d(16*i-14)*x2 + d(16*i-13)*x3 + d(16*i-12)*x4
219  yv2= d(16*i-11)*x1 + d(16*i-10)*x2 + d(16*i- 9)*x3 + d(16*i- 8)*x4
220  yv3= d(16*i- 7)*x1 + d(16*i- 6)*x2 + d(16*i- 5)*x3 + d(16*i- 4)*x4
221  yv4= d(16*i- 3)*x1 + d(16*i- 2)*x2 + d(16*i- 1)*x3 + d(16*i )*x4
222 
223  js= indexl(i-1) + 1
224  je= indexl(i )
225  do j= js, je
226  in = iteml(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 + al(16*j-15)*x1 + al(16*j-14)*x2 + al(16*j-13)*x3 + al(16*j-12)*x4
232  yv2= yv2 + al(16*j-11)*x1 + al(16*j-10)*x2 + al(16*j- 9)*x3 + al(16*j- 8)*x4
233  yv3= yv3 + al(16*j- 7)*x1 + al(16*j- 6)*x2 + al(16*j- 5)*x3 + al(16*j- 4)*x4
234  yv4= yv4 + al(16*j- 3)*x1 + al(16*j- 2)*x2 + al(16*j- 1)*x3 + al(16*j )*x4
235  enddo
236  js= indexu(i-1) + 1
237  je= indexu(i )
238  do j= js, je
239  in = itemu(j)
240  ! if (async_matvec_flg .and. in > N) cycle
241  x1= x(4*in-3)
242  x2= x(4*in-2)
243  x3= x(4*in-1)
244  x4= x(4*in )
245  yv1= yv1 + au(16*j-15)*x1 + au(16*j-14)*x2 + au(16*j-13)*x3 + au(16*j-12)*x4
246  yv2= yv2 + au(16*j-11)*x1 + au(16*j-10)*x2 + au(16*j- 9)*x3 + au(16*j- 8)*x4
247  yv3= yv3 + au(16*j- 7)*x1 + au(16*j- 6)*x2 + au(16*j- 5)*x3 + au(16*j- 4)*x4
248  yv4= yv4 + au(16*j- 3)*x1 + au(16*j- 2)*x2 + au(16*j- 1)*x3 + au(16*j )*x4
249  enddo
250  y(4*i-3)= yv1
251  y(4*i-2)= yv2
252  y(4*i-1)= yv3
253  y(4*i )= yv4
254  enddo
255  enddo
256  !$OMP END PARALLEL
257 
258  !OCL END_CACHE_SUBSECTOR
259  !OCL END_CACHE_SECTOR_SIZE
260 
261  !call stop_collection("loopInMatvec44")
262  !call fapp_stop("loopInMatvec44", 1, 0)
263 
264  end_time = hecmw_wtime()
265  time_ax = time_ax + end_time - start_time
266 
267  endif
268  end subroutine hecmw_matvec_44_inner
269 
270  !C
271  !C***
272  !C*** hecmw_matresid_44
273  !C***
274  !C
275  subroutine hecmw_matresid_44 (hecMESH, hecMAT, X, B, R, time_Ax, COMMtime)
277  implicit none
278  type (hecmwst_local_mesh), intent(in) :: hecmesh
279  type (hecmwst_matrix), intent(in) :: hecmat
280  real(kind=kreal), intent(in) :: x(:), b(:)
281  real(kind=kreal), intent(out) :: r(:)
282  real(kind=kreal), intent(inout) :: time_ax
283  real(kind=kreal), intent(inout), optional :: commtime
284 
285  integer(kind=kint) :: i
286  real(kind=kreal) :: tcomm
287 
288  tcomm = 0.d0
289  call hecmw_matvec_44 (hecmesh, hecmat, x, r, time_ax, tcomm)
290  if (present(commtime)) commtime = commtime + tcomm
291  !$omp parallel default(none),private(i),shared(hecMAT,R,B)
292  !$omp do
293  do i = 1, hecmat%N * 4
294  r(i) = b(i) - r(i)
295  enddo
296  !$omp end do
297  !$omp end parallel
298  end subroutine hecmw_matresid_44
299 
300  !C
301  !C***
302  !C*** hecmw_rel_resid_L2_44
303  !C***
304  !C
305  function hecmw_rel_resid_l2_44 (hecMESH, hecMAT, time_Ax, COMMtime)
308  implicit none
309  real(kind=kreal) :: hecmw_rel_resid_l2_44
310  type ( hecmwst_local_mesh ), intent(in) :: hecmesh
311  type ( hecmwst_matrix ), intent(in) :: hecmat
312  real(kind=kreal), intent(inout) :: time_ax
313  real(kind=kreal), intent(inout), optional :: commtime
314 
315  real(kind=kreal), allocatable :: r(:)
316  real(kind=kreal) :: bnorm2, rnorm2
317  real(kind=kreal) :: tcomm
318 
319  allocate(r(hecmat%NDOF*hecmat%NP))
320 
321  tcomm = 0.d0
322  call hecmw_innerproduct_r(hecmesh, hecmat%NDOF, &
323  hecmat%B, hecmat%B, bnorm2, tcomm)
324  if (bnorm2 == 0.d0) then
325  bnorm2 = 1.d0
326  endif
327  call hecmw_matresid_44(hecmesh, hecmat, hecmat%X, hecmat%B, r, time_ax, tcomm)
328  call hecmw_innerproduct_r(hecmesh, hecmat%NDOF, r, r, rnorm2, tcomm)
329  hecmw_rel_resid_l2_44 = sqrt(rnorm2 / bnorm2)
330 
331  if (present(commtime)) commtime = commtime + tcomm
332 
333  deallocate(r)
334  end function hecmw_rel_resid_l2_44
335 
336  !C
337  !C***
338  !C*** hecmw_Tvec_44
339  !C***
340  !C
341  subroutine hecmw_tvec_44 (hecMESH, X, Y, COMMtime)
343  use m_hecmw_comm_f
344  implicit none
345  type (hecmwst_local_mesh), intent(in) :: hecmesh
346  real(kind=kreal), intent(in) :: x(:)
347  real(kind=kreal), intent(out) :: y(:)
348  real(kind=kreal), intent(inout) :: commtime
349 
350  real(kind=kreal) :: start_time, end_time
351  integer(kind=kint) :: i, j, jj, k, kk
352 
353  start_time= hecmw_wtime()
354  call hecmw_update_r (hecmesh, x, hecmesh%n_node, 4)
355  end_time= hecmw_wtime()
356  commtime = commtime + end_time - start_time
357 
358  !$omp parallel default(none),private(i,k,kk,j,jj),shared(hecMESH,X,Y)
359  !$omp do
360  do i= 1, hecmesh%nn_internal * hecmesh%n_dof
361  y(i)= x(i)
362  enddo
363  !$omp end do
364 
365  !$omp do
366  outer: do i= 1, hecmesh%mpc%n_mpc
367  do j= hecmesh%mpc%mpc_index(i-1) + 1, hecmesh%mpc%mpc_index(i)
368  if (hecmesh%mpc%mpc_dof(j) > 4) cycle outer
369  enddo
370  k = hecmesh%mpc%mpc_index(i-1) + 1
371  kk = 4 * (hecmesh%mpc%mpc_item(k) - 1) + hecmesh%mpc%mpc_dof(k)
372  y(kk) = 0.d0
373  do j= hecmesh%mpc%mpc_index(i-1) + 2, hecmesh%mpc%mpc_index(i)
374  jj = 4 * (hecmesh%mpc%mpc_item(j) - 1) + hecmesh%mpc%mpc_dof(j)
375  y(kk) = y(kk) - hecmesh%mpc%mpc_val(j) * x(jj)
376  enddo
377  enddo outer
378  !$omp end do
379  !$omp end parallel
380 
381  end subroutine hecmw_tvec_44
382 
383  !C
384  !C***
385  !C*** hecmw_Ttvec_44
386  !C***
387  !C
388  subroutine hecmw_ttvec_44 (hecMESH, X, Y, COMMtime)
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  !$omp parallel default(none),private(i,k,kk,j,jj),shared(hecMESH,X,Y)
406  !$omp do
407  do i= 1, hecmesh%nn_internal * hecmesh%n_dof
408  y(i)= x(i)
409  enddo
410  !$omp end do
411 
412  !$omp do
413  outer: do i= 1, hecmesh%mpc%n_mpc
414  do j= hecmesh%mpc%mpc_index(i-1) + 1, hecmesh%mpc%mpc_index(i)
415  if (hecmesh%mpc%mpc_dof(j) > 4) cycle outer
416  enddo
417  k = hecmesh%mpc%mpc_index(i-1) + 1
418  kk = 4 * (hecmesh%mpc%mpc_item(k) - 1) + hecmesh%mpc%mpc_dof(k)
419  y(kk) = 0.d0
420  do j= hecmesh%mpc%mpc_index(i-1) + 2, hecmesh%mpc%mpc_index(i)
421  jj = 4 * (hecmesh%mpc%mpc_item(j) - 1) + hecmesh%mpc%mpc_dof(j)
422  !omp atomic
423  y(jj) = y(jj) - hecmesh%mpc%mpc_val(j) * x(kk)
424  enddo
425  enddo outer
426  !$omp end do
427  !$omp end parallel
428 
429  end subroutine hecmw_ttvec_44
430 
431  !C
432  !C***
433  !C*** hecmw_TtmatTvec_44
434  !C***
435  !C
436  subroutine hecmw_ttmattvec_44 (hecMESH, hecMAT, X, Y, W, time_Ax, COMMtime)
438  implicit none
439  type (hecmwst_local_mesh), intent(in) :: hecmesh
440  type (hecmwst_matrix), intent(in) :: hecmat
441  real(kind=kreal), intent(in) :: x(:)
442  real(kind=kreal), intent(out) :: y(:), w(:)
443  real(kind=kreal), intent(inout) :: time_ax
444  real(kind=kreal), intent(inout) :: commtime
445 
446  call hecmw_tvec_44(hecmesh, x, y, commtime)
447  call hecmw_matvec_44_inner(hecmesh, hecmat, y, w, time_ax, commtime)
448  call hecmw_ttvec_44(hecmesh, w, y, commtime)
449 
450  end subroutine hecmw_ttmattvec_44
451 
452 
453  !C
454  !C***
455  !C*** hecmw_mat_diag_sr_44
456  !C***
457  !C
458  subroutine hecmw_mat_diag_sr_44(hecMESH, hecMAT, COMMtime)
460  use m_hecmw_comm_f
461  implicit none
462  type (hecmwst_local_mesh), intent(in) :: hecmesh
463  type (hecmwst_matrix), intent(inout), target :: hecmat
464  real(kind=kreal), intent(inout), optional :: commtime
465  real(kind=kreal), allocatable :: w(:,:)
466  real(kind=kreal), pointer :: d(:)
467  integer(kind=kint) :: ip
468  real(kind=kreal) :: start_time, end_time
469  allocate(w(4*hecmat%NP,4))
470  d => hecmat%D
471  do ip= 1, hecmat%N
472  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)
473  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)
474  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)
475  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 )
476  enddo
477  start_time= hecmw_wtime()
478  call hecmw_update_r (hecmesh, w(:,1), hecmat%NP, 4)
479  call hecmw_update_r (hecmesh, w(:,2), hecmat%NP, 4)
480  call hecmw_update_r (hecmesh, w(:,3), hecmat%NP, 4)
481  call hecmw_update_r (hecmesh, w(:,4), hecmat%NP, 4)
482  end_time= hecmw_wtime()
483  if (present(commtime)) commtime = commtime + end_time - start_time
484  do ip= hecmat%N+1, hecmat%NP
485  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)
486  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)
487  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)
488  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)
489  enddo
490  deallocate(w)
491  end subroutine hecmw_mat_diag_sr_44
492 
493 end module hecmw_solver_las_44
hecmw_solver_las_44::hecmw_ttvec_44
subroutine, public hecmw_ttvec_44(hecMESH, X, Y, COMMtime)
Definition: hecmw_solver_las_44.f90:389
hecmw_solver_las_44::hecmw_matvec_44_unset_async
subroutine, public hecmw_matvec_44_unset_async
Definition: hecmw_solver_las_44.f90:79
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_44::hecmw_matvec_44_set_async
subroutine, public hecmw_matvec_44_set_async(hecMAT)
Definition: hecmw_solver_las_44.f90:67
hecmw_solver_las_44
Definition: hecmw_solver_las_44.f90:6
hecmw_jad_type::hecmw_jad_is_initialized
integer(kind=kint) function, public hecmw_jad_is_initialized()
Definition: hecmw_jadm.f90:56
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_44::hecmw_rel_resid_l2_44
real(kind=kreal) function, public hecmw_rel_resid_l2_44(hecMESH, hecMAT, time_Ax, COMMtime)
Definition: hecmw_solver_las_44.f90:306
m_hecmw_comm_f
Definition: hecmw_comm_f.F90:6
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
hecmw_solver_las_44::hecmw_mat_diag_sr_44
subroutine, public hecmw_mat_diag_sr_44(hecMESH, hecMAT, COMMtime)
Definition: hecmw_solver_las_44.f90:459
hecmw_solver_las_44::hecmw_matresid_44
subroutine, public hecmw_matresid_44(hecMESH, hecMAT, X, B, R, time_Ax, COMMtime)
Definition: hecmw_solver_las_44.f90:276
hecmw_solver_las_44::hecmw_matvec_44
subroutine, public hecmw_matvec_44(hecMESH, hecMAT, X, Y, time_Ax, COMMtime)
Definition: hecmw_solver_las_44.f90:35
hecmw_solver_las_44::hecmw_ttmattvec_44
subroutine, public hecmw_ttmattvec_44(hecMESH, hecMAT, X, Y, W, time_Ax, COMMtime)
Definition: hecmw_solver_las_44.f90:437
hecmw_matrix_misc
Definition: hecmw_matrix_misc.f90:6
hecmw_solver_misc
Definition: hecmw_solver_misc.f90:6
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_solver_las_44::hecmw_tvec_44
subroutine, public hecmw_tvec_44(hecMESH, X, Y, COMMtime)
Definition: hecmw_solver_las_44.f90:342
hecmw_util::hecmwst_matrix
Definition: hecmw_util_f.F90:444