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

Comparing branches/development/src/brains/SimInfo.cpp (file contents):
Revision 1544 by gezelter, Fri Mar 18 19:31:52 2011 UTC vs.
Revision 1665 by gezelter, Tue Nov 22 20:38:56 2011 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 54 | Line 55
55   #include "math/Vector3.hpp"
56   #include "primitives/Molecule.hpp"
57   #include "primitives/StuntDouble.hpp"
57 #include "UseTheForce/DarkSide/neighborLists_interface.h"
58 #include "UseTheForce/doForces_interface.h"
58   #include "utils/MemoryUtils.hpp"
59   #include "utils/simError.h"
60   #include "selection/SelectionManager.hpp"
61   #include "io/ForceFieldOptions.hpp"
62   #include "UseTheForce/ForceField.hpp"
63   #include "nonbonded/SwitchingFunction.hpp"
65
64   #ifdef IS_MPI
65 < #include "UseTheForce/mpiComponentPlan.h"
66 < #include "UseTheForce/DarkSide/simParallel_interface.h"
69 < #endif
65 > #include <mpi.h>
66 > #endif
67  
68   using namespace std;
69   namespace OpenMD {
# Line 78 | Line 75 | namespace OpenMD {
75      nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
76      nAtoms_(0), nBonds_(0),  nBends_(0), nTorsions_(0), nInversions_(0),
77      nRigidBodies_(0), nIntegrableObjects_(0), nCutoffGroups_(0),
78 <    nConstraints_(0), sman_(NULL), fortranInitialized_(false),
78 >    nConstraints_(0), sman_(NULL), topologyDone_(false),
79      calcBoxDipole_(false), useAtomicVirial_(true) {    
80      
81      MoleculeStamp* molStamp;
# Line 132 | Line 129 | namespace OpenMD {
129      //equal to the total number of atoms minus number of atoms belong to
130      //cutoff group defined in meta-data file plus the number of cutoff
131      //groups defined in meta-data file
135    std::cerr << "nGA = " << nGlobalAtoms_ << "\n";
136    std::cerr << "nCA = " << nCutoffAtoms << "\n";
137    std::cerr << "nG = " << nGroups << "\n";
132  
133      nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
140
141    std::cerr << "nGCG = " << nGlobalCutoffGroups_ << "\n";
134      
135      //every free atom (atom does not belong to rigid bodies) is an
136      //integrable object therefore the total number of integrable objects
# Line 280 | Line 272 | namespace OpenMD {
272      fdf_ = fdf_local;
273   #endif
274      return fdf_;
275 +  }
276 +  
277 +  unsigned int SimInfo::getNLocalCutoffGroups(){
278 +    int nLocalCutoffAtoms = 0;
279 +    Molecule* mol;
280 +    MoleculeIterator mi;
281 +    CutoffGroup* cg;
282 +    Molecule::CutoffGroupIterator ci;
283 +    
284 +    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
285 +      
286 +      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
287 +           cg = mol->nextCutoffGroup(ci)) {
288 +        nLocalCutoffAtoms += cg->getNumAtom();
289 +        
290 +      }        
291 +    }
292 +    
293 +    return nAtoms_ - nLocalCutoffAtoms + nCutoffGroups_;
294    }
295      
296    void SimInfo::calcNdfRaw() {
# Line 687 | Line 698 | namespace OpenMD {
698      Atom* atom;
699      set<AtomType*> atomTypes;
700      
701 <    for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {      
702 <      for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
701 >    for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
702 >      for(atom = mol->beginAtom(ai); atom != NULL;
703 >          atom = mol->nextAtom(ai)) {
704          atomTypes.insert(atom->getAtomType());
705        }      
706      }    
707 <
707 >    
708   #ifdef IS_MPI
709  
710      // loop over the found atom types on this processor, and add their
711      // numerical idents to a vector:
712 <
712 >    
713      vector<int> foundTypes;
714      set<AtomType*>::iterator i;
715      for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
# Line 706 | Line 718 | namespace OpenMD {
718      // count_local holds the number of found types on this processor
719      int count_local = foundTypes.size();
720  
721 <    // count holds the total number of found types on all processors
710 <    // (some will be redundant with the ones found locally):
711 <    int count;
712 <    MPI::COMM_WORLD.Allreduce(&count_local, &count, 1, MPI::INT, MPI::SUM);
721 >    int nproc = MPI::COMM_WORLD.Get_size();
722  
723 <    // create a vector to hold the globally found types, and resize it:
724 <    vector<int> ftGlobal;
725 <    ftGlobal.resize(count);
726 <    vector<int> counts;
723 >    // we need arrays to hold the counts and displacement vectors for
724 >    // all processors
725 >    vector<int> counts(nproc, 0);
726 >    vector<int> disps(nproc, 0);
727  
728 <    int nproc = MPI::COMM_WORLD.Get_size();
729 <    counts.resize(nproc);
730 <    vector<int> disps;
731 <    disps.resize(nproc);
728 >    // fill the counts array
729 >    MPI::COMM_WORLD.Allgather(&count_local, 1, MPI::INT, &counts[0],
730 >                              1, MPI::INT);
731 >  
732 >    // use the processor counts to compute the displacement array
733 >    disps[0] = 0;    
734 >    int totalCount = counts[0];
735 >    for (int iproc = 1; iproc < nproc; iproc++) {
736 >      disps[iproc] = disps[iproc-1] + counts[iproc-1];
737 >      totalCount += counts[iproc];
738 >    }
739  
740 <    // now spray out the foundTypes to all the other processors:
740 >    // we need a (possibly redundant) set of all found types:
741 >    vector<int> ftGlobal(totalCount);
742      
743 +    // now spray out the foundTypes to all the other processors:    
744      MPI::COMM_WORLD.Allgatherv(&foundTypes[0], count_local, MPI::INT,
745 <                               &ftGlobal[0], &counts[0], &disps[0], MPI::INT);
745 >                               &ftGlobal[0], &counts[0], &disps[0],
746 >                               MPI::INT);
747  
748 +    vector<int>::iterator j;
749 +
750      // foundIdents is a stl set, so inserting an already found ident
751      // will have no effect.
752      set<int> foundIdents;
753 <    vector<int>::iterator j;
753 >
754      for (j = ftGlobal.begin(); j != ftGlobal.end(); ++j)
755        foundIdents.insert((*j));
756      
757      // now iterate over the foundIdents and get the actual atom types
758      // that correspond to these:
759      set<int>::iterator it;
760 <    for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
760 >    for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
761        atomTypes.insert( forceField_->getAtomType((*it)) );
762  
763   #endif
764 <    
764 >
765      return atomTypes;        
766    }
767  
# Line 752 | Line 773 | namespace OpenMD {
773        if ( simParams_->getAccumulateBoxDipole() ) {
774          calcBoxDipole_ = true;      
775        }
776 <
776 >    
777      set<AtomType*>::iterator i;
778      set<AtomType*> atomTypes;
779      atomTypes = getSimulatedAtomTypes();    
# Line 765 | Line 786 | namespace OpenMD {
786        usesMetallic |= (*i)->isMetal();
787        usesDirectional |= (*i)->isDirectional();
788      }
789 <
789 >    
790   #ifdef IS_MPI    
791      int temp;
792      temp = usesDirectional;
793      MPI_Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
794 <
794 >    
795      temp = usesMetallic;
796      MPI_Allreduce(&temp, &usesMetallicAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
797 <
797 >    
798      temp = usesElectrostatic;
799      MPI_Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
800 + #else
801 +
802 +    usesDirectionalAtoms_ = usesDirectional;
803 +    usesMetallicAtoms_ = usesMetallic;
804 +    usesElectrostaticAtoms_ = usesElectrostatic;
805 +
806   #endif
807 <    fInfo_.SIM_uses_PBC = usesPeriodicBoundaries_;    
808 <    fInfo_.SIM_uses_DirectionalAtoms = usesDirectionalAtoms_;
809 <    fInfo_.SIM_uses_MetallicAtoms = usesMetallicAtoms_;
810 <    fInfo_.SIM_requires_SkipCorrection = usesElectrostaticAtoms_;
784 <    fInfo_.SIM_requires_SelfCorrection = usesElectrostaticAtoms_;
785 <    fInfo_.SIM_uses_AtomicVirial = usesAtomicVirial_;
807 >    
808 >    requiresPrepair_ = usesMetallicAtoms_ ? true : false;
809 >    requiresSkipCorrection_ = usesElectrostaticAtoms_ ? true : false;
810 >    requiresSelfCorrection_ = usesElectrostaticAtoms_ ? true : false;    
811    }
812  
813 <  void SimInfo::setupFortran() {
814 <    int isError;
815 <    int nExclude, nOneTwo, nOneThree, nOneFour;
816 <    vector<int> fortranGlobalGroupMembership;
813 >
814 >  vector<int> SimInfo::getGlobalAtomIndices() {
815 >    SimInfo::MoleculeIterator mi;
816 >    Molecule* mol;
817 >    Molecule::AtomIterator ai;
818 >    Atom* atom;
819 >
820 >    vector<int> GlobalAtomIndices(getNAtoms(), 0);
821      
822 <    isError = 0;
822 >    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
823 >      
824 >      for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
825 >        GlobalAtomIndices[atom->getLocalIndex()] = atom->getGlobalIndex();
826 >      }
827 >    }
828 >    return GlobalAtomIndices;
829 >  }
830  
831 <    //globalGroupMembership_ is filled by SimCreator    
832 <    for (int i = 0; i < nGlobalAtoms_; i++) {
833 <      fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
831 >
832 >  vector<int> SimInfo::getGlobalGroupIndices() {
833 >    SimInfo::MoleculeIterator mi;
834 >    Molecule* mol;
835 >    Molecule::CutoffGroupIterator ci;
836 >    CutoffGroup* cg;
837 >
838 >    vector<int> GlobalGroupIndices;
839 >    
840 >    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
841 >      
842 >      //local index of cutoff group is trivial, it only depends on the
843 >      //order of travesing
844 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
845 >           cg = mol->nextCutoffGroup(ci)) {
846 >        GlobalGroupIndices.push_back(cg->getGlobalIndex());
847 >      }        
848      }
849 +    return GlobalGroupIndices;
850 +  }
851  
852 +
853 +  void SimInfo::prepareTopology() {
854 +    int nExclude, nOneTwo, nOneThree, nOneFour;
855 +
856      //calculate mass ratio of cutoff group
801    vector<RealType> mfact;
857      SimInfo::MoleculeIterator mi;
858      Molecule* mol;
859      Molecule::CutoffGroupIterator ci;
# Line 807 | Line 862 | namespace OpenMD {
862      Atom* atom;
863      RealType totalMass;
864  
865 <    //to avoid memory reallocation, reserve enough space for mfact
866 <    mfact.reserve(getNCutoffGroups());
865 >    /**
866 >     * The mass factor is the relative mass of an atom to the total
867 >     * mass of the cutoff group it belongs to.  By default, all atoms
868 >     * are their own cutoff groups, and therefore have mass factors of
869 >     * 1.  We need some special handling for massless atoms, which
870 >     * will be treated as carrying the entire mass of the cutoff
871 >     * group.
872 >     */
873 >    massFactors_.clear();
874 >    massFactors_.resize(getNAtoms(), 1.0);
875      
876      for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
877 <      for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
877 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
878 >           cg = mol->nextCutoffGroup(ci)) {
879  
880          totalMass = cg->getMass();
881          for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
882            // Check for massless groups - set mfact to 1 if true
883 <          if (totalMass != 0)
884 <            mfact.push_back(atom->getMass()/totalMass);
883 >          if (totalMass != 0)
884 >            massFactors_[atom->getLocalIndex()] = atom->getMass()/totalMass;
885            else
886 <            mfact.push_back( 1.0 );
886 >            massFactors_[atom->getLocalIndex()] = 1.0;
887          }
888        }      
889      }
# Line 833 | Line 897 | namespace OpenMD {
897          identArray_.push_back(atom->getIdent());
898        }
899      }    
836
837    //fill molMembershipArray
838    //molMembershipArray is filled by SimCreator    
839    vector<int> molMembershipArray(nGlobalAtoms_);
840    for (int i = 0; i < nGlobalAtoms_; i++) {
841      molMembershipArray[i] = globalMolMembership_[i] + 1;
842    }
900      
901 <    //setup fortran simulation
901 >    //scan topology
902  
903      nExclude = excludedInteractions_.getSize();
904      nOneTwo = oneTwoInteractions_.getSize();
# Line 853 | Line 910 | namespace OpenMD {
910      int* oneThreeList = oneThreeInteractions_.getPairList();
911      int* oneFourList = oneFourInteractions_.getPairList();
912  
913 <    setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0],
857 <                   &nExclude, excludeList,
858 <                   &nOneTwo, oneTwoList,
859 <                   &nOneThree, oneThreeList,
860 <                   &nOneFour, oneFourList,
861 <                   &molMembershipArray[0], &mfact[0], &nCutoffGroups_,
862 <                   &fortranGlobalGroupMembership[0], &isError);
863 <    
864 <    if( isError ){
865 <      
866 <      sprintf( painCave.errMsg,
867 <               "There was an error setting the simulation information in fortran.\n" );
868 <      painCave.isFatal = 1;
869 <      painCave.severity = OPENMD_ERROR;
870 <      simError();
871 <    }
872 <    
873 <    
874 <    sprintf( checkPointMsg,
875 <             "succesfully sent the simulation information to fortran.\n");
876 <    
877 <    errorCheckPoint();
878 <    
879 <    // Setup number of neighbors in neighbor list if present
880 <    if (simParams_->haveNeighborListNeighbors()) {
881 <      int nlistNeighbors = simParams_->getNeighborListNeighbors();
882 <      setNeighbors(&nlistNeighbors);
883 <    }
884 <  
885 < #ifdef IS_MPI    
886 <    //SimInfo is responsible for creating localToGlobalAtomIndex and
887 <    //localToGlobalGroupIndex
888 <    vector<int> localToGlobalAtomIndex(getNAtoms(), 0);
889 <    vector<int> localToGlobalCutoffGroupIndex;
890 <    mpiSimData parallelData;
891 <
892 <    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
893 <
894 <      //local index(index in DataStorge) of atom is important
895 <      for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
896 <        localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1;
897 <      }
898 <
899 <      //local index of cutoff group is trivial, it only depends on the order of travesing
900 <      for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
901 <        localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1);
902 <      }        
903 <        
904 <    }
905 <
906 <    //fill up mpiSimData struct
907 <    parallelData.nMolGlobal = getNGlobalMolecules();
908 <    parallelData.nMolLocal = getNMolecules();
909 <    parallelData.nAtomsGlobal = getNGlobalAtoms();
910 <    parallelData.nAtomsLocal = getNAtoms();
911 <    parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
912 <    parallelData.nGroupsLocal = getNCutoffGroups();
913 <    parallelData.myNode = worldRank;
914 <    MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
915 <
916 <    //pass mpiSimData struct and index arrays to fortran
917 <    setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
918 <                    &localToGlobalAtomIndex[0],  &(parallelData.nGroupsLocal),
919 <                    &localToGlobalCutoffGroupIndex[0], &isError);
920 <
921 <    if (isError) {
922 <      sprintf(painCave.errMsg,
923 <              "mpiRefresh errror: fortran didn't like something we gave it.\n");
924 <      painCave.isFatal = 1;
925 <      simError();
926 <    }
927 <
928 <    sprintf(checkPointMsg, " mpiRefresh successful.\n");
929 <    errorCheckPoint();
930 < #endif
931 <
932 <    initFortranFF(&isError);
933 <    if (isError) {
934 <      sprintf(painCave.errMsg,
935 <              "initFortranFF errror: fortran didn't like something we gave it.\n");
936 <      painCave.isFatal = 1;
937 <      simError();
938 <    }
939 <    fortranInitialized_ = true;
913 >    topologyDone_ = true;
914    }
915  
916    void SimInfo::addProperty(GenericData* genData) {

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines