TABLE OF CONTENTS


ABINIT/m_xredistribute [ Modules ]

[ Top ] [ Modules ]

NAME

  m_xredistribute

FUNCTION

  This module contains a function to re-dristribute
  results on different procs.

COPYRIGHT

 Copyright (C) 2008-2024 ABINIT group (MMANCINI)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .

INPUTS

OUTPUT

SOURCE

21 #if defined HAVE_CONFIG_H
22 #include "config.h"
23 #endif
24 
25 #include "abi_common.h"
26 
27 
28 module m_xredistribute
29 
30  use defs_basis
31  use m_errors
32  use m_abicore
33 #if defined HAVE_MPI2
34  use mpi
35 #endif
36 
37  implicit none
38 
39 #if defined HAVE_MPI1
40  include 'mpif.h'
41 #endif
42 
43  private
44 
45  public :: xredistribute         ! Redistribute the work
46 
47 
48  interface xredistribute
49   module procedure xredistribute_mpi_dp
50   module procedure xredistribute_mpi_2d_dp
51  end interface xredistribute
52 
53 
54 CONTAINS  !===========================================================

ABINIT/xredistribute_mpi_2d_dp [ Functions ]

[ Top ] [ Functions ]

NAME

  xredistribute_mpi_2d_dp

FUNCTION

INPUTS

  xval= buffer array
  send_counts= number of sent elements (initial distribution)
  send_displs= postions of values sent by the processor (initial positions)
  rec_counts= number of received elements (final distribution)
  rec_displs= positions of values received by the processors (final position)
  nproc=number of processor
  me=proc me
  spaceComm= MPI communicator

OUTPUT

  ier= exit status, a non-zero value meaning there is an error

SIDE EFFECTS

  recvbuf= received buffer

SOURCE

148 subroutine xredistribute_mpi_2d_dp(xval,send_counts,send_displs,recvbuf,&
149                                    rec_counts,rec_displs,me,nproc,spaceComm,ier)
150 
151 !Arguments-------------------------
152  integer ,intent(in) :: me,nproc
153  real(dp),intent(in) :: xval(:,:)
154  real(dp),intent(inout) :: recvbuf(:,:)
155  integer ,intent(in) :: send_displs(0:nproc-1),send_counts(0:nproc-1)
156  integer ,intent(in) :: rec_displs(0:nproc-1),rec_counts(0:nproc-1)
157  integer ,intent(in) :: spaceComm
158  integer ,intent(out):: ier
159 
160 !Local variables-------------------
161  integer :: size
162  real(dp),allocatable :: totbuff(:)
163  character(500) :: msg
164 ! *********************************************************************
165 
166  ier=0
167 #if defined HAVE_MPI
168  if (spaceComm /= MPI_COMM_SELF .and. spaceComm /= MPI_COMM_NULL) then
169  size = sum(send_counts)
170  if(size /=sum(rec_counts))then
171    msg = 'the total sizes of sent and receved msg are not equal'
172    ABI_ERROR(msg)
173  endif
174 
175  ABI_MALLOC(totbuff,(size))
176    !--join all the vector in to a single one
177    call MPI_ALLGATHERV(xval,send_counts(me),MPI_DOUBLE_PRECISION,totbuff,&
178 &                      send_counts,send_displs,MPI_DOUBLE_PRECISION,spaceComm,ier)
179 
180 
181    !--now distribute the total vector on the procs
182    call MPI_SCATTERV(totbuff,rec_counts,rec_displs,MPI_DOUBLE_PRECISION,&
183 &                    recvbuf,rec_counts(me),MPI_DOUBLE_PRECISION,&
184 &                    0,spaceComm,ier)
185    ABI_FREE(totbuff)
186  end if
187 #endif
188 end subroutine xredistribute_mpi_2d_dp

ABINIT/xredistribute_mpi_dp [ Functions ]

[ Top ] [ Functions ]

NAME

  xredistribute_mpi_dp

FUNCTION

INPUTS

  xval= buffer array
  send_counts= number of sent elements (initial distribution)
  send_displs= postions of values sent by the processor (initial positions)
  rec_counts= number of received elements (final distribution)
  rec_displs= positions of values received by the processors (final position)
  nproc=number of processor
  me=proc me
  spaceComm= MPI communicator

OUTPUT

  ier= exit status, a non-zero value meaning there is an error

SIDE EFFECTS

  recvbuf= received buffer

SOURCE

 80 subroutine xredistribute_mpi_dp(xval,send_counts,send_displs,recvbuf,&
 81                                 rec_counts,rec_displs,me,nproc,spaceComm,ier)
 82 
 83 !Arguments-------------------------
 84  integer ,intent(in) :: me,nproc
 85  real(dp),intent(in) :: xval(:)
 86  real(dp),intent(inout) :: recvbuf(:)
 87  integer ,intent(in) :: send_displs(0:nproc-1),send_counts(0:nproc-1)
 88  integer ,intent(in) :: rec_displs(0:nproc-1),rec_counts(0:nproc-1)
 89  integer ,intent(in) :: spaceComm
 90  integer ,intent(out):: ier
 91 
 92 !Local variables-------------------
 93  integer :: size
 94  real(dp),allocatable :: totbuff(:)
 95  character(500) :: msg
 96 ! *********************************************************************
 97 
 98  ier=0
 99 #if defined HAVE_MPI
100  if (spaceComm /= MPI_COMM_SELF .and. spaceComm /= MPI_COMM_NULL) then
101  size = sum(send_counts)
102  if(size /=sum(rec_counts))then
103    msg = 'the total sizes of sent and receved msg are not equal'
104    ABI_ERROR(msg)
105  endif
106 
107  ABI_MALLOC(totbuff,(size))
108    !--join all the vector in to a single one
109    call MPI_ALLGATHERV(xval,send_counts(me),MPI_DOUBLE_PRECISION,totbuff,&
110 &                      send_counts,send_displs,MPI_DOUBLE_PRECISION,spaceComm,ier)
111 
112 
113    !--now distribute the total vector on the procs
114    call MPI_SCATTERV(totbuff,rec_counts,rec_displs,MPI_DOUBLE_PRECISION,&
115 &                    recvbuf,rec_counts(me),MPI_DOUBLE_PRECISION,&
116 &                    0,spaceComm,ier)
117 
118    ABI_FREE(totbuff)
119  end if
120 #endif
121 end subroutine xredistribute_mpi_dp