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

Comparing:
trunk/src/brains/SimCreator.cpp (file contents), Revision 1241 by gezelter, Fri Apr 25 15:14:47 2008 UTC vs.
branches/development/src/brains/SimCreator.cpp (file contents), Revision 1593 by gezelter, Fri Jul 15 21:35:14 2011 UTC

# Line 6 | Line 6
6   * redistribute this software in source and binary code form, provided
7   * that the following conditions are met:
8   *
9 < * 1. Acknowledgement of the program authors must be made in any
10 < *    publication of scientific results based in part on use of the
11 < *    program.  An acceptable form of acknowledgement is citation of
12 < *    the article in which the program was described (Matthew
13 < *    A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 < *    J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 < *    Parallel Simulation Engine for Molecular Dynamics,"
16 < *    J. Comput. Chem. 26, pp. 252-271 (2005))
17 < *
18 < * 2. Redistributions of source code must retain the above copyright
9 > * 1. Redistributions of source code must retain the above copyright
10   *    notice, this list of conditions and the following disclaimer.
11   *
12 < * 3. Redistributions in binary form must reproduce the above copyright
12 > * 2. Redistributions in binary form must reproduce the above copyright
13   *    notice, this list of conditions and the following disclaimer in the
14   *    documentation and/or other materials provided with the
15   *    distribution.
# Line 37 | Line 28
28   * arising out of the use of or inability to use software, even if the
29   * University of Notre Dame has been advised of the possibility of
30   * such damages.
31 + *
32 + * SUPPORT OPEN SCIENCE!  If you use OpenMD or its source code in your
33 + * research, please cite the appropriate papers when you publish your
34 + * work.  Good starting points are:
35 + *                                                                      
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).                        
40   */
41  
42   /**
# Line 79 | Line 79
79   #include "math/ParallelRandNumGen.hpp"
80   #endif
81  
82 < namespace oopse {
82 > namespace OpenMD {
83    
84    Globals* SimCreator::parseFile(std::istream& rawMetaDataStream, const std::string& filename, int startOfMetaDataBlock ){
85      Globals* simParams = NULL;
# Line 108 | Line 108 | namespace oopse {
108        } else {
109          //get stream size
110          commStatus = MPI_Bcast(&streamSize, 1, MPI_LONG, masterNode, MPI_COMM_WORLD);  
111 <                
111 >
112          char* buf = new char[streamSize];
113          assert(buf);
114                  
# Line 116 | Line 116 | namespace oopse {
116          commStatus = MPI_Bcast(buf, streamSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
117                  
118          ppStream.str(buf);
119 <        delete buf;
119 >        delete [] buf;
120  
121        }
122   #endif            
# Line 145 | Line 145 | namespace oopse {
145        treeParser.initializeASTFactory(factory);
146        treeParser.setASTFactory(&factory);
147        simParams = treeParser.walkTree(parser.getAST());
148
148      }
149  
150        
# Line 215 | Line 214 | namespace oopse {
214        painCave.isFatal = 1;
215        simError();        
216      }
217 <    catch (OOPSEException& e) {
217 >    catch (OpenMDException& e) {
218        sprintf(painCave.errMsg,
219                "%s\n",
220                e.getMessage().c_str());
# Line 235 | Line 234 | namespace oopse {
234    
235    SimInfo*  SimCreator::createSim(const std::string & mdFileName,
236                                    bool loadInitCoords) {
237 <
237 >    
238      const int bufferSize = 65535;
239      char buffer[bufferSize];
240      int lineNo = 0;
# Line 263 | Line 262 | namespace oopse {
262        mdFile_.getline(buffer, bufferSize);
263        ++lineNo;
264        std::string line = trimLeftCopy(buffer);
265 <      i = CaseInsensitiveFind(line, "<OOPSE");
266 <      if (i == string::npos) {
265 >      i = CaseInsensitiveFind(line, "<OpenMD");
266 >      if (static_cast<size_t>(i) == string::npos) {
267 >        // try the older file strings to see if that works:
268 >        i = CaseInsensitiveFind(line, "<OOPSE");
269 >      }
270 >      
271 >      if (static_cast<size_t>(i) == string::npos) {
272 >        // still no luck!
273          sprintf(painCave.errMsg,
274 <                "SimCreator: File: %s is not an OOPSE file!\n",
274 >                "SimCreator: File: %s is not a valid OpenMD file!\n",
275                  mdFileName.c_str());
276          painCave.isFatal = 1;
277          simError();
# Line 330 | Line 335 | namespace oopse {
335      Globals* simParams = parseFile(rawMetaDataStream, mdFileName, metaDataBlockStart+1);
336      
337      //create the force field
338 <    ForceField * ff = ForceFieldFactory::getInstance()
339 <      ->createForceField(simParams->getForceField());
335 <    
338 >    ForceField * ff = ForceFieldFactory::getInstance()->createForceField(simParams->getForceField());
339 >
340      if (ff == NULL) {
341        sprintf(painCave.errMsg,
342                "ForceField Factory can not create %s force field\n",
# Line 365 | Line 369 | namespace oopse {
369      }
370      
371      ff->parse(forcefieldFileName);
368    ff->setFortranForceOptions();
372      //create SimInfo
373      SimInfo * info = new SimInfo(ff, simParams);
374  
# Line 383 | Line 386 | namespace oopse {
386      //create the molecules
387      createMolecules(info);
388      
386    
389      //allocate memory for DataStorage(circular reference, need to
390      //break it)
391      info->setSnapshotManager(new SimSnapshotManager(info));
# Line 395 | Line 397 | namespace oopse {
397      //responsibility to LocalIndexManager.
398      setGlobalIndex(info);
399      
400 <    //Although addExcludePairs is called inside SimInfo's addMolecule
400 >    //Although addInteractionPairs is called inside SimInfo's addMolecule
401      //method, at that point atoms don't have the global index yet
402      //(their global index are all initialized to -1).  Therefore we
403 <    //have to call addExcludePairs explicitly here. A way to work
403 >    //have to call addInteractionPairs explicitly here. A way to work
404      //around is that we can determine the beginning global indices of
405      //atoms before they get created.
406      SimInfo::MoleculeIterator mi;
407      Molecule* mol;
408      for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
409 <      info->addExcludePairs(mol);
409 >      info->addInteractionPairs(mol);
410      }
411      
412      if (loadInitCoords)
413        loadCoordinates(info, mdFileName);    
412    
414      return info;
415    }
416    
# Line 473 | Line 474 | namespace oopse {
474                "\tthe number of molecules.  This will not result in a \n"
475                "\tusable division of atoms for force decomposition.\n"
476                "\tEither try a smaller number of processors, or run the\n"
477 <              "\tsingle-processor version of OOPSE.\n", nProcessors, nGlobalMols);
477 >              "\tsingle-processor version of OpenMD.\n", nProcessors, nGlobalMols);
478        
479        painCave.isFatal = 1;
480        simError();
# Line 707 | Line 708 | namespace oopse {
708      }
709      
710      //fill globalGroupMembership
711 <    std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
711 >    std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), -1);
712      for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {        
713        for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
714          
# Line 717 | Line 718 | namespace oopse {
718          
719        }      
720      }
721 <    
721 >  
722   #ifdef IS_MPI    
723      // Since the globalGroupMembership has been zero filled and we've only
724      // poked values into the atoms we know, we can do an Allreduce
725      // to get the full globalGroupMembership array (We think).
726      // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
727      // docs said we could.
728 <    std::vector<int> tmpGroupMembership(nGlobalAtoms, 0);
728 >    std::vector<int> tmpGroupMembership(info->getNGlobalAtoms(), 0);
729      MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0], nGlobalAtoms,
730 <                  MPI_INT, MPI_SUM, MPI_COMM_WORLD);
730 >                  MPI_INT, MPI_MAX, MPI_COMM_WORLD);
731      info->setGlobalGroupMembership(tmpGroupMembership);
732 +
733 +    cerr << "ggm:\n";
734 +    for (int i = 0; i < tmpGroupMembership.size(); i++)
735 +      cerr << "i = " << i << "\t ggm(i) = " << tmpGroupMembership[i] << "\n";
736 +
737   #else
738      info->setGlobalGroupMembership(globalGroupMembership);
739   #endif
# Line 736 | Line 742 | namespace oopse {
742      std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0);
743      
744      for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
739      
745        for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
746          globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
747        }
748      }
749      
750   #ifdef IS_MPI
751 <    std::vector<int> tmpMolMembership(nGlobalAtoms, 0);
751 >    std::vector<int> tmpMolMembership(info->getNGlobalAtoms(), 0);
752      
753      MPI_Allreduce(&globalMolMembership[0], &tmpMolMembership[0], nGlobalAtoms,
754                    MPI_INT, MPI_SUM, MPI_COMM_WORLD);
# Line 769 | Line 774 | namespace oopse {
774      std::vector<int> numIntegrableObjectsPerMol = nIOPerMol;
775   #endif    
776  
777 < std::vector<int> startingIOIndexForMol(info->getNGlobalMolecules());
778 <
779 < int startingIndex = 0;
780 < for (int i = 0; i < info->getNGlobalMolecules(); i++) {
781 <  startingIOIndexForMol[i] = startingIndex;
782 <  startingIndex += numIntegrableObjectsPerMol[i];
783 < }
784 <
785 < std::vector<StuntDouble*> IOIndexToIntegrableObject(info->getNGlobalIntegrableObjects(), (StuntDouble*)NULL);
786 < for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
777 >    std::vector<int> startingIOIndexForMol(info->getNGlobalMolecules());
778 >    
779 >    int startingIndex = 0;
780 >    for (int i = 0; i < info->getNGlobalMolecules(); i++) {
781 >      startingIOIndexForMol[i] = startingIndex;
782 >      startingIndex += numIntegrableObjectsPerMol[i];
783 >    }
784 >    
785 >    std::vector<StuntDouble*> IOIndexToIntegrableObject(info->getNGlobalIntegrableObjects(), (StuntDouble*)NULL);
786 >    for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
787        int myGlobalIndex = mol->getGlobalIndex();
788        int globalIO = startingIOIndexForMol[myGlobalIndex];
789        for (StuntDouble* integrableObject = mol->beginIntegrableObject(ioi); integrableObject != NULL;
790             integrableObject = mol->nextIntegrableObject(ioi)) {
791 <            integrableObject->setGlobalIntegrableObjectIndex(globalIO);
792 <            IOIndexToIntegrableObject[globalIO] = integrableObject;
793 <            globalIO++;
791 >        integrableObject->setGlobalIntegrableObjectIndex(globalIO);
792 >        IOIndexToIntegrableObject[globalIO] = integrableObject;
793 >        globalIO++;
794        }
795      }
796 <
797 <  info->setIOIndexToIntegrableObject(IOIndexToIntegrableObject);
798 <  
796 >    
797 >    info->setIOIndexToIntegrableObject(IOIndexToIntegrableObject);
798 >    
799    }
800    
801    void SimCreator::loadCoordinates(SimInfo* info, const std::string& mdFileName) {
802      Globals* simParams;
803 +
804      simParams = info->getSimParams();
805      
800    
806      DumpReader reader(info, mdFileName);
807      int nframes = reader.getNFrames();
808 <    
808 >
809      if (nframes > 0) {
810        reader.readFrame(nframes - 1);
811      } else {
# Line 811 | Line 816 | int startingIndex = 0;
816        painCave.isFatal = 1;
817        simError();
818      }
814    
819      //copy the current snapshot to previous snapshot
820      info->getSnapshotManager()->advance();
821    }
822    
823 < } //end namespace oopse
823 > } //end namespace OpenMD
824  
825  

Comparing:
trunk/src/brains/SimCreator.cpp (property svn:keywords), Revision 1241 by gezelter, Fri Apr 25 15:14:47 2008 UTC vs.
branches/development/src/brains/SimCreator.cpp (property svn:keywords), Revision 1593 by gezelter, Fri Jul 15 21:35:14 2011 UTC

# Line 0 | Line 1
1 + Author Id Revision Date

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines