D_Gather_BLACS Subroutine

subroutine D_Gather_BLACS(Mloc, M, descM, unit)

Uses

  • proc~~d_gather_blacs~2~~UsesGraph proc~d_gather_blacs~2 D_Gather_BLACS module~sf_mpi SF_MPI proc~d_gather_blacs~2->module~sf_mpi

Arguments

Type IntentOptional Attributes Name
real(kind=8), intent(in), dimension(:,:) :: Mloc
real(kind=8), intent(inout), dimension(:,:) :: M
integer, intent(in), dimension(9) :: descM
integer, optional :: unit

Calls

proc~~d_gather_blacs~2~~CallsGraph proc~d_gather_blacs~2 D_Gather_BLACS blacs_gridinfo blacs_gridinfo proc~d_gather_blacs~2->blacs_gridinfo dgerv2d dgerv2d proc~d_gather_blacs~2->dgerv2d dgesd2d dgesd2d proc~d_gather_blacs~2->dgesd2d infog2l infog2l proc~d_gather_blacs~2->infog2l interface~bcast_mpi Bcast_MPI proc~d_gather_blacs~2->interface~bcast_mpi

Source Code

subroutine D_Gather_BLACS(Mloc,M,descM,unit)
  USE SF_MPI
  include 'mpif.h'
  real(8),dimension(:,:),intent(in)    :: Mloc
  real(8),dimension(:,:),intent(inout) :: M
  integer,dimension(9),intent(in)      :: descM
  integer,optional                     :: unit
  integer                              :: i,j
  integer                              :: myi,myj
  integer                              :: Ni,Nj
  integer                              :: Nbi,Nbj
  integer                              :: Nrow,Ncol,Qrows
  integer                              :: rankX,rankY
  integer                              :: pNx,pNy
  integer                              :: sendR,sendC
  real(8)                              :: t_start,t_stop
  integer                              :: context
  !
  context = descM(2)
  call blacs_gridinfo( context, pNx, pNy, rankX, rankY)
  !
  Ni  = size(M,1)
  Nj  = size(M,2)
  Nbi = descM(5)
  Nbj = descM(6)
  Qrows=descM(9) ; if(Qrows/=size(Mloc,1))stop "Gather_BLACS error: Qrows != size(Mloc,1)"
  !
  if(rankX==0.AND.rankY==0)call cpu_time(t_start)
  do j=1,Nj,Nbj
     Ncol = Nbj ; if(Nj-j<Nbj-1)Ncol=Nj-j+1
     do i=1,Ni,Nbi
        Nrow = Nbi ; if(Ni-i<Nbi-1)Nrow=Ni-i+1
        call infog2l(i,j,descM, pNx, pNy, rankX, rankY, myi, myj, SendR, SendC)
        if(rankX==SendR .AND. rankY==SendC)then
           call dgesd2d(context,Nrow,Ncol,Mloc(myi,myj),Qrows,0,0)
        endif
        if(rankX==0 .AND. rankY==0)then
           call dgerv2d(context,Nrow,Ncol,M(i,j),Ni,SendR,SendC)
        endif
     enddo
  enddo
  if(rankX==0.AND.rankY==0)call cpu_time(t_stop)
  if(present(unit))then
     if(rankX==0.AND.rankY==0)write(unit,"(A20,F21.12)")"Time Gather :",t_stop-t_start
  endif
  !
  if(rankX==0.AND.rankY==0)call cpu_time(t_start)
  call Bcast_MPI(MPI_COMM_WORLD,M)
  if(rankX==0.AND.rankY==0)call cpu_time(t_stop)
  if(present(unit))then
     if(rankX==0.AND.rankY==0)write(unit,"(A20,F21.12)")"Time Bcast :",t_stop-t_start
  endif
  !
  return
end subroutine D_Gather_BLACS