ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/flucq/FluctuatingChargeObjectiveFunction.cpp
Revision: 1740
Committed: Tue Jun 5 18:01:50 2012 UTC (12 years, 10 months ago) by gezelter
Original Path: branches/development/src/flucq/FluctuatingChargeObjectiveFunction.cpp
File size: 5434 byte(s)
Log Message:
Adding the Flucq constraints and objective function for minimization

File Contents

# Content
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, 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 "flucq/FluctuatingChargeObjectiveFunction.hpp"
44
45 namespace OpenMD{
46
47 FluctuatingChargeObjectiveFunction::FluctuatingChargeObjectiveFunction(SimInfo* info, ForceManager* forceMan, FluctuatingChargeConstraints* fqConstraints)
48 : info_(info), forceMan_(forceMan), fqConstraints_(fqConstraints),
49 thermo(info) {
50 }
51
52 RealType FluctuatingChargeObjectiveFunction::value(const DynamicVector<RealType>& x) {
53 setCoor(x);
54 forceMan_->calcForces();
55 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
56 return curSnapshot->statData[Stats::ELECTROSTATIC_POTENTIAL];
57 }
58
59 void FluctuatingChargeObjectiveFunction::gradient(DynamicVector<RealType>& grad, const DynamicVector<RealType>& x) {
60 int shakeStatus;
61
62 setCoor(x);
63
64 forceMan_->calcForces();
65 fqConstraints_->applyConstraints();
66
67 getGrad(grad);
68 }
69
70 RealType FluctuatingChargeObjectiveFunction::valueAndGradient(DynamicVector<RealType>& grad,
71 const DynamicVector<RealType>& x) {
72
73 setCoor(x);
74
75 forceMan_->calcForces();
76 fqConstraints_->applyConstraints();
77
78 getGrad(grad);
79
80 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
81 return curSnapshot->statData[Stats::ELECTROSTATIC_POTENTIAL];
82 }
83
84 void FluctuatingChargeObjectiveFunction::setCoor(const DynamicVector<RealType> &x) const {
85 SimInfo::MoleculeIterator i;
86 Molecule::FluctuatingChargeIterator j;
87 Molecule* mol;
88 Atom* atom;
89 int index = 0;
90
91 for (mol = info_->beginMolecule(i); mol != NULL;
92 mol = info_->nextMolecule(i)) {
93
94 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
95 atom = mol->nextFluctuatingCharge(j)) {
96
97 atom->setFlucQPos(x[index++]);
98 }
99 }
100 }
101
102 void FluctuatingChargeObjectiveFunction::getGrad(DynamicVector<RealType> &grad) {
103 SimInfo::MoleculeIterator i;
104 Molecule::FluctuatingChargeIterator j;
105 Molecule* mol;
106 Atom* atom;
107
108 int index = 0;
109
110 for (mol = info_->beginMolecule(i); mol != NULL;
111 mol = info_->nextMolecule(i)) {
112
113 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
114 atom = mol->nextFluctuatingCharge(j)) {
115
116 grad[index++] = -atom->getFlucQFrc();
117 }
118 }
119 }
120
121 DynamicVector<RealType> FluctuatingChargeObjectiveFunction::setInitialCoords() {
122 SimInfo::MoleculeIterator i;
123 Molecule::FluctuatingChargeIterator j;
124 Molecule* mol;
125 Atom* atom;
126
127 int index = 0;
128 for (mol = info_->beginMolecule(i); mol != NULL;
129 mol = info_->nextMolecule(i)) {
130
131 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
132 atom = mol->nextFluctuatingCharge(j)) {
133
134 index++;
135 }
136 }
137
138 DynamicVector<RealType> initCoords(index);
139 index = 0;
140 for (mol = info_->beginMolecule(i); mol != NULL;
141 mol = info_->nextMolecule(i)) {
142
143 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
144 atom = mol->nextFluctuatingCharge(j)) {
145
146 initCoords[index++] = atom->getFlucQPos();
147 }
148 }
149 return initCoords;
150 }
151 }

Properties

Name Value
svn:eol-style native