TABLE OF CONTENTS


ABINIT/initmpi_img [ Functions ]

[ Top ] [ Functions ]

NAME

  initmpi_img

FUNCTION

  Initializes the mpi informations for parallelism over images of the cell (npimage>1).

COPYRIGHT

  Copyright (C) 2008-2018 ABINIT group (MT,GG)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .
  For the initials of contributors, see
  ~abinit/doc/developers/contributors.txt.

INPUTS

  dtset <type(dataset_type)>=all input variables in this dataset
  mpi_enreg= informations about MPI parallelization
  option= see below

OUTPUT

  mpi_enreg%my_nimage= number of images of the cell treated by current proc
  ===== if option==1 or option==-1
    mpi_enreg%my_imgtab= indexes of images of the cell treated by current proc
  ===== if option==2 or option==3 or option==-1
    mpi_enreg%comm_cell=Communicator over all processors treating the same image
    mpi_enreg%nproc_cell=size of comm_cell
    mpi_enreg%me_cell=my rank in comm_cell
  ===== if option==3 or option==-1
    mpi_enreg%comm_img=Communicator over all images
    mpi_enreg%nproc_img=size of comm_img
    mpi_enreg%me_img=my rank in comm_img
    mpi_enreg%distrb_img(:)=index of processor treating each image (in comm_img communicator)

PARENTS

      mpi_setup

CHILDREN

      sort_int

SOURCE

 44 #if defined HAVE_CONFIG_H
 45 #include "config.h"
 46 #endif
 47 
 48 #include "abi_common.h"
 49 
 50 subroutine initmpi_img(dtset,mpi_enreg,option)
 51 
 52  use defs_basis
 53  use defs_abitypes
 54  use m_errors
 55  use m_xmpi
 56  use m_profiling_abi
 57  use m_sort
 58 
 59  use m_io_tools,  only: flush_unit
 60 
 61 !This section has been created automatically by the script Abilint (TD).
 62 !Do not modify the following lines by hand.
 63 #undef ABI_FUNC
 64 #define ABI_FUNC 'initmpi_img'
 65 !End of the abilint section
 66 
 67  implicit none
 68 
 69 !Arguments ------------------------------------
 70  integer,intent(in) :: option
 71  type(dataset_type),intent(in) :: dtset
 72  type(MPI_type),intent(inout) :: mpi_enreg
 73 
 74 !Local variables-------------------------------
 75  integer :: imod,irank,iprocmax,iprocmin,jrank
 76  integer :: ndynimage_eff,nimage_eff,nproc_per_image,nrank
 77  logical,parameter :: debug=.false.
 78  integer,allocatable :: ranks(:)
 79  character(len=500) :: msg
 80 
 81 !integer :: group_cell,ierr
 82 
 83 ! ***********************************************************************
 84 
 85  DBG_ENTER("COLL")
 86  if (option/=0) then
 87    mpi_enreg%comm_img=xmpi_comm_self
 88    mpi_enreg%comm_cell=mpi_enreg%comm_world
 89  end if
 90 
 91  if (xmpi_paral==1.and.dtset%npimage>1.and.dtset%optdriver==RUNL_GSTATE) then
 92 
 93 !  Activate flag for parallelization over images
 94    mpi_enreg%paral_img=1
 95 
 96    ndynimage_eff=dtset%ndynimage;if (dtset%ntimimage<=1) ndynimage_eff=0
 97 
 98 !  Print several warnings
 99    if (option==0) then
100      nimage_eff=max(ndynimage_eff,dtset%nimage-ndynimage_eff)
101      if (dtset%npimage>nimage_eff) then
102        write(unit=msg,fmt='(3a,i4,a,i4,4a)') &
103 &       'The number of processors used for the parallelization',ch10,&
104 &       ' over images (npimage=',dtset%npimage,&
105 &       ') is greater than the number of dynamic (or static) images (',nimage_eff,') !',ch10,&
106 &       ' This is unefficient.',ch10
107        MSG_WARNING(msg)
108      end if
109      if (dtset%npimage>mpi_enreg%nproc) then
110        write(unit=msg,fmt='(3a,i6,a,i4,4a)') &
111 &       'The number of processors used for the parallelization',ch10,&
112 &       ' over images (nproc=',mpi_enreg%nproc,&
113 &       ') is smaller than npimage in input file (',dtset%npimage,&
114 &       ')!',ch10,' This is unconsistent.',ch10
115        MSG_ERROR(msg)
116      end if
117      if (mod(nimage_eff,dtset%npimage)/=0) then
118        write(unit=msg,fmt='(3a,i4,a,i4,4a)') &
119 &       'The number of processors used for the parallelization',ch10,&
120 &       ' over images (npimage=',dtset%npimage,&
121 &       ') does not divide the number of dynamic images (',nimage_eff,&
122 &       ') !',ch10,' This is unefficient (charge unbalancing).',ch10
123        MSG_WARNING(msg)
124      end if
125    end if
126 
127 !  # of images treated by current proc
128    nproc_per_image=mpi_enreg%nproc/dtset%npimage
129    iprocmax=nproc_per_image*dtset%npimage-1
130    if (mpi_enreg%me<=iprocmax) then
131      mpi_enreg%my_nimage=(ndynimage_eff/dtset%npimage)+((dtset%nimage-ndynimage_eff)/dtset%npimage)
132      imod=mod(ndynimage_eff,dtset%npimage)-1
133      if (mpi_enreg%me/nproc_per_image<=imod) mpi_enreg%my_nimage=mpi_enreg%my_nimage+1
134      imod=mod((dtset%nimage-ndynimage_eff),dtset%npimage)-1
135      if (mpi_enreg%me/nproc_per_image<=imod) mpi_enreg%my_nimage=mpi_enreg%my_nimage+1
136    else
137      mpi_enreg%my_nimage=0
138    end if
139    if (option==1.or.option==-1) then
140 !    Indexes of images treated by current proc
141      if (mpi_enreg%me<=iprocmax) then
142        ABI_ALLOCATE(mpi_enreg%my_imgtab,(mpi_enreg%my_nimage))
143        nrank=0
144        imod=mpi_enreg%me/nproc_per_image+1;imod=mod(imod,dtset%npimage)
145 !      Dynamic images
146        irank=0
147        do jrank=1,dtset%nimage
148          if (dtset%dynimage(jrank)/=0.and.dtset%ntimimage>1) then
149            irank=irank+1
150            if (mod(irank,dtset%npimage)==imod) then
151              nrank=nrank+1
152              mpi_enreg%my_imgtab(nrank)=jrank
153            end if
154          end if
155        end do
156 !      Static images
157        irank=0
158        do jrank=1,dtset%nimage
159          if (dtset%dynimage(jrank)==0.or.dtset%ntimimage<=1) then
160            irank=irank+1
161            if (mod(irank,dtset%npimage)==imod) then
162              nrank=nrank+1
163              mpi_enreg%my_imgtab(nrank)=jrank
164            end if
165          end if
166        end do
167        if (nrank/=mpi_enreg%my_nimage) then
168          MSG_BUG('Error on nrank !')
169        end if
170 !      Sort images by increasing index (this step is MANDATORY !!)
171        ABI_ALLOCATE(ranks,(nrank))
172        call sort_int(nrank,mpi_enreg%my_imgtab,ranks)
173        ABI_DEALLOCATE(ranks)
174      else
175        ABI_ALLOCATE(mpi_enreg%my_imgtab,(0))
176      end if
177    end if
178    if (option==2.or.option==3.or.option==-1) then
179 !    Communicator over one image
180      if (mpi_enreg%me<=iprocmax) then
181        ABI_ALLOCATE(ranks,(nproc_per_image))
182        iprocmin=(mpi_enreg%me/nproc_per_image)*nproc_per_image
183        ranks=(/((iprocmin+irank-1),irank=1,nproc_per_image)/)
184        mpi_enreg%comm_cell=xmpi_subcomm(mpi_enreg%comm_world,nproc_per_image,ranks)
185        ABI_DEALLOCATE(ranks)
186        mpi_enreg%me_cell=xmpi_comm_rank(mpi_enreg%comm_cell)
187        mpi_enreg%nproc_cell=nproc_per_image
188        if (mpi_enreg%me_cell==0.and.mod(mpi_enreg%me,nproc_per_image)/=0) then
189          MSG_BUG('Error on me_cell !')
190        end if
191      else
192        mpi_enreg%comm_img=xmpi_comm_null
193        mpi_enreg%nproc_cell=0
194        mpi_enreg%me_cell=-1
195      end if
196    end if
197    if (option==3.or.option==-1) then
198 !    Communicator over all images
199      if (mpi_enreg%me<=iprocmax) then
200        ABI_ALLOCATE(ranks,(dtset%npimage))
201        iprocmin=mod(mpi_enreg%me,nproc_per_image)
202        ranks=(/((iprocmin+(irank-1)*nproc_per_image),irank=1,dtset%npimage)/)
203        mpi_enreg%comm_img=xmpi_subcomm(mpi_enreg%comm_world,dtset%npimage,ranks)
204        ABI_DEALLOCATE(ranks)
205        mpi_enreg%me_img=xmpi_comm_rank(mpi_enreg%comm_img)
206        mpi_enreg%nproc_img=dtset%npimage
207        if (iprocmin==0.and.mpi_enreg%me_img==0.and.mpi_enreg%me/=0) then
208          MSG_BUG('Error on me_img!')
209        end if
210        ABI_ALLOCATE(mpi_enreg%distrb_img,(dtset%nimage))
211 !      Dynamic images
212        nrank=0
213        do irank=1,dtset%nimage
214          if (dtset%dynimage(irank)/=0.and.dtset%ntimimage>1) then
215            nrank=nrank+1
216            mpi_enreg%distrb_img(irank)=mod(nrank,dtset%npimage)-1
217            if (mpi_enreg%distrb_img(irank)==-1) mpi_enreg%distrb_img(irank)=dtset%npimage-1
218          end if
219        end do
220 !      Static images
221        nrank=0
222        do irank=1,dtset%nimage
223          if (dtset%dynimage(irank)==0.or.dtset%ntimimage<=1) then
224            nrank=nrank+1
225            mpi_enreg%distrb_img(irank)=mod(nrank,dtset%npimage)-1
226            if (mpi_enreg%distrb_img(irank)==-1) mpi_enreg%distrb_img(irank)=dtset%npimage-1
227          end if
228        end do
229      else
230        mpi_enreg%comm_img=xmpi_comm_null
231        mpi_enreg%nproc_img=0
232        mpi_enreg%me_img=-1
233        ABI_ALLOCATE(mpi_enreg%distrb_img,(0))
234      end if
235    end if
236 
237 !  if (debug) then
238 !  write(200+mpi_enreg%me,*) "=================================="
239 !  write(200+mpi_enreg%me,*) "DEBUGGING STATMENTS IN INITMPI_IMG"
240 !  write(200+mpi_enreg%me,*) "=================================="
241 !  write(200+mpi_enreg%me,*) "option         =",option
242 !  write(200+mpi_enreg%me,*) "MPI_UNDEFINED  =",MPI_UNDEFINED
243 !  write(200+mpi_enreg%me,*) "MPI_IDENT      =",MPI_IDENT
244 !  write(200+mpi_enreg%me,*) "MPI_CONGRUENT  =",MPI_CONGRUENT
245 !  write(200+mpi_enreg%me,*) "MPI_SIMILAR    =",MPI_SIMILAR
246 !  write(200+mpi_enreg%me,*) "MPI_UNEQUAL    =",MPI_UNEQUAL
247 !  write(200+mpi_enreg%me,*) "null_comm      =",MPI_COMM_NULL
248 !  write(200+mpi_enreg%me,*) "self_comm      =",xmpi_comm_self
249 !  write(200+mpi_enreg%me,*) "world_comm     =",mpi_enreg%comm_world
250 !  write(200+mpi_enreg%me,*) "empty_group    =",MPI_GROUP_EMPTY
251 !  write(200+mpi_enreg%me,*) "nimage         =",mpi_enreg%my_nimage
252 !  write(200+mpi_enreg%me,*) "nproc_per_image=",nproc_per_image
253 !  call MPI_COMM_SIZE(mpi_enreg%comm_world,irank,ierr)
254 !  write(200+mpi_enreg%me,*) "Size of world_comm    =",irank
255 !  call MPI_COMM_RANK(mpi_enreg%comm_world,irank,ierr)
256 !  write(200+mpi_enreg%me,*) "My rank in world_comm =",irank
257 !  if (option==1.or.option==-1) then
258 !  write(200+mpi_enreg%me,*) "index_img=",mpi_enreg%my_imgtab(:)
259 !  end if
260 !  if (option==2.or.option==3.or.option==-1) then
261 !  write(200+mpi_enreg%me,*) "nproc_cell  =",mpi_enreg%nproc_cell
262 !  write(200+mpi_enreg%me,*) "me_cell     =",mpi_enreg%me_cell
263 !  call xmpi_comm_group(mpi_enreg%comm_cell,group_cell,ierr)
264 !  write(200+mpi_enreg%me,*) "group_cell  =",group_cell
265 !  write(200+mpi_enreg%me,*) "comm_cell   =",mpi_enreg%comm_cell
266 !  if (group_cell/=MPI_GROUP_EMPTY) then
267 !  call MPI_GROUP_SIZE(group_cell,irank,ierr)
268 !  write(200+mpi_enreg%me,*) "Size of group_cell   =",irank
269 !  call MPI_GROUP_RANK(group_cell,irank,ierr)
270 !  write(200+mpi_enreg%me,*) "My rank in group_cell=",irank
271 !  else
272 !  write(200+mpi_enreg%me,*) "Size of group_cell   =",0
273 !  write(200+mpi_enreg%me,*) "My rank in group_cell=",-1
274 !  end if
275 !  if (mpi_enreg%comm_cell/=MPI_COMM_NULL) then
276 !  call MPI_COMM_SIZE(mpi_enreg%comm_cell,irank,ierr)
277 !  write(200+mpi_enreg%me,*) "Size of comm_cell   =",irank
278 !  call MPI_COMM_RANK(mpi_enreg%comm_cell,irank,ierr)
279 !  write(200+mpi_enreg%me,*) "My rank in comm_cell=",irank
280 !  call MPI_COMM_COMPARE(mpi_enreg%comm_world,mpi_enreg%comm_cell,irank,ierr)
281 !  write(200+mpi_enreg%me,*) "Comparison world_comm/comm_cell=",irank
282 !  call MPI_COMM_COMPARE(xmpi_comm_self,mpi_enreg%comm_cell,irank,ierr)
283 !  write(200+mpi_enreg%me,*) "Comparison self_comm/comm_cell =",irank
284 !  else
285 !  write(200+mpi_enreg%me,*) "Size of comm_cell   =",0
286 !  write(200+mpi_enreg%me,*) "My rank in comm_cell=",-1
287 !  write(200+mpi_enreg%me,*) "Comparison world_comm/comm_cell=",MPI_UNEQUAL
288 !  write(200+mpi_enreg%me,*) "Comparison self_comm/comm_cell =",MPI_UNEQUAL
289 !  end if
290 !  end if
291 !  if (option==3.or.option==-1) then
292 !  write(200+mpi_enreg%me,*) "nproc_img  =",mpi_enreg%nproc_img
293 !  write(200+mpi_enreg%me,*) "me_img     =",mpi_enreg%me_img
294 !  write(200+mpi_enreg%me,*) "img_comm   =",mpi_enreg%comm_img
295 !  if (mpi_enreg%comm_img/=MPI_COMM_NULL) then
296 !  call MPI_COMM_SIZE(mpi_enreg%comm_img,irank,ierr)
297 !  write(200+mpi_enreg%me,*) "Size of img_comm   =",irank
298 !  call MPI_COMM_RANK(mpi_enreg%comm_img,irank,ierr)
299 !  write(200+mpi_enreg%me,*) "My rank in img_comm=",irank
300 !  call MPI_COMM_COMPARE(mpi_enreg%comm_world,mpi_enreg%comm_img,irank,ierr)
301 !  write(200+mpi_enreg%me,*) "Comparison world_comm/img_comm=",irank
302 !  call MPI_COMM_COMPARE(xmpi_comm_self,mpi_enreg%comm_img,irank,ierr)
303 !  write(200+mpi_enreg%me,*) "Comparison self_comm/img_comm =",irank
304 !  else
305 !  write(200+mpi_enreg%me,*) "Size of img_comm   =",0
306 !  write(200+mpi_enreg%me,*) "My rank in img_comm=",-1
307 !  write(200+mpi_enreg%me,*) "Comparison world_comm/img_comm=",MPI_UNEQUAL
308 !  write(200+mpi_enreg%me,*) "Comparison self_comm/img_comm =",MPI_UNEQUAL
309 !  end if
310 !  write(200+mpi_enreg%me,*) "distrb_img=",mpi_enreg%distrb_img(:)
311 !  end if
312 !  write(200+mpi_enreg%me,*)
313 !  call flush_unit(200+mpi_enreg%me)
314 !  if (option==-1) stop
315 !  end if
316 
317  else
318 
319 !  Do not activate flag for parallelization over images
320    mpi_enreg%paral_img=0
321 !  # of images treated by current proc
322    if (dtset%optdriver==RUNL_GSTATE) then
323      mpi_enreg%my_nimage=dtset%nimage
324    else
325      mpi_enreg%my_nimage=1
326    end if
327 !  Indexes of images treated by current proc
328    if (option==1.or.option==-1) then
329      ABI_ALLOCATE(mpi_enreg%my_imgtab,(mpi_enreg%my_nimage))
330      mpi_enreg%my_imgtab=(/(irank,irank=1,mpi_enreg%my_nimage)/)
331    end if
332 !  Communicator over all images
333    if (option==2.or.option==3.or.option==-1) then
334 !    Communicator for one image
335      mpi_enreg%nproc_cell=mpi_enreg%nproc
336      mpi_enreg%me_cell=mpi_enreg%me
337    end if
338    if (option==3.or.option==-1) then
339 !    Communicator over all images
340      mpi_enreg%nproc_img=1
341      mpi_enreg%comm_img=xmpi_comm_self
342      mpi_enreg%me_img=0
343      ABI_ALLOCATE(mpi_enreg%distrb_img,(dtset%nimage))
344      mpi_enreg%distrb_img(:)=0
345    end if
346  end if
347 
348  DBG_EXIT("COLL")
349 
350 end subroutine initmpi_img