CpuNonbondedForce.h 10.6 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

/* 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__

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
        class ComputeDirectTask;
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
      
         --------------------------------------------------------------------------------------- */
      
90
      void setPeriodic(RealVec* periodicBoxVectors);
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
       
      /**---------------------------------------------------------------------------------------
      
         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
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
            
         --------------------------------------------------------------------------------------- */
          
132
      void calculateReciprocalIxn(int numberOfAtoms, float* posq, const std::vector<RealVec>& atomCoordinates,
peastman's avatar
peastman committed
133
                            const std::vector<std::pair<float, float> >& atomParameters, const std::vector<std::set<int> >& exclusions,
134
                            std::vector<RealVec>& forces, double* totalEnergy) const;
135
136
137
138
139
140
      
      /**---------------------------------------------------------------------------------------
      
         Calculate LJ Coulomb pair ixn
      
         @param numberOfAtoms    number of atoms
peastman's avatar
peastman committed
141
         @param posq             atom coordinates and charges
142
         @param atomCoordinates  atom coordinates (periodic boundary conditions not applied)
peastman's avatar
peastman committed
143
         @param atomParameters   atom parameters (sigma/2, 2*sqrt(epsilon))
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
         @param threads          the thread pool to use
149
150
151
      
         --------------------------------------------------------------------------------------- */
          
152
      void calculateDirectIxn(int numberOfAtoms, float* posq, const std::vector<RealVec>& atomCoordinates, const std::vector<std::pair<float, float> >& atomParameters,
153
            const std::vector<std::set<int> >& exclusions, std::vector<AlignedArray<float> >& threadForce, double* totalEnergy, ThreadPool& threads);
154
155
156
157

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

160
protected:
161
162
163
        bool cutoff;
        bool useSwitch;
        bool periodic;
164
        bool triclinic;
165
166
        bool ewald;
        bool pme;
167
        bool tableIsValid;
168
        const CpuNeighborList* neighborList;
169
170
171
        float recipBoxSize[3];
        RealVec periodicBoxVectors[3];
        AlignedArray<fvec4> periodicBoxVec4;
172
173
174
175
176
        float cutoffDistance, switchingDistance;
        float krf, crf;
        float alphaEwald;
        int numRx, numRy, numRz;
        int meshDim[3];
177
        std::vector<float> ewaldScaleTable;
peastman's avatar
peastman committed
178
        float ewaldDX, ewaldDXInv;
179
        std::vector<double> threadEnergy;
180
        // The following variables are used to make information accessible to the individual threads.
peastman's avatar
peastman committed
181
        int numberOfAtoms;
182
        float* posq;
183
        RealVec const* atomCoordinates;
peastman's avatar
peastman committed
184
185
        std::pair<float, float> const* atomParameters;        
        std::set<int> const* exclusions;
186
        std::vector<AlignedArray<float> >* threadForce;
187
        bool includeEnergy;
188
        void* atomicCounter;
189

peastman's avatar
peastman committed
190
191
        static const float TWO_OVER_SQRT_PI;
        static const int NUM_TABLE_POINTS;
peastman's avatar
peastman committed
192
            
193
194
      /**---------------------------------------------------------------------------------------
      
peastman's avatar
peastman committed
195
         Calculate LJ Coulomb pair ixn between two atoms
196
      
peastman's avatar
peastman committed
197
198
         @param atom1            the index of the first atom
         @param atom2            the index of the second atom
199
200
201
202
203
         @param forces           force array (forces added)
         @param totalEnergy      total energy
            
         --------------------------------------------------------------------------------------- */
          
204
      void calculateOneIxn(int atom1, int atom2, float* forces, double* totalEnergy, const fvec4& boxSize, const fvec4& invBoxSize);
peastman's avatar
peastman committed
205
206
207
            
      /**---------------------------------------------------------------------------------------
      
208
         Calculate all the interactions for one atom block.
peastman's avatar
peastman committed
209
      
210
211
212
213
214
215
         @param blockIndex       the index of the atom block
         @param forces           force array (forces added)
         @param totalEnergy      total energy
            
         --------------------------------------------------------------------------------------- */
          
216
      virtual void calculateBlockIxn(int blockIndex, float* forces, double* totalEnergy, const fvec4& boxSize, const fvec4& invBoxSize) = 0;
217
218
219
220
221
222
            
      /**---------------------------------------------------------------------------------------
      
         Calculate all the interactions for one atom block.
      
         @param blockIndex       the index of the atom block
peastman's avatar
peastman committed
223
224
225
226
227
         @param forces           force array (forces added)
         @param totalEnergy      total energy
            
         --------------------------------------------------------------------------------------- */
          
228
      virtual void calculateBlockEwaldIxn(int blockIndex, float* forces, double* totalEnergy, const fvec4& boxSize, const fvec4& invBoxSize) = 0;
peastman's avatar
peastman committed
229

peastman's avatar
peastman committed
230
231
232
233
      /**
       * Compute the displacement and squared distance between two points, optionally using
       * periodic boundary conditions.
       */
234
      void getDeltaR(const fvec4& posI, const fvec4& posJ, fvec4& deltaR, float& r2, bool periodic, const fvec4& boxSize, const fvec4& invBoxSize) const;
235

peastman's avatar
peastman committed
236
237
238
239
      /**
       * Create a lookup table for the scale factor used with Ewald and PME.
       */
      void tabulateEwaldScaleFactor();
240

peastman's avatar
peastman committed
241
      /**
242
       * Compute a fast approximation to erfc(x).
peastman's avatar
peastman committed
243
       */
244
      static float erfcApprox(float x);
245
246
};

247
248
} // namespace OpenMM

249
250
251
// ---------------------------------------------------------------------------------------

#endif // OPENMM_CPU_NONBONDED_FORCE_H__