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

Comparing trunk/src/restraints/RestraintForceManager.cpp (file contents):
Revision 1969 by gezelter, Wed Feb 26 14:14:50 2014 UTC vs.
Revision 2024 by gezelter, Thu Oct 16 19:13:51 2014 UTC

# Line 225 | Line 225 | namespace OpenMD {
225          evaluator.loadScriptString(objectSelection);
226          seleMan.setSelectionSet(evaluator.evaluate());        
227          int selectionCount = seleMan.getSelectionCount();
228 <        
228 >
229 > #ifdef IS_MPI    
230 >        MPI_Allreduce(MPI_IN_PLACE, &selectionCount, 1, MPI_INT, MPI_SUM,
231 >                      MPI_COMM_WORLD);
232 > #endif
233 >                
234          sprintf(painCave.errMsg,
235                  "Restraint Info: The specified restraint objectSelection,\n"
236                  "\t\t%s\n"
237                  "\twill result in %d integrable objects being\n"
238                  "\trestrained.\n", objectSelection.c_str(), selectionCount);
239 +        painCave.severity = OPENMD_INFO;
240          painCave.isFatal = 0;
241          simError();                    
242  
# Line 313 | Line 319 | namespace OpenMD {
319    void RestraintForceManager::calcForces(){
320  
321      ForceManager::calcForces();    
322 <    RealType restPot_local, restPot;
322 >    RealType restPot(0.0);
323  
324 <    restPot_local = doRestraints(1.0);
324 >    restPot = doRestraints(1.0);
325  
326   #ifdef IS_MPI    
327 <    MPI_Allreduce(&restPot_local, &restPot, 1,
328 <                  MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
323 < #else
324 <    restPot = restPot_local;
327 >    MPI_Allreduce(MPI_IN_PLACE, &restPot, 1, MPI_REALTYPE, MPI_SUM,
328 >                  MPI_COMM_WORLD);
329   #endif
330 +
331      currSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot();
327    RealType pot = currSnapshot_->getLongRangePotential();
328    pot += restPot;
329    currSnapshot_->setLongRangePotential(pot);
332      currSnapshot_->setRestraintPotential(restPot);
333  
334      //write out forces and current positions of restrained molecules    
# Line 391 | Line 393 | namespace OpenMD {
393  
394        std::vector<Vector3d> struc;
395        std::vector<Vector3d> forces;
396 +
397        
398        for(sd = (*rm)->beginIntegrableObject(ioi); sd != NULL;
399            sd = (*rm)->nextIntegrableObject(ioi)) {
# Line 409 | Line 412 | namespace OpenMD {
412          index++;
413        }
414        
415 <      unscaledPotential_ += mRest->getUnscaledPotential();
415 >      unscaledPotential_ += mRest->getUnscaledPotential();      
416  
414      restInfo = mRest->getRestraintInfo();
415
417        // only collect data on restraints that we're going to print:
418        if (mRest->getPrintRestraint())
419 +        restInfo = mRest->getRestraintInfo();
420          restInfo_.push_back(restInfo);
421      }
422  
# Line 451 | Line 453 | namespace OpenMD {
453        
454        // phew.  At this point, we should have the pointer to the
455        // correct Object restraint in the variable oRest.
454
456        oRest->setScaleFactor(scalingFactor);
457        
458        Vector3d pos = (*ro)->getPos();
# Line 475 | Line 476 | namespace OpenMD {
476  
477        unscaledPotential_ += oRest->getUnscaledPotential();
478  
478      restInfo = oRest->getRestraintInfo();
479
479        // only collect data on restraints that we're going to print:
480        if (oRest->getPrintRestraint())
481 +        restInfo = oRest->getRestraintInfo();      
482          restInfo_.push_back(restInfo);
483      }
484  

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines