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/base/comm/psb_ihalo.f90

428 lines
13 KiB
Fortran

!!$
!!$ Parallel Sparse BLAS version 2.2
!!$ (C) Copyright 2006/2007/2008
!!$ Salvatore Filippone University of Rome Tor Vergata
!!$ Alfredo Buttari University of Rome Tor Vergata
!!$
!!$ Redistribution and use in source and binary forms, with or without
!!$ modification, are permitted provided that the following conditions
!!$ are met:
!!$ 1. Redistributions of source code must retain the above copyright
!!$ notice, this list of conditions and the following disclaimer.
!!$ 2. Redistributions in binary form must reproduce the above copyright
!!$ notice, this list of conditions, and the following disclaimer in the
!!$ documentation and/or other materials provided with the distribution.
!!$ 3. The name of the PSBLAS group or the names of its contributors may
!!$ not be used to endorse or promote products derived from this
!!$ software without specific written permission.
!!$
!!$ THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
!!$ ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
!!$ TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
!!$ PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE PSBLAS GROUP OR ITS CONTRIBUTORS
!!$ BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
!!$ CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
!!$ SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
!!$ INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
!!$ CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
!!$ ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
!!$ POSSIBILITY OF SUCH DAMAGE.
!!$
!!$
20 years ago
! File: psb_ihalo.f90
!
! Subroutine: psb_ihalom
! This subroutine performs the exchange of the halo elements in a
! distributed dense matrix between all the processes.
20 years ago
!
! Arguments:
20 years ago
! x - integer,dimension(:,:). The local part of the dense matrix.
! desc_a - type(psb_desc_type). The communication descriptor.
! info - integer. Return code
! alpha - real(optional). Scale factor.
20 years ago
! jx - integer(optional). The starting column of the global matrix.
! ik - integer(optional). The number of columns to gather.
! work - real(optional). Work area.
! tran - character(optional). Transpose exchange.
! mode - integer(optional). Communication mode (see Swapdata)
! data - integer Which index list in desc_a should be used
! to retrieve rows, default psb_comm_halo_
! psb_comm_halo_ use halo_index
! psb_comm_ext_ use ext_index
! psb_comm_ovrl_ use ovrl_index
! psb_comm_mov_ use ovr_mst_idx
!
20 years ago
!
subroutine psb_ihalom(x,desc_a,info,alpha,jx,ik,work,tran,mode,data)
use psb_sparse_mod, psb_protect_name => psb_ihalom
20 years ago
use psi_mod
implicit none
integer, intent(inout), target :: x(:,:)
type(psb_desc_type), intent(in) :: desc_a
integer, intent(out) :: info
real(psb_dpk_), intent(in), optional :: alpha
20 years ago
integer, intent(inout), optional, target :: work(:)
integer, intent(in), optional :: mode,jx,ik,data
20 years ago
character, intent(in), optional :: tran
! locals
integer :: ictxt, np, me, &
& err_act, m, n, iix, jjx, ix, ijx, nrow, k, maxk, liwork,&
& imode, err,data_
20 years ago
integer, pointer :: xp(:,:), iwork(:)
17 years ago
character :: tran_
20 years ago
character(len=20) :: name, ch_err
logical :: aliw
20 years ago
name='psb_ihalom'
if(psb_get_errstatus() /= 0) return
20 years ago
info=0
call psb_erractionsave(err_act)
ictxt=psb_cd_get_context(desc_a)
20 years ago
! check on blacs grid
call psb_info(ictxt, me, np)
if (np == -1) then
20 years ago
info = 2010
call psb_errpush(info,name)
goto 9999
endif
ix = 1
if (present(jx)) then
ijx = jx
20 years ago
else
ijx = 1
20 years ago
endif
m = psb_cd_get_global_rows(desc_a)
n = psb_cd_get_global_cols(desc_a)
nrow = psb_cd_get_local_rows(desc_a)
maxk=size(x,2)-ijx+1
20 years ago
if(present(ik)) then
if(ik > maxk) then
k=maxk
else
k=ik
end if
20 years ago
else
k = maxk
20 years ago
end if
if (present(tran)) then
psblas-2.2-maint: base/comm/psb_dhalo.f90 base/comm/psb_ihalo.f90 base/comm/psb_zhalo.f90 base/modules/psb_spmat_type.f90 base/modules/psb_string_mod.f90 base/psblas/psb_dspmm.f90 base/psblas/psb_dspsm.f90 base/psblas/psb_zspmm.f90 base/psblas/psb_zspsm.f90 base/serial/dp/dcoco.f base/serial/dp/dcocr.f base/serial/dp/dcrco.f base/serial/dp/dcrcr.f base/serial/dp/dcrjd.f base/serial/dp/dcsrp1.f base/serial/dp/dcsrrp.f base/serial/dp/djadrp.f base/serial/dp/djadrp1.f base/serial/dp/djdcox.f base/serial/dp/dvtfg.f base/serial/dp/zcoco.f base/serial/dp/zcocr.f base/serial/dp/zcrco.f base/serial/dp/zcrcr.f base/serial/dp/zcrjd.f base/serial/jad/djadsm.f base/serial/psb_cest.f90 base/serial/psb_dcoins.f90 base/serial/psb_dcsprt.f90 base/serial/psb_dfixcoo.f90 base/serial/psb_dipcoo2csc.f90 base/serial/psb_dipcoo2csr.f90 base/serial/psb_dipcsr2coo.f90 base/serial/psb_dnumbmm.f90 base/serial/psb_drwextd.f90 base/serial/psb_dspcnv.f90 base/serial/psb_dspgetrow.f90 base/serial/psb_dspscal.f90 base/serial/psb_dsymbmm.f90 base/serial/psb_dtransp.f90 base/serial/psb_lsame.f90 base/serial/psb_update_mod.f90 base/serial/psb_zcoins.f90 base/serial/psb_zcsprt.f90 base/serial/psb_zfixcoo.f90 base/serial/psb_zipcoo2csc.f90 base/serial/psb_zipcoo2csr.f90 base/serial/psb_zipcsr2coo.f90 base/serial/psb_znumbmm.f90 base/serial/psb_zrwextd.f90 base/serial/psb_zspcnv.f90 base/serial/psb_zspgetrow.f90 base/serial/psb_zspscal.f90 base/serial/psb_zsymbmm.f90 base/serial/psb_ztransc.f90 base/serial/psb_ztransp.f90 base/tools/psb_cdren.f90 base/tools/psb_dsphalo.F90 base/tools/psb_glob_to_loc.f90 base/tools/psb_loc_to_glob.f90 base/tools/psb_zsphalo.F90 krylov/psb_krylov_mod.f90 prec/psb_dbjac_aply.f90 prec/psb_dgprec_aply.f90 prec/psb_dprc_aply.f90 prec/psb_dprecbld.f90 prec/psb_dprecinit.f90 prec/psb_zbjac_aply.f90 prec/psb_zgprec_aply.f90 prec/psb_zprc_aply.f90 prec/psb_zprecbld.f90 prec/psb_zprecinit.f90 util/psb_hbio_mod.f90 util/psb_mat_dist_mod.f90 util/psb_metispart_mod.F90 util/psb_mmio_mod.f90 util/psb_read_mat_mod.f90 Fixed name of TOUPPER and friends with prefix PSB_.
17 years ago
tran_ = psb_toupper(tran)
20 years ago
else
17 years ago
tran_ = 'N'
20 years ago
endif
if (present(data)) then
data_ = data
else
data_ = psb_comm_halo_
endif
20 years ago
if (present(mode)) then
imode = mode
20 years ago
else
imode = IOR(psb_swap_send_,psb_swap_recv_)
20 years ago
endif
! check vector correctness
call psb_chkvect(m,1,size(x,1),ix,ijx,desc_a,info,iix,jjx)
if(info /= 0) then
info=4010
ch_err='psb_chkvect'
call psb_errpush(info,name,a_err=ch_err)
20 years ago
end if
if (iix /= 1) then
info=3040
call psb_errpush(info,name)
20 years ago
end if
err=info
call psb_errcomm(ictxt,err)
if(err /= 0) goto 9999
20 years ago
! we should write an "iscal"
!!$ if(present(alpha)) then
!!$ if(alpha /= 1.d0) then
20 years ago
!!$ do i=0, k-1
!!$ call iscal(nrow,alpha,x(1,jjx+i),1)
!!$ end do
!!$ end if
!!$ end if
liwork=nrow
if (present(work)) then
if(size(work) >= liwork) then
aliw=.false.
iwork => work
else
aliw=.true.
allocate(iwork(liwork),stat=info)
if(info /= 0) then
20 years ago
info=4010
ch_err='psb_realloc'
20 years ago
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
end if
else
aliw=.true.
allocate(iwork(liwork),stat=info)
if(info /= 0) then
info=4010
ch_err='psb_realloc'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
20 years ago
end if
xp => x(iix:size(x,1),jjx:jjx+k-1)
! exchange halo elements
17 years ago
if(tran_ == 'N') then
call psi_swapdata(imode,k,izero,xp,&
& desc_a,iwork,info,data=data_)
17 years ago
else if((tran_ == 'T').or.(tran_ == 'C')) then
call psi_swaptran(imode,k,ione,xp,&
& desc_a,iwork,info)
17 years ago
else
info = 4001
call psb_errpush(info,name,a_err='invalid tran')
goto 9999
20 years ago
end if
if(info /= 0) then
call psb_errpush(4010,name,a_err='PSI_iSwap...')
goto 9999
20 years ago
end if
if (aliw) deallocate(iwork)
nullify(iwork)
20 years ago
call psb_erractionrestore(err_act)
return
9999 continue
call psb_erractionrestore(err_act)
17 years ago
if (err_act == psb_act_abort_) then
call psb_error(ictxt)
return
20 years ago
end if
return
end subroutine psb_ihalom
!!$
!!$ Parallel Sparse BLAS version 2.2
!!$ (C) Copyright 2006/2007/2008
!!$ Salvatore Filippone University of Rome Tor Vergata
!!$ Alfredo Buttari University of Rome Tor Vergata
!!$
!!$ Redistribution and use in source and binary forms, with or without
!!$ modification, are permitted provided that the following conditions
!!$ are met:
!!$ 1. Redistributions of source code must retain the above copyright
!!$ notice, this list of conditions and the following disclaimer.
!!$ 2. Redistributions in binary form must reproduce the above copyright
!!$ notice, this list of conditions, and the following disclaimer in the
!!$ documentation and/or other materials provided with the distribution.
!!$ 3. The name of the PSBLAS group or the names of its contributors may
!!$ not be used to endorse or promote products derived from this
!!$ software without specific written permission.
!!$
!!$ THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
!!$ ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
!!$ TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
!!$ PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE PSBLAS GROUP OR ITS CONTRIBUTORS
!!$ BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
!!$ CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
!!$ SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
!!$ INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
!!$ CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
!!$ ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
!!$ POSSIBILITY OF SUCH DAMAGE.
!!$
!!$
20 years ago
! Subroutine: psb_ihalov
! This subroutine performs the exchange of the halo elements in a
! distributed dense matrix between all the processes.
20 years ago
!
! Arguments:
20 years ago
! x - integer,dimension(:). The local part of the dense matrix.
! desc_a - type(psb_desc_type). The communication descriptor.
! info - integer. Return code
! alpha - real(optional). Scale factor.
! jx - integer(optional). The starting column of the global matrix.
! ik - integer(optional). The number of columns to gather.
! work - real(optional). Work area.
! tran - character(optional). Transpose exchange.
! mode - integer(optional). Communication mode (see Swapdata)
! data - integer Which index list in desc_a should be used
! to retrieve rows, default psb_comm_halo_
! psb_comm_halo_ use halo_index
! psb_comm_ext_ use ext_index
! psb_comm_ovrl_ use ovrl_index
! psb_comm_mov_ use ovr_mst_idx
!
20 years ago
!
subroutine psb_ihalov(x,desc_a,info,alpha,work,tran,mode,data)
use psb_sparse_mod, psb_protect_name => psb_ihalov
20 years ago
use psi_mod
implicit none
integer, intent(inout) :: x(:)
type(psb_desc_type), intent(in) :: desc_a
integer, intent(out) :: info
real(psb_dpk_), intent(in), optional :: alpha
20 years ago
integer, intent(inout), optional, target :: work(:)
integer, intent(in), optional :: mode,data
20 years ago
character, intent(in), optional :: tran
! locals
integer :: ictxt, np, me,&
& err_act, m, n, iix, jjx, ix, ijx, nrow, imode,&
& err, liwork, data_
20 years ago
integer,pointer :: iwork(:)
17 years ago
character :: tran_
20 years ago
character(len=20) :: name, ch_err
logical :: aliw
20 years ago
name='psb_ihalov'
if(psb_get_errstatus() /= 0) return
20 years ago
info=0
call psb_erractionsave(err_act)
ictxt=psb_cd_get_context(desc_a)
20 years ago
! check on blacs grid
call psb_info(ictxt, me, np)
if (np == -1) then
20 years ago
info = 2010
call psb_errpush(info,name)
goto 9999
endif
ix = 1
ijx = 1
m = psb_cd_get_global_rows(desc_a)
n = psb_cd_get_global_cols(desc_a)
nrow = psb_cd_get_local_rows(desc_a)
! ncol = psb_cd_get_local_cols(desc_a)
20 years ago
if (present(tran)) then
psblas-2.2-maint: base/comm/psb_dhalo.f90 base/comm/psb_ihalo.f90 base/comm/psb_zhalo.f90 base/modules/psb_spmat_type.f90 base/modules/psb_string_mod.f90 base/psblas/psb_dspmm.f90 base/psblas/psb_dspsm.f90 base/psblas/psb_zspmm.f90 base/psblas/psb_zspsm.f90 base/serial/dp/dcoco.f base/serial/dp/dcocr.f base/serial/dp/dcrco.f base/serial/dp/dcrcr.f base/serial/dp/dcrjd.f base/serial/dp/dcsrp1.f base/serial/dp/dcsrrp.f base/serial/dp/djadrp.f base/serial/dp/djadrp1.f base/serial/dp/djdcox.f base/serial/dp/dvtfg.f base/serial/dp/zcoco.f base/serial/dp/zcocr.f base/serial/dp/zcrco.f base/serial/dp/zcrcr.f base/serial/dp/zcrjd.f base/serial/jad/djadsm.f base/serial/psb_cest.f90 base/serial/psb_dcoins.f90 base/serial/psb_dcsprt.f90 base/serial/psb_dfixcoo.f90 base/serial/psb_dipcoo2csc.f90 base/serial/psb_dipcoo2csr.f90 base/serial/psb_dipcsr2coo.f90 base/serial/psb_dnumbmm.f90 base/serial/psb_drwextd.f90 base/serial/psb_dspcnv.f90 base/serial/psb_dspgetrow.f90 base/serial/psb_dspscal.f90 base/serial/psb_dsymbmm.f90 base/serial/psb_dtransp.f90 base/serial/psb_lsame.f90 base/serial/psb_update_mod.f90 base/serial/psb_zcoins.f90 base/serial/psb_zcsprt.f90 base/serial/psb_zfixcoo.f90 base/serial/psb_zipcoo2csc.f90 base/serial/psb_zipcoo2csr.f90 base/serial/psb_zipcsr2coo.f90 base/serial/psb_znumbmm.f90 base/serial/psb_zrwextd.f90 base/serial/psb_zspcnv.f90 base/serial/psb_zspgetrow.f90 base/serial/psb_zspscal.f90 base/serial/psb_zsymbmm.f90 base/serial/psb_ztransc.f90 base/serial/psb_ztransp.f90 base/tools/psb_cdren.f90 base/tools/psb_dsphalo.F90 base/tools/psb_glob_to_loc.f90 base/tools/psb_loc_to_glob.f90 base/tools/psb_zsphalo.F90 krylov/psb_krylov_mod.f90 prec/psb_dbjac_aply.f90 prec/psb_dgprec_aply.f90 prec/psb_dprc_aply.f90 prec/psb_dprecbld.f90 prec/psb_dprecinit.f90 prec/psb_zbjac_aply.f90 prec/psb_zgprec_aply.f90 prec/psb_zprc_aply.f90 prec/psb_zprecbld.f90 prec/psb_zprecinit.f90 util/psb_hbio_mod.f90 util/psb_mat_dist_mod.f90 util/psb_metispart_mod.F90 util/psb_mmio_mod.f90 util/psb_read_mat_mod.f90 Fixed name of TOUPPER and friends with prefix PSB_.
17 years ago
tran_ = psb_toupper(tran)
20 years ago
else
17 years ago
tran_ = 'N'
20 years ago
endif
if (present(data)) then
data_ = data
else
data_ = psb_comm_halo_
endif
20 years ago
if (present(mode)) then
imode = mode
20 years ago
else
imode = IOR(psb_swap_send_,psb_swap_recv_)
20 years ago
endif
! check vector correctness
call psb_chkvect(m,1,size(x,1),ix,ijx,desc_a,info,iix,jjx)
if(info /= 0) then
info=4010
ch_err='psb_chkvect'
call psb_errpush(info,name,a_err=ch_err)
20 years ago
end if
if (iix /= 1) then
info=3040
call psb_errpush(info,name)
20 years ago
end if
err=info
call psb_errcomm(ictxt,err)
if(err /= 0) goto 9999
20 years ago
!!$ if(present(alpha)) then
!!$ if(alpha /= 1.d0) then
20 years ago
!!$ call dscal(nrow,alpha,x,1)
!!$ end if
!!$ end if
liwork=nrow
if (present(work)) then
if(size(work) >= liwork) then
aliw=.false.
iwork => work
else
aliw=.true.
allocate(iwork(liwork),stat=info)
if(info /= 0) then
20 years ago
info=4010
ch_err='psb_realloc'
20 years ago
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
end if
else
aliw=.true.
allocate(iwork(liwork),stat=info)
if(info /= 0) then
info=4010
ch_err='psb_realloc'
call psb_errpush(info,name,a_err=ch_err)
goto 9999
end if
20 years ago
end if
! exchange halo elements
17 years ago
if(tran_ == 'N') then
call psi_swapdata(imode,izero,x(iix:size(x)),&
& desc_a,iwork,info,data=data_)
17 years ago
else if((tran_ == 'T').or.(tran_ == 'C')) then
call psi_swaptran(imode,ione,x(iix:size(x)),&
& desc_a,iwork,info)
17 years ago
else
info = 4001
call psb_errpush(info,name,a_err='invalid tran')
goto 9999
20 years ago
end if
if(info /= 0) then
call psb_errpush(4010,name,a_err='PSI_iswapdata')
goto 9999
20 years ago
end if
if (aliw) deallocate(iwork)
nullify(iwork)
20 years ago
call psb_erractionrestore(err_act)
return
9999 continue
call psb_erractionrestore(err_act)
17 years ago
if (err_act == psb_act_abort_) then
call psb_error(ictxt)
return
20 years ago
end if
return
end subroutine psb_ihalov