1 |
|
2 |
|
3 |
!! MPI support for long range forces using force decomposition |
4 |
!! on a square grid of processors. |
5 |
!! Corresponds to mpiSimunation.cpp for C++ |
6 |
!! mpi_module also contains a private interface for mpi f90 routines. |
7 |
!! |
8 |
!! @author Charles F. Vardeman II |
9 |
!! @author Matthew Meineke |
10 |
!! @version $Id: mpiSimulation_module.F90,v 1.5 2003-07-17 19:25:51 chuckv Exp $, $Date: 2003-07-17 19:25:51 $, $Name: not supported by cvs2svn $, $Revision: 1.5 $ |
11 |
|
12 |
module mpiSimulation |
13 |
use definitions |
14 |
#ifdef IS_MPI |
15 |
use mpi |
16 |
implicit none |
17 |
PRIVATE |
18 |
|
19 |
|
20 |
!! PUBLIC Subroutines contained in this module |
21 |
!! gather and scatter are a generic interface |
22 |
!! to gather and scatter routines |
23 |
public :: gather, scatter |
24 |
public :: setupSimParallel |
25 |
public :: replanSimParallel |
26 |
public :: getNcol |
27 |
public :: getNrow |
28 |
public :: isMPISimSet |
29 |
public :: printComponentPlan |
30 |
public :: getMyNode |
31 |
|
32 |
!! PUBLIC Subroutines contained in MPI module |
33 |
public :: mpi_bcast |
34 |
public :: mpi_allreduce |
35 |
public :: mpi_reduce |
36 |
public :: mpi_send |
37 |
public :: mpi_recv |
38 |
public :: mpi_get_processor_name |
39 |
public :: mpi_finalize |
40 |
|
41 |
!! PUBLIC mpi variables |
42 |
public :: mpi_comm_world |
43 |
public :: mpi_character |
44 |
public :: mpi_integer |
45 |
public :: mpi_double_precision |
46 |
public :: mpi_sum |
47 |
public :: mpi_max |
48 |
public :: mpi_status_size |
49 |
public :: mpi_any_source |
50 |
|
51 |
!! Safety logical to prevent access to ComponetPlan until |
52 |
!! set by C++. |
53 |
logical :: ComponentPlanSet = .false. |
54 |
|
55 |
|
56 |
!! generic mpi error declaration. |
57 |
integer,public :: mpi_err |
58 |
|
59 |
|
60 |
|
61 |
!! Include mpiComponentPlan type. mpiComponentPlan is a |
62 |
!! dual header file for both c and fortran. |
63 |
#define __FORTRAN90 |
64 |
#include "mpiComponentPlan.h" |
65 |
|
66 |
|
67 |
|
68 |
!! Tags used during force loop for parallel simulation |
69 |
integer, allocatable, dimension(:) :: tagLocal |
70 |
integer, public, allocatable, dimension(:) :: tagRow |
71 |
integer, public, allocatable, dimension(:) :: tagColumn |
72 |
|
73 |
!! Logical set true if mpiSimulation has been initialized |
74 |
logical :: isSimSet = .false. |
75 |
|
76 |
|
77 |
type (mpiComponentPlan) :: mpiSim |
78 |
|
79 |
!! gs_plan contains plans for gather and scatter routines |
80 |
type, public :: gs_plan |
81 |
private |
82 |
type (mpiComponentPlan), pointer :: gsComponentPlan => NULL() |
83 |
integer :: gsPlanSize !! size of this plan (nDim*nComponents) |
84 |
integer :: globalPlanSize !! size of all components in plan |
85 |
integer, dimension(:), pointer :: displs !! Displacements array for mpi indexed from 0. |
86 |
integer, dimension(:), pointer :: counts !! Counts array for mpi indexed from 0. |
87 |
integer :: myPlanComm !! My communicator for this plan |
88 |
integer :: myPlanRank !! My rank in this plan |
89 |
integer :: planNprocs !! Number of processors in this plan |
90 |
end type gs_plan |
91 |
|
92 |
! plans for different decompositions |
93 |
type (gs_plan), public :: plan_row |
94 |
type (gs_plan), public :: plan_row3d |
95 |
type (gs_plan), public :: plan_col |
96 |
type (gs_plan), public :: plan_col3d |
97 |
type(gs_plan), public :: plan_row_Rotation |
98 |
type(gs_plan), public :: plan_col_Rotation |
99 |
|
100 |
type (mpiComponentPlan), pointer :: simComponentPlan |
101 |
|
102 |
! interface for gather scatter routines |
103 |
!! Generic interface for gather. |
104 |
!! Gathers an local array into row or column array |
105 |
!! Interface provided for integer, double and double |
106 |
!! rank 2 arrays. |
107 |
interface gather |
108 |
module procedure gather_integer |
109 |
module procedure gather_double |
110 |
module procedure gather_double_2d |
111 |
end interface |
112 |
|
113 |
!! Generic interface for scatter. |
114 |
!! Scatters a row or column array, adding componets |
115 |
!! and reducing them to a local nComponent array. |
116 |
!! Interface provided for double and double rank=2 arrays. |
117 |
|
118 |
interface scatter |
119 |
module procedure scatter_double |
120 |
module procedure scatter_double_2d |
121 |
end interface |
122 |
|
123 |
|
124 |
|
125 |
contains |
126 |
|
127 |
!! Sets up mpiComponentPlan with structure passed from C++. |
128 |
subroutine setupSimParallel(thisComponentPlan,ntags,tags,status) |
129 |
! Passed Arguments |
130 |
! integer, intent(inout) :: nDim !! Number of dimensions |
131 |
!! mpiComponentPlan struct from C |
132 |
type (mpiComponentPlan), intent(inout) :: thisComponentPlan |
133 |
!! Number of tags passed, nlocal |
134 |
integer, intent(in) :: ntags |
135 |
!! Result status, 0 = normal, -1 = error |
136 |
integer, intent(out) :: status |
137 |
integer :: localStatus |
138 |
!! Global reference tag for local particles |
139 |
integer, dimension(ntags),intent(inout) :: tags |
140 |
|
141 |
write(*,*) 'mpiSim_mod thinks node', thisComponentPlan%myNode, ' has tags(1) = ', tags(1) |
142 |
|
143 |
|
144 |
|
145 |
status = 0 |
146 |
if (componentPlanSet) then |
147 |
return |
148 |
endif |
149 |
componentPlanSet = .true. |
150 |
|
151 |
!! copy c component plan to fortran |
152 |
mpiSim = thisComponentPlan |
153 |
write(*,*) "Seting up simParallel" |
154 |
|
155 |
call make_Force_Grid(mpiSim,localStatus) |
156 |
if (localStatus /= 0) then |
157 |
write(default_error,*) "Error creating force grid" |
158 |
status = -1 |
159 |
return |
160 |
endif |
161 |
|
162 |
call updateGridComponents(mpiSim,localStatus) |
163 |
if (localStatus /= 0) then |
164 |
write(default_error,*) "Error updating grid components" |
165 |
status = -1 |
166 |
return |
167 |
endif |
168 |
|
169 |
|
170 |
!! initialize gather and scatter plans used in this simulation |
171 |
call plan_gather_scatter(1,mpiSim%myNlocal,& |
172 |
mpiSim,mpiSim%rowComm,plan_row) |
173 |
call plan_gather_scatter(nDim,mpiSim%myNlocal,& |
174 |
mpiSim,mpiSim%rowComm,plan_row3d) |
175 |
call plan_gather_scatter(9,mpiSim%myNlocal,& |
176 |
mpiSim,mpiSim%rowComm,plan_row_Rotation) |
177 |
call plan_gather_scatter(1,mpiSim%myNlocal,& |
178 |
mpiSim,mpiSim%columnComm,plan_col) |
179 |
call plan_gather_scatter(nDim,mpiSim%myNlocal,& |
180 |
mpiSim,mpiSim%columnComm,plan_col3d) |
181 |
call plan_gather_scatter(9,mpiSim%myNlocal,& |
182 |
mpiSim,mpiSim%columnComm,plan_col_Rotation) |
183 |
|
184 |
|
185 |
|
186 |
! Initialize tags |
187 |
call setTags(tags,localStatus) |
188 |
if (localStatus /= 0) then |
189 |
status = -1 |
190 |
return |
191 |
endif |
192 |
isSimSet = .true. |
193 |
|
194 |
! call printComponentPlan(mpiSim,0) |
195 |
end subroutine setupSimParallel |
196 |
|
197 |
subroutine replanSimParallel(thisComponentPlan,status) |
198 |
! Passed Arguments |
199 |
!! mpiComponentPlan struct from C |
200 |
type (mpiComponentPlan), intent(inout) :: thisComponentPlan |
201 |
integer, intent(out) :: status |
202 |
integer :: localStatus |
203 |
integer :: mpierror |
204 |
status = 0 |
205 |
|
206 |
call updateGridComponents(thisComponentPlan,localStatus) |
207 |
if (localStatus /= 0) then |
208 |
status = -1 |
209 |
return |
210 |
endif |
211 |
|
212 |
!! Unplan Gather Scatter plans |
213 |
call unplan_gather_scatter(plan_row) |
214 |
call unplan_gather_scatter(plan_row3d) |
215 |
call unplan_gather_scatter(plan_row_Rotation) |
216 |
call unplan_gather_scatter(plan_col) |
217 |
call unplan_gather_scatter(plan_col3d) |
218 |
call unplan_gather_scatter(plan_col_Rotation) |
219 |
|
220 |
!! initialize gather and scatter plans used in this simulation |
221 |
call plan_gather_scatter(1,thisComponentPlan%myNlocal,& |
222 |
thisComponentPlan,thisComponentPlan%rowComm,plan_row) |
223 |
call plan_gather_scatter(nDim,thisComponentPlan%myNlocal,& |
224 |
thisComponentPlan,thisComponentPlan%rowComm,plan_row3d) |
225 |
call plan_gather_scatter(9,thisComponentPlan%myNlocal,& |
226 |
thisComponentPlan,thisComponentPlan%rowComm,plan_row_Rotation) |
227 |
call plan_gather_scatter(1,thisComponentPlan%myNlocal,& |
228 |
thisComponentPlan,thisComponentPlan%columnComm,plan_col) |
229 |
call plan_gather_scatter(nDim,thisComponentPlan%myNlocal,& |
230 |
thisComponentPlan,thisComponentPlan%rowComm,plan_col3d) |
231 |
call plan_gather_scatter(9,thisComponentPlan%myNlocal,& |
232 |
thisComponentPlan,thisComponentPlan%rowComm,plan_col_Rotation) |
233 |
|
234 |
|
235 |
|
236 |
end subroutine replanSimParallel |
237 |
|
238 |
!! Updates number of row and column components for long range forces. |
239 |
subroutine updateGridComponents(thisComponentPlan,status) |
240 |
type (mpiComponentPlan) :: thisComponentPlan !! mpiComponentPlan |
241 |
|
242 |
!! Status return |
243 |
!! - 0 Success |
244 |
!! - -1 Failure |
245 |
integer, intent(out) :: status |
246 |
integer :: nComponentsLocal |
247 |
integer :: nComponentsRow = 0 |
248 |
integer :: nComponentsColumn = 0 |
249 |
integer :: mpiErrors |
250 |
|
251 |
status = 0 |
252 |
if (.not. componentPlanSet) return |
253 |
|
254 |
if (thisComponentPlan%myNlocal == 0 ) then |
255 |
status = -1 |
256 |
return |
257 |
endif |
258 |
|
259 |
nComponentsLocal = thisComponentPlan%myNlocal |
260 |
|
261 |
write(*,*) "UpdateGridComponents: myNlocal ", nComponentsLocal |
262 |
call mpi_allreduce(nComponentsLocal,nComponentsRow,1,mpi_integer,& |
263 |
mpi_sum,thisComponentPlan%rowComm,mpiErrors) |
264 |
if (mpiErrors /= 0) then |
265 |
status = -1 |
266 |
return |
267 |
endif |
268 |
|
269 |
call mpi_allreduce(nComponentsLocal,nComponentsColumn,1,mpi_integer, & |
270 |
mpi_sum,thisComponentPlan%columnComm,mpiErrors) |
271 |
if (mpiErrors /= 0) then |
272 |
status = -1 |
273 |
return |
274 |
endif |
275 |
|
276 |
thisComponentPlan%nComponentsRow = nComponentsRow |
277 |
thisComponentPlan%nComponentsColumn = nComponentsColumn |
278 |
write(*,*) "UpdateGridComponents: myNRow ",& |
279 |
thisComponentPlan%nComponentsRow |
280 |
write(*,*) "UpdateGridComponents: myNColumn ",& |
281 |
thisComponentPlan%nComponentsColumn |
282 |
|
283 |
end subroutine updateGridComponents |
284 |
|
285 |
|
286 |
!! Creates a square force decomposition of processors into row and column |
287 |
!! communicators. |
288 |
subroutine make_Force_Grid(thisComponentPlan,status) |
289 |
type (mpiComponentPlan) :: thisComponentPlan |
290 |
integer, intent(out) :: status !! status returns -1 if error |
291 |
integer :: nColumnsMax !! Maximum number of columns |
292 |
integer :: nWorldProcessors !! Total number of processors in World comm. |
293 |
integer :: rowIndex !! Row for this processor. |
294 |
integer :: columnIndex !! Column for this processor. |
295 |
integer :: nRows !! Total number of rows. |
296 |
integer :: nColumns !! Total number of columns. |
297 |
integer :: mpiErrors !! MPI errors. |
298 |
integer :: rowCommunicator !! MPI row communicator. |
299 |
integer :: columnCommunicator |
300 |
integer :: myWorldRank |
301 |
integer :: i |
302 |
|
303 |
|
304 |
if (.not. ComponentPlanSet) return |
305 |
status = 0 |
306 |
|
307 |
!! We make a dangerous assumption here that if numberProcessors is |
308 |
!! zero, then we need to get the information from MPI. |
309 |
if (thisComponentPlan%numberProcessors == 0 ) then |
310 |
call mpi_comm_size( MPI_COMM_WORLD, nWorldProcessors,mpiErrors) |
311 |
if ( mpiErrors /= 0 ) then |
312 |
status = -1 |
313 |
return |
314 |
endif |
315 |
call mpi_comm_rank( MPI_COMM_WORLD,myWorldRank,mpiErrors) |
316 |
if ( mpiErrors /= 0 ) then |
317 |
status = -1 |
318 |
return |
319 |
endif |
320 |
|
321 |
else |
322 |
nWorldProcessors = thisComponentPlan%numberProcessors |
323 |
myWorldRank = thisComponentPlan%myNode |
324 |
endif |
325 |
|
326 |
|
327 |
|
328 |
|
329 |
nColumnsMax = nint(sqrt(real(nWorldProcessors,kind=dp))) |
330 |
|
331 |
do i = 1, nColumnsMax |
332 |
if (mod(nWorldProcessors,i) == 0) nColumns = i |
333 |
end do |
334 |
|
335 |
nRows = nWorldProcessors/nColumns |
336 |
|
337 |
rowIndex = myWorldRank/nColumns |
338 |
|
339 |
|
340 |
|
341 |
call mpi_comm_split(mpi_comm_world,rowIndex,0,rowCommunicator,mpiErrors) |
342 |
if ( mpiErrors /= 0 ) then |
343 |
write(default_error,*) "MPI comm split failed at row communicator" |
344 |
status = -1 |
345 |
return |
346 |
endif |
347 |
|
348 |
columnIndex = mod(myWorldRank,nColumns) |
349 |
call mpi_comm_split(mpi_comm_world,columnIndex,0,columnCommunicator,mpiErrors) |
350 |
if ( mpiErrors /= 0 ) then |
351 |
write(default_error,*) "MPI comm split faild at columnCommunicator" |
352 |
status = -1 |
353 |
return |
354 |
endif |
355 |
|
356 |
! Set appropriate components of thisComponentPlan |
357 |
thisComponentPlan%rowComm = rowCommunicator |
358 |
thisComponentPlan%columnComm = columnCommunicator |
359 |
thisComponentPlan%rowIndex = rowIndex |
360 |
thisComponentPlan%columnIndex = columnIndex |
361 |
thisComponentPlan%numberRows = nRows |
362 |
thisComponentPlan%numberColumns = nColumns |
363 |
|
364 |
|
365 |
end subroutine make_Force_Grid |
366 |
|
367 |
|
368 |
!! initalizes a gather scatter plan |
369 |
subroutine plan_gather_scatter( nDim,nComponents,thisComponentPlan, & |
370 |
thisComm, this_plan,status) |
371 |
integer, intent(in) :: nDim !! Number of dimensions for gather scatter plan |
372 |
integer, intent(in) :: nComponents |
373 |
type (mpiComponentPlan), intent(in), target :: thisComponentPlan |
374 |
type (gs_plan), intent(out) :: this_plan !! MPI Component Plan |
375 |
integer, intent(in) :: thisComm !! MPI communicator for this plan |
376 |
|
377 |
integer :: arraySize !! size to allocate plan for |
378 |
integer, intent(out), optional :: status |
379 |
integer :: ierror |
380 |
integer :: i,junk |
381 |
|
382 |
if (present(status)) status = 0 |
383 |
|
384 |
|
385 |
|
386 |
!! Set gsComponetPlan pointer |
387 |
!! to the componet plan we want to use for this gather scatter plan. |
388 |
!! WARNING this could be dangerous since thisComponentPlan was origionally |
389 |
!! allocated in C++ and there is a significant difference between c and |
390 |
!! f95 pointers.... |
391 |
this_plan%gsComponentPlan => thisComponentPlan |
392 |
|
393 |
! Set this plan size for displs array. |
394 |
this_plan%gsPlanSize = nDim * nComponents |
395 |
|
396 |
! Duplicate communicator for this plan |
397 |
call mpi_comm_dup(thisComm,this_plan%myPlanComm,mpi_err) |
398 |
if (mpi_err /= 0) then |
399 |
if (present(status)) status = -1 |
400 |
return |
401 |
end if |
402 |
call mpi_comm_rank(this_plan%myPlanComm,this_plan%myPlanRank,mpi_err) |
403 |
if (mpi_err /= 0) then |
404 |
if (present(status)) status = -1 |
405 |
return |
406 |
end if |
407 |
|
408 |
call mpi_comm_size(this_plan%myPlanComm,this_plan%planNprocs,mpi_err) |
409 |
|
410 |
if (mpi_err /= 0) then |
411 |
if (present(status)) status = -1 |
412 |
return |
413 |
end if |
414 |
|
415 |
!! counts and displacements arrays are indexed from 0 to be compatable |
416 |
!! with MPI arrays. |
417 |
allocate (this_plan%counts(0:this_plan%planNprocs-1),STAT=ierror) |
418 |
if (ierror /= 0) then |
419 |
if (present(status)) status = -1 |
420 |
return |
421 |
end if |
422 |
|
423 |
allocate (this_plan%displs(0:this_plan%planNprocs-1),STAT=ierror) |
424 |
if (ierror /= 0) then |
425 |
if (present(status)) status = -1 |
426 |
return |
427 |
end if |
428 |
|
429 |
!! gather all the local sizes into a size # processors array. |
430 |
call mpi_allgather(this_plan%gsPlanSize,1,mpi_integer,this_plan%counts, & |
431 |
1,mpi_integer,thisComm,mpi_err) |
432 |
|
433 |
if (mpi_err /= 0) then |
434 |
if (present(status)) status = -1 |
435 |
return |
436 |
end if |
437 |
|
438 |
|
439 |
!! figure out the total number of particles in this plan |
440 |
this_plan%globalPlanSize = sum(this_plan%counts) |
441 |
|
442 |
|
443 |
!! initialize plan displacements. |
444 |
this_plan%displs(0) = 0 |
445 |
do i = 1, this_plan%planNprocs - 1,1 |
446 |
this_plan%displs(i) = this_plan%displs(i-1) + this_plan%counts(i-1) |
447 |
end do |
448 |
|
449 |
|
450 |
end subroutine plan_gather_scatter |
451 |
|
452 |
|
453 |
subroutine unplan_gather_scatter(this_plan) |
454 |
type (gs_plan), intent(inout) :: this_plan |
455 |
|
456 |
|
457 |
this_plan%gsComponentPlan => null() |
458 |
call mpi_comm_free(this_plan%myPlanComm,mpi_err) |
459 |
|
460 |
deallocate(this_plan%counts) |
461 |
deallocate(this_plan%displs) |
462 |
|
463 |
end subroutine unplan_gather_scatter |
464 |
|
465 |
subroutine gather_integer( sbuffer, rbuffer, this_plan, status) |
466 |
|
467 |
type (gs_plan), intent(in) :: this_plan |
468 |
integer, dimension(:), intent(in) :: sbuffer |
469 |
integer, dimension(:), intent(in) :: rbuffer |
470 |
integer :: noffset |
471 |
integer, intent(out), optional :: status |
472 |
integer :: i |
473 |
|
474 |
|
475 |
|
476 |
if (present(status)) status = 0 |
477 |
noffset = this_plan%displs(this_plan%myPlanRank) |
478 |
|
479 |
! if (getmyNode() == 1) then |
480 |
! write(*,*) "Node 0 printing allgatherv vars" |
481 |
! write(*,*) "Noffset: ", noffset |
482 |
! write(*,*) "PlanSize: ", this_plan%gsPlanSize |
483 |
! write(*,*) "PlanComm: ", this_plan%myPlanComm |
484 |
! end if |
485 |
|
486 |
call mpi_allgatherv(sbuffer,this_plan%gsPlanSize, mpi_integer, & |
487 |
rbuffer,this_plan%counts,this_plan%displs,mpi_integer, & |
488 |
this_plan%myPlanComm, mpi_err) |
489 |
|
490 |
if (mpi_err /= 0) then |
491 |
if (present(status)) status = -1 |
492 |
endif |
493 |
|
494 |
end subroutine gather_integer |
495 |
|
496 |
subroutine gather_double( sbuffer, rbuffer, this_plan, status) |
497 |
|
498 |
type (gs_plan), intent(in) :: this_plan |
499 |
real( kind = DP ), dimension(:), intent(in) :: sbuffer |
500 |
real( kind = DP ), dimension(:), intent(in) :: rbuffer |
501 |
integer :: noffset |
502 |
integer, intent(out), optional :: status |
503 |
|
504 |
|
505 |
if (present(status)) status = 0 |
506 |
noffset = this_plan%displs(this_plan%myPlanRank) |
507 |
|
508 |
call mpi_allgatherv(sbuffer,this_plan%gsPlanSize, mpi_double_precision, & |
509 |
rbuffer,this_plan%counts,this_plan%displs,mpi_double_precision, & |
510 |
this_plan%myPlanComm, mpi_err) |
511 |
|
512 |
if (mpi_err /= 0) then |
513 |
if (present(status)) status = -1 |
514 |
endif |
515 |
|
516 |
end subroutine gather_double |
517 |
|
518 |
subroutine gather_double_2d( sbuffer, rbuffer, this_plan, status) |
519 |
|
520 |
type (gs_plan), intent(in) :: this_plan |
521 |
real( kind = DP ), dimension(:,:), intent(in) :: sbuffer |
522 |
real( kind = DP ), dimension(:,:), intent(in) :: rbuffer |
523 |
integer :: noffset,i,ierror |
524 |
integer, intent(out), optional :: status |
525 |
|
526 |
external mpi_allgatherv |
527 |
|
528 |
if (present(status)) status = 0 |
529 |
|
530 |
! noffset = this_plan%displs(this_plan%me) |
531 |
|
532 |
call mpi_allgatherv(sbuffer,this_plan%gsPlanSize, mpi_double_precision, & |
533 |
rbuffer,this_plan%counts,this_plan%displs,mpi_double_precision, & |
534 |
this_plan%myPlanComm, mpi_err) |
535 |
|
536 |
if (mpi_err /= 0) then |
537 |
if (present(status)) status = -1 |
538 |
endif |
539 |
|
540 |
end subroutine gather_double_2d |
541 |
|
542 |
subroutine scatter_double( sbuffer, rbuffer, this_plan, status) |
543 |
|
544 |
type (gs_plan), intent(in) :: this_plan |
545 |
real( kind = DP ), dimension(:), intent(in) :: sbuffer |
546 |
real( kind = DP ), dimension(:), intent(in) :: rbuffer |
547 |
integer, intent(out), optional :: status |
548 |
external mpi_reduce_scatter |
549 |
|
550 |
if (present(status)) status = 0 |
551 |
|
552 |
call mpi_reduce_scatter(sbuffer,rbuffer, this_plan%counts, & |
553 |
mpi_double_precision, MPI_SUM, this_plan%myPlanComm, mpi_err) |
554 |
|
555 |
if (mpi_err /= 0) then |
556 |
if (present(status)) status = -1 |
557 |
endif |
558 |
|
559 |
end subroutine scatter_double |
560 |
|
561 |
subroutine scatter_double_2d( sbuffer, rbuffer, this_plan, status) |
562 |
|
563 |
type (gs_plan), intent(in) :: this_plan |
564 |
real( kind = DP ), dimension(:,:), intent(in) :: sbuffer |
565 |
real( kind = DP ), dimension(:,:), intent(in) :: rbuffer |
566 |
integer, intent(out), optional :: status |
567 |
external mpi_reduce_scatter |
568 |
|
569 |
if (present(status)) status = 0 |
570 |
call mpi_reduce_scatter(sbuffer,rbuffer, this_plan%counts, & |
571 |
mpi_double_precision, MPI_SUM, this_plan%myPlanComm, mpi_err) |
572 |
|
573 |
if (mpi_err /= 0) then |
574 |
if (present(status)) status = -1 |
575 |
endif |
576 |
|
577 |
end subroutine scatter_double_2d |
578 |
|
579 |
|
580 |
subroutine setTags(tags,status) |
581 |
integer, dimension(:) :: tags |
582 |
integer :: status |
583 |
|
584 |
integer :: alloc_stat |
585 |
|
586 |
integer :: ncol |
587 |
integer :: nrow |
588 |
|
589 |
status = 0 |
590 |
! allocate row arrays |
591 |
nrow = getNrow(plan_row) |
592 |
ncol = getNcol(plan_col) |
593 |
|
594 |
if (.not. allocated(tagRow)) then |
595 |
allocate(tagRow(nrow),STAT=alloc_stat) |
596 |
if (alloc_stat /= 0 ) then |
597 |
status = -1 |
598 |
return |
599 |
endif |
600 |
else |
601 |
deallocate(tagRow) |
602 |
allocate(tagRow(nrow),STAT=alloc_stat) |
603 |
if (alloc_stat /= 0 ) then |
604 |
status = -1 |
605 |
return |
606 |
endif |
607 |
|
608 |
endif |
609 |
! allocate column arrays |
610 |
if (.not. allocated(tagColumn)) then |
611 |
allocate(tagColumn(ncol),STAT=alloc_stat) |
612 |
if (alloc_stat /= 0 ) then |
613 |
status = -1 |
614 |
return |
615 |
endif |
616 |
else |
617 |
deallocate(tagColumn) |
618 |
allocate(tagColumn(ncol),STAT=alloc_stat) |
619 |
if (alloc_stat /= 0 ) then |
620 |
status = -1 |
621 |
return |
622 |
endif |
623 |
endif |
624 |
|
625 |
call gather(tags,tagRow,plan_row) |
626 |
call gather(tags,tagColumn,plan_col) |
627 |
|
628 |
end subroutine setTags |
629 |
|
630 |
pure function getNcol(thisplan) result(ncol) |
631 |
type (gs_plan), intent(in) :: thisplan |
632 |
integer :: ncol |
633 |
ncol = thisplan%gsComponentPlan%nComponentsColumn |
634 |
end function getNcol |
635 |
|
636 |
pure function getNrow(thisplan) result(nrow) |
637 |
type (gs_plan), intent(in) :: thisplan |
638 |
integer :: nrow |
639 |
nrow = thisplan%gsComponentPlan%nComponentsRow |
640 |
end function getNrow |
641 |
|
642 |
function isMPISimSet() result(isthisSimSet) |
643 |
logical :: isthisSimSet |
644 |
if (isSimSet) then |
645 |
isthisSimSet = .true. |
646 |
else |
647 |
isthisSimSet = .false. |
648 |
endif |
649 |
end function isMPISimSet |
650 |
|
651 |
|
652 |
|
653 |
subroutine printComponentPlan(this_plan,printNode) |
654 |
|
655 |
type (mpiComponentPlan), intent(in) :: this_plan |
656 |
integer, optional :: printNode |
657 |
logical :: print_me = .false. |
658 |
|
659 |
if (present(printNode)) then |
660 |
if (printNode == mpiSim%myNode) print_me = .true. |
661 |
else |
662 |
print_me = .true. |
663 |
endif |
664 |
|
665 |
if (print_me) then |
666 |
write(default_error,*) "SetupSimParallel: writing component plan" |
667 |
|
668 |
write(default_error,*) "nMolGlobal: ", mpiSim%nMolGlobal |
669 |
write(default_error,*) "nAtomsGlobal: ", mpiSim%nAtomsGlobal |
670 |
write(default_error,*) "nBondGlobal: ", mpiSim%nBondsGlobal |
671 |
write(default_error,*) "nTorsionsGlobal: ", mpiSim%nTorsionsGlobal |
672 |
write(default_error,*) "nSRIGlobal: ", mpiSim%nSRIGlobal |
673 |
write(default_error,*) "myNlocal: ", mpiSim%myNlocal |
674 |
write(default_error,*) "myNode: ", mpiSim%myNode |
675 |
write(default_error,*) "numberProcessors: ", mpiSim%numberProcessors |
676 |
write(default_error,*) "rowComm: ", mpiSim%rowComm |
677 |
write(default_error,*) "columnComm: ", mpiSim%columnComm |
678 |
write(default_error,*) "numberRows: ", mpiSim%numberRows |
679 |
write(default_error,*) "numberColumns: ", mpiSim%numberColumns |
680 |
write(default_error,*) "nComponentsRow: ", mpiSim%nComponentsRow |
681 |
write(default_error,*) "nComponentsColumn: ", mpiSim%nComponentsColumn |
682 |
write(default_error,*) "rowIndex: ", mpiSim%rowIndex |
683 |
write(default_error,*) "columnIndex: ", mpiSim%columnIndex |
684 |
endif |
685 |
end subroutine printComponentPlan |
686 |
|
687 |
function getMyNode() result(myNode) |
688 |
integer :: myNode |
689 |
myNode = mpiSim%myNode |
690 |
end function getMyNode |
691 |
|
692 |
#endif // is_mpi |
693 |
end module mpiSimulation |
694 |
|
695 |
|