ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/restraints/ThermoIntegrationForceManager.cpp
Revision: 2115
Committed: Fri Mar 11 00:43:19 2005 UTC (19 years, 4 months ago) by chrisfen
File size: 5989 byte(s)
Log Message:
fixed a bug in MPI restraints

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 namespace oopse {
52
53 ThermoIntegrationForceManager::ThermoIntegrationForceManager(SimInfo* info):
54 ForceManager(info){
55 currSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot();
56 simParam = info_->getSimParams();
57
58 if (simParam->haveThermIntLambda()){
59 tIntLambda_ = simParam->getThermIntLambda();
60 }
61 else{
62 tIntLambda_ = 1.0;
63 sprintf(painCave.errMsg,
64 "ThermoIntegration error: the transformation parameter (lambda) was\n"
65 "\tnot specified. OOPSE will use a default value of %f. To set\n"
66 "\tlambda, use the thermodynamicIntegrationLambda variable.\n",
67 tIntLambda_);
68 painCave.isFatal = 0;
69 simError();
70 }
71
72 if (simParam->haveThermIntK()){
73 tIntK_ = simParam->getThermIntK();
74 }
75 else{
76 tIntK_ = 1.0;
77 sprintf(painCave.errMsg,
78 "ThermoIntegration Warning: the tranformation parameter exponent\n"
79 "\t(k) was not specified. OOPSE will use a default value of %f.\n"
80 "\tTo set k, use the thermodynamicIntegrationK variable.\n",
81 tIntK_);
82 painCave.isFatal = 0;
83 simError();
84 }
85
86 if (simParam->getUseSolidThermInt()) {
87 // build a restraint object
88 restraint_ = new Restraints(info_, tIntLambda_, tIntK_);
89
90 }
91
92 // build the scaling factor used to modulate the forces and torques
93 factor_ = pow(tIntLambda_, tIntK_);
94
95 sprintf(painCave.errMsg,"%f is the factor\n",factor_);
96 painCave.isFatal = 0;
97 simError();
98
99 }
100
101 ThermoIntegrationForceManager::~ThermoIntegrationForceManager(){
102 }
103
104 void ThermoIntegrationForceManager::calcForces(bool needPotential,
105 bool needStress){
106 Snapshot* curSnapshot;
107 SimInfo::MoleculeIterator mi;
108 Molecule* mol;
109 Molecule::IntegrableObjectIterator ii;
110 StuntDouble* integrableObject;
111 Vector3d frc;
112 Vector3d trq;
113 Mat3x3d tempTau;
114
115 // perform the standard calcForces first
116 ForceManager::calcForces(needPotential, needStress);
117
118 curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
119
120 // now scale forces and torques of all the integrableObjects
121
122 for (mol = info_->beginMolecule(mi); mol != NULL;
123 mol = info_->nextMolecule(mi)) {
124 for (integrableObject = mol->beginIntegrableObject(ii);
125 integrableObject != NULL;
126 integrableObject = mol->nextIntegrableObject(ii)) {
127 frc = integrableObject->getFrc();
128 frc *= factor_;
129 integrableObject->setFrc(frc);
130
131 if (integrableObject->isDirectional()){
132 trq = integrableObject->getTrq();
133 trq *= factor_;
134 integrableObject->setTrq(trq);
135 }
136 }
137 }
138
139 // set vraw to be the unmodulated potential
140 lrPot_ = curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL];
141 curSnapshot->statData[Stats::VRAW] = lrPot_;
142
143 // modulate the potential and update the snapshot
144 lrPot_ *= factor_;
145 curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot_;
146
147 // scale the pressure tensor
148 tempTau = curSnapshot->statData.getTau();
149 tempTau *= factor_;
150 curSnapshot->statData.setTau(tempTau);
151
152 // do crystal restraint forces for thermodynamic integration
153 if (simParam->getUseSolidThermInt()) {
154
155 lrPot_ += restraint_->Calc_Restraint_Forces();
156 curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot_;
157
158 vHarm_ = restraint_->getVharm();
159 curSnapshot->statData[Stats::VHARM] = vHarm_;
160 }
161
162 }
163
164 }