ReferenceKernels.h 72.8 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
#ifndef OPENMM_REFERENCEKERNELS_H_
#define OPENMM_REFERENCEKERNELS_H_

/* -------------------------------------------------------------------------- *
 *                                   OpenMM                                   *
 * -------------------------------------------------------------------------- *
 * This is part of the OpenMM molecular simulation toolkit originating from   *
 * Simbios, the NIH National Center for Physics-Based Simulation of           *
 * Biological Structures at Stanford, funded under the NIH Roadmap for        *
 * Medical Research, grant U54 GM072970. See https://simtk.org.               *
 *                                                                            *
12
 * Portions copyright (c) 2008-2023 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
40
#include "SimTKOpenMMRealType.h"
#include "ReferenceNeighborList.h"
41
#include "lepton/CompiledExpression.h"
42
#include "lepton/CustomFunction.h"
43
44
#include <array>
#include <utility>
45

46
47
namespace OpenMM {

48
class ReferenceObc;
49
class ReferenceAndersenThermostat;
50
class ReferenceLangevinMiddleDynamics;
51
52
53
54
class ReferenceCustomBondIxn;
class ReferenceCustomAngleIxn;
class ReferenceCustomTorsionIxn;
class ReferenceCustomExternalIxn;
55
class ReferenceCustomCentroidBondIxn;
56
class ReferenceCustomCompoundBondIxn;
57
class ReferenceCustomCVForce;
58
class ReferenceCustomHbondIxn;
59
class ReferenceCustomManyParticleIxn;
60
class ReferenceGayBerneForce;
61
class ReferenceBrownianDynamics;
62
class ReferenceStochasticDynamics;
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 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
316
317
    /**
     * 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);
318
319
private:
    int numBonds;
320
321
    std::vector<std::vector<int> >bondIndexArray;
    std::vector<std::vector<double> >bondParamArray;
322
    bool usePeriodic;
323
324
};

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

367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
/**
 * 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);
    /**
382
383
384
385
386
387
     * 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
388
     */
389
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
390
391
392
393
394
395
396
    /**
     * 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);
397
398
private:
    int numAngles;
399
400
    std::vector<std::vector<int> >angleIndexArray;
    std::vector<std::vector<double> >angleParamArray;
401
    bool usePeriodic;
402
403
};

404
405
406
407
408
/**
 * 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:
409
    ReferenceCalcCustomAngleForceKernel(std::string name, const Platform& platform) : CalcCustomAngleForceKernel(name, platform), ixn(NULL) {
410
    }
411
    ~ReferenceCalcCustomAngleForceKernel();
412
413
414
415
416
417
418
419
    /**
     * 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);
    /**
420
     * Execute the kernel to calculate the forces and/or energy.
421
     *
422
423
424
425
     * @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
426
     */
427
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
428
429
430
431
432
433
434
    /**
     * 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);
435
436
private:
    int numAngles;
437
    ReferenceCustomAngleIxn* ixn;
438
439
    std::vector<std::vector<int> >angleIndexArray;
    std::vector<std::vector<double> >angleParamArray;
440
    Lepton::CompiledExpression energyExpression, forceExpression;
441
442
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
443
    bool usePeriodic;
444
445
};

446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
/**
 * 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);
    /**
461
462
463
464
465
466
     * 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
467
     */
468
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
469
470
471
472
473
474
475
    /**
     * 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);
476
477
private:
    int numTorsions;
478
479
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
480
    bool usePeriodic;
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
};

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

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

557
558
559
560
561
/**
 * 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:
562
    ReferenceCalcCustomTorsionForceKernel(std::string name, const Platform& platform) : CalcCustomTorsionForceKernel(name, platform), ixn(NULL) {
563
    }
564
    ~ReferenceCalcCustomTorsionForceKernel();
565
566
567
568
569
570
571
572
    /**
     * 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);
    /**
573
     * Execute the kernel to calculate the forces and/or energy.
574
     *
575
576
577
578
     * @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
579
     */
580
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
581
582
583
584
585
586
587
    /**
     * 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);
588
589
private:
    int numTorsions;
590
    ReferenceCustomTorsionIxn* ixn;
591
592
    std::vector<std::vector<int> >torsionIndexArray;
    std::vector<std::vector<double> >torsionParamArray;
593
    Lepton::CompiledExpression energyExpression, forceExpression;
594
595
    std::vector<Lepton::CompiledExpression> energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, energyParamDerivNames;
596
    bool usePeriodic;
597
598
};

599
600
601
602
603
604
605
606
607
608
609
610
611
/**
 * 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
612
     */
613
    void initialize(const System& system, const NonbondedForce& force);
614
    /**
615
616
617
618
619
     * 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
620
     * @param includeReciprocal  true if reciprocal space interactions should be included
621
     * @return the potential energy due to the force
622
     */
623
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy, bool includeDirect, bool includeReciprocal);
624
625
626
627
628
629
630
    /**
     * 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);
631
632
633
634
635
636
637
638
639
    /**
     * 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;
640
641
642
    /**
     * Get the dispersion parameters being used for the dispersion term in LJPME.
     *
643
644
645
646
     * @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
647
     */
648
    void getLJPMEParameters(double& alpha, int& nx, int& ny, int& nz) const;
649
private:
650
    void computeParameters(ContextImpl& context);
Peter Eastman's avatar
Peter Eastman committed
651
    int numParticles, num14;
652
653
    std::vector<std::vector<int> >bonded14IndexArray;
    std::vector<std::vector<double> > particleParamArray, bonded14ParamArray;
654
655
    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
656
    double nonbondedCutoff, switchingDistance, rfDielectric, ewaldAlpha, ewaldDispersionAlpha, dispersionCoefficient;
657
    int kmax[3], gridSize[3], dispersionGridSize[3];
658
    bool useSwitchingFunction, exceptionsArePeriodic;
659
660
661
    std::vector<std::set<int> > exclusions;
    NonbondedMethod nonbondedMethod;
    NeighborList* neighborList;
662
663
};

664
665
666
667
668
/**
 * This kernel is invoked by CustomNonbondedForce to calculate the forces acting on the system.
 */
class ReferenceCalcCustomNonbondedForceKernel : public CalcCustomNonbondedForceKernel {
public:
669
    ReferenceCalcCustomNonbondedForceKernel(std::string name, const Platform& platform) : CalcCustomNonbondedForceKernel(name, platform), forceCopy(NULL) {
670
671
672
673
674
675
676
677
678
679
    }
    ~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);
    /**
680
     * Execute the kernel to calculate the forces and/or energy.
681
     *
682
683
684
685
     * @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
686
     */
687
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
688
689
690
691
692
693
694
    /**
     * 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);
695
private:
696
    void createExpressions(const CustomNonbondedForce& force);
697
    int numParticles;
698
    std::vector<std::vector<double> > particleParamArray;
peastman's avatar
peastman committed
699
    double nonbondedCutoff, switchingDistance, periodicBoxSize[3], longRangeCoefficient;
700
    bool useSwitchingFunction, hasInitializedLongRangeCorrection;
701
    CustomNonbondedForce* forceCopy;
702
    CustomNonbondedForceImpl::LongRangeCorrectionData longRangeCorrectionData;
703
    std::map<std::string, double> globalParamValues;
704
    std::vector<std::set<int> > exclusions;
705
    Lepton::CompiledExpression energyExpression, forceExpression;
706
707
    std::vector<Lepton::CompiledExpression> computedValueExpressions, energyParamDerivExpressions;
    std::vector<std::string> parameterNames, globalParameterNames, computedValueNames, energyParamDerivNames;
708
    std::vector<std::pair<std::set<int>, std::set<int> > > interactionGroups;
709
    std::vector<double> longRangeCoefficientDerivs;
710
    std::map<std::string, int> tabulatedFunctionUpdateCount;
711
712
713
714
    NonbondedMethod nonbondedMethod;
    NeighborList* neighborList;
};

715
/**
716
 * This kernel is invoked by GBSAOBCForce to calculate the forces acting on the system.
717
 */
718
class ReferenceCalcGBSAOBCForceKernel : public CalcGBSAOBCForceKernel {
719
public:
720
    ReferenceCalcGBSAOBCForceKernel(std::string name, const Platform& platform) : CalcGBSAOBCForceKernel(name, platform) {
721
    }
722
    ~ReferenceCalcGBSAOBCForceKernel();
723
    /**
724
     * Initialize the kernel.
725
     * 
726
     * @param system     the System this kernel will be applied to
727
     * @param force      the GBSAOBCForce this kernel will be used for
728
     */
729
    void initialize(const System& system, const GBSAOBCForce& force);
730
    /**
731
732
733
734
735
736
     * 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
737
     */
738
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
739
740
741
742
743
744
745
    /**
     * 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);
746
private:
747
    ReferenceObc* obc;
peastman's avatar
peastman committed
748
    std::vector<double> charges;
749
    bool isPeriodic;
750
751
};

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

806
807
808
809
810
/**
 * 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:
811
    ReferenceCalcCustomExternalForceKernel(std::string name, const Platform& platform) : CalcCustomExternalForceKernel(name, platform), ixn(NULL) {
812
    }
813
    ~ReferenceCalcCustomExternalForceKernel();
814
815
816
817
818
819
820
821
    /**
     * 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);
    /**
822
     * Execute the kernel to calculate the forces and/or energy.
823
     *
824
825
826
827
     * @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
828
     */
829
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
830
831
832
833
834
835
836
    /**
     * Copy changed parameters over to a context.
     *
     * @param context    the context to copy parameters to
     * @param force      the CustomExternalForce to copy the parameters from
     */
    void copyParametersToContext(ContextImpl& context, const CustomExternalForce& force);
837
838
private:
    int numParticles;
839
    ReferenceCustomExternalIxn* ixn;
840
    std::vector<int> particles;
841
    std::vector<std::vector<double> > particleParamArray;
842
    Lepton::CompiledExpression energyExpression, forceExpressionX, forceExpressionY, forceExpressionZ;
843
    std::vector<std::string> parameterNames, globalParameterNames;
peastman's avatar
peastman committed
844
    Vec3* boxVectors;
845
846
};

847
848
849
850
851
/**
 * This kernel is invoked by CustomHbondForce to calculate the forces acting on the system.
 */
class ReferenceCalcCustomHbondForceKernel : public CalcCustomHbondForceKernel {
public:
852
    ReferenceCalcCustomHbondForceKernel(std::string name, const Platform& platform) : CalcCustomHbondForceKernel(name, platform), ixn(NULL) {
853
854
855
856
857
858
859
860
861
862
    }
    ~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);
    /**
863
     * Execute the kernel to calculate the forces and/or energy.
864
     *
865
866
867
868
     * @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
869
     */
870
    double execute(ContextImpl& context, bool includeForces, bool includeEnergy);
871
872
873
874
875
876
877
    /**
     * 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);
878
private:
879
    void createInteraction(const CustomHbondForce& force);
880
    int numDonors, numAcceptors, numParticles;
881
    bool isPeriodic;
882
    std::vector<std::vector<int> > donorParticles, acceptorParticles;
883
    std::vector<std::vector<double> > donorParamArray, acceptorParamArray;
peastman's avatar
peastman committed
884
    double nonbondedCutoff;
885
    ReferenceCustomHbondIxn* ixn;
886
    std::vector<std::set<int> > exclusions;
887
    std::vector<std::string> globalParameterNames;
888
    std::map<std::string, int> tabulatedFunctionUpdateCount;
889
890
};

891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
/**
 * 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:
923
    void createInteraction(const CustomCentroidBondForce& force);
924
    int numBonds, numParticles;
925
926
927
    std::vector<std::vector<int> > bondGroups;
    std::vector<std::vector<int> > groupAtoms;
    std::vector<std::vector<double> > normalizedWeights;
928
    std::vector<std::vector<double> > bondParamArray;
929
    ReferenceCustomCentroidBondIxn* ixn;
930
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
931
    std::map<std::string, int> tabulatedFunctionUpdateCount;
932
    bool usePeriodic;
933
    Vec3* boxVectors;
934
935
};

936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
/**
 * 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);
960
961
962
963
964
965
966
    /**
     * 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);
967
private:
968
    void createInteraction(const CustomCompoundBondForce& force);
969
    int numBonds;
970
    std::vector<std::vector<int> > bondParticles;
971
    std::vector<std::vector<double> > bondParamArray;
972
    ReferenceCustomCompoundBondIxn* ixn;
973
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
974
    std::map<std::string, int> tabulatedFunctionUpdateCount;
975
    bool usePeriodic;
976
    Vec3* boxVectors;
977
978
};

979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
/**
 * 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
1012
    double cutoffDistance;
1013
    std::vector<std::vector<double> > particleParamArray;
1014
1015
    ReferenceCustomManyParticleIxn* ixn;
    std::vector<std::string> globalParameterNames;
1016
    std::map<std::string, int> tabulatedFunctionUpdateCount;
1017
1018
1019
    NonbondedMethod nonbondedMethod;
};

1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
/**
 * 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;
};

1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
/**
 * 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
1067
     * @param innerContext   the context created by the CustomCVForce for computing collective variables
1068
     */
1069
    void initialize(const System& system, const CustomCVForce& force, ContextImpl& innerContext);
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
    /**
     * 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);
1087
1088
1089
1090
1091
1092
1093
    /**
     * 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);
1094
1095
1096
1097
1098
private:
    ReferenceCustomCVForce* ixn;
    std::vector<std::string> globalParameterNames, energyParamDerivNames;
};

1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
/**
 * 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;
};

1134
1135
1136
1137
1138
/**
 * This kernel is invoked by VerletIntegrator to take one time step.
 */
class ReferenceIntegrateVerletStepKernel : public IntegrateVerletStepKernel {
public:
1139
    ReferenceIntegrateVerletStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateVerletStepKernel(name, platform),
1140
        data(data), dynamics(0) {
1141
    }
1142
    ~ReferenceIntegrateVerletStepKernel();
1143
    /**
1144
     * Initialize the kernel.
1145
     * 
1146
1147
     * @param system     the System this kernel will be applied to
     * @param integrator the VerletIntegrator this kernel will be used for
1148
     */
1149
    void initialize(const System& system, const VerletIntegrator& integrator);
1150
1151
1152
    /**
     * Execute the kernel.
     * 
1153
1154
     * @param context    the context in which to execute this kernel
     * @param integrator the VerletIntegrator this kernel is being used for
1155
     */
1156
    void execute(ContextImpl& context, const VerletIntegrator& integrator);
1157
1158
1159
1160
1161
1162
1163
    /**
     * 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);
1164
private:
1165
    ReferencePlatform::PlatformData& data;
1166
    ReferenceVerletDynamics* dynamics;
peastman's avatar
peastman committed
1167
    std::vector<double> masses;
1168
    double prevStepSize;
1169
1170
1171
}
;
/**
1172
 * This kernel is invoked by NoseHooverIntegrator to take one time step.
1173
 */
1174
class ReferenceIntegrateNoseHooverStepKernel : public IntegrateNoseHooverStepKernel {
1175
public:
1176
    ReferenceIntegrateNoseHooverStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateNoseHooverStepKernel(name, platform),
1177
1178
        data(data), dynamics(0) {
    }
1179
    ~ReferenceIntegrateNoseHooverStepKernel();
1180
1181
1182
1183
    /**
     * Initialize the kernel.
     * 
     * @param system     the System this kernel will be applied to
1184
     * @param integrator the NoseHooverIntegrator this kernel will be used for
1185
     */
1186
    void initialize(const System& system, const NoseHooverIntegrator& integrator);
1187
1188
1189
1190
1191
    /**
     * Execute the kernel.
     * 
     * @param context    the context in which to execute this kernel
     * @param integrator the VerletIntegrator this kernel is being used for
1192
1193
     * @param forcesAreValid a reference to the parent integrator's boolean for keeping
     *                       track of the validity of the current forces.
1194
     */
1195
    void execute(ContextImpl& context, const NoseHooverIntegrator& integrator, bool &forcesAreValid);
1196
1197
1198
1199
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
1200
     * @param integrator the NoseHooverIntegrator this kernel is being used for
1201
     */
1202
    double computeKineticEnergy(ContextImpl& context, const NoseHooverIntegrator& integrator);
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
    /**
     * 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);
1238
1239
1240
1241
1242
1243
1244
1245
    /**
     * 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);
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
    /**
     * 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);
1262
1263
private:
    ReferencePlatform::PlatformData& data;
1264
1265
    ReferenceNoseHooverChain* chainPropagator;
    ReferenceNoseHooverDynamics* dynamics;
1266
    std::vector<double> masses;
1267
1268
    std::vector<std::vector<double> > chainPositions;
    std::vector<std::vector<double> > chainVelocities;
1269
    double prevStepSize;
1270
1271
1272
1273
1274
1275
1276
};

/**
 * This kernel is invoked by LangevinIntegrator to take one time step.
 */
class ReferenceIntegrateLangevinStepKernel : public IntegrateLangevinStepKernel {
public:
1277
    ReferenceIntegrateLangevinStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateLangevinStepKernel(name, platform),
1278
        data(data), dynamics(0) {
1279
    }
1280
    ~ReferenceIntegrateLangevinStepKernel();
1281
    /**
Peter Eastman's avatar
Peter Eastman committed
1282
     * Initialize the kernel, setting up the particle masses.
1283
     * 
1284
1285
     * @param system     the System this kernel will be applied to
     * @param integrator the LangevinIntegrator this kernel will be used for
1286
     */
1287
    void initialize(const System& system, const LangevinIntegrator& integrator);
1288
1289
1290
    /**
     * Execute the kernel.
     * 
1291
1292
     * @param context    the context in which to execute this kernel
     * @param integrator the LangevinIntegrator this kernel is being used for
1293
     */
1294
    void execute(ContextImpl& context, const LangevinIntegrator& integrator);
1295
1296
1297
1298
1299
1300
1301
    /**
     * 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);
1302
private:
1303
    ReferencePlatform::PlatformData& data;
1304
    ReferenceStochasticDynamics* dynamics;
peastman's avatar
peastman committed
1305
    std::vector<double> masses;
1306
    double prevTemp, prevFriction, prevStepSize;
1307
1308
};

1309
/**
1310
 * This kernel is invoked by LangevinMiddleIntegrator to take one time step.
1311
 */
1312
class ReferenceIntegrateLangevinMiddleStepKernel : public IntegrateLangevinMiddleStepKernel {
1313
public:
1314
    ReferenceIntegrateLangevinMiddleStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateLangevinMiddleStepKernel(name, platform),
1315
1316
        data(data), dynamics(0) {
    }
1317
    ~ReferenceIntegrateLangevinMiddleStepKernel();
1318
1319
1320
1321
    /**
     * Initialize the kernel, setting up the particle masses.
     * 
     * @param system     the System this kernel will be applied to
1322
     * @param integrator the LangevinMiddleIntegrator this kernel will be used for
1323
     */
1324
    void initialize(const System& system, const LangevinMiddleIntegrator& integrator);
1325
1326
1327
1328
    /**
     * Execute the kernel.
     * 
     * @param context    the context in which to execute this kernel
1329
     * @param integrator the LangevinMiddleIntegrator this kernel is being used for
1330
     */
1331
    void execute(ContextImpl& context, const LangevinMiddleIntegrator& integrator);
1332
1333
1334
1335
    /**
     * Compute the kinetic energy.
     * 
     * @param context    the context in which to execute this kernel
1336
     * @param integrator the LangevinMiddleIntegrator this kernel is being used for
1337
     */
1338
    double computeKineticEnergy(ContextImpl& context, const LangevinMiddleIntegrator& integrator);
1339
1340
private:
    ReferencePlatform::PlatformData& data;
1341
    ReferenceLangevinMiddleDynamics* dynamics;
1342
1343
1344
1345
    std::vector<double> masses;
    double prevTemp, prevFriction, prevStepSize;
};

1346
1347
1348
1349
1350
/**
 * This kernel is invoked by BrownianIntegrator to take one time step.
 */
class ReferenceIntegrateBrownianStepKernel : public IntegrateBrownianStepKernel {
public:
1351
    ReferenceIntegrateBrownianStepKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : IntegrateBrownianStepKernel(name, platform),
1352
        data(data), dynamics(0) {
1353
    }
1354
    ~ReferenceIntegrateBrownianStepKernel();
1355
    /**
1356
     * Initialize the kernel.
1357
     * 
1358
1359
     * @param system     the System this kernel will be applied to
     * @param integrator the BrownianIntegrator this kernel will be used for
1360
     */
1361
    void initialize(const System& system, const BrownianIntegrator& integrator);
1362
1363
1364
    /**
     * Execute the kernel.
     * 
1365
1366
     * @param context    the context in which to execute this kernel
     * @param integrator the BrownianIntegrator this kernel is being used for
1367
     */
1368
    void execute(ContextImpl& context, const BrownianIntegrator& integrator);
1369
1370
1371
1372
1373
1374
1375
    /**
     * 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);
1376
private:
1377
    ReferencePlatform::PlatformData& data;
1378
    ReferenceBrownianDynamics* dynamics;
peastman's avatar
peastman committed
1379
    std::vector<double> masses;
1380
    double prevTemp, prevFriction, prevStepSize;
1381
1382
};

1383
1384
1385
1386
1387
1388
/**
 * 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),
1389
        data(data), dynamics(0) {
1390
1391
1392
1393
1394
1395
    }
    ~ReferenceIntegrateVariableLangevinStepKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
1396
     * @param integrator the VariableLangevinIntegrator this kernel will be used for
1397
1398
1399
1400
1401
1402
     */
    void initialize(const System& system, const VariableLangevinIntegrator& integrator);
    /**
     * Execute the kernel.
     *
     * @param context    the context in which to execute this kernel
1403
     * @param integrator the VariableLangevinIntegrator this kernel is being used for
1404
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1405
     * @return the size of the step that was taken
1406
     */
1407
    double execute(ContextImpl& context, const VariableLangevinIntegrator& integrator, double maxTime);
1408
1409
1410
1411
1412
1413
1414
    /**
     * 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);
1415
1416
1417
private:
    ReferencePlatform::PlatformData& data;
    ReferenceVariableStochasticDynamics* dynamics;
peastman's avatar
peastman committed
1418
    std::vector<double> masses;
1419
1420
1421
    double prevTemp, prevFriction, prevErrorTol;
};

1422
1423
1424
1425
1426
1427
/**
 * 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),
1428
        data(data), dynamics(0) {
1429
1430
1431
1432
1433
1434
    }
    ~ReferenceIntegrateVariableVerletStepKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
1435
     * @param integrator the VariableVerletIntegrator this kernel will be used for
1436
1437
1438
1439
1440
1441
     */
    void initialize(const System& system, const VariableVerletIntegrator& integrator);
    /**
     * Execute the kernel.
     *
     * @param context    the context in which to execute this kernel
1442
     * @param integrator the VariableVerletIntegrator this kernel is being used for
1443
     * @param maxTime    the maximum time beyond which the simulation should not be advanced
1444
     * @return the size of the step that was taken
1445
     */
1446
    double execute(ContextImpl& context, const VariableVerletIntegrator& integrator, double maxTime);
1447
1448
1449
1450
1451
1452
1453
    /**
     * 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);
1454
1455
1456
private:
    ReferencePlatform::PlatformData& data;
    ReferenceVariableVerletDynamics* dynamics;
peastman's avatar
peastman committed
1457
    std::vector<double> masses;
1458
    double prevErrorTol;
1459
1460
};

1461
1462
1463
1464
1465
1466
/**
 * 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),
1467
        data(data), dynamics(0) {
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
    }
    ~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);
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
    /**
     * 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);
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
    /**
     * 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
1532
1533
    std::vector<double> masses, globalValues;
    std::vector<std::vector<OpenMM::Vec3> > perDofValues; 
1534
1535
};

1536
/**
Peter Eastman's avatar
Peter Eastman committed
1537
 * This kernel is invoked by AndersenThermostat at the start of each time step to adjust the particle velocities.
1538
1539
1540
 */
class ReferenceApplyAndersenThermostatKernel : public ApplyAndersenThermostatKernel {
public:
1541
    ReferenceApplyAndersenThermostatKernel(std::string name, const Platform& platform) : ApplyAndersenThermostatKernel(name, platform), thermostat(0) {
1542
    }
1543
    ~ReferenceApplyAndersenThermostatKernel();
1544
    /**
1545
     * Initialize the kernel.
1546
     * 
1547
1548
     * @param system     the System this kernel will be applied to
     * @param thermostat the AndersenThermostat this kernel will be used for
1549
     */
1550
    void initialize(const System& system, const AndersenThermostat& thermostat);
1551
1552
1553
    /**
     * Execute the kernel.
     * 
1554
     * @param context    the context in which to execute this kernel
1555
     */
1556
    void execute(ContextImpl& context);
1557
1558
private:
    ReferenceAndersenThermostat* thermostat;
1559
    std::vector<std::vector<int> > particleGroups;
peastman's avatar
peastman committed
1560
    std::vector<double> masses;
1561
1562
};

1563
1564
1565
1566
1567
/**
 * This kernel is invoked by MonteCarloBarostat to adjust the periodic box volume
 */
class ReferenceApplyMonteCarloBarostatKernel : public ApplyMonteCarloBarostatKernel {
public:
1568
    ReferenceApplyMonteCarloBarostatKernel(std::string name, const Platform& platform) : ApplyMonteCarloBarostatKernel(name, platform), barostat(NULL) {
1569
1570
1571
1572
1573
1574
1575
    }
    ~ReferenceApplyMonteCarloBarostatKernel();
    /**
     * Initialize the kernel.
     *
     * @param system     the System this kernel will be applied to
     * @param barostat   the MonteCarloBarostat this kernel will be used for
1576
     * @param rigidMolecules  whether molecules should be kept rigid while scaling coordinates
1577
     */
1578
    void initialize(const System& system, const Force& barostat, bool rigidMolecules=true);
1579
1580
1581
1582
1583
1584
1585
    /**
     * 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);
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
    /**
     * 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
     */
1598
    void scaleCoordinates(ContextImpl& context, double scaleX, double scaleY, double scaleZ);
1599
    /**
1600
1601
     * Reject the most recent Monte Carlo step, restoring the particle positions to where they were when
     * saveCoordinates() was last called.
1602
1603
1604
1605
1606
     *
     * @param context    the context in which to execute this kernel
     */
    void restoreCoordinates(ContextImpl& context);
private:
1607
    bool rigidMolecules;
1608
1609
1610
    ReferenceMonteCarloBarostat* barostat;
};

1611
1612
1613
1614
1615
/**
 * This kernel is invoked to remove center of mass motion from the system.
 */
class ReferenceRemoveCMMotionKernel : public RemoveCMMotionKernel {
public:
1616
    ReferenceRemoveCMMotionKernel(std::string name, const Platform& platform, ReferencePlatform::PlatformData& data) : RemoveCMMotionKernel(name, platform), data(data) {
1617
1618
    }
    /**
Peter Eastman's avatar
Peter Eastman committed
1619
     * Initialize the kernel, setting up the particle masses.
1620
     * 
1621
1622
     * @param system     the System this kernel will be applied to
     * @param force      the CMMotionRemover this kernel will be used for
1623
     */
1624
    void initialize(const System& system, const CMMotionRemover& force);
1625
1626
1627
    /**
     * Execute the kernel.
     * 
1628
     * @param context    the context in which to execute this kernel
1629
     */
1630
    void execute(ContextImpl& context);
1631
private:
1632
    ReferencePlatform::PlatformData& data;
1633
    std::vector<double> masses;
1634
    int frequency;
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
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
/**
 * 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;
};

1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
/**
 * 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;
};

1712
1713
1714
} // namespace OpenMM

#endif /*OPENMM_REFERENCEKERNELS_H_*/