ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceMatrixDecomposition.cpp
(Generate patch)

Comparing branches/development/src/parallel/ForceMatrixDecomposition.cpp (file contents):
Revision 1721 by gezelter, Thu May 24 14:17:42 2012 UTC vs.
Revision 1761 by gezelter, Fri Jun 22 20:01:37 2012 UTC

# Line 95 | Line 95 | namespace OpenMD {
95      storageLayout_ = sman_->getStorageLayout();
96      ff_ = info_->getForceField();
97      nLocal_ = snap_->getNumberOfAtoms();
98 <    
98 >  
99      nGroups_ = info_->getNLocalCutoffGroups();
100      // gather the information for atomtype IDs (atids):
101      idents = info_->getIdentArray();
# Line 109 | Line 109 | namespace OpenMD {
109      PairList* oneTwo = info_->getOneTwoInteractions();
110      PairList* oneThree = info_->getOneThreeInteractions();
111      PairList* oneFour = info_->getOneFourInteractions();
112 <
112 >    
113 >    if (needVelocities_)
114 >      snap_->cgData.setStorageLayout(DataStorage::dslPosition |
115 >                                     DataStorage::dslVelocity);
116 >    else
117 >      snap_->cgData.setStorageLayout(DataStorage::dslPosition);
118 >    
119   #ifdef IS_MPI
120  
121      MPI::Intracomm row = rowComm.getComm();
# Line 145 | Line 151 | namespace OpenMD {
151      cgRowData.resize(nGroupsInRow_);
152      cgRowData.setStorageLayout(DataStorage::dslPosition);
153      cgColData.resize(nGroupsInCol_);
154 <    cgColData.setStorageLayout(DataStorage::dslPosition);
155 <        
154 >    if (needVelocities_)
155 >      // we only need column velocities if we need them.
156 >      cgColData.setStorageLayout(DataStorage::dslPosition |
157 >                                 DataStorage::dslVelocity);
158 >    else    
159 >      cgColData.setStorageLayout(DataStorage::dslPosition);
160 >      
161      identsRow.resize(nAtomsInRow_);
162      identsCol.resize(nAtomsInCol_);
163      
# Line 164 | Line 175 | namespace OpenMD {
175  
176      pot_row.resize(nAtomsInRow_);
177      pot_col.resize(nAtomsInCol_);
178 +
179 +    expot_row.resize(nAtomsInRow_);
180 +    expot_col.resize(nAtomsInCol_);
181  
182      AtomRowToGlobal.resize(nAtomsInRow_);
183      AtomColToGlobal.resize(nAtomsInCol_);
# Line 450 | Line 464 | namespace OpenMD {
464      }
465    }
466  
453
467    groupCutoffs ForceMatrixDecomposition::getGroupCutoffs(int cg1, int cg2) {
468      int i, j;  
469   #ifdef IS_MPI
# Line 474 | Line 487 | namespace OpenMD {
487    void ForceMatrixDecomposition::zeroWorkArrays() {
488      pairwisePot = 0.0;
489      embeddingPot = 0.0;
490 +    excludedPot = 0.0;
491 +    excludedSelfPot = 0.0;
492  
493   #ifdef IS_MPI
494      if (storageLayout_ & DataStorage::dslForce) {
# Line 492 | Line 507 | namespace OpenMD {
507      fill(pot_col.begin(), pot_col.end(),
508           Vector<RealType, N_INTERACTION_FAMILIES> (0.0));  
509  
510 +    fill(expot_row.begin(), expot_row.end(),
511 +         Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
512 +
513 +    fill(expot_col.begin(), expot_col.end(),
514 +         Vector<RealType, N_INTERACTION_FAMILIES> (0.0));  
515 +
516      if (storageLayout_ & DataStorage::dslParticlePot) {    
517        fill(atomRowData.particlePot.begin(), atomRowData.particlePot.end(),
518             0.0);
# Line 600 | Line 621 | namespace OpenMD {
621      cgPlanVectorColumn->gather(snap_->cgData.position,
622                                 cgColData.position);
623  
624 +
625 +
626 +    if (needVelocities_) {
627 +      // gather up the atomic velocities
628 +      AtomPlanVectorColumn->gather(snap_->atomData.velocity,
629 +                                   atomColData.velocity);
630 +      
631 +      cgPlanVectorColumn->gather(snap_->cgData.velocity,
632 +                                 cgColData.velocity);
633 +    }
634 +
635      
636      // if needed, gather the atomic rotation matrices
637      if (storageLayout_ & DataStorage::dslAmat) {
# Line 759 | Line 791 | namespace OpenMD {
791  
792      vector<potVec> pot_temp(nLocal_,
793                              Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
794 +    vector<potVec> expot_temp(nLocal_,
795 +                              Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
796  
797      // scatter/gather pot_row into the members of my column
798            
799      AtomPlanPotRow->scatter(pot_row, pot_temp);
800 +    AtomPlanPotRow->scatter(expot_row, expot_temp);
801  
802 <    for (int ii = 0;  ii < pot_temp.size(); ii++ )
802 >    for (int ii = 0;  ii < pot_temp.size(); ii++ )
803        pairwisePot += pot_temp[ii];
804 <    
804 >
805 >    for (int ii = 0;  ii < expot_temp.size(); ii++ )
806 >      excludedPot += expot_temp[ii];
807 >        
808 >    if (storageLayout_ & DataStorage::dslParticlePot) {
809 >      // This is the pairwise contribution to the particle pot.  The
810 >      // embedding contribution is added in each of the low level
811 >      // non-bonded routines.  In single processor, this is done in
812 >      // unpackInteractionData, not in collectData.
813 >      for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
814 >        for (int i = 0; i < nLocal_; i++) {
815 >          // factor of two is because the total potential terms are divided
816 >          // by 2 in parallel due to row/ column scatter      
817 >          snap_->atomData.particlePot[i] += 2.0 * pot_temp[i](ii);
818 >        }
819 >      }
820 >    }
821 >
822      fill(pot_temp.begin(), pot_temp.end(),
823           Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
824 +    fill(expot_temp.begin(), expot_temp.end(),
825 +         Vector<RealType, N_INTERACTION_FAMILIES> (0.0));
826        
827      AtomPlanPotColumn->scatter(pot_col, pot_temp);    
828 +    AtomPlanPotColumn->scatter(expot_col, expot_temp);    
829      
830      for (int ii = 0;  ii < pot_temp.size(); ii++ )
831        pairwisePot += pot_temp[ii];    
832 +
833 +    for (int ii = 0;  ii < expot_temp.size(); ii++ )
834 +      excludedPot += expot_temp[ii];    
835 +
836 +    if (storageLayout_ & DataStorage::dslParticlePot) {
837 +      // This is the pairwise contribution to the particle pot.  The
838 +      // embedding contribution is added in each of the low level
839 +      // non-bonded routines.  In single processor, this is done in
840 +      // unpackInteractionData, not in collectData.
841 +      for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
842 +        for (int i = 0; i < nLocal_; i++) {
843 +          // factor of two is because the total potential terms are divided
844 +          // by 2 in parallel due to row/ column scatter      
845 +          snap_->atomData.particlePot[i] += 2.0 * pot_temp[i](ii);
846 +        }
847 +      }
848 +    }
849      
850 +    if (storageLayout_ & DataStorage::dslParticlePot) {
851 +      int npp = snap_->atomData.particlePot.size();
852 +      vector<RealType> ppot_temp(npp, 0.0);
853 +
854 +      // This is the direct or embedding contribution to the particle
855 +      // pot.
856 +      
857 +      AtomPlanRealRow->scatter(atomRowData.particlePot, ppot_temp);
858 +      for (int i = 0; i < npp; i++) {
859 +        snap_->atomData.particlePot[i] += ppot_temp[i];
860 +      }
861 +
862 +      fill(ppot_temp.begin(), ppot_temp.end(), 0.0);
863 +      
864 +      AtomPlanRealColumn->scatter(atomColData.particlePot, ppot_temp);
865 +      for (int i = 0; i < npp; i++) {
866 +        snap_->atomData.particlePot[i] += ppot_temp[i];
867 +      }
868 +    }
869 +
870      for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
871        RealType ploc1 = pairwisePot[ii];
872        RealType ploc2 = 0.0;
# Line 783 | Line 875 | namespace OpenMD {
875      }
876  
877      for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
878 <      RealType ploc1 = embeddingPot[ii];
878 >      RealType ploc1 = excludedPot[ii];
879        RealType ploc2 = 0.0;
880        MPI::COMM_WORLD.Allreduce(&ploc1, &ploc2, 1, MPI::REALTYPE, MPI::SUM);
881 <      embeddingPot[ii] = ploc2;
881 >      excludedPot[ii] = ploc2;
882      }
883  
884 +    // Here be dragons.
885 +    MPI::Intracomm col = colComm.getComm();
886 +
887 +    col.Allreduce(MPI::IN_PLACE,
888 +                  &snap_->frameData.conductiveHeatFlux[0], 3,
889 +                  MPI::REALTYPE, MPI::SUM);
890 +
891 +
892   #endif
893  
894    }
895  
896 +  /**
897 +   * Collects information obtained during the post-pair (and embedding
898 +   * functional) loops onto local data structures.
899 +   */
900 +  void ForceMatrixDecomposition::collectSelfData() {
901 +    snap_ = sman_->getCurrentSnapshot();
902 +    storageLayout_ = sman_->getStorageLayout();
903 +
904 + #ifdef IS_MPI
905 +    for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
906 +      RealType ploc1 = embeddingPot[ii];
907 +      RealType ploc2 = 0.0;
908 +      MPI::COMM_WORLD.Allreduce(&ploc1, &ploc2, 1, MPI::REALTYPE, MPI::SUM);
909 +      embeddingPot[ii] = ploc2;
910 +    }    
911 +    for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
912 +      RealType ploc1 = excludedSelfPot[ii];
913 +      RealType ploc2 = 0.0;
914 +      MPI::COMM_WORLD.Allreduce(&ploc1, &ploc2, 1, MPI::REALTYPE, MPI::SUM);
915 +      excludedSelfPot[ii] = ploc2;
916 +    }    
917 + #endif
918 +    
919 +  }
920 +
921 +
922 +
923    int ForceMatrixDecomposition::getNAtomsInRow() {  
924   #ifdef IS_MPI
925      return nAtomsInRow_;
# Line 831 | Line 958 | namespace OpenMD {
958      
959      snap_->wrapVector(d);
960      return d;    
961 +  }
962 +
963 +  Vector3d ForceMatrixDecomposition::getGroupVelocityColumn(int cg2){
964 + #ifdef IS_MPI
965 +    return cgColData.velocity[cg2];
966 + #else
967 +    return snap_->cgData.velocity[cg2];
968 + #endif
969    }
970  
971 +  Vector3d ForceMatrixDecomposition::getAtomVelocityColumn(int atom2){
972 + #ifdef IS_MPI
973 +    return atomColData.velocity[atom2];
974 + #else
975 +    return snap_->atomData.velocity[atom2];
976 + #endif
977 +  }
978  
979 +
980    Vector3d ForceMatrixDecomposition::getAtomToGroupVectorRow(int atom1, int cg1){
981  
982      Vector3d d;
# Line 899 | Line 1042 | namespace OpenMD {
1042     * We need to exclude some overcounted interactions that result from
1043     * the parallel decomposition.
1044     */
1045 <  bool ForceMatrixDecomposition::skipAtomPair(int atom1, int atom2) {
1046 <    int unique_id_1, unique_id_2;
1045 >  bool ForceMatrixDecomposition::skipAtomPair(int atom1, int atom2, int cg1, int cg2) {
1046 >    int unique_id_1, unique_id_2, group1, group2;
1047          
1048   #ifdef IS_MPI
1049      // in MPI, we have to look up the unique IDs for each atom
1050      unique_id_1 = AtomRowToGlobal[atom1];
1051      unique_id_2 = AtomColToGlobal[atom2];
1052 +    group1 = cgRowToGlobal[cg1];
1053 +    group2 = cgColToGlobal[cg2];
1054   #else
1055      unique_id_1 = AtomLocalToGlobal[atom1];
1056      unique_id_2 = AtomLocalToGlobal[atom2];
1057 +    group1 = cgLocalToGlobal[cg1];
1058 +    group2 = cgLocalToGlobal[cg2];
1059   #endif  
1060  
1061      if (unique_id_1 == unique_id_2) return true;
# Line 919 | Line 1066 | namespace OpenMD {
1066        if ((unique_id_1 + unique_id_2) % 2 == 0) return true;
1067      } else {
1068        if ((unique_id_1 + unique_id_2) % 2 == 1) return true;
1069 +    }
1070 + #endif    
1071 +
1072 + #ifndef IS_MPI
1073 +    if (group1 == group2) {
1074 +      if (unique_id_1 < unique_id_2) return true;
1075      }
1076   #endif
1077      
# Line 1022 | Line 1175 | namespace OpenMD {
1175  
1176   #else
1177      
1025
1026    // cerr << "atoms = " << atom1 << " " << atom2 << "\n";
1027    // cerr << "pos1 = " << snap_->atomData.position[atom1] << "\n";
1028    // cerr << "pos2 = " << snap_->atomData.position[atom2] << "\n";
1029
1178      idat.atypes = make_pair( atypesLocal[atom1], atypesLocal[atom2]);
1031    //idat.atypes = make_pair( ff_->getAtomType(idents[atom1]),
1032    //                         ff_->getAtomType(idents[atom2]) );
1179  
1180      if (storageLayout_ & DataStorage::dslAmat) {
1181        idat.A1 = &(snap_->atomData.aMat[atom1]);
# Line 1084 | Line 1230 | namespace OpenMD {
1230   #ifdef IS_MPI
1231      pot_row[atom1] += RealType(0.5) *  *(idat.pot);
1232      pot_col[atom2] += RealType(0.5) *  *(idat.pot);
1233 +    expot_row[atom1] += RealType(0.5) *  *(idat.excludedPot);
1234 +    expot_col[atom2] += RealType(0.5) *  *(idat.excludedPot);
1235  
1236      atomRowData.force[atom1] += *(idat.f1);
1237      atomColData.force[atom2] -= *(idat.f1);
1238  
1239      if (storageLayout_ & DataStorage::dslFlucQForce) {              
1240 <      atomRowData.flucQFrc[atom1] += *(idat.dVdFQ1);
1241 <      atomColData.flucQFrc[atom2] += *(idat.dVdFQ2);
1240 >      atomRowData.flucQFrc[atom1] -= *(idat.dVdFQ1);
1241 >      atomColData.flucQFrc[atom2] -= *(idat.dVdFQ2);
1242      }
1243  
1244      if (storageLayout_ & DataStorage::dslElectricField) {              
# Line 1098 | Line 1246 | namespace OpenMD {
1246        atomColData.electricField[atom2] += *(idat.eField2);
1247      }
1248  
1101    // should particle pot be done here also?
1249   #else
1250      pairwisePot += *(idat.pot);
1251 +    excludedPot += *(idat.excludedPot);
1252  
1253      snap_->atomData.force[atom1] += *(idat.f1);
1254      snap_->atomData.force[atom2] -= *(idat.f1);
1255  
1256      if (idat.doParticlePot) {
1257 +      // This is the pairwise contribution to the particle pot.  The
1258 +      // embedding contribution is added in each of the low level
1259 +      // non-bonded routines.  In parallel, this calculation is done
1260 +      // in collectData, not in unpackInteractionData.
1261        snap_->atomData.particlePot[atom1] += *(idat.vpair) * *(idat.sw);
1262 <      snap_->atomData.particlePot[atom2] -= *(idat.vpair) * *(idat.sw);
1262 >      snap_->atomData.particlePot[atom2] += *(idat.vpair) * *(idat.sw);
1263      }
1264      
1265      if (storageLayout_ & DataStorage::dslFlucQForce) {              
1266 <      snap_->atomData.flucQFrc[atom1] += *(idat.dVdFQ1);
1266 >      snap_->atomData.flucQFrc[atom1] -= *(idat.dVdFQ1);
1267        snap_->atomData.flucQFrc[atom2] -= *(idat.dVdFQ2);
1268      }
1269  

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines