CpuNonbondedForce.h 10 KB
Newer Older
1
2
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
28

/* Portions copyright (c) 2006-2013 Stanford University and Simbios.
 * 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__

#include "ReferencePairIxn.h"
29
#include <pthread.h>
30
31
32
#include <set>
#include <utility>
#include <vector>
33
#include <smmintrin.h>
34
35
36
// ---------------------------------------------------------------------------------------

class CpuNonbondedForce {
37
38
    public:
        class ThreadData;
39
40
41
42
43
44
45

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

46
       CpuNonbondedForce();
47
48
49
50
51
52
53

      /**---------------------------------------------------------------------------------------
      
         Destructor
      
         --------------------------------------------------------------------------------------- */

54
       ~CpuNonbondedForce();
55
56
57
58
59
60
61
62
63
64
65

      /**---------------------------------------------------------------------------------------
      
         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
      
         --------------------------------------------------------------------------------------- */
      
66
      void setUseCutoff(float distance, const std::vector<std::pair<int, int> >& neighbors, float solventDielectric);
67
68
69
70
71
72
73
74
75

      /**---------------------------------------------------------------------------------------
      
         Set the force to use a switching function on the Lennard-Jones interaction.
      
         @param distance            the switching distance
      
         --------------------------------------------------------------------------------------- */
      
76
      void setUseSwitchingFunction(float distance);
77
78
79
80
81
82
83
84
85
86
87
      
      /**---------------------------------------------------------------------------------------
      
         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.
      
         @param boxSize             the X, Y, and Z widths of the periodic box
      
         --------------------------------------------------------------------------------------- */
      
88
      void setPeriodic(float* periodicBoxSize);
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
       
      /**---------------------------------------------------------------------------------------
      
         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);

     
      /**---------------------------------------------------------------------------------------
      
         Set the force to use Particle-Mesh Ewald (PME) summation.
      
         @param alpha    the Ewald separation parameter
         @param gridSize the dimensions of the mesh
      
         --------------------------------------------------------------------------------------- */
      
      void setUsePME(float alpha, int meshSize[3]);
peastman's avatar
peastman committed
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131

      /**---------------------------------------------------------------------------------------
      
         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))
         @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
            
         --------------------------------------------------------------------------------------- */
          
      void calculateReciprocalIxn(int numberOfAtoms, float* posq, std::vector<OpenMM::RealVec>& atomCoordinates,
                            const std::vector<std::pair<float, float> >& atomParameters, const std::vector<std::set<int> >& exclusions,
132
                            std::vector<OpenMM::RealVec>& forces, float* totalEnergy) const;
133
134
135
136
137
138
      
      /**---------------------------------------------------------------------------------------
      
         Calculate LJ Coulomb pair ixn
      
         @param numberOfAtoms    number of atoms
peastman's avatar
peastman committed
139
140
         @param posq             atom coordinates and charges
         @param atomParameters   atom parameters (sigma/2, 2*sqrt(epsilon))
141
142
143
144
145
146
147
         @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
      
         --------------------------------------------------------------------------------------- */
          
148
149
      void calculateDirectIxn(int numberOfAtoms, float* posq, const std::vector<std::pair<float, float> >& atomParameters,
            const std::vector<std::set<int> >& exclusions, float* forces, float* totalEnergy);
150
151
152
153
154

    /**
     * This routine contains the code executed by each thread.
     */
    void runThread(int index, std::vector<float>& threadForce, double& threadEnergy);
155
156

private:
157
158
159
160
161
162
163
164
165
166
167
168
        bool cutoff;
        bool useSwitch;
        bool periodic;
        bool ewald;
        bool pme;
        const std::vector<std::pair<int, int> >* neighborList;
        float periodicBoxSize[3];
        float cutoffDistance, switchingDistance;
        float krf, crf;
        float alphaEwald;
        int numRx, numRy, numRz;
        int meshDim[3];
peastman's avatar
peastman committed
169
170
        std::vector<float> ewaldScaleX, ewaldScaleY, ewaldScaleDeriv;
        float ewaldDX, ewaldDXInv;
171
172
173
174
175
176
177
        bool isDeleted;
        int numThreads, waitCount;
        std::vector<pthread_t> thread;
        std::vector<ThreadData*> threadData;
        pthread_cond_t startCondition, endCondition;
        pthread_mutex_t lock;
        // The following variables are used to make information accessible to the individual threads.
peastman's avatar
peastman committed
178
        int numberOfAtoms;
179
        float* posq;
peastman's avatar
peastman committed
180
181
        std::pair<float, float> const* atomParameters;        
        std::set<int> const* exclusions;
182
183
        bool includeEnergy;

peastman's avatar
peastman committed
184
185
        static const float TWO_OVER_SQRT_PI;
        static const int NUM_TABLE_POINTS;
peastman's avatar
peastman committed
186
            
187
188
      /**---------------------------------------------------------------------------------------
      
peastman's avatar
peastman committed
189
         Calculate LJ Coulomb pair ixn between two atoms
190
      
peastman's avatar
peastman committed
191
192
         @param atom1            the index of the first atom
         @param atom2            the index of the second atom
193
194
195
196
197
         @param forces           force array (forces added)
         @param totalEnergy      total energy
            
         --------------------------------------------------------------------------------------- */
          
peastman's avatar
peastman committed
198
      void calculateOneIxn(int atom1, int atom2, float* forces, double* totalEnergy, const __m128& boxSize, const __m128& invBoxSize);
peastman's avatar
peastman committed
199
200
201
202
203
204
205
206
207
208
209
210
            
      /**---------------------------------------------------------------------------------------
      
         Calculate LJ Coulomb pair ixn between two atoms
      
         @param atom1            the index of the first atom
         @param atom2            the index of the second atom
         @param forces           force array (forces added)
         @param totalEnergy      total energy
            
         --------------------------------------------------------------------------------------- */
          
peastman's avatar
peastman committed
211
      void calculateOneEwaldIxn(int atom1, int atom2, float* forces, double* totalEnergy, const __m128& boxSize, const __m128& invBoxSize);
peastman's avatar
peastman committed
212

peastman's avatar
peastman committed
213
214
215
216
      /**
       * Compute the displacement and squared distance between two points, optionally using
       * periodic boundary conditions.
       */
peastman's avatar
peastman committed
217
      void getDeltaR(const __m128& posI, const __m128& posJ, __m128& deltaR, float& r2, bool periodic, const __m128& boxSize, const __m128& invBoxSize) const;
218

peastman's avatar
peastman committed
219
220
221
      /**
       * Compute a fast approximation to erfc(x).
       */
peastman's avatar
peastman committed
222
      static float erfcApprox(float x);
peastman's avatar
peastman committed
223
224
225
226
227
228
229
230
231
232

      /**
       * Create a lookup table for the scale factor used with Ewald and PME.
       */
      void tabulateEwaldScaleFactor();
      
      /**
       * Evaluate the scale factor used with Ewald and PME: erfc(alpha*r) + 2*alpha*r*exp(-alpha*alpha*r*r)/sqrt(PI)
       */
      float ewaldScaleFunction(float x);
233
234
235
236
237
};

// ---------------------------------------------------------------------------------------

#endif // OPENMM_CPU_NONBONDED_FORCE_H__