ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/mpiSimulation_module.F90
Revision: 357
Committed: Mon Mar 17 20:54:25 2003 UTC (21 years, 5 months ago) by mmeineke
File size: 21090 byte(s)
Log Message:
working on the fortran init ForceField interface.

File Contents

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