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_dhalo.f90

430 lines
13 KiB
Fortran

!!$
!!$ Parallel Sparse BLAS v2.0
!!$ (C) Copyright 2006 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_dhalo.f90
!
! Subroutine: psb_dhalom
! 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 - real,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
!
20 years ago
!
subroutine psb_dhalom(x,desc_a,info,alpha,jx,ik,work,tran,mode,data)
20 years ago
use psb_descriptor_type
use psb_const_mod
use psi_mod
use psb_check_mod
use psb_realloc_mod
20 years ago
use psb_error_mod
17 years ago
use psb_string_mod
use psb_penv_mod
20 years ago
implicit none
real(kind(1.d0)), intent(inout), target :: x(:,:)
type(psb_desc_type), intent(in) :: desc_a
integer, intent(out) :: info
real(kind(1.d0)), intent(in), optional :: alpha
real(kind(1.d0)), 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, k, maxk, nrow, imode, i,&
& err, liwork,data_
20 years ago
real(kind(1.d0)),pointer :: iwork(:), xp(:,:)
17 years ago
character :: tran_
20 years ago
character(len=20) :: name, ch_err
logical :: aliw
20 years ago
name='psb_dhalom'
if(psb_get_errstatus().ne.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
20 years ago
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.gt.maxk) then
k=maxk
else
k=ik
end if
20 years ago
else
k = maxk
20 years ago
end if
if (present(tran)) then
17 years ago
tran_ = 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)
20 years ago
if(info.ne.0) then
info=4010
ch_err='psb_chkvect'
call psb_errpush(info,name,a_err=ch_err)
20 years ago
end if
if (iix.ne.1) then
info=3040
call psb_errpush(info,name)
20 years ago
end if
err=info
call psb_errcomm(ictxt,err)
20 years ago
if(err.ne.0) goto 9999
if(present(alpha)) then
if(alpha.ne.1.d0) then
do i=0, k-1
call dscal(nrow,alpha,x(1,jjx+i),1)
end do
end if
20 years ago
end if
liwork=nrow
if (present(work)) then
if(size(work).ge.liwork) then
iwork => work
aliw=.false.
else
aliw=.true.
allocate(iwork(liwork),stat=info)
if(info.ne.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.
!!$ write(0,*) 'halom ',liwork
allocate(iwork(liwork),stat=info)
if(info.ne.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
xp => x(iix:size(x,1),jjx:jjx+k-1)
17 years ago
if(tran_ == 'N') then
call psi_swapdata(imode,k,dzero,xp,&
& desc_a,iwork,info,data=data_)
17 years ago
else if((tran_ == 'T').or.(tran_ == 'C')) then
call psi_swaptran(imode,k,done,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.ne.0) then
ch_err='PSI_dSwapdata'
call psb_errpush(4010,name,a_err=ch_err)
goto 9999
20 years ago
end if
if (aliw) deallocate(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_dhalom
!!$
!!$ Parallel Sparse BLAS v2.0
!!$ (C) Copyright 2006 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_dhalov
! This subroutine performs the exchange of the halo elements in a
! distributed dense vector between all the processes.
20 years ago
!
! Arguments:
20 years ago
! x - real,dimension(:). The local part of the dense vector.
! desc_a - type(psb_desc_type). The communication descriptor.
! info - integer. Return code
! alpha - real(optional). Scale factor.
! 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
!
20 years ago
!
subroutine psb_dhalov(x,desc_a,info,alpha,work,tran,mode,data)
20 years ago
use psb_descriptor_type
use psb_const_mod
use psi_mod
use psb_check_mod
use psb_realloc_mod
20 years ago
use psb_error_mod
17 years ago
use psb_string_mod
use psb_penv_mod
20 years ago
implicit none
real(kind(1.d0)), intent(inout) :: x(:)
type(psb_desc_type), intent(in) :: desc_a
integer, intent(out) :: info
real(kind(1.d0)), intent(in), optional :: alpha
real(kind(1.d0)), target, optional :: 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
real(kind(1.d0)),pointer :: iwork(:)
17 years ago
character :: tran_
20 years ago
character(len=20) :: name, ch_err
logical :: aliw
20 years ago
name='psb_dhalov'
if(psb_get_errstatus().ne.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)
20 years ago
if (present(tran)) then
17 years ago
tran_ = 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)
20 years ago
if(info.ne.0) then
info=4010
ch_err='psb_chkvect'
call psb_errpush(info,name,a_err=ch_err)
20 years ago
end if
if (iix.ne.1) then
info=3040
call psb_errpush(info,name)
20 years ago
end if
err=info
call psb_errcomm(ictxt,err)
20 years ago
if(err.ne.0) goto 9999
if(present(alpha)) then
if(alpha.ne.1.d0) then
call dscal(nrow,alpha,x,ione)
end if
20 years ago
end if
liwork=nrow
if (present(work)) then
if(size(work).ge.liwork) then
iwork => work
aliw=.false.
else
aliw=.true.
allocate(iwork(liwork),stat=info)
if(info.ne.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.ne.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
20 years ago
! exchange halo elements
17 years ago
if(tran_ == 'N') then
call psi_swapdata(imode,dzero,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,done,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.ne.0) then
ch_err='PSI_swapdata'
call psb_errpush(4010,name,a_err=ch_err)
goto 9999
20 years ago
end if
if (aliw) deallocate(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_dhalov