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 1767 by gezelter, Fri Jul 6 22:01:58 2012 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/doForces_interface.h"
58 #include "UseTheForce/DarkSide/neighborLists_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"
62 > #include "brains/ForceField.hpp"
63   #include "nonbonded/SwitchingFunction.hpp"
65
66
64   #ifdef IS_MPI
65 < #include "UseTheForce/mpiComponentPlan.h"
66 < #include "UseTheForce/DarkSide/simParallel_interface.h"
70 < #endif
65 > #include <mpi.h>
66 > #endif
67  
68   using namespace std;
69   namespace OpenMD {
# Line 76 | Line 72 | namespace OpenMD {
72      forceField_(ff), simParams_(simParams),
73      ndf_(0), fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
74      nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
75 <    nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
75 >    nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0), nGlobalFluctuatingCharges_(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), nFluctuatingCharges_(0), sman_(NULL), topologyDone_(false),
79      calcBoxDipole_(false), useAtomicVirial_(true) {    
80      
81      MoleculeStamp* molStamp;
# Line 133 | 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
132 +
133      nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
134      
135      //every free atom (atom does not belong to rigid bodies) is an
# Line 228 | Line 225 | namespace OpenMD {
225  
226  
227    void SimInfo::calcNdf() {
228 <    int ndf_local;
228 >    int ndf_local, nfq_local;
229      MoleculeIterator i;
230      vector<StuntDouble*>::iterator j;
231 +    vector<Atom*>::iterator k;
232 +
233      Molecule* mol;
234      StuntDouble* integrableObject;
235 +    Atom* atom;
236  
237      ndf_local = 0;
238 +    nfq_local = 0;
239      
240      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
241        for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
# Line 249 | Line 250 | namespace OpenMD {
250              ndf_local += 3;
251            }
252          }
252            
253        }
254 +      for (atom = mol->beginFluctuatingCharge(k); atom != NULL;
255 +           atom = mol->nextFluctuatingCharge(k)) {
256 +        if (atom->isFluctuatingCharge()) {
257 +          nfq_local++;
258 +        }
259 +      }
260      }
261      
262 +    ndfLocal_ = ndf_local;
263 +
264      // n_constraints is local, so subtract them on each processor
265      ndf_local -= nConstraints_;
266  
267   #ifdef IS_MPI
268      MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
269 +    MPI_Allreduce(&nfq_local,&nGlobalFluctuatingCharges_,1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
270   #else
271      ndf_ = ndf_local;
272 +    nGlobalFluctuatingCharges_ = nfq_local;
273   #endif
274  
275      // nZconstraints_ is global, as are the 3 COM translations for the
# Line 276 | Line 286 | namespace OpenMD {
286   #endif
287      return fdf_;
288    }
289 +  
290 +  unsigned int SimInfo::getNLocalCutoffGroups(){
291 +    int nLocalCutoffAtoms = 0;
292 +    Molecule* mol;
293 +    MoleculeIterator mi;
294 +    CutoffGroup* cg;
295 +    Molecule::CutoffGroupIterator ci;
296      
297 +    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
298 +      
299 +      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
300 +           cg = mol->nextCutoffGroup(ci)) {
301 +        nLocalCutoffAtoms += cg->getNumAtom();
302 +        
303 +      }        
304 +    }
305 +    
306 +    return nAtoms_ - nLocalCutoffAtoms + nCutoffGroups_;
307 +  }
308 +    
309    void SimInfo::calcNdfRaw() {
310      int ndfRaw_local;
311  
# Line 657 | Line 686 | namespace OpenMD {
686    /**
687     * update
688     *
689 <   *  Performs the global checks and variable settings after the objects have been
690 <   *  created.
689 >   *  Performs the global checks and variable settings after the
690 >   *  objects have been created.
691     *
692     */
693 <  void SimInfo::update() {
665 <    
693 >  void SimInfo::update() {  
694      setupSimVariables();
667    setupCutoffs();
668    setupSwitching();
669    setupElectrostatics();
670    setupNeighborlists();
671
672 #ifdef IS_MPI
673    setupFortranParallel();
674 #endif
675    setupFortranSim();
676    fortranInitialized_ = true;
677
695      calcNdf();
696      calcNdfRaw();
697      calcNdfTrans();
698    }
699    
700 +  /**
701 +   * getSimulatedAtomTypes
702 +   *
703 +   * Returns an STL set of AtomType* that are actually present in this
704 +   * simulation.  Must query all processors to assemble this information.
705 +   *
706 +   */
707    set<AtomType*> SimInfo::getSimulatedAtomTypes() {
708      SimInfo::MoleculeIterator mi;
709      Molecule* mol;
# Line 687 | Line 711 | namespace OpenMD {
711      Atom* atom;
712      set<AtomType*> atomTypes;
713      
714 <    for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {      
715 <      for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
714 >    for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
715 >      for(atom = mol->beginAtom(ai); atom != NULL;
716 >          atom = mol->nextAtom(ai)) {
717          atomTypes.insert(atom->getAtomType());
718        }      
719      }    
720 <    return atomTypes;        
721 <  }
720 >    
721 > #ifdef IS_MPI
722  
723 <  /**
724 <   * 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() {
723 >    // loop over the found atom types on this processor, and add their
724 >    // numerical idents to a vector:
725      
726 <    if (simParams_->haveCutoffRadius()) {
727 <      cutoffRadius_ = simParams_->getCutoffRadius();
728 <    } else {      
729 <      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 <    }
726 >    vector<int> foundTypes;
727 >    set<AtomType*>::iterator i;
728 >    for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
729 >      foundTypes.push_back( (*i)->getIdent() );
730  
731 <    InteractionManager::Instance()->setCutoffRadius(cutoffRadius_);
731 >    // count_local holds the number of found types on this processor
732 >    int count_local = foundTypes.size();
733  
734 <    map<string, CutoffMethod> stringToCutoffMethod;
751 <    stringToCutoffMethod["HARD"] = HARD;
752 <    stringToCutoffMethod["SWITCHING_FUNCTION"] = SWITCHING_FUNCTION;
753 <    stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL;    
754 <    stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE;
755 <  
756 <    if (simParams_->haveCutoffMethod()) {
757 <      string cutMeth = toUpperCopy(simParams_->getCutoffMethod());
758 <      map<string, CutoffMethod>::iterator i;
759 <      i = stringToCutoffMethod.find(cutMeth);
760 <      if (i == stringToCutoffMethod.end()) {
761 <        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;        
780 <    }
734 >    int nproc = MPI::COMM_WORLD.Get_size();
735  
736 <    InteractionManager::Instance()->setCutoffMethod(cutoffMethod_);
737 <  }
738 <  
739 <  /**
740 <   * setupSwitching
741 <   *
742 <   * Sets the values of switchingRadius and
743 <   *  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() {
793 <    
794 <    if (simParams_->haveSwitchingRadius()) {
795 <      switchingRadius_ = simParams_->getSwitchingRadius();
796 <      if (switchingRadius_ > cutoffRadius_) {        
797 <        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 <    }          
736 >    // we need arrays to hold the counts and displacement vectors for
737 >    // all processors
738 >    vector<int> counts(nproc, 0);
739 >    vector<int> disps(nproc, 0);
740 >
741 >    // fill the counts array
742 >    MPI::COMM_WORLD.Allgather(&count_local, 1, MPI::INT, &counts[0],
743 >                              1, MPI::INT);
744    
745 <    InteractionManager::Instance()->setSwitchingRadius(switchingRadius_);
745 >    // use the processor counts to compute the displacement array
746 >    disps[0] = 0;    
747 >    int totalCount = counts[0];
748 >    for (int iproc = 1; iproc < nproc; iproc++) {
749 >      disps[iproc] = disps[iproc-1] + counts[iproc-1];
750 >      totalCount += counts[iproc];
751 >    }
752  
753 <    SwitchingFunctionType ft;
753 >    // we need a (possibly redundant) set of all found types:
754 >    vector<int> ftGlobal(totalCount);
755      
756 <    if (simParams_->haveSwitchingFunctionType()) {
757 <      string funcType = simParams_->getSwitchingFunctionType();
758 <      toUpper(funcType);
759 <      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 <    }
756 >    // now spray out the foundTypes to all the other processors:    
757 >    MPI::COMM_WORLD.Allgatherv(&foundTypes[0], count_local, MPI::INT,
758 >                               &ftGlobal[0], &counts[0], &disps[0],
759 >                               MPI::INT);
760  
761 <    InteractionManager::Instance()->setSwitchingFunctionType(ft);
842 <  }
761 >    vector<int>::iterator j;
762  
763 <  /**
764 <   * setupSkinThickness
765 <   *
766 <   *  If the skinThickness was explicitly set, use that value (but check it)
767 <   *  If the skinThickness was not explicitly set: use 1.0 angstroms
768 <   */
769 <  void SimInfo::setupSkinThickness() {    
770 <    if (simParams_->haveSkinThickness()) {
771 <      skinThickness_ = simParams_->getSkinThickness();
772 <    } else {      
773 <      skinThickness_ = 1.0;
774 <      sprintf(painCave.errMsg,
775 <              "SimInfo Warning: No value was set for the skinThickness.\n"
776 <              "\tOpenMD will use a default value of %f Angstroms\n"
777 <              "\tfor this simulation\n", skinThickness_);
778 <      painCave.isFatal = 0;
860 <      simError();
861 <    }            
763 >    // foundIdents is a stl set, so inserting an already found ident
764 >    // will have no effect.
765 >    set<int> foundIdents;
766 >
767 >    for (j = ftGlobal.begin(); j != ftGlobal.end(); ++j)
768 >      foundIdents.insert((*j));
769 >    
770 >    // now iterate over the foundIdents and get the actual atom types
771 >    // that correspond to these:
772 >    set<int>::iterator it;
773 >    for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
774 >      atomTypes.insert( forceField_->getAtomType((*it)) );
775 >
776 > #endif
777 >
778 >    return atomTypes;        
779    }
780  
781 <  void SimInfo::setupSimType() {
781 >  void SimInfo::setupSimVariables() {
782 >    useAtomicVirial_ = simParams_->getUseAtomicVirial();
783 >    // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true
784 >    calcBoxDipole_ = false;
785 >    if ( simParams_->haveAccumulateBoxDipole() )
786 >      if ( simParams_->getAccumulateBoxDipole() ) {
787 >        calcBoxDipole_ = true;      
788 >      }
789 >    
790      set<AtomType*>::iterator i;
791      set<AtomType*> atomTypes;
792 <    atomTypes = getSimulatedAtomTypes();
793 <
794 <    useAtomicVirial_ = simParams_->getUseAtomicVirial();
795 <
796 <    int usesElectrostatic = 0;
872 <    int usesMetallic = 0;
873 <    int usesDirectional = 0;
792 >    atomTypes = getSimulatedAtomTypes();    
793 >    bool usesElectrostatic = false;
794 >    bool usesMetallic = false;
795 >    bool usesDirectional = false;
796 >    bool usesFluctuatingCharges =  false;
797      //loop over all of the atom types
798      for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
799        usesElectrostatic |= (*i)->isElectrostatic();
800        usesMetallic |= (*i)->isMetal();
801        usesDirectional |= (*i)->isDirectional();
802 +      usesFluctuatingCharges |= (*i)->isFluctuatingCharge();
803      }
804  
805 < #ifdef IS_MPI    
806 <    int temp;
805 > #ifdef IS_MPI
806 >    bool temp;
807      temp = usesDirectional;
808 <    MPI_Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
809 <
808 >    MPI::COMM_WORLD.Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI::BOOL,
809 >                              MPI::LOR);
810 >        
811      temp = usesMetallic;
812 <    MPI_Allreduce(&temp, &usesMetallicAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
813 <
812 >    MPI::COMM_WORLD.Allreduce(&temp, &usesMetallicAtoms_, 1, MPI::BOOL,
813 >                              MPI::LOR);
814 >    
815      temp = usesElectrostatic;
816 <    MPI_Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
816 >    MPI::COMM_WORLD.Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI::BOOL,
817 >                              MPI::LOR);
818 >
819 >    temp = usesFluctuatingCharges;
820 >    MPI::COMM_WORLD.Allreduce(&temp, &usesFluctuatingCharges_, 1, MPI::BOOL,
821 >                              MPI::LOR);
822 > #else
823 >
824 >    usesDirectionalAtoms_ = usesDirectional;
825 >    usesMetallicAtoms_ = usesMetallic;
826 >    usesElectrostaticAtoms_ = usesElectrostatic;
827 >    usesFluctuatingCharges_ = usesFluctuatingCharges;
828 >
829   #endif
830 <    fInfo_.SIM_uses_PBC = usesPeriodicBoundaries_;    
831 <    fInfo_.SIM_uses_DirectionalAtoms = usesDirectionalAtoms_;
832 <    fInfo_.SIM_uses_MetallicAtoms = usesMetallicAtoms_;
833 <    fInfo_.SIM_requires_SkipCorrection = usesElectrostaticAtoms_;
896 <    fInfo_.SIM_requires_SelfCorrection = usesElectrostaticAtoms_;
897 <    fInfo_.SIM_uses_AtomicVirial = usesAtomicVirial_;
830 >    
831 >    requiresPrepair_ = usesMetallicAtoms_ ? true : false;
832 >    requiresSkipCorrection_ = usesElectrostaticAtoms_ ? true : false;
833 >    requiresSelfCorrection_ = usesElectrostaticAtoms_ ? true : false;    
834    }
835  
836 <  void SimInfo::setupFortranSim() {
837 <    int isError;
838 <    int nExclude, nOneTwo, nOneThree, nOneFour;
839 <    vector<int> fortranGlobalGroupMembership;
836 >
837 >  vector<int> SimInfo::getGlobalAtomIndices() {
838 >    SimInfo::MoleculeIterator mi;
839 >    Molecule* mol;
840 >    Molecule::AtomIterator ai;
841 >    Atom* atom;
842 >
843 >    vector<int> GlobalAtomIndices(getNAtoms(), 0);
844      
845 <    notifyFortranSkinThickness(&skinThickness_);
845 >    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
846 >      
847 >      for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
848 >        GlobalAtomIndices[atom->getLocalIndex()] = atom->getGlobalIndex();
849 >      }
850 >    }
851 >    return GlobalAtomIndices;
852 >  }
853  
907    int ljsp = cutoffMethod_ == SHIFTED_POTENTIAL ? 1 : 0;
908    int ljsf = cutoffMethod_ == SHIFTED_FORCE ? 1 : 0;
909    notifyFortranCutoffs(&cutoffRadius_, &switchingRadius_, &ljsp, &ljsf);
854  
855 <    isError = 0;
855 >  vector<int> SimInfo::getGlobalGroupIndices() {
856 >    SimInfo::MoleculeIterator mi;
857 >    Molecule* mol;
858 >    Molecule::CutoffGroupIterator ci;
859 >    CutoffGroup* cg;
860  
861 <    //globalGroupMembership_ is filled by SimCreator    
862 <    for (int i = 0; i < nGlobalAtoms_; i++) {
863 <      fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
861 >    vector<int> GlobalGroupIndices;
862 >    
863 >    for (mol = beginMolecule(mi); mol != NULL; mol  = nextMolecule(mi)) {
864 >      
865 >      //local index of cutoff group is trivial, it only depends on the
866 >      //order of travesing
867 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
868 >           cg = mol->nextCutoffGroup(ci)) {
869 >        GlobalGroupIndices.push_back(cg->getGlobalIndex());
870 >      }        
871      }
872 +    return GlobalGroupIndices;
873 +  }
874  
875 +
876 +  void SimInfo::prepareTopology() {
877 +    int nExclude, nOneTwo, nOneThree, nOneFour;
878 +
879      //calculate mass ratio of cutoff group
919    vector<RealType> mfact;
880      SimInfo::MoleculeIterator mi;
881      Molecule* mol;
882      Molecule::CutoffGroupIterator ci;
# Line 925 | Line 885 | namespace OpenMD {
885      Atom* atom;
886      RealType totalMass;
887  
888 <    //to avoid memory reallocation, reserve enough space for mfact
889 <    mfact.reserve(getNCutoffGroups());
888 >    /**
889 >     * The mass factor is the relative mass of an atom to the total
890 >     * mass of the cutoff group it belongs to.  By default, all atoms
891 >     * are their own cutoff groups, and therefore have mass factors of
892 >     * 1.  We need some special handling for massless atoms, which
893 >     * will be treated as carrying the entire mass of the cutoff
894 >     * group.
895 >     */
896 >    massFactors_.clear();
897 >    massFactors_.resize(getNAtoms(), 1.0);
898      
899      for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
900 <      for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
900 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
901 >           cg = mol->nextCutoffGroup(ci)) {
902  
903          totalMass = cg->getMass();
904          for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
905            // Check for massless groups - set mfact to 1 if true
906 <          if (totalMass != 0)
907 <            mfact.push_back(atom->getMass()/totalMass);
906 >          if (totalMass != 0)
907 >            massFactors_[atom->getLocalIndex()] = atom->getMass()/totalMass;
908            else
909 <            mfact.push_back( 1.0 );
909 >            massFactors_[atom->getLocalIndex()] = 1.0;
910          }
911        }      
912      }
913  
914 <    //fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!)
946 <    vector<int> identArray;
914 >    // Build the identArray_
915  
916 <    //to avoid memory reallocation, reserve enough space identArray
917 <    identArray.reserve(getNAtoms());
950 <    
916 >    identArray_.clear();
917 >    identArray_.reserve(getNAtoms());    
918      for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
919        for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
920 <        identArray.push_back(atom->getIdent());
920 >        identArray_.push_back(atom->getIdent());
921        }
922      }    
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    }
923      
924 <    //setup fortran simulation
924 >    //scan topology
925  
926      nExclude = excludedInteractions_.getSize();
927      nOneTwo = oneTwoInteractions_.getSize();
# Line 973 | Line 933 | namespace OpenMD {
933      int* oneThreeList = oneThreeInteractions_.getPairList();
934      int* oneFourList = oneFourInteractions_.getPairList();
935  
936 <    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 <
936 >    topologyDone_ = true;
937    }
938  
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
1063  }
1064
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
939    void SimInfo::addProperty(GenericData* genData) {
940      properties_.addProperty(genData);  
941    }
# Line 1111 | Line 970 | namespace OpenMD {
970      Molecule* mol;
971      RigidBody* rb;
972      Atom* atom;
973 +    CutoffGroup* cg;
974      SimInfo::MoleculeIterator mi;
975      Molecule::RigidBodyIterator rbIter;
976 <    Molecule::AtomIterator atomIter;;
976 >    Molecule::AtomIterator atomIter;
977 >    Molecule::CutoffGroupIterator cgIter;
978  
979      for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
980          
# Line 1124 | Line 985 | namespace OpenMD {
985        for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
986          rb->setSnapshotManager(sman_);
987        }
988 +
989 +      for (cg = mol->beginCutoffGroup(cgIter); cg != NULL; cg = mol->nextCutoffGroup(cgIter)) {
990 +        cg->setSnapshotManager(sman_);
991 +      }
992      }    
993      
994    }
995  
1131  Vector3d SimInfo::getComVel(){
1132    SimInfo::MoleculeIterator i;
1133    Molecule* mol;
996  
1135    Vector3d comVel(0.0);
1136    RealType totalMass = 0.0;
1137    
1138
1139    for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1140      RealType mass = mol->getMass();
1141      totalMass += mass;
1142      comVel += mass * mol->getComVel();
1143    }  
1144
1145 #ifdef IS_MPI
1146    RealType tmpMass = totalMass;
1147    Vector3d tmpComVel(comVel);    
1148    MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1149    MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1150 #endif
1151
1152    comVel /= totalMass;
1153
1154    return comVel;
1155  }
1156
1157  Vector3d SimInfo::getCom(){
1158    SimInfo::MoleculeIterator i;
1159    Molecule* mol;
1160
1161    Vector3d com(0.0);
1162    RealType totalMass = 0.0;
1163    
1164    for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1165      RealType mass = mol->getMass();
1166      totalMass += mass;
1167      com += mass * mol->getCom();
1168    }  
1169
1170 #ifdef IS_MPI
1171    RealType tmpMass = totalMass;
1172    Vector3d tmpCom(com);    
1173    MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1174    MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1175 #endif
1176
1177    com /= totalMass;
1178
1179    return com;
1180
1181  }        
1182
997    ostream& operator <<(ostream& o, SimInfo& info) {
998  
999      return o;
1000    }
1001    
1002 <  
1189 <   /*
1190 <   Returns center of mass and center of mass velocity in one function call.
1191 <   */
1192 <  
1193 <   void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1194 <      SimInfo::MoleculeIterator i;
1195 <      Molecule* mol;
1196 <      
1197 <    
1198 <      RealType totalMass = 0.0;
1199 <    
1200 <
1201 <      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1202 <         RealType mass = mol->getMass();
1203 <         totalMass += mass;
1204 <         com += mass * mol->getCom();
1205 <         comVel += mass * mol->getComVel();          
1206 <      }  
1207 <      
1208 < #ifdef IS_MPI
1209 <      RealType tmpMass = totalMass;
1210 <      Vector3d tmpCom(com);  
1211 <      Vector3d tmpComVel(comVel);
1212 <      MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1213 <      MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1214 <      MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1215 < #endif
1216 <      
1217 <      com /= totalMass;
1218 <      comVel /= totalMass;
1219 <   }        
1220 <  
1221 <   /*
1222 <   Return intertia tensor for entire system and angular momentum Vector.
1223 <
1224 <
1225 <       [  Ixx -Ixy  -Ixz ]
1226 <    J =| -Iyx  Iyy  -Iyz |
1227 <       [ -Izx -Iyz   Izz ]
1228 <    */
1229 <
1230 <   void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1231 <      
1232 <
1233 <      RealType xx = 0.0;
1234 <      RealType yy = 0.0;
1235 <      RealType zz = 0.0;
1236 <      RealType xy = 0.0;
1237 <      RealType xz = 0.0;
1238 <      RealType yz = 0.0;
1239 <      Vector3d com(0.0);
1240 <      Vector3d comVel(0.0);
1241 <      
1242 <      getComAll(com, comVel);
1243 <      
1244 <      SimInfo::MoleculeIterator i;
1245 <      Molecule* mol;
1246 <      
1247 <      Vector3d thisq(0.0);
1248 <      Vector3d thisv(0.0);
1249 <
1250 <      RealType thisMass = 0.0;
1251 <    
1252 <      
1253 <      
1254 <  
1255 <      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1256 <        
1257 <         thisq = mol->getCom()-com;
1258 <         thisv = mol->getComVel()-comVel;
1259 <         thisMass = mol->getMass();
1260 <         // Compute moment of intertia coefficients.
1261 <         xx += thisq[0]*thisq[0]*thisMass;
1262 <         yy += thisq[1]*thisq[1]*thisMass;
1263 <         zz += thisq[2]*thisq[2]*thisMass;
1264 <        
1265 <         // compute products of intertia
1266 <         xy += thisq[0]*thisq[1]*thisMass;
1267 <         xz += thisq[0]*thisq[2]*thisMass;
1268 <         yz += thisq[1]*thisq[2]*thisMass;
1269 <            
1270 <         angularMomentum += cross( thisq, thisv ) * thisMass;
1271 <            
1272 <      }  
1273 <      
1274 <      
1275 <      inertiaTensor(0,0) = yy + zz;
1276 <      inertiaTensor(0,1) = -xy;
1277 <      inertiaTensor(0,2) = -xz;
1278 <      inertiaTensor(1,0) = -xy;
1279 <      inertiaTensor(1,1) = xx + zz;
1280 <      inertiaTensor(1,2) = -yz;
1281 <      inertiaTensor(2,0) = -xz;
1282 <      inertiaTensor(2,1) = -yz;
1283 <      inertiaTensor(2,2) = xx + yy;
1284 <      
1285 < #ifdef IS_MPI
1286 <      Mat3x3d tmpI(inertiaTensor);
1287 <      Vector3d tmpAngMom;
1288 <      MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1289 <      MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1290 < #endif
1291 <              
1292 <      return;
1293 <   }
1294 <
1295 <   //Returns the angular momentum of the system
1296 <   Vector3d SimInfo::getAngularMomentum(){
1297 <      
1298 <      Vector3d com(0.0);
1299 <      Vector3d comVel(0.0);
1300 <      Vector3d angularMomentum(0.0);
1301 <      
1302 <      getComAll(com,comVel);
1303 <      
1304 <      SimInfo::MoleculeIterator i;
1305 <      Molecule* mol;
1306 <      
1307 <      Vector3d thisr(0.0);
1308 <      Vector3d thisp(0.0);
1309 <      
1310 <      RealType thisMass;
1311 <      
1312 <      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {        
1313 <        thisMass = mol->getMass();
1314 <        thisr = mol->getCom()-com;
1315 <        thisp = (mol->getComVel()-comVel)*thisMass;
1316 <        
1317 <        angularMomentum += cross( thisr, thisp );
1318 <        
1319 <      }  
1320 <      
1321 < #ifdef IS_MPI
1322 <      Vector3d tmpAngMom;
1323 <      MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1324 < #endif
1325 <      
1326 <      return angularMomentum;
1327 <   }
1328 <  
1002 >  
1003    StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) {
1004      return IOIndexToIntegrableObject.at(index);
1005    }
# Line 1333 | Line 1007 | namespace OpenMD {
1007    void SimInfo::setIOIndexToIntegrableObject(const vector<StuntDouble*>& v) {
1008      IOIndexToIntegrableObject= v;
1009    }
1336
1337  /* Returns the Volume of the simulation based on a ellipsoid with semi-axes
1338     based on the radius of gyration V=4/3*Pi*R_1*R_2*R_3
1339     where R_i are related to the principle inertia moments R_i = sqrt(C*I_i/N), this reduces to
1340     V = 4/3*Pi*(C/N)^3/2*sqrt(det(I)). See S.E. Baltazar et. al. Comp. Mat. Sci. 37 (2006) 526-536.
1341  */
1342  void SimInfo::getGyrationalVolume(RealType &volume){
1343    Mat3x3d intTensor;
1344    RealType det;
1345    Vector3d dummyAngMom;
1346    RealType sysconstants;
1347    RealType geomCnst;
1348
1349    geomCnst = 3.0/2.0;
1350    /* Get the inertial tensor and angular momentum for free*/
1351    getInertiaTensor(intTensor,dummyAngMom);
1352    
1353    det = intTensor.determinant();
1354    sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1355    volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(det);
1356    return;
1357  }
1358
1359  void SimInfo::getGyrationalVolume(RealType &volume, RealType &detI){
1360    Mat3x3d intTensor;
1361    Vector3d dummyAngMom;
1362    RealType sysconstants;
1363    RealType geomCnst;
1364
1365    geomCnst = 3.0/2.0;
1366    /* Get the inertial tensor and angular momentum for free*/
1367    getInertiaTensor(intTensor,dummyAngMom);
1368    
1369    detI = intTensor.determinant();
1370    sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1371    volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(detI);
1372    return;
1373  }
1010   /*
1011     void SimInfo::setStuntDoubleFromGlobalIndex(vector<StuntDouble*> v) {
1012        assert( v.size() == nAtoms_ + nRigidBodies_);

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines