ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/flucq/FluctuatingChargeObjectiveFunction.cpp
Revision: 1767
Committed: Fri Jul 6 22:01:58 2012 UTC (12 years, 9 months ago) by gezelter
File size: 5641 byte(s)
Log Message:
Various fixes required to compile OpenMD with the MS Visual C++ compiler

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
56 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
57 potVec pot = curSnapshot->getLongRangePotentials();
58 potVec exPot = curSnapshot->getExcludedPotentials();
59
60 return pot[ELECTROSTATIC_FAMILY] + exPot[ELECTROSTATIC_FAMILY];
61 }
62
63 void FluctuatingChargeObjectiveFunction::gradient(DynamicVector<RealType>& grad, const DynamicVector<RealType>& x) {
64
65 setCoor(x);
66
67 forceMan_->calcForces();
68 fqConstraints_->applyConstraints();
69
70 getGrad(grad);
71 }
72
73 RealType FluctuatingChargeObjectiveFunction::valueAndGradient(DynamicVector<RealType>& grad,
74 const DynamicVector<RealType>& x) {
75
76 setCoor(x);
77
78 forceMan_->calcForces();
79 fqConstraints_->applyConstraints();
80
81 getGrad(grad);
82
83 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
84 potVec pot = curSnapshot->getLongRangePotentials();
85 potVec exPot = curSnapshot->getExcludedPotentials();
86
87 return pot[ELECTROSTATIC_FAMILY] + exPot[ELECTROSTATIC_FAMILY];
88 }
89
90 void FluctuatingChargeObjectiveFunction::setCoor(const DynamicVector<RealType> &x) const {
91 SimInfo::MoleculeIterator i;
92 Molecule::FluctuatingChargeIterator j;
93 Molecule* mol;
94 Atom* atom;
95 int index = 0;
96
97 for (mol = info_->beginMolecule(i); mol != NULL;
98 mol = info_->nextMolecule(i)) {
99
100 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
101 atom = mol->nextFluctuatingCharge(j)) {
102
103 atom->setFlucQPos(x[index++]);
104 }
105 }
106 }
107
108 void FluctuatingChargeObjectiveFunction::getGrad(DynamicVector<RealType> &grad) {
109 SimInfo::MoleculeIterator i;
110 Molecule::FluctuatingChargeIterator j;
111 Molecule* mol;
112 Atom* atom;
113
114 int index = 0;
115
116 for (mol = info_->beginMolecule(i); mol != NULL;
117 mol = info_->nextMolecule(i)) {
118
119 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
120 atom = mol->nextFluctuatingCharge(j)) {
121
122 grad[index++] = -atom->getFlucQFrc();
123 }
124 }
125 }
126
127 DynamicVector<RealType> FluctuatingChargeObjectiveFunction::setInitialCoords() {
128 SimInfo::MoleculeIterator i;
129 Molecule::FluctuatingChargeIterator j;
130 Molecule* mol;
131 Atom* atom;
132
133 int index = 0;
134 for (mol = info_->beginMolecule(i); mol != NULL;
135 mol = info_->nextMolecule(i)) {
136
137 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
138 atom = mol->nextFluctuatingCharge(j)) {
139
140 index++;
141 }
142 }
143
144 DynamicVector<RealType> initCoords(index);
145 index = 0;
146 for (mol = info_->beginMolecule(i); mol != NULL;
147 mol = info_->nextMolecule(i)) {
148
149 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
150 atom = mol->nextFluctuatingCharge(j)) {
151
152 initCoords[index++] = atom->getFlucQPos();
153 }
154 }
155 return initCoords;
156 }
157 }

Properties

Name Value
svn:eol-style native