CpuNonbondedForce.h 12 KB
Newer Older
1

peastman's avatar
peastman committed
2
/* Portions copyright (c) 2006-2017 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"
33
34
35
36
37
#include <set>
#include <utility>
#include <vector>
// ---------------------------------------------------------------------------------------

38
39
namespace OpenMM {

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

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

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

55
56
        virtual ~CpuNonbondedForce();
        
57
58
59
60
61
62
63
64
65
66
      /**---------------------------------------------------------------------------------------
      
         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
      
         --------------------------------------------------------------------------------------- */
      
67
      void setUseCutoff(float distance, const CpuNeighborList& neighbors, float solventDielectric);
68
69
70
71
72
73
74
75
76

      /**---------------------------------------------------------------------------------------
      
         Set the force to use a switching function on the Lennard-Jones interaction.
      
         @param distance            the switching distance
      
         --------------------------------------------------------------------------------------- */
      
77
      void setUseSwitchingFunction(float distance);
78
79
80
81
82
83
84
      
      /**---------------------------------------------------------------------------------------
      
         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.
      
85
         @param periodicBoxVectors    the vectors defining the periodic box
86
87
88
      
         --------------------------------------------------------------------------------------- */
      
peastman's avatar
peastman committed
89
      void setPeriodic(Vec3* periodicBoxVectors);
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
       
      /**---------------------------------------------------------------------------------------
      
         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);

     
      /**---------------------------------------------------------------------------------------
106

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

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

112
         --------------------------------------------------------------------------------------- */
113

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

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

         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
127
128
129
130
131
132
133
134
      /**---------------------------------------------------------------------------------------
      
         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))
135
         @param C6Paramrs        C6 parameters for multiplicative representation of dispersion
peastman's avatar
peastman committed
136
137
138
139
140
141
         @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
            
         --------------------------------------------------------------------------------------- */
142

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

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

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

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

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

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

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

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

      /**
       * 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);
279
280
};

281
282
} // namespace OpenMM

283
284
285
// ---------------------------------------------------------------------------------------

#endif // OPENMM_CPU_NONBONDED_FORCE_H__