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

Comparing trunk/src/brains/ForceManager.cpp (file contents):
Revision 1788 by gezelter, Wed Aug 29 20:17:07 2012 UTC vs.
Revision 1993 by gezelter, Tue Apr 29 17:32:31 2014 UTC

# Line 35 | Line 35
35   *                                                                      
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).          
38 > * [3]  Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).          
39   * [4]  Kuang & Gezelter,  J. Chem. Phys. 133, 164101 (2010).
40   * [5]  Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41   */
# Line 44 | Line 44
44   * @file ForceManager.cpp
45   * @author tlin
46   * @date 11/09/2004
47 * @time 10:39am
47   * @version 1.0
48   */
49  
# Line 68 | Line 67 | namespace OpenMD {
67   using namespace std;
68   namespace OpenMD {
69    
70 <  ForceManager::ForceManager(SimInfo * info) : info_(info) {
70 >  ForceManager::ForceManager(SimInfo * info) : info_(info), switcher_(NULL),
71 >                                               initialized_(false) {
72      forceField_ = info_->getForceField();
73      interactionMan_ = new InteractionManager();
74      fDecomp_ = new ForceMatrixDecomposition(info_, interactionMan_);
75 +    thermo = new Thermo(info_);
76    }
77  
78 +  ForceManager::~ForceManager() {
79 +    perturbations_.clear();
80 +    
81 +    delete switcher_;
82 +    delete interactionMan_;
83 +    delete fDecomp_;
84 +    delete thermo;
85 +  }
86 +  
87    /**
88     * setupCutoffs
89     *
# Line 88 | Line 98 | namespace OpenMD {
98     *      simulation for suggested cutoff values (e.g. 2.5 * sigma).
99     *      Use the maximum suggested value that was found.
100     *
101 <   * cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE,
102 <   *                        or SHIFTED_POTENTIAL)
101 >   * cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE, TAYLOR_SHIFTED,
102 >   *                        SHIFTED_POTENTIAL, or EWALD_FULL)
103     *      If cutoffMethod was explicitly set, use that choice.
104     *      If cutoffMethod was not explicitly set, use SHIFTED_FORCE
105     *
# Line 118 | Line 128 | namespace OpenMD {
128      else
129        mdFileVersion = 0;
130    
131 +    // We need the list of simulated atom types to figure out cutoffs
132 +    // as well as long range corrections.
133 +
134 +    set<AtomType*>::iterator i;
135 +    set<AtomType*> atomTypes_;
136 +    atomTypes_ = info_->getSimulatedAtomTypes();
137 +
138      if (simParams_->haveCutoffRadius()) {
139        rCut_ = simParams_->getCutoffRadius();
140      } else {      
# Line 132 | Line 149 | namespace OpenMD {
149          rCut_ = 12.0;
150        } else {
151          RealType thisCut;
152 <        set<AtomType*>::iterator i;
136 <        set<AtomType*> atomTypes;
137 <        atomTypes = info_->getSimulatedAtomTypes();        
138 <        for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
152 >        for (i = atomTypes_.begin(); i != atomTypes_.end(); ++i) {
153            thisCut = interactionMan_->getSuggestedCutoffRadius((*i));
154            rCut_ = max(thisCut, rCut_);
155          }
# Line 157 | Line 171 | namespace OpenMD {
171      stringToCutoffMethod["SWITCHED"] = SWITCHED;
172      stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL;    
173      stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE;
174 +    stringToCutoffMethod["TAYLOR_SHIFTED"] = TAYLOR_SHIFTED;
175 +    stringToCutoffMethod["EWALD_FULL"] = EWALD_FULL;
176    
177      if (simParams_->haveCutoffMethod()) {
178        string cutMeth = toUpperCopy(simParams_->getCutoffMethod());
# Line 166 | Line 182 | namespace OpenMD {
182          sprintf(painCave.errMsg,
183                  "ForceManager::setupCutoffs: Could not find chosen cutoffMethod %s\n"
184                  "\tShould be one of: "
185 <                "HARD, SWITCHED, SHIFTED_POTENTIAL, or SHIFTED_FORCE\n",
185 >                "HARD, SWITCHED, SHIFTED_POTENTIAL, TAYLOR_SHIFTED,\n"
186 >                "\tSHIFTED_FORCE, or EWALD_FULL\n",
187                  cutMeth.c_str());
188          painCave.isFatal = 1;
189          painCave.severity = OPENMD_ERROR;
# Line 210 | Line 227 | namespace OpenMD {
227              cutoffMethod_ = SHIFTED_POTENTIAL;
228            } else if (myMethod == "SHIFTED_FORCE") {
229              cutoffMethod_ = SHIFTED_FORCE;
230 +          } else if (myMethod == "TAYLOR_SHIFTED") {
231 +            cutoffMethod_ = TAYLOR_SHIFTED;
232 +          } else if (myMethod == "EWALD_FULL") {
233 +            cutoffMethod_ = EWALD_FULL;
234            }
235          
236            if (simParams_->haveSwitchingRadius())
237              rSwitch_ = simParams_->getSwitchingRadius();
238  
239 <          if (myMethod == "SHIFTED_POTENTIAL" || myMethod == "SHIFTED_FORCE") {
239 >          if (myMethod == "SHIFTED_POTENTIAL" || myMethod == "SHIFTED_FORCE" ||
240 >              myMethod == "TAYLOR_SHIFTED" || myMethod == "EWALD_FULL") {
241              if (simParams_->haveSwitchingRadius()){
242                sprintf(painCave.errMsg,
243                        "ForceManager::setupCutoffs : DEPRECATED ERROR MESSAGE\n"
# Line 370 | Line 392 | namespace OpenMD {
392      }
393      switcher_->setSwitchType(sft_);
394      switcher_->setSwitch(rSwitch_, rCut_);
373    interactionMan_->setSwitchingRadius(rSwitch_);
395    }
396  
397  
# Line 394 | Line 415 | namespace OpenMD {
415        doParticlePot_ = info_->getSimParams()->getOutputParticlePotential();
416        doHeatFlux_ = info_->getSimParams()->getPrintHeatFlux();
417        if (doHeatFlux_) doParticlePot_ = true;
418 +
419 +      doElectricField_ = info_->getSimParams()->getOutputElectricField();
420 +      doSitePotential_ = info_->getSimParams()->getOutputSitePotential();
421    
422      }
423  
# Line 428 | Line 452 | namespace OpenMD {
452        ElectricField* eField = new ElectricField(info_);
453        perturbations_.push_back(eField);
454      }
431
432    fDecomp_->distributeInitialData();
433
434    initialized_ = true;
455  
456 +    usePeriodicBoundaryConditions_ = info_->getSimParams()->getUsePeriodicBoundaryConditions();
457 +    
458 +    fDecomp_->distributeInitialData();
459 +    
460 +    initialized_ = true;
461 +    
462    }
463 <
463 >  
464    void ForceManager::calcForces() {
465      
466      if (!initialized_) initialize();
467 <
467 >    
468      preCalculation();  
469      shortRangeInteractions();
470      longRangeInteractions();
# Line 612 | Line 638 | namespace OpenMD {
638      // Collect from all nodes.  This should eventually be moved into a
639      // SystemDecomposition, but this is a better place than in
640      // Thermo to do the collection.
641 <    MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &bondPotential, 1, MPI::REALTYPE,
642 <                              MPI::SUM);
643 <    MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &bendPotential, 1, MPI::REALTYPE,
644 <                              MPI::SUM);
645 <    MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &torsionPotential, 1,
646 <                              MPI::REALTYPE, MPI::SUM);
647 <    MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &inversionPotential, 1,
648 <                              MPI::REALTYPE, MPI::SUM);
641 >
642 >    MPI_Allreduce(MPI_IN_PLACE, &bondPotential, 1, MPI_REALTYPE,
643 >                  MPI_SUM, MPI_COMM_WORLD);
644 >    MPI_Allreduce(MPI_IN_PLACE, &bendPotential, 1, MPI_REALTYPE,
645 >                  MPI_SUM, MPI_COMM_WORLD);
646 >    MPI_Allreduce(MPI_IN_PLACE, &torsionPotential, 1,
647 >                  MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
648 >    MPI_Allreduce(MPI_IN_PLACE, &inversionPotential, 1,
649 >                  MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
650   #endif
651  
652      Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
# Line 637 | Line 664 | namespace OpenMD {
664    
665    void ForceManager::longRangeInteractions() {
666  
640
667      Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
668      DataStorage* config = &(curSnapshot->atomData);
669      DataStorage* cgConfig = &(curSnapshot->cgData);
670  
671 +
672      //calculate the center of mass of cutoff group
673  
674      SimInfo::MoleculeIterator mi;
# Line 674 | Line 701 | namespace OpenMD {
701      RealType vij;
702      Vector3d fij, fg, f1;
703      tuple3<RealType, RealType, RealType> cuts;
704 <    RealType rCutSq;
704 >    RealType rCut, rCutSq, rListSq;
705      bool in_switching_region;
706      RealType sw, dswdr, swderiv;
707 <    vector<int> atomListColumn, atomListRow, atomListLocal;
707 >    vector<int> atomListColumn, atomListRow;
708      InteractionData idat;
709      SelfData sdat;
710      RealType mf;
711      RealType vpair;
712      RealType dVdFQ1(0.0);
713      RealType dVdFQ2(0.0);
687    Vector3d eField1(0.0);
688    Vector3d eField2(0.0);
714      potVec longRangePotential(0.0);
715 +    RealType reciprocalPotential(0.0);
716      potVec workPot(0.0);
717      potVec exPot(0.0);
718 +    Vector3d eField1(0.0);
719 +    Vector3d eField2(0.0);
720 +    RealType sPot1(0.0);
721 +    RealType sPot2(0.0);
722 +                  
723      vector<int>::iterator ia, jb;
724  
725      int loopStart, loopEnd;
726 <
726 >    
727 >    idat.rcut = &rCut;
728      idat.vdwMult = &vdwMult;
729      idat.electroMult = &electroMult;
730      idat.pot = &workPot;
# Line 703 | Line 735 | namespace OpenMD {
735      idat.dVdFQ1 = &dVdFQ1;
736      idat.dVdFQ2 = &dVdFQ2;
737      idat.eField1 = &eField1;
738 <    idat.eField2 = &eField2;
738 >    idat.eField2 = &eField2;
739 >    idat.sPot1 = &sPot1;
740 >    idat.sPot2 = &sPot2;
741      idat.f1 = &f1;
742      idat.sw = &sw;
743      idat.shiftedPot = (cutoffMethod_ == SHIFTED_POTENTIAL) ? true : false;
744 <    idat.shiftedForce = (cutoffMethod_ == SHIFTED_FORCE) ? true : false;
744 >    idat.shiftedForce = (cutoffMethod_ == SHIFTED_FORCE || cutoffMethod_ == TAYLOR_SHIFTED) ? true : false;
745      idat.doParticlePot = doParticlePot_;
746 +    idat.doElectricField = doElectricField_;
747 +    idat.doSitePotential = doSitePotential_;
748      sdat.doParticlePot = doParticlePot_;
749      
750      loopEnd = PAIR_LOOP;
# Line 721 | Line 757 | namespace OpenMD {
757      
758        if (iLoop == loopStart) {
759          bool update_nlist = fDecomp_->checkNeighborList();
760 <        if (update_nlist)
761 <          neighborList = fDecomp_->buildNeighborList();
762 <      }            
760 >        if (update_nlist) {
761 >          if (!usePeriodicBoundaryConditions_)
762 >            Mat3x3d bbox = thermo->getBoundingBox();
763 >          fDecomp_->buildNeighborList(neighborList_);
764 >        }
765 >      }
766  
767 <      for (vector<pair<int, int> >::iterator it = neighborList.begin();
768 <             it != neighborList.end(); ++it) {
767 >      for (vector<pair<int, int> >::iterator it = neighborList_.begin();
768 >             it != neighborList_.end(); ++it) {
769                  
770          cg1 = (*it).first;
771          cg2 = (*it).second;
772          
773 <        cuts = fDecomp_->getGroupCutoffs(cg1, cg2);
773 >        fDecomp_->getGroupCutoffs(cg1, cg2, rCut, rCutSq, rListSq);
774  
775          d_grp  = fDecomp_->getIntergroupVector(cg1, cg2);
776  
777 <        curSnapshot->wrapVector(d_grp);        
777 >        // already wrapped in the getIntergroupVector call:
778 >        // curSnapshot->wrapVector(d_grp);        
779          rgrpsq = d_grp.lengthSquare();
740        rCutSq = cuts.second;
780  
781          if (rgrpsq < rCutSq) {
743          idat.rcut = &cuts.first;
782            if (iLoop == PAIR_LOOP) {
783              vij = 0.0;
784 <            fij = V3Zero;
784 >            fij.zero();
785 >            eField1.zero();
786 >            eField2.zero();
787 >            sPot1 = 0.0;
788 >            sPot2 = 0.0;
789            }
790            
791            in_switching_region = switcher_->getSwitch(rgrpsq, sw, dswdr,
# Line 768 | Line 810 | namespace OpenMD {
810                  vpair = 0.0;
811                  workPot = 0.0;
812                  exPot = 0.0;
813 <                f1 = V3Zero;
813 >                f1.zero();
814                  dVdFQ1 = 0.0;
815                  dVdFQ2 = 0.0;
816  
# Line 795 | Line 837 | namespace OpenMD {
837                
838                  r = sqrt( *(idat.r2) );
839                  idat.rij = &r;
840 <              
840 >
841                  if (iLoop == PREPAIR_LOOP) {
842                    interactionMan_->doPrePair(idat);
843                  } else {
# Line 818 | Line 860 | namespace OpenMD {
860                fij += fg;
861  
862                if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
863 <                stressTensor -= outProduct( *(idat.d), fg);
864 <                if (doHeatFlux_)
865 <                  fDecomp_->addToHeatFlux(*(idat.d) * dot(fg, vel2));
866 <                
863 >                if (!fDecomp_->skipAtomPair(atomListRow[0],
864 >                                            atomListColumn[0],
865 >                                            cg1, cg2)) {
866 >                  stressTensor -= outProduct( *(idat.d), fg);
867 >                  if (doHeatFlux_)
868 >                    fDecomp_->addToHeatFlux(*(idat.d) * dot(fg, vel2));
869 >                }                
870                }
871            
872                for (ia = atomListRow.begin();
# Line 888 | Line 933 | namespace OpenMD {
933          }
934        }
935      }
936 <  
936 >    
937      // collects pairwise information
938      fDecomp_->collectData();
939 +    if (cutoffMethod_ == EWALD_FULL) {
940 +      interactionMan_->doReciprocalSpaceSum(reciprocalPotential);
941 +
942 +      curSnapshot->setReciprocalPotential(reciprocalPotential);
943 +    }
944          
945      if (info_->requiresSelfCorrection()) {
946        for (unsigned int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
# Line 912 | Line 962 | namespace OpenMD {
962  
963    }
964  
915  
965    void ForceManager::postCalculation() {
966  
967      vector<Perturbation*>::iterator pi;
# Line 938 | Line 987 | namespace OpenMD {
987      }
988      
989   #ifdef IS_MPI
990 <    MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, stressTensor.getArrayPointer(), 9,
991 <                              MPI::REALTYPE, MPI::SUM);
990 >    MPI_Allreduce(MPI_IN_PLACE, stressTensor.getArrayPointer(), 9,
991 >                  MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
992   #endif
993      curSnapshot->setStressTensor(stressTensor);
994      
995 +    if (info_->getSimParams()->getUseLongRangeCorrections()) {
996 +      /*
997 +      RealType vol = curSnapshot->getVolume();
998 +      RealType Elrc(0.0);
999 +      RealType Wlrc(0.0);
1000 +
1001 +      set<AtomType*>::iterator i;
1002 +      set<AtomType*>::iterator j;
1003 +    
1004 +      RealType n_i, n_j;
1005 +      RealType rho_i, rho_j;
1006 +      pair<RealType, RealType> LRI;
1007 +      
1008 +      for (i = atomTypes_.begin(); i != atomTypes_.end(); ++i) {
1009 +        n_i = RealType(info_->getGlobalCountOfType(*i));
1010 +        rho_i = n_i /  vol;
1011 +        for (j = atomTypes_.begin(); j != atomTypes_.end(); ++j) {
1012 +          n_j = RealType(info_->getGlobalCountOfType(*j));
1013 +          rho_j = n_j / vol;
1014 +          
1015 +          LRI = interactionMan_->getLongRangeIntegrals( (*i), (*j) );
1016 +
1017 +          Elrc += n_i   * rho_j * LRI.first;
1018 +          Wlrc -= rho_i * rho_j * LRI.second;
1019 +        }
1020 +      }
1021 +      Elrc *= 2.0 * NumericConstant::PI;
1022 +      Wlrc *= 2.0 * NumericConstant::PI;
1023 +
1024 +      RealType lrp = curSnapshot->getLongRangePotential();
1025 +      curSnapshot->setLongRangePotential(lrp + Elrc);
1026 +      stressTensor += Wlrc * SquareMatrix3<RealType>::identity();
1027 +      curSnapshot->setStressTensor(stressTensor);
1028 +      */
1029 +    
1030 +    }
1031    }
1032 < } //end namespace OpenMD
1032 > }

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines