CpuCustomNonbondedForce.cpp 13.5 KB
Newer Older
1

2
/* Portions copyright (c) 2009-2016 Stanford University and Simbios.
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
 * Contributors: Peter Eastman
 *
 * 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.
 */

#include <string.h>
#include <sstream>

#include "SimTKOpenMMUtilities.h"
#include "ReferenceForce.h"
#include "CpuCustomNonbondedForce.h"
31
#include "openmm/internal/gmx_atomic.h"
32
33
34
35

using namespace OpenMM;
using namespace std;

36
37
38
39
40
41
42
43
44
class CpuCustomNonbondedForce::ComputeForceTask : public ThreadPool::Task {
public:
    ComputeForceTask(CpuCustomNonbondedForce& owner) : owner(owner) {
    }
    void execute(ThreadPool& threads, int threadIndex) {
        owner.threadComputeForce(threads, threadIndex);
    }
    CpuCustomNonbondedForce& owner;
};
45

46
47
CpuCustomNonbondedForce::ThreadData::ThreadData(const Lepton::CompiledExpression& energyExpression, const Lepton::CompiledExpression& forceExpression, const vector<string>& parameterNames) :
            energyExpression(energyExpression), forceExpression(forceExpression) {
48
49
    energyR = ReferenceForce::getVariablePointer(this->energyExpression, "r");
    forceR = ReferenceForce::getVariablePointer(this->forceExpression, "r");
50
    for (int i = 0; i < (int) parameterNames.size(); i++) {
51
52
        for (int j = 1; j < 3; j++) {
            stringstream name;
53
            name << parameterNames[i] << j;
54
55
56
57
58
59
            energyParticleParams.push_back(ReferenceForce::getVariablePointer(this->energyExpression, name.str()));
            forceParticleParams.push_back(ReferenceForce::getVariablePointer(this->forceExpression, name.str()));
        }
    }
}

60
CpuCustomNonbondedForce::CpuCustomNonbondedForce(const Lepton::CompiledExpression& energyExpression,
61
62
            const Lepton::CompiledExpression& forceExpression, const vector<string>& parameterNames, const vector<set<int> >& exclusions,ThreadPool& threads) :
            cutoff(false), useSwitch(false), periodic(false), paramNames(parameterNames), exclusions(exclusions), threads(threads) {
63
64
    for (int i = 0; i < threads.getNumThreads(); i++)
        threadData.push_back(new ThreadData(energyExpression, forceExpression, parameterNames));
65
66
}

67
68
69
70
CpuCustomNonbondedForce::~CpuCustomNonbondedForce() {
    for (int i = 0; i < (int) threadData.size(); i++)
        delete threadData[i];
}
71

72
void CpuCustomNonbondedForce::setUseCutoff(RealOpenMM distance, const CpuNeighborList& neighbors) {
73
74
75
76
77
78
    cutoff = true;
    cutoffDistance = distance;
    neighborList = &neighbors;
  }

void CpuCustomNonbondedForce::setInteractionGroups(const vector<pair<set<int>, set<int> > >& groups) {
79
80
81
82
83
84
85
86
87
88
89
90
91
    for (int group = 0; group < (int) groups.size(); group++) {
        const set<int>& set1 = groups[group].first;
        const set<int>& set2 = groups[group].second;
        for (set<int>::const_iterator atom1 = set1.begin(); atom1 != set1.end(); ++atom1) {
            for (set<int>::const_iterator atom2 = set2.begin(); atom2 != set2.end(); ++atom2) {
                if (*atom1 == *atom2 || exclusions[*atom1].find(*atom2) != exclusions[*atom1].end())
                    continue; // This is an excluded interaction.
                if (*atom1 > *atom2 && set1.find(*atom2) != set1.end() && set2.find(*atom1) != set2.end())
                    continue; // Both atoms are in both sets, so skip duplicate interactions.
                groupInteractions.push_back(make_pair(*atom1, *atom2));
            }
        }
    }
92
93
}

94
void CpuCustomNonbondedForce::setUseSwitchingFunction(RealOpenMM distance) {
95
96
97
98
    useSwitch = true;
    switchingDistance = distance;
}

99
void CpuCustomNonbondedForce::setPeriodic(RealVec* periodicBoxVectors) {
100
    assert(cutoff);
101
102
103
    assert(periodicBoxVectors[0][0] >= 2.0*cutoffDistance);
    assert(periodicBoxVectors[1][1] >= 2.0*cutoffDistance);
    assert(periodicBoxVectors[2][2] >= 2.0*cutoffDistance);
104
    periodic = true;
105
106
107
108
109
110
111
112
113
114
115
116
117
    this->periodicBoxVectors[0] = periodicBoxVectors[0];
    this->periodicBoxVectors[1] = periodicBoxVectors[1];
    this->periodicBoxVectors[2] = periodicBoxVectors[2];
    recipBoxSize[0] = (float) (1.0/periodicBoxVectors[0][0]);
    recipBoxSize[1] = (float) (1.0/periodicBoxVectors[1][1]);
    recipBoxSize[2] = (float) (1.0/periodicBoxVectors[2][2]);
    periodicBoxVec4.resize(3);
    periodicBoxVec4[0] = fvec4(periodicBoxVectors[0][0], periodicBoxVectors[0][1], periodicBoxVectors[0][2], 0);
    periodicBoxVec4[1] = fvec4(periodicBoxVectors[1][0], periodicBoxVectors[1][1], periodicBoxVectors[1][2], 0);
    periodicBoxVec4[2] = fvec4(periodicBoxVectors[2][0], periodicBoxVectors[2][1], periodicBoxVectors[2][2], 0);
    triclinic = (periodicBoxVectors[0][1] != 0.0 || periodicBoxVectors[0][2] != 0.0 ||
                 periodicBoxVectors[1][0] != 0.0 || periodicBoxVectors[1][2] != 0.0 ||
                 periodicBoxVectors[2][0] != 0.0 || periodicBoxVectors[2][1] != 0.0);
118
}
119
120


