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 1764 by gezelter, Tue Jul 3 18:32:27 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/DarkSide/neighborLists_interface.h"
58 #include "UseTheForce/doForces_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
64   #ifdef IS_MPI
65 < #include "UseTheForce/mpiComponentPlan.h"
66 < #include "UseTheForce/DarkSide/simParallel_interface.h"
69 < #endif
65 > #include <mpi.h>
66 > #endif
67  
68   using namespace std;
69   namespace OpenMD {
# Line 75 | 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 132 | 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
135    std::cerr << "nGA = " << nGlobalAtoms_ << "\n";
136    std::cerr << "nCA = " << nCutoffAtoms << "\n";
137    std::cerr << "nG = " << nGroups << "\n";
132  
133      nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
140
141    std::cerr << "nGCG = " << nGlobalCutoffGroups_ << "\n";
134      
135      //every free atom (atom does not belong to rigid bodies) is an
136      //integrable object therefore the total number of integrable objects
# Line 233 | 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 254 | Line 250 | namespace OpenMD {
250              ndf_local += 3;
251            }
252          }
257            
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 280 | Line 285 | namespace OpenMD {
285      fdf_ = fdf_local;
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() {
# 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 <
720 >    
721   #ifdef IS_MPI
722  
723      // loop over the found atom types on this processor, and add their
724      // numerical idents to a vector:
725 <
725 >    
726      vector<int> foundTypes;
727      set<AtomType*>::iterator i;
728      for (i = atomTypes.begin(); i != atomTypes.end(); ++i)
# Line 706 | Line 731 | namespace OpenMD {
731      // count_local holds the number of found types on this processor
732      int count_local = foundTypes.size();
733  
734 <    // 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);
734 >    int nproc = MPI::COMM_WORLD.Get_size();
735  
736 <    // create a vector to hold the globally found types, and resize it:
737 <    vector<int> ftGlobal;
738 <    ftGlobal.resize(count);
739 <    vector<int> counts;
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 <    int nproc = MPI::COMM_WORLD.Get_size();
742 <    counts.resize(nproc);
743 <    vector<int> disps;
744 <    disps.resize(nproc);
741 >    // fill the counts array
742 >    MPI::COMM_WORLD.Allgather(&count_local, 1, MPI::INT, &counts[0],
743 >                              1, MPI::INT);
744 >  
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 <    // now spray out the foundTypes to all the other processors:
753 >    // we need a (possibly redundant) set of all found types:
754 >    vector<int> ftGlobal(totalCount);
755      
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], MPI::INT);
758 >                               &ftGlobal[0], &counts[0], &disps[0],
759 >                               MPI::INT);
760  
761 +    vector<int>::iterator j;
762 +
763      // foundIdents is a stl set, so inserting an already found ident
764      // will have no effect.
765      set<int> foundIdents;
766 <    vector<int>::iterator j;
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)
773 >    for (it = foundIdents.begin(); it != foundIdents.end(); ++it)
774        atomTypes.insert( forceField_->getAtomType((*it)) );
775  
776   #endif
777 <    
777 >
778      return atomTypes;        
779    }
780  
# Line 752 | Line 786 | namespace OpenMD {
786        if ( simParams_->getAccumulateBoxDipole() ) {
787          calcBoxDipole_ = true;      
788        }
789 <
789 >    
790      set<AtomType*>::iterator i;
791      set<AtomType*> atomTypes;
792      atomTypes = getSimulatedAtomTypes();    
793      int usesElectrostatic = 0;
794      int usesMetallic = 0;
795      int usesDirectional = 0;
796 +    int usesFluctuatingCharges =  0;
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;
807      temp = usesDirectional;
808      MPI_Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
809 <
809 >    
810      temp = usesMetallic;
811      MPI_Allreduce(&temp, &usesMetallicAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);    
812 <
812 >    
813      temp = usesElectrostatic;
814      MPI_Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
815 +
816 +    temp = usesFluctuatingCharges;
817 +    MPI_Allreduce(&temp, &usesFluctuatingCharges_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
818 + #else
819 +
820 +    usesDirectionalAtoms_ = usesDirectional;
821 +    usesMetallicAtoms_ = usesMetallic;
822 +    usesElectrostaticAtoms_ = usesElectrostatic;
823 +    usesFluctuatingCharges_ = usesFluctuatingCharges;
824 +
825   #endif
826 <    fInfo_.SIM_uses_PBC = usesPeriodicBoundaries_;    
827 <    fInfo_.SIM_uses_DirectionalAtoms = usesDirectionalAtoms_;
828 <    fInfo_.SIM_uses_MetallicAtoms = usesMetallicAtoms_;
829 <    fInfo_.SIM_requires_SkipCorrection = usesElectrostaticAtoms_;
784 <    fInfo_.SIM_requires_SelfCorrection = usesElectrostaticAtoms_;
785 <    fInfo_.SIM_uses_AtomicVirial = usesAtomicVirial_;
826 >    
827 >    requiresPrepair_ = usesMetallicAtoms_ ? true : false;
828 >    requiresSkipCorrection_ = usesElectrostaticAtoms_ ? true : false;
829 >    requiresSelfCorrection_ = usesElectrostaticAtoms_ ? true : false;    
830    }
831  
832  
# Line 825 | Line 869 | namespace OpenMD {
869    }
870  
871  
872 <  void SimInfo::setupFortran() {
829 <    int isError;
872 >  void SimInfo::prepareTopology() {
873      int nExclude, nOneTwo, nOneThree, nOneFour;
831    vector<int> fortranGlobalGroupMembership;
832    
833    isError = 0;
874  
835    //globalGroupMembership_ is filled by SimCreator    
836    for (int i = 0; i < nGlobalAtoms_; i++) {
837      fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
838    }
839
875      //calculate mass ratio of cutoff group
841    vector<RealType> mfact;
876      SimInfo::MoleculeIterator mi;
877      Molecule* mol;
878      Molecule::CutoffGroupIterator ci;
# Line 847 | Line 881 | namespace OpenMD {
881      Atom* atom;
882      RealType totalMass;
883  
884 <    //to avoid memory reallocation, reserve enough space for mfact
885 <    mfact.reserve(getNCutoffGroups());
884 >    /**
885 >     * The mass factor is the relative mass of an atom to the total
886 >     * mass of the cutoff group it belongs to.  By default, all atoms
887 >     * are their own cutoff groups, and therefore have mass factors of
888 >     * 1.  We need some special handling for massless atoms, which
889 >     * will be treated as carrying the entire mass of the cutoff
890 >     * group.
891 >     */
892 >    massFactors_.clear();
893 >    massFactors_.resize(getNAtoms(), 1.0);
894      
895      for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {        
896 <      for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
896 >      for (cg = mol->beginCutoffGroup(ci); cg != NULL;
897 >           cg = mol->nextCutoffGroup(ci)) {
898  
899          totalMass = cg->getMass();
900          for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
901            // Check for massless groups - set mfact to 1 if true
902 <          if (totalMass != 0)
903 <            mfact.push_back(atom->getMass()/totalMass);
902 >          if (totalMass != 0)
903 >            massFactors_[atom->getLocalIndex()] = atom->getMass()/totalMass;
904            else
905 <            mfact.push_back( 1.0 );
905 >            massFactors_[atom->getLocalIndex()] = 1.0;
906          }
907        }      
908      }
# Line 873 | Line 916 | namespace OpenMD {
916          identArray_.push_back(atom->getIdent());
917        }
918      }    
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    }
919      
920 <    //setup fortran simulation
920 >    //scan topology
921  
922      nExclude = excludedInteractions_.getSize();
923      nOneTwo = oneTwoInteractions_.getSize();
# Line 893 | Line 929 | namespace OpenMD {
929      int* oneThreeList = oneThreeInteractions_.getPairList();
930      int* oneFourList = oneFourInteractions_.getPairList();
931  
932 <    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;
932 >    topologyDone_ = true;
933    }
934  
935    void SimInfo::addProperty(GenericData* genData) {
# Line 1018 | Line 989 | namespace OpenMD {
989      
990    }
991  
1021  Vector3d SimInfo::getComVel(){
1022    SimInfo::MoleculeIterator i;
1023    Molecule* mol;
992  
1025    Vector3d comVel(0.0);
1026    RealType totalMass = 0.0;
1027    
1028
1029    for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1030      RealType mass = mol->getMass();
1031      totalMass += mass;
1032      comVel += mass * mol->getComVel();
1033    }  
1034
1035 #ifdef IS_MPI
1036    RealType tmpMass = totalMass;
1037    Vector3d tmpComVel(comVel);    
1038    MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1039    MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1040 #endif
1041
1042    comVel /= totalMass;
1043
1044    return comVel;
1045  }
1046
1047  Vector3d SimInfo::getCom(){
1048    SimInfo::MoleculeIterator i;
1049    Molecule* mol;
1050
1051    Vector3d com(0.0);
1052    RealType totalMass = 0.0;
1053    
1054    for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1055      RealType mass = mol->getMass();
1056      totalMass += mass;
1057      com += mass * mol->getCom();
1058    }  
1059
1060 #ifdef IS_MPI
1061    RealType tmpMass = totalMass;
1062    Vector3d tmpCom(com);    
1063    MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1064    MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1065 #endif
1066
1067    com /= totalMass;
1068
1069    return com;
1070
1071  }        
1072
993    ostream& operator <<(ostream& o, SimInfo& info) {
994  
995      return o;
996    }
997    
998 <  
1079 <   /*
1080 <   Returns center of mass and center of mass velocity in one function call.
1081 <   */
1082 <  
1083 <   void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1084 <      SimInfo::MoleculeIterator i;
1085 <      Molecule* mol;
1086 <      
1087 <    
1088 <      RealType totalMass = 0.0;
1089 <    
1090 <
1091 <      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1092 <         RealType mass = mol->getMass();
1093 <         totalMass += mass;
1094 <         com += mass * mol->getCom();
1095 <         comVel += mass * mol->getComVel();          
1096 <      }  
1097 <      
1098 < #ifdef IS_MPI
1099 <      RealType tmpMass = totalMass;
1100 <      Vector3d tmpCom(com);  
1101 <      Vector3d tmpComVel(comVel);
1102 <      MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1103 <      MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1104 <      MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1105 < #endif
1106 <      
1107 <      com /= totalMass;
1108 <      comVel /= totalMass;
1109 <   }        
1110 <  
1111 <   /*
1112 <   Return intertia tensor for entire system and angular momentum Vector.
1113 <
1114 <
1115 <       [  Ixx -Ixy  -Ixz ]
1116 <    J =| -Iyx  Iyy  -Iyz |
1117 <       [ -Izx -Iyz   Izz ]
1118 <    */
1119 <
1120 <   void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1121 <      
1122 <
1123 <      RealType xx = 0.0;
1124 <      RealType yy = 0.0;
1125 <      RealType zz = 0.0;
1126 <      RealType xy = 0.0;
1127 <      RealType xz = 0.0;
1128 <      RealType yz = 0.0;
1129 <      Vector3d com(0.0);
1130 <      Vector3d comVel(0.0);
1131 <      
1132 <      getComAll(com, comVel);
1133 <      
1134 <      SimInfo::MoleculeIterator i;
1135 <      Molecule* mol;
1136 <      
1137 <      Vector3d thisq(0.0);
1138 <      Vector3d thisv(0.0);
1139 <
1140 <      RealType thisMass = 0.0;
1141 <    
1142 <      
1143 <      
1144 <  
1145 <      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1146 <        
1147 <         thisq = mol->getCom()-com;
1148 <         thisv = mol->getComVel()-comVel;
1149 <         thisMass = mol->getMass();
1150 <         // Compute moment of intertia coefficients.
1151 <         xx += thisq[0]*thisq[0]*thisMass;
1152 <         yy += thisq[1]*thisq[1]*thisMass;
1153 <         zz += thisq[2]*thisq[2]*thisMass;
1154 <        
1155 <         // compute products of intertia
1156 <         xy += thisq[0]*thisq[1]*thisMass;
1157 <         xz += thisq[0]*thisq[2]*thisMass;
1158 <         yz += thisq[1]*thisq[2]*thisMass;
1159 <            
1160 <         angularMomentum += cross( thisq, thisv ) * thisMass;
1161 <            
1162 <      }  
1163 <      
1164 <      
1165 <      inertiaTensor(0,0) = yy + zz;
1166 <      inertiaTensor(0,1) = -xy;
1167 <      inertiaTensor(0,2) = -xz;
1168 <      inertiaTensor(1,0) = -xy;
1169 <      inertiaTensor(1,1) = xx + zz;
1170 <      inertiaTensor(1,2) = -yz;
1171 <      inertiaTensor(2,0) = -xz;
1172 <      inertiaTensor(2,1) = -yz;
1173 <      inertiaTensor(2,2) = xx + yy;
1174 <      
1175 < #ifdef IS_MPI
1176 <      Mat3x3d tmpI(inertiaTensor);
1177 <      Vector3d tmpAngMom;
1178 <      MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1179 <      MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1180 < #endif
1181 <              
1182 <      return;
1183 <   }
1184 <
1185 <   //Returns the angular momentum of the system
1186 <   Vector3d SimInfo::getAngularMomentum(){
1187 <      
1188 <      Vector3d com(0.0);
1189 <      Vector3d comVel(0.0);
1190 <      Vector3d angularMomentum(0.0);
1191 <      
1192 <      getComAll(com,comVel);
1193 <      
1194 <      SimInfo::MoleculeIterator i;
1195 <      Molecule* mol;
1196 <      
1197 <      Vector3d thisr(0.0);
1198 <      Vector3d thisp(0.0);
1199 <      
1200 <      RealType thisMass;
1201 <      
1202 <      for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {        
1203 <        thisMass = mol->getMass();
1204 <        thisr = mol->getCom()-com;
1205 <        thisp = (mol->getComVel()-comVel)*thisMass;
1206 <        
1207 <        angularMomentum += cross( thisr, thisp );
1208 <        
1209 <      }  
1210 <      
1211 < #ifdef IS_MPI
1212 <      Vector3d tmpAngMom;
1213 <      MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1214 < #endif
1215 <      
1216 <      return angularMomentum;
1217 <   }
1218 <  
998 >  
999    StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) {
1000      return IOIndexToIntegrableObject.at(index);
1001    }
# Line 1223 | Line 1003 | namespace OpenMD {
1003    void SimInfo::setIOIndexToIntegrableObject(const vector<StuntDouble*>& v) {
1004      IOIndexToIntegrableObject= v;
1005    }
1226
1227  /* Returns the Volume of the simulation based on a ellipsoid with semi-axes
1228     based on the radius of gyration V=4/3*Pi*R_1*R_2*R_3
1229     where R_i are related to the principle inertia moments R_i = sqrt(C*I_i/N), this reduces to
1230     V = 4/3*Pi*(C/N)^3/2*sqrt(det(I)). See S.E. Baltazar et. al. Comp. Mat. Sci. 37 (2006) 526-536.
1231  */
1232  void SimInfo::getGyrationalVolume(RealType &volume){
1233    Mat3x3d intTensor;
1234    RealType det;
1235    Vector3d dummyAngMom;
1236    RealType sysconstants;
1237    RealType geomCnst;
1238
1239    geomCnst = 3.0/2.0;
1240    /* Get the inertial tensor and angular momentum for free*/
1241    getInertiaTensor(intTensor,dummyAngMom);
1242    
1243    det = intTensor.determinant();
1244    sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1245    volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(det);
1246    return;
1247  }
1248
1249  void SimInfo::getGyrationalVolume(RealType &volume, RealType &detI){
1250    Mat3x3d intTensor;
1251    Vector3d dummyAngMom;
1252    RealType sysconstants;
1253    RealType geomCnst;
1254
1255    geomCnst = 3.0/2.0;
1256    /* Get the inertial tensor and angular momentum for free*/
1257    getInertiaTensor(intTensor,dummyAngMom);
1258    
1259    detI = intTensor.determinant();
1260    sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1261    volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(detI);
1262    return;
1263  }
1006   /*
1007     void SimInfo::setStuntDoubleFromGlobalIndex(vector<StuntDouble*> v) {
1008        assert( v.size() == nAtoms_ + nRigidBodies_);

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines