Add detailed measurements.

non-diag
sfilippone 10 months ago
parent ae7fad95d4
commit e9d1238b43

@ -83,6 +83,9 @@ subroutine psb_dspmv_vect(alpha,a,x,beta,y,desc_a,info,&
character(len=20) :: name, ch_err character(len=20) :: name, ch_err
logical :: aliw, doswap_ logical :: aliw, doswap_
integer(psb_ipk_) :: debug_level, debug_unit integer(psb_ipk_) :: debug_level, debug_unit
logical, parameter :: do_timings=.true.
integer(psb_ipk_), save :: mv_phase1=-1, mv_phase2=-1, mv_phase3=-1, mv_phase4=-1
integer(psb_ipk_), save :: mv_phase11=-1, mv_phase12=-1
name='psb_dspmv' name='psb_dspmv'
info=psb_success_ info=psb_success_
@ -130,6 +133,19 @@ subroutine psb_dspmv_vect(alpha,a,x,beta,y,desc_a,info,&
call psb_errpush(info,name) call psb_errpush(info,name)
goto 9999 goto 9999
end if end if
if ((do_timings).and.(mv_phase1==-1)) &
& mv_phase1 = psb_get_timer_idx("SPMM: and send ")
if ((do_timings).and.(mv_phase2==-1)) &
& mv_phase2 = psb_get_timer_idx("SPMM: and cmp ad")
if ((do_timings).and.(mv_phase3==-1)) &
& mv_phase3 = psb_get_timer_idx("SPMM: and rcv")
if ((do_timings).and.(mv_phase4==-1)) &
& mv_phase4 = psb_get_timer_idx("SPMM: and cmp and")
if ((do_timings).and.(mv_phase11==-1)) &
& mv_phase11 = psb_get_timer_idx("SPMM: noand exch ")
if ((do_timings).and.(mv_phase12==-1)) &
& mv_phase12 = psb_get_timer_idx("SPMM: noand cmp")
m = desc_a%get_global_rows() m = desc_a%get_global_rows()
n = desc_a%get_global_cols() n = desc_a%get_global_cols()
@ -184,18 +200,22 @@ subroutine psb_dspmv_vect(alpha,a,x,beta,y,desc_a,info,&
logical, parameter :: do_timings=.true. logical, parameter :: do_timings=.true.
real(psb_dpk_) :: t1, t2, t3, t4, t5 real(psb_dpk_) :: t1, t2, t3, t4, t5
if (do_timings) call psb_barrier(ctxt) if (do_timings) call psb_barrier(ctxt)
if (do_timings) t1= psb_wtime() if (do_timings) call psb_tic(mv_phase1)
if (doswap_) call psi_swapdata(psb_swap_send_,& if (doswap_) call psi_swapdata(psb_swap_send_,&
& dzero,x%v,desc_a,iwork,info,data=psb_comm_halo_) & dzero,x%v,desc_a,iwork,info,data=psb_comm_halo_)
if (do_timings) t2= psb_wtime() if (do_timings) call psb_toc(mv_phase1)
if (do_timings) call psb_tic(mv_phase2)
call a%ad%spmm(alpha,x%v,beta,y%v,info) call a%ad%spmm(alpha,x%v,beta,y%v,info)
if (do_timings) t3= psb_wtime() if (do_timings) call psb_toc(mv_phase2)
if (do_timings) call psb_tic(mv_phase3)
if (doswap_) call psi_swapdata(psb_swap_recv_,& if (doswap_) call psi_swapdata(psb_swap_recv_,&
& dzero,x%v,desc_a,iwork,info,data=psb_comm_halo_) & dzero,x%v,desc_a,iwork,info,data=psb_comm_halo_)
if (do_timings) call psb_toc(mv_phase3)
if (do_timings) call psb_tic(mv_phase4)
if (do_timings) t4= psb_wtime() if (do_timings) t4= psb_wtime()
call a%and%spmm(alpha,x%v,done,y%v,info) call a%and%spmm(alpha,x%v,done,y%v,info)
if (do_timings) t5= psb_wtime() if (do_timings) call psb_toc(mv_phase4)
if (do_timings) write(0,*) me,' SPMM:',t2-t1,t3-t2,t4-t3,t5-t4
end block end block
else else
@ -203,15 +223,16 @@ subroutine psb_dspmv_vect(alpha,a,x,beta,y,desc_a,info,&
logical, parameter :: do_timings=.true. logical, parameter :: do_timings=.true.
real(psb_dpk_) :: t1, t2, t3, t4, t5 real(psb_dpk_) :: t1, t2, t3, t4, t5
if (do_timings) call psb_barrier(ctxt) if (do_timings) call psb_barrier(ctxt)
if (do_timings) t1= psb_wtime()
if (do_timings) call psb_tic(mv_phase11)
if (doswap_) then if (doswap_) then
call psi_swapdata(ior(psb_swap_send_,psb_swap_recv_),& call psi_swapdata(ior(psb_swap_send_,psb_swap_recv_),&
& dzero,x%v,desc_a,iwork,info,data=psb_comm_halo_) & dzero,x%v,desc_a,iwork,info,data=psb_comm_halo_)
end if end if
if (do_timings) t2= psb_wtime() if (do_timings) call psb_toc(mv_phase11)
if (do_timings) call psb_tic(mv_phase12)
call psb_csmm(alpha,a,x,beta,y,info) call psb_csmm(alpha,a,x,beta,y,info)
if (do_timings) t3= psb_wtime() if (do_timings) call psb_toc(mv_phase12)
if (do_timings) write(0,*) me,' SPMM:',t2-t1,t3-t2
end block end block
end if end if

@ -142,7 +142,7 @@ contains
! the rhs. ! the rhs.
! !
subroutine psb_d_gen_pde3d(ctxt,idim,a,bv,xv,desc_a,afmt,info,& subroutine psb_d_gen_pde3d(ctxt,idim,a,bv,xv,desc_a,afmt,info,&
& f,amold,vmold,imold,partition,nrl,iv) & f,amold,vmold,imold,partition,nrl,iv,tnd)
use psb_base_mod use psb_base_mod
use psb_util_mod use psb_util_mod
! !
@ -173,7 +173,7 @@ contains
class(psb_d_base_vect_type), optional :: vmold class(psb_d_base_vect_type), optional :: vmold
class(psb_i_base_vect_type), optional :: imold class(psb_i_base_vect_type), optional :: imold
integer(psb_ipk_), optional :: partition, nrl,iv(:) integer(psb_ipk_), optional :: partition, nrl,iv(:)
logical, optional :: tnd
! Local variables. ! Local variables.
integer(psb_ipk_), parameter :: nb=20 integer(psb_ipk_), parameter :: nb=20
@ -202,6 +202,7 @@ contains
real(psb_dpk_) :: t0, t1, t2, t3, tasb, talc, ttot, tgen, tcdasb real(psb_dpk_) :: t0, t1, t2, t3, tasb, talc, ttot, tgen, tcdasb
integer(psb_ipk_) :: err_act integer(psb_ipk_) :: err_act
procedure(d_func_3d), pointer :: f_ procedure(d_func_3d), pointer :: f_
logical :: tnd_
character(len=20) :: name, ch_err,tmpfmt character(len=20) :: name, ch_err,tmpfmt
info = psb_success_ info = psb_success_
@ -495,9 +496,9 @@ contains
t1 = psb_wtime() t1 = psb_wtime()
if (info == psb_success_) then if (info == psb_success_) then
if (present(amold)) then if (present(amold)) then
call psb_spasb(a,desc_a,info,mold=amold) call psb_spasb(a,desc_a,info,mold=amold,bld_and=tnd)
else else
call psb_spasb(a,desc_a,info,afmt=afmt) call psb_spasb(a,desc_a,info,afmt=afmt,bld_and=tnd)
end if end if
end if end if
call psb_barrier(ctxt) call psb_barrier(ctxt)
@ -549,13 +550,14 @@ program pdgenspmv
use psb_base_mod use psb_base_mod
use psb_util_mod use psb_util_mod
use psb_d_pde3d_mod use psb_d_pde3d_mod
implicit none implicit none
! input parameters ! input parameters
character(len=20) :: kmethd, ptype character(len=20) :: kmethd, ptype
character(len=5) :: afmt character(len=5) :: afmt
integer(psb_ipk_) :: idim integer(psb_ipk_) :: idim
logical :: tnd
! miscellaneous ! miscellaneous
real(psb_dpk_), parameter :: one = done real(psb_dpk_), parameter :: one = done
real(psb_dpk_) :: t1, t2, tprec, flops, tflops, tt1, tt2, bdwdth real(psb_dpk_) :: t1, t2, tprec, flops, tflops, tt1, tt2, bdwdth
@ -606,14 +608,14 @@ program pdgenspmv
! !
! get parameters ! get parameters
! !
call get_parms(ctxt,afmt,idim) call get_parms(ctxt,afmt,idim,tnd)
call psb_init_timers()
! !
! allocate and fill in the coefficient matrix, rhs and initial guess ! allocate and fill in the coefficient matrix, rhs and initial guess
! !
call psb_barrier(ctxt) call psb_barrier(ctxt)
t1 = psb_wtime() t1 = psb_wtime()
call psb_gen_pde3d(ctxt,idim,a,bv,xv,desc_a,afmt,info) call psb_gen_pde3d(ctxt,idim,a,bv,xv,desc_a,afmt,info,tnd=tnd)
call psb_barrier(ctxt) call psb_barrier(ctxt)
t2 = psb_wtime() - t1 t2 = psb_wtime() - t1
if(info /= psb_success_) then if(info /= psb_success_) then
@ -694,7 +696,7 @@ program pdgenspmv
write(psb_out_unit,'("Total memory occupation for DESC_A: ",i12)')descsize write(psb_out_unit,'("Total memory occupation for DESC_A: ",i12)')descsize
end if end if
call psb_print_timers(ctxt)
! !
! cleanup storage and exit ! cleanup storage and exit
@ -721,10 +723,11 @@ contains
! !
! get iteration parameters from standard input ! get iteration parameters from standard input
! !
subroutine get_parms(ctxt,afmt,idim) subroutine get_parms(ctxt,afmt,idim,tnd)
type(psb_ctxt_type) :: ctxt type(psb_ctxt_type) :: ctxt
character(len=*) :: afmt character(len=*) :: afmt
integer(psb_ipk_) :: idim integer(psb_ipk_) :: idim
logical :: tnd
integer(psb_ipk_) :: np, iam integer(psb_ipk_) :: np, iam
integer(psb_ipk_) :: intbuf(10), ip integer(psb_ipk_) :: intbuf(10), ip
@ -733,9 +736,11 @@ contains
if (iam == 0) then if (iam == 0) then
read(psb_inp_unit,*) afmt read(psb_inp_unit,*) afmt
read(psb_inp_unit,*) idim read(psb_inp_unit,*) idim
read(psb_inp_unit,*) tnd
endif endif
call psb_bcast(ctxt,afmt) call psb_bcast(ctxt,afmt)
call psb_bcast(ctxt,idim) call psb_bcast(ctxt,idim)
call psb_bcast(ctxt,tnd)
if (iam == 0) then if (iam == 0) then
write(psb_out_unit,'("Testing matrix : ell1")') write(psb_out_unit,'("Testing matrix : ell1")')
@ -743,6 +748,8 @@ contains
write(psb_out_unit,'("Number of processors : ",i0)')np write(psb_out_unit,'("Number of processors : ",i0)')np
write(psb_out_unit,'("Data distribution : BLOCK")') write(psb_out_unit,'("Data distribution : BLOCK")')
write(psb_out_unit,'(" ")') write(psb_out_unit,'(" ")')
write(psb_out_unit,'("Storage format ",a)') afmt
write(psb_out_unit,'("Testing overlap ND ",l8)') tnd
end if end if
return return

@ -868,7 +868,7 @@ program psb_d_pde3d
call psb_errpush(info,name,a_err=ch_err) call psb_errpush(info,name,a_err=ch_err)
goto 9999 goto 9999
end if end if
call psb_print_timers(ctxt)
call psb_exit(ctxt) call psb_exit(ctxt)
stop stop

@ -5,7 +5,7 @@ CSR Storage format for matrix A: CSR COO
200 Domain size (acutal system is this**3 (pde3d) or **2 (pde2d) ) 200 Domain size (acutal system is this**3 (pde3d) or **2 (pde2d) )
3 Partition: 1 BLOCK 3 3D 3 Partition: 1 BLOCK 3 3D
2 Stopping criterion 1 2 2 Stopping criterion 1 2
0008 MAXIT 0200 MAXIT
10 ITRACE 10 ITRACE
002 IRST restart for RGMRES and BiCGSTABL 002 IRST restart for RGMRES and BiCGSTABL
INVK Block Solver ILU,ILUT,INVK,AINVT,AORTH INVK Block Solver ILU,ILUT,INVK,AINVT,AORTH

Loading…
Cancel
Save