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 1550 by gezelter, Wed Apr 27 21:49:59 2011 UTC vs.
Revision 1601 by gezelter, Thu Aug 4 20:04:35 2011 UTC

# Line 71 | Line 71 | namespace OpenMD {
71      nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
72      nAtoms_(0), nBonds_(0),  nBends_(0), nTorsions_(0), nInversions_(0),
73      nRigidBodies_(0), nIntegrableObjects_(0), nCutoffGroups_(0),
74 <    nConstraints_(0), sman_(NULL), fortranInitialized_(false),
74 >    nConstraints_(0), sman_(NULL), topologyDone_(false),
75      calcBoxDipole_(false), useAtomicVirial_(true) {    
76      
77      MoleculeStamp* molStamp;
# Line 125 | Line 125 | namespace OpenMD {
125      //equal to the total number of atoms minus number of atoms belong to
126      //cutoff group defined in meta-data file plus the number of cutoff
127      //groups defined in meta-data file
128    std::cerr << "nGA = " << nGlobalAtoms_ << "\n";
129    std::cerr << "nCA = " << nCutoffAtoms << "\n";
130    std::cerr << "nG = " << nGroups << "\n";
128  
129      nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
133
134    std::cerr << "nGCG = " << nGlobalCutoffGroups_ << "\n";
130      
131      //every free atom (atom does not belong to rigid bodies) is an
132      //integrable object therefore the total number of integrable objects
# Line 274 | Line 269 | namespace OpenMD {
269   #endif
270      return fdf_;
271    }
272 +  
273 +  unsigned int SimInfo::getNLocalCutoffGroups(){
274 +    int nLocalCutoffAtoms = 0;
275 +    Molecule* mol;
276 +    MoleculeIterator mi;
277 +    CutoffGroup* cg;
278 +    Molecule::CutoffGroupIterator ci;
279      
280 +    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
281 +      
282 +      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
283 +           cg = mol->nextCutoffGroup(ci)) {
284 +        nLocalCutoffAtoms += cg->getNumAtom();
285 +        
286 +      }        
287 +    }
288 +    
289 +    return nAtoms_ - nLocalCutoffAtoms + nCutoffGroups_;
290 +  }
291 +    
292    void SimInfo::calcNdfRaw() {
293      int ndfRaw_local;
294  
# Line 680 | Line 694 | namespace OpenMD {
694      Atom* atom;
695      set<AtomType*> atomTypes;
696      
697 <    for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {      
698 <      for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
697 >    for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
698 >      for(atom = mol->beginAtom(ai); atom != NULL;
699 >          atom = mol->nextAtom(ai)) {
700          atomTypes.insert(atom->getAtomType());
701        }      
702      }    
703 <
703 >    
704   #ifdef IS_MPI
705  
706      // loop over the found atom types on this processor, and add their
707      // numerical idents to a vector:
708 <
708 >    
709      vector<int> foundTypes;
710      set<AtomType*>::iterator i;
711      for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
# Line 699 | Line 714 | namespace OpenMD {
714      // count_local holds the number of found types on this processor
715      int count_local = foundTypes.size();
716  
717 <    // count holds the total number of found types on all processors
703 <    // (some will be redundant with the ones found locally):
704 <    int count;
705 <    MPI::COMM_WORLD.Allreduce(&count_local, &count, 1, MPI::INT, MPI::SUM);
717 >    int nproc = MPI::COMM_WORLD.Get_size();
718  
719 <    // create a vector to hold the globally found types, and resize it:
720 <    vector<int> ftGlobal;
721 <    ftGlobal.resize(count);
722 <    vector<int> counts;
719 >    // we need arrays to hold the counts and displacement vectors for
720 >    // all processors
721 >    vector<int> counts(nproc, 0);
722 >    vector<int> disps(nproc, 0);
723  
724 <    int nproc = MPI::COMM_WORLD.Get_size();
725 <    counts.resize(nproc);
726 <    vector<int> disps;
727 <    disps.resize(nproc);
724 >    // fill the counts array
725 >    MPI::COMM_WORLD.Allgather(&count_local, 1, MPI::INT, &counts[0],
726 >                              1, MPI::INT);
727 >  
728 >    // use the processor counts to compute the displacement array
729 >    disps[0] = 0;    
730 >    int totalCount = counts[0];
731 >    for (int iproc = 1; iproc < nproc; iproc++) {
732 >      disps[iproc] = disps[iproc-1] + counts[iproc-1];
733 >      totalCount += counts[iproc];
734 >    }
735  
736 <    // now spray out the foundTypes to all the other processors:
736 >    // we need a (possibly redundant) set of all found types:
737 >    vector<int> ftGlobal(totalCount);
738      
739 +    // now spray out the foundTypes to all the other processors:    
740      MPI::COMM_WORLD.Allgatherv(&foundTypes[0], count_local, MPI::INT,
741 <                               &ftGlobal[0], &counts[0], &disps[0], MPI::INT);
741 >                               &ftGlobal[0], &counts[0], &disps[0],
742 >                               MPI::INT);
743  
744 +    vector<int>::iterator j;
745 +
746      // foundIdents is a stl set, so inserting an already found ident
747      // will have no effect.
748      set<int> foundIdents;
749 <    vector<int>::iterator j;
749 >
750      for (j = ftGlobal.begin(); j != ftGlobal.end(); ++j)
751        foundIdents.insert((*j));
752      
753      // now iterate over the foundIdents and get the actual atom types
754      // that correspond to these:
755      set<int>::iterator it;
756 <    for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
756 >    for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
757        atomTypes.insert( forceField_->getAtomType((*it)) );
758  
759   #endif
760 <    
760 >
761      return atomTypes;        
762    }
763  
# Line 745 | Line 769 | namespace OpenMD {
769        if ( simParams_->getAccumulateBoxDipole() ) {
770          calcBoxDipole_ = true;      
771        }
772 <
772 >    
773      set<AtomType*>::iterator i;
774      set<AtomType*> atomTypes;
775      atomTypes = getSimulatedAtomTypes();    
# Line 758 | Line 782 | namespace OpenMD {
782        usesMetallic |= (*i)->isMetal();
783        usesDirectional |= (*i)->isDirectional();
784      }
785 <
785 >    
786   #ifdef IS_MPI    
787      int temp;
788      temp = usesDirectional;
789      MPI_Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
790 <
790 >    
791      temp = usesMetallic;
792      MPI_Allreduce(&temp, &usesMetallicAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
793 <
793 >    
794      temp = usesElectrostatic;
795      MPI_Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
796 + #else
797 +
798 +    usesDirectionalAtoms_ = usesDirectional;
799 +    usesMetallicAtoms_ = usesMetallic;
800 +    usesElectrostaticAtoms_ = usesElectrostatic;
801 +
802   #endif
803 <    fInfo_.SIM_uses_PBC = usesPeriodicBoundaries_;    
804 <    fInfo_.SIM_uses_DirectionalAtoms = usesDirectionalAtoms_;
805 <    fInfo_.SIM_uses_MetallicAtoms = usesMetallicAtoms_;
806 <    fInfo_.SIM_requires_SkipCorrection = usesElectrostaticAtoms_;
777 <    fInfo_.SIM_requires_SelfCorrection = usesElectrostaticAtoms_;
778 <    fInfo_.SIM_uses_AtomicVirial = usesAtomicVirial_;
803 >    
804 >    requiresPrepair_ = usesMetallicAtoms_ ? true : false;
805 >    requiresSkipCorrection_ = usesElectrostaticAtoms_ ? true : false;
806 >    requiresSelfCorrection_ = usesElectrostaticAtoms_ ? true : false;    
807    }
808  
809  
# Line 818 | Line 846 | namespace OpenMD {
846    }
847  
848  
849 <  void SimInfo::setupFortran() {
822 <    int isError;
849 >  void SimInfo::prepareTopology() {
850      int nExclude, nOneTwo, nOneThree, nOneFour;
824    vector<int> fortranGlobalGroupMembership;
825    
826    isError = 0;
851  
828    //globalGroupMembership_ is filled by SimCreator    
829    for (int i = 0; i < nGlobalAtoms_; i++) {
830      fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
831    }
832
852      //calculate mass ratio of cutoff group
834    vector<RealType> mfact;
853      SimInfo::MoleculeIterator mi;
854      Molecule* mol;
855      Molecule::CutoffGroupIterator ci;
# Line 840 | Line 858 | namespace OpenMD {
858      Atom* atom;
859      RealType totalMass;
860  
861 <    //to avoid memory reallocation, reserve enough space for mfact
862 <    mfact.reserve(getNCutoffGroups());
861 >    /**
862 >     * The mass factor is the relative mass of an atom to the total
863 >     * mass of the cutoff group it belongs to.  By default, all atoms
864 >     * are their own cutoff groups, and therefore have mass factors of
865 >     * 1.  We need some special handling for massless atoms, which
866 >     * will be treated as carrying the entire mass of the cutoff
867 >     * group.
868 >     */
869 >    massFactors_.clear();
870 >    massFactors_.resize(getNAtoms(), 1.0);
871      
872      for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
873 <      for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
873 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
874 >           cg = mol->nextCutoffGroup(ci)) {
875  
876          totalMass = cg->getMass();
877          for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
878            // Check for massless groups - set mfact to 1 if true
879 <          if (totalMass != 0)
880 <            mfact.push_back(atom->getMass()/totalMass);
879 >          if (totalMass != 0)
880 >            massFactors_[atom->getLocalIndex()] = atom->getMass()/totalMass;
881            else
882 <            mfact.push_back( 1.0 );
882 >            massFactors_[atom->getLocalIndex()] = 1.0;
883          }
884        }      
885      }
# Line 866 | Line 893 | namespace OpenMD {
893          identArray_.push_back(atom->getIdent());
894        }
895      }    
869
870    //fill molMembershipArray
871    //molMembershipArray is filled by SimCreator    
872    vector<int> molMembershipArray(nGlobalAtoms_);
873    for (int i = 0; i < nGlobalAtoms_; i++) {
874      molMembershipArray[i] = globalMolMembership_[i] + 1;
875    }
896      
897 <    //setup fortran simulation
897 >    //scan topology
898  
899      nExclude = excludedInteractions_.getSize();
900      nOneTwo = oneTwoInteractions_.getSize();
# Line 886 | Line 906 | namespace OpenMD {
906      int* oneThreeList = oneThreeInteractions_.getPairList();
907      int* oneFourList = oneFourInteractions_.getPairList();
908  
909 <    //setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray_[0],
890 <    //               &nExclude, excludeList,
891 <    //               &nOneTwo, oneTwoList,
892 <    //               &nOneThree, oneThreeList,
893 <    //               &nOneFour, oneFourList,
894 <    //               &molMembershipArray[0], &mfact[0], &nCutoffGroups_,
895 <    //               &fortranGlobalGroupMembership[0], &isError);
896 <    
897 <    // if( isError ){
898 <    //  
899 <    //  sprintf( painCave.errMsg,
900 <    //         "There was an error setting the simulation information in fortran.\n" );
901 <    //  painCave.isFatal = 1;
902 <    //  painCave.severity = OPENMD_ERROR;
903 <    //  simError();
904 <    //}
905 <    
906 <    
907 <    // sprintf( checkPointMsg,
908 <    //          "succesfully sent the simulation information to fortran.\n");
909 <    
910 <    // errorCheckPoint();
911 <    
912 <    // Setup number of neighbors in neighbor list if present
913 <    //if (simParams_->haveNeighborListNeighbors()) {
914 <    //  int nlistNeighbors = simParams_->getNeighborListNeighbors();
915 <    //  setNeighbors(&nlistNeighbors);
916 <    //}
917 <  
918 < #ifdef IS_MPI    
919 <    // mpiSimData parallelData;
920 <
921 <    //fill up mpiSimData struct
922 <    // parallelData.nMolGlobal = getNGlobalMolecules();
923 <    // parallelData.nMolLocal = getNMolecules();
924 <    // parallelData.nAtomsGlobal = getNGlobalAtoms();
925 <    // parallelData.nAtomsLocal = getNAtoms();
926 <    // parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
927 <    // parallelData.nGroupsLocal = getNCutoffGroups();
928 <    // parallelData.myNode = worldRank;
929 <    // MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
930 <
931 <    //pass mpiSimData struct and index arrays to fortran
932 <    //setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
933 <    //                &localToGlobalAtomIndex[0],  &(parallelData.nGroupsLocal),
934 <    //                &localToGlobalCutoffGroupIndex[0], &isError);
935 <
936 <    // if (isError) {
937 <    //   sprintf(painCave.errMsg,
938 <    //           "mpiRefresh errror: fortran didn't like something we gave it.\n");
939 <    //   painCave.isFatal = 1;
940 <    //   simError();
941 <    // }
942 <
943 <    // sprintf(checkPointMsg, " mpiRefresh successful.\n");
944 <    // errorCheckPoint();
945 < #endif
946 <
947 <    // initFortranFF(&isError);
948 <    // if (isError) {
949 <    //   sprintf(painCave.errMsg,
950 <    //           "initFortranFF errror: fortran didn't like something we gave it.\n");
951 <    //   painCave.isFatal = 1;
952 <    //   simError();
953 <    // }
954 <    // fortranInitialized_ = true;
909 >    topologyDone_ = true;
910    }
911  
912    void SimInfo::addProperty(GenericData* genData) {

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines