ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/mpiSimulation.cpp
(Generate patch)

Comparing trunk/OOPSE/libmdtools/mpiSimulation.cpp (file contents):
Revision 432 by chuckv, Thu Mar 27 23:33:40 2003 UTC vs.
Revision 829 by gezelter, Tue Oct 28 16:03:37 2003 UTC

# Line 1 | Line 1
1   #ifdef IS_MPI
2   #include <iostream>
3 < #include <cstdlib>
4 < #include <cstring>
5 < #include <cmath>
3 > #include <stdlib.h>
4 > #include <string.h>
5 > #include <math.h>
6   #include <mpi.h>
7 #include <mpi++.h>
7  
8   #include "mpiSimulation.hpp"
9   #include "simError.h"
10   #include "fortranWrappers.hpp"
11   #include "randomSPRNG.hpp"
12  
14 #define BASE_SEED 123456789
15
13   mpiSimulation* mpiSim;
14  
15   mpiSimulation::mpiSimulation(SimInfo* the_entryPlug)
# Line 20 | Line 17 | mpiSimulation::mpiSimulation(SimInfo* the_entryPlug)
17    entryPlug = the_entryPlug;
18    mpiPlug = new mpiSimData;
19    
20 <  mpiPlug->numberProcessors = MPI::COMM_WORLD.Get_size();
20 >  MPI_Comm_size(MPI_COMM_WORLD, &(mpiPlug->numberProcessors) );
21    mpiPlug->myNode = worldRank;
22  
23    MolToProcMap = new int[entryPlug->n_mol];
# Line 60 | Line 57 | int* mpiSimulation::divideLabor( void ){
57    int old_atoms, add_atoms, new_atoms;
58  
59    int nTarget;
60 <  int molIndex, atomIndex, compIndex, compStart;
60 >  int molIndex, atomIndex;
61    int done;
65  int nLocal, molLocal;
62    int i, j, loops, which_proc, nmol_local, natoms_local;
63    int nmol_global, natoms_global;
64 <  int local_index, index;
65 <  int smallDiff, bigDiff;
70 <  int baseSeed = BASE_SEED;
64 >  int local_index;
65 >  int baseSeed = entryPlug->getSeed();
66  
72  int testSum;
73
67    nComponents = entryPlug->nComponents;
68    compStamps = entryPlug->compStamps;
69    componentsNmol = entryPlug->componentsNmol;
# Line 83 | Line 76 | int* mpiSimulation::divideLabor( void ){
76    mpiPlug->nSRIGlobal = entryPlug->n_SRI;
77    mpiPlug->nMolGlobal = entryPlug->n_mol;
78  
79 +  
80    myRandom = new randomSPRNG( baseSeed );
81  
82 <  a = (double)mpiPlug->nMolGlobal / (double)mpiPlug->nAtomsGlobal;
82 >  a = 3.0 * (double)mpiPlug->nMolGlobal / (double)mpiPlug->nAtomsGlobal;
83  
84    // Initialize things that we'll send out later:
85    for (i = 0; i < mpiPlug->numberProcessors; i++ ) {
# Line 136 | Line 130 | int* mpiSimulation::divideLabor( void ){
130          add_atoms = compStamps[MolComponentType[i]]->getNAtoms();
131          new_atoms = old_atoms + add_atoms;
132  
139        // If the processor already had too many atoms, just skip this
140        // processor and try again.
141
133          // If we've been through this loop too many times, we need
134          // to just give up and assign the molecule to this processor
135          // and be done with it.
# Line 161 | Line 152 | int* mpiSimulation::divideLabor( void ){
152            done = 1;
153            continue;
154          }
164
165        if (old_atoms >= nTarget) continue;
155      
156          // If we can add this molecule to this processor without sending
157          // it above nTarget, then go ahead and do it:
# Line 179 | Line 168 | int* mpiSimulation::divideLabor( void ){
168          }
169  
170  
171 <        // The only situation left is where old_atoms < nTarget, but
172 <        // new_atoms > nTarget.   We want to accept this with some
173 <        // probability that dies off the farther we are from nTarget
171 >        // The only situation left is when new_atoms > nTarget.  We
172 >        // want to accept this with some probability that dies off the
173 >        // farther we are from nTarget
174  
175          // roughly:  x = new_atoms - nTarget
176          //           Pacc(x) = exp(- a * x)
177 <        // where a = 1 / (average atoms per molecule)
177 >        // where a = penalty / (average atoms per molecule)
178  
179          x = (double) (new_atoms - nTarget);
180          y = myRandom->getRandom();
181 <        
182 <        if (exp(- a * x) > y) {
181 >      
182 >        if (y < exp(- a * x)) {
183            MolToProcMap[i] = which_proc;
184            AtomsPerProc[which_proc] += add_atoms;
185            for (j = 0 ; j < add_atoms; j++ ) {
# Line 208 | Line 197 | int* mpiSimulation::divideLabor( void ){
197  
198      // Spray out this nonsense to all other processors:
199  
200 <    MPI::COMM_WORLD.Bcast(MolToProcMap, mpiPlug->nMolGlobal,
201 <                          MPI_INT, 0);
200 >    MPI_Bcast(MolToProcMap, mpiPlug->nMolGlobal,
201 >              MPI_INT, 0, MPI_COMM_WORLD);
202  
203 <    MPI::COMM_WORLD.Bcast(AtomToProcMap, mpiPlug->nAtomsGlobal,
204 <                          MPI_INT, 0);
203 >    MPI_Bcast(AtomToProcMap, mpiPlug->nAtomsGlobal,
204 >              MPI_INT, 0, MPI_COMM_WORLD);
205  
206 <    MPI::COMM_WORLD.Bcast(MolComponentType, mpiPlug->nMolGlobal,
207 <                          MPI_INT, 0);
206 >    MPI_Bcast(MolComponentType, mpiPlug->nMolGlobal,
207 >              MPI_INT, 0, MPI_COMM_WORLD);
208  
209 <    MPI::COMM_WORLD.Bcast(AtomsPerProc, mpiPlug->numberProcessors,
210 <                          MPI_INT, 0);    
209 >    MPI_Bcast(AtomsPerProc, mpiPlug->numberProcessors,
210 >              MPI_INT, 0, MPI_COMM_WORLD);    
211    } else {
212  
213      // Listen to your marching orders from processor 0:
214      
215 <    MPI::COMM_WORLD.Bcast(MolToProcMap, mpiPlug->nMolGlobal,
216 <                          MPI_INT, 0);
215 >    MPI_Bcast(MolToProcMap, mpiPlug->nMolGlobal,
216 >              MPI_INT, 0, MPI_COMM_WORLD);
217      
218 <    MPI::COMM_WORLD.Bcast(AtomToProcMap, mpiPlug->nAtomsGlobal,
219 <                          MPI_INT, 0);
218 >    MPI_Bcast(AtomToProcMap, mpiPlug->nAtomsGlobal,
219 >              MPI_INT, 0, MPI_COMM_WORLD);
220  
221 <    MPI::COMM_WORLD.Bcast(MolComponentType, mpiPlug->nMolGlobal,
222 <                          MPI_INT, 0);
221 >    MPI_Bcast(MolComponentType, mpiPlug->nMolGlobal,
222 >              MPI_INT, 0, MPI_COMM_WORLD);
223      
224 <    MPI::COMM_WORLD.Bcast(AtomsPerProc, mpiPlug->numberProcessors,
225 <                          MPI_INT, 0);
224 >    MPI_Bcast(AtomsPerProc, mpiPlug->numberProcessors,
225 >              MPI_INT, 0, MPI_COMM_WORLD);
226  
227  
228    }
# Line 255 | Line 244 | int* mpiSimulation::divideLabor( void ){
244      }
245    }
246  
247 <  std::cerr << "proc = " << mpiPlug->myNode << " atoms = " << natoms_local << "\n";
248 <
249 <  MPI::COMM_WORLD.Allreduce(&nmol_local,&nmol_global,1,MPI_INT,MPI_SUM);
250 <  MPI::COMM_WORLD.Allreduce(&natoms_local,&natoms_global,1,MPI_INT,MPI_SUM);
247 >  MPI_Allreduce(&nmol_local,&nmol_global,1,MPI_INT,MPI_SUM,
248 >                MPI_COMM_WORLD);
249 >  MPI_Allreduce(&natoms_local,&natoms_global,1,MPI_INT,
250 >                MPI_SUM, MPI_COMM_WORLD);
251    
252    if( nmol_global != entryPlug->n_mol ){
253      sprintf( painCave.errMsg,
# Line 303 | Line 292 | void mpiSimulation::mpiRefresh( void ){
292    int isError, i;
293    int *globalIndex = new int[mpiPlug->myNlocal];
294  
295 <  for(i=0; i<mpiPlug->myNlocal; i++) globalIndex[i] = entryPlug->atoms[i]->getGlobalIndex();
295 >  // Fortran indexing needs to be increased by 1 in order to get the 2 languages to
296 >  // not barf
297  
298 +  for(i=0; i<mpiPlug->myNlocal; i++) globalIndex[i] = entryPlug->atoms[i]->getGlobalIndex()+1;
299 +
300    
301    isError = 0;
302    setFsimParallel( mpiPlug, &(entryPlug->n_atoms), globalIndex, &isError );

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines