ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/restraints/ThermoIntegrationForceManager.cpp
Revision: 2103
Committed: Thu Mar 10 16:15:13 2005 UTC (19 years, 4 months ago) by chrisfen
File size: 5748 byte(s)
Log Message:
fixing 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 "primitives/Molecule.hpp"
46 #include "utils/simError.h"
47 #include "utils/OOPSEConstant.hpp"
48 #include "utils/StringUtils.hpp"
49
50 namespace oopse {
51
52 ThermoIntegrationForceManager::ThermoIntegrationForceManager(SimInfo* info):
53 ForceManager(info){
54 currSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot();
55 simParam = info_->getSimParams();
56
57 if (simParam->haveThermIntLambda()){
58 tIntLambda_ = simParam->getThermIntLambda();
59 }
60 else{
61 tIntLambda_ = 1.0;
62 sprintf(painCave.errMsg,
63 "ThermoIntegration error: the transformation parameter (lambda) was\n"
64 "\tnot specified. OOPSE will use a default value of %f. To set\n"
65 "\tlambda, use the thermodynamicIntegrationLambda variable.\n",
66 tIntLambda_);
67 painCave.isFatal = 0;
68 simError();
69 }
70
71 if (simParam->haveThermIntK()){
72 tIntK_ = simParam->getThermIntK();
73 }
74 else{
75 tIntK_ = 1.0;
76 sprintf(painCave.errMsg,
77 "ThermoIntegration Warning: the tranformation parameter exponent\n"
78 "\t(k) was not specified. OOPSE will use a default value of %f.\n"
79 "\tTo set k, use the thermodynamicIntegrationK variable.\n",
80 tIntK_);
81 painCave.isFatal = 0;
82 simError();
83 }
84
85 if (simParam->getUseSolidThermInt()) {
86 // build a restraint object
87 restraint_ = new Restraints(info_, tIntLambda_, tIntK_);
88
89 }
90
91 // build the scaling factor used to modulate the forces and torques
92 factor_ = pow(tIntLambda_, tIntK_);
93
94 printf("%f is the factor\n",factor_);
95
96 }
97
98 ThermoIntegrationForceManager::~ThermoIntegrationForceManager(){
99 }
100
101 void ThermoIntegrationForceManager::calcForces(bool needPotential,
102 bool needStress){
103 Snapshot* curSnapshot;
104 SimInfo::MoleculeIterator mi;
105 Molecule* mol;
106 Molecule::IntegrableObjectIterator ii;
107 StuntDouble* integrableObject;
108 Vector3d frc;
109 Vector3d trq;
110
111 // perform the standard calcForces first
112 ForceManager::calcForces(needPotential, needStress);
113
114 curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
115
116 // now scale forces and torques of all the integrableObjects
117
118 for (mol = info_->beginMolecule(mi); mol != NULL;
119 mol = info_->nextMolecule(mi)) {
120 for (integrableObject = mol->beginIntegrableObject(ii);
121 integrableObject != NULL;
122 integrableObject = mol->nextIntegrableObject(ii)) {
123 frc = integrableObject->getFrc();
124 frc *= factor_;
125 integrableObject->setFrc(frc);
126
127 if (integrableObject->isDirectional()){
128 trq = integrableObject->getTrq();
129 trq *= factor_;
130 integrableObject->setTrq(trq);
131 }
132 }
133
134 // set vraw to be the unmodulated potential
135 lrPot_ = curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL];
136 curSnapshot->statData[Stats::VRAW] = lrPot_;
137
138 // modulate the potential and update the snapshot
139 lrPot_ *= factor_;
140 curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot_;
141 }
142
143 // do crystal restraint forces for thermodynamic integration
144 if (simParam->getUseSolidThermInt()) {
145
146 lrPot_ += restraint_->Calc_Restraint_Forces();
147 curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot_;
148
149 vHarm_ = restraint_->getVharm();
150 curSnapshot->statData[Stats::VHARM] = vHarm_;
151 }
152
153 }
154
155 }