ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/restraints/ThermoIntegrationForceManager.cpp
Revision: 2990
Committed: Thu Aug 31 20:12:12 2006 UTC (17 years, 9 months ago) by chrisfen
File size: 7274 byte(s)
Log Message:
converted samples over to the OOPSE-4 file format

Fixed an MPI bug in thermodynamic integration

File Contents

# Content
1 /*
2 * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3 *
4 * The University of Notre Dame grants you ("Licensee") a
5 * non-exclusive, royalty free, license to use, modify and
6 * redistribute this software in source and binary code form, provided
7 * that the following conditions are met:
8 *
9 * 1. Acknowledgement of the program authors must be made in any
10 * publication of scientific results based in part on use of the
11 * program. An acceptable form of acknowledgement is citation of
12 * the article in which the program was described (Matthew
13 * A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 * J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 * Parallel Simulation Engine for Molecular Dynamics,"
16 * J. Comput. Chem. 26, pp. 252-271 (2005))
17 *
18 * 2. Redistributions of source code must retain the above copyright
19 * notice, this list of conditions and the following disclaimer.
20 *
21 * 3. Redistributions in binary form must reproduce the above copyright
22 * notice, this list of conditions and the following disclaimer in the
23 * documentation and/or other materials provided with the
24 * distribution.
25 *
26 * This software is provided "AS IS," without a warranty of any
27 * kind. All express or implied conditions, representations and
28 * warranties, including any implied warranty of merchantability,
29 * fitness for a particular purpose or non-infringement, are hereby
30 * excluded. The University of Notre Dame and its licensors shall not
31 * be liable for any damages suffered by licensee as a result of
32 * using, modifying or distributing the software or its
33 * derivatives. In no event will the University of Notre Dame or its
34 * licensors be liable for any lost revenue, profit or data, or for
35 * direct, indirect, special, consequential, incidental or punitive
36 * damages, however caused and regardless of the theory of liability,
37 * arising out of the use of or inability to use software, even if the
38 * University of Notre Dame has been advised of the possibility of
39 * such damages.
40 */
41
42 #include <cmath>
43 #include "restraints/ThermoIntegrationForceManager.hpp"
44 #include "integrators/Integrator.hpp"
45 #include "math/SquareMatrix3.hpp"
46 #include "primitives/Molecule.hpp"
47 #include "utils/simError.h"
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):
59 ForceManager(info){
60 currSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot();
61 simParam = info_->getSimParams();
62
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\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->haveThermodynamicIntegrationK()){
79 tIntK_ = simParam->getThermodynamicIntegrationK();
80 }
81 else{
82 tIntK_ = 1.0;
83 sprintf(painCave.errMsg,
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();
91 }
92
93 if (simParam->getUseSolidThermInt()) {
94 // build a restraint object
95 restraint_ = new Restraints(info_, tIntLambda_, tIntK_);
96
97 }
98
99 // build the scaling factor used to modulate the forces and torques
100 factor_ = pow(tIntLambda_, tIntK_);
101
102 }
103
104 ThermoIntegrationForceManager::~ThermoIntegrationForceManager(){
105 }
106
107 void ThermoIntegrationForceManager::calcForces(bool needPotential,
108 bool needStress){
109 Snapshot* curSnapshot;
110 SimInfo::MoleculeIterator mi;
111 Molecule* mol;
112 Molecule::IntegrableObjectIterator ii;
113 StuntDouble* integrableObject;
114 Vector3d frc;
115 Vector3d trq;
116 Mat3x3d tempTau;
117
118 // perform the standard calcForces first
119 ForceManager::calcForces(needPotential, needStress);
120
121 curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
122
123 // now scale forces and torques of all the integrableObjects
124
125 for (mol = info_->beginMolecule(mi); mol != NULL;
126 mol = info_->nextMolecule(mi)) {
127 for (integrableObject = mol->beginIntegrableObject(ii);
128 integrableObject != NULL;
129 integrableObject = mol->nextIntegrableObject(ii)) {
130 frc = integrableObject->getFrc();
131 frc *= factor_;
132 integrableObject->setFrc(frc);
133
134 if (integrableObject->isDirectional()){
135 trq = integrableObject->getTrq();
136 trq *= factor_;
137 integrableObject->setTrq(trq);
138 }
139 }
140 }
141
142 // set vraw to be the unmodulated potential
143 lrPot_ = curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL];
144 curSnapshot->statData[Stats::VRAW] = lrPot_;
145
146 // modulate the potential and update the snapshot
147 lrPot_ *= factor_;
148 curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot_;
149
150 // scale the pressure tensor
151 tempTau = curSnapshot->statData.getTau();
152 tempTau *= factor_;
153 curSnapshot->statData.setTau(tempTau);
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();
160 curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot_;
161
162 vHarm_ = restraint_->getVharm();
163 curSnapshot->statData[Stats::VHARM] = vHarm_;
164 }
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
172 // do the MPI crystal restraint forces for each processor
173 if (simParam->getUseSolidThermInt()) {
174 tempLRPot = restraint_->Calc_Restraint_Forces();
175 tempVHarm = restraint_->getVharm();
176 }
177
178 // master receives and accumulates the restraint info
179 if (worldRank == 0) {
180 for(int i = 0 ; i < nproc; ++i) {
181 if (i == worldRank) {
182 lrPot_ += tempLRPot;
183 vHarm_ += tempVHarm;
184 } else {
185 MPI_Recv(&tempLRPot, 1, MPI_REALTYPE, i,
186 TAKE_THIS_TAG_REAL, MPI_COMM_WORLD, &ierr);
187 MPI_Recv(&tempVHarm, 1, MPI_REALTYPE, i,
188 TAKE_THIS_TAG_REAL, MPI_COMM_WORLD, &ierr);
189 lrPot_ += tempLRPot;
190 vHarm_ += tempVHarm;
191 }
192 }
193
194 // give the final values to stats
195 curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot_;
196 curSnapshot->statData[Stats::VHARM] = vHarm_;
197
198 } else {
199 // pack up and send the appropriate info to the master node
200 for(int j = 1; j < nproc; ++j) {
201 if (worldRank == j) {
202
203 MPI_Send(&tempLRPot, 1, MPI_REALTYPE, 0,
204 TAKE_THIS_TAG_REAL, MPI_COMM_WORLD);
205 MPI_Send(&tempVHarm, 1, MPI_REALTYPE, 0,
206 TAKE_THIS_TAG_REAL, MPI_COMM_WORLD);
207 }
208 }
209 }
210 #endif //is_mpi
211 }
212
213 }