CpuKernels.h 31.1 KB
Newer Older
1
2
3
4
5
6
#ifndef OPENMM_CPUKERNELS_H_
#define OPENMM_CPUKERNELS_H_

/* -------------------------------------------------------------------------- *
 *                                   OpenMM                                   *
 * -------------------------------------------------------------------------- *
Evan Pretti's avatar
Evan Pretti committed
7
8
 * This is part of the OpenMM molecular simulation toolkit.                   *
 * See https://openmm.org/development.                                        *
9
 *                                                                            *
10
 * Portions copyright (c) 2013-2025 Stanford University and the Authors.      *
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
 * 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.                                     *
 * -------------------------------------------------------------------------- */

33
#include "CpuBondForce.h"
34
#include "CpuConstantPotentialForce.h"
35
#include "CpuCustomGBForce.h"
36
#include "CpuCustomManyParticleForce.h"
37
#include "CpuCustomNonbondedForce.h"
38
#include "CpuGayBerneForce.h"
39
#include "CpuGBSAOBCForce.h"
40
#include "CpuLangevinMiddleDynamics.h"
Evan Pretti's avatar
Evan Pretti committed
41
#include "CpuLCPOForce.h"
42
#include "CpuNeighborList.h"
43
#include "CpuNonbondedForce.h"
44
#include "CpuPlatform.h"
45
#include "ReferenceKernels.h"
46
47
#include "openmm/kernels.h"
#include "openmm/System.h"
48
#include "openmm/internal/CustomNonbondedForceImpl.h"
49
50
#include <array>
#include <tuple>
51
52
53

namespace OpenMM {

54
55
56
57
58
59
60
61
62
63
/**
 * This kernel is invoked at the beginning and end of force and energy computations.  It gives the
 * Platform a chance to clear buffers and do other initialization at the beginning, and to do any
 * necessary work at the end to determine the final results.
 */
class CpuCalcForcesAndEnergyKernel : public CalcForcesAndEnergyKernel {
public:
    CpuCalcForcesAndEnergyKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data, ContextImpl& context);
    /**
     * Initialize the kernel.
Evan Pretti's avatar
Evan Pretti committed
64
     *
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
     * @param system     the System this kernel will be applied to
     */
    void initialize(const System& system);
    /**
     * This is called at the beginning of each force/energy computation, before calcForcesAndEnergy() has been called on
     * any ForceImpl.
     *
     * @param context       the context in which to execute this kernel
     * @param includeForce  true if forces should be computed
     * @param includeEnergy true if potential energy should be computed
     * @param groups        a set of bit flags for which force groups to include
     */
    void beginComputation(ContextImpl& context, bool includeForce, bool includeEnergy, int groups);
    /**
     * This is called at the end of each force/energy computation, after calcForcesAndEnergy() has been called on
     * every ForceImpl.
     *
     * @param context       the context in which to execute this kernel
     * @param includeForce  true if forces should be computed
     * @param includeEnergy true if potential energy should be computed
     * @param groups        a set of bit flags for which force groups to include
86
87
     * @param valid         the method may set this to false to indicate the results are invalid and the force/energy
     *                      calculation should be repeated
88
89
90
91
     * @return the potential energy of the system.  This value is added to all values returned by ForceImpls'
     * calcForcesAndEnergy() methods.  That is, each force kernel may <i>either</i> return its contribution to the
     * energy directly, <i>or</i> add it to an internal buffer so that it will be included here.
     */
92
    double finishComputation(ContextImpl& context, bool includeForce, bool includeEnergy, int groups, bool& valid);
93
94
95
private:
    CpuPlatform::PlatformData& data;
    Kernel referenceKernel;
peastman's avatar
peastman committed
96
    std::vector<Vec3> lastPositions;
97
98
};

99
100
101
102
103
104
105
106
107
108
109
/**
 * This kernel provides methods for setting and retrieving various state data: time, positions,
 * velocities, and forces.
 */
class CpuUpdateStateDataKernel : public ReferenceUpdateStateDataKernel {
public:
    CpuUpdateStateDataKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data, ReferencePlatform::PlatformData& refdata) :
            ReferenceUpdateStateDataKernel(name, platform, refdata), data(data) {
    }
    /**
     * Create a checkpoint recording the current state of the Context.
Evan Pretti's avatar
Evan Pretti committed
110
     *
111
112
113
114
115
     * @param stream    an output stream the checkpoint data should be written to
     */
    void createCheckpoint(ContextImpl& context, std::ostream& stream);
    /**
     * Load a checkpoint that was written by createCheckpoint().
Evan Pretti's avatar
Evan Pretti committed
116
     *
117
118
119
120
121
122
123
     * @param stream    an input stream the checkpoint data should be read from
     */
    void loadCheckpoint(ContextImpl& context, std::istream& stream);
private:
    CpuPlatform::PlatformData& data;
};

peastman's avatar
peastman committed
124
125
126
127
128
129
/**
 * This kernel is invoked by HarmonicAngleForce to calculate the forces acting on the system and the energy of the system.
 */
class CpuCalcHarmonicAngleForceKernel : public CalcHarmonicAngleForceKernel {
public:
    CpuCalcHarmonicAngleForceKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data) :
130
            CalcHarmonicAngleForceKernel(name, platform), data(data), usePeriodic(false) {
peastman's avatar
peastman committed
131
132
133
    }
    /**
     * Initialize the kernel.
Evan Pretti's avatar
Evan Pretti committed
134
     *
peastman's avatar
peastman committed
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
     * @param system     the System this kernel will be applied to
     * @param force      the HarmonicAngleForce this kernel will be used for
     */
    void initialize(const System& system, const HarmonicAngleForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
151
152
     * @param firstAngle the index of the first bond whose parameters might have changed
     * @param lastAngle  the index of the last bond whose parameters might have changed
peastman's avatar
peastman committed
153
     */
154
    void copyParametersToContext(ContextImpl& context, const HarmonicAngleForce& force, int firstAngle, int lastAngle);
peastman's avatar
peastman committed
155
156
157
private:
    CpuPlatform::PlatformData& data;
    int numAngles;
158
159
    std::vector<std::vector<int> > angleIndexArray;
    std::vector<std::vector<double> > angleParamArray;
peastman's avatar
peastman committed
160
    CpuBondForce bondForce;
161
    bool usePeriodic;
peastman's avatar
peastman committed
162
163
};

164
165
166
167
168
169
/**
 * This kernel is invoked by PeriodicTorsionForce to calculate the forces acting on the system and the energy of the system.
 */
class CpuCalcPeriodicTorsionForceKernel : public CalcPeriodicTorsionForceKernel {
public:
    CpuCalcPeriodicTorsionForceKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data) :
170
            CalcPeriodicTorsionForceKernel(name, platform), data(data), usePeriodic(false) {
171
172
173
    }
    /**
     * Initialize the kernel.
Evan Pretti's avatar
Evan Pretti committed
174
     *
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
     * @param system     the System this kernel will be applied to
     * @param force      the PeriodicTorsionForce this kernel will be used for
     */
    void initialize(const System& system, const PeriodicTorsionForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
191
192
193
194
     * @param context      the context to copy parameters to
     * @param force        the PeriodicTorsionForce to copy the parameters from
     * @param firstTorsion the index of the first torsion whose parameters might have changed
     * @param lastTorsion  the index of the last torsion whose parameters might have changed
195
     */
196
    void copyParametersToContext(ContextImpl& context, const PeriodicTorsionForce& force, int firstTorsion, int lastTorsion);
197
198
199
private:
    CpuPlatform::PlatformData& data;
    int numTorsions;
200
201
    std::vector<std::vector<int> > torsionIndexArray;
    std::vector<std::vector<double> > torsionParamArray;
202
    CpuBondForce bondForce;
203
    bool usePeriodic;
204
205
206
207
208
209
210
211
};

/**
 * This kernel is invoked by RBTorsionForce to calculate the forces acting on the system and the energy of the system.
 */
class CpuCalcRBTorsionForceKernel : public CalcRBTorsionForceKernel {
public:
    CpuCalcRBTorsionForceKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data) :
212
            CalcRBTorsionForceKernel(name, platform), data(data), usePeriodic(false) {
213
214
215
    }
    /**
     * Initialize the kernel.
Evan Pretti's avatar
Evan Pretti committed
216
     *
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
     * @param system     the System this kernel will be applied to
     * @param force      the RBTorsionForce this kernel will be used for
     */
    void initialize(const System& system, const RBTorsionForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the RBTorsionForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const RBTorsionForce& force);
private:
    CpuPlatform::PlatformData& data;
    int numTorsions;
240
241
    std::vector<std::vector<int> > torsionIndexArray;
    std::vector<std::vector<double> > torsionParamArray;
242
    CpuBondForce bondForce;
243
    bool usePeriodic;
244
245
};

246
247
248
249
250
/**
 * This kernel is invoked by NonbondedForce to calculate the forces acting on the system.
 */
class CpuCalcNonbondedForceKernel : public CalcNonbondedForceKernel {
public:
251
    CpuCalcNonbondedForceKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data);
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
    ~CpuCalcNonbondedForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the NonbondedForce this kernel will be used for
     */
    void initialize(const System& system, const NonbondedForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @param includeDirect  true if direct space interactions should be included
     * @param includeReciprocal  true if reciprocal space interactions should be included
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy, bool includeDirect, bool includeReciprocal);
    /**
     * Copy changed parameters over to a context.
     *
274
275
276
277
278
279
     * @param context        the context to copy parameters to
     * @param force          the NonbondedForce to copy the parameters from
     * @param firstParticle  the index of the first particle whose parameters might have changed
     * @param lastParticle   the index of the last particle whose parameters might have changed
     * @param firstException the index of the first exception whose parameters might have changed
     * @param lastException  the index of the last exception whose parameters might have changed
280
     */
281
    void copyParametersToContext(ContextImpl& context, const NonbondedForce& force, int firstParticle, int lastParticle, int firstException, int lastException);
282
283
    /**
     * Get the parameters being used for PME.
284
     *
285
286
287
288
289
290
     * @param alpha   the separation parameter
     * @param nx      the number of grid points along the X axis
     * @param ny      the number of grid points along the Y axis
     * @param nz      the number of grid points along the Z axis
     */
    void getPMEParameters(double& alpha, int& nx, int& ny, int& nz) const;
291
292
293
    /**
     * Get the parameters being used for the dispersion term in LJPME.
     *
294
295
296
297
     * @param alpha   the separation parameter
     * @param nx      the number of grid points along the X axis
     * @param ny      the number of grid points along the Y axis
     * @param nz      the number of grid points along the Z axis
298
     */
299
    void getLJPMEParameters(double& alpha, int& nx, int& ny, int& nz) const;
300
private:
301
    class PmeIO;
302
    void computeParameters(ContextImpl& context, bool offsetsOnly);
303
    CpuPlatform::PlatformData& data;
304
    int numParticles, num14, chargePosqIndex, ljPosqIndex;
305
306
    std::vector<std::vector<int> > bonded14IndexArray;
    std::vector<std::vector<double> > bonded14ParamArray;
307
    std::map<int, int> nb14Index;
308
    double nonbondedCutoff, switchingDistance, rfDielectric, ewaldAlpha, ewaldDispersionAlpha, ewaldSelfEnergy, dispersionCoefficient, totalCharge;
309
    int kmax[3], gridSize[3], dispersionGridSize[3];
310
    bool useSwitchingFunction, exceptionsArePeriodic, useOptimizedPme, hasInitializedPme, hasInitializedDispersionPme, hasParticleOffsets, hasExceptionOffsets;
311
    std::vector<std::set<int> > exclusions;
peastman's avatar
peastman committed
312
    std::vector<std::pair<float, float> > particleParams;
313
    std::vector<float> C6params;
314
    std::vector<float> charges;
315
316
317
    std::vector<std::array<double, 3> > baseParticleParams, baseExceptionParams;
    std::vector<std::vector<std::tuple<double, double, double, int> > > particleParamOffsets, exceptionParamOffsets;
    std::vector<std::string> paramNames;
318
    std::vector<double> paramValues;
319
    NonbondedMethod nonbondedMethod;
320
    CpuNonbondedForce* nonbonded;
321
    Kernel optimizedPme, optimizedDispersionPme;
peastman's avatar
peastman committed
322
    CpuBondForce bondForce;
323
324
};

325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
/**
 * This kernel is invoked by ConstantPotentialForce to calculate the forces acting on the system.
 */
class CpuCalcConstantPotentialForceKernel : public CalcConstantPotentialForceKernel {
public:
    CpuCalcConstantPotentialForceKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data);
    ~CpuCalcConstantPotentialForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the ConstantPotentialForce this kernel will be used for
     */
    void initialize(const System& system, const ConstantPotentialForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
     * @param context        the context to copy parameters to
     * @param force          the ConstantPotentialForce to copy the parameters from
     * @param firstParticle  the index of the first particle whose parameters might have changed
     * @param lastParticle   the index of the last particle whose parameters might have changed
     * @param firstException the index of the first exception whose parameters might have changed
     * @param lastException  the index of the last exception whose parameters might have changed
     * @param firstElectrode the index of the first electrode whose parameters might have changed
     * @param lastElectrode  the index of the last electrode whose parameters might have changed
     */
    void copyParametersToContext(ContextImpl& context, const ConstantPotentialForce& force, int firstParticle, int lastParticle, int firstException, int lastException, int firstElectrode, int lastElectrode);
    /**
     * Get the parameters being used for PME.
     *
     * @param alpha   the separation parameter
     * @param nx      the number of grid points along the X axis
     * @param ny      the number of grid points along the Y axis
     * @param nz      the number of grid points along the Z axis
     */
    void getPMEParameters(double& alpha, int& nx, int& ny, int& nz) const;
    /**
     * Get the charges on all particles.
     *
     * @param context       the context to copy parameters to
     * @param[out] charges  a vector to populate with particle charges
     */
    void getCharges(ContextImpl& context, std::vector<double>& charges);
private:
    void checkBoxSize(const Vec3* boxVectors);
    void ensurePmeInitialized(ContextImpl& context);
private:
    CpuPlatform::PlatformData& data;
    int numParticles, num14, numElectrodeParticles, chargePosqIndex;
    std::vector<double> setCharges;
    std::vector<float> charges;
    std::vector<std::vector<double> > bonded14ParamArray;
    std::vector<std::vector<int> > bonded14IndexArray;
    std::map<int, int> nb14Index;
    std::vector<std::set<int> > exclusions;
    std::vector<int> sysToElec, elecToSys, sysElec, elecElec;
    std::vector<std::array<double, 3> > electrodeParams;
    double nonbondedCutoff, ewaldAlpha, cgErrorTol, chargeTarget;
    int gridSize[3];
    bool exceptionsArePeriodic, hasInitializedPme, useChargeConstraint;
    Vec3 externalField;
    CpuConstantPotentialForce* constantPotential;
    CpuConstantPotentialSolver* solver;
    CpuBondForce bondForce;
    Kernel pmeKernel;
};


402
403
404
405
406
/**
 * This kernel is invoked by CustomNonbondedForce to calculate the forces acting on the system.
 */
class CpuCalcCustomNonbondedForceKernel : public CalcCustomNonbondedForceKernel {
public:
407
    CpuCalcCustomNonbondedForceKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data);
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
    ~CpuCalcCustomNonbondedForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomNonbondedForce this kernel will be used for
     */
    void initialize(const System& system, const CustomNonbondedForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomNonbondedForce to copy the parameters from
430
431
     * @param firstParticle  the index of the first particle whose parameters might have changed
     * @param lastParticle   the index of the last particle whose parameters might have changed
432
     */
433
    void copyParametersToContext(ContextImpl& context, const CustomNonbondedForce& force, int firstParticle, int lastParticle);
434
435
private:
    void createInteraction(const CustomNonbondedForce& force);
436
437
    CpuPlatform::PlatformData& data;
    int numParticles;
438
    std::vector<std::vector<double> > particleParamArray;
439
440
441
    double nonbondedCutoff, switchingDistance, periodicBoxSize[3], longRangeCoefficient;
    bool useSwitchingFunction, hasInitializedLongRangeCorrection;
    CustomNonbondedForce* forceCopy;
442
    CustomNonbondedForceImpl::LongRangeCorrectionData longRangeCorrectionData;
443
444
    std::map<std::string, double> globalParamValues;
    std::vector<std::set<int> > exclusions;
445
    std::vector<std::string> parameterNames, globalParameterNames, computedValueNames, energyParamDerivNames;
446
    std::vector<std::pair<std::set<int>, std::set<int> > > interactionGroups;
447
    std::vector<double> longRangeCoefficientDerivs;
448
    std::map<std::string, int> tabulatedFunctionUpdateCount;
449
    NonbondedMethod nonbondedMethod;
450
    CpuCustomNonbondedForce* nonbonded;
451
452
};

453
454
455
456
457
458
459
460
461
462
463
/**
 * This kernel is invoked by GBSAOBCForce to calculate the forces acting on the system.
 */
class CpuCalcGBSAOBCForceKernel : public CalcGBSAOBCForceKernel {
public:
    CpuCalcGBSAOBCForceKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data) : CalcGBSAOBCForceKernel(name, platform),
            data(data) {
    }
    ~CpuCalcGBSAOBCForceKernel();
    /**
     * Initialize the kernel.
Evan Pretti's avatar
Evan Pretti committed
464
     *
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
     * @param system     the System this kernel will be applied to
     * @param force      the GBSAOBCForce this kernel will be used for
     */
    void initialize(const System& system, const GBSAOBCForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the GBSAOBCForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const GBSAOBCForce& force);
private:
    CpuPlatform::PlatformData& data;
487
    int posqIndex;
488
    std::vector<std::pair<float, float> > particleParams;
489
    std::vector<float> charges;
490
491
492
    CpuGBSAOBCForce obc;
};

493
494
495
496
497
498
/**
 * This kernel is invoked by CustomGBForce to calculate the forces acting on the system.
 */
class CpuCalcCustomGBForceKernel : public CalcCustomGBForceKernel {
public:
    CpuCalcCustomGBForceKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data) :
Peter Eastman's avatar
Peter Eastman committed
499
            CalcCustomGBForceKernel(name, platform), data(data), ixn(NULL), neighborList(NULL) {
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
    }
    ~CpuCalcCustomGBForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomGBForce this kernel will be used for
     */
    void initialize(const System& system, const CustomGBForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomGBForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CustomGBForce& force);
private:
526
    void createInteraction(const CustomGBForce& force);
527
528
529
    CpuPlatform::PlatformData& data;
    int numParticles;
    bool isPeriodic;
530
    std::vector<std::vector<double> > particleParamArray;
peastman's avatar
peastman committed
531
    double nonbondedCutoff;
532
    CpuCustomGBForce* ixn;
Peter Eastman's avatar
Peter Eastman committed
533
    CpuNeighborList* neighborList;
534
    std::vector<std::set<int> > exclusions;
535
    std::vector<std::string> particleParameterNames, globalParameterNames, energyParamDerivNames, valueNames;
536
537
    std::vector<OpenMM::CustomGBForce::ComputationType> valueTypes;
    std::vector<OpenMM::CustomGBForce::ComputationType> energyTypes;
538
    std::map<std::string, int> tabulatedFunctionUpdateCount;
539
540
541
    NonbondedMethod nonbondedMethod;
};

542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
/**
 * This kernel is invoked by CustomManyParticleForce to calculate the forces acting on the system and the energy of the system.
 */
class CpuCalcCustomManyParticleForceKernel : public CalcCustomManyParticleForceKernel {
public:
    CpuCalcCustomManyParticleForceKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data) : CalcCustomManyParticleForceKernel(name, platform),
            data(data), ixn(NULL) {
    }
    ~CpuCalcCustomManyParticleForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomManyParticleForce this kernel will be used for
     */
    void initialize(const System& system, const CustomManyParticleForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomManyParticleForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CustomManyParticleForce& force);
private:
    CpuPlatform::PlatformData& data;
    int numParticles;
peastman's avatar
peastman committed
577
    double cutoffDistance;
578
    std::vector<std::vector<double> > particleParamArray;
579
580
    CpuCustomManyParticleForce* ixn;
    std::vector<std::string> globalParameterNames;
581
    std::map<std::string, int> tabulatedFunctionUpdateCount;
582
583
584
    NonbondedMethod nonbondedMethod;
};

585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
/**
 * This kernel is invoked by GayBerneForce to calculate the forces acting on the system.
 */
class CpuCalcGayBerneForceKernel : public CalcGayBerneForceKernel {
public:
    CpuCalcGayBerneForceKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data) : CalcGayBerneForceKernel(name, platform),
            data(data), ixn(NULL) {
    }
    ~CpuCalcGayBerneForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the GayBerneForce this kernel will be used for
     */
    void initialize(const System& system, const GayBerneForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the GayBerneForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const GayBerneForce& force);
private:
    CpuPlatform::PlatformData& data;
    CpuGayBerneForce* ixn;
};

Evan Pretti's avatar
Evan Pretti committed
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
/**
 * This kernel is invoked by LCPOForce to calculate the forces acting on the system and the energy of the system.
 */
class CpuCalcLCPOForceKernel : public CalcLCPOForceKernel {
public:
    CpuCalcLCPOForceKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data) : CalcLCPOForceKernel(name, platform),
            data(data), ixn(NULL) {
    }
    ~CpuCalcLCPOForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the LCPOForce this kernel will be used for
     */
    void initialize(const System& system, const LCPOForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
     * @param context        the context to copy parameters to
     * @param force          the LCPOForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const LCPOForce& force);
private:
    CpuPlatform::PlatformData& data;
    CpuLCPOForce* ixn;
    bool doInteraction;
    double oneBodyEnergy;
    std::vector<int> activeParticles;
    std::vector<fvec4> parameters;
};

662
/**
663
 * This kernel is invoked by LangevinMiddleIntegrator to take one time step.
664
 */
665
class CpuIntegrateLangevinMiddleStepKernel : public IntegrateLangevinMiddleStepKernel {
666
public:
667
    CpuIntegrateLangevinMiddleStepKernel(std::string name, const Platform& platform, CpuPlatform::PlatformData& data) : IntegrateLangevinMiddleStepKernel(name, platform),
668
669
            data(data), dynamics(0) {
    }
670
    ~CpuIntegrateLangevinMiddleStepKernel();
671
672
    /**
     * Initialize the kernel, setting up the particle masses.
Evan Pretti's avatar
Evan Pretti committed
673
     *
674
     * @param system     the System this kernel will be applied to
675
     * @param integrator the LangevinMiddleIntegrator this kernel will be used for
676
     */
677
    void initialize(const System& system, const LangevinMiddleIntegrator& integrator);
678
679
    /**
     * Execute the kernel.
Evan Pretti's avatar
Evan Pretti committed
680
     *
681
     * @param context    the context in which to execute this kernel
682
     * @param integrator the LangevinMiddleIntegrator this kernel is being used for
683
     */
684
    void execute(ContextImpl& context, const LangevinMiddleIntegrator& integrator);
685
686
    /**
     * Compute the kinetic energy.
Evan Pretti's avatar
Evan Pretti committed
687
     *
688
     * @param context    the context in which to execute this kernel
689
     * @param integrator the LangevinMiddleIntegrator this kernel is being used for
690
     */
691
    double computeKineticEnergy(ContextImpl& context, const LangevinMiddleIntegrator& integrator);
692
693
private:
    CpuPlatform::PlatformData& data;
694
    CpuLangevinMiddleDynamics* dynamics;
695
696
697
698
    std::vector<double> masses;
    double prevTemp, prevFriction, prevStepSize;
};

699
700
701
702
} // namespace OpenMM

#endif /*OPENMM_CPUKERNELS_H_*/