You cannot select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.
psblas3/test/unitTest/test_psb_swapdata.pf

1632 lines
32 KiB
Plaintext

module test_psb_swapdata
use pfunit_mod
use psb_base_mod
implicit none
include 'mpif.h'
interface
subroutine psi_sswap_xchg_v(iictxt,iicomm,flag,beta,y,xchg,info)
use psi_mod, psb_protect_name => psi_sswap_xchg_v
use psb_error_mod
use psb_realloc_mod
use psb_desc_mod
use psb_penv_mod
use psb_s_base_vect_mod
use iso_fortran_env
implicit none
integer(psb_ipk_), intent(in) :: iictxt,iicomm,flag
integer(psb_ipk_), intent(out) :: info
real(psb_spk_) :: y(:)
real(psb_spk_) :: beta
class(psb_xch_idx_type), intent(inout) :: xchg
! locals
integer(psb_mpik_) :: ictxt, icomm, np, me,&
& proc_to_comm, p2ptag, iret
integer(psb_ipk_) :: nesd, nerv,&
& err_act, i, idx_pt, totsnd_, totrcv_,p1,p2,isz,rp1,rp2,&
& snd_pt, rcv_pt, pnti, n, ip, img, nxch, myself
integer :: count
real(psb_spk_), allocatable, save :: buffer(:)[:], sndbuf(:)
type(event_type), allocatable, save :: ufg(:)[:]
type(event_type), allocatable, save :: clear[:]
integer, save :: last_clear_count = 0
logical :: swap_mpi, swap_sync, swap_send, swap_recv,&
& albf,do_send,do_recv
integer(psb_ipk_) :: ierr(5)
character(len=20) :: name
end subroutine psi_sswap_xchg_v
subroutine psi_sswap_xchg_m(iictxt,iicomm,flag,m,beta,y,xchg,info)
use psi_mod, psb_protect_name => psi_sswap_xchg_m
use psb_error_mod
use psb_realloc_mod
use psb_desc_mod
use psb_penv_mod
use psb_s_base_vect_mod
use iso_fortran_env
implicit none
integer(psb_ipk_), intent(in) :: iictxt,iicomm,flag, m
integer(psb_ipk_), intent(out) :: info
real(psb_spk_) :: y(:,:)
real(psb_spk_) :: beta
class(psb_xch_idx_type), intent(inout) :: xchg
! locals
integer(psb_mpik_) :: ictxt, icomm, np, me,&
& proc_to_comm, p2ptag, iret
integer(psb_ipk_) :: nesd, nerv,&
& err_act, i, idx_pt, totsnd_, totrcv_,p1,p2,isz,rp1,rp2,&
& snd_pt, rcv_pt, pnti, n, ip, img, nxch, myself
integer :: count
real(psb_spk_), allocatable, save :: buffer(:)[:], sndbuf(:)
type(event_type), allocatable, save :: ufg(:)[:]
type(event_type), allocatable, save :: clear[:]
integer, save :: last_clear_count = 0
logical :: swap_mpi, swap_sync, swap_send, swap_recv,&
& albf,do_send,do_recv
integer(psb_ipk_) :: ierr(5)
character(len=20) :: name
end subroutine psi_sswap_xchg_m
subroutine psi_dswap_xchg_v(iictxt,iicomm,flag,beta,y,xchg,info)
use psi_mod, psb_protect_name => psi_dswap_xchg_v
use psb_error_mod
use psb_realloc_mod
use psb_desc_mod
use psb_penv_mod
use psb_d_base_vect_mod
use iso_fortran_env
implicit none
integer(psb_ipk_), intent(in) :: iictxt,iicomm,flag
integer(psb_ipk_), intent(out) :: info
real(psb_dpk_) :: y(:)
real(psb_dpk_) :: beta
class(psb_xch_idx_type), intent(inout) :: xchg
! locals
integer(psb_mpik_) :: ictxt, icomm, np, me,&
& proc_to_comm, p2ptag, iret
integer(psb_ipk_) :: nesd, nerv,&
& err_act, i, idx_pt, totsnd_, totrcv_,p1,p2,isz,rp1,rp2,&
& snd_pt, rcv_pt, pnti, n, ip, img, nxch, myself
integer :: count
real(psb_dpk_), allocatable, save :: buffer(:)[:], sndbuf(:)
type(event_type), allocatable, save :: ufg(:)[:]
type(event_type), allocatable, save :: clear[:]
integer, save :: last_clear_count = 0
logical :: swap_mpi, swap_sync, swap_send, swap_recv,&
& albf,do_send,do_recv
integer(psb_ipk_) :: ierr(5)
character(len=20) :: name
end subroutine psi_dswap_xchg_v
end interface
interface
subroutine psi_dswap_xchg_m(iictxt,iicomm,flag,m,beta,y,xchg,info)
use psi_mod, psb_protect_name => psi_dswap_xchg_m
use psb_error_mod
use psb_realloc_mod
use psb_desc_mod
use psb_penv_mod
use psb_d_base_vect_mod
use iso_fortran_env
implicit none
integer(psb_ipk_), intent(in) :: iictxt,iicomm,flag,m
integer(psb_ipk_), intent(out) :: info
real(psb_dpk_) :: y(:,:)
real(psb_dpk_) :: beta
class(psb_xch_idx_type), intent(inout) :: xchg
! locals
integer(psb_mpik_) :: ictxt, icomm, np, me,&
& proc_to_comm, p2ptag, iret
integer(psb_ipk_) :: nesd, nerv,&
& err_act, i, idx_pt, totsnd_, totrcv_,p1,p2,isz,rp1,rp2,&
& snd_pt, rcv_pt, pnti, n, ip, img, nxch, myself
integer :: count
real(psb_dpk_), allocatable, save :: buffer(:)[:], sndbuf(:)
type(event_type), allocatable, save :: ufg(:)[:]
type(event_type), allocatable, save :: clear[:]
integer, save :: last_clear_count = 0
logical :: swap_mpi, swap_sync, swap_send, swap_recv,&
& albf,do_send,do_recv
integer(psb_ipk_) :: ierr(5)
character(len=20) :: name
end subroutine psi_dswap_xchg_m
end interface
contains
@test(nimgs=[std])
subroutine test_psb_sswapdatav_2imgs(this)
implicit none
Class(CafTestMethod), intent(inout) :: this
integer :: msg, me, i=0, np, j, info
integer, parameter :: nrows=6
integer :: icontxt, mid, true
integer, allocatable :: vg(:), ia(:)
real(psb_spk_), allocatable :: val(:)
real(psb_spk_), allocatable :: y(:), check(:)
class(psb_xch_idx_type), pointer :: xchg
integer(psb_ipk_) :: iictxt, icomm, flag
type(psb_desc_type):: desc_a
type(psb_s_vect_type) :: x
np = this%getNumImages()
if (np < 2) then
print*,'You need at least 2 processes to run this test.'
return
endif
call psb_init(icontxt,np,MPI_COMM_WORLD)
!call psb_info(icontxt, me, np)
me = this_image()
!Allocate vectors
allocate(vg(nrows))
allocate(ia(nrows))
allocate(val(nrows))
allocate(y(nrows))
i = 0
do j=1,size(vg,1)
vg(j)= i
i = i+1
if (i==np) then
i=0
endif
enddo
!Use only 2 processes
!Assuming nrows is a multiple of 2 so mid is an integer
!Distribute equally to the two processes
mid=nrows/2
do i=1, mid
vg(i)=0
enddo
do i=mid+1, nrows
vg(i)=1
enddo
do i=1,size(ia,1)
ia(i)=i
enddo
do i=1,mid
val(i)=1.
enddo
do i=mid + 1,nrows
val(i)=2.
enddo
call psb_cdall(icontxt,desc_a,info, vg=vg)
if ( me == 1) call psb_cdins(nz=1,ia=(/mid/),ja=(/mid+1/), desc_a=desc_a, info=info)
if ( me == 2) call psb_cdins(nz=1,ia=(/mid+1/),ja=(/mid/), desc_a=desc_a, info=info)
call psb_cdasb(desc_a, info)
call psb_geall(x,desc_a,info)
call psb_geins(m=nrows, irw=ia, val=val,x=x, desc_a=desc_a, info=info)
call psb_geasb(x,desc_a,info)
call psb_barrier(icontxt)
y = x%get_vect()
!Let's modify x, so we need to update halo indices
if ((me == 1).or.(me == 2)) then
y(mid +1)=y(mid+1) + 2.0
endif
call psb_barrier(icontxt)
! END OF SETUP
iictxt = desc_a%get_context()
icomm = desc_a%get_mpic()
call desc_a%get_list(psb_comm_halo_,xchg,info)
flag = IOR(psb_swap_send_, psb_swap_recv_)
call psi_sswap_xchg_v(iictxt,icomm,flag,0.0,y,xchg,info)
!GETTING BACK X
call psb_barrier(icontxt)
!Let's build the expected solution
if (allocated(check)) deallocate(check)
if ((me == 1).or.(me==2)) then
allocate(check(mid+1))
else
allocate(check(1))
endif
if (me == 1 ) then
check(1:mid)=1.0
check(mid + 1)=2.0
else if (me == 2) then
check(1:mid)=2.0
check(mid + 1)=1.0
else
check(1)=0.0
endif
!call psb_barrier(icontxt)
if ((me==1).or.(me==2)) then
true = 1
else
true=0
endif
@assertEqual(real(true*check),real(true*y))
deallocate(vg,ia,val,y,check)
call psb_gefree(x, desc_a, info)
call psb_cdfree(desc_a, info)
call psb_exit(icontxt)
end subroutine test_psb_sswapdatav_2imgs
@test(nimgs=[std])
subroutine test_psb_swapdatam_2imgs(this)
implicit none
Class(CafTestMethod), intent(inout) :: this
integer :: msg, me, i=0, np, j, info
integer, parameter :: nrows=6
integer :: icontxt, mid, true
integer, allocatable :: vg(:), ia(:)
real(psb_dpk_), allocatable :: val(:)
real(psb_dpk_), allocatable :: y(:,:), check(:), v(:)
class(psb_xch_idx_type), pointer :: xchg
integer(psb_ipk_) :: iictxt, icomm, flag
type(psb_desc_type):: desc_a
type(psb_d_vect_type) :: x
np = this%getNumImages()
if (np < 2) then
print*,'You need at least 2 processes to run this test.'
return
endif
call psb_init(icontxt,np,MPI_COMM_WORLD)
!call psb_info(icontxt, me, np)
me = this_image()
!Allocate vectors
allocate(vg(nrows))
allocate(ia(nrows))
allocate(val(nrows))
allocate(v(nrows))
i = 0
do j=1,size(vg,1)
vg(j)= i
i = i+1
if (i==np) then
i=0
endif
enddo
!Use only 2 processes
!Assuming nrows is a multiple of 2 so mid is an integer
!Distribute equally to the two processes
mid=nrows/2
do i=1, mid
vg(i)=0
enddo
do i=mid+1, nrows
vg(i)=1
enddo
do i=1,size(ia,1)
ia(i)=i
enddo
do i=1,mid
val(i)=1.
enddo
do i=mid + 1,nrows
val(i)=2.
enddo
call psb_cdall(icontxt,desc_a,info, vg=vg)
if ( me == 1) call psb_cdins(nz=1,ia=(/mid/),ja=(/mid+1/), desc_a=desc_a, info=info)
if ( me == 2) call psb_cdins(nz=1,ia=(/mid+1/),ja=(/mid/), desc_a=desc_a, info=info)
call psb_cdasb(desc_a, info)
call psb_geall(x,desc_a,info)
call psb_geins(m=nrows, irw=ia, val=val,x=x, desc_a=desc_a, info=info)
call psb_geasb(x,desc_a,info)
call psb_barrier(icontxt)
v = x%get_vect()
allocate(y(size(v,1),1))
y(:,1)=v
!Let's modify x, so we need to update halo indices
if ((me == 1).or.(me == 2)) then
y(mid +1,1)=y(mid+1,1) + 2.0d0
endif
call psb_barrier(icontxt)
! END OF SETUP
iictxt = desc_a%get_context()
icomm = desc_a%get_mpic()
call desc_a%get_list(psb_comm_halo_,xchg,info)
flag = IOR(psb_swap_send_, psb_swap_recv_)
call psi_dswap_xchg_m(iictxt,icomm,flag,1,0.0d0,y,xchg,info)
!GETTING BACK X
call psb_barrier(icontxt)
!Let's build the expected solution
if (allocated(check)) deallocate(check)
if ((me == 1).or.(me==2)) then
allocate(check(mid+1))
else
allocate(check(1))
endif
if (me == 1 ) then
check(1:mid)=1.0d0
check(mid + 1)=2.0d0
else if (me == 2) then
check(1:mid)=2.0d0
check(mid + 1)=1.0d0
else
check(1)=0.0d0
endif
!call psb_barrier(icontxt)
if ((me==1).or.(me==2)) then
true = 1
else
true=0
endif
@assertEqual(real(true*check),real(true*y(:,1)))
deallocate(vg,ia,val,y,v,check)
call psb_gefree(x, desc_a, info)
call psb_cdfree(desc_a, info)
call psb_exit(icontxt)
end subroutine test_psb_swapdatam_2imgs
@test(nimgs=[std])
subroutine test_psb_swapdatav_2imgs(this)
implicit none
Class(CafTestMethod), intent(inout) :: this
integer :: msg, me, i=0, np, j, info
integer, parameter :: nrows=6
integer :: icontxt, mid, true
integer, allocatable :: vg(:), ia(:)
real(psb_dpk_), allocatable :: val(:)
real(psb_dpk_), allocatable :: y(:), check(:)
class(psb_xch_idx_type), pointer :: xchg
integer(psb_ipk_) :: iictxt, icomm, flag
type(psb_desc_type):: desc_a
type(psb_d_vect_type) :: x
np = this%getNumImages()
if (np < 2) then
print*,'You need at least 2 processes to run this test.'
return
endif
call psb_init(icontxt,np,MPI_COMM_WORLD)
!call psb_info(icontxt, me, np)
me = this_image()
!Allocate vectors
allocate(vg(nrows))
allocate(ia(nrows))
allocate(val(nrows))
allocate(y(nrows))
i = 0
do j=1,size(vg,1)
vg(j)= i
i = i+1
if (i==np) then
i=0
endif
enddo
!Use only 2 processes
!Assuming nrows is a multiple of 2 so mid is an integer
!Distribute equally to the two processes
mid=nrows/2
do i=1, mid
vg(i)=0
enddo
do i=mid+1, nrows
vg(i)=1
enddo
do i=1,size(ia,1)
ia(i)=i
enddo
do i=1,mid
val(i)=1.
enddo
do i=mid + 1,nrows
val(i)=2.
enddo
call psb_cdall(icontxt,desc_a,info, vg=vg)
if ( me == 1) call psb_cdins(nz=1,ia=(/mid/),ja=(/mid+1/), desc_a=desc_a, info=info)
if ( me == 2) call psb_cdins(nz=1,ia=(/mid+1/),ja=(/mid/), desc_a=desc_a, info=info)
call psb_cdasb(desc_a, info)
call psb_geall(x,desc_a,info)
call psb_geins(m=nrows, irw=ia, val=val,x=x, desc_a=desc_a, info=info)
call psb_geasb(x,desc_a,info)
call psb_barrier(icontxt)
y = x%get_vect()
!Let's modify x, so we need to update halo indices
if ((me == 1).or.(me == 2)) then
y(mid +1)=y(mid+1) + 2.0d0
endif
call psb_barrier(icontxt)
! END OF SETUP
iictxt = desc_a%get_context()
icomm = desc_a%get_mpic()
call desc_a%get_list(psb_comm_halo_,xchg,info)
flag = IOR(psb_swap_send_, psb_swap_recv_)
call psi_dswap_xchg_v(iictxt,icomm,flag,0.0d0,y,xchg,info)
!GETTING BACK X
call psb_barrier(icontxt)
!Let's build the expected solution
if (allocated(check)) deallocate(check)
if ((me == 1).or.(me==2)) then
allocate(check(mid+1))
else
allocate(check(1))
endif
if (me == 1 ) then
check(1:mid)=1.0d0
check(mid + 1)=2.0d0
else if (me == 2) then
check(1:mid)=2.0d0
check(mid + 1)=1.0d0
else
check(1)=0.0d0
endif
!call psb_barrier(icontxt)
if ((me==1).or.(me==2)) then
true = 1
else
true=0
endif
@assertEqual(real(true*check),real(true*y))
deallocate(vg,ia,val,y,check)
call psb_gefree(x, desc_a, info)
call psb_cdfree(desc_a, info)
call psb_exit(icontxt)
end subroutine test_psb_swapdatav_2imgs
@test(nimgs=[std])
subroutine test_psb_swapdata_2imgs(this)
implicit none
Class(CafTestMethod), intent(inout) :: this
integer :: msg, me, i=0, np, j, info
integer, parameter :: nrows=6
integer :: icontxt, mid, true
integer, allocatable :: vg(:), ia(:)
real(psb_dpk_), allocatable :: val(:)
real(psb_dpk_), allocatable :: v(:), check(:)
class(psb_xch_idx_type), pointer :: xchg
integer(psb_ipk_) :: iictxt, icomm, flag
type(psb_desc_type):: desc_a
type(psb_d_vect_type) :: x
np = this%getNumImages()
if (np < 2) then
print*,'You need at least 2 processes to run this test.'
return
endif
call psb_init(icontxt,np,MPI_COMM_WORLD)
!call psb_info(icontxt, me, np)
me = this_image()
!Allocate vectors
allocate(vg(nrows))
allocate(ia(nrows))
allocate(val(nrows))
allocate(v(nrows))
i = 0
do j=1,size(vg,1)
vg(j)= i
i = i+1
if (i==np) then
i=0
endif
enddo
!Use only 2 processes
!Assuming nrows is a multiple of 2 so mid is an integer
!Distribute equally to the two processes
mid=nrows/2
do i=1, mid
vg(i)=0
enddo
do i=mid+1, nrows
vg(i)=1
enddo
do i=1,size(ia,1)
ia(i)=i
enddo
do i=1,mid
val(i)=1.
enddo
do i=mid + 1,nrows
val(i)=2.
enddo
call psb_cdall(icontxt,desc_a,info, vg=vg)
if ( me == 1) call psb_cdins(nz=1,ia=(/mid/),ja=(/mid+1/), desc_a=desc_a, info=info)
if ( me == 2) call psb_cdins(nz=1,ia=(/mid+1/),ja=(/mid/), desc_a=desc_a, info=info)
call psb_cdasb(desc_a, info)
call psb_geall(x,desc_a,info)
call psb_geins(m=nrows, irw=ia, val=val,x=x, desc_a=desc_a, info=info)
call psb_geasb(x,desc_a,info)
call psb_barrier(icontxt)
v = x%get_vect()
!Let's modify x, so we need to update halo indices
if ((me == 1).or.(me == 2)) then
x%v%v(mid +1)=x%v%v(mid+1) + 2.0d0
endif
call psb_barrier(icontxt)
! END OF SETUP
iictxt = desc_a%get_context()
icomm = desc_a%get_mpic()
call desc_a%get_list(psb_comm_halo_,xchg,info)
flag = IOR(psb_swap_send_, psb_swap_recv_)
call psi_dswap_xchg_vect(iictxt,icomm,flag,0.0d0,x%v,xchg,info)
!GETTING BACK X
call psb_barrier(icontxt)
v = x%get_vect()
!Let's build the expected solution
if (allocated(check)) deallocate(check)
if ((me == 1).or.(me==2)) then
allocate(check(mid+1))
else
allocate(check(1))
endif
if (me == 1 ) then
check(1:mid)=1.0d0
check(mid + 1)=2.0d0
else if (me == 2) then
check(1:mid)=2.0d0
check(mid + 1)=1.0d0
else
check(1)=0.0d0
endif
!call psb_barrier(icontxt)
if ((me==1).or.(me==2)) then
true = 1
else
true=0
endif
@assertEqual(real(true*check),real(true*v))
deallocate(vg,ia,val,v,check)
call psb_gefree(x, desc_a, info)
call psb_cdfree(desc_a, info)
call psb_exit(icontxt)
end subroutine test_psb_swapdata_2imgs
@test(nimgs=[std])
subroutine test_psb_swapdata_4imgs(this)
implicit none
Class(CafTestMethod), intent(inout) :: this
integer :: msg, me, i=0, np, j, info, nz
integer, parameter :: nrows = 8
integer :: icontxt, mid, true
integer, allocatable :: vg(:), ia(:), ja(:), irw(:)
real(psb_dpk_), allocatable :: val(:)
real(psb_dpk_), allocatable :: v(:), check(:)
integer(psb_ipk_) :: iictxt, icomm, flag
class(psb_xch_idx_type), pointer :: xchg
type(psb_desc_type):: desc_a
type(psb_d_vect_type) :: x
np = this%getNumImages()
if (np < 4) then
print*,'You need at least 4 processes to run this test.'
return
endif
call psb_init(icontxt,np,MPI_COMM_WORLD)
me = this_image()
!Allocate vectors
allocate(vg(nrows))
allocate(val(nrows))
allocate(v(nrows))
!Use only 4 processes
!Assuming nrows is a multiple of 4 so mid is an integer
!Distribute equally to the two processes
mid=nrows/4
do i=1, mid
vg(i)=0
enddo
do i=mid+1, 2*mid
vg(i)=1
enddo
do i=2*mid + 1, 3*mid
vg(i)=2
enddo
do i=3*mid+1, nrows
vg(i)=3
enddo
if (me == 1) nz = 5
if (me == 2) nz = 7
if (me == 3) nz = 7
if (me == 4) nz = 5
if (me > 4) nz = 0
allocate(ia(nz),ja(nz))
if (me == 1) then
ia(1)=2
ja(1)=1
ia(2)=1
ja(2)=2
ia(3)=2
ja(3)=3
ia(4)=1
ja(4)=4
ia(5)=2
ja(5)=5
endif
if (me == 2) then
ia(1)=4
ja(1)=1
ia(2)=3
ja(2)=2
ia(3)=4
ja(3)=3
ia(4)=3
ja(4)=4
ia(5)=4
ja(5)=5
ia(6)=3
ja(6)=6
ia(7)=4
ja(7)=6
endif
if (me == 3) then
ia(1)=5
ja(1)=2
ia(2)=6
ja(2)=3
ia(3)=5
ja(3)=4
ia(4)=6
ja(4)=5
ia(5)=5
ja(5)=6
ia(6)=6
ja(6)=7
ia(7)=5
ja(7)=8
endif
if (me == 4) then
ia(1)=7
ja(1)=4
ia(2)=8
ja(2)=5
ia(3)=7
ja(3)=6
ia(4)=8
ja(4)=7
ia(5)=7
ja(5)=8
endif
do i=1,mid
val(i)=1.
enddo
do i= mid + 1, 2*mid
val(i)=2.
enddo
do i=2*mid + 1, 3*mid
val(i)=3.
enddo
do i=3*mid + 1, nrows
val(i)=4.
enddo
call psb_cdall(icontxt,desc_a,info, vg=vg)
call psb_cdins(nz=nz,ia=ia,ja=ja, desc_a=desc_a, info=info)
call psb_cdasb(desc_a, info)
allocate(irw(nrows))
do i=1,nrows
irw(i)=i
enddo
call psb_geall(x,desc_a,info)
call psb_geins(m=nrows, irw=irw, val=val,x=x, desc_a=desc_a, info=info)
call psb_geasb(x,desc_a,info)
if (me==1) nz = 5
if (me==2) nz = 6
if (me==3) nz = 7
if (me==4) nz = 5
if (me > 4) nz = 1
allocate (check(nz))
if (me == 1) then
check(1)=2
check(2)=2
check(3)=8
check(4)=8
check(5)=18
endif
if (me == 2) then
check(1)=8
check(2)=8
check(3)=2
check(4)=2
check(5)=18
check(6)=18
endif
if (me == 3) then
check(1)=18
check(2)=18
check(3)=1
check(4)=8
check(5)=8
check(6)=32
check(7)=32
endif
if (me == 4) then
check(1)=32
check(2)=32
check(3)=8
check(4)=18
check(5)=18
endif
! END OF SETUP
call psb_barrier(icontxt)
!We can do something better here
x%v%v = x%v%v*2*me
call psb_barrier(icontxt)
iictxt = desc_a%get_context()
icomm = desc_a%get_mpic()
call desc_a%get_list(psb_comm_halo_,xchg,info)
flag = IOR(psb_swap_send_, psb_swap_recv_)
call psi_dswap_xchg_vect(iictxt,icomm,flag,0.0d0,x%v,xchg,info)
call psb_barrier(icontxt)
v = x%get_vect()
!call psb_barrier(icontxt)
if ((me==1).or.(me==2)) then
true = 1
else
true=0
endif
@assertEqual(real(true*check),real(true*v))
deallocate(vg,ia,val,v,check)
call psb_gefree(x, desc_a, info)
call psb_cdfree(desc_a, info)
call psb_exit(icontxt)
end subroutine test_psb_swapdata_4imgs
@test(nimgs=[std])
subroutine test_psb_swapdata_8imgs(this)
implicit none
Class(CafTestMethod), intent(inout) :: this
integer :: msg, me, i=0, np, j, info, nz
integer, parameter :: nrows = 8
integer :: icontxt, true
integer, allocatable :: vg(:), ia(:), ja(:), irw(:)
real(psb_dpk_), allocatable :: val(:)
real(psb_dpk_), allocatable :: v(:), check(:)
integer(psb_ipk_) :: iictxt, icomm, flag
class(psb_xch_idx_type), pointer :: xchg
type(psb_desc_type):: desc_a
type(psb_d_vect_type) :: x
np = this%getNumImages()
if (np < 8) then
print*,'You need at least 8 processes to run this test.'
return
endif
call psb_init(icontxt,np,MPI_COMM_WORLD)
me = this_image()
!Allocate vectors
allocate(vg(nrows))
allocate(val(nrows))
allocate(v(nrows))
!Use only 8 processes
!Each process has a line
do i=1, nrows
vg(i)=i-1
enddo
if (me == 1) nz = 1
if (me == 2) nz = 2
if (me == 3) nz = 3
if ((me >= 4).and.(me <= 8)) nz = 2
if (me > 8) nz = 0
allocate(ia(nz),ja(nz))
if (me == 1) then
ia(1)=1
ja(1)=6
endif
if (me == 2) then
ia(1)=2
ja(1)=1
ia(2)=2
ja(2)=7
endif
if (me == 3) then
ia(1)=3
ja(1)=1
ia(2)=3
ja(2)=2
ia(3)=3
ja(3)=8
endif
if (me == 4) then
ia(1)=4
ja(1)=2
ia(2)=4
ja(2)=3
endif
if (me == 5) then
ia(1)=5
ja(1)=3
ia(2)=5
ja(2)=4
endif
if (me == 6) then
ia(1)=6
ja(1)=4
ia(2)=6
ja(2)=5
endif
if (me == 7) then
ia(1)=7
ja(1)=5
ia(2)=7
ja(2)=6
endif
if (me == 8) then
ia(1)=8
ja(1)=6
ia(2)=8
ja(2)=7
endif
do i=1,nrows
val(i)=me
enddo
call psb_cdall(icontxt,desc_a,info, vg=vg)
call psb_cdins(nz=nz,ia=ia,ja=ja, desc_a=desc_a, info=info)
call psb_cdasb(desc_a, info)
allocate(irw(nrows))
do i=1,nrows
irw(i)=i
enddo
call psb_geall(x,desc_a,info)
call psb_geins(m=nrows, irw=irw, val=val,x=x, desc_a=desc_a, info=info)
call psb_geasb(x,desc_a,info)
if (me == 1) nz = 2
if (me == 2) nz = 3
if (me == 3) nz = 4
if ((me >= 4).and.(me <= 8)) nz = 3
if (me > 8) nz = 1
if (allocated(check)) deallocate(check)
allocate (check(nz))
if (me == 1) then
check(1)=2
check(2)=12
endif
if (me == 2) then
check(1)=4
check(2)=2
check(3)=14
endif
if (me == 3) then
check(1)=6
check(2)=2
check(3)=4
check(4)=16
endif
if (me == 4) then
check(1)=8
check(2)=4
check(3)=6
endif
if (me == 5) then
check(1)=10
check(2)=6
check(3)=8
endif
if (me == 6) then
check(1)=12
check(2)=8
check(3)=10
endif
if (me == 7) then
check(1)=14
check(2)=10
check(3)=12
endif
if (me == 8) then
check(1)=16
check(2)=12
check(3)=14
endif
if (me > 8) then
check(1)=0
endif
! END OF SETUP
call psb_barrier(icontxt)
!We can do something better here
x%v%v = x%v%v + me
call psb_barrier(icontxt)
iictxt = desc_a%get_context()
icomm = desc_a%get_mpic()
call desc_a%get_list(psb_comm_halo_,xchg,info)
flag = IOR(psb_swap_send_, psb_swap_recv_)
call psi_dswap_xchg_vect(iictxt,icomm,flag,0.0d0,x%v,xchg,info)
call psb_barrier(icontxt)
v = x%get_vect()
!call psb_barrier(icontxt)
if ((me==1).or.(me==2)) then
true = 1
else
true=0
endif
@assertEqual(real(true*check),real(true*v))
deallocate(vg,ia,val,v,check)
call psb_gefree(x, desc_a, info)
call psb_cdfree(desc_a, info)
call psb_exit(icontxt)
end subroutine test_psb_swapdata_8imgs
@test(nimgs=[std])
subroutine test_psb_swapdata_8imgs_b(this)
implicit none
Class(CafTestMethod), intent(inout) :: this
integer :: msg, me, i=0, np, j, info, nz
integer, parameter :: nrows = 8
integer :: icontxt, true
integer, allocatable :: vg(:), ia(:), ja(:), irw(:)
real(psb_dpk_), allocatable :: val(:)
real(psb_dpk_), allocatable :: v(:), check(:)
integer(psb_ipk_) :: iictxt, icomm, flag
class(psb_xch_idx_type), pointer :: xchg
type(psb_desc_type):: desc_a
type(psb_d_vect_type) :: x
np = this%getNumImages()
if (np < 8) then
print*,'You need at least 8 processes to run this test.'
return
endif
call psb_init(icontxt,np,MPI_COMM_WORLD)
me = this_image()
!Allocate vectors
allocate(vg(nrows))
allocate(val(nrows))
allocate(v(nrows))
!Use only 8 processes
!Each process has a line
do i=1, nrows
vg(i)=i-1
enddo
if (me == 1) nz = 7
if ((me >= 2).and.(me <= 8)) nz = 1
if (me > 8) nz = 0
allocate(ia(nz),ja(nz))
if (me == 1) then
ia(1)=1
ja(1)=2
ia(2)=1
ja(2)=3
ia(3)=1
ja(3)=4
ia(4)=1
ja(4)=5
ia(5)=1
ja(5)=6
ia(6)=1
ja(6)=7
ia(7)=1
ja(7)=8
endif
if (me == 2) then
ia(1)=2
ja(1)=1
endif
if (me == 3) then
ia(1)=3
ja(1)=1
endif
if (me == 4) then
ia(1)=4
ja(1)=1
endif
if (me == 5) then
ia(1)=5
ja(1)=1
endif
if (me == 6) then
ia(1)=6
ja(1)=1
endif
if (me == 7) then
ia(1)=7
ja(1)=1
endif
if (me == 8) then
ia(1)=8
ja(1)=1
endif
do i=1,nrows
val(i)=me
enddo
call psb_cdall(icontxt,desc_a,info, vg=vg)
call psb_cdins(nz=nz,ia=ia,ja=ja, desc_a=desc_a, info=info)
call psb_cdasb(desc_a, info)
allocate(irw(nrows))
do i=1,nrows
irw(i)=i
enddo
call psb_geall(x,desc_a,info)
call psb_geins(m=nrows, irw=irw, val=val,x=x, desc_a=desc_a, info=info)
call psb_geasb(x,desc_a,info)
if (me == 1) nz = 8
if ((me >= 2).and.(me <= 8)) nz = 2
if (me > 8) nz = 1
if (allocated(check)) deallocate(check)
allocate (check(nz))
if (me == 1) then
check(1)=2
check(2)=4
check(3)=6
check(4)=8
check(5)=10
check(6)=12
check(7)=14
check(8)=16
endif
if (me == 2) then
check(1)=4
check(2)=2
endif
if (me == 3) then
check(1)=6
check(2)=2
endif
if (me == 4) then
check(1)=8
check(2)=2
endif
if (me == 5) then
check(1)=10
check(2)=2
endif
if (me == 6) then
check(1)=12
check(2)=2
endif
if (me == 7) then
check(1)=14
check(2)=2
endif
if (me == 8) then
check(1)=16
check(2)=2
endif
if (me > 8) then
check(1)=0
endif
! END OF SETUP
call psb_barrier(icontxt)
!We can do something better here
x%v%v = x%v%v + me
call psb_barrier(icontxt)
iictxt = desc_a%get_context()
icomm = desc_a%get_mpic()
call desc_a%get_list(psb_comm_halo_,xchg,info)
flag = IOR(psb_swap_send_, psb_swap_recv_)
call psi_dswap_xchg_vect(iictxt,icomm,flag,0.0d0,x%v,xchg,info)
call psb_barrier(icontxt)
v = x%get_vect()
!call psb_barrier(icontxt)
if ((me==1).or.(me==2)) then
true = 1
else
true=0
endif
@assertEqual(real(true*check),real(true*v))
deallocate(vg,ia,val,v,check)
call psb_gefree(x, desc_a, info)
call psb_cdfree(desc_a, info)
call psb_exit(icontxt)
end subroutine test_psb_swapdata_8imgs_b
@test(nimgs=[std])
subroutine test_psb_swapdatatran_2imgs(this)
implicit none
Class(CafTestMethod), intent(inout) :: this
integer :: msg, me, i=0, np, j, info
integer, parameter :: nrows=6
integer :: icontxt, mid, true
integer, allocatable :: vg(:), ia(:)
real(psb_dpk_), allocatable :: val(:)
real(psb_dpk_), allocatable :: v(:), check(:)
class(psb_xch_idx_type), pointer :: xchg
integer(psb_ipk_) :: iictxt, icomm, flag
type(psb_desc_type):: desc_a
type(psb_d_vect_type) :: x
np = this%getNumImages()
if (np < 2) then
print*,'You need at least 2 processes to run this test.'
return
endif
call psb_init(icontxt,np,MPI_COMM_WORLD)
!call psb_info(icontxt, me, np)
me = this_image()
!Allocate vectors
allocate(vg(nrows))
allocate(ia(nrows))
allocate(val(nrows))
allocate(v(nrows))
i = 0
do j=1,size(vg,1)
vg(j)= i
i = i+1
if (i==np) then
i=0
endif
enddo
!Use only 2 processes
!Assuming nrows is a multiple of 2 so mid is an integer
!Distribute equally to the two processes
mid=nrows/2
do i=1, mid
vg(i)=0
enddo
do i=mid+1, nrows
vg(i)=1
enddo
do i=1,size(ia,1)
ia(i)=i
enddo
do i=1,mid
val(i)=1.
enddo
do i=mid + 1,nrows
val(i)=2.
enddo
call psb_cdall(icontxt,desc_a,info, vg=vg)
if ( me == 1) call psb_cdins(nz=1,ia=(/mid/),ja=(/mid+1/), desc_a=desc_a, info=info)
if ( me == 2) call psb_cdins(nz=1,ia=(/mid+1/),ja=(/mid/), desc_a=desc_a, info=info)
call psb_cdasb(desc_a, info)
call psb_geall(x,desc_a,info)
call psb_geins(m=nrows, irw=ia, val=val,x=x, desc_a=desc_a, info=info)
call psb_geasb(x,desc_a,info)
call psb_barrier(icontxt)
!v = x%get_vect()
!Let's modify x, so we need to update halo indices
if ((me == 1).or.(me == 2)) then
call psb_geaxpby(dble(me),x,0.0d0,x,desc_a,info)
!x%v%v(mid +1)=x%v%v(mid+1) + 2.0d0
endif
call psb_barrier(icontxt)
! END OF SETUP
v = x%get_vect()
iictxt = desc_a%get_context()
icomm = desc_a%get_mpic()
call desc_a%get_list(psb_comm_halo_,xchg,info)
flag = IOR(psb_swap_send_, psb_swap_recv_)
sync all
call psi_dswaptran_xchg_vect(iictxt,icomm,flag,0.0d0,x%v,xchg,info)
!GETTING BACK X
call psb_barrier(icontxt)
v = x%get_vect()
sync all
!Let's build the expected solution
if ((me == 1).or.(me==2)) then
allocate(check(mid+1))
else
allocate(check(1))
endif
if (me == 1 )then
check(1:mid)=1.0d0
check(mid + 1)=2.0d0
else if (me == 2) then
check(1)=2.0d0
check(mid-1:mid)=4.0d0
check(mid + 1)=1.0d0
else
check(1)=0.0d0
endif
!call psb_barrier(icontxt)
if ((me==1).or.(me==2)) then
true = 1
else
true=0
endif
@assertEqual(real(true*check),real(true*v))
deallocate(vg,ia,val,v,check)
call psb_gefree(x, desc_a, info)
call psb_cdfree(desc_a, info)
call psb_exit(icontxt)
end subroutine test_psb_swapdatatran_2imgs
@test(nimgs=[std])
subroutine test_psb_sswapdatam_2imgs(this)
implicit none
Class(CafTestMethod), intent(inout) :: this
integer :: msg, me, i=0, np, j, info
integer, parameter :: nrows=6
integer :: icontxt, mid, true
integer, allocatable :: vg(:), ia(:)
real(psb_spk_), allocatable :: val(:)
real(psb_spk_), allocatable :: y(:,:), check(:), v(:)
class(psb_xch_idx_type), pointer :: xchg
integer(psb_ipk_) :: iictxt, icomm, flag
type(psb_desc_type):: desc_a
type(psb_s_vect_type) :: x
np = this%getNumImages()
if (np < 2) then
print*,'You need at least 2 processes to run this test.'
return
endif
call psb_init(icontxt,np,MPI_COMM_WORLD)
!call psb_info(icontxt, me, np)
me = this_image()
!Allocate vectors
allocate(vg(nrows))
allocate(ia(nrows))
allocate(val(nrows))
allocate(v(nrows))
i = 0
do j=1,size(vg,1)
vg(j)= i
i = i+1
if (i==np) then
i=0
endif
enddo
!Use only 2 processes
!Assuming nrows is a multiple of 2 so mid is an integer
!Distribute equally to the two processes
mid=nrows/2
do i=1, mid
vg(i)=0
enddo
do i=mid+1, nrows
vg(i)=1
enddo
do i=1,size(ia,1)
ia(i)=i
enddo
do i=1,mid
val(i)=1.
enddo
do i=mid + 1,nrows
val(i)=2.
enddo
call psb_cdall(icontxt,desc_a,info, vg=vg)
if ( me == 1) call psb_cdins(nz=1,ia=(/mid/),ja=(/mid+1/), desc_a=desc_a, info=info)
if ( me == 2) call psb_cdins(nz=1,ia=(/mid+1/),ja=(/mid/), desc_a=desc_a, info=info)
call psb_cdasb(desc_a, info)
call psb_geall(x,desc_a,info)
call psb_geins(m=nrows, irw=ia, val=val,x=x, desc_a=desc_a, info=info)
call psb_geasb(x,desc_a,info)
call psb_barrier(icontxt)
v = x%get_vect()
allocate(y(size(v,1),1))
y(:,1)=v
!Let's modify x, so we need to update halo indices
if ((me == 1).or.(me == 2)) then
y(mid +1,1)=y(mid+1,1) + 2.0
endif
call psb_barrier(icontxt)
! END OF SETUP
iictxt = desc_a%get_context()
icomm = desc_a%get_mpic()
call desc_a%get_list(psb_comm_halo_,xchg,info)
flag = IOR(psb_swap_send_, psb_swap_recv_)
print*,'size of y', size(y,1), size(y,2)
call psi_sswap_xchg_m(iictxt,icomm,flag,1,0.0,y,xchg,info)
!GETTING BACK X
call psb_barrier(icontxt)
!Let's build the expected solution
if (allocated(check)) deallocate(check)
if ((me == 1).or.(me==2)) then
allocate(check(mid+1))
else
allocate(check(1))
endif
if (me == 1 ) then
check(1:mid)=1.0
check(mid + 1)=2.0
else if (me == 2) then
check(1:mid)=2.0
check(mid + 1)=1.0
else
check(1)=0.0
endif
!call psb_barrier(icontxt)
if ((me==1).or.(me==2)) then
true = 1
else
true=0
endif
@assertEqual(real(true*check),real(true*y(:,1)))
deallocate(vg,ia,val,y,v,check)
call psb_gefree(x, desc_a, info)
call psb_cdfree(desc_a, info)
call psb_exit(icontxt)
end subroutine test_psb_sswapdatam_2imgs
end module test_psb_swapdata