ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/restraints/ThermoIntegrationForceManager.cpp
Revision: 998
Committed: Mon Jul 3 13:18:43 2006 UTC (18 years, 10 months ago) by chrisfen
Original Path: trunk/src/restraints/ThermoIntegrationForceManager.cpp
File size: 5843 byte(s)
Log Message:
Added simulation box dipole moment accumulation for the purposes of calculating dielectric constants

File Contents

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