| 1 |
/* |
| 2 |
* Copyright (c) 2012 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. Redistributions of source code must retain the above copyright |
| 10 |
* notice, this list of conditions and the following disclaimer. |
| 11 |
* |
| 12 |
* 2. Redistributions in binary form must reproduce the above copyright |
| 13 |
* notice, this list of conditions and the following disclaimer in the |
| 14 |
* documentation and/or other materials provided with the |
| 15 |
* distribution. |
| 16 |
* |
| 17 |
* This software is provided "AS IS," without a warranty of any |
| 18 |
* kind. All express or implied conditions, representations and |
| 19 |
* warranties, including any implied warranty of merchantability, |
| 20 |
* fitness for a particular purpose or non-infringement, are hereby |
| 21 |
* excluded. The University of Notre Dame and its licensors shall not |
| 22 |
* be liable for any damages suffered by licensee as a result of |
| 23 |
* using, modifying or distributing the software or its |
| 24 |
* derivatives. In no event will the University of Notre Dame or its |
| 25 |
* licensors be liable for any lost revenue, profit or data, or for |
| 26 |
* direct, indirect, special, consequential, incidental or punitive |
| 27 |
* damages, however caused and regardless of the theory of liability, |
| 28 |
* arising out of the use of or inability to use software, even if the |
| 29 |
* University of Notre Dame has been advised of the possibility of |
| 30 |
* such damages. |
| 31 |
* |
| 32 |
* SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your |
| 33 |
* research, please cite the appropriate papers when you publish your |
| 34 |
* work. Good starting points are: |
| 35 |
* |
| 36 |
* [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). |
| 37 |
* [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). |
| 38 |
* [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008). |
| 39 |
* [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010). |
| 40 |
* [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011). |
| 41 |
*/ |
| 42 |
|
| 43 |
#include "optimization/PotentialEnergyObjectiveFunction.hpp" |
| 44 |
|
| 45 |
namespace OpenMD{ |
| 46 |
|
| 47 |
PotentialEnergyObjectiveFunction::PotentialEnergyObjectiveFunction(SimInfo* info, ForceManager* forceMan) |
| 48 |
: info_(info), forceMan_(forceMan), thermo(info) { |
| 49 |
shake_ = new Shake(info_); |
| 50 |
} |
| 51 |
|
| 52 |
|
| 53 |
|
| 54 |
RealType PotentialEnergyObjectiveFunction::value(const DynamicVector<RealType>& x) { |
| 55 |
setCoor(x); |
| 56 |
shake_->constraintR(); |
| 57 |
forceMan_->calcForces(); |
| 58 |
shake_->constraintF(); |
| 59 |
return thermo.getPotential(); |
| 60 |
} |
| 61 |
|
| 62 |
void PotentialEnergyObjectiveFunction::gradient(DynamicVector<RealType>& grad, const DynamicVector<RealType>& x) { |
| 63 |
|
| 64 |
setCoor(x); |
| 65 |
shake_->constraintR(); |
| 66 |
forceMan_->calcForces(); |
| 67 |
shake_->constraintF(); |
| 68 |
getGrad(grad); |
| 69 |
} |
| 70 |
|
| 71 |
RealType PotentialEnergyObjectiveFunction::valueAndGradient(DynamicVector<RealType>& grad, |
| 72 |
const DynamicVector<RealType>& x) { |
| 73 |
setCoor(x); |
| 74 |
shake_->constraintR(); |
| 75 |
forceMan_->calcForces(); |
| 76 |
shake_->constraintF(); |
| 77 |
getGrad(grad); |
| 78 |
return thermo.getPotential(); |
| 79 |
} |
| 80 |
|
| 81 |
void PotentialEnergyObjectiveFunction::setCoor(const DynamicVector<RealType> &x) const { |
| 82 |
Vector3d position; |
| 83 |
Vector3d eulerAngle; |
| 84 |
SimInfo::MoleculeIterator i; |
| 85 |
Molecule::IntegrableObjectIterator j; |
| 86 |
Molecule* mol; |
| 87 |
StuntDouble* sd; |
| 88 |
int index = 0; |
| 89 |
|
| 90 |
for (mol = info_->beginMolecule(i); mol != NULL; |
| 91 |
mol = info_->nextMolecule(i)) { |
| 92 |
|
| 93 |
for (sd = mol->beginIntegrableObject(j); sd != NULL; |
| 94 |
sd = mol->nextIntegrableObject(j)) { |
| 95 |
|
| 96 |
position[0] = x[index++]; |
| 97 |
position[1] = x[index++]; |
| 98 |
position[2] = x[index++]; |
| 99 |
|
| 100 |
sd->setPos(position); |
| 101 |
|
| 102 |
if (sd->isDirectional()) { |
| 103 |
eulerAngle[0] = x[index++]; |
| 104 |
eulerAngle[1] = x[index++]; |
| 105 |
eulerAngle[2] = x[index++]; |
| 106 |
|
| 107 |
sd->setEuler(eulerAngle); |
| 108 |
|
| 109 |
if (sd->isRigidBody()) { |
| 110 |
RigidBody* rb = static_cast<RigidBody*>(sd); |
| 111 |
rb->updateAtoms(); |
| 112 |
} |
| 113 |
|
| 114 |
} |
| 115 |
} |
| 116 |
} |
| 117 |
} |
| 118 |
|
| 119 |
void PotentialEnergyObjectiveFunction::getGrad(DynamicVector<RealType> &grad) { |
| 120 |
SimInfo::MoleculeIterator i; |
| 121 |
Molecule::IntegrableObjectIterator j; |
| 122 |
Molecule* mol; |
| 123 |
StuntDouble* sd; |
| 124 |
std::vector<RealType> myGrad; |
| 125 |
|
| 126 |
int index = 0; |
| 127 |
|
| 128 |
for (mol = info_->beginMolecule(i); mol != NULL; |
| 129 |
mol = info_->nextMolecule(i)) { |
| 130 |
|
| 131 |
for (sd = mol->beginIntegrableObject(j); sd != NULL; |
| 132 |
sd = mol->nextIntegrableObject(j)) { |
| 133 |
|
| 134 |
myGrad = sd->getGrad(); |
| 135 |
|
| 136 |
for (size_t k = 0; k < myGrad.size(); ++k) { |
| 137 |
grad[index++] = myGrad[k]; |
| 138 |
} |
| 139 |
|
| 140 |
} |
| 141 |
} |
| 142 |
} |
| 143 |
|
| 144 |
DynamicVector<RealType> PotentialEnergyObjectiveFunction::setInitialCoords() { |
| 145 |
SimInfo::MoleculeIterator i; |
| 146 |
Molecule::IntegrableObjectIterator j; |
| 147 |
Molecule* mol; |
| 148 |
StuntDouble* sd; |
| 149 |
|
| 150 |
Vector3d pos; |
| 151 |
Vector3d eulerAngle; |
| 152 |
|
| 153 |
DynamicVector<RealType> xinit(info_->getNdfLocal(), 0.0); |
| 154 |
|
| 155 |
int index = 0; |
| 156 |
|
| 157 |
for (mol = info_->beginMolecule(i); mol != NULL; |
| 158 |
mol = info_->nextMolecule(i)) { |
| 159 |
|
| 160 |
for (sd = mol->beginIntegrableObject(j); sd != NULL; |
| 161 |
sd = mol->nextIntegrableObject(j)) { |
| 162 |
|
| 163 |
pos = sd->getPos(); |
| 164 |
|
| 165 |
xinit[index++] = pos[0]; |
| 166 |
xinit[index++] = pos[1]; |
| 167 |
xinit[index++] = pos[2]; |
| 168 |
|
| 169 |
if (sd->isDirectional()) { |
| 170 |
eulerAngle = sd->getEuler(); |
| 171 |
xinit[index++] = eulerAngle[0]; |
| 172 |
xinit[index++] = eulerAngle[1]; |
| 173 |
xinit[index++] = eulerAngle[2]; |
| 174 |
} |
| 175 |
|
| 176 |
} |
| 177 |
} |
| 178 |
return xinit; |
| 179 |
} |
| 180 |
} |