ReferenceKernels.h 60.8 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
#ifndef OPENMM_REFERENCEKERNELS_H_
#define OPENMM_REFERENCEKERNELS_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
#include "ReferencePlatform.h"
36
#include "openmm/kernels.h"
37
38
#include "SimTKOpenMMRealType.h"
#include "ReferenceNeighborList.h"
39
#include "lepton/CompiledExpression.h"
40
#include "lepton/CustomFunction.h"
41
42
#include <array>
#include <utility>
43

44
45
namespace OpenMM {

46
class ReferenceObc;
47
class ReferenceAndersenThermostat;
48
49
50
51
class ReferenceCustomBondIxn;
class ReferenceCustomAngleIxn;
class ReferenceCustomTorsionIxn;
class ReferenceCustomExternalIxn;
52
class ReferenceCustomCentroidBondIxn;
53
class ReferenceCustomCompoundBondIxn;
54
class ReferenceCustomCVForce;
55
class ReferenceCustomHbondIxn;
56
class ReferenceCustomManyParticleIxn;
57
class ReferenceGayBerneForce;
58
class ReferenceBrownianDynamics;
59
class ReferenceStochasticDynamics;
60
class ReferenceConstraintAlgorithm;
61
class ReferenceMonteCarloBarostat;
62
class ReferenceVariableStochasticDynamics;
63
class ReferenceVariableVerletDynamics;
64
class ReferenceVerletDynamics;
65
class ReferenceCustomDynamics;
66

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

111
/**
112
113
 * This kernel provides methods for setting and retrieving various state data: time, positions,
 * velocities, and forces.
114
 */
115
class ReferenceUpdateStateDataKernel : public UpdateStateDataKernel {
116
public:
117
    ReferenceUpdateStateDataKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : UpdateStateDataKernel(name, platform), data(data) {
118
119
120
121
122
123
124
125
126
127
128
129
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     */
    void initialize(const System& system);
    /**
     * Get the current time (in picoseconds).
     *
     * @param context    the context in which to execute this kernel
     */
130
    double getTime(const ContextImpl& context) const;
131
132
133
134
135
    /**
     * Set the current time (in picoseconds).
     *
     * @param context    the context in which to execute this kernel
     */
136
    void setTime(ContextImpl& context, double time);
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
    /**
     * Get the positions of all particles.
     *
     * @param positions  on exit, this contains the particle positions
     */
    void getPositions(ContextImpl& context, std::vector<Vec3>& positions);
    /**
     * Set the positions of all particles.
     *
     * @param positions  a vector containg the particle positions
     */
    void setPositions(ContextImpl& context, const std::vector<Vec3>& positions);
    /**
     * Get the velocities of all particles.
     *
     * @param velocities  on exit, this contains the particle velocities
     */
    void getVelocities(ContextImpl& context, std::vector<Vec3>& velocities);
    /**
     * Set the velocities of all particles.
     *
     * @param velocities  a vector containg the particle velocities
     */
    void setVelocities(ContextImpl& context, const std::vector<Vec3>& velocities);
    /**
     * Get the current forces on all particles.
     *
     * @param forces  on exit, this contains the forces
     */
    void getForces(ContextImpl& context, std::vector<Vec3>& forces);
167
168
169
170
171
172
    /**
     * Get the current derivatives of the energy with respect to context parameters.
     *
     * @param derivs  on exit, this contains the derivatives
     */
    void getEnergyParameterDerivatives(ContextImpl& context, std::map<std::string, double>& derivs);
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
    /**
     * 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
     */
    void getPeriodicBoxVectors(ContextImpl& context, Vec3& a, Vec3& b, Vec3& c) const;
    /**
     * 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
     */
188
    void setPeriodicBoxVectors(ContextImpl& context, const Vec3& a, const Vec3& b, const Vec3& c);
Peter Eastman's avatar
Peter Eastman committed
189
190
191
192
193
194
195
196
197
198
199
200
    /**
     * Create a checkpoint recording the current state of the Context.
     * 
     * @param stream    an output stream the checkpoint data should be written to
     */
    void createCheckpoint(ContextImpl& context, std::ostream& stream);
    /**
     * Load a checkpoint that was written by createCheckpoint().
     * 
     * @param stream    an input stream the checkpoint data should be read from
     */
    void loadCheckpoint(ContextImpl& context, std::istream& stream);
201
202
203
204
private:
    ReferencePlatform::PlatformData& data;
};

205
206
207
208
209
210
/**
 * This kernel modifies the positions of particles to enforce distance constraints.
 */
class ReferenceApplyConstraintsKernel : public ApplyConstraintsKernel {
public:
    ReferenceApplyConstraintsKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) :
211
            ApplyConstraintsKernel(name, platform), data(data) {
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
    }
    ~ReferenceApplyConstraintsKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     */
    void initialize(const System& system);
    /**
     * 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.
     */
    void apply(ContextImpl& context, double tol);
227
228
229
230
231
232
233
    /**
     * 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.
     */
    void applyToVelocities(ContextImpl& context, double tol);
234
235
private:
    ReferencePlatform::PlatformData& data;
peastman's avatar
peastman committed
236
237
    std::vector<double> masses;
    std::vector<double> inverseMasses;
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 ReferenceVirtualSitesKernel : public VirtualSitesKernel {
public:
    ReferenceVirtualSitesKernel(std::string name, const Platform& platform) : VirtualSitesKernel(name, platform) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     */
    void initialize(const System& system);
    /**
     * Compute the virtual site locations.
     *
     * @param context    the context in which to execute this kernel
     */
    void computePositions(ContextImpl& context);
};

261
/**
262
 * This kernel is invoked by HarmonicBondForce to calculate the forces acting on the system and the energy of the system.
263
 */
264
class ReferenceCalcHarmonicBondForceKernel : public CalcHarmonicBondForceKernel {
265
public:
266
    ReferenceCalcHarmonicBondForceKernel(std::string name, const Platform& platform) : CalcHarmonicBondForceKernel(name, platform) {
267
268
    }
    /**
269
     * Initialize the kernel.
270
     * 
271
     * @param system     the System this kernel will be applied to
272
273
274
275
     * @param force      the HarmonicBondForce this kernel will be used for
     */
    void initialize(const System& system, const HarmonicBondForce& force);
    /**
276
277
278
279
280
281
     * 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
282
     */
283
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
284
285
286
287
288
289
290
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the HarmonicBondForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const HarmonicBondForce& force);
291
292
private:
    int numBonds;
293
294
    std::vector<std::vector<int> >bondIndexArray;
    std::vector<std::vector<double> >bondParamArray;
295
    bool usePeriodic;
296
297
};

298
299
300
301
302
/**
 * This kernel is invoked by CustomBondForce to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcCustomBondForceKernel : public CalcCustomBondForceKernel {
public:
303
    ReferenceCalcCustomBondForceKernel(std::string name, const Platform& platform) : CalcCustomBondForceKernel(name, platform), ixn(NULL) {
304
    }
305
    ~ReferenceCalcCustomBondForceKernel();
306
307
308
309
310
311
312
313
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomBondForce this kernel will be used for
     */
    void initialize(const System& system, const CustomBondForce& force);
    /**
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
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
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
     */
    void copyParametersToContext(ContextImpl& context, const CustomBondForce& force);
329
330
private:
    int numBonds;
331
    ReferenceCustomBondIxn* ixn;
332
333
    std::vector<std::vector<int> >bondIndexArray;
    std::vector<std::vector<double> >bondParamArray;
334
    Lepton::CompiledExpression energyExpression, forceExpression;
335
336
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
337
    bool usePeriodic;
338
339
};

340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
/**
 * This kernel is invoked by HarmonicAngleForce to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcHarmonicAngleForceKernel : public CalcHarmonicAngleForceKernel {
public:
    ReferenceCalcHarmonicAngleForceKernel(std::string name, const Platform& platform) : CalcHarmonicAngleForceKernel(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
     */
    void initialize(const System& system, const HarmonicAngleForce& force);
    /**
355
356
357
358
359
360
     * 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
361
     */
362
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
363
364
365
366
367
368
369
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the HarmonicAngleForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const HarmonicAngleForce& force);
370
371
private:
    int numAngles;
372
373
    std::vector<std::vector<int> >angleIndexArray;
    std::vector<std::vector<double> >angleParamArray;
374
    bool usePeriodic;
375
376
};

377
378
379
380
381
/**
 * This kernel is invoked by CustomAngleForce to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcCustomAngleForceKernel : public CalcCustomAngleForceKernel {
public:
382
    ReferenceCalcCustomAngleForceKernel(std::string name, const Platform& platform) : CalcCustomAngleForceKernel(name, platform), ixn(NULL) {
383
    }
384
    ~ReferenceCalcCustomAngleForceKernel();
385
386
387
388
389
390
391
392
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomAngleForce this kernel will be used for
     */
    void initialize(const System& system, const CustomAngleForce& force);
    /**
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
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
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
     */
    void copyParametersToContext(ContextImpl& context, const CustomAngleForce& force);
408
409
private:
    int numAngles;
410
    ReferenceCustomAngleIxn* ixn;
411
412
    std::vector<std::vector<int> >angleIndexArray;
    std::vector<std::vector<double> >angleParamArray;
413
    Lepton::CompiledExpression energyExpression, forceExpression;
414
415
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
416
    bool usePeriodic;
417
418
};

419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
/**
 * This kernel is invoked by PeriodicTorsionForce to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcPeriodicTorsionForceKernel : public CalcPeriodicTorsionForceKernel {
public:
    ReferenceCalcPeriodicTorsionForceKernel(std::string name, const Platform& platform) : CalcPeriodicTorsionForceKernel(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
     */
    void initialize(const System& system, const PeriodicTorsionForce& force);
    /**
434
435
436
437
438
439
     * 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
440
     */
441
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
442
443
444
445
446
447
448
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the PeriodicTorsionForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const PeriodicTorsionForce& force);
449
450
private:
    int numTorsions;
451
452
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
453
    bool usePeriodic;
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
};