121
void CpuCustomNonbondedForce::calculatePairIxn(int numberOfAtoms, float* posq, vector<RealVec>& atomCoordinates, RealOpenMM** atomParameters,
122
                                             RealOpenMM* fixedParameters, const map<string, double>& globalParameters,
123
                                             vector<AlignedArray<float> >& threadForce, bool includeForce, bool includeEnergy, double& totalEnergy) {
124
125
126
127
128
129
    // Record the parameters for the threads.
    
    this->numberOfAtoms = numberOfAtoms;
    this->posq = posq;
    this->atomCoordinates = &atomCoordinates[0];
    this->atomParameters = atomParameters;
130
    this->globalParameters = &globalParameters;
131
    this->threadForce = &threadForce;
132
133
    this->includeForce = includeForce;
    this->includeEnergy = includeEnergy;
134
135
136
137
138
    threadEnergy.resize(threads.getNumThreads());
    gmx_atomic_t counter;
    gmx_atomic_set(&counter, 0);
    this->atomicCounter = &counter;
    
139
    // Signal the threads to start running and wait for them to finish.
140
    
141
142
143
    ComputeForceTask task(*this);
    threads.execute(task);
    threads.waitForThreads();
144
    
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
    // Combine the energies from all the threads.
    
    if (includeEnergy) {
        int numThreads = threads.getNumThreads();
        for (int i = 0; i < numThreads; i++)
            totalEnergy += threadEnergy[i];
    }
}

void CpuCustomNonbondedForce::threadComputeForce(ThreadPool& threads, int threadIndex) {
    // Compute this thread's subset of interactions.

    int numThreads = threads.getNumThreads();
    threadEnergy[threadIndex] = 0;
    double& energy = threadEnergy[threadIndex];
    float* forces = &(*threadForce)[threadIndex][0];
    ThreadData& data = *threadData[threadIndex];
    for (map<string, double>::const_iterator iter = globalParameters->begin(); iter != globalParameters->end(); ++iter) {
        ReferenceForce::setVariable(ReferenceForce::getVariablePointer(data.energyExpression, iter->first), iter->second);
        ReferenceForce::setVariable(ReferenceForce::getVariablePointer(data.forceExpression, iter->first), iter->second);
165
    }
166
167
    fvec4 boxSize(periodicBoxVectors[0][0], periodicBoxVectors[1][1], periodicBoxVectors[2][2], 0);
    fvec4 invBoxSize(recipBoxSize[0], recipBoxSize[1], recipBoxSize[2], 0);
168
    if (groupInteractions.size() > 0) {
169
170
        // The user has specified interaction groups, so compute only the requested interactions.
        
171
172
173
174
175
176
177
178
179
180
181
        while (true) {
            int i = gmx_atomic_fetch_add(reinterpret_cast<gmx_atomic_t*>(atomicCounter), 1);
            if (i >= groupInteractions.size())
                break;
            int atom1 = groupInteractions[i].first;
            int atom2 = groupInteractions[i].second;
            for (int j = 0; j < (int) paramNames.size(); j++) {
                ReferenceForce::setVariable(data.energyParticleParams[j*2], atomParameters[atom1][j]);
                ReferenceForce::setVariable(data.energyParticleParams[j*2+1], atomParameters[atom2][j]);
                ReferenceForce::setVariable(data.forceParticleParams[j*2], atomParameters[atom1][j]);
                ReferenceForce::setVariable(data.forceParticleParams[j*2+1], atomParameters[atom2][j]);
182
            }
183
            calculateOneIxn(atom1, atom2, data, forces, energy, boxSize, invBoxSize);
184
185
186
187
        }
    }
    else if (cutoff) {
        // We are using a cutoff, so get the interactions from the neighbor list.
188
189
190
191
192

        while (true) {
            int blockIndex = gmx_atomic_fetch_add(reinterpret_cast<gmx_atomic_t*>(atomicCounter), 1);
            if (blockIndex >= neighborList->getNumBlocks())
                break;
193
194
            const int blockSize = neighborList->getBlockSize();
            const int* blockAtom = &neighborList->getSortedAtoms()[blockSize*blockIndex];
195
196
197
198
199
            const vector<int>& neighbors = neighborList->getBlockNeighbors(blockIndex);
            const vector<char>& exclusions = neighborList->getBlockExclusions(blockIndex);
            for (int i = 0; i < (int) neighbors.size(); i++) {
                int first = neighbors[i];
                for (int j = 0; j < (int) paramNames.size(); j++) {
200
201
                    ReferenceForce::setVariable(data.energyParticleParams[j*2], atomParameters[first][j]);
                    ReferenceForce::setVariable(data.forceParticleParams[j*2], atomParameters[first][j]);
202
                }
203
                for (int k = 0; k < blockSize; k++) {
204
205
206
                    if ((exclusions[i] & (1<<k)) == 0) {
                        int second = blockAtom[k];
                        for (int j = 0; j < (int) paramNames.size(); j++) {
207
208
                            ReferenceForce::setVariable(data.energyParticleParams[j*2+1], atomParameters[second][j]);
                            ReferenceForce::setVariable(data.forceParticleParams[j*2+1], atomParameters[second][j]);
209
                        }
210
                        calculateOneIxn(first, second, data, forces, energy, boxSize, invBoxSize);
211
212
213
214
215
216
217
218
                    }
                }
            }
        }
    }
    else {
        // Every particle interacts with every other one.
        
219
220
221
222
        while (true) {
            int ii = gmx_atomic_fetch_add(reinterpret_cast<gmx_atomic_t*>(atomicCounter), 1);
            if (ii >= numberOfAtoms)
                break;
223
224
225
            for (int jj = ii+1; jj < numberOfAtoms; jj++) {
                if (exclusions[jj].find(ii) == exclusions[jj].end()) {
                    for (int j = 0; j < (int) paramNames.size(); j++) {
226
227
228
229
                        ReferenceForce::setVariable(data.energyParticleParams[j*2], atomParameters[ii][j]);
                        ReferenceForce::setVariable(data.energyParticleParams[j*2+1], atomParameters[jj][j]);
                        ReferenceForce::setVariable(data.forceParticleParams[j*2], atomParameters[ii][j]);
                        ReferenceForce::setVariable(data.forceParticleParams[j*2+1], atomParameters[jj][j]);
230
                    }
231
                    calculateOneIxn(ii, jj, data, forces, energy, boxSize, invBoxSize);
232
233
234
235
236
237
                }
            }
        }
    }
}

