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 1549 by gezelter, Wed Apr 27 18:38:15 2011 UTC vs.
Revision 1597 by gezelter, Tue Jul 26 15:49:24 2011 UTC

# Line 54 | Line 54
54   #include "math/Vector3.hpp"
55   #include "primitives/Molecule.hpp"
56   #include "primitives/StuntDouble.hpp"
57 #include "UseTheForce/DarkSide/neighborLists_interface.h"
58 #include "UseTheForce/doForces_interface.h"
57   #include "utils/MemoryUtils.hpp"
58   #include "utils/simError.h"
59   #include "selection/SelectionManager.hpp"
# Line 63 | Line 61
61   #include "UseTheForce/ForceField.hpp"
62   #include "nonbonded/SwitchingFunction.hpp"
63  
66 #ifdef IS_MPI
67 #include "UseTheForce/mpiComponentPlan.h"
68 #include "UseTheForce/DarkSide/simParallel_interface.h"
69 #endif
70
64   using namespace std;
65   namespace OpenMD {
66    
# Line 78 | 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 132 | 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
135    std::cerr << "nGA = " << nGlobalAtoms_ << "\n";
136    std::cerr << "nCA = " << nCutoffAtoms << "\n";
137    std::cerr << "nG = " << nGroups << "\n";
128  
129      nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
140
141    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 281 | 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 687 | 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 706 | 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
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);
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;
723 <
724 <    int nproc = MPI::COMM_WORLD.Get_size();
725 <    counts.resize(nproc);
726 <    vector<int> disps;
727 <    disps.resize(nproc);
728 <
729 <    // now spray out the foundTypes to all the other processors:
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 >    // 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 >    // 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 752 | 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 765 | 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_;
784 <    fInfo_.SIM_requires_SelfCorrection = usesElectrostaticAtoms_;
785 <    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 798 | Line 819 | namespace OpenMD {
819        
820        for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
821          GlobalAtomIndices[atom->getLocalIndex()] = atom->getGlobalIndex();
822 +        cerr << "LI = " << atom->getLocalIndex() << "GAI = " << GlobalAtomIndices[atom->getLocalIndex()] << "\n";
823        }
824      }
825      return GlobalAtomIndices;
# Line 819 | Line 841 | namespace OpenMD {
841        for (cg = mol->beginCutoffGroup(ci); cg != NULL;
842             cg = mol->nextCutoffGroup(ci)) {
843          GlobalGroupIndices.push_back(cg->getGlobalIndex());
844 +        cerr << "LI, GGI = " << GlobalGroupIndices.size() << " " << cg->getGlobalIndex() << "\n";
845        }        
846      }
847      return GlobalGroupIndices;
848    }
849  
850  
851 <  void SimInfo::setupFortran() {
829 <    int isError;
851 >  void SimInfo::prepareTopology() {
852      int nExclude, nOneTwo, nOneThree, nOneFour;
831    vector<int> fortranGlobalGroupMembership;
832    
833    isError = 0;
853  
835    //globalGroupMembership_ is filled by SimCreator    
836    for (int i = 0; i < nGlobalAtoms_; i++) {
837      fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
838    }
839
854      //calculate mass ratio of cutoff group
841    vector<RealType> mfact;
855      SimInfo::MoleculeIterator mi;
856      Molecule* mol;
857      Molecule::CutoffGroupIterator ci;
# Line 847 | Line 860 | namespace OpenMD {
860      Atom* atom;
861      RealType totalMass;
862  
863 <    //to avoid memory reallocation, reserve enough space for mfact
864 <    mfact.reserve(getNCutoffGroups());
863 >    /**
864 >     * The mass factor is the relative mass of an atom to the total
865 >     * mass of the cutoff group it belongs to.  By default, all atoms
866 >     * are their own cutoff groups, and therefore have mass factors of
867 >     * 1.  We need some special handling for massless atoms, which
868 >     * will be treated as carrying the entire mass of the cutoff
869 >     * group.
870 >     */
871 >    massFactors_.clear();
872 >    massFactors_.resize(getNAtoms(), 1.0);
873      
874      for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
875 <      for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
875 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
876 >           cg = mol->nextCutoffGroup(ci)) {
877  
878          totalMass = cg->getMass();
879          for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
880            // Check for massless groups - set mfact to 1 if true
881 <          if (totalMass != 0)
882 <            mfact.push_back(atom->getMass()/totalMass);
881 >          if (totalMass != 0)
882 >            massFactors_[atom->getLocalIndex()] = atom->getMass()/totalMass;
883            else
884 <            mfact.push_back( 1.0 );
884 >            massFactors_[atom->getLocalIndex()] = 1.0;
885          }
886        }      
887      }
# Line 873 | Line 895 | namespace OpenMD {
895          identArray_.push_back(atom->getIdent());
896        }
897      }    
876
877    //fill molMembershipArray
878    //molMembershipArray is filled by SimCreator    
879    vector<int> molMembershipArray(nGlobalAtoms_);
880    for (int i = 0; i < nGlobalAtoms_; i++) {
881      molMembershipArray[i] = globalMolMembership_[i] + 1;
882    }
898      
899 <    //setup fortran simulation
899 >    //scan topology
900  
901      nExclude = excludedInteractions_.getSize();
902      nOneTwo = oneTwoInteractions_.getSize();
# Line 893 | Line 908 | namespace OpenMD {
908      int* oneThreeList = oneThreeInteractions_.getPairList();
909      int* oneFourList = oneFourInteractions_.getPairList();
910  
911 <    setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray_[0],
897 <                   &nExclude, excludeList,
898 <                   &nOneTwo, oneTwoList,
899 <                   &nOneThree, oneThreeList,
900 <                   &nOneFour, oneFourList,
901 <                   &molMembershipArray[0], &mfact[0], &nCutoffGroups_,
902 <                   &fortranGlobalGroupMembership[0], &isError);
903 <    
904 <    if( isError ){
905 <      
906 <      sprintf( painCave.errMsg,
907 <               "There was an error setting the simulation information in fortran.\n" );
908 <      painCave.isFatal = 1;
909 <      painCave.severity = OPENMD_ERROR;
910 <      simError();
911 <    }
912 <    
913 <    
914 <    sprintf( checkPointMsg,
915 <             "succesfully sent the simulation information to fortran.\n");
916 <    
917 <    errorCheckPoint();
918 <    
919 <    // Setup number of neighbors in neighbor list if present
920 <    if (simParams_->haveNeighborListNeighbors()) {
921 <      int nlistNeighbors = simParams_->getNeighborListNeighbors();
922 <      setNeighbors(&nlistNeighbors);
923 <    }
924 <  
925 < #ifdef IS_MPI    
926 <    mpiSimData parallelData;
927 <
928 <    //fill up mpiSimData struct
929 <    parallelData.nMolGlobal = getNGlobalMolecules();
930 <    parallelData.nMolLocal = getNMolecules();
931 <    parallelData.nAtomsGlobal = getNGlobalAtoms();
932 <    parallelData.nAtomsLocal = getNAtoms();
933 <    parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
934 <    parallelData.nGroupsLocal = getNCutoffGroups();
935 <    parallelData.myNode = worldRank;
936 <    MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
937 <
938 <    //pass mpiSimData struct and index arrays to fortran
939 <    //setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
940 <    //                &localToGlobalAtomIndex[0],  &(parallelData.nGroupsLocal),
941 <    //                &localToGlobalCutoffGroupIndex[0], &isError);
942 <
943 <    if (isError) {
944 <      sprintf(painCave.errMsg,
945 <              "mpiRefresh errror: fortran didn't like something we gave it.\n");
946 <      painCave.isFatal = 1;
947 <      simError();
948 <    }
949 <
950 <    sprintf(checkPointMsg, " mpiRefresh successful.\n");
951 <    errorCheckPoint();
952 < #endif
953 <
954 <    initFortranFF(&isError);
955 <    if (isError) {
956 <      sprintf(painCave.errMsg,
957 <              "initFortranFF errror: fortran didn't like something we gave it.\n");
958 <      painCave.isFatal = 1;
959 <      simError();
960 <    }
961 <    fortranInitialized_ = true;
911 >    topologyDone_ = true;
912    }
913  
914    void SimInfo::addProperty(GenericData* genData) {

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines