ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimCreator.cpp
(Generate patch)

Comparing branches/development/src/brains/SimCreator.cpp (file contents):
Revision 1577 by gezelter, Wed Jun 8 20:26:56 2011 UTC vs.
Revision 1803 by gezelter, Wed Oct 3 14:20:07 2012 UTC

# Line 36 | Line 36
36   * [1]  Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).            
37   * [2]  Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).          
38   * [3]  Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).          
39 < * [4]  Vardeman & Gezelter, in progress (2009).                        
39 > * [4]  Kuang & Gezelter,  J. Chem. Phys. 133, 164101 (2010).
40 > * [5]  Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41   */
42  
43   /**
# Line 55 | Line 56
56   #include "brains/SimCreator.hpp"
57   #include "brains/SimSnapshotManager.hpp"
58   #include "io/DumpReader.hpp"
59 < #include "UseTheForce/ForceFieldFactory.hpp"
59 > #include "brains/ForceField.hpp"
60   #include "utils/simError.h"
61   #include "utils/StringUtils.hpp"
62   #include "math/SeqRandNumGen.hpp"
# Line 75 | Line 76
76   #include "antlr/NoViableAltForCharException.hpp"
77   #include "antlr/NoViableAltException.hpp"
78  
79 + #include "types/DirectionalAdapter.hpp"
80 + #include "types/MultipoleAdapter.hpp"
81 + #include "types/EAMAdapter.hpp"
82 + #include "types/SuttonChenAdapter.hpp"
83 + #include "types/PolarizableAdapter.hpp"
84 + #include "types/FixedChargeAdapter.hpp"
85 + #include "types/FluctuatingChargeAdapter.hpp"
86 +
87   #ifdef IS_MPI
88 + #include "mpi.h"
89   #include "math/ParallelRandNumGen.hpp"
90   #endif
91  
92   namespace OpenMD {
93    
94 <  Globals* SimCreator::parseFile(std::istream& rawMetaDataStream, const std::string& filename, int startOfMetaDataBlock ){
94 >  Globals* SimCreator::parseFile(std::istream& rawMetaDataStream, const std::string& filename, int mdFileVersion, int startOfMetaDataBlock ){
95      Globals* simParams = NULL;
96      try {
97  
# Line 92 | Line 102 | namespace OpenMD {
102        const int masterNode = 0;
103        int commStatus;
104        if (worldRank == masterNode) {
105 < #endif
106 <                
105 >        commStatus = MPI_Bcast(&mdFileVersion, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
106 > #endif                
107          SimplePreprocessor preprocessor;
108          preprocessor.preprocess(rawMetaDataStream, filename, startOfMetaDataBlock, ppStream);
109                  
# Line 106 | Line 116 | namespace OpenMD {
116              
117                  
118        } else {
119 +
120 +        commStatus = MPI_Bcast(&mdFileVersion, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
121 +
122          //get stream size
123          commStatus = MPI_Bcast(&streamSize, 1, MPI_LONG, masterNode, MPI_COMM_WORLD);  
124  
# Line 229 | Line 242 | namespace OpenMD {
242        simError();
243      }
244  
245 +    simParams->setMDfileVersion(mdFileVersion);
246      return simParams;
247    }
248    
# Line 242 | Line 256 | namespace OpenMD {
256      int metaDataBlockStart = -1;
257      int metaDataBlockEnd = -1;
258      int i;
259 <    int mdOffset;
259 >    streamoff mdOffset;
260 >    int mdFileVersion;
261  
262 +
263   #ifdef IS_MPI            
264      const int masterNode = 0;
265      if (worldRank == masterNode) {
266   #endif
267  
268 <      std::ifstream mdFile_(mdFileName.c_str());
268 >      std::ifstream mdFile_;
269 >      mdFile_.open(mdFileName.c_str(), ifstream::in | ifstream::binary);
270        
271        if (mdFile_.fail()) {
272          sprintf(painCave.errMsg,
# Line 276 | Line 293 | namespace OpenMD {
293          painCave.isFatal = 1;
294          simError();
295        }
296 +      
297 +      // found the correct opening string, now try to get the file
298 +      // format version number.
299  
300 +      StringTokenizer tokenizer(line, "=<> \t\n\r");
301 +      std::string fileType = tokenizer.nextToken();
302 +      toUpper(fileType);
303 +
304 +      mdFileVersion = 0;
305 +
306 +      if (fileType == "OPENMD") {
307 +        while (tokenizer.hasMoreTokens()) {
308 +          std::string token(tokenizer.nextToken());
309 +          toUpper(token);
310 +          if (token == "VERSION") {
311 +            mdFileVersion = tokenizer.nextTokenAsInt();
312 +            break;
313 +          }
314 +        }
315 +      }
316 +            
317        //scan through the input stream and find MetaData tag        
318        while(mdFile_.getline(buffer, bufferSize)) {
319          ++lineNo;
# Line 332 | Line 369 | namespace OpenMD {
369      std::stringstream rawMetaDataStream(mdRawData);
370  
371      //parse meta-data file
372 <    Globals* simParams = parseFile(rawMetaDataStream, mdFileName, metaDataBlockStart+1);
372 >    Globals* simParams = parseFile(rawMetaDataStream, mdFileName, mdFileVersion,
373 >                                   metaDataBlockStart + 1);
374      
375      //create the force field
376 <    ForceField * ff = ForceFieldFactory::getInstance()->createForceField(simParams->getForceField());
376 >    ForceField * ff = new ForceField(simParams->getForceField());
377  
378      if (ff == NULL) {
379        sprintf(painCave.errMsg,
# Line 386 | Line 424 | namespace OpenMD {
424      //create the molecules
425      createMolecules(info);
426      
427 +    //find the storage layout
428 +
429 +    int storageLayout = computeStorageLayout(info);
430 +
431      //allocate memory for DataStorage(circular reference, need to
432      //break it)
433 <    info->setSnapshotManager(new SimSnapshotManager(info));
433 >    info->setSnapshotManager(new SimSnapshotManager(info, storageLayout));
434      
435      //set the global index of atoms, rigidbodies and cutoffgroups
436      //(only need to be set once, the global index will never change
# Line 465 | Line 507 | namespace OpenMD {
507      int nGlobalMols = info->getNGlobalMolecules();
508      std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error condition:
509      
510 <    MPI_Comm_size(MPI_COMM_WORLD, &nProcessors);
510 >    nProcessors = MPI::COMM_WORLD.Get_size();
511      
512      if (nProcessors > nGlobalMols) {
513        sprintf(painCave.errMsg,
# Line 503 | Line 545 | namespace OpenMD {
545        nTarget = (int)(precast + 0.5);
546        
547        for(i = 0; i < nGlobalMols; i++) {
548 +
549          done = 0;
550          loops = 0;
551          
# Line 527 | Line 570 | namespace OpenMD {
570            // and be done with it.
571            
572            if (loops > 100) {
573 +
574              sprintf(painCave.errMsg,
575 <                    "I've tried 100 times to assign molecule %d to a "
576 <                    " processor, but can't find a good spot.\n"
577 <                    "I'm assigning it at random to processor %d.\n",
578 <                    i, which_proc);
579 <            
575 >                    "There have been 100 attempts to assign molecule %d to an\n"
576 >                    "\tunderworked processor, but there's no good place to\n"
577 >                    "\tleave it.  OpenMD is assigning it at random to processor %d.\n",
578 >                    i, which_proc);
579 >          
580              painCave.isFatal = 0;
581 +            painCave.severity = OPENMD_INFO;
582              simError();
583              
584              molToProcMap[i] = which_proc;
# Line 578 | Line 623 | namespace OpenMD {
623        }
624        
625        delete myRandom;
626 <      
626 >
627        // Spray out this nonsense to all other processors:
628 <      
584 <      MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
628 >      MPI::COMM_WORLD.Bcast(&molToProcMap[0], nGlobalMols, MPI::INT, 0);
629      } else {
630        
631        // Listen to your marching orders from processor 0:
632 <      
633 <      MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
632 >      MPI::COMM_WORLD.Bcast(&molToProcMap[0], nGlobalMols, MPI::INT, 0);
633 >
634      }
635      
636      info->setMolToProcMap(molToProcMap);
# Line 609 | Line 653 | namespace OpenMD {
653   #endif
654          
655          stampId = info->getMoleculeStampId(i);
656 <        Molecule * mol = molCreator.createMolecule(info->getForceField(), info->getMoleculeStamp(stampId),
657 <                                                   stampId, i, info->getLocalIndexManager());
656 >        Molecule * mol = molCreator.createMolecule(info->getForceField(),
657 >                                                   info->getMoleculeStamp(stampId),
658 >                                                   stampId, i,
659 >                                                   info->getLocalIndexManager());
660          
661          info->addMolecule(mol);
662          
# Line 623 | Line 669 | namespace OpenMD {
669      } //end for(int i=0)  
670    }
671      
672 +  int SimCreator::computeStorageLayout(SimInfo* info) {
673 +
674 +    Globals* simParams = info->getSimParams();
675 +    int nRigidBodies = info->getNGlobalRigidBodies();
676 +    set<AtomType*> atomTypes = info->getSimulatedAtomTypes();
677 +    set<AtomType*>::iterator i;
678 +    bool hasDirectionalAtoms = false;
679 +    bool hasFixedCharge = false;
680 +    bool hasDipoles = false;    
681 +    bool hasQuadrupoles = false;    
682 +    bool hasPolarizable = false;    
683 +    bool hasFluctuatingCharge = false;    
684 +    bool hasMetallic = false;
685 +    int storageLayout = 0;
686 +    storageLayout |= DataStorage::dslPosition;
687 +    storageLayout |= DataStorage::dslVelocity;
688 +    storageLayout |= DataStorage::dslForce;
689 +
690 +    for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
691 +
692 +      DirectionalAdapter da = DirectionalAdapter( (*i) );
693 +      MultipoleAdapter ma = MultipoleAdapter( (*i) );
694 +      EAMAdapter ea = EAMAdapter( (*i) );
695 +      SuttonChenAdapter sca = SuttonChenAdapter( (*i) );
696 +      PolarizableAdapter pa = PolarizableAdapter( (*i) );
697 +      FixedChargeAdapter fca = FixedChargeAdapter( (*i) );
698 +      FluctuatingChargeAdapter fqa = FluctuatingChargeAdapter( (*i) );
699 +
700 +      if (da.isDirectional()){
701 +        hasDirectionalAtoms = true;
702 +      }
703 +      if (ma.isDipole()){
704 +        hasDipoles = true;
705 +      }
706 +      if (ma.isQuadrupole()){
707 +        hasQuadrupoles = true;
708 +      }
709 +      if (ea.isEAM() || sca.isSuttonChen()){
710 +        hasMetallic = true;
711 +      }
712 +      if ( fca.isFixedCharge() ){
713 +        hasFixedCharge = true;
714 +      }
715 +      if ( fqa.isFluctuatingCharge() ){
716 +        hasFluctuatingCharge = true;
717 +      }
718 +      if ( pa.isPolarizable() ){
719 +        hasPolarizable = true;
720 +      }
721 +    }
722 +    
723 +    if (nRigidBodies > 0 || hasDirectionalAtoms) {
724 +      storageLayout |= DataStorage::dslAmat;
725 +      if(storageLayout & DataStorage::dslVelocity) {
726 +        storageLayout |= DataStorage::dslAngularMomentum;
727 +      }
728 +      if (storageLayout & DataStorage::dslForce) {
729 +        storageLayout |= DataStorage::dslTorque;
730 +      }
731 +    }
732 +    if (hasDipoles) {
733 +      storageLayout |= DataStorage::dslDipole;
734 +    }
735 +    if (hasQuadrupoles) {
736 +      storageLayout |= DataStorage::dslQuadrupole;
737 +    }
738 +    if (hasFixedCharge || hasFluctuatingCharge) {
739 +      storageLayout |= DataStorage::dslSkippedCharge;
740 +    }
741 +    if (hasMetallic) {
742 +      storageLayout |= DataStorage::dslDensity;
743 +      storageLayout |= DataStorage::dslFunctional;
744 +      storageLayout |= DataStorage::dslFunctionalDerivative;
745 +    }
746 +    if (hasPolarizable) {
747 +      storageLayout |= DataStorage::dslElectricField;
748 +    }
749 +    if (hasFluctuatingCharge){
750 +      storageLayout |= DataStorage::dslFlucQPosition;
751 +      if(storageLayout & DataStorage::dslVelocity) {
752 +        storageLayout |= DataStorage::dslFlucQVelocity;
753 +      }
754 +      if (storageLayout & DataStorage::dslForce) {
755 +        storageLayout |= DataStorage::dslFlucQForce;
756 +      }
757 +    }
758 +    
759 +    // if the user has asked for them, make sure we've got the memory for the
760 +    // objects defined.
761 +
762 +    if (simParams->getOutputParticlePotential()) {
763 +      storageLayout |= DataStorage::dslParticlePot;
764 +    }
765 +
766 +    if (simParams->havePrintHeatFlux()) {
767 +      if (simParams->getPrintHeatFlux()) {
768 +        storageLayout |= DataStorage::dslParticlePot;
769 +      }
770 +    }
771 +
772 +    if (simParams->getOutputElectricField()) {
773 +      storageLayout |= DataStorage::dslElectricField;
774 +    }
775 +
776 +    if (simParams->getOutputFluctuatingCharges()) {
777 +      storageLayout |= DataStorage::dslFlucQPosition;
778 +      storageLayout |= DataStorage::dslFlucQVelocity;
779 +      storageLayout |= DataStorage::dslFlucQForce;
780 +    }
781 +
782 +    return storageLayout;
783 +  }
784 +
785    void SimCreator::setGlobalIndex(SimInfo *info) {
786      SimInfo::MoleculeIterator mi;
787      Molecule::AtomIterator ai;
# Line 637 | Line 796 | namespace OpenMD {
796      int beginRigidBodyIndex;
797      int beginCutoffGroupIndex;
798      int nGlobalAtoms = info->getNGlobalAtoms();
799 <
641 <    /**@todo fixme */
642 < #ifndef IS_MPI
799 >    int nGlobalRigidBodies = info->getNGlobalRigidBodies();
800      
801      beginAtomIndex = 0;
645    beginRigidBodyIndex = 0;
646    beginCutoffGroupIndex = 0;
647    
648 #else
649    
650    int nproc;
651    int myNode;
652    
653    myNode = worldRank;
654    MPI_Comm_size(MPI_COMM_WORLD, &nproc);
655    
656    std::vector < int > tmpAtomsInProc(nproc, 0);
657    std::vector < int > tmpRigidBodiesInProc(nproc, 0);
658    std::vector < int > tmpCutoffGroupsInProc(nproc, 0);
659    std::vector < int > NumAtomsInProc(nproc, 0);
660    std::vector < int > NumRigidBodiesInProc(nproc, 0);
661    std::vector < int > NumCutoffGroupsInProc(nproc, 0);
662    
663    tmpAtomsInProc[myNode] = info->getNAtoms();
664    tmpRigidBodiesInProc[myNode] = info->getNRigidBodies();
665    tmpCutoffGroupsInProc[myNode] = info->getNCutoffGroups();
666    
667    //do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups
668    MPI_Allreduce(&tmpAtomsInProc[0], &NumAtomsInProc[0], nproc, MPI_INT,
669                  MPI_SUM, MPI_COMM_WORLD);
670    MPI_Allreduce(&tmpRigidBodiesInProc[0], &NumRigidBodiesInProc[0], nproc,
671                  MPI_INT, MPI_SUM, MPI_COMM_WORLD);
672    MPI_Allreduce(&tmpCutoffGroupsInProc[0], &NumCutoffGroupsInProc[0], nproc,
673                  MPI_INT, MPI_SUM, MPI_COMM_WORLD);
674    
675    beginAtomIndex = 0;
676    beginRigidBodyIndex = 0;
677    beginCutoffGroupIndex = 0;
678    
679    for(int i = 0; i < myNode; i++) {
680      beginAtomIndex += NumAtomsInProc[i];
681      beginRigidBodyIndex += NumRigidBodiesInProc[i];
682      beginCutoffGroupIndex += NumCutoffGroupsInProc[i];
683    }
684    
685 #endif
686    
802      //rigidbody's index begins right after atom's
803 <    beginRigidBodyIndex += info->getNGlobalAtoms();
804 <    
805 <    for(mol = info->beginMolecule(mi); mol != NULL;
806 <        mol = info->nextMolecule(mi)) {
803 >    beginRigidBodyIndex = info->getNGlobalAtoms();
804 >    beginCutoffGroupIndex = 0;
805 >
806 >    for(int i = 0; i < info->getNGlobalMolecules(); i++) {
807        
808 <      //local index(index in DataStorge) of atom is important
809 <      for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
810 <        atom->setGlobalIndex(beginAtomIndex++);
808 > #ifdef IS_MPI      
809 >      if (info->getMolToProc(i) == worldRank) {
810 > #endif        
811 >        // stuff to do if I own this molecule
812 >        mol = info->getMoleculeByGlobalIndex(i);
813 >
814 >        //local index(index in DataStorge) of atom is important
815 >        for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
816 >          atom->setGlobalIndex(beginAtomIndex++);
817 >        }
818 >        
819 >        for(rb = mol->beginRigidBody(ri); rb != NULL;
820 >            rb = mol->nextRigidBody(ri)) {
821 >          rb->setGlobalIndex(beginRigidBodyIndex++);
822 >        }
823 >        
824 >        //local index of cutoff group is trivial, it only depends on
825 >        //the order of travesing
826 >        for(cg = mol->beginCutoffGroup(ci); cg != NULL;
827 >            cg = mol->nextCutoffGroup(ci)) {
828 >          cg->setGlobalIndex(beginCutoffGroupIndex++);
829 >        }        
830 >        
831 > #ifdef IS_MPI        
832 >      }  else {
833 >
834 >        // stuff to do if I don't own this molecule
835 >        
836 >        int stampId = info->getMoleculeStampId(i);
837 >        MoleculeStamp* stamp = info->getMoleculeStamp(stampId);
838 >
839 >        beginAtomIndex += stamp->getNAtoms();
840 >        beginRigidBodyIndex += stamp->getNRigidBodies();
841 >        beginCutoffGroupIndex += stamp->getNCutoffGroups() + stamp->getNFreeAtoms();
842        }
843 <      
844 <      for(rb = mol->beginRigidBody(ri); rb != NULL;
845 <          rb = mol->nextRigidBody(ri)) {
846 <        rb->setGlobalIndex(beginRigidBodyIndex++);
701 <      }
702 <      
703 <      //local index of cutoff group is trivial, it only depends on the order of travesing
704 <      for(cg = mol->beginCutoffGroup(ci); cg != NULL;
705 <          cg = mol->nextCutoffGroup(ci)) {
706 <        cg->setGlobalIndex(beginCutoffGroupIndex++);
707 <      }
708 <    }
709 <    
843 > #endif          
844 >
845 >    } //end for(int i=0)  
846 >
847      //fill globalGroupMembership
848      std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
849      for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {        
# Line 726 | Line 863 | namespace OpenMD {
863      // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
864      // docs said we could.
865      std::vector<int> tmpGroupMembership(info->getNGlobalAtoms(), 0);
866 <    MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0], nGlobalAtoms,
867 <                  MPI_INT, MPI_SUM, MPI_COMM_WORLD);
866 >    MPI::COMM_WORLD.Allreduce(&globalGroupMembership[0],
867 >                              &tmpGroupMembership[0], nGlobalAtoms,
868 >                              MPI::INT, MPI::SUM);
869      info->setGlobalGroupMembership(tmpGroupMembership);
870   #else
871      info->setGlobalGroupMembership(globalGroupMembership);
872   #endif
873      
874      //fill molMembership
875 <    std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0);
875 >    std::vector<int> globalMolMembership(info->getNGlobalAtoms() +
876 >                                         info->getNGlobalRigidBodies(), 0);
877      
878 <    for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
878 >    for(mol = info->beginMolecule(mi); mol != NULL;
879 >        mol = info->nextMolecule(mi)) {
880        for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
881          globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
882        }
883 +      for (rb = mol->beginRigidBody(ri); rb != NULL;
884 +           rb = mol->nextRigidBody(ri)) {
885 +        globalMolMembership[rb->getGlobalIndex()] = mol->getGlobalIndex();
886 +      }
887      }
888      
889   #ifdef IS_MPI
890 <    std::vector<int> tmpMolMembership(info->getNGlobalAtoms(), 0);
890 >    std::vector<int> tmpMolMembership(info->getNGlobalAtoms() +
891 >                                      info->getNGlobalRigidBodies(), 0);
892 >    MPI::COMM_WORLD.Allreduce(&globalMolMembership[0], &tmpMolMembership[0],
893 >                              nGlobalAtoms + nGlobalRigidBodies,
894 >                              MPI::INT, MPI::SUM);
895      
748    MPI_Allreduce(&globalMolMembership[0], &tmpMolMembership[0], nGlobalAtoms,
749                  MPI_INT, MPI_SUM, MPI_COMM_WORLD);
750    
896      info->setGlobalMolMembership(tmpMolMembership);
897   #else
898      info->setGlobalMolMembership(globalMolMembership);
# Line 757 | Line 902 | namespace OpenMD {
902      // here the molecules are listed by their global indices.
903  
904      std::vector<int> nIOPerMol(info->getNGlobalMolecules(), 0);
905 <    for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
905 >    for (mol = info->beginMolecule(mi); mol != NULL;
906 >         mol = info->nextMolecule(mi)) {
907        nIOPerMol[mol->getGlobalIndex()] = mol->getNIntegrableObjects();      
908      }
909      
910   #ifdef IS_MPI
911      std::vector<int> numIntegrableObjectsPerMol(info->getNGlobalMolecules(), 0);
912 <    MPI_Allreduce(&nIOPerMol[0], &numIntegrableObjectsPerMol[0],
913 <                  info->getNGlobalMolecules(), MPI_INT, MPI_SUM, MPI_COMM_WORLD);
912 >    MPI::COMM_WORLD.Allreduce(&nIOPerMol[0], &numIntegrableObjectsPerMol[0],
913 >                              info->getNGlobalMolecules(), MPI::INT, MPI::SUM);
914   #else
915      std::vector<int> numIntegrableObjectsPerMol = nIOPerMol;
916   #endif    
# Line 778 | Line 924 | namespace OpenMD {
924      }
925      
926      std::vector<StuntDouble*> IOIndexToIntegrableObject(info->getNGlobalIntegrableObjects(), (StuntDouble*)NULL);
927 <    for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
927 >    for (mol = info->beginMolecule(mi); mol != NULL;
928 >         mol = info->nextMolecule(mi)) {
929        int myGlobalIndex = mol->getGlobalIndex();
930        int globalIO = startingIOIndexForMol[myGlobalIndex];
931 <      for (StuntDouble* integrableObject = mol->beginIntegrableObject(ioi); integrableObject != NULL;
932 <           integrableObject = mol->nextIntegrableObject(ioi)) {
933 <        integrableObject->setGlobalIntegrableObjectIndex(globalIO);
934 <        IOIndexToIntegrableObject[globalIO] = integrableObject;
931 >      for (StuntDouble* sd = mol->beginIntegrableObject(ioi); sd != NULL;
932 >           sd = mol->nextIntegrableObject(ioi)) {
933 >        sd->setGlobalIntegrableObjectIndex(globalIO);
934 >        IOIndexToIntegrableObject[globalIO] = sd;
935          globalIO++;
936        }
937      }
938 <    
938 >      
939      info->setIOIndexToIntegrableObject(IOIndexToIntegrableObject);
940      
941    }

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines