"vscode:/vscode.git/clone" did not exist on "ac37b2e76cbedb4c233a01809268124be9e4eb0c"
ReferenceLJCoulombIxn.h 9.18 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
41
42
43
      const OpenMM::NeighborList* neighborList;
      RealOpenMM periodicBoxSize[3];
      RealOpenMM cutoffDistance;
      RealOpenMM krf, crf;
44
45
      RealOpenMM alphaEwald;
      int numRx, numRy, numRz;
46
47
48
49
50
51

      // parameter indices

      static const int SigIndex = 0;
      static const int EpsIndex = 1;
      static const int   QIndex = 2;
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
            
      /**---------------------------------------------------------------------------------------
      
         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
      
         @return ReferenceForce::DefaultReturn
            
         --------------------------------------------------------------------------------------- */
          
      int calculateOneIxn( int atom1, int atom2, RealOpenMM** atomCoordinates,
                            RealOpenMM** atomParameters, RealOpenMM** forces,
                            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
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
      /**---------------------------------------------------------------------------------------
      
         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
      
         @return ReferenceForce::DefaultReturn
      
         --------------------------------------------------------------------------------------- */
      
      int setUseCutoff( RealOpenMM distance, const OpenMM::NeighborList& neighbors, RealOpenMM solventDielectric );
      
      /**---------------------------------------------------------------------------------------
      
         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
      
         @return ReferenceForce::DefaultReturn
      
         --------------------------------------------------------------------------------------- */
      
      int setPeriodic( RealOpenMM* boxSize );      
119
120
121
      
      /**---------------------------------------------------------------------------------------
      
122
         Set the force to use Ewald summation.
123
      
124
125
126
127
         @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
128
129
130
      
         --------------------------------------------------------------------------------------- */
      
131
      void setUseEwald(RealOpenMM alpha, int kmaxx, int kmaxy, int kmaxz);
132

133

134
135
136
137
138
139
140
      /**---------------------------------------------------------------------------------------
      
         Calculate parameters for LJ 1-4 ixn
      
         @param c6               c6
         @param c12              c12
         @param q1               q1 charge atom
141
         @param epsfacSqrt       epsfacSqrt (what is this?)
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
         @param parameters       output parameters:
                                    parameter[SigIndex]  = sqrt(c6*c6/c12)
                                    parameter[EpsIndex]  = 0.5*( (c12/c6)**1/6 )
                                    parameter[QIndex]    = epsfactorSqrt*q1
      
         @return ReferenceForce::DefaultReturn
      
         --------------------------------------------------------------------------------------- */
      
      int getDerivedParameters( RealOpenMM c6, RealOpenMM c12, RealOpenMM q1, 
                                RealOpenMM epsfacSqrt,
                                RealOpenMM* parameters ) const;
      
      /**---------------------------------------------------------------------------------------
      
         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
      
         @return ReferenceForce::DefaultReturn
            
         --------------------------------------------------------------------------------------- */
          
      int calculatePairIxn( int numberOfAtoms, RealOpenMM** atomCoordinates,
                            RealOpenMM** atomParameters, int** exclusions,
                            RealOpenMM* fixedParameters, RealOpenMM** forces,
                            RealOpenMM* energyByAtom, RealOpenMM* totalEnergy ) const;
179
180

private:
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
      /**---------------------------------------------------------------------------------------
      
         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
196

197
198
199
200
201
202
203
204
         @return ReferenceForce::DefaultReturn
            
         --------------------------------------------------------------------------------------- */
          
      int calculateEwaldIxn( int numberOfAtoms, RealOpenMM** atomCoordinates,
                            RealOpenMM** atomParameters, int** exclusions,
                            RealOpenMM* fixedParameters, RealOpenMM** forces,
                            RealOpenMM* energyByAtom, RealOpenMM* totalEnergy ) const;
205
206
207
208
209
210
      
};

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

#endif // __ReferenceLJCoulombIxn_H__