ReferenceKernels.h 73.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-2024 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
#include "openmm/internal/CustomCPPForceImpl.h"
38
#include "openmm/internal/CustomNonbondedForceImpl.h"
39
#include "openmm/internal/windowsExport.h"
40
41
#include "SimTKOpenMMRealType.h"
#include "ReferenceNeighborList.h"
42
#include "lepton/CompiledExpression.h"
43
#include "lepton/CustomFunction.h"
44
45
#include <array>
#include <utility>
46

47
48
namespace OpenMM {

49
class ReferenceObc;
50
class ReferenceAndersenThermostat;
51
class ReferenceLangevinMiddleDynamics;
52
53
54
55
class ReferenceCustomBondIxn;
class ReferenceCustomAngleIxn;
class ReferenceCustomTorsionIxn;
class ReferenceCustomExternalIxn;
56
class ReferenceCustomCentroidBondIxn;
57
class ReferenceCustomCompoundBondIxn;
58
class ReferenceCustomCVForce;
59
class ReferenceCustomHbondIxn;
60
class ReferenceCustomManyParticleIxn;
61
class ReferenceGayBerneForce;
62
class ReferenceBrownianDynamics;
63
class ReferenceConstraintAlgorithm;
64
class ReferenceNoseHooverChain;
65
class ReferenceMonteCarloBarostat;
66
class ReferenceNoseHooverDynamics;
67
class ReferenceVariableStochasticDynamics;
68
class ReferenceVariableVerletDynamics;
69
class ReferenceVerletDynamics;
70
class ReferenceCustomDynamics;
71

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

116
/**
117
118
 * This kernel provides methods for setting and retrieving various state data: time, positions,
 * velocities, and forces.
119
 */
120
class OPENMM_EXPORT ReferenceUpdateStateDataKernel : public UpdateStateDataKernel {
121
public:
122
    ReferenceUpdateStateDataKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : UpdateStateDataKernel(name, platform), data(data) {
123
124
125
126
127
128
129
130
131
132
133
134
    }
    /**
     * 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
     */
135
    double getTime(const ContextImpl& context) const;
136
137
138
139
140
    /**
     * Set the current time (in picoseconds).
     *
     * @param context    the context in which to execute this kernel
     */
141
    void setTime(ContextImpl& context, double time);
142
143
144
145
146
147
148
149
150
151
152
153
    /**
     * Get the current step count
     *
     * @param context    the context in which to execute this kernel
     */
    long long getStepCount(const ContextImpl& context) const;
    /**
     * Set the current step count
     *
     * @param context    the context in which to execute this kernel
     */
    void setStepCount(const ContextImpl& context, long long count);
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
    /**
     * 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);
178
179
180
181
182
183
184
185
186
    /**
     * Compute velocities, shifted in time to account for a leapfrog integrator.  The shift
     * is based on the most recently computed forces.
     * 
     * @param context     the context in which to execute this kernel
     * @param timeShift   the amount by which to shift the velocities in time
     * @param velocities  the shifted velocities are returned in this
     */
    void computeShiftedVelocities(ContextImpl& context, double timeShift, std::vector<Vec3>& velocities);
187
188
189
190
191
192
    /**
     * Get the current forces on all particles.
     *
     * @param forces  on exit, this contains the forces
     */
    void getForces(ContextImpl& context, std::vector<Vec3>& forces);
193
194
195
196
197
198
    /**
     * 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);
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
    /**
     * 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
     */
214
    void setPeriodicBoxVectors(ContextImpl& context, const Vec3& a, const Vec3& b, const Vec3& c);
Peter Eastman's avatar
Peter Eastman committed
215
216
217
218
219
220
221
222
223
224
225
226
    /**
     * 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);
227
228
private:
    ReferencePlatform::PlatformData& data;
229
    std::vector<double> masses;
230
231
};

232
233
234
235
236
237
/**
 * 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) :
238
            ApplyConstraintsKernel(name, platform), data(data) {
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
    }
    ~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);
254
255
256
257
258
259
260
    /**
     * 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);
261
262
private:
    ReferencePlatform::PlatformData& data;
peastman's avatar
peastman committed
263
264
    std::vector<double> masses;
    std::vector<double> inverseMasses;
265
266
};

267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
/**
 * 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);
};

288
/**
289
 * This kernel is invoked by HarmonicBondForce to calculate the forces acting on the system and the energy of the system.
290
 */
291
class ReferenceCalcHarmonicBondForceKernel : public CalcHarmonicBondForceKernel {
292
public:
293
    ReferenceCalcHarmonicBondForceKernel(std::string name, const Platform& platform) : CalcHarmonicBondForceKernel(name, platform) {
294
295
    }
    /**
296
     * Initialize the kernel.
297
     * 
298
     * @param system     the System this kernel will be applied to
299
300
301
302
     * @param force      the HarmonicBondForce this kernel will be used for
     */
    void initialize(const System& system, const HarmonicBondForce& force);
    /**
303
304
305
306
307
308
     * 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
309
     */
310
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
311
312
313
314
315
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the HarmonicBondForce to copy the parameters from
316
317
     * @param firstBond  the index of the first bond whose parameters might have changed
     * @param lastBond   the index of the last bond whose parameters might have changed
318
     */
319
    void copyParametersToContext(ContextImpl& context, const HarmonicBondForce& force, int firstBond, int lastBond);
320
321
private:
    int numBonds;
322
323
    std::vector<std::vector<int> >bondIndexArray;
    std::vector<std::vector<double> >bondParamArray;
324
    bool usePeriodic;
325
326
};

327
328
329
330
331
/**
 * 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:
332
    ReferenceCalcCustomBondForceKernel(std::string name, const Platform& platform) : CalcCustomBondForceKernel(name, platform), ixn(NULL) {
333
    }
334
    ~ReferenceCalcCustomBondForceKernel();
335
336
337
338
339
340
341
342
    /**
     * 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);
    /**
343
     * Execute the kernel to calculate the forces and/or energy.
344
     *
345
346
347
348
     * @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
349
     */
350
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
351
352
353
    /**
     * Copy changed parameters over to a context.
     *
354
355
     * @param firstBond  the index of the first bond whose parameters might have changed
     * @param lastBond   the index of the last bond whose parameters might have changed
356
     */
357
    void copyParametersToContext(ContextImpl& context, const CustomBondForce& force, int firstBond, int lastBond);
358
359
private:
    int numBonds;
360
    ReferenceCustomBondIxn* ixn;
361
362
    std::vector<std::vector<int> >bondIndexArray;
    std::vector<std::vector<double> >bondParamArray;
363
    Lepton::CompiledExpression energyExpression, forceExpression;
364
365
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
366
    bool usePeriodic;
367
368
};

369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
/**
 * 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);
    /**
384
385
386
387
388
389
     * 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
390
     */
391
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
392
393
394
    /**
     * Copy changed parameters over to a context.
     *
395
396
     * @param firstAngle the index of the first bond whose parameters might have changed
     * @param lastAngle  the index of the last bond whose parameters might have changed
397
     */
398
    void copyParametersToContext(ContextImpl& context, const HarmonicAngleForce& force, int firstAngle, int lastAngle);
399
400
private:
    int numAngles;
401
402
    std::vector<std::vector<int> >angleIndexArray;
    std::vector<std::vector<double> >angleParamArray;
403
    bool usePeriodic;
404
405
};

406
407
408
409
410
/**
 * 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:
411
    ReferenceCalcCustomAngleForceKernel(std::string name, const Platform& platform) : CalcCustomAngleForceKernel(name, platform), ixn(NULL) {
412
    }
413
    ~ReferenceCalcCustomAngleForceKernel();
414
415
416
417
418
419
420
421
    /**
     * 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);
    /**
422
     * Execute the kernel to calculate the forces and/or energy.
423
     *
424
425
426
427
     * @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
428
     */
429
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
430
431
432
    /**
     * Copy changed parameters over to a context.
     *
433
434
     * @param firstAngle the index of the first bond whose parameters might have changed
     * @param lastAngle  the index of the last bond whose parameters might have changed
435
     */
436
    void copyParametersToContext(ContextImpl& context, const CustomAngleForce& force, int firstAngle, int lastAngle);
437
438
private:
    int numAngles;
439
    ReferenceCustomAngleIxn* ixn;
440
441
    std::vector<std::vector<int> >angleIndexArray;
    std::vector<std::vector<double> >angleParamArray;
442
    Lepton::CompiledExpression energyExpression, forceExpression;
443
444
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
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 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);
    /**
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
    /**
     * Copy changed parameters over to a context.
     *
474
475
476
477
     * @param context      the context to copy parameters to
     * @param force        the PeriodicTorsionForce to copy the parameters from
     * @param firstTorsion the index of the first torsion whose parameters might have changed
     * @param lastTorsion  the index of the last torsion whose parameters might have changed
478
     */
479
    void copyParametersToContext(ContextImpl& context, const PeriodicTorsionForce& force, int firstTorsion, int lastTorsion);
480
481
private:
    int numTorsions;
482
483
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
484
    bool usePeriodic;
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
};

/**
 * 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);
    /**
502
503
504
505
506
507
     * 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
508
     */
509
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
510
511
512
513
514
515
516
    /**
     * 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);
517
518
private:
    int numTorsions;
519
520
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
521
    bool usePeriodic;
522
523
};

524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
/**
 * 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);
    /**
539
     * Execute the kernel to calculate the forces and/or energy.
540
     *
541
542
543
544
     * @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
545
     */
546
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
547
548
549
550
551
552
553
    /**
     * 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);
554
private:
peastman's avatar
peastman committed
555
    std::vector<std::vector<std::vector<double> > > coeff;
556
557
    std::vector<int> torsionMaps;
    std::vector<std::vector<int> > torsionIndices;
558
    bool usePeriodic;
559
560
};

561
562
563
564
565
/**
 * 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:
566
    ReferenceCalcCustomTorsionForceKernel(std::string name, const Platform& platform) : CalcCustomTorsionForceKernel(name, platform), ixn(NULL) {
567
    }
568
    ~ReferenceCalcCustomTorsionForceKernel();
569
570
571
572
573
574
575
576
    /**
     * 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);
    /**
577
     * Execute the kernel to calculate the forces and/or energy.
578
     *
579
580
581
582
     * @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
583
     */
584
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
585
586
587
    /**
     * Copy changed parameters over to a context.
     *
588
589
590
591
     * @param context      the context to copy parameters to
     * @param force        the CustomTorsionForce to copy the parameters from
     * @param firstTorsion the index of the first torsion whose parameters might have changed
     * @param lastTorsion  the index of the last torsion whose parameters might have changed
592
     */
593
    void copyParametersToContext(ContextImpl& context, const CustomTorsionForce& force, int firstTorsion, int lastTorsion);
594
595
private:
    int numTorsions;
596
    ReferenceCustomTorsionIxn* ixn;
597
598
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
599
    Lepton::CompiledExpression energyExpression, forceExpression;
600
601
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
602
    bool usePeriodic;
603
604
};

605
606
607
608
609
610
611
612
613
614
615
616
617
/**
 * 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
618
     */
619
    void initialize(const System& system, const NonbondedForce& force);
620
    /**
621
622
623
624
625
     * 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
626
     * @param includeReciprocal  true if reciprocal space interactions should be included
627
     * @return the potential energy due to the force
628
     */
629
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy, bool includeDirect, bool includeReciprocal);
630
631
632
    /**
     * Copy changed parameters over to a context.
     *
633
634
635
636
637
638
     * @param context        the context to copy parameters to
     * @param force          the NonbondedForce to copy the parameters from
     * @param firstParticle  the index of the first particle whose parameters might have changed
     * @param lastParticle   the index of the last particle whose parameters might have changed
     * @param firstException the index of the first exception whose parameters might have changed
     * @param lastException  the index of the last exception whose parameters might have changed
639
     */
640
    void copyParametersToContext(ContextImpl& context, const NonbondedForce& force, int firstParticle, int lastParticle, int firstException, int lastException);
641
642
643
644
645
646
647
648
649
    /**
     * 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;
650
651
652
    /**
     * Get the dispersion parameters being used for the dispersion term in LJPME.
     *
653
654
655
656
     * @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
657
     */
658
    void getLJPMEParameters(double& alpha, int& nx, int& ny, int& nz) const;
659
private:
660
    void computeParameters(ContextImpl& context);
Peter Eastman's avatar
Peter Eastman committed
661
    int numParticles, num14;
662
663
    std::vector<std::vector<int> >bonded14IndexArray;
    std::vector<std::vector<double> > particleParamArray, bonded14ParamArray;
664
665
    std::vector<std::array<double, 3> > baseParticleParams, baseExceptionParams;
    std::map<std::pair<std::string, int>, std::array<double, 3> > particleParamOffsets, exceptionParamOffsets;
666
    std::map<int, int> nb14Index;
peastman's avatar
peastman committed
667
    double nonbondedCutoff, switchingDistance, rfDielectric, ewaldAlpha, ewaldDispersionAlpha, dispersionCoefficient;
668
    int kmax[3], gridSize[3], dispersionGridSize[3];
669
    bool useSwitchingFunction, exceptionsArePeriodic;
670
671
672
    std::vector<std::set<int> > exclusions;
    NonbondedMethod nonbondedMethod;
    NeighborList* neighborList;
673
674
};

675
676
677
678
679
/**
 * This kernel is invoked by CustomNonbondedForce to calculate the forces acting on the system.
 */
class ReferenceCalcCustomNonbondedForceKernel : public CalcCustomNonbondedForceKernel {
public:
680
    ReferenceCalcCustomNonbondedForceKernel(std::string name, const Platform& platform) : CalcCustomNonbondedForceKernel(name, platform), forceCopy(NULL) {
681
682
683
684
685
686
687
688
689
690
    }
    ~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);
    /**
691
     * Execute the kernel to calculate the forces and/or energy.
692
     *
693
694
695
696
     * @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
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomNonbondedForce to copy the parameters from
704
705
     * @param firstParticle  the index of the first particle whose parameters might have changed
     * @param lastParticle   the index of the last particle whose parameters might have changed
706
     */
707
    void copyParametersToContext(ContextImpl& context, const CustomNonbondedForce& force, int firstParticle, int lastParticle);
708
private:
709
    void createExpressions(const CustomNonbondedForce& force);
710
    int numParticles;
711
    std::vector<std::vector<double> > particleParamArray;
peastman's avatar
peastman committed
712
    double nonbondedCutoff, switchingDistance, periodicBoxSize[3], longRangeCoefficient;
713
    bool useSwitchingFunction, hasInitializedLongRangeCorrection;
714
    CustomNonbondedForce* forceCopy;
715
    CustomNonbondedForceImpl::LongRangeCorrectionData longRangeCorrectionData;
716
    std::map<std::string, double> globalParamValues;
717
    std::vector<std::set<int> > exclusions;
718
    Lepton::CompiledExpression energyExpression, forceExpression;
719
720
    std::vector<Lepton::CompiledExpression> computedValueExpressions, energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, computedValueNames, energyParamDerivNames;
721
    std::vector<std::pair<std::set<int>, std::set<int> > > interactionGroups;
722
    std::vector<double> longRangeCoefficientDerivs;
723
    std::map<std::string, int> tabulatedFunctionUpdateCount;
724
725
726
727
    NonbondedMethod nonbondedMethod;
    NeighborList* neighborList;
};

728
/**
729
 * This kernel is invoked by GBSAOBCForce to calculate the forces acting on the system.
730
 */
731
class ReferenceCalcGBSAOBCForceKernel : public CalcGBSAOBCForceKernel {
732
public:
733
    ReferenceCalcGBSAOBCForceKernel(std::string name, const Platform& platform) : CalcGBSAOBCForceKernel(name, platform) {
734
    }
735
    ~ReferenceCalcGBSAOBCForceKernel();
736
    /**
737
     * Initialize the kernel.
738
     * 
739
     * @param system     the System this kernel will be applied to
740
     * @param force      the GBSAOBCForce this kernel will be used for
741
     */
742
    void initialize(const System& system, const GBSAOBCForce& force);
743
    /**
744
745
746
747
748
749
     * 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
750
     */
751
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
752
753
754
755
756
757
758
    /**
     * 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);
759
private:
760
    ReferenceObc* obc;
peastman's avatar
peastman committed
761
    std::vector<double> charges;
762
    bool isPeriodic;
763
764
};

765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
/**
 * 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);
    /**
781
     * Execute the kernel to calculate the forces and/or energy.
782
     *
783
784
785
786
     * @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
787
     */
788
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
789
790
791
792
793
794
795
    /**
     * 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);
796
private:
797
    void createExpressions(const CustomGBForce& force);
798
    int numParticles;
799
    bool isPeriodic;
800
    std::vector<std::vector<double> > particleParamArray;
peastman's avatar
peastman committed
801
    double nonbondedCutoff;
802
    std::vector<std::set<int> > exclusions;
803
    std::vector<std::string> particleParameterNames, globalParameterNames, energyParamDerivNames, valueNames;
804
805
806
    std::vector<Lepton::CompiledExpression> valueExpressions;
    std::vector<std::vector<Lepton::CompiledExpression> > valueDerivExpressions;
    std::vector<std::vector<Lepton::CompiledExpression> > valueGradientExpressions;
807
    std::vector<std::vector<Lepton::CompiledExpression> > valueParamDerivExpressions;
808
    std::vector<OpenMM::CustomGBForce::ComputationType> valueTypes;
809
810
811
    std::vector<Lepton::CompiledExpression> energyExpressions;
    std::vector<std::vector<Lepton::CompiledExpression> > energyDerivExpressions;
    std::vector<std::vector<Lepton::CompiledExpression> > energyGradientExpressions;
812
    std::vector<std::vector<Lepton::CompiledExpression> > energyParamDerivExpressions;
813
    std::vector<OpenMM::CustomGBForce::ComputationType> energyTypes;
814
    std::map<std::string, int> tabulatedFunctionUpdateCount;
815
816
817
818
    NonbondedMethod nonbondedMethod;
    NeighborList* neighborList;
};

819
820
821
822
823
/**
 * 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:
824
    ReferenceCalcCustomExternalForceKernel(std::string name, const Platform& platform) : CalcCustomExternalForceKernel(name, platform), ixn(NULL) {
825
    }
826
    ~ReferenceCalcCustomExternalForceKernel();
827
828
829
830
831
832
833
834
    /**
     * 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);
    /**
835
     * Execute the kernel to calculate the forces and/or energy.
836
     *
837
838
839
840
     * @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
841
     */
842
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
843
844
845
    /**
     * Copy changed parameters over to a context.
     *
846
847
848
849
     * @param context        the context to copy parameters to
     * @param force          the CustomExternalForce to copy the parameters from
     * @param firstParticle  the index of the first particle whose parameters might have changed
     * @param lastParticle   the index of the last particle whose parameters might have changed
850
     */
851
    void copyParametersToContext(ContextImpl& context, const CustomExternalForce& force, int firstParticle, int lastParticle);
852
853
private:
    int numParticles;
854
    ReferenceCustomExternalIxn* ixn;
855
    std::vector<int> particles;
856
    std::vector<std::vector<double> > particleParamArray;
857
    Lepton::CompiledExpression energyExpression, forceExpressionX, forceExpressionY, forceExpressionZ;
858
    std::vector<std::string> parameterNames, globalParameterNames;
peastman's avatar
peastman committed
859
    Vec3* boxVectors;
860
861
};

862
863
864
865
866
/**
 * This kernel is invoked by CustomHbondForce to calculate the forces acting on the system.
 */
class ReferenceCalcCustomHbondForceKernel : public CalcCustomHbondForceKernel {
public:
867
    ReferenceCalcCustomHbondForceKernel(std::string name, const Platform& platform) : CalcCustomHbondForceKernel(name, platform), ixn(NULL) {
868
869
870
871
872
873
874
875
876
877
    }
    ~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);
    /**
878
     * Execute the kernel to calculate the forces and/or energy.
879
     *
880
881
882
883
     * @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
884
     */
885
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
886
887
888
889
890
891
892
    /**
     * 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);
893
private:
894
    void createInteraction(const CustomHbondForce& force);
895
    int numDonors, numAcceptors, numParticles;
896
    bool isPeriodic;
897
    std::vector<std::vector<int> > donorParticles, acceptorParticles;
898
    std::vector<std::vector<double> > donorParamArray, acceptorParamArray;
peastman's avatar
peastman committed
899
    double nonbondedCutoff;
900
    ReferenceCustomHbondIxn* ixn;
901
    std::vector<std::set<int> > exclusions;
902
    std::vector<std::string> globalParameterNames;
903
    std::map<std::string, int> tabulatedFunctionUpdateCount;
904
905
};

906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
/**
 * 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:
938
    void createInteraction(const CustomCentroidBondForce& force);
939
    int numBonds, numParticles;
940
941
942
    std::vector<std::vector<int> > bondGroups;
    std::vector<std::vector<int> > groupAtoms;
    std::vector<std::vector<double> > normalizedWeights;
943
    std::vector<std::vector<double> > bondParamArray;
944
    ReferenceCustomCentroidBondIxn* ixn;
945
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
946
    std::map<std::string, int> tabulatedFunctionUpdateCount;
947
    bool usePeriodic;
948
    Vec3* boxVectors;
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
/**
 * 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);
975
976
977
978
979
980
981
    /**
     * 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);
982
private:
983
    void createInteraction(const CustomCompoundBondForce& force);
984
    int numBonds;
985
    std::vector<std::vector<int> > bondParticles;
986
    std::vector<std::vector<double> > bondParamArray;
987
    ReferenceCustomCompoundBondIxn* ixn;
988
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
989
    std::map<std::string, int> tabulatedFunctionUpdateCount;
990
    bool usePeriodic;
991
    Vec3* boxVectors;
992
993
};

994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
/**
 * This kernel is invoked by 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
1027
    double cutoffDistance;
1028
    std::vector<std::vector<double> > particleParamArray;
1029
1030
    ReferenceCustomManyParticleIxn* ixn;
    std::vector<std::string> globalParameterNames;
1031
    std::map<std::string, int> tabulatedFunctionUpdateCount;
1032
1033
1034
    NonbondedMethod nonbondedMethod;
};

1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
/**
 * 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;
};

1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
/**
 * 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
1082
     * @param innerContext   the context created by the CustomCVForce for computing collective variables
1083
     */
1084
    void initialize(const System& system, const CustomCVForce& force, ContextImpl& innerContext);
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
    /**
     * 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);
1102
1103
1104
1105
1106
1107
1108
    /**
     * 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);
1109
1110
1111
1112
1113
private:
    ReferenceCustomCVForce* ixn;
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
};

1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
/**
 * 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;
};

1149
1150
1151
1152
1153
/**
 * This kernel is invoked by VerletIntegrator to take one time step.
 */
class ReferenceIntegrateVerletStepKernel : public IntegrateVerletStepKernel {
public:
1154
    ReferenceIntegrateVerletStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateVerletStepKernel(name, platform),
1155
        data(data), dynamics(0) {
1156
    }
1157
    ~ReferenceIntegrateVerletStepKernel();
1158
    /**
1159
     * Initialize the kernel.
1160
     * 
1161
1162
     * @param system     the System this kernel will be applied to
     * @param integrator the VerletIntegrator this kernel will be used for
1163
     */
1164
    void initialize(const System& system, const VerletIntegrator& integrator);
1165
1166
1167
    /**
     * Execute the kernel.
     * 
1168
1169
     * @param context    the context in which to execute this kernel
     * @param integrator the VerletIntegrator this kernel is being used for
1170
     */
1171
    void execute(ContextImpl& context, const VerletIntegrator& integrator);
1172
1173
1174
1175
1176
1177
1178
    /**
     * 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);
1179
private:
1180
    ReferencePlatform::PlatformData& data;
1181
    ReferenceVerletDynamics* dynamics;
peastman's avatar
peastman committed
1182
    std::vector<double> masses;
1183
    double prevStepSize;
1184
1185
1186
}
;
/**
1187
 * This kernel is invoked by NoseHooverIntegrator to take one time step.
1188
 */
1189
class ReferenceIntegrateNoseHooverStepKernel : public IntegrateNoseHooverStepKernel {
1190
public:
1191
    ReferenceIntegrateNoseHooverStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateNoseHooverStepKernel(name, platform),
1192
1193
        data(data), dynamics(0) {
    }
1194
    ~ReferenceIntegrateNoseHooverStepKernel();
1195
1196
1197
1198
    /**
     * Initialize the kernel.
     * 
     * @param system     the System this kernel will be applied to
1199
     * @param integrator the NoseHooverIntegrator this kernel will be used for
1200
     */
1201
    void initialize(const System& system, const NoseHooverIntegrator& integrator);
1202
1203
1204
1205
1206
    /**
     * Execute the kernel.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the VerletIntegrator this kernel is being used for
1207
1208
     * @param forcesAreValid a reference to the parent integrator's boolean for keeping
     *                       track of the validity of the current forces.
1209
     */
1210
    void execute(ContextImpl& context, const NoseHooverIntegrator& integrator, bool &forcesAreValid);
1211
1212
1213
1214
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
1215
     * @param integrator the NoseHooverIntegrator this kernel is being used for
1216
     */
1217
    double computeKineticEnergy(ContextImpl& context, const NoseHooverIntegrator& integrator);
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
    /**
     * Execute the kernel that propagates the Nose Hoover chain and determines the velocity scale factor.
     * 
     * @param context  the context in which to execute this kernel
     * @param noseHooverChain the object describing the chain to be propagated.
     * @param kineticEnergy the {center of mass, relative} kineticEnergies of the particles being thermostated by this chain.
     * @param timeStep the time step used by the integrator.
     * @return the velocity scale factor to apply to the particles associated with this heat bath.
     */
    std::pair<double, double> propagateChain(ContextImpl& context, const NoseHooverChain &noseHooverChain, std::pair<double, double> kineticEnergy, double timeStep);
    /**
     * 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.
     */
    double computeHeatBathEnergy(ContextImpl& context, const NoseHooverChain &noseHooverChain);
    /**
     * 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.
     * @param downloadValue whether the computed value should be downloaded and returned.
     */
    std::pair<double, double> computeMaskedKineticEnergy(ContextImpl& context, const NoseHooverChain &noseHooverChain, bool downloadValue);
    /**
     * 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.
     * @param scaleFactor the multiplicative factor by which {absolute, relative} velocities are scaled.
     */
    void scaleVelocities(ContextImpl& context, const NoseHooverChain &noseHooverChain, std::pair<double, double> scaleFactor);
1253
1254
1255
1256
1257
1258
1259
1260
    /**
     * Write the chain states to a checkpoint.
     */
    void createCheckpoint(ContextImpl& context, std::ostream& stream) const;
    /**
     * Load the chain states from a checkpoint.
     */
    void loadCheckpoint(ContextImpl& context, std::istream& stream);
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
    /**
     * Get the internal states of all chains.
     * 
     * @param context       the context for which to get the states
     * @param positions     element [i][j] contains the position of bead j for chain i
     * @param velocities    element [i][j] contains the velocity of bead j for chain i
     */
    void getChainStates(ContextImpl& context, std::vector<std::vector<double> >& positions, std::vector<std::vector<double> >& velocities) const;
    /**
     * Set the internal states of all chains.
     * 
     * @param context       the context for which to get the states
     * @param positions     element [i][j] contains the position of bead j for chain i
     * @param velocities    element [i][j] contains the velocity of bead j for chain i
     */
    void setChainStates(ContextImpl& context, const std::vector<std::vector<double> >& positions, const std::vector<std::vector<double> >& velocities);
1277
1278
private:
    ReferencePlatform::PlatformData& data;
1279
1280
    ReferenceNoseHooverChain* chainPropagator;
    ReferenceNoseHooverDynamics* dynamics;
1281
    std::vector<double> masses;
1282
1283
    std::vector<std::vector<double> > chainPositions;
    std::vector<std::vector<double> > chainVelocities;
1284
    double prevStepSize;
1285
1286
};

1287
/**
1288
 * This kernel is invoked by LangevinMiddleIntegrator to take one time step.
1289
 */
1290
class ReferenceIntegrateLangevinMiddleStepKernel : public IntegrateLangevinMiddleStepKernel {
1291
public:
1292
    ReferenceIntegrateLangevinMiddleStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateLangevinMiddleStepKernel(name, platform),
1293
1294
        data(data), dynamics(0) {
    }
1295
    ~ReferenceIntegrateLangevinMiddleStepKernel();
1296
1297
1298
1299
    /**
     * Initialize the kernel, setting up the particle masses.
     * 
     * @param system     the System this kernel will be applied to
1300
     * @param integrator the LangevinMiddleIntegrator this kernel will be used for
1301
     */
1302
    void initialize(const System& system, const LangevinMiddleIntegrator& integrator);
1303
1304
1305
1306
    /**
     * Execute the kernel.
     * 
     * @param context    the context in which to execute this kernel
1307
     * @param integrator the LangevinMiddleIntegrator this kernel is being used for
1308
     */
1309
    void execute(ContextImpl& context, const LangevinMiddleIntegrator& integrator);
1310
1311
1312
1313
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
1314
     * @param integrator the LangevinMiddleIntegrator this kernel is being used for
1315
     */
1316
    double computeKineticEnergy(ContextImpl& context, const LangevinMiddleIntegrator& integrator);
1317
1318
private:
    ReferencePlatform::PlatformData& data;
1319
    ReferenceLangevinMiddleDynamics* dynamics;
1320
1321
1322
1323
    std::vector<double> masses;
    double prevTemp, prevFriction, prevStepSize;
};

1324
1325
1326
1327
1328
/**
 * This kernel is invoked by BrownianIntegrator to take one time step.
 */
class ReferenceIntegrateBrownianStepKernel : public IntegrateBrownianStepKernel {
public:
1329
    ReferenceIntegrateBrownianStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateBrownianStepKernel(name, platform),
1330
        data(data), dynamics(0) {
1331
    }
1332
    ~ReferenceIntegrateBrownianStepKernel();
1333
    /**
1334
     * Initialize the kernel.
1335
     * 
1336
1337
     * @param system     the System this kernel will be applied to
     * @param integrator the BrownianIntegrator this kernel will be used for
1338
     */
1339
    void initialize(const System& system, const BrownianIntegrator& integrator);
1340
1341
1342
    /**
     * Execute the kernel.
     * 
1343
1344
     * @param context    the context in which to execute this kernel
     * @param integrator the BrownianIntegrator this kernel is being used for
1345
     */
1346
    void execute(ContextImpl& context, const BrownianIntegrator& integrator);
1347
1348
1349
1350
1351
1352
1353
    /**
     * 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);
1354
private:
1355
    ReferencePlatform::PlatformData& data;
1356
    ReferenceBrownianDynamics* dynamics;
peastman's avatar
peastman committed
1357
    std::vector<double> masses;
1358
    double prevTemp, prevFriction, prevStepSize;
1359
1360
};

1361
1362
1363
1364
1365
1366
/**
 * 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),
1367
        data(data), dynamics(0) {
1368
1369
1370
1371
1372
1373
    }
    ~ReferenceIntegrateVariableLangevinStepKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
1374
     * @param integrator the VariableLangevinIntegrator this kernel will be used for
1375
1376
1377
1378
1379
1380
     */
    void initialize(const System& system, const VariableLangevinIntegrator& integrator);
    /**
     * Execute the kernel.
     *
     * @param context    the context in which to execute this kernel
1381
     * @param integrator the VariableLangevinIntegrator this kernel is being used for
1382
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1383
     * @return the size of the step that was taken
1384
     */
1385
    double execute(ContextImpl& context, const VariableLangevinIntegrator& integrator, double maxTime);
1386
1387
1388
1389
1390
1391
1392
    /**
     * 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);
1393
1394
1395
private:
    ReferencePlatform::PlatformData& data;
    ReferenceVariableStochasticDynamics* dynamics;
peastman's avatar
peastman committed
1396
    std::vector<double> masses;
1397
1398
1399
    double prevTemp, prevFriction, prevErrorTol;
};

1400
1401
1402
1403
1404
1405
/**
 * 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),
1406
        data(data), dynamics(0) {
1407
1408
1409
1410
1411
1412
    }
    ~ReferenceIntegrateVariableVerletStepKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
1413
     * @param integrator the VariableVerletIntegrator this kernel will be used for
1414
1415
1416
1417
1418
1419
     */
    void initialize(const System& system, const VariableVerletIntegrator& integrator);
    /**
     * Execute the kernel.
     *
     * @param context    the context in which to execute this kernel
1420
     * @param integrator the VariableVerletIntegrator this kernel is being used for
1421
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1422
     * @return the size of the step that was taken
1423
     */
1424
    double execute(ContextImpl& context, const VariableVerletIntegrator& integrator, double maxTime);
1425
1426
1427
1428
1429
1430
1431
    /**
     * 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);
1432
1433
1434
private:
    ReferencePlatform::PlatformData& data;
    ReferenceVariableVerletDynamics* dynamics;
peastman's avatar
peastman committed
1435
    std::vector<double> masses;
1436
    double prevErrorTol;
1437
1438
};

1439
1440
1441
1442
1443
1444
/**
 * 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),
1445
        data(data), dynamics(0) {
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
    }
    ~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);
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
    /**
     * 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);
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
    /**
     * 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
1510
1511
    std::vector<double> masses, globalValues;
    std::vector<std::vector<OpenMM::Vec3> > perDofValues; 
1512
1513
};

1514
/**
Peter Eastman's avatar
Peter Eastman committed
1515
 * This kernel is invoked by AndersenThermostat at the start of each time step to adjust the particle velocities.
1516
1517
1518
 */
class ReferenceApplyAndersenThermostatKernel : public ApplyAndersenThermostatKernel {
public:
1519
    ReferenceApplyAndersenThermostatKernel(std::string name, const Platform& platform) : ApplyAndersenThermostatKernel(name, platform), thermostat(0) {
1520
    }
1521
    ~ReferenceApplyAndersenThermostatKernel();
1522
    /**
1523
     * Initialize the kernel.
1524
     * 
1525
1526
     * @param system     the System this kernel will be applied to
     * @param thermostat the AndersenThermostat this kernel will be used for
1527
     */
1528
    void initialize(const System& system, const AndersenThermostat& thermostat);
1529
1530
1531
    /**
     * Execute the kernel.
     * 
1532
     * @param context    the context in which to execute this kernel
1533
     */
1534
    void execute(ContextImpl& context);
1535
1536
private:
    ReferenceAndersenThermostat* thermostat;
1537
    std::vector<std::vector<int> > particleGroups;
peastman's avatar
peastman committed
1538
    std::vector<double> masses;
1539
1540
};

1541
1542
1543
1544
1545
/**
 * This kernel is invoked by MonteCarloBarostat to adjust the periodic box volume
 */
class ReferenceApplyMonteCarloBarostatKernel : public ApplyMonteCarloBarostatKernel {
public:
1546
    ReferenceApplyMonteCarloBarostatKernel(std::string name, const Platform& platform) : ApplyMonteCarloBarostatKernel(name, platform), barostat(NULL) {
1547
1548
1549
1550
1551
1552
1553
    }
    ~ReferenceApplyMonteCarloBarostatKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param barostat   the MonteCarloBarostat this kernel will be used for
1554
     * @param rigidMolecules  whether molecules should be kept rigid while scaling coordinates
1555
     */
1556
    void initialize(const System& system, const Force& barostat, bool rigidMolecules=true);
1557
1558
1559
1560
1561
1562
1563
    /**
     * Save the coordinates before attempting a Monte Carlo step.  This allows us to restore them
     * if the step is rejected.
     *
     * @param context    the context in which to execute this kernel
     */
    void saveCoordinates(ContextImpl& context);
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
    /**
     * 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
     */
1576
    void scaleCoordinates(ContextImpl& context, double scaleX, double scaleY, double scaleZ);
1577
    /**
1578
1579
     * Reject the most recent Monte Carlo step, restoring the particle positions to where they were when
     * saveCoordinates() was last called.
1580
1581
1582
1583
1584
     *
     * @param context    the context in which to execute this kernel
     */
    void restoreCoordinates(ContextImpl& context);
private:
1585
    bool rigidMolecules;
1586
1587
1588
    ReferenceMonteCarloBarostat* barostat;
};

1589
1590
1591
1592
1593
/**
 * This kernel is invoked to remove center of mass motion from the system.
 */
class ReferenceRemoveCMMotionKernel : public RemoveCMMotionKernel {
public:
1594
    ReferenceRemoveCMMotionKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : RemoveCMMotionKernel(name, platform), data(data) {
1595
1596
    }
    /**
Peter Eastman's avatar
Peter Eastman committed
1597
     * Initialize the kernel, setting up the particle masses.
1598
     * 
1599
1600
     * @param system     the System this kernel will be applied to
     * @param force      the CMMotionRemover this kernel will be used for
1601
     */
1602
    void initialize(const System& system, const CMMotionRemover& force);
1603
1604
1605
    /**
     * Execute the kernel.
     * 
1606
     * @param context    the context in which to execute this kernel
1607
     */
1608
    void execute(ContextImpl& context);
1609
private:
1610
    ReferencePlatform::PlatformData& data;
1611
    std::vector<double> masses;
1612
    int frequency;
1613
1614
};

1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
/**
 * This kernel is invoked by ATMForce to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcATMForceKernel : public CalcATMForceKernel {
public:
    ReferenceCalcATMForceKernel(std::string name, const Platform& platform) : CalcATMForceKernel(name, platform) {
    }
    /**
     * Initialize the kernel.
     * 
     * @param system     the System this kernel will be applied to
     * @param force      the ATMForce this kernel will be used for
     */
    void initialize(const System& system, const ATMForce& force);
    /**
     * Scale the forces from the inner contexts and apply them to the main context.
     *
     * @param context        the context in which to execute this kernel
     * @param innerContext1  the first inner context
     * @param innerContext2  the second inner context
     * @param dEdu0          the derivative of the final energy with respect to the first inner context's energy
     * @param dEdu1          the derivative of the final energy with respect to the second inner context's energy
     * @param energyParamDerivs  derivatives of the final energy with respect to global parameters
     */
    void applyForces(ContextImpl& context, ContextImpl& innerContext0, ContextImpl& innerContext1,
                     double dEdu0, double dEdu1, const std::map<std::string, double>& energyParamDerivs);
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the ATMForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const ATMForce& force);
    /**
     * Copy state information to the inner contexts.
     *
     * @param context        the context in which to execute this kernel
     * @param innerContext1  the first context created by the ATMForce for computing displaced energy
     * @param innerContext2  the second context created by the ATMForce for computing displaced energy
     */
    void copyState(ContextImpl& context, ContextImpl& innerContext0, ContextImpl& innerContext1);
private:
    int numParticles;
    std::vector<Vec3> displ1;
    std::vector<Vec3> displ0;
};

1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
/**
 * This kernel is invoked by CustomCPPForceImpl to calculate the forces acting on the system and the energy of the system.
 */
class ReferenceCalcCustomCPPForceKernel : public CalcCustomCPPForceKernel {
public:
    ReferenceCalcCustomCPPForceKernel(std::string name, const Platform& platform) : CalcCustomCPPForceKernel(name, platform), force(NULL) {
    }
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param force      the CustomCPPForceImpl this kernel will be used for
     */
    void initialize(const System& system, CustomCPPForceImpl& 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);
private:
    CustomCPPForceImpl* force;
    std::vector<Vec3> forces;
};

1690
1691
1692
} // namespace OpenMM

#endif /*OPENMM_REFERENCEKERNELS_H_*/