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