-
Notifications
You must be signed in to change notification settings - Fork 0
/
m_paral_atom.F90
375 lines (321 loc) · 9.7 KB
/
m_paral_atom.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
!!****m* ABINIT/m_paral_atom
!! NAME
!! m_paral_atom
!!
!! FUNCTION
!! This module provides routines and methods used to manage the parallelisation/distribution
!! of PAW data over atomic sites
!!
!! COPYRIGHT
!! Copyright (C) 2012-2022 ABINIT group (MT, MD)
!! 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.
!!
!! NOTES
!! FOR DEVELOPPERS: in order to preserve the portability of libPAW library,
!! please consult ~abinit/src/??_libpaw/libpaw-coding-rules.txt
!!
!! SOURCE
#include "libpaw.h"
MODULE m_paral_atom
USE_DEFS
USE_MSG_HANDLING
USE_MPI_WRAPPERS
USE_MEMORY_PROFILING
implicit none
private
!public procedures.
public :: get_my_natom
public :: get_my_atmtab
public :: free_my_atmtab
public :: get_proc_atmtab
public :: get_atm_proc
!!***
CONTAINS
!!***
!----------------------------------------------------------------------
!!****f* m_paral_atom/get_my_natom
!! NAME
!! get_my_natom
!!
!! FUNCTION
!! Given the total number of atoms, return the number of atoms treated by current process
!!
!! INPUTS
!! comm_atom=communicator over atoms
!! natom=total number of atoms
!!
!! OUTPUT
!! my_natom=number of atoms treated by current process
!!
!! SOURCE
subroutine get_my_natom(comm_atom,my_natom,natom)
!Arguments ---------------------------------------------
!scalars
integer,intent(in) :: comm_atom,natom
integer,intent(out) :: my_natom
!arrays
!Local variables ---------------------------------------
!scalars
integer :: me,nproc
!arrays
! *************************************************************************
my_natom=natom
if (comm_atom/=xpaw_mpi_comm_self.and.comm_atom/=xpaw_mpi_comm_null) then
nproc=xpaw_mpi_comm_size(comm_atom)
me=xpaw_mpi_comm_rank(comm_atom)
my_natom=natom/nproc
if (me<=(mod(natom,nproc)-1)) my_natom=natom/nproc + 1
endif
end subroutine get_my_natom
!!***
!----------------------------------------------------------------------
!!****f* m_paral_atom/get_my_atmtab
!! NAME
!! get_my_atmtab
!!
!! FUNCTION
!! Given the total number of atoms and a MPI communicator return a table
!! containing the indexes of atoms treated by current processor.
!!
!! INPUTS
!! comm_atom=communicator over atoms
!! my_natom_ref= --optional-- a reference value for the local number of atoms
!! (just for checking purposes)
!! natom=total number of atoms
!!
!! OUTPUT
!!
!! SIDE EFFECTS
!! my_atmtab(:)=indexes of atoms treated by current process
!! nothing is done if my_atmtab(:) is already associated to a target
!! my_atmtab_allocated=true if my_atmtab is allocated
!! paral_atom=flag controlling parallelism over atoms
!!
!! SOURCE
subroutine get_my_atmtab(comm_atom,my_atmtab,my_atmtab_allocated,paral_atom,natom,&
& my_natom_ref) ! optional argument
!Arguments ---------------------------------------------
!scalars
integer,intent(in) :: comm_atom,natom
integer,intent(in),optional :: my_natom_ref
logical,intent(inout) :: my_atmtab_allocated,paral_atom
!arrays
integer,pointer :: my_atmtab(:)
!Local variables ---------------------------------------
!scalars
integer :: iatom,me,my_natom,natom_bef,nmod,nproc
character(len=500) :: msg
!arrays
! *************************************************************************
my_atmtab_allocated=.false.
if (.not.paral_atom) return
if (comm_atom==xpaw_mpi_comm_self.or.comm_atom==xpaw_mpi_comm_null) paral_atom=.false.
if (paral_atom) then
nproc=xpaw_mpi_comm_size(comm_atom)
paral_atom=(nproc>1)
if (paral_atom) then
if (.not.associated(my_atmtab)) then
! Get local number of atoms
me=xpaw_mpi_comm_rank(comm_atom)
my_natom=natom/nproc
if (me<=(mod(natom,nproc)-1)) my_natom=natom/nproc + 1
! Get table of atoms
if (my_natom>0) then
LIBPAW_POINTER_ALLOCATE(my_atmtab,(my_natom))
my_atmtab_allocated=.true.
if (my_natom==natom) then
my_atmtab(1:my_natom)=(/(iatom,iatom=1,natom)/)
else
! The atoms are distributed contigously by egal part
! (the rest is distributed on all the procs)
nmod=mod(natom,nproc)
if (me<=(nmod-1)) then
natom_bef=me*(natom/nproc)+me
else
natom_bef=me*(natom/nproc)+nmod
endif
do iatom=1,my_natom
my_atmtab(iatom)=iatom+natom_bef
enddo
end if
end if
else
my_natom=size(my_atmtab)
end if
if (present(my_natom_ref).and.(my_natom>0)) then
if (my_natom_ref/=size(my_atmtab)) then
msg='my_atmtab should have a size equal to my_natom !'
LIBPAW_BUG(msg)
end if
end if
end if
endif
end subroutine get_my_atmtab
!!***
!----------------------------------------------------------------------
!!****f* m_paral_atom/free_my_atmtab
!! NAME
!! free_my_atmtab
!!
!! FUNCTION
!! Cleanly deallocate a table of atom indexes (my_atmtab)
!!
!! INPUTS
!!
!! OUTPUT
!!
!! SIDE EFFECTS
!! my_atmtab_allocated=true if my_atmtab is allocated
!! my_atmtab(:)=indexes of atoms treated by current process
!! nothing is done if my_atmtab(:) is already associated to a target
!!
!! SOURCE
subroutine free_my_atmtab(my_atmtab,my_atmtab_allocated)
!Arguments ---------------------------------------------
!scalars
logical,intent(inout) :: my_atmtab_allocated
!arrays
integer,pointer :: my_atmtab(:)
!Local variables ---------------------------------------
!scalars
!arrays
! *************************************************************************
if (my_atmtab_allocated) then
LIBPAW_POINTER_DEALLOCATE(my_atmtab)
nullify(my_atmtab)
my_atmtab_allocated=.false.
end if
end subroutine free_my_atmtab
!!***
!----------------------------------------------------------------------
!!****f* m_paral_atom/get_proc_atmtab
!! NAME
!! get_proc_atmtab
!!
!! FUNCTION
!! Given the total number of atoms and the size of a communicator,
!! return a table containing the indexes of atoms treated by a processor (with index iproc)
!!
!! INPUTS
!! comm_atom_size= size of communicator (over atoms)
!! iproc= rank of the processor
!! natom= total number of atoms
!!
!! OUTPUT
!! atmtab(natom_out)= indexes of atoms treated by process iproc
!! natom_out= number of atoms treated by process iproc
!!
!! NOTES
!! In case of modification of the distribution of atom over proc,
!! get_atmtab must be modified accordingly
!!
!! SOURCE
subroutine get_proc_atmtab(iproc,atmtab,natom_out,natom,comm_atom_size)
!Arguments ---------------------------------------------
!scalars
integer, intent(in) :: comm_atom_size,iproc,natom
!arrays
integer,intent(out) :: natom_out
integer,allocatable, intent(out):: atmtab(:)
!Local variables ---------------------------------------
!scalars
integer :: iatom,natom_bef,nmod,nproc
!arrays
! *************************************************************************
nproc=comm_atom_size
natom_out=natom/nproc ; if (iproc<=(mod(natom,nproc)-1)) natom_out=natom/nproc+1
! Get table of atoms
if (natom_out>0) then
LIBPAW_ALLOCATE(atmtab,(natom_out))
! The atoms are distributed contigously by egal part
! The rest is distributed on all the procs
! (see get_my_atmtab)
nmod=mod(natom,nproc)
if (iproc<=(nmod-1)) then
natom_bef=iproc*(natom/nproc)+iproc
else
natom_bef=iproc*(natom/nproc)+nmod
end if
do iatom=1,natom_out
atmtab(iatom)=iatom+natom_bef
end do
else
natom_out=0
LIBPAW_ALLOCATE(atmtab,(0))
end if
end subroutine get_proc_atmtab
!!***
!----------------------------------------------------------------------
!!****f* m_paral_atom/get_atm_proc
!! NAME
!! get_atm_proc
!!
!! FUNCTION
!! Given a list of atoms and a MPI communicator size, return a table
!! containing the corresponding processor indexes.
!!
!! COPYRIGHT
!! Copyright (C) 2012-2022 ABINIT group (MD)
!! 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
!! atom_list(:)= index of atoms
!! nproc=size of communicator over atoms
!! natom=total number of atoms
!!
!! OUTPUT
!! proc_list(:) = index of procs
!!
!! NOTES
!! The atoms are distributed contigously by egal part; the rest is distributed
!! on all the procs (see get_my_atmtab).
!! In case of modification of the distribution of atom over proc, this routine
!! must be modified accordingly.
!!
!! SOURCE
subroutine get_atm_proc(atom_list,natom,nproc,proc_list)
!Arguments ---------------------------------------------
!scalars
integer, intent(in) :: natom,nproc
!arrays
integer, intent(in) :: atom_list(:)
integer, intent(out) :: proc_list(:)
!Local variables ---------------------------------------
!scalars
integer :: nb_atom,dn,dn1,iatom,natomlim,iatm,jproclim,nmod
!arrays
! *************************************************************************
nmod=mod(natom,nproc);nb_atom=size(atom_list)
if (nmod==0) then
dn=natom/nproc
do iatm =1, nb_atom
iatom=atom_list(iatm)
proc_list(iatm)=(iatom-1)/dn
end do
else
dn=natom/nproc
dn1=natom/nproc + 1
! Under (jproclim+1), 1 atome by proc is added
! The rest nmod is distributed among jproclim+1 first procs
jproclim=nmod -1
natomlim=dn1*(jproclim+1)
do iatm=1,nb_atom
iatom=atom_list(iatm)
if (iatom<=natomlim) then
proc_list(iatm)=(iatom -1 )/dn1
else
proc_list(iatm)=jproclim + 1 + (iatom - 1 -(natomlim))/dn
end if
enddo
end if
end subroutine get_atm_proc
!!***
!----------------------------------------------------------------------
END MODULE m_paral_atom
!!***