ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/optimization/PotentialEnergyObjectiveFunction.cpp
Revision: 1750
Committed: Thu Jun 7 12:53:46 2012 UTC (12 years, 10 months ago) by gezelter
File size: 6258 byte(s)
Log Message:
Fixing some bugs in optimization, fixing status functions so that they
dump correctly (although some things are deferred until the Stats is
accumulator-based).

File Contents

# User Rev Content
1 gezelter 1741 /*
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, 24107 (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 gezelter 1748 : info_(info), forceMan_(forceMan), thermo(info) {
49     shake_ = new Shake(info_);
50 gezelter 1741 }
51    
52    
53    
54     RealType PotentialEnergyObjectiveFunction::value(const DynamicVector<RealType>& x) {
55     setCoor(x);
56 gezelter 1748 shake_->constraintR();
57 gezelter 1741 forceMan_->calcForces();
58 gezelter 1748 shake_->constraintF();
59 gezelter 1741 return thermo.getPotential();
60     }
61    
62     void PotentialEnergyObjectiveFunction::gradient(DynamicVector<RealType>& grad, const DynamicVector<RealType>& x) {
63    
64 gezelter 1748 setCoor(x);
65     shake_->constraintR();
66     forceMan_->calcForces();
67     shake_->constraintF();
68     getGrad(grad);
69 gezelter 1741 }
70    
71     RealType PotentialEnergyObjectiveFunction::valueAndGradient(DynamicVector<RealType>& grad,
72     const DynamicVector<RealType>& x) {
73 gezelter 1748 setCoor(x);
74     shake_->constraintR();
75     forceMan_->calcForces();
76     shake_->constraintF();
77 gezelter 1741 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* integrableObject;
88     int index = 0;
89    
90     for (mol = info_->beginMolecule(i); mol != NULL;
91     mol = info_->nextMolecule(i)) {
92     for (integrableObject = mol->beginIntegrableObject(j);
93     integrableObject != NULL;
94     integrableObject = mol->nextIntegrableObject(j)) {
95    
96     position[0] = x[index++];
97     position[1] = x[index++];
98     position[2] = x[index++];
99    
100     integrableObject->setPos(position);
101    
102     if (integrableObject->isDirectional()) {
103     eulerAngle[0] = x[index++];
104     eulerAngle[1] = x[index++];
105     eulerAngle[2] = x[index++];
106    
107     integrableObject->setEuler(eulerAngle);
108 gezelter 1748
109     if (integrableObject->isRigidBody()) {
110     RigidBody* rb = static_cast<RigidBody*>(integrableObject);
111     rb->updateAtoms();
112     }
113     }
114 gezelter 1741 }
115     }
116     }
117    
118     void PotentialEnergyObjectiveFunction::getGrad(DynamicVector<RealType> &grad) {
119     SimInfo::MoleculeIterator i;
120     Molecule::IntegrableObjectIterator j;
121     Molecule* mol;
122     StuntDouble* integrableObject;
123     std::vector<RealType> myGrad;
124    
125     int index = 0;
126    
127     for (mol = info_->beginMolecule(i); mol != NULL;
128     mol = info_->nextMolecule(i)) {
129     for (integrableObject = mol->beginIntegrableObject(j);
130     integrableObject != NULL;
131     integrableObject = mol->nextIntegrableObject(j)) {
132     myGrad = integrableObject->getGrad();
133 gezelter 1750
134     for (size_t k = 0; k < myGrad.size(); ++k) {
135 gezelter 1741 grad[index++] = myGrad[k];
136     }
137     }
138     }
139     }
140    
141     DynamicVector<RealType> PotentialEnergyObjectiveFunction::setInitialCoords() {
142     SimInfo::MoleculeIterator i;
143     Molecule::IntegrableObjectIterator j;
144     Molecule* mol;
145     StuntDouble* integrableObject;
146    
147     Vector3d pos;
148     Vector3d eulerAngle;
149    
150     DynamicVector<RealType> xinit(info_->getNdfLocal(), 0.0);
151    
152     int index = 0;
153    
154     for (mol = info_->beginMolecule(i); mol != NULL;
155     mol = info_->nextMolecule(i)) {
156     for (integrableObject = mol->beginIntegrableObject(j);
157     integrableObject != NULL;
158     integrableObject = mol->nextIntegrableObject(j)) {
159    
160     pos = integrableObject->getPos();
161    
162     xinit[index++] = pos[0];
163     xinit[index++] = pos[1];
164     xinit[index++] = pos[2];
165    
166     if (integrableObject->isDirectional()) {
167     eulerAngle = integrableObject->getEuler();
168     xinit[index++] = eulerAngle[0];
169     xinit[index++] = eulerAngle[1];
170     xinit[index++] = eulerAngle[2];
171     }
172     }
173     }
174     return xinit;
175     }
176     }

Properties

Name Value
svn:eol-style native