32     type (hecmwST_matrix), 
target :: hecMAT
 
   33     type (hecmwST_local_mesh) :: hecMESH
 
   35     integer(kind=kint) :: error
 
   36     integer(kind=kint) :: ITER, METHOD, PRECOND, NSET, METHOD2
 
   37     integer(kind=kint) :: iterPREmax
 
   38     integer(kind=kint) :: ITERlog, TIMElog
 
   39     real(kind=
kreal)   :: resid, sigma_diag, thresh, filter, resid2
 
   40     real(kind=
kreal)   :: time_setup, time_comm, time_sol, tr
 
   41     real(kind=
kreal)   :: time_ax, time_precond
 
   43     integer(kind=kint) :: NREST
 
   44     real(kind=
kreal)   :: sigma
 
   46     integer(kind=kint) :: totalmpc, MPC_METHOD
 
   47     integer(kind=kint) :: auto_sigma_diag
 
   67     if (sigma_diag.lt.0.d0) 
then 
   78     totalmpc = hecmesh%mpc%n_mpc
 
   81     if (totalmpc > 0 .and. mpc_method == 2) 
then 
  106           hecmat%symmetric = .true.
 
  107           call hecmw_solve_cg( hecmesh, hecmat, iter, resid, error, time_setup, time_sol, time_comm )
 
  109           hecmat%symmetric = .false.
 
  112           hecmat%symmetric = .false.
 
  113           call hecmw_solve_gmres( hecmesh,hecmat, iter, resid, error, time_setup, time_sol, time_comm )
 
  115           hecmat%symmetric = .false.
 
  116           call hecmw_solve_gpbicg( hecmesh,hecmat, iter, resid, error, time_setup, time_sol, time_comm )
 
  118           hecmat%symmetric = .false.
 
  119           call hecmw_solve_gmresr( hecmesh,hecmat, iter, resid, error, time_setup, time_sol, time_comm )
 
  121           hecmat%symmetric = .false.
 
  122           call hecmw_solve_gmresren( hecmesh,hecmat, iter, resid, error, time_setup, time_sol, time_comm )
 
  131         if ((precond>=10 .and. precond<20) .and. auto_sigma_diag==1 .and. sigma_diag<2.d0) 
then 
  132           sigma_diag = sigma_diag + 0.1
 
  133           if (hecmesh%my_rank.eq.0) 
write(*,*) 
'Increasing SIGMA_DIAG to', sigma_diag
 
  135         elseif (method==1 .and. method2>1) 
then 
  136           if (auto_sigma_diag.eq.1) sigma_diag = 1.0
 
  151     if (hecmesh%my_rank.eq.0 .and. (iterlog.eq.1 .or. timelog.ge.1)) 
then 
  152       write(*,
"(a,1pe12.5)")
'### Relative residual =', resid2
 
  160     if (totalmpc > 0 .and. mpc_method == 2) 
then 
  167     if (hecmesh%my_rank.eq.0 .and. timelog.ge.1) 
then 
  168       tr= (time_sol-time_comm)/(time_sol+1.d-24)*100.d0
 
  169       write (*,
'(/a)')          
'### summary of linear solver' 
  170       write (*,
'(i10,a, 1pe16.6)')      iter, 
' iterations  ', resid
 
  171       write (*,
'(a, 1pe16.6 )') 
'    set-up time      : ', time_setup
 
  172       write (*,
'(a, 1pe16.6 )') 
'    solver time      : ', time_sol
 
  173       write (*,
'(a, 1pe16.6 )') 
'    solver/comm time : ', time_comm
 
  174       write (*,
'(a, 1pe16.6 )') 
'    solver/matvec    : ', time_ax
 
  175       write (*,
'(a, 1pe16.6 )') 
'    solver/precond   : ', time_precond
 
  177         write (*,
'(a, 1pe16.6 )') 
'    solver/1 iter    : ', time_sol / iter
 
  178       write (*,
'(a, 1pe16.6/)') 
'    work ratio (%)   : ', tr
 
  191     type (hecmwST_local_mesh) :: hecMESH
 
  193     integer(kint) :: N, NP, NDOF, NPU, NPL, NZ
 
  194     integer(kint) :: base, i, count_Ax
 
  195     real(kreal) :: time_Ax, size_matrix, flop_matrix, size_vector, memory_size, tmp, num
 
  196     real(kreal) :: t_max, t_min, t_avg, t_sd
 
  197     character(2) :: SI(0:6) = [
'  ',
' K',
' M',
' G',
' T',
' P',
' E']
 
  204     npu = hecmat%indexU(n)
 
  205     npl = hecmat%indexL(n)
 
  208     size_matrix = kreal*nz*ndof**2 &
 
  211     size_vector = kreal*n*ndof &
 
  213     memory_size = size_matrix + size_vector
 
  214     flop_matrix = 2.0d0*nz*ndof**2
 
  221     i = int(log(num) / log(dble(base)))
 
  222     tmp = 1.0d0/base**i * num
 
  223     if(hecmesh%my_rank == 0)
then 
  224       write (*,
"(a,f11.3,a,a)") 
"memory amount of coef. matrix: ", tmp, si(i),
"B" 
  227     num = count_ax*memory_size/time_ax
 
  228     i = int(log(num) / log(dble(base)))
 
  229     tmp = 1.0d0/base**i * num
 
  231     if(hecmesh%my_rank == 0)
then 
  232       write(*,
"(a,f11.3,a,a)") 
"matvec memory band width     : ", tmp, si(i),
"B/s" 
  233       write(*,
"(a,f11.3)") 
'  Max     :',t_max
 
  234       write(*,
"(a,f11.3)") 
'  Min     :',t_min
 
  235       write(*,
"(a,f11.3)") 
'  Avg     :',t_avg
 
  236       write(*,
"(a,f11.3)") 
'  Std Dev :',t_sd
 
  239     num = count_ax*flop_matrix/time_ax
 
  240     i = int(log(num) / log(dble(base)))
 
  241     tmp = 1.0d0/base**i * num
 
  243     if(hecmesh%my_rank == 0)
then 
  244       write(*,
"(a,f11.3,a,a)") 
"matvec FLOPs                 : ", tmp, si(i),
"FLOPs" 
  245       write(*,
"(a,f11.3)") 
'  Max     :',t_max
 
  246       write(*,
"(a,f11.3)") 
'  Min     :',t_min
 
  247       write(*,
"(a,f11.3)") 
'  Avg     :',t_avg
 
  248       write(*,
"(a,f11.3)") 
'  Std Dev :',t_sd
 
  258     type (hecmwST_local_mesh) :: hecMESH
 
  259     type (hecmwST_matrix), 
target :: hecMAT
 
  260     integer (kind=kint)::PRECOND,iterPREmax,i,j,error
 
  267       do j = 1, hecmat%NDOF
 
  268         if (dabs(hecmat%D(hecmat%NDOF*hecmat%NDOF*(i-1)+(j-1)*(hecmat%NDOF+1)+1)).eq.0.d0) 
then 
  275     if (error.ne.0 .and. (precond.lt.10 .and. iterpremax.gt.0)) 
then 
  289     real(kind=
kreal),    
dimension(1) :: rhs
 
  290     integer (kind=kint)::precond,iterpremax,i,j,error
 
  302       do j = 1, hecmat%NDOF
 
  303         rhs(1)=rhs(1) + hecmat%B(hecmat%NDOF*(i-1)+j)**2
 
  306     if (hecmesh%mpc%n_mpc > 0) 
then 
  307       do i= 1, hecmesh%mpc%n_mpc
 
  308         rhs(1)= rhs(1) + hecmesh%mpc%mpc_const(i)**2
 
  313     if (rhs(1).eq.0.d0) 
then 
  328     type (hecmwST_local_mesh) :: hecMESH
 
  329     type (hecmwST_matrix), 
target :: hecMAT
 
  330     integer(kind=kint) :: METHOD
 
  331     integer(kind=kint) :: ITER, PRECOND, NSET, iterPREmax, NREST,NBFGS
 
  332     integer(kind=kint) :: ITERlog, TIMElog
 
  334     character(len=30) :: msg_precond
 
  335     character(len=30) :: msg_method
 
  351         msg_method=
"BiCGSTAB" 
  360            msg_method=
"SUP-GMRESR" 
  363         msg_method=
"GMRESR-EN" 
  365         msg_method=
"Unlabeled" 
  375         msg_precond=
"DirectMUMPS" 
  377         write(msg_precond,
"(a,i0,a)") 
"ILU(",precond-10,
")" 
  383         msg_precond=
"Unlabeled" 
  385     if (hecmesh%my_rank.eq.0 .and. (iterlog.eq.1 .or. timelog.ge.1)) 
then 
  386       write (*,
'(a,i0,a,i0,a,a,a,a,a,i0)') 
'### ',hecmat%NDOF,
'x',hecmat%NDOF,
' BLOCK ', &
 
  387         &   trim(msg_method),
", ",trim(msg_precond),
", ", iterpremax
 
subroutine, public hecmw_mat_dump(hecMAT, hecMESH)
subroutine, public hecmw_mat_dump_solution(hecMAT)
real(kind=kreal) function, public hecmw_mat_get_sigma_diag(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_iterpremax(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_nrest(hecMAT)
subroutine, public hecmw_mat_set_flag_diverged(hecMAT, flag_diverged)
real(kind=kreal) function, public hecmw_mat_get_resid(hecMAT)
subroutine, public hecmw_mat_set_flag_converged(hecMAT, flag_converged)
integer(kind=kint) function, public hecmw_mat_get_nbfgs(hecMAT)
subroutine, public hecmw_mat_recycle_precond_setting(hecMAT)
subroutine, public hecmw_mat_set_sigma_diag(hecMAT, sigma_diag)
integer(kind=kint) function, public hecmw_mat_get_iterlog(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_timelog(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_method2(hecMAT)
real(kind=kreal) function, public hecmw_mat_get_thresh(hecMAT)
subroutine, public hecmw_mat_set_flag_mpcmatvec(hecMAT, flag_mpcmatvec)
integer(kind=kint) function, public hecmw_mat_get_method(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_precond(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_nset(hecMAT)
real(kind=kreal) function, public hecmw_mat_get_filter(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_mpc_method(hecMAT)
real(kind=kreal) function, public hecmw_mat_get_sigma(hecMAT)
integer(kind=kint) function, public hecmw_mat_get_iter(hecMAT)
subroutine, public hecmw_precond_clear_timer
real(kind=kreal) function, public hecmw_precond_get_timer()
subroutine hecmw_solve_bicgstab(hecMESH, hecMAT, ITER, RESID, error, Tset, Tsol, Tcomm)
subroutine, public hecmw_solve_cg(hecMESH, hecMAT, ITER, RESID, error, Tset, Tsol, Tcomm)
subroutine, public hecmw_solve_gmres(hecMESH, hecMAT, ITER, RESID, error, Tset, Tsol, Tcomm)
subroutine, public hecmw_solve_gmresr(hecMESH, hecMAT, ITER, RESID, error, Tset, Tsol, Tcomm)
subroutine, public hecmw_solve_gmresren(hecMESH, hecMAT, ITER, RESID, error, Tset, Tsol, Tcomm)
subroutine, public hecmw_solve_gpbicg(hecMESH, hecMAT, ITER, RESID, error, Tset, Tsol, Tcomm)
subroutine hecmw_solve_iterative_printmsg(hecMESH, hecMAT, METHOD)
subroutine hecmw_output_flops(hecMESH, hecMAT, count_Ax, time_Ax)
logical function hecmw_solve_check_zerorhs(hecMESH, hecMAT)
subroutine hecmw_solve_check_zerodiag(hecMESH, hecMAT)
subroutine hecmw_solve_iterative(hecMESH, hecMAT)
subroutine, public hecmw_matvec_clear_timer
subroutine, public hecmw_matvec_set_async(hecMAT)
real(kind=kreal) function, public hecmw_rel_resid_l2(hecMESH, hecMAT, COMMtime)
real(kind=kreal) function, public hecmw_matvec_get_timer()
subroutine, public hecmw_matvec_unset_async
subroutine hecmw_time_statistics(hecMESH, time, t_max, t_min, t_avg, t_sd)
integer(kind=kint), parameter hecmw_sum
integer(kind=4), parameter kint
integer(kind=kint), parameter hecmw_max
integer(kind=4), parameter kreal
subroutine hecmw_allreduce_i1(hecMESH, s, ntag)
subroutine hecmw_allreduce_r(hecMESH, val, n, ntag)
subroutine hecmw_allreduce_r1(hecMESH, s, ntag)
integer(kind=kint), parameter hecmw_solver_error_diverge_pc
integer(kind=kint), parameter hecmw_solver_error_diverge_nan
integer(kind=kint), parameter hecmw_solver_error_incons_pc
integer(kind=kint), parameter hecmw_solver_error_zero_diag
integer(kind=kint), parameter hecmw_solver_error_zero_rhs
subroutine hecmw_solve_error(hecMESH, IFLAG)
integer(kind=kint), parameter hecmw_solver_error_diverge_mat