| 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 | } | 
| 96 |  | 
| 97 | ThermoIntegrationForceManager::~ThermoIntegrationForceManager(){ | 
| 98 | } | 
| 99 |  | 
| 100 | void ThermoIntegrationForceManager::calcForces(bool needPotential, | 
| 101 | bool needStress){ | 
| 102 | Snapshot* curSnapshot; | 
| 103 | SimInfo::MoleculeIterator mi; | 
| 104 | Molecule* mol; | 
| 105 | Molecule::IntegrableObjectIterator ii; | 
| 106 | StuntDouble* integrableObject; | 
| 107 | Vector3d frc; | 
| 108 | Vector3d trq; | 
| 109 | Mat3x3d tempTau; | 
| 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 |  | 
| 135 | // set vraw to be the unmodulated potential | 
| 136 | lrPot_ = curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL]; | 
| 137 | curSnapshot->statData[Stats::VRAW] = lrPot_; | 
| 138 |  | 
| 139 | // modulate the potential and update the snapshot | 
| 140 | lrPot_ *= factor_; | 
| 141 | curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot_; | 
| 142 |  | 
| 143 | // scale the pressure tensor | 
| 144 | tempTau = curSnapshot->statData.getTau(); | 
| 145 | tempTau *= factor_; | 
| 146 | curSnapshot->statData.setTau(tempTau); | 
| 147 |  | 
| 148 | // do crystal restraint forces for thermodynamic integration | 
| 149 | if (simParam->getUseSolidThermInt()) { | 
| 150 |  | 
| 151 | lrPot_ += restraint_->Calc_Restraint_Forces(); | 
| 152 | curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot_; | 
| 153 |  | 
| 154 | vHarm_ = restraint_->getVharm(); | 
| 155 | curSnapshot->statData[Stats::VHARM] = vHarm_; | 
| 156 | } | 
| 157 |  | 
| 158 | } | 
| 159 |  | 
| 160 | } |