ReferenceLJCoulombIxn.h 8.54 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 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 __ReferenceLJCoulombIxn_H__
#define __ReferenceLJCoulombIxn_H__

#include "ReferencePairIxn.h"
29
#include "ReferenceNeighborList.h"
30
31
32
33
34
35

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

class ReferenceLJCoulombIxn : public ReferencePairIxn {

   private:
36
37
38
       
      bool cutoff;
      bool periodic;
39
      bool ewald;
40
      bool pme;
41
42
43
44
      const OpenMM::NeighborList* neighborList;
      RealOpenMM periodicBoxSize[3];
      RealOpenMM cutoffDistance;
      RealOpenMM krf, crf;
45
46
      RealOpenMM alphaEwald;
      int numRx, numRy, numRz;
47
      int meshDim[3];
48
49
50
51
52
53

      // parameter indices

      static const int SigIndex = 0;
      static const int EpsIndex = 1;
      static const int   QIndex = 2;
54
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 energyByAtom     atom energy
         @param totalEnergy      total energy
            
         --------------------------------------------------------------------------------------- */
          
69
      void calculateOneIxn( int atom1, int atom2, std::vector<OpenMM::RealVec>& atomCoordinates,
70
                            RealOpenMM** atomParameters, std::vector<OpenMM::RealVec>& forces,
71
72
                            RealOpenMM* energyByAtom, RealOpenMM* totalEnergy ) const;

73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91

   public:

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

       ReferenceLJCoulombIxn( );

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

       ~ReferenceLJCoulombIxn( );

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
      
         --------------------------------------------------------------------------------------- */
      
102
      void setUseCutoff( RealOpenMM distance, const OpenMM::NeighborList& neighbors, RealOpenMM solventDielectric );
103
104
105
106
107
108
109
110
111
112
113
      
      /**---------------------------------------------------------------------------------------
      
         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
      
         --------------------------------------------------------------------------------------- */
      
114
      void setPeriodic( OpenMM::RealVec& boxSize );
115
       
116
117
      /**---------------------------------------------------------------------------------------
      
118
         Set the force to use Ewald summation.
119
      
120
121
122
123
         @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
124
125
126
      
         --------------------------------------------------------------------------------------- */
      
127
      void setUseEwald(RealOpenMM alpha, int kmaxx, int kmaxy, int kmaxz);
128

129
130
131
132
133
     
      /**---------------------------------------------------------------------------------------
      
         Set the force to use Particle-Mesh Ewald (PME) summation.
      
134
135
         @param alpha    the Ewald separation parameter
         @param gridSize the dimensions of the mesh
136
137
138
      
         --------------------------------------------------------------------------------------- */
      
139
      void setUsePME(RealOpenMM alpha, int meshSize[3]);
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
      
      /**---------------------------------------------------------------------------------------
      
         Calculate LJ Coulomb pair ixn
      
         @param numberOfAtoms    number of atoms
         @param atomCoordinates  atom coordinates
         @param atomParameters   atom parameters (charges, c6, c12, ...)     atomParameters[atomIndex][paramterIndex]
         @param exclusions       atom exclusion indices                      exclusions[atomIndex][atomToExcludeIndex]
                                 exclusions[atomIndex][0] = number of exclusions
                                 exclusions[atomIndex][1-no.] = atom indices of atoms to excluded from
                                 interacting w/ atom atomIndex
         @param fixedParameters  non atom parameters (not currently used)
         @param forces           force array (forces added)
         @param energyByAtom     atom energy
         @param totalEnergy      total energy
      
         --------------------------------------------------------------------------------------- */
          
159
      void calculatePairIxn( int numberOfAtoms, std::vector<OpenMM::RealVec>& atomCoordinates,
160
                            RealOpenMM** atomParameters, int** exclusions,
161
                            RealOpenMM* fixedParameters, std::vector<OpenMM::RealVec>& forces,
162
                            RealOpenMM* energyByAtom, RealOpenMM* totalEnergy ) const;
163
164

private:
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
      /**---------------------------------------------------------------------------------------
      
         Calculate Ewald ixn
      
         @param numberOfAtoms    number of atoms
         @param atomCoordinates  atom coordinates
         @param atomParameters   atom parameters (charges, c6, c12, ...)     atomParameters[atomIndex][paramterIndex]
         @param exclusions       atom exclusion indices                      exclusions[atomIndex][atomToExcludeIndex]
                                 exclusions[atomIndex][0] = number of exclusions
                                 exclusions[atomIndex][1-no.] = atom indices of atoms to excluded from
                                 interacting w/ atom atomIndex
         @param fixedParameters  non atom parameters (not currently used)
         @param forces           force array (forces added)
         @param energyByAtom     atom energy
         @param totalEnergy      total energy
180
181
182
            
         --------------------------------------------------------------------------------------- */
          
183
      void calculateEwaldIxn( int numberOfAtoms, std::vector<OpenMM::RealVec>& atomCoordinates,
184
                            RealOpenMM** atomParameters, int** exclusions,
185
                            RealOpenMM* fixedParameters, std::vector<OpenMM::RealVec>& forces,
186
                            RealOpenMM* energyByAtom, RealOpenMM* totalEnergy ) const;
187
188
189
190
191
};

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

#endif // __ReferenceLJCoulombIxn_H__