! ! 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_d_collective_mod use psi_penv_mod interface psb_max module procedure psb_dmaxs, psb_dmaxv, psb_dmaxm, & & psb_dmaxs_ec, psb_dmaxv_ec, psb_dmaxm_ec end interface interface psb_min module procedure psb_dmins, psb_dminv, psb_dminm, & & psb_dmins_ec, psb_dminv_ec, psb_dminm_ec end interface psb_min interface psb_nrm2 module procedure psb_d_nrm2s, psb_d_nrm2v, & & psb_d_nrm2s_ec, psb_d_nrm2v_ec end interface psb_nrm2 interface psb_sum module procedure psb_dsums, psb_dsumv, psb_dsumm, & & psb_dsums_ec, psb_dsumv_ec, psb_dsumm_ec end interface interface psb_amx module procedure psb_damxs, psb_damxv, psb_damxm, & & psb_damxs_ec, psb_damxv_ec, psb_damxm_ec end interface interface psb_amn module procedure psb_damns, psb_damnv, psb_damnm, & & psb_damns_ec, psb_damnv_ec, psb_damnm_ec end interface interface psb_bcast module procedure psb_dbcasts, psb_dbcastv, psb_dbcastm, & & psb_dbcasts_ec, psb_dbcastv_ec, psb_dbcastm_ec end interface psb_bcast interface psb_scan_sum module procedure psb_dscan_sums, psb_dscan_sumv end interface psb_scan_sum interface psb_exscan_sum module procedure psb_dexscan_sums, psb_dexscan_sumv end interface psb_exscan_sum contains ! !!!!!!!!!!!!!!!!!!!!!! ! ! Reduction operations ! ! !!!!!!!!!!!!!!!!!!!!!! ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!! ! ! MAX ! ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!! subroutine psb_dmaxs(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_dpk_), intent(inout) :: dat integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_) :: 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_dpk_,mpi_max,ictxt,info) dat = dat_ else call mpi_reduce(dat,dat_,1,psb_mpi_r_dpk_,mpi_max,root_,ictxt,info) if (iam == root_) dat = dat_ endif #endif end subroutine psb_dmaxs subroutine psb_dmaxv(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_dpk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_), 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_dpk_,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_dpk_,mpi_max,root_,ictxt,info) else call psb_realloc(1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_dpk_,mpi_max,root_,ictxt,info) end if endif #endif end subroutine psb_dmaxv subroutine psb_dmaxm(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_dpk_), intent(inout) :: dat(:,:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_), 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_dpk_,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_dpk_,mpi_max,root_,ictxt,info) else call psb_realloc(1,1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_dpk_,mpi_max,root_,ictxt,info) end if endif #endif end subroutine psb_dmaxm subroutine psb_dmaxs_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_dmaxs_ec subroutine psb_dmaxv_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_dmaxv_ec subroutine psb_dmaxm_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_dmaxm_ec ! ! MIN: Minimum Value ! subroutine psb_dmins(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_dpk_), intent(inout) :: dat integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_) :: 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_dpk_,mpi_min,ictxt,info) dat = dat_ else call mpi_reduce(dat,dat_,1,psb_mpi_r_dpk_,mpi_min,root_,ictxt,info) if (iam == root_) dat = dat_ endif #endif end subroutine psb_dmins subroutine psb_dminv(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_dpk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_), 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_dpk_,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_dpk_,mpi_min,root_,ictxt,info) else call psb_realloc(1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_dpk_,mpi_min,root_,ictxt,info) end if endif #endif end subroutine psb_dminv subroutine psb_dminm(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_dpk_), intent(inout) :: dat(:,:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_), 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_dpk_,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_dpk_,mpi_min,root_,ictxt,info) else call psb_realloc(1,1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_dpk_,mpi_min,root_,ictxt,info) end if endif #endif end subroutine psb_dminm subroutine psb_dmins_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_dmins_ec subroutine psb_dminv_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_dminv_ec subroutine psb_dminm_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_dminm_ec ! !!!!!!!!!!!! ! ! Norm 2, only for reals ! ! !!!!!!!!!!!! subroutine psb_d_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_dpk_), intent(inout) :: dat integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_) :: 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_dpk_,mpi_dnrm2_op,ictxt,info) dat = dat_ else call mpi_reduce(dat,dat_,1,psb_mpi_r_dpk_,mpi_dnrm2_op,root_,ictxt,info) if (iam == root_) dat = dat_ endif #endif end subroutine psb_d_nrm2s subroutine psb_d_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_dpk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_), 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_dpk_,& & mpi_dnrm2_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_dpk_,& & mpi_dnrm2_op,root_,ictxt,info) else call psb_realloc(1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_dpk_,& & mpi_dnrm2_op,root_,ictxt,info) end if endif #endif end subroutine psb_d_nrm2v subroutine psb_d_nrm2s_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_d_nrm2s_ec subroutine psb_d_nrm2v_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_d_nrm2v_ec ! ! SUM ! subroutine psb_dsums(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_dpk_), intent(inout) :: dat integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_) :: 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_dpk_,mpi_sum,ictxt,info) dat = dat_ else call mpi_reduce(dat,dat_,1,psb_mpi_r_dpk_,mpi_sum,root_,ictxt,info) if (iam == root_) dat = dat_ endif #endif end subroutine psb_dsums subroutine psb_dsumv(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_dpk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_), 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_dpk_,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_dpk_,mpi_sum,root_,ictxt,info) else call psb_realloc(1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_dpk_,mpi_sum,root_,ictxt,info) end if endif #endif end subroutine psb_dsumv subroutine psb_dsumm(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_dpk_), intent(inout) :: dat(:,:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_), 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_dpk_,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_dpk_,mpi_sum,root_,ictxt,info) else call psb_realloc(1,1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_dpk_,mpi_sum,root_,ictxt,info) end if endif #endif end subroutine psb_dsumm subroutine psb_dsums_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_dsums_ec subroutine psb_dsumv_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_dsumv_ec subroutine psb_dsumm_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_dsumm_ec ! ! AMX: Maximum Absolute Value ! subroutine psb_damxs(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_dpk_), intent(inout) :: dat integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_) :: 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_dpk_,mpi_damx_op,ictxt,info) dat = dat_ else call mpi_reduce(dat,dat_,1,psb_mpi_r_dpk_,mpi_damx_op,root_,ictxt,info) if (iam == root_) dat = dat_ endif #endif end subroutine psb_damxs subroutine psb_damxv(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_dpk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_), 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_dpk_,mpi_damx_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_dpk_,mpi_damx_op,root_,ictxt,info) else call psb_realloc(1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_dpk_,mpi_damx_op,root_,ictxt,info) end if endif #endif end subroutine psb_damxv subroutine psb_damxm(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_dpk_), intent(inout) :: dat(:,:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_), 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_dpk_,mpi_damx_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_dpk_,mpi_damx_op,root_,ictxt,info) else call psb_realloc(1,1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_dpk_,mpi_damx_op,root_,ictxt,info) end if endif #endif end subroutine psb_damxm subroutine psb_damxs_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_damxs_ec subroutine psb_damxv_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_damxv_ec subroutine psb_damxm_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_damxm_ec ! ! AMN: Minimum Absolute Value ! subroutine psb_damns(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_dpk_), intent(inout) :: dat integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_) :: 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_dpk_,mpi_damn_op,ictxt,info) dat = dat_ else call mpi_reduce(dat,dat_,1,psb_mpi_r_dpk_,mpi_damn_op,root_,ictxt,info) if (iam == root_) dat = dat_ endif #endif end subroutine psb_damns subroutine psb_damnv(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_dpk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_), 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_dpk_,mpi_damn_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_dpk_,mpi_damn_op,root_,ictxt,info) else call psb_realloc(1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_dpk_,mpi_damn_op,root_,ictxt,info) end if endif #endif end subroutine psb_damnv subroutine psb_damnm(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_dpk_), intent(inout) :: dat(:,:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_), 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_dpk_,mpi_damn_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_dpk_,mpi_damn_op,root_,ictxt,info) else call psb_realloc(1,1,dat_,iinfo) call mpi_reduce(dat,dat_,size(dat),psb_mpi_r_dpk_,mpi_damn_op,root_,ictxt,info) end if endif #endif end subroutine psb_damnm subroutine psb_damns_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_damns_ec subroutine psb_damnv_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_damnv_ec subroutine psb_damnm_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_damnm_ec ! ! BCAST Broadcast ! subroutine psb_dbcasts(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_dpk_), 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_dpk_,root_,ictxt,info) #endif end subroutine psb_dbcasts subroutine psb_dbcastv(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_dpk_), 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_dpk_,root_,ictxt,info) #endif end subroutine psb_dbcastv subroutine psb_dbcastm(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_dpk_), 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_dpk_,root_,ictxt,info) #endif end subroutine psb_dbcastm subroutine psb_dbcasts_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_dbcasts_ec subroutine psb_dbcastv_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_dbcastv_ec subroutine psb_dbcastm_ec(ictxt,dat,root) implicit none integer(psb_epk_), intent(in) :: ictxt real(psb_dpk_), 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_dbcastm_ec ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!! ! ! SCAN ! ! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!! subroutine psb_dscan_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_dpk_), intent(inout) :: dat real(psb_dpk_) :: 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_dpk_,mpi_sum,ictxt,info) dat = dat_ #endif end subroutine psb_dscan_sums subroutine psb_dexscan_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_dpk_), intent(inout) :: dat real(psb_dpk_) :: 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_dpk_,mpi_sum,ictxt,info) dat = dat_ #else dat = dzero #endif end subroutine psb_dexscan_sums subroutine psb_dscan_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_dpk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_), 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_dpk_,mpi_sum,ictxt,info) #endif end subroutine psb_dscan_sumv subroutine psb_dexscan_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_dpk_), intent(inout) :: dat(:) integer(psb_mpk_), intent(in), optional :: root integer(psb_mpk_) :: root_ real(psb_dpk_), 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_dpk_,mpi_sum,ictxt,info) #else dat = dzero #endif end subroutine psb_dexscan_sumv end module psi_d_collective_mod