NonbondedForce.cpp 10.8 KB
Newer Older
1
2
3
4
5
6
7
8
/* -------------------------------------------------------------------------- *
 *                                   OpenMM                                   *
 * -------------------------------------------------------------------------- *
 * This is part of the OpenMM molecular simulation toolkit originating from   *
 * Simbios, the NIH National Center for Physics-Based Simulation of           *
 * Biological Structures at Stanford, funded under the NIH Roadmap for        *
 * Medical Research, grant U54 GM072970. See https://simtk.org.               *
 *                                                                            *
9
 * Portions copyright (c) 2008-2016 Stanford University and the Authors.      *
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
 * Authors: Peter Eastman                                                     *
 * Contributors:                                                              *
 *                                                                            *
 * Permission is hereby granted, free of charge, to any person obtaining a    *
 * copy of this software and associated documentation files (the "Software"), *
 * to deal in the Software without restriction, including without limitation  *
 * the rights to use, copy, modify, merge, publish, distribute, sublicense,   *
 * and/or sell copies of the Software, and to permit persons to whom the      *
 * Software is furnished to do so, subject to the following conditions:       *
 *                                                                            *
 * The above copyright notice and this permission notice shall be included in *
 * all copies or substantial portions of the Software.                        *
 *                                                                            *
 * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR *
 * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,   *
 * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL    *
 * THE AUTHORS, CONTRIBUTORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM,    *
 * DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR      *
 * OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE  *
 * USE OR OTHER DEALINGS IN THE SOFTWARE.                                     *
 * -------------------------------------------------------------------------- */

32
33
34
#include "openmm/Force.h"
#include "openmm/OpenMMException.h"
#include "openmm/NonbondedForce.h"
35
#include "openmm/internal/AssertionUtilities.h"
36
#include "openmm/internal/NonbondedForceImpl.h"
37
#include <cmath>
38
39
#include <map>
#include <sstream>
40
#include <utility>
41
42

using namespace OpenMM;
43
using std::map;
44
45
using std::pair;
using std::set;
46
47
using std::string;
using std::stringstream;
48
using std::vector;
49

50
NonbondedForce::NonbondedForce() : nonbondedMethod(NoCutoff), cutoffDistance(1.0), switchingDistance(-1.0), rfDielectric(78.3),
51
        ewaldErrorTol(5e-4), alpha(0.0), useSwitchingFunction(false), useDispersionCorrection(true), recipForceGroup(-1), nx(0), ny(0), nz(0) {
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
}

NonbondedForce::NonbondedMethod NonbondedForce::getNonbondedMethod() const {
    return nonbondedMethod;
}

void NonbondedForce::setNonbondedMethod(NonbondedMethod method) {
    nonbondedMethod = method;
}

double NonbondedForce::getCutoffDistance() const {
    return cutoffDistance;
}

void NonbondedForce::setCutoffDistance(double distance) {
    cutoffDistance = distance;
}

70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
bool NonbondedForce::getUseSwitchingFunction() const {
    return useSwitchingFunction;
}

void NonbondedForce::setUseSwitchingFunction(bool use) {
    useSwitchingFunction = use;
}

double NonbondedForce::getSwitchingDistance() const {
    return switchingDistance;
}

void NonbondedForce::setSwitchingDistance(double distance) {
    switchingDistance = distance;
}

86
87
88
89
90
91
92
93
double NonbondedForce::getReactionFieldDielectric() const {
    return rfDielectric;
}

void NonbondedForce::setReactionFieldDielectric(double dielectric) {
    rfDielectric = dielectric;
}

94
95
96
97
double NonbondedForce::getEwaldErrorTolerance() const {
    return ewaldErrorTol;
}

98
void NonbondedForce::setEwaldErrorTolerance(double tol) {
99
100
101
    ewaldErrorTol = tol;
}

102
103
104
105
106
107
108
void NonbondedForce::getPMEParameters(double& alpha, int& nx, int& ny, int& nz) const {
    alpha = this->alpha;
    nx = this->nx;
    ny = this->ny;
    nz = this->nz;
}

109
110
111
112
113
114
115
void NonbondedForce::getLJPMEParameters(double& dalpha, int& dnx, int& dny, int& dnz) const {
    dalpha = this->dalpha;
    dnx = this->dnx;
    dny = this->dny;
    dnz = this->dnz;
}

116
117
118
119
120
121
122
void NonbondedForce::setPMEParameters(double alpha, int nx, int ny, int nz) {
    this->alpha = alpha;
    this->nx = nx;
    this->ny = ny;
    this->nz = nz;
}

123
124
125
126
127
128
129
void NonbondedForce::setLJPMEParameters(double dalpha, int dnx, int dny, int dnz) {
    this->dalpha = dalpha;
    this->dnx = dnx;
    this->dny = dny;
    this->dnz = dnz;
}

130
131
132
133
void NonbondedForce::getPMEParametersInContext(const Context& context, double& alpha, int& nx, int& ny, int& nz) const {
    dynamic_cast<const NonbondedForceImpl&>(getImplInContext(context)).getPMEParameters(alpha, nx, ny, nz);
}

134
135
136
137
void NonbondedForce::getLJPMEParametersInContext(const Context& context, double& dalpha, int& dnx, int& dny, int& dnz) const {
    dynamic_cast<const NonbondedForceImpl&>(getImplInContext(context)).getLJPMEParameters(dalpha, dnx, dny, dnz);
}

138
int NonbondedForce::addParticle(double charge, double sigma, double epsilon) {
139
    particles.push_back(ParticleInfo(charge, sigma, epsilon));
140
    return particles.size()-1;
141
142
}

143
void NonbondedForce::getParticleParameters(int index, double& charge, double& sigma, double& epsilon) const {
144
    ASSERT_VALID_INDEX(index, particles);
Peter Eastman's avatar
Peter Eastman committed
145
    charge = particles[index].charge;
146
147
    sigma = particles[index].sigma;
    epsilon = particles[index].epsilon;
148
149
}

150
void NonbondedForce::setParticleParameters(int index, double charge, double sigma, double epsilon) {
151
    ASSERT_VALID_INDEX(index, particles);
Peter Eastman's avatar
Peter Eastman committed
152
    particles[index].charge = charge;
153
154
    particles[index].sigma = sigma;
    particles[index].epsilon = epsilon;
155
156
}

157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
int NonbondedForce::addException(int particle1, int particle2, double chargeProd, double sigma, double epsilon, bool replace) {
    map<pair<int, int>, int>::iterator iter = exceptionMap.find(pair<int, int>(particle1, particle2));
    int newIndex;
    if (iter == exceptionMap.end())
        iter = exceptionMap.find(pair<int, int>(particle2, particle1));
    if (iter != exceptionMap.end()) {
        if (!replace) {
            stringstream msg;
            msg << "NonbondedForce: There is already an exception for particles ";
            msg << particle1;
            msg << " and ";
            msg << particle2;
            throw OpenMMException(msg.str());
        }
        exceptions[iter->second] = ExceptionInfo(particle1, particle2, chargeProd, sigma, epsilon);
        newIndex = iter->second;
        exceptionMap.erase(iter->first);
    }
    else {
        exceptions.push_back(ExceptionInfo(particle1, particle2, chargeProd, sigma, epsilon));
        newIndex = exceptions.size()-1;
    }
    exceptionMap[pair<int, int>(particle1, particle2)] = newIndex;
    return newIndex;
181
182
}
void NonbondedForce::getExceptionParameters(int index, int& particle1, int& particle2, double& chargeProd, double& sigma, double& epsilon) const {
183
    ASSERT_VALID_INDEX(index, exceptions);
184
185
186
187
188
    particle1 = exceptions[index].particle1;
    particle2 = exceptions[index].particle2;
    chargeProd = exceptions[index].chargeProd;
    sigma = exceptions[index].sigma;
    epsilon = exceptions[index].epsilon;
189
190
}

191
void NonbondedForce::setExceptionParameters(int index, int particle1, int particle2, double chargeProd, double sigma, double epsilon) {
192
    ASSERT_VALID_INDEX(index, exceptions);
193
194
195
196
197
    exceptions[index].particle1 = particle1;
    exceptions[index].particle2 = particle2;
    exceptions[index].chargeProd = chargeProd;
    exceptions[index].sigma = sigma;
    exceptions[index].epsilon = epsilon;
198
199
}

200
ForceImpl* NonbondedForce::createImpl() const {
201
202
    return new NonbondedForceImpl(*this);
}
203
204

void NonbondedForce::createExceptionsFromBonds(const vector<pair<int, int> >& bonds, double coulomb14Scale, double lj14Scale) {
205
206
207
    for (int i = 0; i < (int) bonds.size(); ++i)
        if (bonds[i].first < 0 || bonds[i].second < 0 || bonds[i].first >= particles.size() || bonds[i].second >= particles.size())
            throw OpenMMException("createExceptionsFromBonds: Illegal particle index in list of bonds");
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253

    // Find particles separated by 1, 2, or 3 bonds.

    vector<set<int> > exclusions(particles.size());
    vector<set<int> > bonded12(exclusions.size());
    for (int i = 0; i < (int) bonds.size(); ++i) {
        bonded12[bonds[i].first].insert(bonds[i].second);
        bonded12[bonds[i].second].insert(bonds[i].first);
    }
    for (int i = 0; i < (int) exclusions.size(); ++i)
        addExclusionsToSet(bonded12, exclusions[i], i, i, 2);

    // Find particles separated by 1 or 2 bonds and create the exceptions.

    for (int i = 0; i < (int) exclusions.size(); ++i) {
        set<int> bonded13;
        addExclusionsToSet(bonded12, bonded13, i, i, 1);
        for (set<int>::const_iterator iter = exclusions[i].begin(); iter != exclusions[i].end(); ++iter)
            if (*iter < i) {
                if (bonded13.find(*iter) == bonded13.end()) {
                    // This is a 1-4 interaction.

                    const ParticleInfo& particle1 = particles[*iter];
                    const ParticleInfo& particle2 = particles[i];
                    const double chargeProd = coulomb14Scale*particle1.charge*particle2.charge;
                    const double sigma = 0.5*(particle1.sigma+particle2.sigma);
                    const double epsilon = lj14Scale*std::sqrt(particle1.epsilon*particle2.epsilon);
                    addException(*iter, i, chargeProd, sigma, epsilon);
                }
                else {
                    // This interaction should be completely excluded.

                    addException(*iter, i, 0.0, 1.0, 0.0);
                }
            }
    }
}

void NonbondedForce::addExclusionsToSet(const vector<set<int> >& bonded12, set<int>& exclusions, int baseParticle, int fromParticle, int currentLevel) const {
    for (set<int>::const_iterator iter = bonded12[fromParticle].begin(); iter != bonded12[fromParticle].end(); ++iter) {
        if (*iter != baseParticle)
            exclusions.insert(*iter);
        if (currentLevel > 0)
            addExclusionsToSet(bonded12, exclusions, baseParticle, *iter, currentLevel-1);
    }
}
254
255
256
257
258
259

int NonbondedForce::getReciprocalSpaceForceGroup() const {
    return recipForceGroup;
}

void NonbondedForce::setReciprocalSpaceForceGroup(int group) {
260
261
    if (group < -1 || group > 31)
        throw OpenMMException("Force group must be between -1 and 31");
262
263
    recipForceGroup = group;
}
264
265
266
267

void NonbondedForce::updateParametersInContext(Context& context) {
    dynamic_cast<NonbondedForceImpl&>(getImplInContext(context)).updateParametersInContext(getContextImpl(context));
}