ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/flucq/FluctuatingChargeObjectiveFunction.cpp
Revision: 1761
Committed: Fri Jun 22 20:01:37 2012 UTC (12 years, 10 months ago) by gezelter
File size: 5662 byte(s)
Log Message:
fixes for fluctuating charges

File Contents

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

Properties

Name Value
svn:eol-style native