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 1586 by gezelter, Tue Jun 21 06:34: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 273 | Line 268 | namespace OpenMD {
268      fdf_ = fdf_local;
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() {
# Line 745 | Line 759 | namespace OpenMD {
759        if ( simParams_->getAccumulateBoxDipole() ) {
760          calcBoxDipole_ = true;      
761        }
762 <
762 >    
763      set<AtomType*>::iterator i;
764      set<AtomType*> atomTypes;
765      atomTypes = getSimulatedAtomTypes();    
# Line 758 | Line 772 | namespace OpenMD {
772        usesMetallic |= (*i)->isMetal();
773        usesDirectional |= (*i)->isDirectional();
774      }
775 <
775 >    
776   #ifdef IS_MPI    
777      int temp;
778      temp = usesDirectional;
779      MPI_Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
780 <
780 >    
781      temp = usesMetallic;
782      MPI_Allreduce(&temp, &usesMetallicAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
783 <
783 >    
784      temp = usesElectrostatic;
785      MPI_Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
786 + #else
787 +
788 +    usesDirectionalAtoms_ = usesDirectional;
789 +    usesMetallicAtoms_ = usesMetallic;
790 +    usesElectrostaticAtoms_ = usesElectrostatic;
791 +
792   #endif
793 <    fInfo_.SIM_uses_PBC = usesPeriodicBoundaries_;    
794 <    fInfo_.SIM_uses_DirectionalAtoms = usesDirectionalAtoms_;
795 <    fInfo_.SIM_uses_MetallicAtoms = usesMetallicAtoms_;
796 <    fInfo_.SIM_requires_SkipCorrection = usesElectrostaticAtoms_;
777 <    fInfo_.SIM_requires_SelfCorrection = usesElectrostaticAtoms_;
778 <    fInfo_.SIM_uses_AtomicVirial = usesAtomicVirial_;
793 >    
794 >    requiresPrepair_ = usesMetallicAtoms_ ? true : false;
795 >    requiresSkipCorrection_ = usesElectrostaticAtoms_ ? true : false;
796 >    requiresSelfCorrection_ = usesElectrostaticAtoms_ ? true : false;    
797    }
798  
799  
# Line 818 | Line 836 | namespace OpenMD {
836    }
837  
838  
839 <  void SimInfo::setupFortran() {
822 <    int isError;
839 >  void SimInfo::prepareTopology() {
840      int nExclude, nOneTwo, nOneThree, nOneFour;
824    vector<int> fortranGlobalGroupMembership;
825    
826    isError = 0;
841  
828    //globalGroupMembership_ is filled by SimCreator    
829    for (int i = 0; i < nGlobalAtoms_; i++) {
830      fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
831    }
832
842      //calculate mass ratio of cutoff group
834    vector<RealType> mfact;
843      SimInfo::MoleculeIterator mi;
844      Molecule* mol;
845      Molecule::CutoffGroupIterator ci;
# Line 840 | Line 848 | namespace OpenMD {
848      Atom* atom;
849      RealType totalMass;
850  
851 <    //to avoid memory reallocation, reserve enough space for mfact
852 <    mfact.reserve(getNCutoffGroups());
851 >    /**
852 >     * The mass factor is the relative mass of an atom to the total
853 >     * mass of the cutoff group it belongs to.  By default, all atoms
854 >     * are their own cutoff groups, and therefore have mass factors of
855 >     * 1.  We need some special handling for massless atoms, which
856 >     * will be treated as carrying the entire mass of the cutoff
857 >     * group.
858 >     */
859 >    massFactors_.clear();
860 >    massFactors_.resize(getNAtoms(), 1.0);
861      
862 +    cerr << "mfs in si = " << massFactors_.size() << "\n";
863      for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
864 <      for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
864 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
865 >           cg = mol->nextCutoffGroup(ci)) {
866  
867          totalMass = cg->getMass();
868          for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
869            // Check for massless groups - set mfact to 1 if true
870 <          if (totalMass != 0)
871 <            mfact.push_back(atom->getMass()/totalMass);
870 >          if (totalMass != 0)
871 >            massFactors_[atom->getLocalIndex()] = atom->getMass()/totalMass;
872            else
873 <            mfact.push_back( 1.0 );
873 >            massFactors_[atom->getLocalIndex()] = 1.0;
874          }
875        }      
876      }
# Line 866 | Line 884 | namespace OpenMD {
884          identArray_.push_back(atom->getIdent());
885        }
886      }    
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    }
887      
888 <    //setup fortran simulation
888 >    //scan topology
889  
890      nExclude = excludedInteractions_.getSize();
891      nOneTwo = oneTwoInteractions_.getSize();
# Line 886 | Line 897 | namespace OpenMD {
897      int* oneThreeList = oneThreeInteractions_.getPairList();
898      int* oneFourList = oneFourInteractions_.getPairList();
899  
900 <    //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;
900 >    topologyDone_ = true;
901    }
902  
903    void SimInfo::addProperty(GenericData* genData) {

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines