ReferenceLJCoulombIxn.h 9.25 KB
Newer Older
1

2
/* Portions copyright (c) 2006-2020 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
28
 * 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 __ReferenceLJCoulombIxn_H__
#define __ReferenceLJCoulombIxn_H__

#include "ReferencePairIxn.h"
29
#include "ReferenceNeighborList.h"
30

31
namespace OpenMM {
32

33
class ReferenceLJCoulombIxn {
34
35

   private:
36
37
       
      bool cutoff;
38
      bool useSwitch;
39
      bool periodic, periodicExceptions;
40
      bool ewald;
41
      bool pme, ljpme;
42
      const OpenMM::NeighborList* neighborList;
peastman's avatar
peastman committed
43
44
45
      OpenMM::Vec3 periodicBoxVectors[3];
      double cutoffDistance, switchingDistance;
      double krf, crf;
peastman's avatar
peastman committed
46
      double alphaEwald, alphaDispersionEwald;
47
      int numRx, numRy, numRz;
48
      int meshDim[3], dispersionMeshDim[3];
49
50
51
52
53
54

      // parameter indices

      static const int SigIndex = 0;
      static const int EpsIndex = 1;
      static const int   QIndex = 2;
55
56
57
58
59
60
61
62
63
64
65
66
67
68
            
      /**---------------------------------------------------------------------------------------
      
         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 atomCoordinates  atom coordinates
         @param atomParameters   atom parameters (charges, c6, c12, ...)     atomParameters[atomIndex][paramterIndex]
         @param forces           force array (forces added)
         @param totalEnergy      total energy
            
         --------------------------------------------------------------------------------------- */
          
peastman's avatar
peastman committed
69
      void calculateOneIxn(int atom1, int atom2, std::vector<OpenMM::Vec3>& atomCoordinates,
70
71
                           std::vector<std::vector<double> >& atomParameters, std::vector<OpenMM::Vec3>& forces,
                           double* totalEnergy) const;
72

73
74
75
76
77
78
79
80
81

   public:

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

82
       ReferenceLJCoulombIxn();
83
84
85
86
87
88
89

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

90
       ~ReferenceLJCoulombIxn();
91

92
93
94
95
96
97
98
99
100
101
      /**---------------------------------------------------------------------------------------
      
         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
      
         --------------------------------------------------------------------------------------- */
      
peastman's avatar
peastman committed
102
      void setUseCutoff(double distance, const OpenMM::NeighborList& neighbors, double solventDielectric);
103
104
105
106
107
108
109
110
111

      /**---------------------------------------------------------------------------------------
      
         Set the force to use a switching function on the Lennard-Jones interaction.
      
         @param distance            the switching distance
      
         --------------------------------------------------------------------------------------- */
      
peastman's avatar
peastman committed
112
      void setUseSwitchingFunction(double distance);
113
114
115
116
117
118
119
      
      /**---------------------------------------------------------------------------------------
      
         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.
      
120
         @param vectors    the vectors defining the periodic box
121
122
123
      
         --------------------------------------------------------------------------------------- */
      
peastman's avatar
peastman committed
124
      void setPeriodic(OpenMM::Vec3* vectors);
125
       
126
127
      /**---------------------------------------------------------------------------------------
      
128
         Set the force to use Ewald summation.
129
      
130
131
132
133
         @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
134
135
136
      
         --------------------------------------------------------------------------------------- */
      
peastman's avatar
peastman committed
137
      void setUseEwald(double alpha, int kmaxx, int kmaxy, int kmaxz);
138

139
140
     
      /**---------------------------------------------------------------------------------------
141

142
         Set the force to use Particle-Mesh Ewald (PME) summation.
143

144
145
         @param alpha    the Ewald separation parameter
         @param gridSize the dimensions of the mesh
146

147
148
         --------------------------------------------------------------------------------------- */
      
peastman's avatar
peastman committed
149
      void setUsePME(double alpha, int meshSize[3]);
150
      
151
152
153
154
155
156
157
158
159
      /**---------------------------------------------------------------------------------------

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

         @param dalpha    the dispersion Ewald separation parameter
         @param dgridSize the dimensions of the dispersion mesh

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

peastman's avatar
peastman committed
160
      void setUseLJPME(double dalpha, int dmeshSize[3]);
161
162
163
164
165
166
167
168
      
      /**---------------------------------------------------------------------------------------

         Set whether exceptions use periodic boundary conditions.

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

      void setPeriodicExceptions(bool periodic);
169

170
171
172
173
174
175
176
      /**---------------------------------------------------------------------------------------
      
         Calculate LJ Coulomb pair ixn
      
         @param numberOfAtoms    number of atoms
         @param atomCoordinates  atom coordinates
         @param atomParameters   atom parameters (charges, c6, c12, ...)     atomParameters[atomIndex][paramterIndex]
177
178
         @param exclusions       atom exclusion indices
                                 exclusions[atomIndex] contains the list of exclusions for that atom
179
180
         @param forces           force array (forces added)
         @param totalEnergy      total energy
181
182
         @param includeDirect      true if direct space interactions should be included
         @param includeReciprocal  true if reciprocal space interactions should be included
183
184
185
      
         --------------------------------------------------------------------------------------- */
          
peastman's avatar
peastman committed
186
      void calculatePairIxn(int numberOfAtoms, std::vector<OpenMM::Vec3>& atomCoordinates,
187
188
                            std::vector<std::vector<double> >& atomParameters, std::vector<std::set<int> >& exclusions,
                            std::vector<OpenMM::Vec3>& forces, double* totalEnergy, bool includeDirect, bool includeReciprocal) const;
189
190

private:
191
192
193
194
195
196
197
      /**---------------------------------------------------------------------------------------
      
         Calculate Ewald ixn
      
         @param numberOfAtoms    number of atoms
         @param atomCoordinates  atom coordinates
         @param atomParameters   atom parameters (charges, c6, c12, ...)     atomParameters[atomIndex][paramterIndex]
198
199
         @param exclusions       atom exclusion indices
                                 exclusions[atomIndex] contains the list of exclusions for that atom
200
201
         @param forces           force array (forces added)
         @param totalEnergy      total energy
202
203
         @param includeDirect      true if direct space interactions should be included
         @param includeReciprocal  true if reciprocal space interactions should be included
204
205
206
            
         --------------------------------------------------------------------------------------- */
          
peastman's avatar
peastman committed
207
      void calculateEwaldIxn(int numberOfAtoms, std::vector<OpenMM::Vec3>& atomCoordinates,
208
209
                             std::vector<std::vector<double> >& atomParameters, std::vector<std::set<int> >& exclusions,
                             std::vector<OpenMM::Vec3>& forces, double* totalEnergy, bool includeDirect, bool includeReciprocal) const;
210
211
};

212
} // namespace OpenMM
213
214

#endif // __ReferenceLJCoulombIxn_H__