ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/flucq/FluctuatingChargeObjectiveFunction.cpp
Revision: 1927
Committed: Wed Aug 14 20:19:19 2013 UTC (11 years, 8 months ago) by gezelter
File size: 5753 byte(s)
Log Message:
FlucRho/EAM initial commit

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

Properties

Name Value
svn:eol-style native