/**
 * This kernel is invoked by RBTorsionForce to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcRBTorsionForceKernel : public CalcRBTorsionForceKernel {
public:
    ReferenceCalcRBTorsionForceKernel(std::string name, const Platform& platform) : CalcRBTorsionForceKernel(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
     */
    void initialize(const System& system, const RBTorsionForce& force);
    /**
471
472
473
474
475
476
     * 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
477
     */
478
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
479
480
481
482
483
484
485
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the RBTorsionForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const RBTorsionForce& force);
486
487
private:
    int numTorsions;
488
489
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
490
    bool usePeriodic;
491
492
};

493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
/**
 * This kernel is invoked by CMAPTorsionForce to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcCMAPTorsionForceKernel : public CalcCMAPTorsionForceKernel {
public:
    ReferenceCalcCMAPTorsionForceKernel(std::string name, const Platform& platform) : CalcCMAPTorsionForceKernel(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
     */
    void initialize(const System& system, const CMAPTorsionForce& force);
    /**
508
     * Execute the kernel to calculate the forces and/or energy.
509
     *
510
511
512
513
     * @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
514
     */
515
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
516
517
518
519
520
521
522
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CMAPTorsionForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CMAPTorsionForce& force);
523
private:
peastman's avatar
peastman committed
524
    std::vector<std::vector<std::vector<double> > > coeff;
525
526
    std::vector<int> torsionMaps;
    std::vector<std::vector<int> > torsionIndices;
527
    bool usePeriodic;
528
529
};

530
531
532
533
534
/**
 * This kernel is invoked by CustomTorsionForce to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcCustomTorsionForceKernel : public CalcCustomTorsionForceKernel {
public:
535
    ReferenceCalcCustomTorsionForceKernel(std::string name, const Platform& platform) : CalcCustomTorsionForceKernel(name, platform), ixn(NULL) {
536
    }
537
    ~ReferenceCalcCustomTorsionForceKernel();
538
539
540
541
542
543
544
545
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomTorsionForce this kernel will be used for
     */
    void initialize(const System& system, const CustomTorsionForce& force);
    /**
546
     * Execute the kernel to calculate the forces and/or energy.
547
     *
548
549
550
551
     * @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
552
     */
553
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
554
555
556
557
558
559
560
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomTorsionForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CustomTorsionForce& force);
561
562
private:
    int numTorsions;
563
    ReferenceCustomTorsionIxn* ixn;
564
565
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
566
    Lepton::CompiledExpression energyExpression, forceExpression;
567
568
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
569
    bool usePeriodic;
570
571
};

572
573
574
575
576
577
578
579
580
581
582
583
584
/**
 * This kernel is invoked by NonbondedForce to calculate the forces acting on the system.
 */
class ReferenceCalcNonbondedForceKernel : public CalcNonbondedForceKernel {
public:
    ReferenceCalcNonbondedForceKernel(std::string name, const Platform& platform) : CalcNonbondedForceKernel(name, platform) {
    }
    ~ReferenceCalcNonbondedForceKernel();
    /**
     * Initialize the kernel.
     * 
     * @param system     the System this kernel will be applied to
     * @param force      the NonbondedForce this kernel will be used for
585
     */
586
    void initialize(const System& system, const NonbondedForce& force);
587
    /**
588
589
590
591
592
     * 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
593
     * @param includeReciprocal  true if reciprocal space interactions should be included
594
     * @return the potential energy due to the force
595
     */
596
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy, bool includeDirect, bool includeReciprocal);
597
598
599
600
601
602
603
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the NonbondedForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const NonbondedForce& force);
604
605
606
607
608
609
610
611
612
    /**
     * Get the parameters being used for PME.
     * 
     * @param alpha   the separation parameter
     * @param nx      the number of grid points along the X axis
     * @param ny      the number of grid points along the Y axis
     * @param nz      the number of grid points along the Z axis
     */
    void getPMEParameters(double& alpha, int& nx, int& ny, int& nz) const;
613
614
615
    /**
     * Get the dispersion parameters being used for the dispersion term in LJPME.
     *
616
617
618
619
     * @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
620
     */
621
    void getLJPMEParameters(double& alpha, int& nx, int& ny, int& nz) const;
622
private:
623
    void computeParameters(ContextImpl& context);
Peter Eastman's avatar
Peter Eastman committed
624
    int numParticles, num14;
625
626
    std::vector<std::vector<int> >bonded14IndexArray;
    std::vector<std::vector<double> > particleParamArray, bonded14ParamArray;
627
628
    std::vector<std::array<double, 3> > baseParticleParams, baseExceptionParams;
    std::map<std::pair<std::string, int>, std::array<double, 3> > particleParamOffsets, exceptionParamOffsets;
peastman's avatar
peastman committed
629
    double nonbondedCutoff, switchingDistance, rfDielectric, ewaldAlpha, ewaldDispersionAlpha, dispersionCoefficient;
630
    int kmax[3], gridSize[3], dispersionGridSize[3];
631
    bool useSwitchingFunction;
632
633
634
    std::vector<std::set<int> > exclusions;
    NonbondedMethod nonbondedMethod;
    NeighborList* neighborList;
635
636
};

637
638
639
640
641
/**
 * This kernel is invoked by CustomNonbondedForce to calculate the forces acting on the system.
 */
class ReferenceCalcCustomNonbondedForceKernel : public CalcCustomNonbondedForceKernel {
public:
642
    ReferenceCalcCustomNonbondedForceKernel(std::string name, const Platform& platform) : CalcCustomNonbondedForceKernel(name, platform), forceCopy(NULL) {
643
644
645
646
647
648
649
650
651
652
    }
    ~ReferenceCalcCustomNonbondedForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomNonbondedForce this kernel will be used for
     */
    void initialize(const System& system, const CustomNonbondedForce& force);
    /**
653
     * Execute the kernel to calculate the forces and/or energy.
654
     *
655
656
657
658
     * @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
659
     */
660
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
661
662
663
664
665
666
667
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomNonbondedForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CustomNonbondedForce& force);
668
private:
669
    int numParticles;
670
    std::vector<std::vector<double> > particleParamArray;
peastman's avatar
peastman committed
671
    double nonbondedCutoff, switchingDistance, periodicBoxSize[3], longRangeCoefficient;
672
    bool useSwitchingFunction, hasInitializedLongRangeCorrection;
673
674
    CustomNonbondedForce* forceCopy;
    std::map<std::string, double> globalParamValues;
675
    std::vector<std::set<int> > exclusions;
676
    Lepton::CompiledExpression energyExpression, forceExpression;
677
678
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
679
    std::vector<std::pair<std::set<int>, std::set<int> > > interactionGroups;
680
    std::vector<double> longRangeCoefficientDerivs;
681
682
683
684
    NonbondedMethod nonbondedMethod;
    NeighborList* neighborList;
};

685
/**
686
 * This kernel is invoked by GBSAOBCForce to calculate the forces acting on the system.
687
 */
688
class ReferenceCalcGBSAOBCForceKernel : public CalcGBSAOBCForceKernel {
689
public:
690
    ReferenceCalcGBSAOBCForceKernel(std::string name, const Platform& platform) : CalcGBSAOBCForceKernel(name, platform) {
691
    }
692
    ~ReferenceCalcGBSAOBCForceKernel();
693
    /**
694
     * Initialize the kernel.
695
     * 
696
     * @param system     the System this kernel will be applied to
697
     * @param force      the GBSAOBCForce this kernel will be used for
698
     */
699
    void initialize(const System& system, const GBSAOBCForce& force);
700
    /**
701
702
703
704
705
706
     * 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
707
     */
708
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
709
710
711
712
713
714
715
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the GBSAOBCForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const GBSAOBCForce& force);
716
private:
717
    ReferenceObc* obc;
peastman's avatar
peastman committed
718
    std::vector<double> charges;
719
    bool isPeriodic;
720
721
};

722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
/**
 * This kernel is invoked by CustomGBForce to calculate the forces acting on the system.
 */
class ReferenceCalcCustomGBForceKernel : public CalcCustomGBForceKernel {
public:
    ReferenceCalcCustomGBForceKernel(std::string name, const Platform& platform) : CalcCustomGBForceKernel(name, platform) {
    }
    ~ReferenceCalcCustomGBForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomGBForce this kernel will be used for
     */
    void initialize(const System& system, const CustomGBForce& force);
    /**
738
     * Execute the kernel to calculate the forces and/or energy.
739
     *
740
741
742
743
     * @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
744
     */
745
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
746
747
748
749
750
751
752
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomGBForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CustomGBForce& force);
753
754
private:
    int numParticles;
755
    bool isPeriodic;
756
    std::vector<std::vector<double> > particleParamArray;
peastman's avatar
peastman committed
757
    double nonbondedCutoff;
758
    std::vector<std::set<int> > exclusions;
759
    std::vector<std::string> particleParameterNames, globalParameterNames, energyParamDerivNames, valueNames;
760
761
762
    std::vector<Lepton::CompiledExpression> valueExpressions;
    std::vector<std::vector<Lepton::CompiledExpression> > valueDerivExpressions;
    std::vector<std::vector<Lepton::CompiledExpression> > valueGradientExpressions;
763
    std::vector<std::vector<Lepton::CompiledExpression> > valueParamDerivExpressions;
764
    std::vector<OpenMM::CustomGBForce::ComputationType> valueTypes;
765
766
767
    std::vector<Lepton::CompiledExpression> energyExpressions;
    std::vector<std::vector<Lepton::CompiledExpression> > energyDerivExpressions;
    std::vector<std::vector<Lepton::CompiledExpression> > energyGradientExpressions;
768
    std::vector<std::vector<Lepton::CompiledExpression> > energyParamDerivExpressions;
769
770
771
772
773
    std::vector<OpenMM::CustomGBForce::ComputationType> energyTypes;
    NonbondedMethod nonbondedMethod;
    NeighborList* neighborList;
};

774
775
776
777
778
/**
 * This kernel is invoked by CustomExternalForce to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcCustomExternalForceKernel : public CalcCustomExternalForceKernel {
public:
779
    ReferenceCalcCustomExternalForceKernel(std::string name, const Platform& platform) : CalcCustomExternalForceKernel(name, platform), ixn(NULL) {
780
    }
781
    ~ReferenceCalcCustomExternalForceKernel();
782
783
784
785
786
787
788
789
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomExternalForce this kernel will be used for
     */
    void initialize(const System& system, const CustomExternalForce& force);
    /**
790
     * Execute the kernel to calculate the forces and/or energy.
791
     *
792
793
794
795
     * @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
796
     */
797
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
798
799
800
801
802
803
804
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomExternalForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CustomExternalForce& force);
805
private:
806
    class PeriodicDistanceFunction;
807
    int numParticles;
808
    ReferenceCustomExternalIxn* ixn;
809
    std::vector<int> particles;
810
    std::vector<std::vector<double> > particleParamArray;
811
    Lepton::CompiledExpression energyExpression, forceExpressionX, forceExpressionY, forceExpressionZ;
812
    std::vector<std::string> parameterNames, globalParameterNames;
peastman's avatar
peastman committed
813
    Vec3* boxVectors;
814
815
816
817
};

class ReferenceCalcCustomExternalForceKernel::PeriodicDistanceFunction : public Lepton::CustomFunction {
public:
peastman's avatar
peastman committed
818
819
    Vec3** boxVectorHandle;
    PeriodicDistanceFunction(Vec3** boxVectorHandle);
820
821
822
823
    int getNumArguments() const;
    double evaluate(const double* arguments) const;
    double evaluateDerivative(const double* arguments, const int* derivOrder) const;
    Lepton::CustomFunction* clone() const;
824
825
};

826
827
828
829
830
/**
 * This kernel is invoked by CustomHbondForce to calculate the forces acting on the system.
 */
class ReferenceCalcCustomHbondForceKernel : public CalcCustomHbondForceKernel {
public:
831
    ReferenceCalcCustomHbondForceKernel(std::string name, const Platform& platform) : CalcCustomHbondForceKernel(name, platform), ixn(NULL) {
832
833
834
835
836
837
838
839
840
841
    }
    ~ReferenceCalcCustomHbondForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomHbondForce this kernel will be used for
     */
    void initialize(const System& system, const CustomHbondForce& force);
    /**
842
     * Execute the kernel to calculate the forces and/or energy.
843
     *
844
845
846
847
     * @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
848
     */
849
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
850
851
852
853
854
855
856
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomHbondForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CustomHbondForce& force);
857
858
private:
    int numDonors, numAcceptors, numParticles;
859
    bool isPeriodic;
860
    std::vector<std::vector<double> > donorParamArray, acceptorParamArray;
peastman's avatar
peastman committed
861
    double nonbondedCutoff;
862
    ReferenceCustomHbondIxn* ixn;
863
    std::vector<std::set<int> > exclusions;
864
    std::vector<std::string> globalParameterNames;
865
866
};

867
868
869
870
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
/**
 * This kernel is invoked by CustomCentroidBondForce to calculate the forces acting on the system.
 */
class ReferenceCalcCustomCentroidBondForceKernel : public CalcCustomCentroidBondForceKernel {
public:
    ReferenceCalcCustomCentroidBondForceKernel(std::string name, const Platform& platform) : CalcCustomCentroidBondForceKernel(name, platform), ixn(NULL) {
    }
    ~ReferenceCalcCustomCentroidBondForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomCentroidBondForce this kernel will be used for
     */
    void initialize(const System& system, const CustomCentroidBondForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomCentroidBondForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CustomCentroidBondForce& force);
private:
    int numBonds, numParticles;
900
    std::vector<std::vector<double> > bondParamArray;
901
    ReferenceCustomCentroidBondIxn* ixn;
902
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
903
    bool usePeriodic;
904
905
};

906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
/**
 * This kernel is invoked by CustomCompoundBondForce to calculate the forces acting on the system.
 */
class ReferenceCalcCustomCompoundBondForceKernel : public CalcCustomCompoundBondForceKernel {
public:
    ReferenceCalcCustomCompoundBondForceKernel(std::string name, const Platform& platform) : CalcCustomCompoundBondForceKernel(name, platform), ixn(NULL) {
    }
    ~ReferenceCalcCustomCompoundBondForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomCompoundBondForce this kernel will be used for
     */
    void initialize(const System& system, const CustomCompoundBondForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
930
931
932
933
934
935
936
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomCompoundBondForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CustomCompoundBondForce& force);
937
private:
938
    int numBonds;
939
    std::vector<std::vector<double> > bondParamArray;
940
    ReferenceCustomCompoundBondIxn* ixn;
941
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
942
    bool usePeriodic;
943
944
};

945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
/**
 * This kernel is invoked by CustomManyParticleForce to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcCustomManyParticleForceKernel : public CalcCustomManyParticleForceKernel {
public:
    ReferenceCalcCustomManyParticleForceKernel(std::string name, const Platform& platform) : CalcCustomManyParticleForceKernel(name, platform), ixn(NULL) {
    }
    ~ReferenceCalcCustomManyParticleForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomManyParticleForce this kernel will be used for
     */
    void initialize(const System& system, const CustomManyParticleForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomManyParticleForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CustomManyParticleForce& force);
private:
    int numParticles;
peastman's avatar
peastman committed
978
    double cutoffDistance;
979
    std::vector<std::vector<double> > particleParamArray;
980
981
982
983
984
    ReferenceCustomManyParticleIxn* ixn;
    std::vector<std::string> globalParameterNames;
    NonbondedMethod nonbondedMethod;
};

985
986
987
988
989
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
/**
 * This kernel is invoked by GayBerneForce to calculate the forces acting on the system.
 */
class ReferenceCalcGayBerneForceKernel : public CalcGayBerneForceKernel {
public:
    ReferenceCalcGayBerneForceKernel(std::string name, const Platform& platform) : CalcGayBerneForceKernel(name, platform), ixn(NULL) {
    }
    ~ReferenceCalcGayBerneForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the GayBerneForce this kernel will be used for
     */
    void initialize(const System& system, const GayBerneForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the GayBerneForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const GayBerneForce& force);
private:
    ReferenceGayBerneForce* ixn;
};

1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
/**
 * This kernel is invoked by CustomCVForce to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcCustomCVForceKernel : public CalcCustomCVForceKernel {
public:
    ReferenceCalcCustomCVForceKernel(std::string name, const Platform& platform) : CalcCustomCVForceKernel(name, platform), ixn(NULL) {
    }
    ~ReferenceCalcCustomCVForceKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomCVForce this kernel will be used for
1032
     * @param innerContext   the context created by the CustomCVForce for computing collective variables
1033
     */
1034
    void initialize(const System& system, const CustomCVForce& force, ContextImpl& innerContext);
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
    /**
     * 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
     */
    double execute(ContextImpl& context, ContextImpl& innerContext, bool includeForces, bool includeEnergy);
    /**
     * 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
     */
    void copyState(ContextImpl& context, ContextImpl& innerContext);
1052
1053
1054
1055
1056
1057
1058
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomCVForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CustomCVForce& force);
1059
1060
1061
1062
1063
private:
    ReferenceCustomCVForce* ixn;
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
};

1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
/**
 * This kernel is invoked by RMSDForce to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcRMSDForceKernel : public CalcRMSDForceKernel {
public:
    ReferenceCalcRMSDForceKernel(std::string name, const Platform& platform) : CalcRMSDForceKernel(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
     */
    void initialize(const System& system, const RMSDForce& force);
    /**
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
     */
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the RMSDForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const RMSDForce& force);
private:
    std::vector<Vec3> referencePos;
    std::vector<int> particles;
};

1099
1100
1101
1102
1103
/**
 * This kernel is invoked by VerletIntegrator to take one time step.
 */
class ReferenceIntegrateVerletStepKernel : public IntegrateVerletStepKernel {
public:
1104
    ReferenceIntegrateVerletStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateVerletStepKernel(name, platform),
1105
        data(data), dynamics(0) {
1106
    }
1107
    ~ReferenceIntegrateVerletStepKernel();
1108
    /**
1109
     * Initialize the kernel.
1110
     * 
1111
1112
     * @param system     the System this kernel will be applied to
     * @param integrator the VerletIntegrator this kernel will be used for
1113
     */
1114
    void initialize(const System& system, const VerletIntegrator& integrator);
1115
1116
1117
    /**
     * Execute the kernel.
     * 
1118
1119
     * @param context    the context in which to execute this kernel
     * @param integrator the VerletIntegrator this kernel is being used for
1120
     */
1121
    void execute(ContextImpl& context, const VerletIntegrator& integrator);
1122
1123
1124
1125
1126
1127
1128
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the VerletIntegrator this kernel is being used for
     */
    double computeKineticEnergy(ContextImpl& context, const VerletIntegrator& integrator);
1129
private:
1130
    ReferencePlatform::PlatformData& data;
1131
    ReferenceVerletDynamics* dynamics;
peastman's avatar
peastman committed
1132
    std::vector<double> masses;
1133
    double prevStepSize;
1134
1135
1136
1137
1138
1139
1140
};

/**
 * This kernel is invoked by LangevinIntegrator to take one time step.
 */
class ReferenceIntegrateLangevinStepKernel : public IntegrateLangevinStepKernel {
public:
1141
    ReferenceIntegrateLangevinStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateLangevinStepKernel(name, platform),
1142
        data(data), dynamics(0) {
1143
    }
1144
    ~ReferenceIntegrateLangevinStepKernel();
1145
    /**
Peter Eastman's avatar
Peter Eastman committed
1146
     * Initialize the kernel, setting up the particle masses.
1147
     * 
1148
1149
     * @param system     the System this kernel will be applied to
     * @param integrator the LangevinIntegrator this kernel will be used for
1150
     */
1151
    void initialize(const System& system, const LangevinIntegrator& integrator);
1152
1153
1154
    /**
     * Execute the kernel.
     * 
1155
1156
     * @param context    the context in which to execute this kernel
     * @param integrator the LangevinIntegrator this kernel is being used for
1157
     */
1158
    void execute(ContextImpl& context, const LangevinIntegrator& integrator);
1159
1160
1161
1162
1163
1164
1165
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the LangevinIntegrator this kernel is being used for
     */
    double computeKineticEnergy(ContextImpl& context, const LangevinIntegrator& integrator);
1166
private:
1167
    ReferencePlatform::PlatformData& data;
1168
    ReferenceStochasticDynamics* dynamics;
peastman's avatar
peastman committed
1169
    std::vector<double> masses;
1170
    double prevTemp, prevFriction, prevStepSize;
1171
1172
1173
1174
1175
1176
1177
};

/**
 * This kernel is invoked by BrownianIntegrator to take one time step.
 */
class ReferenceIntegrateBrownianStepKernel : public IntegrateBrownianStepKernel {
public:
1178
    ReferenceIntegrateBrownianStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateBrownianStepKernel(name, platform),
1179
        data(data), dynamics(0) {
1180
    }
1181
    ~ReferenceIntegrateBrownianStepKernel();
1182
    /**
1183
     * Initialize the kernel.
1184
     * 
1185
1186
     * @param system     the System this kernel will be applied to
     * @param integrator the BrownianIntegrator this kernel will be used for
1187
     */
1188
    void initialize(const System& system, const BrownianIntegrator& integrator);
1189
1190
1191
    /**
     * Execute the kernel.
     * 
1192
1193
     * @param context    the context in which to execute this kernel
     * @param integrator the BrownianIntegrator this kernel is being used for
1194
     */
1195
    void execute(ContextImpl& context, const BrownianIntegrator& integrator);
1196
1197
1198
1199
1200
1201
1202
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the BrownianIntegrator this kernel is being used for
     */
    double computeKineticEnergy(ContextImpl& context, const BrownianIntegrator& integrator);
1203
private:
1204
    ReferencePlatform::PlatformData& data;
1205
    ReferenceBrownianDynamics* dynamics;
peastman's avatar
peastman committed
1206
    std::vector<double> masses;
1207
    double prevTemp, prevFriction, prevStepSize;
1208
1209
};

1210
1211
1212
1213
1214
1215
/**
 * This kernel is invoked by VariableLangevinIntegrator to take one time step.
 */
class ReferenceIntegrateVariableLangevinStepKernel : public IntegrateVariableLangevinStepKernel {
public:
    ReferenceIntegrateVariableLangevinStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateVariableLangevinStepKernel(name, platform),
1216
        data(data), dynamics(0) {
1217
1218
1219
1220
1221
1222
    }
    ~ReferenceIntegrateVariableLangevinStepKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
1223
     * @param integrator the VariableLangevinIntegrator this kernel will be used for
1224
1225
1226
1227
1228
1229
     */
    void initialize(const System& system, const VariableLangevinIntegrator& integrator);
    /**
     * Execute the kernel.
     *
     * @param context    the context in which to execute this kernel
1230
     * @param integrator the VariableLangevinIntegrator this kernel is being used for
1231
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1232
     * @return the size of the step that was taken
1233
     */
1234
    double execute(ContextImpl& context, const VariableLangevinIntegrator& integrator, double maxTime);
1235
1236
1237
1238
1239
1240
1241
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the VariableLangevinIntegrator this kernel is being used for
     */
    double computeKineticEnergy(ContextImpl& context, const VariableLangevinIntegrator& integrator);
1242
1243
1244
private:
    ReferencePlatform::PlatformData& data;
    ReferenceVariableStochasticDynamics* dynamics;
peastman's avatar
peastman committed
1245
    std::vector<double> masses;
1246
1247
1248
    double prevTemp, prevFriction, prevErrorTol;
};

1249
1250
1251
1252
1253
1254
/**
 * This kernel is invoked by VariableVerletIntegrator to take one time step.
 */
class ReferenceIntegrateVariableVerletStepKernel : public IntegrateVariableVerletStepKernel {
public:
    ReferenceIntegrateVariableVerletStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateVariableVerletStepKernel(name, platform),
1255
        data(data), dynamics(0) {
1256
1257
1258
1259
1260
1261
    }
    ~ReferenceIntegrateVariableVerletStepKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
1262
     * @param integrator the VariableVerletIntegrator this kernel will be used for
1263
1264
1265
1266
1267
1268
     */
    void initialize(const System& system, const VariableVerletIntegrator& integrator);
    /**
     * Execute the kernel.
     *
     * @param context    the context in which to execute this kernel
1269
     * @param integrator the VariableVerletIntegrator this kernel is being used for
1270
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1271
     * @return the size of the step that was taken
1272
     */
1273
    double execute(ContextImpl& context, const VariableVerletIntegrator& integrator, double maxTime);
1274
1275
1276
1277
1278
1279
1280
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the VariableVerletIntegrator this kernel is being used for
     */
    double computeKineticEnergy(ContextImpl& context, const VariableVerletIntegrator& integrator);
1281
1282
1283
private:
    ReferencePlatform::PlatformData& data;
    ReferenceVariableVerletDynamics* dynamics;
peastman's avatar
peastman committed
1284
    std::vector<double> masses;
1285
    double prevErrorTol;
1286
1287
};

1288
1289
1290
1291
1292
1293
/**
 * This kernel is invoked by CustomIntegrator to take one time step.
 */
class ReferenceIntegrateCustomStepKernel : public IntegrateCustomStepKernel {
public:
    ReferenceIntegrateCustomStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateCustomStepKernel(name, platform),
1294
        data(data), dynamics(0) {
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
    }
    ~ReferenceIntegrateCustomStepKernel();
    /**
     * Initialize the kernel.
     * 
     * @param system     the System this kernel will be applied to
     * @param integrator the CustomIntegrator this kernel will be used for
     */
    void initialize(const System& system, const CustomIntegrator& integrator);
    /**
     * 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.
     */
    void execute(ContextImpl& context, CustomIntegrator& integrator, bool& forcesAreValid);
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
    /**
     * 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.
     */
    double computeKineticEnergy(ContextImpl& context, CustomIntegrator& integrator, bool& forcesAreValid);
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
    /**
     * 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
     */
    void getGlobalVariables(ContextImpl& context, std::vector<double>& values) const;
    /**
     * Set the values of all global variables.
     *
     * @param context   the context in which to execute this kernel
     * @param values    a vector containing the values
     */
    void setGlobalVariables(ContextImpl& context, const std::vector<double>& values);
    /**
     * 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
     */
    void getPerDofVariable(ContextImpl& context, int variable, std::vector<Vec3>& values) const;
    /**
     * 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
     */
    void setPerDofVariable(ContextImpl& context, int variable, const std::vector<Vec3>& values);
private:
    ReferencePlatform::PlatformData& data;
    ReferenceCustomDynamics* dynamics;
peastman's avatar
peastman committed
1359
1360
    std::vector<double> masses, globalValues;
    std::vector<std::vector<OpenMM::Vec3> > perDofValues; 
1361
1362
};

1363
/**
Peter Eastman's avatar
Peter Eastman committed
1364
 * This kernel is invoked by AndersenThermostat at the start of each time step to adjust the particle velocities.
1365
1366
1367
 */
class ReferenceApplyAndersenThermostatKernel : public ApplyAndersenThermostatKernel {
public:
1368
    ReferenceApplyAndersenThermostatKernel(std::string name, const Platform& platform) : ApplyAndersenThermostatKernel(name, platform), thermostat(0) {
1369
    }
1370
    ~ReferenceApplyAndersenThermostatKernel();
1371
    /**
1372
     * Initialize the kernel.
1373
     * 
1374
1375
     * @param system     the System this kernel will be applied to
     * @param thermostat the AndersenThermostat this kernel will be used for
1376
     */
1377
    void initialize(const System& system, const AndersenThermostat& thermostat);
1378
1379
1380
    /**
     * Execute the kernel.
     * 
1381
     * @param context    the context in which to execute this kernel
1382
     */
1383
    void execute(ContextImpl& context);
1384
1385
private:
    ReferenceAndersenThermostat* thermostat;
1386
    std::vector<std::vector<int> > particleGroups;
peastman's avatar
peastman committed
1387
    std::vector<double> masses;
1388
1389
};

1390
1391
1392
1393
1394
/**
 * This kernel is invoked by MonteCarloBarostat to adjust the periodic box volume
 */
class ReferenceApplyMonteCarloBarostatKernel : public ApplyMonteCarloBarostatKernel {
public:
1395
    ReferenceApplyMonteCarloBarostatKernel(std::string name, const Platform& platform) : ApplyMonteCarloBarostatKernel(name, platform), barostat(NULL) {
1396
1397
1398
1399
1400
1401
1402
1403
    }
    ~ReferenceApplyMonteCarloBarostatKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param barostat   the MonteCarloBarostat this kernel will be used for
     */
1404
    void initialize(const System& system, const Force& barostat);
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
    /**
     * 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
     */
1417
    void scaleCoordinates(ContextImpl& context, double scaleX, double scaleY, double scaleZ);
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
    /**
     * Reject the most recent Monte Carlo step, restoring the particle positions to where they were before
     * scaleCoordinates() was last called.
     *
     * @param context    the context in which to execute this kernel
     */
    void restoreCoordinates(ContextImpl& context);
private:
    ReferenceMonteCarloBarostat* barostat;
};

1429
1430
1431
1432
1433
/**
 * This kernel is invoked to remove center of mass motion from the system.
 */
class ReferenceRemoveCMMotionKernel : public RemoveCMMotionKernel {
public:
1434
    ReferenceRemoveCMMotionKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : RemoveCMMotionKernel(name, platform), data(data) {
1435
1436
    }
    /**
Peter Eastman's avatar
Peter Eastman committed
1437
     * Initialize the kernel, setting up the particle masses.
1438
     * 
1439
1440
     * @param system     the System this kernel will be applied to
     * @param force      the CMMotionRemover this kernel will be used for
1441
     */
1442
    void initialize(const System& system, const CMMotionRemover& force);
1443
1444
1445
    /**
     * Execute the kernel.
     * 
1446
     * @param context    the context in which to execute this kernel
1447
     */
1448
    void execute(ContextImpl& context);
1449
private:
1450
    ReferencePlatform::PlatformData& data;
1451
    std::vector<double> masses;
1452
    int frequency;
1453
1454
};

1455
1456
1457
} // namespace OpenMM

#endif /*OPENMM_REFERENCEKERNELS_H_*/