ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/optimization/PotentialEnergyObjectiveFunction.cpp
Revision: 1741
Committed: Tue Jun 5 18:02:44 2012 UTC (12 years, 10 months ago) by gezelter
File size: 5912 byte(s)
Log Message:
Adding initial import of optimization library

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

Properties

Name Value
svn:eol-style native