--- trunk/src/restraints/ThermoIntegrationForceManager.cpp 2005/04/15 22:04:00 507 +++ trunk/src/restraints/ThermoIntegrationForceManager.cpp 2006/08/31 22:34:49 1029 @@ -48,6 +48,11 @@ #include "utils/OOPSEConstant.hpp" #include "utils/StringUtils.hpp" +#ifdef IS_MPI +#include +#define TAKE_THIS_TAG_REAL 2 +#endif //is_mpi + namespace oopse { ThermoIntegrationForceManager::ThermoIntegrationForceManager(SimInfo* info): @@ -55,29 +60,31 @@ namespace oopse { currSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot(); simParam = info_->getSimParams(); - if (simParam->haveThermIntLambda()){ - tIntLambda_ = simParam->getThermIntLambda(); + if (simParam->haveThermodynamicIntegrationLambda()){ + tIntLambda_ = simParam->getThermodynamicIntegrationLambda(); } else{ tIntLambda_ = 1.0; sprintf(painCave.errMsg, - "ThermoIntegration error: the transformation parameter (lambda) was\n" - "\tnot specified. OOPSE will use a default value of %f. To set\n" - "\tlambda, use the thermodynamicIntegrationLambda variable.\n", + "ThermoIntegration error: the transformation parameter\n" + "\t(lambda) was not specified. OOPSE will use a default\n" + "\tvalue of %f. To set lambda, use the \n" + "\tthermodynamicIntegrationLambda variable.\n", tIntLambda_); painCave.isFatal = 0; simError(); } - if (simParam->haveThermIntK()){ - tIntK_ = simParam->getThermIntK(); + if (simParam->haveThermodynamicIntegrationK()){ + tIntK_ = simParam->getThermodynamicIntegrationK(); } else{ tIntK_ = 1.0; sprintf(painCave.errMsg, - "ThermoIntegration Warning: the tranformation parameter exponent\n" - "\t(k) was not specified. OOPSE will use a default value of %f.\n" - "\tTo set k, use the thermodynamicIntegrationK variable.\n", + "ThermoIntegration Warning: the tranformation parameter\n" + "\texponent (k) was not specified. OOPSE will use a default\n" + "\tvalue of %f. To set k, use the thermodynamicIntegrationK\n" + "\tvariable.\n", tIntK_); painCave.isFatal = 0; simError(); @@ -144,8 +151,9 @@ namespace oopse { tempTau = curSnapshot->statData.getTau(); tempTau *= factor_; curSnapshot->statData.setTau(tempTau); - - // do crystal restraint forces for thermodynamic integration +#ifndef IS_MPI + // do the single processor crystal restraint forces for + // thermodynamic integration if (simParam->getUseSolidThermInt()) { lrPot_ += restraint_->Calc_Restraint_Forces(); @@ -154,7 +162,53 @@ namespace oopse { vHarm_ = restraint_->getVharm(); curSnapshot->statData[Stats::VHARM] = vHarm_; } - +#else + double tempLRPot = 0.0; + double tempVHarm = 0.0; + MPI_Status ierr; + int nproc; + MPI_Comm_size(MPI_COMM_WORLD, &nproc); + vHarm_ = 0.0; + + // do the MPI crystal restraint forces for each processor + if (simParam->getUseSolidThermInt()) { + tempLRPot = restraint_->Calc_Restraint_Forces(); + tempVHarm = restraint_->getVharm(); + } + + // master receives and accumulates the restraint info + if (worldRank == 0) { + for(int i = 0 ; i < nproc; ++i) { + if (i == worldRank) { + lrPot_ += tempLRPot; + vHarm_ += tempVHarm; + } else { + MPI_Recv(&tempLRPot, 1, MPI_REALTYPE, i, + TAKE_THIS_TAG_REAL, MPI_COMM_WORLD, &ierr); + MPI_Recv(&tempVHarm, 1, MPI_REALTYPE, i, + TAKE_THIS_TAG_REAL, MPI_COMM_WORLD, &ierr); + lrPot_ += tempLRPot; + vHarm_ += tempVHarm; + } + } + + // give the final values to stats + curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot_; + curSnapshot->statData[Stats::VHARM] = vHarm_; + + } else { + // pack up and send the appropriate info to the master node + for(int j = 1; j < nproc; ++j) { + if (worldRank == j) { + + MPI_Send(&tempLRPot, 1, MPI_REALTYPE, 0, + TAKE_THIS_TAG_REAL, MPI_COMM_WORLD); + MPI_Send(&tempVHarm, 1, MPI_REALTYPE, 0, + TAKE_THIS_TAG_REAL, MPI_COMM_WORLD); + } + } + } +#endif //is_mpi } }