1 |
!! |
2 |
!! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
!! |
4 |
!! The University of Notre Dame grants you ("Licensee") a |
5 |
!! non-exclusive, royalty free, license to use, modify and |
6 |
!! redistribute this software in source and binary code form, provided |
7 |
!! that the following conditions are met: |
8 |
!! |
9 |
!! 1. Acknowledgement of the program authors must be made in any |
10 |
!! publication of scientific results based in part on use of the |
11 |
!! program. An acceptable form of acknowledgement is citation of |
12 |
!! the article in which the program was described (Matthew |
13 |
!! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
!! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
!! Parallel Simulation Engine for Molecular Dynamics," |
16 |
!! J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
!! |
18 |
!! 2. Redistributions of source code must retain the above copyright |
19 |
!! notice, this list of conditions and the following disclaimer. |
20 |
!! |
21 |
!! 3. Redistributions in binary form must reproduce the above copyright |
22 |
!! notice, this list of conditions and the following disclaimer in the |
23 |
!! documentation and/or other materials provided with the |
24 |
!! distribution. |
25 |
!! |
26 |
!! This software is provided "AS IS," without a warranty of any |
27 |
!! kind. All express or implied conditions, representations and |
28 |
!! warranties, including any implied warranty of merchantability, |
29 |
!! fitness for a particular purpose or non-infringement, are hereby |
30 |
!! excluded. The University of Notre Dame and its licensors shall not |
31 |
!! be liable for any damages suffered by licensee as a result of |
32 |
!! using, modifying or distributing the software or its |
33 |
!! derivatives. In no event will the University of Notre Dame or its |
34 |
!! licensors be liable for any lost revenue, profit or data, or for |
35 |
!! direct, indirect, special, consequential, incidental or punitive |
36 |
!! damages, however caused and regardless of the theory of liability, |
37 |
!! arising out of the use of or inability to use software, even if the |
38 |
!! University of Notre Dame has been advised of the possibility of |
39 |
!! such damages. |
40 |
!! |
41 |
|
42 |
|
43 |
!! MPI support for long range forces using force decomposition |
44 |
!! on a square grid of processors. |
45 |
!! Corresponds to mpiSimunation.cpp for C++ |
46 |
!! mpi_module also contains a private interface for mpi f90 routines. |
47 |
!! |
48 |
!! @author Charles F. Vardeman II |
49 |
!! @author Matthew Meineke |
50 |
!! @version $Id: simParallel.F90,v 1.5 2005-09-07 22:23:20 chuckv Exp $, $Date: 2005-09-07 22:23:20 $, $Name: not supported by cvs2svn $, $Revision: 1.5 $ |
51 |
|
52 |
module mpiSimulation |
53 |
use definitions |
54 |
#ifdef IS_MPI |
55 |
use oopseMPI |
56 |
implicit none |
57 |
PRIVATE |
58 |
|
59 |
|
60 |
!! PUBLIC Subroutines contained in this module |
61 |
!! gather and scatter are a generic interface |
62 |
!! to gather and scatter routines |
63 |
public :: gather, scatter |
64 |
public :: setupSimParallel |
65 |
public :: replanSimParallel |
66 |
public :: getNatomsInCol |
67 |
public :: getNatomsInRow |
68 |
public :: getNgroupsInCol |
69 |
public :: getNgroupsInRow |
70 |
public :: isMPISimSet |
71 |
public :: printComponentPlan |
72 |
public :: getMyNode |
73 |
|
74 |
!! PUBLIC Subroutines contained in MPI module |
75 |
public :: mpi_bcast |
76 |
public :: mpi_allreduce |
77 |
! public :: mpi_reduce |
78 |
public :: mpi_send |
79 |
public :: mpi_recv |
80 |
public :: mpi_get_processor_name |
81 |
public :: mpi_finalize |
82 |
|
83 |
!! PUBLIC mpi variables |
84 |
public :: mpi_comm_world |
85 |
public :: mpi_character |
86 |
public :: mpi_integer |
87 |
public :: mpi_lor |
88 |
public :: mpi_logical |
89 |
public :: mpi_double_precision |
90 |
public :: mpi_sum |
91 |
public :: mpi_max |
92 |
public :: mpi_status_size |
93 |
public :: mpi_any_source |
94 |
|
95 |
|
96 |
|
97 |
!! Safety logical to prevent access to ComponetPlan until |
98 |
!! set by C++. |
99 |
logical, save :: ComponentPlanSet = .false. |
100 |
|
101 |
|
102 |
!! generic mpi error declaration. |
103 |
integer, public :: mpi_err |
104 |
|
105 |
#ifdef PROFILE |
106 |
public :: printCommTime |
107 |
public :: getCommTime |
108 |
real,save :: commTime = 0.0 |
109 |
real :: commTimeInitial,commTimeFinal |
110 |
#endif |
111 |
|
112 |
!! Include mpiComponentPlan type. mpiComponentPlan is a |
113 |
!! dual header file for both c and fortran. |
114 |
#define __FORTRAN90 |
115 |
#include "UseTheForce/mpiComponentPlan.h" |
116 |
|
117 |
|
118 |
!! Tags used during force loop for parallel simulation |
119 |
integer, public, allocatable, dimension(:) :: AtomLocalToGlobal |
120 |
integer, public, allocatable, dimension(:) :: AtomRowToGlobal |
121 |
integer, public, allocatable, dimension(:) :: AtomColToGlobal |
122 |
integer, public, allocatable, dimension(:) :: GroupLocalToGlobal |
123 |
integer, public, allocatable, dimension(:) :: GroupRowToGlobal |
124 |
integer, public, allocatable, dimension(:) :: GroupColToGlobal |
125 |
|
126 |
!! Logical set true if mpiSimulation has been initialized |
127 |
logical, save :: isSimSet = .false. |
128 |
|
129 |
|
130 |
type (mpiComponentPlan), save :: mpiSim |
131 |
|
132 |
!! gs_plan contains plans for gather and scatter routines |
133 |
type, public :: gs_plan |
134 |
private |
135 |
type (mpiComponentPlan), pointer :: gsComponentPlan => NULL() |
136 |
integer :: gsPlanSize !! size of this plan (nDim*nComponents) |
137 |
integer :: globalPlanSize !! size of all components in plan |
138 |
integer, dimension(:), pointer :: displs !! Displacements array for mpi indexed from 0. |
139 |
integer, dimension(:), pointer :: counts !! Counts array for mpi indexed from 0. |
140 |
integer :: myPlanComm !! My communicator for this plan |
141 |
integer :: myPlanRank !! My rank in this plan |
142 |
integer :: planNprocs !! Number of processors in this plan |
143 |
end type gs_plan |
144 |
|
145 |
! plans for different decompositions |
146 |
type (gs_plan), public, save :: plan_atom_row |
147 |
type (gs_plan), public, save :: plan_atom_row_3d |
148 |
type (gs_plan), public, save :: plan_atom_col |
149 |
type (gs_plan), public, save :: plan_atom_col_3d |
150 |
type (gs_plan), public, save :: plan_atom_row_Rotation |
151 |
type (gs_plan), public, save :: plan_atom_col_Rotation |
152 |
type (gs_plan), public, save :: plan_group_row |
153 |
type (gs_plan), public, save :: plan_group_col |
154 |
type (gs_plan), public, save :: plan_group_row_3d |
155 |
type (gs_plan), public, save :: plan_group_col_3d |
156 |
|
157 |
type (mpiComponentPlan), pointer :: simComponentPlan |
158 |
|
159 |
! interface for gather scatter routines |
160 |
!! Generic interface for gather. |
161 |
!! Gathers an local array into row or column array |
162 |
!! Interface provided for integer, double and double |
163 |
!! rank 2 arrays. |
164 |
interface gather |
165 |
module procedure gather_integer |
166 |
module procedure gather_double |
167 |
module procedure gather_double_2d |
168 |
end interface |
169 |
|
170 |
!! Generic interface for scatter. |
171 |
!! Scatters a row or column array, adding componets |
172 |
!! and reducing them to a local nComponent array. |
173 |
!! Interface provided for double and double rank=2 arrays. |
174 |
|
175 |
interface scatter |
176 |
module procedure scatter_double |
177 |
module procedure scatter_double_2d |
178 |
end interface |
179 |
|
180 |
|
181 |
|
182 |
contains |
183 |
|
184 |
!! Sets up mpiComponentPlan with structure passed from C++. |
185 |
subroutine setupSimParallel(thisComponentPlan, nAtomTags, atomTags, & |
186 |
nGroupTags, groupTags, status) |
187 |
!! Passed Arguments |
188 |
!! mpiComponentPlan struct from C |
189 |
type (mpiComponentPlan), intent(inout) :: thisComponentPlan |
190 |
!! Number of tags passed |
191 |
integer, intent(in) :: nAtomTags, nGroupTags |
192 |
!! Result status, 0 = normal, -1 = error |
193 |
integer, intent(out) :: status |
194 |
integer :: localStatus |
195 |
!! Global reference tag for local particles |
196 |
integer, dimension(nAtomTags), intent(inout) :: atomTags |
197 |
integer, dimension(nGroupTags), intent(inout) :: groupTags |
198 |
|
199 |
!write(*,*) 'mpiSim_mod thinks node', thisComponentPlan%myNode, & |
200 |
! ' has atomTags(1) = ', atomTags(1) |
201 |
|
202 |
status = 0 |
203 |
if (componentPlanSet) then |
204 |
return |
205 |
endif |
206 |
componentPlanSet = .true. |
207 |
|
208 |
!! copy c component plan to fortran |
209 |
mpiSim = thisComponentPlan |
210 |
!write(*,*) "Seting up simParallel" |
211 |
|
212 |
call make_Force_Grid(mpiSim, localStatus) |
213 |
if (localStatus /= 0) then |
214 |
write(default_error,*) "Error creating force grid" |
215 |
status = -1 |
216 |
return |
217 |
endif |
218 |
|
219 |
call updateGridComponents(mpiSim, localStatus) |
220 |
if (localStatus /= 0) then |
221 |
write(default_error,*) "Error updating grid components" |
222 |
status = -1 |
223 |
return |
224 |
endif |
225 |
|
226 |
!! initialize gather and scatter plans used in this simulation |
227 |
call plan_gather_scatter(1, mpiSim%nAtomsLocal, & |
228 |
mpiSim, mpiSim%rowComm, plan_atom_row) |
229 |
call plan_gather_scatter(nDim, mpiSim%nAtomsLocal, & |
230 |
mpiSim, mpiSim%rowComm, plan_atom_row_3d) |
231 |
call plan_gather_scatter(9, mpiSim%nAtomsLocal, & |
232 |
mpiSim, mpiSim%rowComm, plan_atom_row_Rotation) |
233 |
call plan_gather_scatter(1, mpiSim%nGroupsLocal, & |
234 |
mpiSim, mpiSim%rowComm, plan_group_row) |
235 |
call plan_gather_scatter(nDim, mpiSim%nGroupsLocal, & |
236 |
mpiSim, mpiSim%rowComm, plan_group_row_3d) |
237 |
|
238 |
call plan_gather_scatter(1, mpiSim%nAtomsLocal, & |
239 |
mpiSim, mpiSim%columnComm, plan_atom_col) |
240 |
call plan_gather_scatter(nDim, mpiSim%nAtomsLocal, & |
241 |
mpiSim, mpiSim%columnComm, plan_atom_col_3d) |
242 |
call plan_gather_scatter(9, mpiSim%nAtomsLocal, & |
243 |
mpiSim, mpiSim%columnComm, plan_atom_col_Rotation) |
244 |
call plan_gather_scatter(1, mpiSim%nGroupsLocal, & |
245 |
mpiSim, mpiSim%columnComm, plan_group_col) |
246 |
call plan_gather_scatter(nDim, mpiSim%nGroupsLocal, & |
247 |
mpiSim, mpiSim%columnComm, plan_group_col_3d) |
248 |
|
249 |
! Initialize tags |
250 |
|
251 |
call setAtomTags(atomTags,localStatus) |
252 |
if (localStatus /= 0) then |
253 |
status = -1 |
254 |
return |
255 |
endif |
256 |
|
257 |
|
258 |
call setGroupTags(groupTags,localStatus) |
259 |
if (localStatus /= 0) then |
260 |
status = -1 |
261 |
return |
262 |
endif |
263 |
|
264 |
isSimSet = .true. |
265 |
|
266 |
! call printComponentPlan(mpiSim,0) |
267 |
end subroutine setupSimParallel |
268 |
|
269 |
subroutine replanSimParallel(thisComponentPlan,status) |
270 |
! Passed Arguments |
271 |
!! mpiComponentPlan struct from C |
272 |
type (mpiComponentPlan), intent(inout) :: thisComponentPlan |
273 |
integer, intent(out) :: status |
274 |
integer :: localStatus |
275 |
integer :: mpierror |
276 |
status = 0 |
277 |
|
278 |
call updateGridComponents(thisComponentPlan,localStatus) |
279 |
if (localStatus /= 0) then |
280 |
status = -1 |
281 |
return |
282 |
endif |
283 |
|
284 |
!! Unplan Gather Scatter plans |
285 |
call unplan_gather_scatter(plan_atom_row) |
286 |
call unplan_gather_scatter(plan_atom_row_3d) |
287 |
call unplan_gather_scatter(plan_atom_row_Rotation) |
288 |
call unplan_gather_scatter(plan_group_row) |
289 |
call unplan_gather_scatter(plan_group_row_3d) |
290 |
|
291 |
call unplan_gather_scatter(plan_atom_col) |
292 |
call unplan_gather_scatter(plan_atom_col_3d) |
293 |
call unplan_gather_scatter(plan_atom_col_Rotation) |
294 |
call unplan_gather_scatter(plan_group_col) |
295 |
call unplan_gather_scatter(plan_group_col_3d) |
296 |
|
297 |
!! initialize gather and scatter plans used in this simulation |
298 |
call plan_gather_scatter(1, mpiSim%nAtomsLocal, & |
299 |
mpiSim, mpiSim%rowComm, plan_atom_row) |
300 |
call plan_gather_scatter(nDim, mpiSim%nAtomsLocal, & |
301 |
mpiSim, mpiSim%rowComm, plan_atom_row_3d) |
302 |
call plan_gather_scatter(9, mpiSim%nAtomsLocal, & |
303 |
mpiSim, mpiSim%rowComm, plan_atom_row_Rotation) |
304 |
call plan_gather_scatter(1, mpiSim%nGroupsLocal, & |
305 |
mpiSim, mpiSim%rowComm, plan_group_row) |
306 |
call plan_gather_scatter(nDim, mpiSim%nGroupsLocal, & |
307 |
mpiSim, mpiSim%rowComm, plan_group_row_3d) |
308 |
|
309 |
call plan_gather_scatter(1, mpiSim%nAtomsLocal, & |
310 |
mpiSim, mpiSim%columnComm, plan_atom_col) |
311 |
call plan_gather_scatter(nDim, mpiSim%nAtomsLocal, & |
312 |
mpiSim, mpiSim%columnComm, plan_atom_col_3d) |
313 |
call plan_gather_scatter(9, mpiSim%nAtomsLocal, & |
314 |
mpiSim, mpiSim%columnComm, plan_atom_col_Rotation) |
315 |
call plan_gather_scatter(1, mpiSim%nGroupsLocal, & |
316 |
mpiSim, mpiSim%columnComm, plan_group_col) |
317 |
call plan_gather_scatter(nDim, mpiSim%nGroupsLocal, & |
318 |
mpiSim, mpiSim%columnComm, plan_group_col_3d) |
319 |
|
320 |
end subroutine replanSimParallel |
321 |
|
322 |
!! Updates number of row and column components for long range forces. |
323 |
subroutine updateGridComponents(thisComponentPlan, status) |
324 |
type (mpiComponentPlan) :: thisComponentPlan !! mpiComponentPlan |
325 |
|
326 |
!! Status return |
327 |
!! - 0 Success |
328 |
!! - -1 Failure |
329 |
integer, intent(out) :: status |
330 |
integer :: nAtomsLocal |
331 |
integer :: nAtomsInRow = 0 |
332 |
integer :: nAtomsInColumn = 0 |
333 |
integer :: nGroupsLocal |
334 |
integer :: nGroupsInRow = 0 |
335 |
integer :: nGroupsInColumn = 0 |
336 |
integer :: mpiErrors |
337 |
|
338 |
status = 0 |
339 |
if (.not. componentPlanSet) return |
340 |
|
341 |
if (thisComponentPlan%nAtomsLocal == 0) then |
342 |
status = -1 |
343 |
return |
344 |
endif |
345 |
if (thisComponentPlan%nGroupsLocal == 0) then |
346 |
write(*,*) 'tcp%ngl = ', thisComponentPlan%nGroupsLocal |
347 |
status = -1 |
348 |
return |
349 |
endif |
350 |
|
351 |
nAtomsLocal = thisComponentPlan%nAtomsLocal |
352 |
nGroupsLocal = thisComponentPlan%nGroupsLocal |
353 |
|
354 |
call mpi_allreduce(nAtomsLocal, nAtomsInRow, 1, mpi_integer, & |
355 |
mpi_sum, thisComponentPlan%rowComm, mpiErrors) |
356 |
if (mpiErrors /= 0) then |
357 |
status = -1 |
358 |
return |
359 |
endif |
360 |
|
361 |
call mpi_allreduce(nAtomsLocal, nAtomsInColumn, 1, mpi_integer, & |
362 |
mpi_sum, thisComponentPlan%columnComm, mpiErrors) |
363 |
if (mpiErrors /= 0) then |
364 |
status = -1 |
365 |
return |
366 |
endif |
367 |
|
368 |
call mpi_allreduce(nGroupsLocal, nGroupsInRow, 1, mpi_integer, & |
369 |
mpi_sum, thisComponentPlan%rowComm, mpiErrors) |
370 |
if (mpiErrors /= 0) then |
371 |
status = -1 |
372 |
return |
373 |
endif |
374 |
|
375 |
call mpi_allreduce(nGroupsLocal, nGroupsInColumn, 1, mpi_integer, & |
376 |
mpi_sum, thisComponentPlan%columnComm, mpiErrors) |
377 |
if (mpiErrors /= 0) then |
378 |
status = -1 |
379 |
return |
380 |
endif |
381 |
|
382 |
thisComponentPlan%nAtomsInRow = nAtomsInRow |
383 |
thisComponentPlan%nAtomsInColumn = nAtomsInColumn |
384 |
thisComponentPlan%nGroupsInRow = nGroupsInRow |
385 |
thisComponentPlan%nGroupsInColumn = nGroupsInColumn |
386 |
|
387 |
end subroutine updateGridComponents |
388 |
|
389 |
|
390 |
!! Creates a square force decomposition of processors into row and column |
391 |
!! communicators. |
392 |
subroutine make_Force_Grid(thisComponentPlan,status) |
393 |
type (mpiComponentPlan) :: thisComponentPlan |
394 |
integer, intent(out) :: status !! status returns -1 if error |
395 |
integer :: nColumnsMax !! Maximum number of columns |
396 |
integer :: nWorldProcessors !! Total number of processors in World comm. |
397 |
integer :: rowIndex !! Row for this processor. |
398 |
integer :: columnIndex !! Column for this processor. |
399 |
integer :: nRows !! Total number of rows. |
400 |
integer :: nColumns !! Total number of columns. |
401 |
integer :: mpiErrors !! MPI errors. |
402 |
integer :: rowCommunicator !! MPI row communicator. |
403 |
integer :: columnCommunicator |
404 |
integer :: myWorldRank |
405 |
integer :: i |
406 |
|
407 |
|
408 |
if (.not. ComponentPlanSet) return |
409 |
status = 0 |
410 |
|
411 |
!! We make a dangerous assumption here that if numberProcessors is |
412 |
!! zero, then we need to get the information from MPI. |
413 |
if (thisComponentPlan%nProcessors == 0 ) then |
414 |
call mpi_comm_size( MPI_COMM_WORLD, nWorldProcessors,mpiErrors) |
415 |
if ( mpiErrors /= 0 ) then |
416 |
status = -1 |
417 |
return |
418 |
endif |
419 |
call mpi_comm_rank( MPI_COMM_WORLD,myWorldRank,mpiErrors) |
420 |
if ( mpiErrors /= 0 ) then |
421 |
status = -1 |
422 |
return |
423 |
endif |
424 |
|
425 |
else |
426 |
nWorldProcessors = thisComponentPlan%nProcessors |
427 |
myWorldRank = thisComponentPlan%myNode |
428 |
endif |
429 |
|
430 |
nColumnsMax = nint(sqrt(real(nWorldProcessors,kind=dp))) |
431 |
|
432 |
do i = 1, nColumnsMax |
433 |
if (mod(nWorldProcessors,i) == 0) nColumns = i |
434 |
end do |
435 |
|
436 |
nRows = nWorldProcessors/nColumns |
437 |
|
438 |
rowIndex = myWorldRank/nColumns |
439 |
|
440 |
call mpi_comm_split(mpi_comm_world,rowIndex,0,rowCommunicator,mpiErrors) |
441 |
if ( mpiErrors /= 0 ) then |
442 |
write(default_error,*) "MPI comm split failed at row communicator" |
443 |
status = -1 |
444 |
return |
445 |
endif |
446 |
|
447 |
columnIndex = mod(myWorldRank,nColumns) |
448 |
call mpi_comm_split(mpi_comm_world,columnIndex,0,columnCommunicator,mpiErrors) |
449 |
if ( mpiErrors /= 0 ) then |
450 |
write(default_error,*) "MPI comm split faild at columnCommunicator" |
451 |
status = -1 |
452 |
return |
453 |
endif |
454 |
|
455 |
! Set appropriate components of thisComponentPlan |
456 |
thisComponentPlan%rowComm = rowCommunicator |
457 |
thisComponentPlan%columnComm = columnCommunicator |
458 |
thisComponentPlan%rowIndex = rowIndex |
459 |
thisComponentPlan%columnIndex = columnIndex |
460 |
thisComponentPlan%nRows = nRows |
461 |
thisComponentPlan%nColumns = nColumns |
462 |
|
463 |
end subroutine make_Force_Grid |
464 |
|
465 |
!! initalizes a gather scatter plan |
466 |
subroutine plan_gather_scatter( nDim, nObjects, thisComponentPlan, & |
467 |
thisComm, this_plan, status) |
468 |
integer, intent(in) :: nDim !! Number of dimensions for gather scatter plan |
469 |
integer, intent(in) :: nObjects |
470 |
type (mpiComponentPlan), intent(in), target :: thisComponentPlan |
471 |
type (gs_plan), intent(out) :: this_plan !! MPI Component Plan |
472 |
integer, intent(in) :: thisComm !! MPI communicator for this plan |
473 |
|
474 |
integer :: arraySize !! size to allocate plan for |
475 |
integer, intent(out), optional :: status |
476 |
integer :: ierror |
477 |
integer :: i,junk |
478 |
|
479 |
if (present(status)) status = 0 |
480 |
|
481 |
!! Set gsComponentPlan pointer |
482 |
!! to the componet plan we want to use for this gather scatter plan. |
483 |
!! WARNING this could be dangerous since thisComponentPlan was origionally |
484 |
!! allocated in C++ and there is a significant difference between c and |
485 |
!! f95 pointers.... |
486 |
this_plan%gsComponentPlan => thisComponentPlan |
487 |
|
488 |
! Set this plan size for displs array. |
489 |
this_plan%gsPlanSize = nDim * nObjects |
490 |
|
491 |
! Duplicate communicator for this plan |
492 |
call mpi_comm_dup(thisComm, this_plan%myPlanComm, mpi_err) |
493 |
if (mpi_err /= 0) then |
494 |
if (present(status)) status = -1 |
495 |
return |
496 |
end if |
497 |
call mpi_comm_rank(this_plan%myPlanComm, this_plan%myPlanRank, mpi_err) |
498 |
if (mpi_err /= 0) then |
499 |
if (present(status)) status = -1 |
500 |
return |
501 |
end if |
502 |
|
503 |
call mpi_comm_size(this_plan%myPlanComm, this_plan%planNprocs, mpi_err) |
504 |
|
505 |
if (mpi_err /= 0) then |
506 |
if (present(status)) status = -1 |
507 |
return |
508 |
end if |
509 |
|
510 |
!! counts and displacements arrays are indexed from 0 to be compatable |
511 |
!! with MPI arrays. |
512 |
allocate (this_plan%counts(0:this_plan%planNprocs-1),STAT=ierror) |
513 |
if (ierror /= 0) then |
514 |
if (present(status)) status = -1 |
515 |
return |
516 |
end if |
517 |
|
518 |
allocate (this_plan%displs(0:this_plan%planNprocs-1),STAT=ierror) |
519 |
if (ierror /= 0) then |
520 |
if (present(status)) status = -1 |
521 |
return |
522 |
end if |
523 |
|
524 |
!! gather all the local sizes into a size # processors array. |
525 |
call mpi_allgather(this_plan%gsPlanSize,1,mpi_integer,this_plan%counts, & |
526 |
1,mpi_integer,thisComm,mpi_err) |
527 |
|
528 |
if (mpi_err /= 0) then |
529 |
if (present(status)) status = -1 |
530 |
return |
531 |
end if |
532 |
|
533 |
!! figure out the total number of particles in this plan |
534 |
this_plan%globalPlanSize = sum(this_plan%counts) |
535 |
|
536 |
!! initialize plan displacements. |
537 |
this_plan%displs(0) = 0 |
538 |
do i = 1, this_plan%planNprocs - 1,1 |
539 |
this_plan%displs(i) = this_plan%displs(i-1) + this_plan%counts(i-1) |
540 |
end do |
541 |
end subroutine plan_gather_scatter |
542 |
|
543 |
subroutine unplan_gather_scatter(this_plan) |
544 |
type (gs_plan), intent(inout) :: this_plan |
545 |
|
546 |
this_plan%gsComponentPlan => null() |
547 |
call mpi_comm_free(this_plan%myPlanComm,mpi_err) |
548 |
|
549 |
deallocate(this_plan%counts) |
550 |
deallocate(this_plan%displs) |
551 |
|
552 |
end subroutine unplan_gather_scatter |
553 |
|
554 |
subroutine gather_integer( sbuffer, rbuffer, this_plan, status) |
555 |
|
556 |
type (gs_plan), intent(inout) :: this_plan |
557 |
integer, dimension(:), intent(inout) :: sbuffer |
558 |
integer, dimension(:), intent(inout) :: rbuffer |
559 |
integer :: noffset |
560 |
integer, intent(out), optional :: status |
561 |
integer :: i |
562 |
|
563 |
if (present(status)) status = 0 |
564 |
noffset = this_plan%displs(this_plan%myPlanRank) |
565 |
|
566 |
! if (getmyNode() == 1) then |
567 |
! write(*,*) "Node 0 printing allgatherv vars" |
568 |
! write(*,*) "Noffset: ", noffset |
569 |
! write(*,*) "PlanSize: ", this_plan%gsPlanSize |
570 |
! write(*,*) "PlanComm: ", this_plan%myPlanComm |
571 |
! end if |
572 |
|
573 |
call mpi_allgatherv(sbuffer,this_plan%gsPlanSize, mpi_integer, & |
574 |
rbuffer,this_plan%counts,this_plan%displs,mpi_integer, & |
575 |
this_plan%myPlanComm, mpi_err) |
576 |
|
577 |
if (mpi_err /= 0) then |
578 |
if (present(status)) status = -1 |
579 |
endif |
580 |
|
581 |
end subroutine gather_integer |
582 |
|
583 |
subroutine gather_double( sbuffer, rbuffer, this_plan, status) |
584 |
|
585 |
type (gs_plan), intent(in) :: this_plan |
586 |
real( kind = DP ), dimension(:), intent(inout) :: sbuffer |
587 |
real( kind = DP ), dimension(:), intent(inout) :: rbuffer |
588 |
integer :: noffset |
589 |
integer, intent(out), optional :: status |
590 |
|
591 |
|
592 |
if (present(status)) status = 0 |
593 |
noffset = this_plan%displs(this_plan%myPlanRank) |
594 |
#ifdef PROFILE |
595 |
call cpu_time(commTimeInitial) |
596 |
#endif |
597 |
call mpi_allgatherv(sbuffer,this_plan%gsPlanSize, mpi_double_precision, & |
598 |
rbuffer,this_plan%counts,this_plan%displs,mpi_double_precision, & |
599 |
this_plan%myPlanComm, mpi_err) |
600 |
#ifdef PROFILE |
601 |
call cpu_time(commTimeFinal) |
602 |
commTime = commTime + commTimeFinal - commTimeInitial |
603 |
#endif |
604 |
|
605 |
if (mpi_err /= 0) then |
606 |
if (present(status)) status = -1 |
607 |
endif |
608 |
|
609 |
end subroutine gather_double |
610 |
|
611 |
subroutine gather_double_2d( sbuffer, rbuffer, this_plan, status) |
612 |
|
613 |
type (gs_plan), intent(in) :: this_plan |
614 |
real( kind = DP ), dimension(:,:), intent(inout) :: sbuffer |
615 |
real( kind = DP ), dimension(:,:), intent(inout) :: rbuffer |
616 |
integer :: noffset,i,ierror |
617 |
integer, intent(out), optional :: status |
618 |
|
619 |
external mpi_allgatherv |
620 |
|
621 |
if (present(status)) status = 0 |
622 |
|
623 |
! noffset = this_plan%displs(this_plan%me) |
624 |
#ifdef PROFILE |
625 |
call cpu_time(commTimeInitial) |
626 |
#endif |
627 |
|
628 |
call mpi_allgatherv(sbuffer,this_plan%gsPlanSize, mpi_double_precision, & |
629 |
rbuffer,this_plan%counts,this_plan%displs,mpi_double_precision, & |
630 |
this_plan%myPlanComm, mpi_err) |
631 |
|
632 |
#ifdef PROFILE |
633 |
call cpu_time(commTimeFinal) |
634 |
commTime = commTime + commTimeFinal - commTimeInitial |
635 |
#endif |
636 |
|
637 |
if (mpi_err /= 0) then |
638 |
if (present(status)) status = -1 |
639 |
endif |
640 |
|
641 |
end subroutine gather_double_2d |
642 |
|
643 |
subroutine scatter_double( sbuffer, rbuffer, this_plan, status) |
644 |
|
645 |
type (gs_plan), intent(in) :: this_plan |
646 |
real( kind = DP ), dimension(:), intent(inout) :: sbuffer |
647 |
real( kind = DP ), dimension(:), intent(inout) :: rbuffer |
648 |
integer, intent(out), optional :: status |
649 |
external mpi_reduce_scatter |
650 |
|
651 |
if (present(status)) status = 0 |
652 |
|
653 |
#ifdef PROFILE |
654 |
call cpu_time(commTimeInitial) |
655 |
#endif |
656 |
call mpi_reduce_scatter(sbuffer,rbuffer, this_plan%counts, & |
657 |
mpi_double_precision, MPI_SUM, this_plan%myPlanComm, mpi_err) |
658 |
#ifdef PROFILE |
659 |
call cpu_time(commTimeFinal) |
660 |
commTime = commTime + commTimeFinal - commTimeInitial |
661 |
#endif |
662 |
|
663 |
if (mpi_err /= 0) then |
664 |
if (present(status)) status = -1 |
665 |
endif |
666 |
|
667 |
end subroutine scatter_double |
668 |
|
669 |
subroutine scatter_double_2d( sbuffer, rbuffer, this_plan, status) |
670 |
|
671 |
type (gs_plan), intent(in) :: this_plan |
672 |
real( kind = DP ), dimension(:,:), intent(inout) :: sbuffer |
673 |
real( kind = DP ), dimension(:,:), intent(inout) :: rbuffer |
674 |
integer, intent(out), optional :: status |
675 |
external mpi_reduce_scatter |
676 |
|
677 |
if (present(status)) status = 0 |
678 |
#ifdef PROFILE |
679 |
call cpu_time(commTimeInitial) |
680 |
#endif |
681 |
|
682 |
call mpi_reduce_scatter(sbuffer,rbuffer, this_plan%counts, & |
683 |
mpi_double_precision, MPI_SUM, this_plan%myPlanComm, mpi_err) |
684 |
#ifdef PROFILE |
685 |
call cpu_time(commTimeFinal) |
686 |
commTime = commTime + commTimeFinal - commTimeInitial |
687 |
#endif |
688 |
|
689 |
if (mpi_err /= 0) then |
690 |
if (present(status)) status = -1 |
691 |
endif |
692 |
|
693 |
end subroutine scatter_double_2d |
694 |
|
695 |
subroutine setAtomTags(tags, status) |
696 |
integer, dimension(:) :: tags |
697 |
integer :: status |
698 |
|
699 |
integer :: alloc_stat |
700 |
|
701 |
integer :: nAtomsInCol |
702 |
integer :: nAtomsInRow |
703 |
|
704 |
status = 0 |
705 |
! allocate row arrays |
706 |
nAtomsInRow = getNatomsInRow(plan_atom_row) |
707 |
nAtomsInCol = getNatomsInCol(plan_atom_col) |
708 |
|
709 |
if(.not. allocated(AtomLocalToGlobal)) then |
710 |
allocate(AtomLocalToGlobal(size(tags)),STAT=alloc_stat) |
711 |
if (alloc_stat /= 0 ) then |
712 |
status = -1 |
713 |
return |
714 |
endif |
715 |
else |
716 |
deallocate(AtomLocalToGlobal) |
717 |
allocate(AtomLocalToGlobal(size(tags)),STAT=alloc_stat) |
718 |
if (alloc_stat /= 0 ) then |
719 |
status = -1 |
720 |
return |
721 |
endif |
722 |
|
723 |
endif |
724 |
|
725 |
AtomLocalToGlobal = tags |
726 |
|
727 |
if (.not. allocated(AtomRowToGlobal)) then |
728 |
allocate(AtomRowToGlobal(nAtomsInRow),STAT=alloc_stat) |
729 |
if (alloc_stat /= 0 ) then |
730 |
status = -1 |
731 |
return |
732 |
endif |
733 |
else |
734 |
deallocate(AtomRowToGlobal) |
735 |
allocate(AtomRowToGlobal(nAtomsInRow),STAT=alloc_stat) |
736 |
if (alloc_stat /= 0 ) then |
737 |
status = -1 |
738 |
return |
739 |
endif |
740 |
|
741 |
endif |
742 |
! allocate column arrays |
743 |
if (.not. allocated(AtomColToGlobal)) then |
744 |
allocate(AtomColToGlobal(nAtomsInCol),STAT=alloc_stat) |
745 |
if (alloc_stat /= 0 ) then |
746 |
status = -1 |
747 |
return |
748 |
endif |
749 |
else |
750 |
deallocate(AtomColToGlobal) |
751 |
allocate(AtomColToGlobal(nAtomsInCol),STAT=alloc_stat) |
752 |
if (alloc_stat /= 0 ) then |
753 |
status = -1 |
754 |
return |
755 |
endif |
756 |
endif |
757 |
|
758 |
call gather(tags, AtomRowToGlobal, plan_atom_row) |
759 |
call gather(tags, AtomColToGlobal, plan_atom_col) |
760 |
|
761 |
end subroutine setAtomTags |
762 |
|
763 |
subroutine setGroupTags(tags, status) |
764 |
integer, dimension(:) :: tags |
765 |
integer :: status |
766 |
|
767 |
integer :: alloc_stat |
768 |
|
769 |
integer :: nGroupsInCol |
770 |
integer :: nGroupsInRow |
771 |
|
772 |
status = 0 |
773 |
|
774 |
nGroupsInRow = getNgroupsInRow(plan_group_row) |
775 |
nGroupsInCol = getNgroupsInCol(plan_group_col) |
776 |
|
777 |
if(allocated(GroupLocalToGlobal)) then |
778 |
deallocate(GroupLocalToGlobal) |
779 |
endif |
780 |
allocate(GroupLocalToGlobal(size(tags)),STAT=alloc_stat) |
781 |
if (alloc_stat /= 0 ) then |
782 |
status = -1 |
783 |
return |
784 |
endif |
785 |
|
786 |
GroupLocalToGlobal = tags |
787 |
|
788 |
if(allocated(GroupRowToGlobal)) then |
789 |
deallocate(GroupRowToGlobal) |
790 |
endif |
791 |
allocate(GroupRowToGlobal(nGroupsInRow),STAT=alloc_stat) |
792 |
if (alloc_stat /= 0 ) then |
793 |
status = -1 |
794 |
return |
795 |
endif |
796 |
|
797 |
if(allocated(GroupColToGlobal)) then |
798 |
deallocate(GroupColToGlobal) |
799 |
endif |
800 |
allocate(GroupColToGlobal(nGroupsInCol),STAT=alloc_stat) |
801 |
if (alloc_stat /= 0 ) then |
802 |
status = -1 |
803 |
return |
804 |
endif |
805 |
|
806 |
call gather(tags, GroupRowToGlobal, plan_group_row) |
807 |
call gather(tags, GroupColToGlobal, plan_group_col) |
808 |
|
809 |
end subroutine setGroupTags |
810 |
|
811 |
function getNatomsInCol(thisplan) result(nInCol) |
812 |
type (gs_plan), intent(in) :: thisplan |
813 |
integer :: nInCol |
814 |
nInCol = thisplan%gsComponentPlan%nAtomsInColumn |
815 |
end function getNatomsInCol |
816 |
|
817 |
function getNatomsInRow(thisplan) result(nInRow) |
818 |
type (gs_plan), intent(in) :: thisplan |
819 |
integer :: nInRow |
820 |
nInRow = thisplan%gsComponentPlan%nAtomsInRow |
821 |
end function getNatomsInRow |
822 |
|
823 |
function getNgroupsInCol(thisplan) result(nInCol) |
824 |
type (gs_plan), intent(in) :: thisplan |
825 |
integer :: nInCol |
826 |
nInCol = thisplan%gsComponentPlan%nGroupsInColumn |
827 |
end function getNgroupsInCol |
828 |
|
829 |
function getNgroupsInRow(thisplan) result(nInRow) |
830 |
type (gs_plan), intent(in) :: thisplan |
831 |
integer :: nInRow |
832 |
nInRow = thisplan%gsComponentPlan%nGroupsInRow |
833 |
end function getNgroupsInRow |
834 |
|
835 |
function isMPISimSet() result(isthisSimSet) |
836 |
logical :: isthisSimSet |
837 |
if (isSimSet) then |
838 |
isthisSimSet = .true. |
839 |
else |
840 |
isthisSimSet = .false. |
841 |
endif |
842 |
end function isMPISimSet |
843 |
|
844 |
subroutine printComponentPlan(this_plan,printNode) |
845 |
|
846 |
type (mpiComponentPlan), intent(in) :: this_plan |
847 |
integer, optional :: printNode |
848 |
logical :: print_me = .false. |
849 |
|
850 |
if (present(printNode)) then |
851 |
if (printNode == mpiSim%myNode) print_me = .true. |
852 |
else |
853 |
print_me = .true. |
854 |
endif |
855 |
|
856 |
if (print_me) then |
857 |
write(default_error,*) "SetupSimParallel: writing component plan" |
858 |
|
859 |
write(default_error,*) "nMolGlobal: ", mpiSim%nMolGlobal |
860 |
write(default_error,*) "nAtomsGlobal: ", mpiSim%nAtomsGlobal |
861 |
write(default_error,*) "nAtomsLocal: ", mpiSim%nAtomsLocal |
862 |
write(default_error,*) "myNode: ", mpiSim%myNode |
863 |
write(default_error,*) "nProcessors: ", mpiSim%nProcessors |
864 |
write(default_error,*) "rowComm: ", mpiSim%rowComm |
865 |
write(default_error,*) "columnComm: ", mpiSim%columnComm |
866 |
write(default_error,*) "nRows: ", mpiSim%nRows |
867 |
write(default_error,*) "nColumns: ", mpiSim%nColumns |
868 |
write(default_error,*) "nAtomsInRow: ", mpiSim%nAtomsInRow |
869 |
write(default_error,*) "nAtomsInColumn: ", mpiSim%nAtomsInColumn |
870 |
write(default_error,*) "rowIndex: ", mpiSim%rowIndex |
871 |
write(default_error,*) "columnIndex: ", mpiSim%columnIndex |
872 |
endif |
873 |
end subroutine printComponentPlan |
874 |
|
875 |
function getMyNode() result(myNode) |
876 |
integer :: myNode |
877 |
myNode = mpiSim%myNode |
878 |
end function getMyNode |
879 |
|
880 |
#ifdef PROFILE |
881 |
subroutine printCommTime() |
882 |
write(*,*) "MPI communication time is: ", commTime |
883 |
end subroutine printCommTime |
884 |
|
885 |
function getCommTime() result(comm_time) |
886 |
real :: comm_time |
887 |
comm_time = commTime |
888 |
end function getCommTime |
889 |
|
890 |
#endif |
891 |
|
892 |
#endif // is_mpi |
893 |
end module mpiSimulation |
894 |
|