! ! Parallel Sparse BLAS version 3.5 ! (C) Copyright 2006-2018 ! 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: ! 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. ! ! module psi_s_collective_mod use psi_penv_mod interface psb_max module procedure psb_smaxs, psb_smaxv, psb_smaxm, & & psb_smaxs_ec, psb_smaxv_ec, psb_smaxm_ec end interface interface psb_min module procedure psb_smins, psb_sminv, psb_sminm, & & psb_smins_ec, psb_sminv_ec, psb_sminm_ec end interface psb_min interface psb_nrm2 module procedure psb_s_nrm2s, psb_s_nrm2v, & & psb_s_nrm2s_ec, psb_s_nrm2v_ec end interface psb_nrm2 interface psb_sum module procedure psb_ssums, psb_ssumv, psb_ssumm, & & psb_ssums_ec, psb_ssumv_ec, psb_ssumm_ec end interface interface psb_amx module procedure psb_samxs, psb_samxv, psb_samxm, & & psb_samxs_ec, psb_samxv_ec, psb_samxm_ec end interface interface psb_amn module procedure psb_samns, psb_samnv, psb_samnm, & & psb_samns_ec, psb_samnv_ec, psb_samnm_ec end interface interface psb_bcast module procedure psb_sbcasts, psb_sbcastv, psb_sbcastm, & & psb_sbcasts_ec, psb_sbcastv_ec, psb_sbcastm_ec end interface psb_bcast interface psb_scan_sum module procedure psb_sscan_sums, psb_sscan_sumv end interface psb_scan_sum interface psb_exscan_sum module procedure psb_sexscan_sums, psb_sexscan_sumv end interface psb_exscan_sum contains ! !!!!!!!!!!!!!!!!!!!!!! ! ! Reduction operations ! ! !!!!!!!!!!!!!!!!!!!!!! ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!! ! ! MAX ! ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!! subroutine psb_smaxs(ictxt,dat,root) #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_) :: dat_ integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif if (root_ == -1) then call mpi_allreduce(dat,dat_,1,psb_mpi_r_spk_,mpi_max,ictxt,info) dat = dat_ else call mpi_reduce(dat,dat_,1,psb_mpi_r_spk_,mpi_max,root_,ictxt,info) if (iam == root_) dat = dat_ endif #endif end subroutine psb_smaxs subroutine psb_smaxv(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_), allocatable :: dat_(:) integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif 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_r_spk_,mpi_max,ictxt,info) else if (iam == root_) then call psb_realloc(size(dat),dat_,iinfo) dat_ = dat call mpi_reduce(dat_,dat,size(dat),psb_mpi_r_spk_,mpi_max,root_,ictxt,info) else call psb_realloc(1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_spk_,mpi_max,root_,ictxt,info) end if endif #endif end subroutine psb_smaxv subroutine psb_smaxm(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:,:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_), allocatable :: dat_(:,:) integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif 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_r_spk_,mpi_max,ictxt,info) else 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_r_spk_,mpi_max,root_,ictxt,info) else call psb_realloc(1,1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_spk_,mpi_max,root_,ictxt,info) end if endif #endif end subroutine psb_smaxm subroutine psb_smaxs_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_max(ictxt_,dat,root_) else call psb_max(ictxt_,dat) end if end subroutine psb_smaxs_ec subroutine psb_smaxv_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_max(ictxt_,dat,root_) else call psb_max(ictxt_,dat) end if end subroutine psb_smaxv_ec subroutine psb_smaxm_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:,:) integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_max(ictxt_,dat,root_) else call psb_max(ictxt_,dat) end if end subroutine psb_smaxm_ec ! ! MIN: Minimum Value ! subroutine psb_smins(ictxt,dat,root) #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_) :: dat_ integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif if (root_ == -1) then call mpi_allreduce(dat,dat_,1,psb_mpi_r_spk_,mpi_min,ictxt,info) dat = dat_ else call mpi_reduce(dat,dat_,1,psb_mpi_r_spk_,mpi_min,root_,ictxt,info) if (iam == root_) dat = dat_ endif #endif end subroutine psb_smins subroutine psb_sminv(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_), allocatable :: dat_(:) integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif 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_r_spk_,mpi_min,ictxt,info) else if (iam == root_) then call psb_realloc(size(dat),dat_,iinfo) dat_ = dat call mpi_reduce(dat_,dat,size(dat),psb_mpi_r_spk_,mpi_min,root_,ictxt,info) else call psb_realloc(1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_spk_,mpi_min,root_,ictxt,info) end if endif #endif end subroutine psb_sminv subroutine psb_sminm(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:,:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_), allocatable :: dat_(:,:) integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif 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_r_spk_,mpi_min,ictxt,info) else 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_r_spk_,mpi_min,root_,ictxt,info) else call psb_realloc(1,1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_spk_,mpi_min,root_,ictxt,info) end if endif #endif end subroutine psb_sminm subroutine psb_smins_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_min(ictxt_,dat,root_) else call psb_min(ictxt_,dat) end if end subroutine psb_smins_ec subroutine psb_sminv_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_min(ictxt_,dat,root_) else call psb_min(ictxt_,dat) end if end subroutine psb_sminv_ec subroutine psb_sminm_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:,:) integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_min(ictxt_,dat,root_) else call psb_min(ictxt_,dat) end if end subroutine psb_sminm_ec ! !!!!!!!!!!!! ! ! Norm 2, only for reals ! ! !!!!!!!!!!!! subroutine psb_s_nrm2s(ictxt,dat,root) #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_) :: dat_ integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif if (root_ == -1) then call mpi_allreduce(dat,dat_,1,psb_mpi_r_spk_,mpi_snrm2_op,ictxt,info) dat = dat_ else call mpi_reduce(dat,dat_,1,psb_mpi_r_spk_,mpi_snrm2_op,root_,ictxt,info) if (iam == root_) dat = dat_ endif #endif end subroutine psb_s_nrm2s subroutine psb_s_nrm2v(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_), allocatable :: dat_(:) integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif 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_r_spk_,& & mpi_snrm2_op,ictxt,info) else if (iam == root_) then call psb_realloc(size(dat),dat_,iinfo) dat_ = dat call mpi_reduce(dat_,dat,size(dat),psb_mpi_r_spk_,& & mpi_snrm2_op,root_,ictxt,info) else call psb_realloc(1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_spk_,& & mpi_snrm2_op,root_,ictxt,info) end if endif #endif end subroutine psb_s_nrm2v subroutine psb_s_nrm2s_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_nrm2(ictxt_,dat,root_) else call psb_nrm2(ictxt_,dat) end if end subroutine psb_s_nrm2s_ec subroutine psb_s_nrm2v_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_nrm2(ictxt_,dat,root_) else call psb_nrm2(ictxt_,dat) end if end subroutine psb_s_nrm2v_ec ! ! SUM ! subroutine psb_ssums(ictxt,dat,root) #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_) :: dat_ integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif if (root_ == -1) then call mpi_allreduce(dat,dat_,1,psb_mpi_r_spk_,mpi_sum,ictxt,info) dat = dat_ else call mpi_reduce(dat,dat_,1,psb_mpi_r_spk_,mpi_sum,root_,ictxt,info) if (iam == root_) dat = dat_ endif #endif end subroutine psb_ssums subroutine psb_ssumv(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_), allocatable :: dat_(:) integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif 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_r_spk_,mpi_sum,ictxt,info) else if (iam == root_) then call psb_realloc(size(dat),dat_,iinfo) dat_ = dat call mpi_reduce(dat_,dat,size(dat),psb_mpi_r_spk_,mpi_sum,root_,ictxt,info) else call psb_realloc(1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_spk_,mpi_sum,root_,ictxt,info) end if endif #endif end subroutine psb_ssumv subroutine psb_ssumm(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:,:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_), allocatable :: dat_(:,:) integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif 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_r_spk_,mpi_sum,ictxt,info) else 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_r_spk_,mpi_sum,root_,ictxt,info) else call psb_realloc(1,1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_spk_,mpi_sum,root_,ictxt,info) end if endif #endif end subroutine psb_ssumm subroutine psb_ssums_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_sum(ictxt_,dat,root_) else call psb_sum(ictxt_,dat) end if end subroutine psb_ssums_ec subroutine psb_ssumv_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_sum(ictxt_,dat,root_) else call psb_sum(ictxt_,dat) end if end subroutine psb_ssumv_ec subroutine psb_ssumm_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:,:) integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_sum(ictxt_,dat,root_) else call psb_sum(ictxt_,dat) end if end subroutine psb_ssumm_ec ! ! AMX: Maximum Absolute Value ! subroutine psb_samxs(ictxt,dat,root) #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_) :: dat_ integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif if (root_ == -1) then call mpi_allreduce(dat,dat_,1,psb_mpi_r_spk_,mpi_samx_op,ictxt,info) dat = dat_ else call mpi_reduce(dat,dat_,1,psb_mpi_r_spk_,mpi_samx_op,root_,ictxt,info) if (iam == root_) dat = dat_ endif #endif end subroutine psb_samxs subroutine psb_samxv(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_), allocatable :: dat_(:) integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif 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_r_spk_,mpi_samx_op,ictxt,info) else if (iam == root_) then call psb_realloc(size(dat),dat_,iinfo) dat_ = dat call mpi_reduce(dat_,dat,size(dat),psb_mpi_r_spk_,mpi_samx_op,root_,ictxt,info) else call psb_realloc(1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_spk_,mpi_samx_op,root_,ictxt,info) end if endif #endif end subroutine psb_samxv subroutine psb_samxm(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:,:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_), allocatable :: dat_(:,:) integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif 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_r_spk_,mpi_samx_op,ictxt,info) else 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_r_spk_,mpi_samx_op,root_,ictxt,info) else call psb_realloc(1,1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_spk_,mpi_samx_op,root_,ictxt,info) end if endif #endif end subroutine psb_samxm subroutine psb_samxs_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_amx(ictxt_,dat,root_) else call psb_amx(ictxt_,dat) end if end subroutine psb_samxs_ec subroutine psb_samxv_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_amx(ictxt_,dat,root_) else call psb_amx(ictxt_,dat) end if end subroutine psb_samxv_ec subroutine psb_samxm_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:,:) integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_amx(ictxt_,dat,root_) else call psb_amx(ictxt_,dat) end if end subroutine psb_samxm_ec ! ! AMN: Minimum Absolute Value ! subroutine psb_samns(ictxt,dat,root) #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_) :: dat_ integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif if (root_ == -1) then call mpi_allreduce(dat,dat_,1,psb_mpi_r_spk_,mpi_samn_op,ictxt,info) dat = dat_ else call mpi_reduce(dat,dat_,1,psb_mpi_r_spk_,mpi_samn_op,root_,ictxt,info) if (iam == root_) dat = dat_ endif #endif end subroutine psb_samns subroutine psb_samnv(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_), allocatable :: dat_(:) integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif 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_r_spk_,mpi_samn_op,ictxt,info) else if (iam == root_) then call psb_realloc(size(dat),dat_,iinfo) dat_ = dat call mpi_reduce(dat_,dat,size(dat),psb_mpi_r_spk_,mpi_samn_op,root_,ictxt,info) else call psb_realloc(1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_spk_,mpi_samn_op,root_,ictxt,info) end if endif #endif end subroutine psb_samnv subroutine psb_samnm(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:,:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_), allocatable :: dat_(:,:) integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = -1 endif 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_r_spk_,mpi_samn_op,ictxt,info) else 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_r_spk_,mpi_samn_op,root_,ictxt,info) else call psb_realloc(1,1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_spk_,mpi_samn_op,root_,ictxt,info) end if endif #endif end subroutine psb_samnm subroutine psb_samns_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_amn(ictxt_,dat,root_) else call psb_amn(ictxt_,dat) end if end subroutine psb_samns_ec subroutine psb_samnv_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_amn(ictxt_,dat,root_) else call psb_amn(ictxt_,dat) end if end subroutine psb_samnv_ec subroutine psb_samnm_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:,:) integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_amn(ictxt_,dat,root_) else call psb_amn(ictxt_,dat) end if end subroutine psb_samnm_ec ! ! BCAST Broadcast ! subroutine psb_sbcasts(ictxt,dat,root) #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = psb_root_ endif call mpi_bcast(dat,1,psb_mpi_r_spk_,root_,ictxt,info) #endif end subroutine psb_sbcasts subroutine psb_sbcastv(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = psb_root_ endif call mpi_bcast(dat,size(dat),psb_mpi_r_spk_,root_,ictxt,info) #endif end subroutine psb_sbcastv subroutine psb_sbcastm(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:,:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) if (present(root)) then root_ = root else root_ = psb_root_ endif call mpi_bcast(dat,size(dat),psb_mpi_r_spk_,root_,ictxt,info) #endif end subroutine psb_sbcastm subroutine psb_sbcasts_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_bcast(ictxt_,dat,root_) else call psb_bcast(ictxt_,dat) end if end subroutine psb_sbcasts_ec subroutine psb_sbcastv_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_bcast(ictxt_,dat,root_) else call psb_bcast(ictxt_,dat) end if end subroutine psb_sbcastv_ec subroutine psb_sbcastm_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:,:) integer(psb_epk_), intent(in), optional :: root integer(psb_mpk_) :: ictxt_, root_ ictxt_ = ictxt if (present(root)) then root_ = root call psb_bcast(ictxt_,dat,root_) else call psb_bcast(ictxt_,dat) end if end subroutine psb_sbcastm_ec ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!! ! ! SCAN ! ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!! subroutine psb_sscan_sums(ictxt,dat) #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat real(psb_spk_) :: dat_ integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) call mpi_scan(dat,dat_,1,psb_mpi_r_spk_,mpi_sum,ictxt,info) dat = dat_ #endif end subroutine psb_sscan_sums subroutine psb_sexscan_sums(ictxt,dat) #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat real(psb_spk_) :: dat_ integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) call mpi_exscan(dat,dat_,1,psb_mpi_r_spk_,mpi_sum,ictxt,info) dat = dat_ #else dat = szero #endif end subroutine psb_sexscan_sums subroutine psb_sscan_sumv(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_), allocatable :: dat_(:) integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) call psb_realloc(size(dat),dat_,iinfo) dat_ = dat if (iinfo == psb_success_) & & call mpi_scan(dat,dat_,size(dat),psb_mpi_r_spk_,mpi_sum,ictxt,info) #endif end subroutine psb_sscan_sumv subroutine psb_sexscan_sumv(ictxt,dat,root) use psb_realloc_mod #ifdef MPI_MOD use mpi #endif implicit none #ifdef MPI_H include 'mpif.h' #endif integer(psb_mpk_), intent(in) :: ictxt real(psb_spk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_spk_), allocatable :: dat_(:) integer(psb_mpk_) :: iam, np, info integer(psb_ipk_) :: iinfo #if !defined(SERIAL_MPI) call psb_info(ictxt,iam,np) call psb_realloc(size(dat),dat_,iinfo) dat_ = dat if (iinfo == psb_success_) & & call mpi_exscan(dat,dat_,size(dat),psb_mpi_r_spk_,mpi_sum,ictxt,info) #else dat = szero #endif end subroutine psb_sexscan_sumv end module psi_s_collective_mod