ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/mpiSimulation.cpp
Revision: 441
Committed: Tue Apr 1 16:50:14 2003 UTC (21 years, 3 months ago) by chuckv
File size: 8560 byte(s)
Log Message:
more bug fixes....

File Contents

# User Rev Content
1 mmeineke 377 #ifdef IS_MPI
2 chuckv 432 #include <iostream>
3 mmeineke 377 #include <cstdlib>
4     #include <cstring>
5 gezelter 416 #include <cmath>
6 mmeineke 377 #include <mpi.h>
7     #include <mpi++.h>
8    
9     #include "mpiSimulation.hpp"
10     #include "simError.h"
11     #include "fortranWrappers.hpp"
12 gezelter 405 #include "randomSPRNG.hpp"
13 mmeineke 377
14 gezelter 416 #define BASE_SEED 123456789
15 mmeineke 377
16     mpiSimulation* mpiSim;
17    
18     mpiSimulation::mpiSimulation(SimInfo* the_entryPlug)
19     {
20     entryPlug = the_entryPlug;
21     mpiPlug = new mpiSimData;
22    
23     mpiPlug->numberProcessors = MPI::COMM_WORLD.Get_size();
24     mpiPlug->myNode = worldRank;
25 gezelter 405
26     MolToProcMap = new int[entryPlug->n_mol];
27     MolComponentType = new int[entryPlug->n_mol];
28     AtomToProcMap = new int[entryPlug->n_atoms];
29    
30 mmeineke 377 mpiSim = this;
31     wrapMeSimParallel( this );
32     }
33    
34    
35     mpiSimulation::~mpiSimulation(){
36    
37 mmeineke 418 delete[] MolToProcMap;
38     delete[] MolComponentType;
39     delete[] AtomToProcMap;
40    
41 mmeineke 377 delete mpiPlug;
42     // perhaps we should let fortran know the party is over.
43    
44     }
45    
46     int* mpiSimulation::divideLabor( void ){
47    
48     int* globalIndex;
49    
50     int nComponents;
51     MoleculeStamp** compStamps;
52 gezelter 416 randomSPRNG *myRandom;
53 mmeineke 377 int* componentsNmol;
54 gezelter 405 int* AtomsPerProc;
55 mmeineke 377
56     double numerator;
57     double denominator;
58     double precast;
59 gezelter 405 double x, y, a;
60     int old_atoms, add_atoms, new_atoms;
61 mmeineke 377
62     int nTarget;
63     int molIndex, atomIndex, compIndex, compStart;
64     int done;
65     int nLocal, molLocal;
66 gezelter 416 int i, j, loops, which_proc, nmol_local, natoms_local;
67     int nmol_global, natoms_global;
68     int local_index, index;
69 mmeineke 377 int smallDiff, bigDiff;
70 gezelter 416 int baseSeed = BASE_SEED;
71 mmeineke 377
72     int testSum;
73    
74     nComponents = entryPlug->nComponents;
75     compStamps = entryPlug->compStamps;
76     componentsNmol = entryPlug->componentsNmol;
77 gezelter 405 AtomsPerProc = new int[mpiPlug->numberProcessors];
78    
79 mmeineke 377 mpiPlug->nAtomsGlobal = entryPlug->n_atoms;
80     mpiPlug->nBondsGlobal = entryPlug->n_bonds;
81     mpiPlug->nBendsGlobal = entryPlug->n_bends;
82     mpiPlug->nTorsionsGlobal = entryPlug->n_torsions;
83     mpiPlug->nSRIGlobal = entryPlug->n_SRI;
84     mpiPlug->nMolGlobal = entryPlug->n_mol;
85    
86 gezelter 416 myRandom = new randomSPRNG( baseSeed );
87 chuckv 404
88 chuckv 434 a = 3.0 * (double)mpiPlug->nMolGlobal / (double)mpiPlug->nAtomsGlobal;
89 chuckv 404
90 gezelter 405 // Initialize things that we'll send out later:
91     for (i = 0; i < mpiPlug->numberProcessors; i++ ) {
92     AtomsPerProc[i] = 0;
93     }
94     for (i = 0; i < mpiPlug->nMolGlobal; i++ ) {
95     // default to an error condition:
96     MolToProcMap[i] = -1;
97     MolComponentType[i] = -1;
98     }
99     for (i = 0; i < mpiPlug->nAtomsGlobal; i++ ) {
100     // default to an error condition:
101     AtomToProcMap[i] = -1;
102     }
103    
104     if (mpiPlug->myNode == 0) {
105     numerator = (double) entryPlug->n_atoms;
106     denominator = (double) mpiPlug->numberProcessors;
107     precast = numerator / denominator;
108     nTarget = (int)( precast + 0.5 );
109 chuckv 404
110 gezelter 405 // Build the array of molecule component types first
111     molIndex = 0;
112     for (i=0; i < nComponents; i++) {
113     for (j=0; j < componentsNmol[i]; j++) {
114     MolComponentType[molIndex] = i;
115     molIndex++;
116     }
117     }
118 chuckv 404
119 gezelter 405 atomIndex = 0;
120 chuckv 404
121 gezelter 405 for (i = 0; i < molIndex; i++ ) {
122 chuckv 404
123 gezelter 405 done = 0;
124     loops = 0;
125 chuckv 404
126 gezelter 405 while( !done ){
127     loops++;
128    
129     // Pick a processor at random
130 chuckv 404
131 gezelter 416 which_proc = (int) (myRandom->getRandom() * mpiPlug->numberProcessors);
132 chuckv 404
133 gezelter 405 // How many atoms does this processor have?
134    
135     old_atoms = AtomsPerProc[which_proc];
136 chuckv 432 add_atoms = compStamps[MolComponentType[i]]->getNAtoms();
137     new_atoms = old_atoms + add_atoms;
138 chuckv 404
139 gezelter 405 // If we've been through this loop too many times, we need
140     // to just give up and assign the molecule to this processor
141     // and be done with it.
142    
143     if (loops > 100) {
144     sprintf( painCave.errMsg,
145     "I've tried 100 times to assign molecule %d to a "
146     " processor, but can't find a good spot.\n"
147     "I'm assigning it at random to processor %d.\n",
148     i, which_proc);
149     painCave.isFatal = 0;
150     simError();
151    
152     MolToProcMap[i] = which_proc;
153     AtomsPerProc[which_proc] += add_atoms;
154     for (j = 0 ; j < add_atoms; j++ ) {
155 mmeineke 422 AtomToProcMap[atomIndex] = which_proc;
156     atomIndex++;
157 gezelter 405 }
158     done = 1;
159     continue;
160     }
161 chuckv 432
162     // If we can add this molecule to this processor without sending
163     // it above nTarget, then go ahead and do it:
164    
165     if (new_atoms <= nTarget) {
166     MolToProcMap[i] = which_proc;
167     AtomsPerProc[which_proc] += add_atoms;
168     for (j = 0 ; j < add_atoms; j++ ) {
169     AtomToProcMap[atomIndex] = which_proc;
170     atomIndex++;
171     }
172     done = 1;
173     continue;
174     }
175    
176    
177 chuckv 434 // The only situation left is when new_atoms > nTarget. We
178     // want to accept this with some probability that dies off the
179     // farther we are from nTarget
180 gezelter 405
181     // roughly: x = new_atoms - nTarget
182     // Pacc(x) = exp(- a * x)
183 chuckv 434 // where a = penalty / (average atoms per molecule)
184 gezelter 405
185     x = (double) (new_atoms - nTarget);
186 gezelter 416 y = myRandom->getRandom();
187 chuckv 434
188     if (y < exp(- a * x)) {
189 gezelter 405 MolToProcMap[i] = which_proc;
190     AtomsPerProc[which_proc] += add_atoms;
191     for (j = 0 ; j < add_atoms; j++ ) {
192 mmeineke 422 AtomToProcMap[atomIndex] = which_proc;
193     atomIndex++;
194     }
195 gezelter 405 done = 1;
196     continue;
197     } else {
198     continue;
199     }
200    
201 mmeineke 377 }
202     }
203 gezelter 405
204     // Spray out this nonsense to all other processors:
205    
206 mmeineke 418 MPI::COMM_WORLD.Bcast(MolToProcMap, mpiPlug->nMolGlobal,
207 gezelter 405 MPI_INT, 0);
208    
209 mmeineke 418 MPI::COMM_WORLD.Bcast(AtomToProcMap, mpiPlug->nAtomsGlobal,
210 gezelter 405 MPI_INT, 0);
211    
212 mmeineke 418 MPI::COMM_WORLD.Bcast(MolComponentType, mpiPlug->nMolGlobal,
213 gezelter 405 MPI_INT, 0);
214    
215 mmeineke 418 MPI::COMM_WORLD.Bcast(AtomsPerProc, mpiPlug->numberProcessors,
216 gezelter 405 MPI_INT, 0);
217     } else {
218    
219     // Listen to your marching orders from processor 0:
220 mmeineke 377
221 mmeineke 418 MPI::COMM_WORLD.Bcast(MolToProcMap, mpiPlug->nMolGlobal,
222 gezelter 405 MPI_INT, 0);
223 mmeineke 377
224 mmeineke 418 MPI::COMM_WORLD.Bcast(AtomToProcMap, mpiPlug->nAtomsGlobal,
225 gezelter 405 MPI_INT, 0);
226    
227 mmeineke 418 MPI::COMM_WORLD.Bcast(MolComponentType, mpiPlug->nMolGlobal,
228 gezelter 405 MPI_INT, 0);
229    
230 mmeineke 418 MPI::COMM_WORLD.Bcast(AtomsPerProc, mpiPlug->numberProcessors,
231 gezelter 405 MPI_INT, 0);
232 chuckv 432
233    
234 mmeineke 377 }
235    
236    
237 gezelter 405 // Let's all check for sanity:
238    
239     nmol_local = 0;
240     for (i = 0 ; i < mpiPlug->nMolGlobal; i++ ) {
241     if (MolToProcMap[i] == mpiPlug->myNode) {
242     nmol_local++;
243     }
244 mmeineke 377 }
245    
246 gezelter 405 natoms_local = 0;
247     for (i = 0; i < mpiPlug->nAtomsGlobal; i++) {
248     if (AtomToProcMap[i] == mpiPlug->myNode) {
249     natoms_local++;
250     }
251     }
252 mmeineke 377
253 gezelter 405 MPI::COMM_WORLD.Allreduce(&nmol_local,&nmol_global,1,MPI_INT,MPI_SUM);
254     MPI::COMM_WORLD.Allreduce(&natoms_local,&natoms_global,1,MPI_INT,MPI_SUM);
255 mmeineke 377
256 gezelter 405 if( nmol_global != entryPlug->n_mol ){
257     sprintf( painCave.errMsg,
258     "The sum of all nmol_local, %d, did not equal the "
259     "total number of molecules, %d.\n",
260     nmol_global, entryPlug->n_mol );
261     painCave.isFatal = 1;
262     simError();
263 mmeineke 377 }
264 gezelter 405
265     if( natoms_global != entryPlug->n_atoms ){
266     sprintf( painCave.errMsg,
267     "The sum of all natoms_local, %d, did not equal the "
268     "total number of atoms, %d.\n",
269     natoms_global, entryPlug->n_atoms );
270     painCave.isFatal = 1;
271     simError();
272     }
273 mmeineke 377
274     sprintf( checkPointMsg,
275     "Successfully divided the molecules among the processors.\n" );
276     MPIcheckPoint();
277    
278 gezelter 405 mpiPlug->myNMol = nmol_local;
279     mpiPlug->myNlocal = natoms_local;
280 mmeineke 377
281     globalIndex = new int[mpiPlug->myNlocal];
282 gezelter 405 local_index = 0;
283     for (i = 0; i < mpiPlug->nAtomsGlobal; i++) {
284     if (AtomToProcMap[i] == mpiPlug->myNode) {
285 mmeineke 422 globalIndex[local_index] = i;
286 chuckv 406 local_index++;
287 gezelter 405 }
288 mmeineke 377 }
289 chuckv 432
290 gezelter 416 return globalIndex;
291 mmeineke 377 }
292    
293    
294     void mpiSimulation::mpiRefresh( void ){
295    
296     int isError, i;
297     int *globalIndex = new int[mpiPlug->myNlocal];
298    
299 chuckv 441 // Fortran indexing needs to be increased by 1 in order to get the 2 languages to
300     // not barf
301 mmeineke 377
302 chuckv 441 for(i=0; i<mpiPlug->myNlocal; i++) globalIndex[i] = entryPlug->atoms[i]->getGlobalIndex()+1;
303    
304 mmeineke 377
305     isError = 0;
306     setFsimParallel( mpiPlug, &(entryPlug->n_atoms), globalIndex, &isError );
307     if( isError ){
308    
309     sprintf( painCave.errMsg,
310     "mpiRefresh errror: fortran didn't like something we gave it.\n" );
311     painCave.isFatal = 1;
312     simError();
313     }
314    
315     delete[] globalIndex;
316    
317     sprintf( checkPointMsg,
318     " mpiRefresh successful.\n" );
319     MPIcheckPoint();
320     }
321    
322    
323     #endif // is_mpi