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 1532 by gezelter, Wed Dec 29 19:59:21 2010 UTC vs.
Revision 1601 by gezelter, Thu Aug 4 20:04:35 2011 UTC

# Line 54 | Line 54
54   #include "math/Vector3.hpp"
55   #include "primitives/Molecule.hpp"
56   #include "primitives/StuntDouble.hpp"
57 #include "UseTheForce/doForces_interface.h"
58 #include "UseTheForce/DarkSide/neighborLists_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
67 #ifdef IS_MPI
68 #include "UseTheForce/mpiComponentPlan.h"
69 #include "UseTheForce/DarkSide/simParallel_interface.h"
70 #endif
71
64   using namespace std;
65   namespace OpenMD {
66    
# Line 79 | 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 133 | 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 +
129      nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
130      
131      //every free atom (atom does not belong to rigid bodies) is an
# Line 275 | 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 657 | Line 669 | namespace OpenMD {
669    /**
670     * update
671     *
672 <   *  Performs the global checks and variable settings after the objects have been
673 <   *  created.
672 >   *  Performs the global checks and variable settings after the
673 >   *  objects have been created.
674     *
675     */
676 <  void SimInfo::update() {
665 <    
676 >  void SimInfo::update() {  
677      setupSimVariables();
667    setupCutoffs();
668    setupSwitching();
669    setupElectrostatics();
670    setupNeighborlists();
671
672 #ifdef IS_MPI
673    setupFortranParallel();
674 #endif
675    setupFortranSim();
676    fortranInitialized_ = true;
677
678      calcNdf();
679      calcNdfRaw();
680      calcNdfTrans();
681    }
682    
683 +  /**
684 +   * getSimulatedAtomTypes
685 +   *
686 +   * Returns an STL set of AtomType* that are actually present in this
687 +   * simulation.  Must query all processors to assemble this information.
688 +   *
689 +   */
690    set<AtomType*> SimInfo::getSimulatedAtomTypes() {
691      SimInfo::MoleculeIterator mi;
692      Molecule* mol;
# 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 <    return atomTypes;        
704 <  }
703 >    
704 > #ifdef IS_MPI
705  
706 <  /**
707 <   * setupCutoffs
700 <   *
701 <   * Sets the values of cutoffRadius and cutoffMethod
702 <   *
703 <   * cutoffRadius : realType
704 <   *  If the cutoffRadius was explicitly set, use that value.
705 <   *  If the cutoffRadius was not explicitly set:
706 <   *      Are there electrostatic atoms?  Use 12.0 Angstroms.
707 <   *      No electrostatic atoms?  Poll the atom types present in the
708 <   *      simulation for suggested cutoff values (e.g. 2.5 * sigma).
709 <   *      Use the maximum suggested value that was found.
710 <   *
711 <   * cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE, SHIFTED_POTENTIAL)
712 <   *      If cutoffMethod was explicitly set, use that choice.
713 <   *      If cutoffMethod was not explicitly set, use SHIFTED_FORCE
714 <   */
715 <  void SimInfo::setupCutoffs() {
706 >    // loop over the found atom types on this processor, and add their
707 >    // numerical idents to a vector:
708      
709 <    if (simParams_->haveCutoffRadius()) {
710 <      cutoffRadius_ = simParams_->getCutoffRadius();
711 <    } else {      
712 <      if (usesElectrostaticAtoms_) {
721 <        sprintf(painCave.errMsg,
722 <                "SimInfo: No value was set for the cutoffRadius.\n"
723 <                "\tOpenMD will use a default value of 12.0 angstroms"
724 <                "\tfor the cutoffRadius.\n");
725 <        painCave.isFatal = 0;
726 <        painCave.severity = OPENMD_INFO;
727 <        simError();
728 <        cutoffRadius_ = 12.0;
729 <      } else {
730 <        RealType thisCut;
731 <        set<AtomType*>::iterator i;
732 <        set<AtomType*> atomTypes;
733 <        atomTypes = getSimulatedAtomTypes();        
734 <        for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
735 <          thisCut = InteractionManager::Instance()->getSuggestedCutoffRadius((*i));
736 <          cutoffRadius_ = max(thisCut, cutoffRadius_);
737 <        }
738 <        sprintf(painCave.errMsg,
739 <                "SimInfo: No value was set for the cutoffRadius.\n"
740 <                "\tOpenMD will use %lf angstroms.\n",
741 <                cutoffRadius_);
742 <        painCave.isFatal = 0;
743 <        painCave.severity = OPENMD_INFO;
744 <        simError();
745 <      }            
746 <    }
709 >    vector<int> foundTypes;
710 >    set<AtomType*>::iterator i;
711 >    for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
712 >      foundTypes.push_back( (*i)->getIdent() );
713  
714 <    InteractionManager::Instance()->setCutoffRadius(cutoffRadius_);
714 >    // count_local holds the number of found types on this processor
715 >    int count_local = foundTypes.size();
716  
717 <    map<string, CutoffMethod> stringToCutoffMethod;
718 <    stringToCutoffMethod["HARD"] = HARD;
719 <    stringToCutoffMethod["SWITCHING_FUNCTION"] = SWITCHING_FUNCTION;
720 <    stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL;    
721 <    stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE;
717 >    int nproc = MPI::COMM_WORLD.Get_size();
718 >
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 <    if (simParams_->haveCutoffMethod()) {
729 <      string cutMeth = toUpperCopy(simParams_->getCutoffMethod());
730 <      map<string, CutoffMethod>::iterator i;
731 <      i = stringToCutoffMethod.find(cutMeth);
732 <      if (i == stringToCutoffMethod.end()) {
733 <        sprintf(painCave.errMsg,
762 <                "SimInfo: Could not find chosen cutoffMethod %s\n"
763 <                "\tShould be one of: "
764 <                "HARD, SWITCHING_FUNCTION, SHIFTED_POTENTIAL, or SHIFTED_FORCE\n",
765 <                cutMeth.c_str());
766 <        painCave.isFatal = 1;
767 <        painCave.severity = OPENMD_ERROR;
768 <        simError();
769 <      } else {
770 <        cutoffMethod_ = i->second;
771 <      }
772 <    } else {
773 <      sprintf(painCave.errMsg,
774 <              "SimInfo: No value was set for the cutoffMethod.\n"
775 <              "\tOpenMD will use SHIFTED_FORCE.\n");
776 <        painCave.isFatal = 0;
777 <        painCave.severity = OPENMD_INFO;
778 <        simError();
779 <        cutoffMethod_ = SHIFTED_FORCE;        
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 <    InteractionManager::Instance()->setCutoffMethod(cutoffMethod_);
737 <  }
784 <  
785 <  /**
786 <   * setupSwitching
787 <   *
788 <   * Sets the values of switchingRadius and
789 <   *  If the switchingRadius was explicitly set, use that value (but check it)
790 <   *  If the switchingRadius was not explicitly set: use 0.85 * cutoffRadius_
791 <   */
792 <  void SimInfo::setupSwitching() {
736 >    // we need a (possibly redundant) set of all found types:
737 >    vector<int> ftGlobal(totalCount);
738      
739 <    if (simParams_->haveSwitchingRadius()) {
740 <      switchingRadius_ = simParams_->getSwitchingRadius();
741 <      if (switchingRadius_ > cutoffRadius_) {        
742 <        sprintf(painCave.errMsg,
798 <                "SimInfo: switchingRadius (%f) is larger than cutoffRadius(%f)\n",
799 <                switchingRadius_, cutoffRadius_);
800 <        painCave.isFatal = 1;
801 <        painCave.severity = OPENMD_ERROR;
802 <        simError();
803 <      }
804 <    } else {      
805 <      switchingRadius_ = 0.85 * cutoffRadius_;
806 <      sprintf(painCave.errMsg,
807 <              "SimInfo: No value was set for the switchingRadius.\n"
808 <              "\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n"
809 <              "\tswitchingRadius = %f. for this simulation\n", switchingRadius_);
810 <      painCave.isFatal = 0;
811 <      painCave.severity = OPENMD_WARNING;
812 <      simError();
813 <    }          
814 <  
815 <    InteractionManager::Instance()->setSwitchingRadius(switchingRadius_);
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],
742 >                               MPI::INT);
743  
744 <    SwitchingFunctionType ft;
818 <    
819 <    if (simParams_->haveSwitchingFunctionType()) {
820 <      string funcType = simParams_->getSwitchingFunctionType();
821 <      toUpper(funcType);
822 <      if (funcType == "CUBIC") {
823 <        ft = cubic;
824 <      } else {
825 <        if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
826 <          ft = fifth_order_poly;
827 <        } else {
828 <          // throw error        
829 <          sprintf( painCave.errMsg,
830 <                   "SimInfo : Unknown switchingFunctionType. (Input file specified %s .)\n"
831 <                   "\tswitchingFunctionType must be one of: "
832 <                   "\"cubic\" or \"fifth_order_polynomial\".",
833 <                   funcType.c_str() );
834 <          painCave.isFatal = 1;
835 <          painCave.severity = OPENMD_ERROR;
836 <          simError();
837 <        }          
838 <      }
839 <    }
744 >    vector<int>::iterator j;
745  
746 <    InteractionManager::Instance()->setSwitchingFunctionType(ft);
747 <  }
746 >    // foundIdents is a stl set, so inserting an already found ident
747 >    // will have no effect.
748 >    set<int> foundIdents;
749  
750 <  /**
751 <   * setupSkinThickness
752 <   *
753 <   *  If the skinThickness was explicitly set, use that value (but check it)
754 <   *  If the skinThickness was not explicitly set: use 1.0 angstroms
755 <   */
756 <  void SimInfo::setupSkinThickness() {    
757 <    if (simParams_->haveSkinThickness()) {
758 <      skinThickness_ = simParams_->getSkinThickness();
759 <    } else {      
760 <      skinThickness_ = 1.0;
761 <      sprintf(painCave.errMsg,
856 <              "SimInfo Warning: No value was set for the skinThickness.\n"
857 <              "\tOpenMD will use a default value of %f Angstroms\n"
858 <              "\tfor this simulation\n", skinThickness_);
859 <      painCave.isFatal = 0;
860 <      simError();
861 <    }            
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)
757 >      atomTypes.insert( forceField_->getAtomType((*it)) );
758 >
759 > #endif
760 >
761 >    return atomTypes;        
762    }
763  
764 <  void SimInfo::setupSimType() {
764 >  void SimInfo::setupSimVariables() {
765 >    useAtomicVirial_ = simParams_->getUseAtomicVirial();
766 >    // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true
767 >    calcBoxDipole_ = false;
768 >    if ( simParams_->haveAccumulateBoxDipole() )
769 >      if ( simParams_->getAccumulateBoxDipole() ) {
770 >        calcBoxDipole_ = true;      
771 >      }
772 >    
773      set<AtomType*>::iterator i;
774      set<AtomType*> atomTypes;
775 <    atomTypes = getSimulatedAtomTypes();
868 <
869 <    useAtomicVirial_ = simParams_->getUseAtomicVirial();
870 <
775 >    atomTypes = getSimulatedAtomTypes();    
776      int usesElectrostatic = 0;
777      int usesMetallic = 0;
778      int usesDirectional = 0;
# Line 877 | 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_;
896 <    fInfo_.SIM_requires_SelfCorrection = usesElectrostaticAtoms_;
897 <    fInfo_.SIM_uses_AtomicVirial = usesAtomicVirial_;
803 >    
804 >    requiresPrepair_ = usesMetallicAtoms_ ? true : false;
805 >    requiresSkipCorrection_ = usesElectrostaticAtoms_ ? true : false;
806 >    requiresSelfCorrection_ = usesElectrostaticAtoms_ ? true : false;    
807    }
808  
809 <  void SimInfo::setupFortranSim() {
810 <    int isError;
811 <    int nExclude, nOneTwo, nOneThree, nOneFour;
812 <    vector<int> fortranGlobalGroupMembership;
809 >
810 >  vector<int> SimInfo::getGlobalAtomIndices() {
811 >    SimInfo::MoleculeIterator mi;
812 >    Molecule* mol;
813 >    Molecule::AtomIterator ai;
814 >    Atom* atom;
815 >
816 >    vector<int> GlobalAtomIndices(getNAtoms(), 0);
817      
818 <    notifyFortranSkinThickness(&skinThickness_);
818 >    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
819 >      
820 >      for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
821 >        GlobalAtomIndices[atom->getLocalIndex()] = atom->getGlobalIndex();
822 >      }
823 >    }
824 >    return GlobalAtomIndices;
825 >  }
826  
907    int ljsp = cutoffMethod_ == SHIFTED_POTENTIAL ? 1 : 0;
908    int ljsf = cutoffMethod_ == SHIFTED_FORCE ? 1 : 0;
909    notifyFortranCutoffs(&cutoffRadius_, &switchingRadius_, &ljsp, &ljsf);
827  
828 <    isError = 0;
828 >  vector<int> SimInfo::getGlobalGroupIndices() {
829 >    SimInfo::MoleculeIterator mi;
830 >    Molecule* mol;
831 >    Molecule::CutoffGroupIterator ci;
832 >    CutoffGroup* cg;
833  
834 <    //globalGroupMembership_ is filled by SimCreator    
835 <    for (int i = 0; i < nGlobalAtoms_; i++) {
836 <      fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
834 >    vector<int> GlobalGroupIndices;
835 >    
836 >    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
837 >      
838 >      //local index of cutoff group is trivial, it only depends on the
839 >      //order of travesing
840 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
841 >           cg = mol->nextCutoffGroup(ci)) {
842 >        GlobalGroupIndices.push_back(cg->getGlobalIndex());
843 >      }        
844      }
845 +    return GlobalGroupIndices;
846 +  }
847  
848 +
849 +  void SimInfo::prepareTopology() {
850 +    int nExclude, nOneTwo, nOneThree, nOneFour;
851 +
852      //calculate mass ratio of cutoff group
919    vector<RealType> mfact;
853      SimInfo::MoleculeIterator mi;
854      Molecule* mol;
855      Molecule::CutoffGroupIterator ci;
# Line 925 | 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      }
886  
887 <    //fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!)
946 <    vector<int> identArray;
887 >    // Build the identArray_
888  
889 <    //to avoid memory reallocation, reserve enough space identArray
890 <    identArray.reserve(getNAtoms());
950 <    
889 >    identArray_.clear();
890 >    identArray_.reserve(getNAtoms());    
891      for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
892        for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
893 <        identArray.push_back(atom->getIdent());
893 >        identArray_.push_back(atom->getIdent());
894        }
895      }    
956
957    //fill molMembershipArray
958    //molMembershipArray is filled by SimCreator    
959    vector<int> molMembershipArray(nGlobalAtoms_);
960    for (int i = 0; i < nGlobalAtoms_; i++) {
961      molMembershipArray[i] = globalMolMembership_[i] + 1;
962    }
896      
897 <    //setup fortran simulation
897 >    //scan topology
898  
899      nExclude = excludedInteractions_.getSize();
900      nOneTwo = oneTwoInteractions_.getSize();
# Line 973 | Line 906 | namespace OpenMD {
906      int* oneThreeList = oneThreeInteractions_.getPairList();
907      int* oneFourList = oneFourInteractions_.getPairList();
908  
909 <    setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0],
977 <                   &nExclude, excludeList,
978 <                   &nOneTwo, oneTwoList,
979 <                   &nOneThree, oneThreeList,
980 <                   &nOneFour, oneFourList,
981 <                   &molMembershipArray[0], &mfact[0], &nCutoffGroups_,
982 <                   &fortranGlobalGroupMembership[0], &isError);
983 <    
984 <    if( isError ){
985 <      
986 <      sprintf( painCave.errMsg,
987 <               "There was an error setting the simulation information in fortran.\n" );
988 <      painCave.isFatal = 1;
989 <      painCave.severity = OPENMD_ERROR;
990 <      simError();
991 <    }
992 <    
993 <    
994 <    sprintf( checkPointMsg,
995 <             "succesfully sent the simulation information to fortran.\n");
996 <    
997 <    errorCheckPoint();
998 <    
999 <    // Setup number of neighbors in neighbor list if present
1000 <    if (simParams_->haveNeighborListNeighbors()) {
1001 <      int nlistNeighbors = simParams_->getNeighborListNeighbors();
1002 <      setNeighbors(&nlistNeighbors);
1003 <    }
1004 <  
1005 <
1006 <  }
1007 <
1008 <
1009 <  void SimInfo::setupFortranParallel() {
1010 < #ifdef IS_MPI    
1011 <    //SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex
1012 <    vector<int> localToGlobalAtomIndex(getNAtoms(), 0);
1013 <    vector<int> localToGlobalCutoffGroupIndex;
1014 <    SimInfo::MoleculeIterator mi;
1015 <    Molecule::AtomIterator ai;
1016 <    Molecule::CutoffGroupIterator ci;
1017 <    Molecule* mol;
1018 <    Atom* atom;
1019 <    CutoffGroup* cg;
1020 <    mpiSimData parallelData;
1021 <    int isError;
1022 <
1023 <    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
1024 <
1025 <      //local index(index in DataStorge) of atom is important
1026 <      for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
1027 <        localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1;
1028 <      }
1029 <
1030 <      //local index of cutoff group is trivial, it only depends on the order of travesing
1031 <      for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
1032 <        localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1);
1033 <      }        
1034 <        
1035 <    }
1036 <
1037 <    //fill up mpiSimData struct
1038 <    parallelData.nMolGlobal = getNGlobalMolecules();
1039 <    parallelData.nMolLocal = getNMolecules();
1040 <    parallelData.nAtomsGlobal = getNGlobalAtoms();
1041 <    parallelData.nAtomsLocal = getNAtoms();
1042 <    parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
1043 <    parallelData.nGroupsLocal = getNCutoffGroups();
1044 <    parallelData.myNode = worldRank;
1045 <    MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
1046 <
1047 <    //pass mpiSimData struct and index arrays to fortran
1048 <    setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
1049 <                    &localToGlobalAtomIndex[0],  &(parallelData.nGroupsLocal),
1050 <                    &localToGlobalCutoffGroupIndex[0], &isError);
1051 <
1052 <    if (isError) {
1053 <      sprintf(painCave.errMsg,
1054 <              "mpiRefresh errror: fortran didn't like something we gave it.\n");
1055 <      painCave.isFatal = 1;
1056 <      simError();
1057 <    }
1058 <
1059 <    sprintf(checkPointMsg, " mpiRefresh successful.\n");
1060 <    errorCheckPoint();
1061 <
1062 < #endif
909 >    topologyDone_ = true;
910    }
911  
1065
1066  void SimInfo::setupSwitchingFunction() {    
1067
1068  }
1069
1070  void SimInfo::setupAccumulateBoxDipole() {    
1071
1072    // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true
1073    if ( simParams_->haveAccumulateBoxDipole() )
1074      if ( simParams_->getAccumulateBoxDipole() ) {
1075        calcBoxDipole_ = true;
1076      }
1077
1078  }
1079
912    void SimInfo::addProperty(GenericData* genData) {
913      properties_.addProperty(genData);  
914    }
# Line 1111 | Line 943 | namespace OpenMD {
943      Molecule* mol;
944      RigidBody* rb;
945      Atom* atom;
946 +    CutoffGroup* cg;
947      SimInfo::MoleculeIterator mi;
948      Molecule::RigidBodyIterator rbIter;
949 <    Molecule::AtomIterator atomIter;;
949 >    Molecule::AtomIterator atomIter;
950 >    Molecule::CutoffGroupIterator cgIter;
951  
952      for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
953          
# Line 1124 | Line 958 | namespace OpenMD {
958        for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
959          rb->setSnapshotManager(sman_);
960        }
961 +
962 +      for (cg = mol->beginCutoffGroup(cgIter); cg != NULL; cg = mol->nextCutoffGroup(cgIter)) {
963 +        cg->setSnapshotManager(sman_);
964 +      }
965      }    
966      
967    }

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines