kernels.h 59.5 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
#ifndef OPENMM_KERNELS_H_
#define OPENMM_KERNELS_H_

/* -------------------------------------------------------------------------- *
 *                                   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.               *
 *                                                                            *
12
 * Portions copyright (c) 2008-2018 Stanford University and the Authors.      *
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
 * 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.                                     *
 * -------------------------------------------------------------------------- */

35
36
#include "openmm/AndersenThermostat.h"
#include "openmm/BrownianIntegrator.h"
37
#include "openmm/CMAPTorsionForce.h"
38
#include "openmm/CMMotionRemover.h"
39
#include "openmm/CustomAngleForce.h"
40
#include "openmm/CustomBondForce.h"
41
#include "openmm/CustomCentroidBondForce.h"
peastman's avatar
peastman committed
42
#include "openmm/CustomCVForce.h"
43
#include "openmm/CustomCompoundBondForce.h"
44
#include "openmm/CustomExternalForce.h"
45
#include "openmm/CustomGBForce.h"
46
#include "openmm/CustomHbondForce.h"
47
#include "openmm/CustomIntegrator.h"
48
#include "openmm/CustomNonbondedForce.h"
49
#include "openmm/CustomManyParticleForce.h"
50
#include "openmm/CustomTorsionForce.h"
peastman's avatar
peastman committed
51
#include "openmm/GayBerneForce.h"
52
53
54
55
56
#include "openmm/GBSAOBCForce.h"
#include "openmm/HarmonicAngleForce.h"
#include "openmm/HarmonicBondForce.h"
#include "openmm/KernelImpl.h"
#include "openmm/LangevinIntegrator.h"
57
#include "openmm/MonteCarloBarostat.h"
58
59
#include "openmm/PeriodicTorsionForce.h"
#include "openmm/RBTorsionForce.h"
60
#include "openmm/RMSDForce.h"
61
62
#include "openmm/NonbondedForce.h"
#include "openmm/System.h"
63
64
#include "openmm/VariableLangevinIntegrator.h"
#include "openmm/VariableVerletIntegrator.h"
65
#include "openmm/VerletIntegrator.h"
Peter Eastman's avatar
Peter Eastman committed
66
#include <iosfwd>
67
68
69
70
71
72
#include <set>
#include <string>
#include <vector>

namespace OpenMM {

73
/**
74
75
76
 * 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.
77
 */
78
class CalcForcesAndEnergyKernel : public KernelImpl {
79
80
public:
    static std::string Name() {
81
        return "CalcForcesAndEnergyKernel";
82
    }
83
    CalcForcesAndEnergyKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
84
85
86
87
88
89
90
91
    }
    /**
     * Initialize the kernel.
     * 
     * @param system     the System this kernel will be applied to
     */
    virtual void initialize(const System& system) = 0;
    /**
92
     * This is called at the beginning of each force/energy computation, before calcForcesAndEnergy() has been called on
93
     * any ForceImpl.
94
     * 
95
96
97
     * @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
98
     * @param groups        a set of bit flags for which force groups to include
99
     */
100
    virtual void beginComputation(ContextImpl& context, bool includeForce, bool includeEnergy, int groups) = 0;
101
    /**
102
     * This is called at the end of each force/energy computation, after calcForcesAndEnergy() has been called on
103
104
     * every ForceImpl.
     *
105
106
107
     * @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
108
     * @param groups        a set of bit flags for which force groups to include
109
110
     * @param valid         the method may set this to false to indicate the results are invalid and the force/energy
     *                      calculation should be repeated
111
     * @return the potential energy of the system.  This value is added to all values returned by ForceImpls'
112
     * calcForcesAndEnergy() methods.  That is, each force kernel may <i>either</i> return its contribution to the
113
114
     * energy directly, <i>or</i> add it to an internal buffer so that it will be included here.
     */
115
    virtual double finishComputation(ContextImpl& context, bool includeForce, bool includeEnergy, int groups, bool& valid) = 0;
116
117
};

118
/**
119
120
 * This kernel provides methods for setting and retrieving various state data: time, positions,
 * velocities, and forces.
121
 */
122
class UpdateStateDataKernel : public KernelImpl {
123
124
125
126
public:
    static std::string Name() {
        return "UpdateTime";
    }
127
    UpdateStateDataKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
128
129
130
131
132
133
134
135
136
137
138
139
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     */
    virtual void initialize(const System& system) = 0;
    /**
     * Get the current time (in picoseconds).
     *
     * @param context    the context in which to execute this kernel
     */
140
    virtual double getTime(const ContextImpl& context) const = 0;
141
142
143
144
145
146
    /**
     * Set the current time (in picoseconds).
     *
     * @param context    the context in which to execute this kernel
     * @param time       the time
     */
147
    virtual void setTime(ContextImpl& context, double time) = 0;
148
149
150
151
152
153
154
155
156
    /**
     * Get the positions of all particles.
     *
     * @param positions  on exit, this contains the particle positions
     */
    virtual void getPositions(ContextImpl& context, std::vector<Vec3>& positions) = 0;
    /**
     * Set the positions of all particles.
     *
157
     * @param positions  a vector containing the particle positions
158
159
160
161
162
163
164
165
166
167
168
     */
    virtual void setPositions(ContextImpl& context, const std::vector<Vec3>& positions) = 0;
    /**
     * Get the velocities of all particles.
     *
     * @param velocities  on exit, this contains the particle velocities
     */
    virtual void getVelocities(ContextImpl& context, std::vector<Vec3>& velocities) = 0;
    /**
     * Set the velocities of all particles.
     *
169
     * @param velocities  a vector containing the particle velocities
170
171
172
173
174
175
176
177
     */
    virtual void setVelocities(ContextImpl& context, const std::vector<Vec3>& velocities) = 0;
    /**
     * Get the current forces on all particles.
     *
     * @param forces  on exit, this contains the forces
     */
    virtual void getForces(ContextImpl& context, std::vector<Vec3>& forces) = 0;
178
179
180
181
182
183
    /**
     * Get the current derivatives of the energy with respect to context parameters.
     *
     * @param derivs  on exit, this contains the derivatives
     */
    virtual void getEnergyParameterDerivatives(ContextImpl& context, std::map<std::string, double>& derivs) = 0;
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
    /**
     * Get the current periodic box vectors.
     *
     * @param a      on exit, this contains the vector defining the first edge of the periodic box
     * @param b      on exit, this contains the vector defining the second edge of the periodic box
     * @param c      on exit, this contains the vector defining the third edge of the periodic box
     */
    virtual void getPeriodicBoxVectors(ContextImpl& context, Vec3& a, Vec3& b, Vec3& c) const = 0;
    /**
     * Set the current periodic box vectors.
     *
     * @param a      the vector defining the first edge of the periodic box
     * @param b      the vector defining the second edge of the periodic box
     * @param c      the vector defining the third edge of the periodic box
     */
199
    virtual void setPeriodicBoxVectors(ContextImpl& context, const Vec3& a, const Vec3& b, const Vec3& c) = 0;
Peter Eastman's avatar
Peter Eastman committed
200
201
202
203
204
205
206
207
208
209
210
211
    /**
     * Create a checkpoint recording the current state of the Context.
     * 
     * @param stream    an output stream the checkpoint data should be written to
     */
    virtual void createCheckpoint(ContextImpl& context, std::ostream& stream) = 0;
    /**
     * Load a checkpoint that was written by createCheckpoint().
     * 
     * @param stream    an input stream the checkpoint data should be read from
     */
    virtual void loadCheckpoint(ContextImpl& context, std::istream& stream) = 0;
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
/**
 * This kernel modifies the positions of particles to enforce distance constraints.
 */
class ApplyConstraintsKernel : public KernelImpl {
public:
    static std::string Name() {
        return "ApplyConstraints";
    }
    ApplyConstraintsKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     */
    virtual void initialize(const System& system) = 0;
    /**
     * Update particle positions to enforce constraints.
     *
     * @param context    the context in which to execute this kernel
     * @param tol        the distance tolerance within which constraints must be satisfied.
     */
    virtual void apply(ContextImpl& context, double tol) = 0;
237
238
239
240
241
242
243
    /**
     * Update particle velocities to enforce constraints.
     *
     * @param context    the context in which to execute this kernel
     * @param tol        the velocity tolerance within which constraints must be satisfied.
     */
    virtual void applyToVelocities(ContextImpl& context, double tol) = 0;
244
245
};

246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
/**
 * This kernel recomputes the positions of virtual sites.
 */
class VirtualSitesKernel : public KernelImpl {
public:
    static std::string Name() {
        return "VirtualSites";
    }
    VirtualSitesKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     */
    virtual void initialize(const System& system) = 0;
    /**
     * Compute the virtual site locations.
     *
     * @param context    the context in which to execute this kernel
     */
    virtual void computePositions(ContextImpl& context) = 0;
};

270
/**
271
 * This kernel is invoked by HarmonicBondForce to calculate the forces acting on the system and the energy of the system.
272
 */
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
class CalcHarmonicBondForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcHarmonicBondForce";
    }
    CalcHarmonicBondForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     * 
     * @param system     the System this kernel will be applied to
     * @param force      the HarmonicBondForce this kernel will be used for
     */
    virtual void initialize(const System& system, const HarmonicBondForce& force) = 0;
    /**
288
     * Execute the kernel to calculate the forces and/or energy.
289
     * 
290
291
292
293
     * @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
294
     */
295
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
296
297
298
299
300
301
302
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the HarmonicBondForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const HarmonicBondForce& force) = 0;
303
304
};

305
306
307
308
309
310
311
312
313
314
315
316
317
318
/**
 * This kernel is invoked by CustomBondForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcCustomBondForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcCustomBondForce";
    }
    CalcCustomBondForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
319
     * @param force      the CustomBondForce this kernel will be used for
320
321
322
     */
    virtual void initialize(const System& system, const CustomBondForce& force) = 0;
    /**
323
     * Execute the kernel to calculate the forces and/or energy.
324
     *
325
326
327
328
     * @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
329
     */
330
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
331
332
333
334
335
336
337
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomBondForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const CustomBondForce& force) = 0;
338
339
};

340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
/**
 * This kernel is invoked by HarmonicAngleForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcHarmonicAngleForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcHarmonicAngleForce";
    }
    CalcHarmonicAngleForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     * 
     * @param system     the System this kernel will be applied to
     * @param force      the HarmonicAngleForce this kernel will be used for
     */
    virtual void initialize(const System& system, const HarmonicAngleForce& force) = 0;
    /**
358
359
360
361
362
363
     * 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
364
     */
365
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
366
367
368
369
370
371
372
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the HarmonicAngleForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const HarmonicAngleForce& force) = 0;
373
374
};

375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
/**
 * This kernel is invoked by CustomAngleForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcCustomAngleForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcCustomAngleForce";
    }
    CalcCustomAngleForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomAngleForce this kernel will be used for
     */
    virtual void initialize(const System& system, const CustomAngleForce& force) = 0;
    /**
393
     * Execute the kernel to calculate the forces and/or energy.
394
     *
395
396
397
398
     * @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
399
     */
400
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
401
402
403
404
405
406
407
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomAngleForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const CustomAngleForce& force) = 0;
408
409
};

410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
/**
 * This kernel is invoked by PeriodicTorsionForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcPeriodicTorsionForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcPeriodicTorsionForce";
    }
    CalcPeriodicTorsionForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     * 
     * @param system     the System this kernel will be applied to
     * @param force      the PeriodicTorsionForce this kernel will be used for
     */
    virtual void initialize(const System& system, const PeriodicTorsionForce& force) = 0;
    /**
428
429
430
431
432
433
     * 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
434
     */
435
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
436
437
438
439
440
441
442
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the PeriodicTorsionForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const PeriodicTorsionForce& force) = 0;
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
};

/**
 * This kernel is invoked by RBTorsionForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcRBTorsionForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcRBTorsionForce";
    }
    CalcRBTorsionForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     * 
     * @param system     the System this kernel will be applied to
     * @param force      the RBTorsionForce this kernel will be used for
     */
    virtual void initialize(const System& system, const RBTorsionForce& force) = 0;
    /**
463
464
465
466
467
468
     * 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
469
     */
470
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
471
472
473
474
475
476
477
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the RBTorsionForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const RBTorsionForce& force) = 0;
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
};

/**
 * This kernel is invoked by CMAPTorsionForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcCMAPTorsionForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcCMAPTorsionForce";
    }
    CalcCMAPTorsionForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CMAPTorsionForce this kernel will be used for
     */
    virtual void initialize(const System& system, const CMAPTorsionForce& force) = 0;
    /**
498
     * Execute the kernel to calculate the forces and/or energy.
499
     *
500
501
502
503
     * @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
504
     */
505
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
506
507
508
509
510
511
512
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CMAPTorsionForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const CMAPTorsionForce& force) = 0;
513
514
};

515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
/**
 * This kernel is invoked by CustomTorsionForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcCustomTorsionForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcCustomTorsionForce";
    }
    CalcCustomTorsionForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomTorsionForce this kernel will be used for
     */
    virtual void initialize(const System& system, const CustomTorsionForce& force) = 0;
    /**
533
     * Execute the kernel to calculate the forces and/or energy.
534
     *
535
536
537
538
     * @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
539
     */
540
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
541
542
543
544
545
546
547
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomTorsionForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const CustomTorsionForce& force) = 0;
548
549
};

550
551
552
553
/**
 * This kernel is invoked by NonbondedForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcNonbondedForceKernel : public KernelImpl {
554
public:
555
556
557
    enum NonbondedMethod {
        NoCutoff = 0,
        CutoffNonPeriodic = 1,
558
        CutoffPeriodic = 2,
559
        Ewald = 3,
560
561
        PME = 4,
        LJPME = 5
562
    };
563
    static std::string Name() {
564
        return "CalcNonbondedForce";
565
    }
566
    CalcNonbondedForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
567
568
    }
    /**
569
     * Initialize the kernel.
570
     * 
571
     * @param system     the System this kernel will be applied to
572
     * @param force      the NonbondedForce this kernel will be used for
573
     */
574
    virtual void initialize(const System& system, const NonbondedForce& force) = 0;
575
    /**
576
577
578
579
580
     * 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
581
582
     * @param includeDirect  true if direct space interactions should be included
     * @param includeReciprocal  true if reciprocal space interactions should be included
583
     * @return the potential energy due to the force
584
     */
585
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy, bool includeDirect, bool includeReciprocal) = 0;
586
587
588
589
590
591
592
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the NonbondedForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const NonbondedForce& force) = 0;
593
594
    /**
     * Get the parameters being used for PME.
595
     *
596
597
598
599
600
601
     * @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
     */
    virtual void getPMEParameters(double& alpha, int& nx, int& ny, int& nz) const = 0;
602
603
604
    /**
     * Get the parameters being used for the dispersion terms in LJPME.
     *
605
606
607
608
     * @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
609
     */
610
    virtual void getLJPMEParameters(double& alpha, int& nx, int& ny, int& nz) const = 0;
611
612
};

613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
/**
 * This kernel is invoked by CustomNonbondedForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcCustomNonbondedForceKernel : public KernelImpl {
public:
    enum NonbondedMethod {
        NoCutoff = 0,
        CutoffNonPeriodic = 1,
        CutoffPeriodic = 2
    };
    static std::string Name() {
        return "CalcCustomNonbondedForce";
    }
    CalcCustomNonbondedForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomNonbondedForce this kernel will be used for
     */
    virtual void initialize(const System& system, const CustomNonbondedForce& force) = 0;
    /**
636
     * Execute the kernel to calculate the forces and/or energy.
637
     *
638
639
640
641
     * @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
642
     */
643
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
644
645
646
647
648
649
650
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomNonbondedForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const CustomNonbondedForce& force) = 0;
651
652
};

653
/**
654
 * This kernel is invoked by GBSAOBCForce to calculate the forces acting on the system and the energy of the system.
655
 */
656
class CalcGBSAOBCForceKernel : public KernelImpl {
657
658
public:
    static std::string Name() {
Mark Friedrichs's avatar
Mark Friedrichs committed
659
        return "CalcGBSAOBCForce";
660
    }
661
    CalcGBSAOBCForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
662
663
    }
    /**
664
     * Initialize the kernel.
665
     * 
666
     * @param system     the System this kernel will be applied to
667
     * @param force      the GBSAOBCForce this kernel will be used for
668
     */
669
    virtual void initialize(const System& system, const GBSAOBCForce& force) = 0;
670
    /**
671
672
673
674
675
676
     * 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
677
     */
678
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
679
680
681
682
683
684
685
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the GBSAOBCForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const GBSAOBCForce& force) = 0;
686
687
};

688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
/**
 * This kernel is invoked by CustomGBForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcCustomGBForceKernel : public KernelImpl {
public:
    enum NonbondedMethod {
        NoCutoff = 0,
        CutoffNonPeriodic = 1,
        CutoffPeriodic = 2
    };
    static std::string Name() {
        return "CalcCustomGBForce";
    }
    CalcCustomGBForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomGBForce this kernel will be used for
     */
    virtual void initialize(const System& system, const CustomGBForce& force) = 0;
    /**
711
     * Execute the kernel to calculate the forces and/or energy.
712
     *
713
714
715
716
     * @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
717
     */
718
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
719
720
721
722
723
724
725
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomGBForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const CustomGBForce& force) = 0;
726
727
};

728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
/**
 * This kernel is invoked by CustomExternalForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcCustomExternalForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcCustomExternalForce";
    }
    CalcCustomExternalForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomExternalForce this kernel will be used for
     */
    virtual void initialize(const System& system, const CustomExternalForce& force) = 0;
    /**
746
     * Execute the kernel to calculate the forces and/or energy.
747
     *
748
749
750
751
     * @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
752
     */
753
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
754
755
756
757
758
759
760
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomExternalForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const CustomExternalForce& force) = 0;
761
762
};

763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
/**
 * This kernel is invoked by CustomHbondForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcCustomHbondForceKernel : public KernelImpl {
public:
    enum NonbondedMethod {
        NoCutoff = 0,
        CutoffNonPeriodic = 1,
        CutoffPeriodic = 2
    };
    static std::string Name() {
        return "CalcCustomHbondForce";
    }
    CalcCustomHbondForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomHbondForce this kernel will be used for
     */
    virtual void initialize(const System& system, const CustomHbondForce& force) = 0;
    /**
786
     * Execute the kernel to calculate the forces and/or energy.
787
     *
788
789
790
791
     * @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
792
     */
793
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
794
795
796
797
798
799
800
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomHbondForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const CustomHbondForce& force) = 0;
801
802
};

803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
/**
 * This kernel is invoked by CustomCentroidBondForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcCustomCentroidBondForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcCustomCentroidBondForce";
    }
    CalcCustomCentroidBondForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomCentroidBondForce this kernel will be used for
     */
    virtual void initialize(const System& system, const CustomCentroidBondForce& force) = 0;
    /**
     * 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
     */
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomCentroidBondForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const CustomCentroidBondForce& force) = 0;
};

838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
/**
 * This kernel is invoked by CustomCompoundBondForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcCustomCompoundBondForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcCustomCompoundBondForce";
    }
    CalcCustomCompoundBondForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomCompoundBondForce this kernel will be used for
     */
    virtual void initialize(const System& system, const CustomCompoundBondForce& force) = 0;
    /**
     * 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
     */
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
864
865
866
867
868
869
870
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomCompoundBondForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const CustomCompoundBondForce& force) = 0;
871
872
};

873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
/**
 * This kernel is invoked by CustomManyParticleForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcCustomManyParticleForceKernel : public KernelImpl {
public:
    enum NonbondedMethod {
        NoCutoff = 0,
        CutoffNonPeriodic = 1,
        CutoffPeriodic = 2
    };
    static std::string Name() {
        return "CalcCustomManyParticleForce";
    }
    CalcCustomManyParticleForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomManyParticleForce this kernel will be used for
     */
    virtual void initialize(const System& system, const CustomManyParticleForce& force) = 0;
    /**
     * 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
     */
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomManyParticleForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const CustomManyParticleForce& force) = 0;
};

peastman's avatar
peastman committed
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
/**
 * This kernel is invoked by GayBerneForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcGayBerneForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcGayBerneForce";
    }
    CalcGayBerneForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the GayBerneForce this kernel will be used for
     */
    virtual void initialize(const System& system, const GayBerneForce& force) = 0;
    /**
     * 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
     */
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the GayBerneForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const GayBerneForce& force) = 0;
};

peastman's avatar
peastman committed
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
/**
 * This kernel is invoked by CustomCVForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcCustomCVForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcCustomCVForce";
    }
    CalcCustomCVForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomCVForce this kernel will be used for
963
     * @param innerContext   the context created by the CustomCVForce for computing collective variables
peastman's avatar
peastman committed
964
     */
965
    virtual void initialize(const System& system, const CustomCVForce& force, ContextImpl& innerContext) = 0;
peastman's avatar
peastman committed
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param innerContext   the context created by the CustomCVForce for computing collective variables
     * @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
     */
    virtual double execute(ContextImpl& context, ContextImpl& innerContext, bool includeForces, bool includeEnergy) = 0;
    /**
     * Copy state information to the inner context.
     *
     * @param context        the context in which to execute this kernel
     * @param innerContext   the context created by the CustomCVForce for computing collective variables
     */
    virtual void copyState(ContextImpl& context, ContextImpl& innerContext) = 0;
983
984
985
986
987
988
989
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomCVForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const CustomCVForce& force) = 0;
peastman's avatar
peastman committed
990
991
};

992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
/**
 * This kernel is invoked by RMSDForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcRMSDForceKernel : public KernelImpl {
public:
    static std::string Name() {
        return "CalcRMSDForce";
    }
    CalcRMSDForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the RMSDForce this kernel will be used for
     */
    virtual void initialize(const System& system, const RMSDForce& force) = 0;
    /**
     * 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
     */
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the RMSDForce to copy the parameters from
     */
    virtual void copyParametersToContext(ContextImpl& context, const RMSDForce& force) = 0;
};

1027
1028
1029
1030
1031
1032
1033
1034
/**
 * This kernel is invoked by VerletIntegrator to take one time step.
 */
class IntegrateVerletStepKernel : public KernelImpl {
public:
    static std::string Name() {
        return "IntegrateVerletStep";
    }
1035
    IntegrateVerletStepKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
1036
1037
    }
    /**
1038
     * Initialize the kernel.
1039
     * 
1040
1041
     * @param system     the System this kernel will be applied to
     * @param integrator the VerletIntegrator this kernel will be used for
1042
     */
1043
    virtual void initialize(const System& system, const VerletIntegrator& integrator) = 0;
1044
1045
1046
    /**
     * Execute the kernel.
     * 
1047
1048
     * @param context    the context in which to execute this kernel
     * @param integrator the VerletIntegrator this kernel is being used for
1049
     */
1050
    virtual void execute(ContextImpl& context, const VerletIntegrator& integrator) = 0;
1051
1052
1053
1054
1055
1056
1057
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the VerletIntegrator this kernel is being used for
     */
    virtual double computeKineticEnergy(ContextImpl& context, const VerletIntegrator& integrator) = 0;
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
};

/**
 * This kernel is invoked by LangevinIntegrator to take one time step.
 */
class IntegrateLangevinStepKernel : public KernelImpl {
public:
    static std::string Name() {
        return "IntegrateLangevinStep";
    }
1068
    IntegrateLangevinStepKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
1069
1070
    }
    /**
1071
     * Initialize the kernel.
1072
     * 
1073
1074
     * @param system     the System this kernel will be applied to
     * @param integrator the LangevinIntegrator this kernel will be used for
1075
     */
1076
    virtual void initialize(const System& system, const LangevinIntegrator& integrator) = 0;
1077
1078
1079
    /**
     * Execute the kernel.
     * 
1080
1081
     * @param context    the context in which to execute this kernel
     * @param integrator the LangevinIntegrator this kernel is being used for
1082
     */
1083
    virtual void execute(ContextImpl& context, const LangevinIntegrator& integrator) = 0;
1084
1085
1086
1087
1088
1089
1090
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the LangevinIntegrator this kernel is being used for
     */
    virtual double computeKineticEnergy(ContextImpl& context, const LangevinIntegrator& integrator) = 0;
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
};

/**
 * This kernel is invoked by BrownianIntegrator to take one time step.
 */
class IntegrateBrownianStepKernel : public KernelImpl {
public:
    static std::string Name() {
        return "IntegrateBrownianStep";
    }
1101
    IntegrateBrownianStepKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
1102
1103
    }
    /**
1104
     * Initialize the kernel.
1105
     * 
1106
1107
     * @param system     the System this kernel will be applied to
     * @param integrator the BrownianIntegrator this kernel will be used for
1108
     */
1109
    virtual void initialize(const System& system, const BrownianIntegrator& integrator) = 0;
1110
1111
1112
    /**
     * Execute the kernel.
     * 
1113
1114
     * @param context    the context in which to execute this kernel
     * @param integrator the BrownianIntegrator this kernel is being used for
1115
     */
1116
    virtual void execute(ContextImpl& context, const BrownianIntegrator& integrator) = 0;
1117
1118
1119
1120
1121
1122
1123
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the BrownianIntegrator this kernel is being used for
     */
    virtual double computeKineticEnergy(ContextImpl& context, const BrownianIntegrator& integrator) = 0;
1124
1125
};

1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
/**
 * This kernel is invoked by VariableLangevinIntegrator to take one time step.
 */
class IntegrateVariableLangevinStepKernel : public KernelImpl {
public:
    static std::string Name() {
        return "IntegrateVariableLangevinStep";
    }
    IntegrateVariableLangevinStepKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param integrator the VariableLangevinIntegrator this kernel will be used for
     */
    virtual void initialize(const System& system, const VariableLangevinIntegrator& integrator) = 0;
    /**
     * Execute the kernel.
     *
     * @param context    the context in which to execute this kernel
1147
     * @param integrator the VariableLangevinIntegrator this kernel is being used for
1148
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1149
     * @return the size of the step that was taken
1150
     */
1151
    virtual double execute(ContextImpl& context, const VariableLangevinIntegrator& integrator, double maxTime) = 0;
1152
1153
1154
1155
1156
1157
1158
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the VariableLangevinIntegrator this kernel is being used for
     */
    virtual double computeKineticEnergy(ContextImpl& context, const VariableLangevinIntegrator& integrator) = 0;
1159
1160
};

1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
/**
 * This kernel is invoked by VariableVerletIntegrator to take one time step.
 */
class IntegrateVariableVerletStepKernel : public KernelImpl {
public:
    static std::string Name() {
        return "IntegrateVariableVerletStep";
    }
    IntegrateVariableVerletStepKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param integrator the VariableVerletIntegrator this kernel will be used for
     */
    virtual void initialize(const System& system, const VariableVerletIntegrator& integrator) = 0;
    /**
     * Execute the kernel.
     *
     * @param context    the context in which to execute this kernel
1182
     * @param integrator the VariableVerletIntegrator this kernel is being used for
1183
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1184
     * @return the size of the step that was taken
1185
     */
1186
    virtual double execute(ContextImpl& context, const VariableVerletIntegrator& integrator, double maxTime) = 0;
1187
1188
1189
1190
1191
1192
1193
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the VariableVerletIntegrator this kernel is being used for
     */
    virtual double computeKineticEnergy(ContextImpl& context, const VariableVerletIntegrator& integrator) = 0;
1194
1195
};

1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
/**
 * This kernel is invoked by CustomIntegrator to take one time step.
 */
class IntegrateCustomStepKernel : public KernelImpl {
public:
    static std::string Name() {
        return "IntegrateCustomStep";
    }
    IntegrateCustomStepKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     * 
     * @param system     the System this kernel will be applied to
     * @param integrator the CustomIntegrator this kernel will be used for
     */
    virtual void initialize(const System& system, const CustomIntegrator& integrator) = 0;
    /**
     * Execute the kernel.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the CustomIntegrator this kernel is being used for
     * @param forcesAreValid if the context has been modified since the last time step, this will be
     *                       false to show that cached forces are invalid and must be recalculated.
     *                       On exit, this should specify whether the cached forces are valid at the
     *                       end of the step.
     */
    virtual void execute(ContextImpl& context, CustomIntegrator& integrator, bool& forcesAreValid) = 0;
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the CustomIntegrator this kernel is being used for
     * @param forcesAreValid if the context has been modified since the last time step, this will be
     *                       false to show that cached forces are invalid and must be recalculated.
     *                       On exit, this should specify whether the cached forces are valid at the
     *                       end of the step.
     */
    virtual double computeKineticEnergy(ContextImpl& context, CustomIntegrator& integrator, bool& forcesAreValid) = 0;
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
    /**
     * Get the values of all global variables.
     *
     * @param context   the context in which to execute this kernel
     * @param values    on exit, this contains the values
     */
    virtual void getGlobalVariables(ContextImpl& context, std::vector<double>& values) const = 0;
    /**
     * Set the values of all global variables.
     *
     * @param context   the context in which to execute this kernel
     * @param values    a vector containing the values
     */
    virtual void setGlobalVariables(ContextImpl& context, const std::vector<double>& values) = 0;
    /**
     * Get the values of a per-DOF variable.
     *
     * @param context   the context in which to execute this kernel
     * @param variable  the index of the variable to get
     * @param values    on exit, this contains the values
     */
    virtual void getPerDofVariable(ContextImpl& context, int variable, std::vector<Vec3>& values) const = 0;
    /**
     * Set the values of a per-DOF variable.
     *
     * @param context   the context in which to execute this kernel
     * @param variable  the index of the variable to get
     * @param values    a vector containing the values
     */
    virtual void setPerDofVariable(ContextImpl& context, int variable, const std::vector<Vec3>& values) = 0;
};

1267
/**
Peter Eastman's avatar
Peter Eastman committed
1268
 * This kernel is invoked by AndersenThermostat at the start of each time step to adjust the particle velocities.
1269
1270
1271
1272
1273
1274
 */
class ApplyAndersenThermostatKernel : public KernelImpl {
public:
    static std::string Name() {
        return "ApplyAndersenThermostat";
    }
1275
    ApplyAndersenThermostatKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
1276
1277
    }
    /**
1278
     * Initialize the kernel.
1279
     * 
1280
1281
     * @param system     the System this kernel will be applied to
     * @param thermostat the AndersenThermostat this kernel will be used for
1282
     */
1283
    virtual void initialize(const System& system, const AndersenThermostat& thermostat) = 0;
1284
1285
1286
    /**
     * Execute the kernel.
     * 
1287
     * @param context    the context in which to execute this kernel
1288
     */
1289
    virtual void execute(ContextImpl& context) = 0;
1290
1291
};

1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
/**
 * This kernel is invoked by MonteCarloBarostat to adjust the periodic box volume
 */
class ApplyMonteCarloBarostatKernel : public KernelImpl {
public:
    static std::string Name() {
        return "ApplyMonteCarloBarostat";
    }
    ApplyMonteCarloBarostatKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param barostat   the MonteCarloBarostat this kernel will be used for
     */
1308
    virtual void initialize(const System& system, const Force& barostat) = 0;
Lee-Ping's avatar
Lee-Ping committed
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
    /**
     * Attempt a Monte Carlo step, scaling particle positions (or cluster centers) by a specified value.
     * This version scales the x, y, and z positions independently.
     * This is called BEFORE the periodic box size is modified.  It should begin by translating each particle
     * or cluster into the first periodic box, so that coordinates will still be correct after the box size
     * is changed.
     *
     * @param context    the context in which to execute this kernel
     * @param scaleX     the scale factor by which to multiply particle x-coordinate
     * @param scaleY     the scale factor by which to multiply particle y-coordinate
     * @param scaleZ     the scale factor by which to multiply particle z-coordinate
     */
    virtual void scaleCoordinates(ContextImpl& context, double scaleX, double scaleY, double scaleZ) = 0;
    /**
     * Reject the most recent Monte Carlo step, restoring the particle positions to where they were before
1324
1325
1326
1327
1328
1329
1330
     * scaleCoordinates() was last called.
     *
     * @param context    the context in which to execute this kernel
     */
    virtual void restoreCoordinates(ContextImpl& context) = 0;
};

1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
/**
 * This kernel is invoked to remove center of mass motion from the system.
 */
class RemoveCMMotionKernel : public KernelImpl {
public:
    static std::string Name() {
        return "RemoveCMMotion";
    }
    RemoveCMMotionKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
1342
     * Initialize the kernel.
1343
     * 
1344
1345
     * @param system     the System this kernel will be applied to
     * @param force      the CMMotionRemover this kernel will be used for
1346
     */
1347
    virtual void initialize(const System& system, const CMMotionRemover& force) = 0;
1348
1349
1350
    /**
     * Execute the kernel.
     * 
1351
     * @param context    the context in which to execute this kernel
1352
     */
1353
    virtual void execute(ContextImpl& context) = 0;
1354
1355
};

1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
/**
 * This kernel performs the reciprocal space calculation for PME.  In most cases, this
 * calculation is done directly by CalcNonbondedForceKernel so this kernel is unneeded.
 * In some cases it may want to outsource the work to a different kernel.  In particular,
 * GPU based platforms sometimes use a CPU based implementation provided by a separate
 * plugin.
 */
class CalcPmeReciprocalForceKernel : public KernelImpl {
public:
    class IO;
    static std::string Name() {
        return "CalcPmeReciprocalForce";
    }
    CalcPmeReciprocalForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     * 
     * @param gridx        the x size of the PME grid
     * @param gridy        the y size of the PME grid
     * @param gridz        the z size of the PME grid
     * @param numParticles the number of particles in the system
     * @param alpha        the Ewald blending parameter
1379
     * @param deterministic whether it should attempt to make the resulting forces deterministic
1380
     */
1381
    virtual void initialize(int gridx, int gridy, int gridz, int numParticles, double alpha, bool deterministic) = 0;
1382
1383
    /**
     * Begin computing the force and energy.
peastman's avatar
peastman committed
1384
1385
1386
1387
     *
     * @param io                  an object that coordinates data transfer
     * @param periodicBoxVectors  the vectors defining the periodic box (measured in nm)
     * @param includeEnergy       true if potential energy should be computed
1388
     */
peastman's avatar
peastman committed
1389
    virtual void beginComputation(IO& io, const Vec3* periodicBoxVectors, bool includeEnergy) = 0;
1390
1391
1392
1393
1394
1395
1396
    /**
     * Finish computing the force and energy.
     * 
     * @param io   an object that coordinates data transfer
     * @return the potential energy due to the PME reciprocal space interactions
     */
    virtual double finishComputation(IO& io) = 0;
1397
1398
1399
1400
1401
1402
1403
1404
1405
    /**
     * 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
     */
    virtual void getPMEParameters(double& alpha, int& nx, int& ny, int& nz) const = 0;
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
};

/**
 * Any class that uses CalcPmeReciprocalForceKernel should create an implementation of this
 * class, then pass it to the kernel to manage communication with it.
 */
class CalcPmeReciprocalForceKernel::IO {
public:
    /**
     * Get a pointer to the atom charges and positions.  This array should contain four
     * elements for each atom: x, y, z, and q in that order.
     */
    virtual float* getPosq() = 0;
    /**
     * Record the forces calculated by the kernel.
     * 
     * @param force    an array containing four elements for each atom.  The first three
     *                 are the x, y, and z components of the force, while the fourth element
     *                 should be ignored.
     */
    virtual void setForce(float* force) = 0;
};


1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
/**
 * This kernel performs the dispersion reciprocal space calculation for LJPME.  In most cases, this
 * calculation is done directly by CalcNonbondedForceKernel so this kernel is unneeded.
 * In some cases it may want to outsource the work to a different kernel.  In particular,
 * GPU based platforms sometimes use a CPU based implementation provided by a separate
 * plugin.
 */
class CalcDispersionPmeReciprocalForceKernel : public KernelImpl {
public:
    class IO;
    static std::string Name() {
        return "CalcDispersionPmeReciprocalForce";
    }
    CalcDispersionPmeReciprocalForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
    }
    /**
     * Initialize the kernel.
     * 
     * @param gridx        the x size of the PME grid
     * @param gridy        the y size of the PME grid
     * @param gridz        the z size of the PME grid
     * @param numParticles the number of particles in the system
     * @param alpha        the Ewald blending parameter
1453
     * @param deterministic whether it should attempt to make the resulting forces deterministic
1454
     */
1455
    virtual void initialize(int gridx, int gridy, int gridz, int numParticles, double alpha, bool deterministic) = 0;
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
    /**
     * Begin computing the force and energy.
     *
     * @param io                  an object that coordinates data transfer
     * @param periodicBoxVectors  the vectors defining the periodic box (measured in nm)
     * @param includeEnergy       true if potential energy should be computed
     */
    virtual void beginComputation(IO& io, const Vec3* periodicBoxVectors, bool includeEnergy) = 0;
    /**
     * Finish computing the force and energy.
     * 
     * @param io   an object that coordinates data transfer
     * @return the potential energy due to the PME reciprocal space interactions
     */
    virtual double finishComputation(IO& io) = 0;
    /**
     * 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
     */
    virtual void getPMEParameters(double& alpha, int& nx, int& ny, int& nz) const = 0;
};

1482
1483
1484
} // namespace OpenMM

#endif /*OPENMM_KERNELS_H_*/