ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/simParallel.F90
Revision: 3111
Committed: Mon Jan 8 21:29:50 2007 UTC (18 years, 3 months ago) by chuckv
File size: 30864 byte(s)
Log Message:
Added more error reporting to simParallel.

File Contents

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