kernels.h 50.6 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-2015 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/CustomCompoundBondForce.h"
42
#include "openmm/CustomExternalForce.h"
43
#include "openmm/CustomGBForce.h"
44
#include "openmm/CustomHbondForce.h"
45
#include "openmm/CustomIntegrator.h"
46
#include "openmm/CustomNonbondedForce.h"
47
#include "openmm/CustomManyParticleForce.h"
48
#include "openmm/CustomTorsionForce.h"
49
50
51
52
53
54
#include "openmm/GBSAOBCForce.h"
#include "openmm/GBVIForce.h"
#include "openmm/HarmonicAngleForce.h"
#include "openmm/HarmonicBondForce.h"
#include "openmm/KernelImpl.h"
#include "openmm/LangevinIntegrator.h"
55
#include "openmm/MonteCarloBarostat.h"
56
57
58
59
#include "openmm/PeriodicTorsionForce.h"
#include "openmm/RBTorsionForce.h"
#include "openmm/NonbondedForce.h"
#include "openmm/System.h"
60
61
#include "openmm/VariableLangevinIntegrator.h"
#include "openmm/VariableVerletIntegrator.h"
62
#include "openmm/VerletIntegrator.h"
Peter Eastman's avatar
Peter Eastman committed
63
#include <iosfwd>
64
65
66
67
68
69
#include <set>
#include <string>
#include <vector>

namespace OpenMM {

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

115
/**
116
117
 * This kernel provides methods for setting and retrieving various state data: time, positions,
 * velocities, and forces.
118
 */
119
class UpdateStateDataKernel : public KernelImpl {
120
121
122
123
public:
    static std::string Name() {
        return "UpdateTime";
    }
124
    UpdateStateDataKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
125
126
127
128
129
130
131
132
133
134
135
136
    }
    /**
     * 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
     */
137
    virtual double getTime(const ContextImpl& context) const = 0;
138
139
140
141
142
143
    /**
     * Set the current time (in picoseconds).
     *
     * @param context    the context in which to execute this kernel
     * @param time       the time
     */
144
    virtual void setTime(ContextImpl& context, double time) = 0;
145
146
147
148
149
150
151
152
153
    /**
     * 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.
     *
154
     * @param positions  a vector containing the particle positions
155
156
157
158
159
160
161
162
163
164
165
     */
    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.
     *
166
     * @param velocities  a vector containing the particle velocities
167
168
169
170
171
172
173
174
     */
    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;
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
    /**
     * 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
     */
    virtual void setPeriodicBoxVectors(ContextImpl& context, const Vec3& a, const Vec3& b, const Vec3& c) const = 0;
Peter Eastman's avatar
Peter Eastman committed
191
192
193
194
195
196
197
198
199
200
201
202
    /**
     * 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;
203
204
};

205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
/**
 * 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;
228
229
230
231
232
233
234
    /**
     * 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;
235
236
};

237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
/**
 * 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;
};

261
/**
262
 * This kernel is invoked by HarmonicBondForce to calculate the forces acting on the system and the energy of the system.
263
 */
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
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;
    /**
279
     * Execute the kernel to calculate the forces and/or energy.
280
     * 
281
282
283
284
     * @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
285
     */
286
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
287
288
289
290
291
292
293
    /**
     * 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;
294
295
};

296
297
298
299
300
301
302
303
304
305
306
307
308
309
/**
 * 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
310
     * @param force      the CustomBondForce this kernel will be used for
311
312
313
     */
    virtual void initialize(const System& system, const CustomBondForce& force) = 0;
    /**
314
     * Execute the kernel to calculate the forces and/or energy.
315
     *
316
317
318
319
     * @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
320
     */
321
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
322
323
324
325
326
327
328
    /**
     * 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;
329
330
};

331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
/**
 * 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;
    /**
349
350
351
352
353
354
     * 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
355
     */
356
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
357
358
359
360
361
362
363
    /**
     * 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;
364
365
};

366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
/**
 * 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;
    /**
384
     * Execute the kernel to calculate the forces and/or energy.
385
     *
386
387
388
389
     * @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
390
     */
391
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
392
393
394
395
396
397
398
    /**
     * 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;
399
400
};

401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
/**
 * 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;
    /**
419
420
421
422
423
424
     * 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
425
     */
426
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
427
428
429
430
431
432
433
    /**
     * 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;
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
};

/**
 * 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;
    /**
454
455
456
457
458
459
     * 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
460
     */
461
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
462
463
464
465
466
467
468
    /**
     * 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;
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
};

/**
 * 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;
    /**
489
     * Execute the kernel to calculate the forces and/or energy.
490
     *
491
492
493
494
     * @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
495
     */
496
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
497
498
499
500
501
502
503
    /**
     * 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;
504
505
};

506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
/**
 * 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;
    /**
524
     * Execute the kernel to calculate the forces and/or energy.
525
     *
526
527
528
529
     * @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
530
     */
531
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
532
533
534
535
536
537
538
    /**
     * 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;
539
540
};

541
542
543
544
/**
 * This kernel is invoked by NonbondedForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcNonbondedForceKernel : public KernelImpl {
545
public:
546
547
548
    enum NonbondedMethod {
        NoCutoff = 0,
        CutoffNonPeriodic = 1,
549
        CutoffPeriodic = 2,
550
551
        Ewald = 3,
        PME = 4
552
    };
553
    static std::string Name() {
554
        return "CalcNonbondedForce";
555
    }
556
    CalcNonbondedForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
557
558
    }
    /**
559
     * Initialize the kernel.
560
     * 
561
     * @param system     the System this kernel will be applied to
562
     * @param force      the NonbondedForce this kernel will be used for
563
     */
564
    virtual void initialize(const System& system, const NonbondedForce& force) = 0;
565
    /**
566
567
568
569
570
     * 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
571
572
     * @param includeDirect  true if direct space interactions should be included
     * @param includeReciprocal  true if reciprocal space interactions should be included
573
     * @return the potential energy due to the force
574
     */
575
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy, bool includeDirect, bool includeReciprocal) = 0;
576
577
578
579
580
581
582
    /**
     * 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;
583
584
};

585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
/**
 * 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;
    /**
608
     * Execute the kernel to calculate the forces and/or energy.
609
     *
610
611
612
613
     * @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
614
     */
615
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
616
617
618
619
620
621
622
    /**
     * 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;
623
624
};

625
/**
626
 * This kernel is invoked by GBSAOBCForce to calculate the forces acting on the system and the energy of the system.
627
 */
628
class CalcGBSAOBCForceKernel : public KernelImpl {
629
630
public:
    static std::string Name() {
Mark Friedrichs's avatar
Mark Friedrichs committed
631
        return "CalcGBSAOBCForce";
632
    }
633
    CalcGBSAOBCForceKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
634
635
    }
    /**
636
     * Initialize the kernel.
637
     * 
638
     * @param system     the System this kernel will be applied to
639
     * @param force      the GBSAOBCForce this kernel will be used for
640
     */
641
    virtual void initialize(const System& system, const GBSAOBCForce& force) = 0;
642
    /**
643
644
645
646
647
648
     * 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
649
     */
650
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
651
652
653
654
655
656
657
    /**
     * 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;
658
659
};

Mark Friedrichs's avatar
Mark Friedrichs committed
660
661
662
663
664
665
/**
 * This kernel is invoked by GBVIForce to calculate the forces acting on the system and the energy of the system.
 */
class CalcGBVIForceKernel : public KernelImpl {
public:
    static std::string Name() {
Mark Friedrichs's avatar
Mark Friedrichs committed
666
        return "CalcGBVIForce";
Mark Friedrichs's avatar
Mark Friedrichs committed
667
668
669
670
671
672
673
674
675
676
677
678
    }
    CalcGBVIForceKernel(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 GBVIForce this kernel will be used for
     * @param scaledRadii scaled radii
     */
    virtual void initialize(const System& system, const GBVIForce& force, const std::vector<double>& scaledRadii) = 0;
    /**
679
680
681
682
683
684
     * 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
Mark Friedrichs's avatar
Mark Friedrichs committed
685
     */
686
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
Mark Friedrichs's avatar
Mark Friedrichs committed
687
688
};

689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
/**
 * 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;
    /**
712
     * Execute the kernel to calculate the forces and/or energy.
713
     *
714
715
716
717
     * @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
718
     */
719
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
720
721
722
723
724
725
726
    /**
     * 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;
727
728
};

729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
/**
 * 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;
    /**
747
     * Execute the kernel to calculate the forces and/or energy.
748
     *
749
750
751
752
     * @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
753
     */
754
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
755
756
757
758
759
760
761
    /**
     * 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;
762
763
};

764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
/**
 * 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;
    /**
787
     * Execute the kernel to calculate the forces and/or energy.
788
     *
789
790
791
792
     * @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
793
     */
794
    virtual double execute(ContextImpl& context, bool includeForces, bool includeEnergy) = 0;
795
796
797
798
799
800
801
    /**
     * 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;
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
/**
 * 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;
830
831
832
833
834
835
836
    /**
     * 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;
837
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
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
/**
 * 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;
};

879
880
881
882
883
884
885
886
/**
 * This kernel is invoked by VerletIntegrator to take one time step.
 */
class IntegrateVerletStepKernel : public KernelImpl {
public:
    static std::string Name() {
        return "IntegrateVerletStep";
    }
887
    IntegrateVerletStepKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
888
889
    }
    /**
890
     * Initialize the kernel.
891
     * 
892
893
     * @param system     the System this kernel will be applied to
     * @param integrator the VerletIntegrator this kernel will be used for
894
     */
895
    virtual void initialize(const System& system, const VerletIntegrator& integrator) = 0;
896
897
898
    /**
     * Execute the kernel.
     * 
899
900
     * @param context    the context in which to execute this kernel
     * @param integrator the VerletIntegrator this kernel is being used for
901
     */
902
    virtual void execute(ContextImpl& context, const VerletIntegrator& integrator) = 0;
903
904
905
906
907
908
909
    /**
     * 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;
910
911
912
913
914
915
916
917
918
919
};

/**
 * This kernel is invoked by LangevinIntegrator to take one time step.
 */
class IntegrateLangevinStepKernel : public KernelImpl {
public:
    static std::string Name() {
        return "IntegrateLangevinStep";
    }
920
    IntegrateLangevinStepKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
921
922
    }
    /**
923
     * Initialize the kernel.
924
     * 
925
926
     * @param system     the System this kernel will be applied to
     * @param integrator the LangevinIntegrator this kernel will be used for
927
     */
928
    virtual void initialize(const System& system, const LangevinIntegrator& integrator) = 0;
929
930
931
    /**
     * Execute the kernel.
     * 
932
933
     * @param context    the context in which to execute this kernel
     * @param integrator the LangevinIntegrator this kernel is being used for
934
     */
935
    virtual void execute(ContextImpl& context, const LangevinIntegrator& integrator) = 0;
936
937
938
939
940
941
942
    /**
     * 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;
943
944
945
946
947
948
949
950
951
952
};

/**
 * This kernel is invoked by BrownianIntegrator to take one time step.
 */
class IntegrateBrownianStepKernel : public KernelImpl {
public:
    static std::string Name() {
        return "IntegrateBrownianStep";
    }
953
    IntegrateBrownianStepKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
954
955
    }
    /**
956
     * Initialize the kernel.
957
     * 
958
959
     * @param system     the System this kernel will be applied to
     * @param integrator the BrownianIntegrator this kernel will be used for
960
     */
961
    virtual void initialize(const System& system, const BrownianIntegrator& integrator) = 0;
962
963
964
    /**
     * Execute the kernel.
     * 
965
966
     * @param context    the context in which to execute this kernel
     * @param integrator the BrownianIntegrator this kernel is being used for
967
     */
968
    virtual void execute(ContextImpl& context, const BrownianIntegrator& integrator) = 0;
969
970
971
972
973
974
975
    /**
     * 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;
976
977
};

978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
/**
 * 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
999
     * @param integrator the VariableLangevinIntegrator this kernel is being used for
1000
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1001
     * @return the size of the step that was taken
1002
     */
1003
    virtual double execute(ContextImpl& context, const VariableLangevinIntegrator& integrator, double maxTime) = 0;
1004
1005
1006
1007
1008
1009
1010
    /**
     * 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;
1011
1012
};

1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
/**
 * 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
1034
     * @param integrator the VariableVerletIntegrator this kernel is being used for
1035
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1036
     * @return the size of the step that was taken
1037
     */
1038
    virtual double execute(ContextImpl& context, const VariableVerletIntegrator& integrator, double maxTime) = 0;
1039
1040
1041
1042
1043
1044
1045
    /**
     * 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;
1046
1047
};

1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
/**
 * 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;
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
    /**
     * 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;
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
    /**
     * 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;
};

1119
/**
Peter Eastman's avatar
Peter Eastman committed
1120
 * This kernel is invoked by AndersenThermostat at the start of each time step to adjust the particle velocities.
1121
1122
1123
1124
1125
1126
 */
class ApplyAndersenThermostatKernel : public KernelImpl {
public:
    static std::string Name() {
        return "ApplyAndersenThermostat";
    }
1127
    ApplyAndersenThermostatKernel(std::string name, const Platform& platform) : KernelImpl(name, platform) {
1128
1129
    }
    /**
1130
     * Initialize the kernel.
1131
     * 
1132
1133
     * @param system     the System this kernel will be applied to
     * @param thermostat the AndersenThermostat this kernel will be used for
1134
     */
1135
    virtual void initialize(const System& system, const AndersenThermostat& thermostat) = 0;
1136
1137
1138
    /**
     * Execute the kernel.
     * 
1139
     * @param context    the context in which to execute this kernel
1140
     */
1141
    virtual void execute(ContextImpl& context) = 0;
1142
1143
};

1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
/**
 * 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
     */
1160
    virtual void initialize(const System& system, const Force& barostat) = 0;
Lee-Ping's avatar
Lee-Ping committed
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
    /**
     * 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
1176
1177
1178
1179
1180
1181
1182
     * scaleCoordinates() was last called.
     *
     * @param context    the context in which to execute this kernel
     */
    virtual void restoreCoordinates(ContextImpl& context) = 0;
};

1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
/**
 * 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) {
    }
    /**
1194
     * Initialize the kernel.
1195
     * 
1196
1197
     * @param system     the System this kernel will be applied to
     * @param force      the CMMotionRemover this kernel will be used for
1198
     */
1199
    virtual void initialize(const System& system, const CMMotionRemover& force) = 0;
1200
1201
1202
    /**
     * Execute the kernel.
     * 
1203
     * @param context    the context in which to execute this kernel
1204
     */
1205
    virtual void execute(ContextImpl& context) = 0;
1206
1207
};

1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
/**
 * 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
     */
    virtual void initialize(int gridx, int gridy, int gridz, int numParticles, double alpha) = 0;
    /**
     * Begin computing the force and energy.
peastman's avatar
peastman committed
1235
1236
1237
1238
     *
     * @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
1239
     */
peastman's avatar
peastman committed
1240
    virtual void beginComputation(IO& io, const Vec3* periodicBoxVectors, bool includeEnergy) = 0;
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
1267
1268
1269
1270
1271
    /**
     * 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;
};

/**
 * 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;
};


1272
1273
1274
} // namespace OpenMM

#endif /*OPENMM_KERNELS_H_*/