ReferenceKernels.h 60 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
class ReferenceCustomCentroidBondIxn;
49
class ReferenceCustomCompoundBondIxn;
50
class ReferenceCustomCVForce;
51
class ReferenceCustomHbondIxn;
52
class ReferenceCustomManyParticleIxn;
53
class ReferenceGayBerneForce;
54
class ReferenceBrownianDynamics;
55
class ReferenceStochasticDynamics;
56
class ReferenceConstraintAlgorithm;
57
class ReferenceMonteCarloBarostat;
58
class ReferenceVariableStochasticDynamics;
59
class ReferenceVariableVerletDynamics;
60
class ReferenceVerletDynamics;
61
class ReferenceCustomDynamics;
62

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

107
/**
108
109
 * This kernel provides methods for setting and retrieving various state data: time, positions,
 * velocities, and forces.
110
 */
111
class ReferenceUpdateStateDataKernel : public UpdateStateDataKernel {
112
public:
113
    ReferenceUpdateStateDataKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : UpdateStateDataKernel(name, platform), data(data) {
114
115
116
117
118
119
120
121
122
123
124
125
    }
    /**
     * 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
     */
126
    double getTime(const ContextImpl& context) const;
127
128
129
130
131
    /**
     * Set the current time (in picoseconds).
     *
     * @param context    the context in which to execute this kernel
     */
132
    void setTime(ContextImpl& context, double time);
133
134
135
136
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
    /**
     * 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);
163
164
165
166
167
168
    /**
     * 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);
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
    /**
     * 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
     */
184
    void setPeriodicBoxVectors(ContextImpl& context, const Vec3& a, const Vec3& b, const Vec3& c);
Peter Eastman's avatar
Peter Eastman committed
185
186
187
188
189
190
191
192
193
194
195
196
    /**
     * 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);
197
198
199
200
private:
    ReferencePlatform::PlatformData& data;
};

201
202
203
204
205
206
/**
 * 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) :
207
            ApplyConstraintsKernel(name, platform), data(data) {
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
    }
    ~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);
223
224
225
226
227
228
229
    /**
     * 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);
230
231
private:
    ReferencePlatform::PlatformData& data;
peastman's avatar
peastman committed
232
233
    std::vector<double> masses;
    std::vector<double> inverseMasses;
234
235
};

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

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

294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
/**
 * 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:
    ReferenceCalcCustomBondForceKernel(std::string name, const Platform& platform) : CalcCustomBondForceKernel(name, platform) {
    }
    /**
     * 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);
    /**
309
     * Execute the kernel to calculate the forces and/or energy.
310
     *
311
312
313
314
     * @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
315
     */
316
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
317
318
319
320
321
322
323
    /**
     * 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);
324
325
private:
    int numBonds;
326
327
    std::vector<std::vector<int> >bondIndexArray;
    std::vector<std::vector<double> >bondParamArray;
328
    Lepton::CompiledExpression energyExpression, forceExpression;
329
330
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
331
    bool usePeriodic;
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 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);
    /**
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
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
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
     */
    void copyParametersToContext(ContextImpl& context, const HarmonicAngleForce& force);
364
365
private:
    int numAngles;
366
367
    std::vector<std::vector<int> >angleIndexArray;
    std::vector<std::vector<double> >angleParamArray;
368
    bool usePeriodic;
369
370
};

371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
/**
 * 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:
    ReferenceCalcCustomAngleForceKernel(std::string name, const Platform& platform) : CalcCustomAngleForceKernel(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
     */
    void initialize(const System& system, const CustomAngleForce& force);
    /**
386
     * Execute the kernel to calculate the forces and/or energy.
387
     *
388
389
390
391
     * @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
392
     */
393
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
394
395
396
397
398
399
400
    /**
     * 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);
401
402
private:
    int numAngles;
403
404
    std::vector<std::vector<int> >angleIndexArray;
    std::vector<std::vector<double> >angleParamArray;
405
    Lepton::CompiledExpression energyExpression, forceExpression;
406
407
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
408
    bool usePeriodic;
409
410
};

411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
/**
 * 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);
    /**
426
427
428
429
430
431
     * 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
432
     */
433
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
434
435
436
437
438
439
440
    /**
     * 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);
441
442
private:
    int numTorsions;
443
444
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
445
    bool usePeriodic;
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
};

/**
 * This kernel is invoked by RBTorsionForce to calculate the forces acting on the system and the energy of the system.
 */
class 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);
    /**
463
464
465
466
467
468
     * Execute the kernel to calculate the forces and/or energy.
     *
     * @param context        the context in which to execute this kernel
     * @param includeForces  true if forces should be calculated
     * @param includeEnergy  true if the energy should be calculated
     * @return the potential energy due to the force
469
     */
470
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
471
472
473
474
475
476
477
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the RBTorsionForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const RBTorsionForce& force);
478
479
private:
    int numTorsions;
480
481
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
482
    bool usePeriodic;
483
484
};

485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
/**
 * 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);
    /**
500
     * Execute the kernel to calculate the forces and/or energy.
501
     *
502
503
504
505
     * @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
506
     */
507
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
508
509
510
511
512
513
514
    /**
     * 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);
515
private:
peastman's avatar
peastman committed
516
    std::vector<std::vector<std::vector<double> > > coeff;
517
518
    std::vector<int> torsionMaps;
    std::vector<std::vector<int> > torsionIndices;
519
    bool usePeriodic;
520
521
};

522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
/**
 * 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:
    ReferenceCalcCustomTorsionForceKernel(std::string name, const Platform& platform) : CalcCustomTorsionForceKernel(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
     */
    void initialize(const System& system, const CustomTorsionForce& force);
    /**
537
     * Execute the kernel to calculate the forces and/or energy.
538
     *
539
540
541
542
     * @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
543
     */
544
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
545
546
547
548
549
550
551
    /**
     * 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);
552
553
private:
    int numTorsions;
554
555
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
556
    Lepton::CompiledExpression energyExpression, forceExpression;
557
558
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
559
    bool usePeriodic;
560
561
};

562
563
564
565
566
567
568
569
570
571
572
573
574
/**
 * 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
575
     */
576
    void initialize(const System& system, const NonbondedForce& force);
577
    /**
578
579
580
581
582
     * 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
583
     * @param includeReciprocal  true if reciprocal space interactions should be included
584
     * @return the potential energy due to the force
585
     */
586
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy, bool includeDirect, bool includeReciprocal);
587
588
589
590
591
592
593
    /**
     * 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);
594
595
596
597
598
599
600
601
602
    /**
     * 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;
603
604
605
    /**
     * Get the dispersion parameters being used for the dispersion term in LJPME.
     *
606
607
608
609
     * @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
610
     */
611
    void getLJPMEParameters(double& alpha, int& nx, int& ny, int& nz) const;
612
private:
613
    void computeParameters(ContextImpl& context);
Peter Eastman's avatar
Peter Eastman committed
614
    int numParticles, num14;
615
616
    std::vector<std::vector<int> >bonded14IndexArray;
    std::vector<std::vector<double> > particleParamArray, bonded14ParamArray;
617
618
    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
619
    double nonbondedCutoff, switchingDistance, rfDielectric, ewaldAlpha, ewaldDispersionAlpha, dispersionCoefficient;
620
    int kmax[3], gridSize[3], dispersionGridSize[3];
621
    bool useSwitchingFunction;
622
623
624
    std::vector<std::set<int> > exclusions;
    NonbondedMethod nonbondedMethod;
    NeighborList* neighborList;
625
626
};

627
628
629
630
631
/**
 * This kernel is invoked by CustomNonbondedForce to calculate the forces acting on the system.
 */
class ReferenceCalcCustomNonbondedForceKernel : public CalcCustomNonbondedForceKernel {
public:
632
    ReferenceCalcCustomNonbondedForceKernel(std::string name, const Platform& platform) : CalcCustomNonbondedForceKernel(name, platform), forceCopy(NULL) {
633
634
635
636
637
638
639
640
641
642
    }
    ~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);
    /**
643
     * Execute the kernel to calculate the forces and/or energy.
644
     *
645
646
647
648
     * @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
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
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 CustomNonbondedForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CustomNonbondedForce& force);
658
private:
659
    int numParticles;
660
    std::vector<std::vector<double> > particleParamArray;
peastman's avatar
peastman committed
661
    double nonbondedCutoff, switchingDistance, periodicBoxSize[3], longRangeCoefficient;
662
    bool useSwitchingFunction, hasInitializedLongRangeCorrection;
663
664
    CustomNonbondedForce* forceCopy;
    std::map<std::string, double> globalParamValues;
665
    std::vector<std::set<int> > exclusions;
666
    Lepton::CompiledExpression energyExpression, forceExpression;
667
668
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
669
    std::vector<std::pair<std::set<int>, std::set<int> > > interactionGroups;
670
    std::vector<double> longRangeCoefficientDerivs;
671
672
673
674
    NonbondedMethod nonbondedMethod;
    NeighborList* neighborList;
};

675
/**
676
 * This kernel is invoked by GBSAOBCForce to calculate the forces acting on the system.
677
 */
678
class ReferenceCalcGBSAOBCForceKernel : public CalcGBSAOBCForceKernel {
679
public:
680
    ReferenceCalcGBSAOBCForceKernel(std::string name, const Platform& platform) : CalcGBSAOBCForceKernel(name, platform) {
681
    }
682
    ~ReferenceCalcGBSAOBCForceKernel();
683
    /**
684
     * Initialize the kernel.
685
     * 
686
     * @param system     the System this kernel will be applied to
687
     * @param force      the GBSAOBCForce this kernel will be used for
688
     */
689
    void initialize(const System& system, const GBSAOBCForce& force);
690
    /**
691
692
693
694
695
696
     * 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
697
     */
698
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
699
700
701
702
703
704
705
    /**
     * 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);
706
private:
707
    ReferenceObc* obc;
peastman's avatar
peastman committed
708
    std::vector<double> charges;
709
    bool isPeriodic;
710
711
};

712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
/**
 * 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);
    /**
728
     * Execute the kernel to calculate the forces and/or energy.
729
     *
730
731
732
733
     * @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
734
     */
735
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
736
737
738
739
740
741
742
    /**
     * 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);
743
744
private:
    int numParticles;
745
    bool isPeriodic;
746
    std::vector<std::vector<double> > particleParamArray;
peastman's avatar
peastman committed
747
    double nonbondedCutoff;
748
    std::vector<std::set<int> > exclusions;
749
    std::vector<std::string> particleParameterNames, globalParameterNames, energyParamDerivNames, valueNames;
750
751
752
    std::vector<Lepton::CompiledExpression> valueExpressions;
    std::vector<std::vector<Lepton::CompiledExpression> > valueDerivExpressions;
    std::vector<std::vector<Lepton::CompiledExpression> > valueGradientExpressions;
753
    std::vector<std::vector<Lepton::CompiledExpression> > valueParamDerivExpressions;
754
    std::vector<OpenMM::CustomGBForce::ComputationType> valueTypes;
755
756
757
    std::vector<Lepton::CompiledExpression> energyExpressions;
    std::vector<std::vector<Lepton::CompiledExpression> > energyDerivExpressions;
    std::vector<std::vector<Lepton::CompiledExpression> > energyGradientExpressions;
758
    std::vector<std::vector<Lepton::CompiledExpression> > energyParamDerivExpressions;
759
760
761
762
763
    std::vector<OpenMM::CustomGBForce::ComputationType> energyTypes;
    NonbondedMethod nonbondedMethod;
    NeighborList* neighborList;
};

764
765
766
767
768
769
770
771
772
773
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:
    ReferenceCalcCustomExternalForceKernel(std::string name, const Platform& platform) : CalcCustomExternalForceKernel(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
     */
    void initialize(const System& system, const CustomExternalForce& force);
    /**
779
     * Execute the kernel to calculate the forces and/or energy.
780
     *
781
782
783
784
     * @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
785
     */
786
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
787
788
789
790
791
792
793
    /**
     * 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);
794
private:
795
    class PeriodicDistanceFunction;
796
797
    int numParticles;
    std::vector<int> particles;
798
    std::vector<std::vector<double> > particleParamArray;
799
    Lepton::CompiledExpression energyExpression, forceExpressionX, forceExpressionY, forceExpressionZ;
800
    std::vector<std::string> parameterNames, globalParameterNames;
peastman's avatar
peastman committed
801
    Vec3* boxVectors;
802
803
804
805
};

class ReferenceCalcCustomExternalForceKernel::PeriodicDistanceFunction : public Lepton::CustomFunction {
public:
peastman's avatar
peastman committed
806
807
    Vec3** boxVectorHandle;
    PeriodicDistanceFunction(Vec3** boxVectorHandle);
808
809
810
811
    int getNumArguments() const;
    double evaluate(const double* arguments) const;
    double evaluateDerivative(const double* arguments, const int* derivOrder) const;
    Lepton::CustomFunction* clone() const;
812
813
};

814
815
816
817
818
/**
 * This kernel is invoked by CustomHbondForce to calculate the forces acting on the system.
 */
class ReferenceCalcCustomHbondForceKernel : public CalcCustomHbondForceKernel {
public:
819
    ReferenceCalcCustomHbondForceKernel(std::string name, const Platform& platform) : CalcCustomHbondForceKernel(name, platform), ixn(NULL) {
820
821
822
823
824
825
826
827
828
829
    }
    ~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);
    /**
830
     * Execute the kernel to calculate the forces and/or energy.
831
     *
832
833
834
835
     * @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
836
     */
837
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
838
839
840
841
842
843
844
    /**
     * 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);
845
846
private:
    int numDonors, numAcceptors, numParticles;
847
    bool isPeriodic;
848
    std::vector<std::vector<double> > donorParamArray, acceptorParamArray;
peastman's avatar
peastman committed
849
    double nonbondedCutoff;
850
    ReferenceCustomHbondIxn* ixn;
851
    std::vector<std::set<int> > exclusions;
852
    std::vector<std::string> globalParameterNames;
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
879
880
881
882
883
884
885
886
887
/**
 * 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;
888
    std::vector<std::vector<double> > bondParamArray;
889
    ReferenceCustomCentroidBondIxn* ixn;
890
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
891
    bool usePeriodic;
892
893
};

894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
/**
 * 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);
918
919
920
921
922
923
924
    /**
     * 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);
925
private:
926
    int numBonds;
927
    std::vector<std::vector<double> > bondParamArray;
928
    ReferenceCustomCompoundBondIxn* ixn;
929
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
930
    bool usePeriodic;
931
932
};

933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
/**
 * 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
966
    double cutoffDistance;
967
    std::vector<std::vector<double> > particleParamArray;
968
969
970
971
972
    ReferenceCustomManyParticleIxn* ixn;
    std::vector<std::string> globalParameterNames;
    NonbondedMethod nonbondedMethod;
};

973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
/**
 * 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;
};

1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
/**
 * 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
1020
     * @param innerContext   the context created by the CustomCVForce for computing collective variables
1021
     */
1022
    void initialize(const System& system, const CustomCVForce& force, ContextImpl& innerContext);
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
    /**
     * 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);
private:
    ReferenceCustomCVForce* ixn;
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
};

1045
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
1076
1077
1078
1079
/**
 * 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;
};

1080
1081
1082
1083
1084
/**
 * This kernel is invoked by VerletIntegrator to take one time step.
 */
class ReferenceIntegrateVerletStepKernel : public IntegrateVerletStepKernel {
public:
1085
    ReferenceIntegrateVerletStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateVerletStepKernel(name, platform),
1086
        data(data), dynamics(0) {
1087
    }
1088
    ~ReferenceIntegrateVerletStepKernel();
1089
    /**
1090
     * Initialize the kernel.
1091
     * 
1092
1093
     * @param system     the System this kernel will be applied to
     * @param integrator the VerletIntegrator this kernel will be used for
1094
     */
1095
    void initialize(const System& system, const VerletIntegrator& integrator);
1096
1097
1098
    /**
     * Execute the kernel.
     * 
1099
1100
     * @param context    the context in which to execute this kernel
     * @param integrator the VerletIntegrator this kernel is being used for
1101
     */
1102
    void execute(ContextImpl& context, const VerletIntegrator& integrator);
1103
1104
1105
1106
1107
1108
1109
    /**
     * 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);
1110
private:
1111
    ReferencePlatform::PlatformData& data;
1112
    ReferenceVerletDynamics* dynamics;
peastman's avatar
peastman committed
1113
    std::vector<double> masses;
1114
    double prevStepSize;
1115
1116
1117
1118
1119
1120
1121
};

/**
 * This kernel is invoked by LangevinIntegrator to take one time step.
 */
class ReferenceIntegrateLangevinStepKernel : public IntegrateLangevinStepKernel {
public:
1122
    ReferenceIntegrateLangevinStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateLangevinStepKernel(name, platform),
1123
        data(data), dynamics(0) {
1124
    }
1125
    ~ReferenceIntegrateLangevinStepKernel();
1126
    /**
Peter Eastman's avatar
Peter Eastman committed
1127
     * Initialize the kernel, setting up the particle masses.
1128
     * 
1129
1130
     * @param system     the System this kernel will be applied to
     * @param integrator the LangevinIntegrator this kernel will be used for
1131
     */
1132
    void initialize(const System& system, const LangevinIntegrator& integrator);
1133
1134
1135
    /**
     * Execute the kernel.
     * 
1136
1137
     * @param context    the context in which to execute this kernel
     * @param integrator the LangevinIntegrator this kernel is being used for
1138
     */
1139
    void execute(ContextImpl& context, const LangevinIntegrator& integrator);
1140
1141
1142
1143
1144
1145
1146
    /**
     * 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);
1147
private:
1148
    ReferencePlatform::PlatformData& data;
1149
    ReferenceStochasticDynamics* dynamics;
peastman's avatar
peastman committed
1150
    std::vector<double> masses;
1151
    double prevTemp, prevFriction, prevStepSize;
1152
1153
1154
1155
1156
1157
1158
};

/**
 * This kernel is invoked by BrownianIntegrator to take one time step.
 */
class ReferenceIntegrateBrownianStepKernel : public IntegrateBrownianStepKernel {
public:
1159
    ReferenceIntegrateBrownianStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateBrownianStepKernel(name, platform),
1160
        data(data), dynamics(0) {
1161
    }
1162
    ~ReferenceIntegrateBrownianStepKernel();
1163
    /**
1164
     * Initialize the kernel.
1165
     * 
1166
1167
     * @param system     the System this kernel will be applied to
     * @param integrator the BrownianIntegrator this kernel will be used for
1168
     */
1169
    void initialize(const System& system, const BrownianIntegrator& integrator);
1170
1171
1172
    /**
     * Execute the kernel.
     * 
1173
1174
     * @param context    the context in which to execute this kernel
     * @param integrator the BrownianIntegrator this kernel is being used for
1175
     */
1176
    void execute(ContextImpl& context, const BrownianIntegrator& integrator);
1177
1178
1179
1180
1181
1182
1183
    /**
     * 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);
1184
private:
1185
    ReferencePlatform::PlatformData& data;
1186
    ReferenceBrownianDynamics* dynamics;
peastman's avatar
peastman committed
1187
    std::vector<double> masses;
1188
    double prevTemp, prevFriction, prevStepSize;
1189
1190
};

1191
1192
1193
1194
1195
1196
/**
 * 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),
1197
        data(data), dynamics(0) {
1198
1199
1200
1201
1202
1203
    }
    ~ReferenceIntegrateVariableLangevinStepKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
1204
     * @param integrator the VariableLangevinIntegrator this kernel will be used for
1205
1206
1207
1208
1209
1210
     */
    void initialize(const System& system, const VariableLangevinIntegrator& integrator);
    /**
     * Execute the kernel.
     *
     * @param context    the context in which to execute this kernel
1211
     * @param integrator the VariableLangevinIntegrator this kernel is being used for
1212
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1213
     * @return the size of the step that was taken
1214
     */
1215
    double execute(ContextImpl& context, const VariableLangevinIntegrator& integrator, double maxTime);
1216
1217
1218
1219
1220
1221
1222
    /**
     * 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);
1223
1224
1225
private:
    ReferencePlatform::PlatformData& data;
    ReferenceVariableStochasticDynamics* dynamics;
peastman's avatar
peastman committed
1226
    std::vector<double> masses;
1227
1228
1229
    double prevTemp, prevFriction, prevErrorTol;
};

1230
1231
1232
1233
1234
1235
/**
 * 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),
1236
        data(data), dynamics(0) {
1237
1238
1239
1240
1241
1242
    }
    ~ReferenceIntegrateVariableVerletStepKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
1243
     * @param integrator the VariableVerletIntegrator this kernel will be used for
1244
1245
1246
1247
1248
1249
     */
    void initialize(const System& system, const VariableVerletIntegrator& integrator);
    /**
     * Execute the kernel.
     *
     * @param context    the context in which to execute this kernel
1250
     * @param integrator the VariableVerletIntegrator this kernel is being used for
1251
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1252
     * @return the size of the step that was taken
1253
     */
1254
    double execute(ContextImpl& context, const VariableVerletIntegrator& integrator, double maxTime);
1255
1256
1257
1258
1259
1260
1261
    /**
     * 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);
1262
1263
1264
private:
    ReferencePlatform::PlatformData& data;
    ReferenceVariableVerletDynamics* dynamics;
peastman's avatar
peastman committed
1265
    std::vector<double> masses;
1266
    double prevErrorTol;
1267
1268
};

1269
1270
1271
1272
1273
1274
/**
 * 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),
1275
        data(data), dynamics(0) {
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
    }
    ~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);
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
    /**
     * 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);
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
    /**
     * 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
1340
1341
    std::vector<double> masses, globalValues;
    std::vector<std::vector<OpenMM::Vec3> > perDofValues; 
1342
1343
};

1344
/**
Peter Eastman's avatar
Peter Eastman committed
1345
 * This kernel is invoked by AndersenThermostat at the start of each time step to adjust the particle velocities.
1346
1347
1348
 */
class ReferenceApplyAndersenThermostatKernel : public ApplyAndersenThermostatKernel {
public:
1349
    ReferenceApplyAndersenThermostatKernel(std::string name, const Platform& platform) : ApplyAndersenThermostatKernel(name, platform), thermostat(0) {
1350
    }
1351
    ~ReferenceApplyAndersenThermostatKernel();
1352
    /**
1353
     * Initialize the kernel.
1354
     * 
1355
1356
     * @param system     the System this kernel will be applied to
     * @param thermostat the AndersenThermostat this kernel will be used for
1357
     */
1358
    void initialize(const System& system, const AndersenThermostat& thermostat);
1359
1360
1361
    /**
     * Execute the kernel.
     * 
1362
     * @param context    the context in which to execute this kernel
1363
     */
1364
    void execute(ContextImpl& context);
1365
1366
private:
    ReferenceAndersenThermostat* thermostat;
1367
    std::vector<std::vector<int> > particleGroups;
peastman's avatar
peastman committed
1368
    std::vector<double> masses;
1369
1370
};

1371
1372
1373
1374
1375
/**
 * This kernel is invoked by MonteCarloBarostat to adjust the periodic box volume
 */
class ReferenceApplyMonteCarloBarostatKernel : public ApplyMonteCarloBarostatKernel {
public:
1376
    ReferenceApplyMonteCarloBarostatKernel(std::string name, const Platform& platform) : ApplyMonteCarloBarostatKernel(name, platform), barostat(NULL) {
1377
1378
1379
1380
1381
1382
1383
1384
    }
    ~ReferenceApplyMonteCarloBarostatKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param barostat   the MonteCarloBarostat this kernel will be used for
     */
1385
    void initialize(const System& system, const Force& barostat);
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
    /**
     * 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
     */
1398
    void scaleCoordinates(ContextImpl& context, double scaleX, double scaleY, double scaleZ);
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
    /**
     * 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;
};

1410
1411
1412
1413
1414
/**
 * This kernel is invoked to remove center of mass motion from the system.
 */
class ReferenceRemoveCMMotionKernel : public RemoveCMMotionKernel {
public:
1415
    ReferenceRemoveCMMotionKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : RemoveCMMotionKernel(name, platform), data(data) {
1416
1417
    }
    /**
Peter Eastman's avatar
Peter Eastman committed
1418
     * Initialize the kernel, setting up the particle masses.
1419
     * 
1420
1421
     * @param system     the System this kernel will be applied to
     * @param force      the CMMotionRemover this kernel will be used for
1422
     */
1423
    void initialize(const System& system, const CMMotionRemover& force);
1424
1425
1426
    /**
     * Execute the kernel.
     * 
1427
     * @param context    the context in which to execute this kernel
1428
     */
1429
    void execute(ContextImpl& context);
1430
private:
1431
    ReferencePlatform::PlatformData& data;
1432
    std::vector<double> masses;
1433
    int frequency;
1434
1435
};

1436
1437
1438
} // namespace OpenMM

#endif /*OPENMM_REFERENCEKERNELS_H_*/