CpuNonbondedForce.h 12 KB
Newer Older
1

peastman's avatar
peastman committed
2
/* Portions copyright (c) 2006-2018 Stanford University and Simbios.
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
 * Contributors: Pande Group
 *
 * 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.
 */

#ifndef OPENMM_CPU_NONBONDED_FORCE_H__
#define OPENMM_CPU_NONBONDED_FORCE_H__

28
#include "AlignedArray.h"
29
#include "CpuNeighborList.h"
30
#include "ReferencePairIxn.h"
31
#include "openmm/internal/ThreadPool.h"
32
#include "openmm/internal/vectorize.h"
peastman's avatar
peastman committed
33
#include <atomic>
34
35
36
37
38
#include <set>
#include <utility>
#include <vector>
// ---------------------------------------------------------------------------------------

39
40
namespace OpenMM {

41
class CpuNonbondedForce {
42
    public:
43
44
45
46
47
48
49

      /**---------------------------------------------------------------------------------------
      
         Constructor
      
         --------------------------------------------------------------------------------------- */

50
       CpuNonbondedForce();
51
52
53
54
       
        /**
         * Virtual destructor.
         */
55

56
57
        virtual ~CpuNonbondedForce();
        
58
59
60
61
62
63
64
65
66
67
      /**---------------------------------------------------------------------------------------
      
         Set the force to use a cutoff.
      
         @param distance            the cutoff distance
         @param neighbors           the neighbor list to use
         @param solventDielectric   the dielectric constant of the bulk solvent
      
         --------------------------------------------------------------------------------------- */
      
68
      void setUseCutoff(float distance, const CpuNeighborList& neighbors, float solventDielectric);
69
70
71
72
73
74
75
76
77

      /**---------------------------------------------------------------------------------------
      
         Set the force to use a switching function on the Lennard-Jones interaction.
      
         @param distance            the switching distance
      
         --------------------------------------------------------------------------------------- */
      
78
      void setUseSwitchingFunction(float distance);
79
80
81
82
83
84
85
      
      /**---------------------------------------------------------------------------------------
      
         Set the force to use periodic boundary conditions.  This requires that a cutoff has
         already been set, and the smallest side of the periodic box is at least twice the cutoff
         distance.
      
86
         @param periodicBoxVectors    the vectors defining the periodic box
87
88
89
      
         --------------------------------------------------------------------------------------- */
      
peastman's avatar
peastman committed
90
      void setPeriodic(Vec3* periodicBoxVectors);
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
       
      /**---------------------------------------------------------------------------------------
      
         Set the force to use Ewald summation.
      
         @param alpha  the Ewald separation parameter
         @param kmaxx  the largest wave vector in the x direction
         @param kmaxy  the largest wave vector in the y direction
         @param kmaxz  the largest wave vector in the z direction
      
         --------------------------------------------------------------------------------------- */
      
      void setUseEwald(float alpha, int kmaxx, int kmaxy, int kmaxz);

     
      /**---------------------------------------------------------------------------------------
107

108
         Set the force to use Particle-Mesh Ewald (PME) summation.
109

110
111
         @param alpha    the Ewald separation parameter
         @param gridSize the dimensions of the mesh
112

113
         --------------------------------------------------------------------------------------- */
114

115
      void setUsePME(float alpha, int meshSize[3]);
peastman's avatar
peastman committed
116

117
118
119
120
121
122
123
124
125
126
127
      /**---------------------------------------------------------------------------------------

         Set the force to use Particle-Mesh Ewald (PME) summation for dispersion.

         @param alpha    the Ewald separation parameter
         @param gridSize the dimensions of the mesh

         --------------------------------------------------------------------------------------- */

      void setUseLJPME(float alpha, int meshSize[3]);

peastman's avatar
peastman committed
128
129
130
131
132
133
134
135
      /**---------------------------------------------------------------------------------------
      
         Calculate Ewald ixn
      
         @param numberOfAtoms    number of atoms
         @param posq             atom coordinates and charges
         @param atomCoordinates  atom coordinates (in format needed by PME)
         @param atomParameters   atom parameters (sigma/2, 2*sqrt(epsilon))
136
         @param C6Paramrs        C6 parameters for multiplicative representation of dispersion
peastman's avatar
peastman committed
137
138
139
140
141
142
         @param exclusions       atom exclusion indices
                                 exclusions[atomIndex] contains the list of exclusions for that atom
         @param forces           force array (forces added)
         @param totalEnergy      total energy
            
         --------------------------------------------------------------------------------------- */
143

peastman's avatar
peastman committed
144
      void calculateReciprocalIxn(int numberOfAtoms, float* posq, const std::vector<Vec3>& atomCoordinates,
145
                                  const std::vector<std::pair<float, float> >& atomParameters, const std::vector<float> &C6params,
peastman's avatar
peastman committed
146
                                  const std::vector<std::set<int> >& exclusions, std::vector<Vec3>& forces, double* totalEnergy) const;
147
148
149
150
151
152
      
      /**---------------------------------------------------------------------------------------
      
         Calculate LJ Coulomb pair ixn
      
         @param numberOfAtoms    number of atoms
peastman's avatar
peastman committed
153
         @param posq             atom coordinates and charges
154
         @param atomCoordinates  atom coordinates (periodic boundary conditions not applied)
peastman's avatar
peastman committed
155
         @param atomParameters   atom parameters (sigma/2, 2*sqrt(epsilon))
156
157
158
159
         @param exclusions       atom exclusion indices
                                 exclusions[atomIndex] contains the list of exclusions for that atom
         @param forces           force array (forces added)
         @param totalEnergy      total energy
160
         @param threads          the thread pool to use
161
162
163
      
         --------------------------------------------------------------------------------------- */
          
peastman's avatar
peastman committed
164
      void calculateDirectIxn(int numberOfAtoms, float* posq, const std::vector<Vec3>& atomCoordinates, const std::vector<std::pair<float, float> >& atomParameters,
165
            const std::vector<float>& C6params, const std::vector<std::set<int> >& exclusions, std::vector<AlignedArray<float> >& threadForce, double* totalEnergy, ThreadPool& threads);
166
167
168
169

    /**
     * This routine contains the code executed by each thread.
     */
170
    void threadComputeDirect(ThreadPool& threads, int threadIndex);
171

172
protected:
173
174
175
        bool cutoff;
        bool useSwitch;
        bool periodic;
176
        bool triclinic;
177
        bool ewald;
178
179
        bool ljpme, pme;
        bool tableIsValid, expTableIsValid;
180
        const CpuNeighborList* neighborList;
181
        float recipBoxSize[3];
peastman's avatar
peastman committed
182
        Vec3 periodicBoxVectors[3];
183
        AlignedArray<fvec4> periodicBoxVec4;
184
185
        float cutoffDistance, switchingDistance;
        float krf, crf;
186
        float alphaEwald, alphaDispersionEwald;
187
        int numRx, numRy, numRz;
188
        int meshDim[3], dispersionMeshDim[3];
189
        std::vector<float> erfcTable, ewaldScaleTable;
190
191
        std::vector<float> exptermsTable, dExptermsTable;
        float ewaldDX, ewaldDXInv, erfcDXInv, exptermsDX, exptermsDXInv;
192
        std::vector<double> threadEnergy;
193
        // The following variables are used to make information accessible to the individual threads.
peastman's avatar
peastman committed
194
        int numberOfAtoms;
195
        float* posq;
peastman's avatar
peastman committed
196
        Vec3 const* atomCoordinates;
peastman's avatar
peastman committed
197
        std::pair<float, float> const* atomParameters;        
198
        float const *C6params;
peastman's avatar
peastman committed
199
        std::set<int> const* exclusions;
200
        std::vector<AlignedArray<float> >* threadForce;
201
        bool includeEnergy;
202
203
        float inverseRcut6;
        float inverseRcut6Expterm;
peastman's avatar
peastman committed
204
        std::atomic<int> atomicCounter;
205

peastman's avatar
peastman committed
206
207
        static const float TWO_OVER_SQRT_PI;
        static const int NUM_TABLE_POINTS;
208

209
210
      /**---------------------------------------------------------------------------------------
      
peastman's avatar
peastman committed
211
         Calculate LJ Coulomb pair ixn between two atoms
212
      
peastman's avatar
peastman committed
213
214
         @param atom1            the index of the first atom
         @param atom2            the index of the second atom
215
216
217
218
219
         @param forces           force array (forces added)
         @param totalEnergy      total energy
            
         --------------------------------------------------------------------------------------- */
          
220
      void calculateOneIxn(int atom1, int atom2, float* forces, double* totalEnergy, const fvec4& boxSize, const fvec4& invBoxSize);
peastman's avatar
peastman committed
221
222
223
            
      /**---------------------------------------------------------------------------------------
      
224
         Calculate all the interactions for one atom block.
peastman's avatar
peastman committed
225
      
226
227
228
229
230
231
         @param blockIndex       the index of the atom block
         @param forces           force array (forces added)
         @param totalEnergy      total energy
            
         --------------------------------------------------------------------------------------- */
          
232
      virtual void calculateBlockIxn(int blockIndex, float* forces, double* totalEnergy, const fvec4& boxSize, const fvec4& invBoxSize) = 0;
233
234
235
236
237
238
            
      /**---------------------------------------------------------------------------------------
      
         Calculate all the interactions for one atom block.
      
         @param blockIndex       the index of the atom block
peastman's avatar
peastman committed
239
240
241
242
243
         @param forces           force array (forces added)
         @param totalEnergy      total energy
            
         --------------------------------------------------------------------------------------- */
          
244
      virtual void calculateBlockEwaldIxn(int blockIndex, float* forces, double* totalEnergy, const fvec4& boxSize, const fvec4& invBoxSize) = 0;
peastman's avatar
peastman committed
245

peastman's avatar
peastman committed
246
247
248
249
      /**
       * Compute the displacement and squared distance between two points, optionally using
       * periodic boundary conditions.
       */
250
      void getDeltaR(const fvec4& posI, const fvec4& posJ, fvec4& deltaR, float& r2, bool periodic, const fvec4& boxSize, const fvec4& invBoxSize) const;
251

peastman's avatar
peastman committed
252
253
254
255
      /**
       * Create a lookup table for the scale factor used with Ewald and PME.
       */
      void tabulateEwaldScaleFactor();
256

257
258
259
260
261
      /**
       * Create a lookup table for the scale factor used with dispersion PME.
       */
      void tabulateExpTerms();

peastman's avatar
peastman committed
262
      /**
263
       * Compute a fast approximation to erfc(x).
peastman's avatar
peastman committed
264
       */
265
      float erfcApprox(float x);
266
267
268
269
270
271
272
273
274
275
276
277
278
279

      /**
       * Compute a fast approximation to (1.0 - EXP(-dar^2) * (1.0 + dar^2 + 0.5*dar^4))
       * where dar = (dispersionAlpha * R)
       * needed for LJPME energies.
       */
      float exptermsApprox(float R);

      /**
       * Compute a fast approximation to (1.0 - EXP(-dar^2) * (1.0 + dar^2 + 0.5*dar^4 + dar^6/6.0))
       * where dar = (dispersionAlpha * R)
       * needed for LJPME forces.
       */
      float dExptermsApprox(float R);
280
281
};

282
283
} // namespace OpenMM

284
285
286
// ---------------------------------------------------------------------------------------

#endif // OPENMM_CPU_NONBONDED_FORCE_H__