ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/perturbations/UniformGradient.cpp
Revision: 2026
Committed: Wed Oct 22 12:23:59 2014 UTC (10 years, 6 months ago) by gezelter
File size: 6652 byte(s)
Log Message:
Starting to add support for UniformGradient. 
Changed Vector3d input type to a more general std::vector<RealType> input.  This change alters RNEMD and UniformField inputs.

File Contents

# User Rev Content
1 gezelter 2026 /*
2     * Copyright (c) 2014 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 "perturbations/UniformGradient.hpp"
44     #include "types/FixedChargeAdapter.hpp"
45     #include "types/FluctuatingChargeAdapter.hpp"
46     #include "types/MultipoleAdapter.hpp"
47     #include "primitives/Molecule.hpp"
48     #include "nonbonded/NonBondedInteraction.hpp"
49     #include "utils/PhysicalConstants.hpp"
50    
51     namespace OpenMD {
52    
53     UniformGradient::UniformGradient(SimInfo* info) : info_(info),
54     doUniformGradient(false),
55     doParticlePot(false),
56     initialized(false) {
57     simParams = info_->getSimParams();
58     }
59    
60     void UniformGradient::initialize() {
61     if (simParams->haveUniformGradient()) {
62     doUniformGradient = true;
63     std::vector<RealType> pv = simParams->getUniformGradient();
64     if (pv.size() != 5) {
65     sprintf(painCave.errMsg,
66     "UniformGradient: Incorrect number of parameters specified.\n"
67     "\tthere should be 5 parameters, but %lu were specified.\n", pv.size());
68     painCave.isFatal = 1;
69     simError();
70     }
71     pars_.a = pv[0];
72     pars_.b = pv[1];
73     pars_.c = pv[2];
74     pars_.alpha = pv[3];
75     pars_.beta = pv[4];
76    
77     Grad_(0,0) = pars_.alpha;
78     Grad_(0,1) = pars_.a;
79     Grad_(0,2) = pars_.b;
80     Grad_(1,0) = Grad_(0,1);
81     Grad_(1,1) = pars_.beta;
82     Grad_(1,2) = pars_.c;
83     Grad_(2,0) = Grad_(0,2);
84     Grad_(2,1) = Grad_(1,2);
85     Grad_(2,2) = - (Grad_(0,0) + Grad_(1,1));
86     }
87     int storageLayout_ = info_->getSnapshotManager()->getStorageLayout();
88     if (storageLayout_ & DataStorage::dslParticlePot) doParticlePot = true;
89     initialized = true;
90     }
91    
92     void UniformGradient::applyPerturbation() {
93    
94     if (!initialized) initialize();
95    
96     SimInfo::MoleculeIterator i;
97     Molecule::AtomIterator j;
98     Molecule* mol;
99     Atom* atom;
100     AtomType* atype;
101     potVec longRangePotential(0.0);
102    
103     RealType C;
104     Vector3d D;
105     Mat3x3d Q;
106    
107     RealType U;
108     RealType fPot;
109     Vector3d t;
110     Vector3d f;
111    
112     Vector3d r;
113     Vector3d EF;
114    
115     bool isCharge;
116    
117     if (doUniformGradient) {
118    
119     U = 0.0;
120     fPot = 0.0;
121    
122     for (mol = info_->beginMolecule(i); mol != NULL;
123     mol = info_->nextMolecule(i)) {
124    
125     for (atom = mol->beginAtom(j); atom != NULL;
126     atom = mol->nextAtom(j)) {
127    
128     isCharge = false;
129     C = 0.0;
130    
131     atype = atom->getAtomType();
132    
133     r = atom->getPos();
134     EF = Grad_ * r;
135    
136     if (atype->isElectrostatic()) {
137     atom->addElectricField(EF * PhysicalConstants::chargeFieldConvert);
138     }
139    
140     FixedChargeAdapter fca = FixedChargeAdapter(atype);
141     if ( fca.isFixedCharge() ) {
142     isCharge = true;
143     C = fca.getCharge();
144     }
145    
146     FluctuatingChargeAdapter fqa = FluctuatingChargeAdapter(atype);
147     if ( fqa.isFluctuatingCharge() ) {
148     isCharge = true;
149     C += atom->getFlucQPos();
150     atom->addFlucQFrc( dot(r, EF)
151     * PhysicalConstants::chargeFieldConvert );
152     }
153    
154     if (isCharge) {
155     f = EF * C * PhysicalConstants::chargeFieldConvert;
156     atom->addFrc(f);
157    
158     U = -dot(r, f);
159     if (doParticlePot) {
160     atom->addParticlePot(U);
161     }
162     fPot += U;
163     }
164    
165     MultipoleAdapter ma = MultipoleAdapter(atype);
166     if (ma.isDipole() ) {
167     D = atom->getDipole() * PhysicalConstants::dipoleFieldConvert;
168    
169     f = D * Grad_;
170     atom->addFrc(f);
171    
172     t = cross(D, EF);
173     atom->addTrq(t);
174    
175     U = -dot(D, EF);
176     if (doParticlePot) {
177     atom->addParticlePot(U);
178     }
179     fPot += U;
180     }
181    
182     if (ma.isQuadrupole() ) {
183     Q = atom->getQuadrupole() * PhysicalConstants::dipoleFieldConvert;
184    
185     t = 2.0 * mCross(Q, Grad_);
186     atom->addTrq(t);
187    
188     U = -doubleDot(Q, Grad_);
189     if (doParticlePot) {
190     atom->addParticlePot(U);
191     }
192     fPot += U;
193     }
194     }
195     }
196    
197     #ifdef IS_MPI
198     MPI_Allreduce(MPI_IN_PLACE, &fPot, 1, MPI_REALTYPE,
199     MPI_SUM, MPI_COMM_WORLD);
200     #endif
201    
202     Snapshot* snap = info_->getSnapshotManager()->getCurrentSnapshot();
203     longRangePotential = snap->getLongRangePotentials();
204     longRangePotential[ELECTROSTATIC_FAMILY] += fPot;
205     snap->setLongRangePotential(longRangePotential);
206     }
207     }
208     }

Properties

Name Value
svn:executable *