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

Comparing trunk/src/restraints/ThermoIntegrationForceManager.cpp (file contents):
Revision 507 by gezelter, Fri Apr 15 22:04:00 2005 UTC vs.
Revision 1029 by chrisfen, Thu Aug 31 22:34:49 2006 UTC

# Line 48 | Line 48
48   #include "utils/OOPSEConstant.hpp"
49   #include "utils/StringUtils.hpp"
50  
51 + #ifdef IS_MPI
52 + #include <mpi.h>
53 + #define TAKE_THIS_TAG_REAL 2
54 + #endif //is_mpi
55 +
56   namespace oopse {
57    
58    ThermoIntegrationForceManager::ThermoIntegrationForceManager(SimInfo* info):
# Line 55 | Line 60 | namespace oopse {
60        currSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot();
61        simParam = info_->getSimParams();
62      
63 <      if (simParam->haveThermIntLambda()){
64 <        tIntLambda_ = simParam->getThermIntLambda();
63 >      if (simParam->haveThermodynamicIntegrationLambda()){
64 >        tIntLambda_ = simParam->getThermodynamicIntegrationLambda();
65        }
66        else{
67          tIntLambda_ = 1.0;
68          sprintf(painCave.errMsg,
69 <                "ThermoIntegration error: the transformation parameter (lambda) was\n"
70 <                "\tnot specified. OOPSE will use a default value of %f. To set\n"
71 <                "\tlambda, use the thermodynamicIntegrationLambda variable.\n",
69 >                "ThermoIntegration error: the transformation parameter\n"
70 >                "\t(lambda) was not specified. OOPSE will use a default\n"
71 >                "\tvalue of %f. To set lambda, use the \n"
72 >                "\tthermodynamicIntegrationLambda variable.\n",
73                  tIntLambda_);
74          painCave.isFatal = 0;
75          simError();
76        }
77      
78 <      if (simParam->haveThermIntK()){
79 <        tIntK_ = simParam->getThermIntK();
78 >      if (simParam->haveThermodynamicIntegrationK()){
79 >        tIntK_ = simParam->getThermodynamicIntegrationK();
80        }
81        else{
82          tIntK_ = 1.0;
83          sprintf(painCave.errMsg,
84 <                "ThermoIntegration Warning: the tranformation parameter exponent\n"
85 <                "\t(k) was not specified. OOPSE will use a default value of %f.\n"
86 <                "\tTo set k, use the thermodynamicIntegrationK variable.\n",
84 >                "ThermoIntegration Warning: the tranformation parameter\n"
85 >                "\texponent (k) was not specified. OOPSE will use a default\n"
86 >                "\tvalue of %f. To set k, use the thermodynamicIntegrationK\n"
87 >                "\tvariable.\n",
88                  tIntK_);
89          painCave.isFatal = 0;
90          simError();      
# Line 144 | Line 151 | namespace oopse {
151      tempTau = curSnapshot->statData.getTau();
152      tempTau *= factor_;
153      curSnapshot->statData.setTau(tempTau);
154 <
155 <    // do crystal restraint forces for thermodynamic integration
154 > #ifndef IS_MPI
155 >    // do the single processor crystal restraint forces for
156 >    // thermodynamic integration
157      if (simParam->getUseSolidThermInt()) {
158        
159        lrPot_ += restraint_->Calc_Restraint_Forces();
# Line 154 | Line 162 | namespace oopse {
162        vHarm_ = restraint_->getVharm();
163        curSnapshot->statData[Stats::VHARM] = vHarm_;
164      }
165 <    
165 > #else
166 >    double tempLRPot = 0.0;
167 >    double tempVHarm = 0.0;
168 >    MPI_Status ierr;
169 >    int nproc;
170 >    MPI_Comm_size(MPI_COMM_WORLD, &nproc);
171 >    vHarm_ = 0.0;
172 >
173 >    // do the MPI crystal restraint forces for each processor
174 >    if (simParam->getUseSolidThermInt()) {
175 >      tempLRPot = restraint_->Calc_Restraint_Forces();
176 >      tempVHarm = restraint_->getVharm();
177 >    }
178 >
179 >    // master receives and accumulates the restraint info
180 >    if (worldRank == 0) {
181 >      for(int i = 0 ; i < nproc; ++i) {
182 >        if (i == worldRank) {
183 >          lrPot_ += tempLRPot;
184 >          vHarm_ += tempVHarm;
185 >        } else {
186 >          MPI_Recv(&tempLRPot, 1, MPI_REALTYPE, i,
187 >                   TAKE_THIS_TAG_REAL, MPI_COMM_WORLD, &ierr);
188 >          MPI_Recv(&tempVHarm, 1, MPI_REALTYPE, i,
189 >                   TAKE_THIS_TAG_REAL, MPI_COMM_WORLD, &ierr);
190 >          lrPot_ += tempLRPot;
191 >          vHarm_ += tempVHarm;
192 >        }
193 >      }
194 >
195 >      // give the final values to stats
196 >      curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot_;
197 >      curSnapshot->statData[Stats::VHARM] = vHarm_;
198 >
199 >    } else {
200 >      // pack up and send the appropriate info to the master node
201 >      for(int j = 1; j < nproc; ++j) {
202 >        if (worldRank == j) {
203 >
204 >          MPI_Send(&tempLRPot, 1, MPI_REALTYPE, 0,
205 >                   TAKE_THIS_TAG_REAL, MPI_COMM_WORLD);
206 >          MPI_Send(&tempVHarm, 1, MPI_REALTYPE, 0,
207 >                   TAKE_THIS_TAG_REAL, MPI_COMM_WORLD);
208 >        }
209 >      }
210 >    }
211 > #endif //is_mpi
212    }
213    
214   }

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines