Added routine for logical all reduce and applied in psb_mask

merge-paraggr-newops
Cirdans-Home 5 years ago
parent 17e24bdcf0
commit 7f42d63275

@ -1,9 +1,9 @@
!
!
! Parallel Sparse BLAS version 3.5
! (C) Copyright 2006-2018
! Salvatore Filippone
! Alfredo Buttari
!
! Salvatore Filippone
! Alfredo Buttari
!
! Redistribution and use in source and binary forms, with or without
! modification, are permitted provided that the following conditions
! are met:
@ -15,7 +15,7 @@
! 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
@ -27,8 +27,8 @@
! 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.
!
!
!
!
module psi_collective_mod
use psi_penv_mod
use psi_m_collective_mod
@ -42,10 +42,10 @@ module psi_collective_mod
module procedure psb_hbcasts, psb_hbcastv,&
& psb_hbcasts_ec, psb_hbcastv_ec,&
& psb_lbcasts, psb_lbcastv, &
& psb_lbcasts_ec, psb_lbcastv_ec
& psb_lbcasts_ec, psb_lbcastv_ec
end interface psb_bcast
#if defined(SHORT_INTEGERS)
interface psb_sum
module procedure psb_i2sums, psb_i2sumv, psb_i2summ, &
@ -55,12 +55,12 @@ module psi_collective_mod
contains
subroutine psb_hbcasts(ictxt,dat,root,length)
#ifdef MPI_MOD
use mpi
#endif
implicit none
implicit none
#ifdef MPI_H
include 'mpif.h'
#endif
@ -85,7 +85,7 @@ contains
call psb_info(ictxt,iam,np)
call mpi_bcast(dat,length_,MPI_CHARACTER,root_,ictxt,info)
#endif
#endif
end subroutine psb_hbcasts
@ -93,7 +93,7 @@ contains
#ifdef MPI_MOD
use mpi
#endif
implicit none
implicit none
#ifdef MPI_H
include 'mpif.h'
#endif
@ -110,24 +110,24 @@ contains
root_ = psb_root_
endif
length_ = len(dat)
size_ = size(dat)
size_ = size(dat)
call psb_info(ictxt,iam,np)
call mpi_bcast(dat,length_*size_,MPI_CHARACTER,root_,ictxt,info)
#endif
#endif
end subroutine psb_hbcastv
subroutine psb_hbcasts_ec(ictxt,dat,root)
implicit none
implicit none
integer(psb_epk_), intent(in) :: ictxt
character(len=*), intent(inout) :: dat
integer(psb_epk_), intent(in), optional :: root
integer(psb_mpk_) :: ictxt_, root_
ictxt_ = ictxt
if (present(root)) then
if (present(root)) then
root_ = root
call psb_bcast(ictxt_,dat,root_)
else
@ -136,14 +136,14 @@ contains
end subroutine psb_hbcasts_ec
subroutine psb_hbcastv_ec(ictxt,dat,root)
implicit none
implicit none
integer(psb_epk_), intent(in) :: ictxt
character(len=*), intent(inout) :: dat(:)
integer(psb_epk_), intent(in), optional :: root
integer(psb_mpk_) :: ictxt_, root_
ictxt_ = ictxt
if (present(root)) then
if (present(root)) then
root_ = root
call psb_bcast(ictxt_,dat,root_)
else
@ -152,12 +152,12 @@ contains
end subroutine psb_hbcastv_ec
subroutine psb_lbcasts(ictxt,dat,root)
#ifdef MPI_MOD
use mpi
#endif
implicit none
implicit none
#ifdef MPI_H
include 'mpif.h'
#endif
@ -176,16 +176,41 @@ contains
call psb_info(ictxt,iam,np)
call mpi_bcast(dat,1,MPI_LOGICAL,root_,ictxt,info)
#endif
#endif
end subroutine psb_lbcasts
subroutine psb_lallreduceand(ictxt,dat,rec)
#ifdef MPI_MOD
use mpi
#endif
implicit none
#ifdef MPI_H
include 'mpif.h'
#endif
integer(psb_mpk_), intent(in) :: ictxt
logical, intent(inout) :: dat
logical, intent(inout), optional :: rec
integer(psb_mpk_) :: iam, np, info
#if !defined(SERIAL_MPI)
call psb_info(ictxt,iam,np)
if (present(rec)) then
call mpi_allreduce(dat,rec,1,MPI_LOGICAL,MPI_LAND,ictxt,info)
else
call mpi_allreduce(MPI_IN_PLACE,dat,1,MPI_LOGICAL,MPI_LAND,ictxt,info)
endif
#endif
end subroutine psb_lallreduceand
subroutine psb_lbcastv(ictxt,dat,root)
#ifdef MPI_MOD
use mpi
#endif
implicit none
implicit none
#ifdef MPI_H
include 'mpif.h'
#endif
@ -204,20 +229,20 @@ contains
call psb_info(ictxt,iam,np)
call mpi_bcast(dat,size(dat),MPI_LOGICAL,root_,ictxt,info)
#endif
#endif
end subroutine psb_lbcastv
subroutine psb_lbcasts_ec(ictxt,dat,root)
implicit none
implicit none
integer(psb_epk_), intent(in) :: ictxt
logical, intent(inout) :: dat
integer(psb_epk_), intent(in), optional :: root
integer(psb_mpk_) :: ictxt_, root_
ictxt_ = ictxt
if (present(root)) then
if (present(root)) then
root_ = root
call psb_bcast(ictxt_,dat,root_)
else
@ -226,14 +251,14 @@ contains
end subroutine psb_lbcasts_ec
subroutine psb_lbcastv_ec(ictxt,dat,root)
implicit none
implicit none
integer(psb_epk_), intent(in) :: ictxt
logical, intent(inout) :: dat(:)
integer(psb_epk_), intent(in), optional :: root
integer(psb_mpk_) :: ictxt_, root_
ictxt_ = ictxt
if (present(root)) then
if (present(root)) then
root_ = root
call psb_bcast(ictxt_,dat,root_)
else
@ -242,14 +267,14 @@ contains
end subroutine psb_lbcastv_ec
#if defined(SHORT_INTEGERS)
subroutine psb_i2sums(ictxt,dat,root)
#ifdef MPI_MOD
use mpi
#endif
implicit none
implicit none
#ifdef MPI_H
include 'mpif.h'
#endif
@ -265,12 +290,12 @@ contains
call psb_info(ictxt,iam,np)
if (present(root)) then
if (present(root)) then
root_ = root
else
root_ = -1
endif
if (root_ == -1) then
if (root_ == -1) then
call mpi_allreduce(dat,dat_,1,psb_mpi_i2pk_,mpi_sum,ictxt,info)
dat = dat_
else
@ -278,7 +303,7 @@ contains
if (iam == root_) dat = dat_
endif
#endif
#endif
end subroutine psb_i2sums
subroutine psb_i2sumv(ictxt,dat,root)
@ -286,7 +311,7 @@ contains
#ifdef MPI_MOD
use mpi
#endif
implicit none
implicit none
#ifdef MPI_H
include 'mpif.h'
#endif
@ -302,18 +327,18 @@ contains
call psb_info(ictxt,iam,np)
if (present(root)) then
if (present(root)) then
root_ = root
else
root_ = -1
endif
if (root_ == -1) then
if (root_ == -1) then
call psb_realloc(size(dat),dat_,iinfo)
dat_=dat
if (iinfo == psb_success_) call mpi_allreduce(dat_,dat,size(dat),&
& psb_mpi_i2pk_,mpi_sum,ictxt,info)
else
if (iam == root_) then
if (iam == root_) then
call psb_realloc(size(dat),dat_,iinfo)
dat_=dat
call mpi_reduce(dat_,dat,size(dat),psb_mpi_i2pk_,mpi_sum,root_,ictxt,info)
@ -321,7 +346,7 @@ contains
call mpi_reduce(dat,dat_,size(dat),psb_mpi_i2pk_,mpi_sum,root_,ictxt,info)
end if
endif
#endif
#endif
end subroutine psb_i2sumv
subroutine psb_i2summ(ictxt,dat,root)
@ -329,7 +354,7 @@ contains
#ifdef MPI_MOD
use mpi
#endif
implicit none
implicit none
#ifdef MPI_H
include 'mpif.h'
#endif
@ -345,18 +370,18 @@ contains
#if !defined(SERIAL_MPI)
call psb_info(ictxt,iam,np)
if (present(root)) then
if (present(root)) then
root_ = root
else
root_ = -1
endif
if (root_ == -1) then
if (root_ == -1) then
call psb_realloc(size(dat,1),size(dat,2),dat_,iinfo)
dat_=dat
if (iinfo == psb_success_) call mpi_allreduce(dat_,dat,size(dat),&
& psb_mpi_i2pk_,mpi_sum,ictxt,info)
else
if (iam == root_) then
if (iam == root_) then
call psb_realloc(size(dat,1),size(dat,2),dat_,iinfo)
dat_=dat
call mpi_reduce(dat_,dat,size(dat),psb_mpi_i2pk_,mpi_sum,root_,ictxt,info)
@ -364,18 +389,18 @@ contains
call mpi_reduce(dat,dat_,size(dat),psb_mpi_i2pk_,mpi_sum,root_,ictxt,info)
end if
endif
#endif
#endif
end subroutine psb_i2summ
subroutine psb_i2sums_ec(ictxt,dat,root)
implicit none
implicit none
integer(psb_epk_), intent(in) :: ictxt
integer(psb_i2pk_), intent(inout) :: dat
integer(psb_epk_), intent(in), optional :: root
integer(psb_mpk_) :: ictxt_, root_
ictxt_ = ictxt
if (present(root)) then
if (present(root)) then
root_ = root
call psb_sum(ictxt_,dat,root_)
else
@ -384,14 +409,14 @@ contains
end subroutine psb_i2sums_ec
subroutine psb_i2sumv_ec(ictxt,dat,root)
implicit none
implicit none
integer(psb_epk_), intent(in) :: ictxt
integer(psb_i2pk_), intent(inout) :: dat(:)
integer(psb_epk_), intent(in), optional :: root
integer(psb_mpk_) :: ictxt_, root_
ictxt_ = ictxt
if (present(root)) then
if (present(root)) then
root_ = root
call psb_sum(ictxt_,dat,root_)
else
@ -400,14 +425,14 @@ contains
end subroutine psb_i2sumv_ec
subroutine psb_i2summ_ec(ictxt,dat,root)
implicit none
implicit none
integer(psb_epk_), intent(in) :: ictxt
integer(psb_i2pk_), intent(inout) :: dat(:,:)
integer(psb_epk_), intent(in), optional :: root
integer(psb_mpk_) :: ictxt_, root_
ictxt_ = ictxt
if (present(root)) then
if (present(root)) then
root_ = root
call psb_sum(ictxt_,dat,root_)
else

@ -214,6 +214,8 @@ subroutine psb_dmask_vect(c,x,m,t,desc_a,info)
call m%mask(c,x,t,info)
end if
call psb_lallreduceand(ictxt,t)
call psb_erractionrestore(err_act)
return

@ -214,6 +214,8 @@ subroutine psb_smask_vect(c,x,m,t,desc_a,info)
call m%mask(c,x,t,info)
end if
call psb_lallreduceand(ictxt,t)
call psb_erractionrestore(err_act)
return

Loading…
Cancel
Save