238
239
240
void CpuCustomNonbondedForce::calculateOneIxn(int ii, int jj, ThreadData& data, 
        float* forces, double& totalEnergy, const fvec4& boxSize, const fvec4& invBoxSize) {
    // Get deltaR, R2, and R between 2 atoms
241
242
243
244
245
246
247
248
249
250
251
252

    fvec4 deltaR;
    fvec4 posI(posq+4*ii);
    fvec4 posJ(posq+4*jj);
    float r2;
    getDeltaR(posI, posJ, deltaR, r2, boxSize, invBoxSize);
    if (cutoff && r2 >= cutoffDistance*cutoffDistance)
        return;
    float r = sqrtf(r2);

    // accumulate forces

253
254
255
256
    ReferenceForce::setVariable(data.energyR, r);
    ReferenceForce::setVariable(data.forceR, r);
    double dEdR = (includeForce ? data.forceExpression.evaluate()/r : 0.0);
    double energy = (includeEnergy ? data.energyExpression.evaluate() : 0.0);
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
    if (useSwitch) {
        if (r > switchingDistance) {
            RealOpenMM t = (r-switchingDistance)/(cutoffDistance-switchingDistance);
            RealOpenMM switchValue = 1+t*t*t*(-10+t*(15-t*6));
            RealOpenMM switchDeriv = t*t*(-30+t*(60-t*30))/(cutoffDistance-switchingDistance);
            dEdR = switchValue*dEdR + energy*switchDeriv/r;
            energy *= switchValue;
        }
    }
    fvec4 result = deltaR*dEdR;
    (fvec4(forces+4*ii)+result).store(forces+4*ii);
    (fvec4(forces+4*jj)-result).store(forces+4*jj);

    // accumulate energies

272
    totalEnergy += energy;
273
274
275
276
277
}

void CpuCustomNonbondedForce::getDeltaR(const fvec4& posI, const fvec4& posJ, fvec4& deltaR, float& r2, const fvec4& boxSize, const fvec4& invBoxSize) const {
    deltaR = posJ-posI;
    if (periodic) {
278
279
280
281
282
283
284
285
286
        if (triclinic) {
            deltaR -= periodicBoxVec4[2]*floorf(deltaR[2]*recipBoxSize[2]+0.5f);
            deltaR -= periodicBoxVec4[1]*floorf(deltaR[1]*recipBoxSize[1]+0.5f);
            deltaR -= periodicBoxVec4[0]*floorf(deltaR[0]*recipBoxSize[0]+0.5f);
        }
        else {
            fvec4 base = round(deltaR*invBoxSize)*boxSize;
            deltaR = deltaR-base;
        }
287
288
289
    }
    r2 = dot3(deltaR, deltaR);
}