ReferenceKernels.h 65.2 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 ReferenceNoseHooverChain;
62
class ReferenceMonteCarloBarostat;
63
class ReferenceVelocityVerletDynamics;
64
class ReferenceVariableStochasticDynamics;
65
class ReferenceVariableVerletDynamics;
66
class ReferenceVerletDynamics;
67
class ReferenceCustomDynamics;
68

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

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

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

242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
/**
 * 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);
};

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

300
301
302
303
304
/**
 * 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:
305
    ReferenceCalcCustomBondForceKernel(std::string name, const Platform& platform) : CalcCustomBondForceKernel(name, platform), ixn(NULL) {
306
    }
307
    ~ReferenceCalcCustomBondForceKernel();
308
309
310
311
312
313
314
315
    /**
     * 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);
    /**
316
     * Execute the kernel to calculate the forces and/or energy.
317
     *
318
319
320
321
     * @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
322
     */
323
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
324
325
326
327
328
329
330
    /**
     * 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);
331
332
private:
    int numBonds;
333
    ReferenceCustomBondIxn* ixn;
334
335
    std::vector<std::vector<int> >bondIndexArray;
    std::vector<std::vector<double> >bondParamArray;
336
    Lepton::CompiledExpression energyExpression, forceExpression;
337
338
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
339
    bool usePeriodic;
340
341
};

342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
/**
 * 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);
    /**
357
358
359
360
361
362
     * 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
363
     */
364
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
365
366
367
368
369
370
371
    /**
     * 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);
372
373
private:
    int numAngles;
374
375
    std::vector<std::vector<int> >angleIndexArray;
    std::vector<std::vector<double> >angleParamArray;
376
    bool usePeriodic;
377
378
};

379
380
381
382
383
/**
 * This kernel is invoked by CustomAngleForce to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcCustomAngleForceKernel : public CalcCustomAngleForceKernel {
public:
384
    ReferenceCalcCustomAngleForceKernel(std::string name, const Platform& platform) : CalcCustomAngleForceKernel(name, platform), ixn(NULL) {
385
    }
386
    ~ReferenceCalcCustomAngleForceKernel();
387
388
389
390
391
392
393
394
    /**
     * 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);
    /**
395
     * Execute the kernel to calculate the forces and/or energy.
396
     *
397
398
399
400
     * @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
401
     */
402
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
403
404
405
406
407
408
409
    /**
     * 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);
410
411
private:
    int numAngles;
412
    ReferenceCustomAngleIxn* ixn;
413
414
    std::vector<std::vector<int> >angleIndexArray;
    std::vector<std::vector<double> >angleParamArray;
415
    Lepton::CompiledExpression energyExpression, forceExpression;
416
417
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
418
    bool usePeriodic;
419
420
};

421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
/**
 * 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);
    /**
436
437
438
439
440
441
     * 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
442
     */
443
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
444
445
446
447
448
449
450
    /**
     * 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);
451
452
private:
    int numTorsions;
453
454
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
455
    bool usePeriodic;
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
};

/**
 * 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);
    /**
473
474
475
476
477
478
     * 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
479
     */
480
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
481
482
483
484
485
486
487
    /**
     * 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);
488
489
private:
    int numTorsions;
490
491
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
492
    bool usePeriodic;
493
494
};

495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
/**
 * 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);
    /**
510
     * Execute the kernel to calculate the forces and/or energy.
511
     *
512
513
514
515
     * @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
516
     */
517
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
518
519
520
521
522
523
524
    /**
     * 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);
525
private:
peastman's avatar
peastman committed
526
    std::vector<std::vector<std::vector<double> > > coeff;
527
528
    std::vector<int> torsionMaps;
    std::vector<std::vector<int> > torsionIndices;
529
    bool usePeriodic;
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:
537
    ReferenceCalcCustomTorsionForceKernel(std::string name, const Platform& platform) : CalcCustomTorsionForceKernel(name, platform), ixn(NULL) {
538
    }
539
    ~ReferenceCalcCustomTorsionForceKernel();
540
541
542
543
544
545
546
547
    /**
     * 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);
    /**
548
     * Execute the kernel to calculate the forces and/or energy.
549
     *
550
551
552
553
     * @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
554
     */
555
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
556
557
558
559
560
561
562
    /**
     * 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);
563
564
private:
    int numTorsions;
565
    ReferenceCustomTorsionIxn* ixn;
566
567
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
568
    Lepton::CompiledExpression energyExpression, forceExpression;
569
570
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
571
    bool usePeriodic;
572
573
};

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

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

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

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

776
777
778
779
780
/**
 * 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:
781
    ReferenceCalcCustomExternalForceKernel(std::string name, const Platform& platform) : CalcCustomExternalForceKernel(name, platform), ixn(NULL) {
782
    }
783
    ~ReferenceCalcCustomExternalForceKernel();
784
785
786
787
788
789
790
791
    /**
     * 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);
    /**
792
     * Execute the kernel to calculate the forces and/or energy.
793
     *
794
795
796
797
     * @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
798
     */
799
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
800
801
802
803
804
805
806
    /**
     * 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);
807
private:
808
    class PeriodicDistanceFunction;
809
    int numParticles;
810
    ReferenceCustomExternalIxn* ixn;
811
    std::vector<int> particles;
812
    std::vector<std::vector<double> > particleParamArray;
813
    Lepton::CompiledExpression energyExpression, forceExpressionX, forceExpressionY, forceExpressionZ;
814
    std::vector<std::string> parameterNames, globalParameterNames;
peastman's avatar
peastman committed
815
    Vec3* boxVectors;
816
817
818
819
};

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

828
829
830
831
832
/**
 * This kernel is invoked by CustomHbondForce to calculate the forces acting on the system.
 */
class ReferenceCalcCustomHbondForceKernel : public CalcCustomHbondForceKernel {
public:
833
    ReferenceCalcCustomHbondForceKernel(std::string name, const Platform& platform) : CalcCustomHbondForceKernel(name, platform), ixn(NULL) {
834
835
836
837
838
839
840
841
842
843
    }
    ~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);
    /**
844
     * Execute the kernel to calculate the forces and/or energy.
845
     *
846
847
848
849
     * @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
850
     */
851
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
852
853
854
855
856
857
858
    /**
     * 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);
859
860
private:
    int numDonors, numAcceptors, numParticles;
861
    bool isPeriodic;
862
    std::vector<std::vector<double> > donorParamArray, acceptorParamArray;
peastman's avatar
peastman committed
863
    double nonbondedCutoff;
864
    ReferenceCustomHbondIxn* ixn;
865
    std::vector<std::set<int> > exclusions;
866
    std::vector<std::string> globalParameterNames;
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
900
901
/**
 * 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;
902
    std::vector<std::vector<double> > bondParamArray;
903
    ReferenceCustomCentroidBondIxn* ixn;
904
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
905
    bool usePeriodic;
906
907
};

908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
/**
 * 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);
932
933
934
935
936
937
938
    /**
     * 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);
939
private:
940
    int numBonds;
941
    std::vector<std::vector<double> > bondParamArray;
942
    ReferenceCustomCompoundBondIxn* ixn;
943
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
944
    bool usePeriodic;
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
978
979
/**
 * 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
980
    double cutoffDistance;
981
    std::vector<std::vector<double> > particleParamArray;
982
983
984
985
986
    ReferenceCustomManyParticleIxn* ixn;
    std::vector<std::string> globalParameterNames;
    NonbondedMethod nonbondedMethod;
};

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

1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
/**
 * 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
1034
     * @param innerContext   the context created by the CustomCVForce for computing collective variables
1035
     */
1036
    void initialize(const System& system, const CustomCVForce& force, ContextImpl& innerContext);
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
    /**
     * 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);
1054
1055
1056
1057
1058
1059
1060
    /**
     * 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);
1061
1062
1063
1064
1065
private:
    ReferenceCustomCVForce* ixn;
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
};

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

1101
1102
1103
1104
1105
/**
 * This kernel is invoked by VerletIntegrator to take one time step.
 */
class ReferenceIntegrateVerletStepKernel : public IntegrateVerletStepKernel {
public:
1106
    ReferenceIntegrateVerletStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateVerletStepKernel(name, platform),
1107
        data(data), dynamics(0) {
1108
    }
1109
    ~ReferenceIntegrateVerletStepKernel();
1110
    /**
1111
     * Initialize the kernel.
1112
     * 
1113
1114
     * @param system     the System this kernel will be applied to
     * @param integrator the VerletIntegrator this kernel will be used for
1115
     */
1116
    void initialize(const System& system, const VerletIntegrator& integrator);
1117
1118
1119
    /**
     * Execute the kernel.
     * 
1120
1121
     * @param context    the context in which to execute this kernel
     * @param integrator the VerletIntegrator this kernel is being used for
1122
     */
1123
    void execute(ContextImpl& context, const VerletIntegrator& integrator);
1124
1125
1126
1127
1128
1129
1130
    /**
     * 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);
1131
private:
1132
    ReferencePlatform::PlatformData& data;
1133
    ReferenceVerletDynamics* dynamics;
peastman's avatar
peastman committed
1134
    std::vector<double> masses;
1135
    double prevStepSize;
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
}
;
/**
 * This kernel is invoked by VelocityVerletIntegrator to take one time step.
 */
class ReferenceIntegrateVelocityVerletStepKernel : public IntegrateVelocityVerletStepKernel {
public:
    ReferenceIntegrateVelocityVerletStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateVelocityVerletStepKernel(name, platform),
        data(data), dynamics(0) {
    }
    ~ReferenceIntegrateVelocityVerletStepKernel();
    /**
     * Initialize the kernel.
     * 
     * @param system     the System this kernel will be applied to
     * @param integrator the VelocityVerletIntegrator this kernel will be used for
     */
1153
    void initialize(const System& system, const NoseHooverIntegrator& integrator);
1154
1155
1156
1157
1158
    /**
     * Execute the kernel.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the VerletIntegrator this kernel is being used for
1159
1160
     * @param forcesAreValid a reference to the parent integrator's boolean for keeping
     *                       track of the validity of the current forces.
1161
     */
1162
    void execute(ContextImpl& context, const NoseHooverIntegrator& integrator, bool &forcesAreValid);
1163
1164
1165
1166
1167
1168
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the VelocityVerletIntegrator this kernel is being used for
     */
1169
    double computeKineticEnergy(ContextImpl& context, const NoseHooverIntegrator& integrator);
1170
1171
1172
1173
1174
private:
    ReferencePlatform::PlatformData& data;
    ReferenceVelocityVerletDynamics* dynamics;
    std::vector<double> masses;
    double prevStepSize;
1175
1176
1177
1178
1179
1180
1181
};

/**
 * This kernel is invoked by LangevinIntegrator to take one time step.
 */
class ReferenceIntegrateLangevinStepKernel : public IntegrateLangevinStepKernel {
public:
1182
    ReferenceIntegrateLangevinStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateLangevinStepKernel(name, platform),
1183
        data(data), dynamics(0) {
1184
    }
1185
    ~ReferenceIntegrateLangevinStepKernel();
1186
    /**
Peter Eastman's avatar
Peter Eastman committed
1187
     * Initialize the kernel, setting up the particle masses.
1188
     * 
1189
1190
     * @param system     the System this kernel will be applied to
     * @param integrator the LangevinIntegrator this kernel will be used for
1191
     */
1192
    void initialize(const System& system, const LangevinIntegrator& integrator);
1193
1194
1195
    /**
     * Execute the kernel.
     * 
1196
1197
     * @param context    the context in which to execute this kernel
     * @param integrator the LangevinIntegrator this kernel is being used for
1198
     */
1199
    void execute(ContextImpl& context, const LangevinIntegrator& integrator);
1200
1201
1202
1203
1204
1205
1206
    /**
     * 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);
1207
private:
1208
    ReferencePlatform::PlatformData& data;
1209
    ReferenceStochasticDynamics* dynamics;
peastman's avatar
peastman committed
1210
    std::vector<double> masses;
1211
    double prevTemp, prevFriction, prevStepSize;
1212
1213
1214
1215
1216
1217
1218
};

/**
 * This kernel is invoked by BrownianIntegrator to take one time step.
 */
class ReferenceIntegrateBrownianStepKernel : public IntegrateBrownianStepKernel {
public:
1219
    ReferenceIntegrateBrownianStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateBrownianStepKernel(name, platform),
1220
        data(data), dynamics(0) {
1221
    }
1222
    ~ReferenceIntegrateBrownianStepKernel();
1223
    /**
1224
     * Initialize the kernel.
1225
     * 
1226
1227
     * @param system     the System this kernel will be applied to
     * @param integrator the BrownianIntegrator this kernel will be used for
1228
     */
1229
    void initialize(const System& system, const BrownianIntegrator& integrator);
1230
1231
1232
    /**
     * Execute the kernel.
     * 
1233
1234
     * @param context    the context in which to execute this kernel
     * @param integrator the BrownianIntegrator this kernel is being used for
1235
     */
1236
    void execute(ContextImpl& context, const BrownianIntegrator& integrator);
1237
1238
1239
1240
1241
1242
1243
    /**
     * 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);
1244
private:
1245
    ReferencePlatform::PlatformData& data;
1246
    ReferenceBrownianDynamics* dynamics;
peastman's avatar
peastman committed
1247
    std::vector<double> masses;
1248
    double prevTemp, prevFriction, prevStepSize;
1249
1250
};

1251
1252
1253
1254
1255
1256
/**
 * 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),
1257
        data(data), dynamics(0) {
1258
1259
1260
1261
1262
1263
    }
    ~ReferenceIntegrateVariableLangevinStepKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
1264
     * @param integrator the VariableLangevinIntegrator this kernel will be used for
1265
1266
1267
1268
1269
1270
     */
    void initialize(const System& system, const VariableLangevinIntegrator& integrator);
    /**
     * Execute the kernel.
     *
     * @param context    the context in which to execute this kernel
1271
     * @param integrator the VariableLangevinIntegrator this kernel is being used for
1272
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1273
     * @return the size of the step that was taken
1274
     */
1275
    double execute(ContextImpl& context, const VariableLangevinIntegrator& integrator, double maxTime);
1276
1277
1278
1279
1280
1281
1282
    /**
     * 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);
1283
1284
1285
private:
    ReferencePlatform::PlatformData& data;
    ReferenceVariableStochasticDynamics* dynamics;
peastman's avatar
peastman committed
1286
    std::vector<double> masses;
1287
1288
1289
    double prevTemp, prevFriction, prevErrorTol;
};

1290
1291
1292
1293
1294
1295
/**
 * 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),
1296
        data(data), dynamics(0) {
1297
1298
1299
1300
1301
1302
    }
    ~ReferenceIntegrateVariableVerletStepKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
1303
     * @param integrator the VariableVerletIntegrator this kernel will be used for
1304
1305
1306
1307
1308
1309
     */
    void initialize(const System& system, const VariableVerletIntegrator& integrator);
    /**
     * Execute the kernel.
     *
     * @param context    the context in which to execute this kernel
1310
     * @param integrator the VariableVerletIntegrator this kernel is being used for
1311
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1312
     * @return the size of the step that was taken
1313
     */
1314
    double execute(ContextImpl& context, const VariableVerletIntegrator& integrator, double maxTime);
1315
1316
1317
1318
1319
1320
1321
    /**
     * 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);
1322
1323
1324
private:
    ReferencePlatform::PlatformData& data;
    ReferenceVariableVerletDynamics* dynamics;
peastman's avatar
peastman committed
1325
    std::vector<double> masses;
1326
    double prevErrorTol;
1327
1328
};

1329
1330
1331
1332
1333
1334
/**
 * 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),
1335
        data(data), dynamics(0) {
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
    }
    ~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);
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
    /**
     * 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);
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
    /**
     * 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
1400
1401
    std::vector<double> masses, globalValues;
    std::vector<std::vector<OpenMM::Vec3> > perDofValues; 
1402
1403
};

1404
/**
Peter Eastman's avatar
Peter Eastman committed
1405
 * This kernel is invoked by AndersenThermostat at the start of each time step to adjust the particle velocities.
1406
1407
1408
 */
class ReferenceApplyAndersenThermostatKernel : public ApplyAndersenThermostatKernel {
public:
1409
    ReferenceApplyAndersenThermostatKernel(std::string name, const Platform& platform) : ApplyAndersenThermostatKernel(name, platform), thermostat(0) {
1410
    }
1411
    ~ReferenceApplyAndersenThermostatKernel();
1412
    /**
1413
     * Initialize the kernel.
1414
     * 
1415
1416
     * @param system     the System this kernel will be applied to
     * @param thermostat the AndersenThermostat this kernel will be used for
1417
     */
1418
    void initialize(const System& system, const AndersenThermostat& thermostat);
1419
1420
1421
    /**
     * Execute the kernel.
     * 
1422
     * @param context    the context in which to execute this kernel
1423
     */
1424
    void execute(ContextImpl& context);
1425
1426
private:
    ReferenceAndersenThermostat* thermostat;
1427
    std::vector<std::vector<int> > particleGroups;
peastman's avatar
peastman committed
1428
    std::vector<double> masses;
1429
1430
};

1431
1432
1433
1434
/**
 * This kernel is invoked by NoseHooverChain at the start of each time step to adjust the thermostat
 * and update the associated particle velocities.
 */
1435
class ReferenceNoseHooverChainKernel : public NoseHooverChainKernel {
1436
public:
1437
    ReferenceNoseHooverChainKernel(std::string name, const Platform& platform) : NoseHooverChainKernel(name, platform), chainPropagator(0) {
1438
    }
1439
    ~ReferenceNoseHooverChainKernel();
1440
1441
1442
    /**
     * Initialize the kernel.
     */
1443
    virtual void initialize();
1444
    /**
1445
     * Execute the kernel that propagates the Nose Hoover chain and determines the velocity scale factor.
1446
     * 
1447
1448
     * @param context  the context in which to execute this kernel
     * @param noseHooverChain the object describing the chain to be propagated.
1449
     * @param kineticEnergy the {absolute, relative} kinetic energies of the particles being thermostated by this chain.
1450
     * @param timeStep the time step used by the integrator.
1451
     * @return the velocity scale factors to apply to the {absolute, relative} motion of particles associated with this heat bath.
1452
     */
1453
    virtual std::pair<double, double> propagateChain(ContextImpl& context, const NoseHooverChain &nhc, std::pair<double, double> kineticEnergy, double timeStep);
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
    /**
     * Execute the kernal that computes the total (kinetic + potential) heat bath energy.
     *
     * @param context the context in which to execute this kernel
     * @param noseHooverChain the chain whose energy is to be determined.
     * @return the total heat bath energy.
     */
    virtual double computeHeatBathEnergy(ContextImpl& context, const NoseHooverChain &nhc);
    /**
     * Execute the kernel that computes the kinetic energy for a subset of atoms,
     * or the relative kinetic energy of Drude particles with respect to their parent atoms
     *
     * @param context the context in which to execute this kernel
     * @param noseHooverChain the chain whose energy is to be determined.
1468
1469
     * @param downloadValue whether the computed value should be downloaded and returned.
     *
1470
     */
1471
     virtual std::pair<double, double> computeMaskedKineticEnergy(ContextImpl& context, const NoseHooverChain &noseHooverChain, bool downloadValue);
1472
1473
1474
1475
1476
1477

    /**
     * Execute the kernel that scales the velocities of particles associated with a nose hoover chain
     *
     * @param context the context in which to execute this kernel
     * @param noseHooverChain the chain whose energy is to be determined.
1478
     * @param scaleFactor the multiplicative factor by which {absolute, relative} velocities are scaled.
1479
     */
1480
    virtual void scaleVelocities(ContextImpl& context, const NoseHooverChain &noseHooverChain, std::pair<double, double> scaleFactor);
1481
1482


1483
1484
1485
1486
private:
    ReferenceNoseHooverChain* chainPropagator;
};

1487
1488
1489
1490
1491
/**
 * This kernel is invoked by MonteCarloBarostat to adjust the periodic box volume
 */
class ReferenceApplyMonteCarloBarostatKernel : public ApplyMonteCarloBarostatKernel {
public:
1492
    ReferenceApplyMonteCarloBarostatKernel(std::string name, const Platform& platform) : ApplyMonteCarloBarostatKernel(name, platform), barostat(NULL) {
1493
1494
1495
1496
1497
1498
1499
1500
    }
    ~ReferenceApplyMonteCarloBarostatKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param barostat   the MonteCarloBarostat this kernel will be used for
     */
1501
    void initialize(const System& system, const Force& barostat);
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
    /**
     * 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
     */
1514
    void scaleCoordinates(ContextImpl& context, double scaleX, double scaleY, double scaleZ);
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
    /**
     * 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;
};

1526
1527
1528
1529
1530
/**
 * This kernel is invoked to remove center of mass motion from the system.
 */
class ReferenceRemoveCMMotionKernel : public RemoveCMMotionKernel {
public:
1531
    ReferenceRemoveCMMotionKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : RemoveCMMotionKernel(name, platform), data(data) {
1532
1533
    }
    /**
Peter Eastman's avatar
Peter Eastman committed
1534
     * Initialize the kernel, setting up the particle masses.
1535
     * 
1536
1537
     * @param system     the System this kernel will be applied to
     * @param force      the CMMotionRemover this kernel will be used for
1538
     */
1539
    void initialize(const System& system, const CMMotionRemover& force);
1540
1541
1542
    /**
     * Execute the kernel.
     * 
1543
     * @param context    the context in which to execute this kernel
1544
     */
1545
    void execute(ContextImpl& context);
1546
private:
1547
    ReferencePlatform::PlatformData& data;
1548
    std::vector<double> masses;
1549
    int frequency;
1550
1551
};

1552
1553
1554
} // namespace OpenMM

#endif /*OPENMM_REFERENCEKERNELS_H_*/