-
Notifications
You must be signed in to change notification settings - Fork 1
/
ComputeDihedrals.C
290 lines (243 loc) · 9.74 KB
/
ComputeDihedrals.C
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
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
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
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
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
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
/**
*** Copyright (c) 1995, 1996, 1997, 1998, 1999, 2000 by
*** The Board of Trustees of the University of Illinois.
*** All rights reserved.
**/
#include "InfoStream.h"
#include "ComputeDihedrals.h"
#include "Molecule.h"
#include "Parameters.h"
#include "Node.h"
#include "ReductionMgr.h"
#include "Lattice.h"
#include "PressureProfile.h"
#include "Debug.h"
// static initialization
int DihedralElem::pressureProfileSlabs = 0;
int DihedralElem::pressureProfileAtomTypes = 1;
BigReal DihedralElem::pressureProfileThickness = 0;
BigReal DihedralElem::pressureProfileMin = 0;
void DihedralElem::getMoleculePointers
(Molecule* mol, int* count, int32*** byatom, Dihedral** structarray)
{
#ifdef MEM_OPT_VERSION
NAMD_die("Should not be called in DihedralElem::getMoleculePointers in memory optimized version!");
#else
*count = mol->numDihedrals;
*byatom = mol->dihedralsByAtom;
*structarray = mol->dihedrals;
#endif
}
void DihedralElem::getParameterPointers(Parameters *p, const DihedralValue **v) {
*v = p->dihedral_array;
}
void DihedralElem::computeForce(BigReal *reduction,
BigReal *pressureProfileData)
{
DebugM(3, "::computeForce() localIndex = " << localIndex[0] << " "
<< localIndex[1] << " " << localIndex[2] << std::endl);
// Calculate the vectors between atoms
const Position & pos0 = p[0]->x[localIndex[0]].position;
const Lattice & lattice = p[0]->p->lattice;
const Position & pos1 = p[1]->x[localIndex[1]].position;
const Vector r12 = lattice.delta(pos0,pos1);
const Position & pos2 = p[2]->x[localIndex[2]].position;
const Vector r23 = lattice.delta(pos1,pos2);
const Position & pos3 = p[3]->x[localIndex[3]].position;
const Vector r34 = lattice.delta(pos2,pos3);
// Calculate the cross products and distances
Vector A = cross(r12,r23);
register BigReal rAinv = A.rlength();
Vector B = cross(r23,r34);
register BigReal rBinv = B.rlength();
Vector C = cross(r23,A);
register BigReal rCinv = C.rlength();
// Calculate the sin and cos
BigReal cos_phi = (A*B)*(rAinv*rBinv);
BigReal sin_phi = (C*B)*(rCinv*rBinv);
BigReal phi= -atan2(sin_phi,cos_phi);
BigReal K=0; // energy
BigReal K1=0; // force
// get the dihedral information
int multiplicity = value->multiplicity;
// Loop through the multiple parameter sets for this
// bond. We will only loop more than once if this
// has multiple parameter sets from Charmm22
for (int mult_num=0; mult_num<multiplicity; mult_num++)
{
/* get angle information */
Real k = value->values[mult_num].k * scale;
Real delta = value->values[mult_num].delta;
int n = value->values[mult_num].n;
// Calculate the energy
if (n)
{
// Periodicity is greater than 0, so use cos form
K += k*(1+cos(n*phi - delta));
K1 += -n*k*sin(n*phi - delta);
}
else
{
// Periodicity is 0, so just use the harmonic form
BigReal diff = phi-delta;
if (diff < -PI) diff += TWOPI;
else if (diff > PI) diff -= TWOPI;
K += k*diff*diff;
K1 += 2.0*k*diff;
}
} /* for multiplicity */
Force f1,f2,f3;
// Normalize B
//rB = 1.0/rB;
B *= rBinv;
// Next, we want to calculate the forces. In order
// to do that, we first need to figure out whether the
// sin or cos form will be more stable. For this,
// just look at the value of phi
if (fabs(sin_phi) > 0.1)
{
// use the sin version to avoid 1/cos terms
// Normalize A
A *= rAinv;
Vector dcosdA;
Vector dcosdB;
dcosdA.x = rAinv*(cos_phi*A.x-B.x);
dcosdA.y = rAinv*(cos_phi*A.y-B.y);
dcosdA.z = rAinv*(cos_phi*A.z-B.z);
dcosdB.x = rBinv*(cos_phi*B.x-A.x);
dcosdB.y = rBinv*(cos_phi*B.y-A.y);
dcosdB.z = rBinv*(cos_phi*B.z-A.z);
K1 = K1/sin_phi;
f1.x = K1*(r23.y*dcosdA.z - r23.z*dcosdA.y);
f1.y = K1*(r23.z*dcosdA.x - r23.x*dcosdA.z);
f1.z = K1*(r23.x*dcosdA.y - r23.y*dcosdA.x);
f3.x = K1*(r23.z*dcosdB.y - r23.y*dcosdB.z);
f3.y = K1*(r23.x*dcosdB.z - r23.z*dcosdB.x);
f3.z = K1*(r23.y*dcosdB.x - r23.x*dcosdB.y);
f2.x = K1*(r12.z*dcosdA.y - r12.y*dcosdA.z
+ r34.y*dcosdB.z - r34.z*dcosdB.y);
f2.y = K1*(r12.x*dcosdA.z - r12.z*dcosdA.x
+ r34.z*dcosdB.x - r34.x*dcosdB.z);
f2.z = K1*(r12.y*dcosdA.x - r12.x*dcosdA.y
+ r34.x*dcosdB.y - r34.y*dcosdB.x);
}
else
{
// This angle is closer to 0 or 180 than it is to
// 90, so use the cos version to avoid 1/sin terms
// Normalize C
// rC = 1.0/rC;
C *= rCinv;
Vector dsindC;
Vector dsindB;
dsindC.x = rCinv*(sin_phi*C.x-B.x);
dsindC.y = rCinv*(sin_phi*C.y-B.y);
dsindC.z = rCinv*(sin_phi*C.z-B.z);
dsindB.x = rBinv*(sin_phi*B.x-C.x);
dsindB.y = rBinv*(sin_phi*B.y-C.y);
dsindB.z = rBinv*(sin_phi*B.z-C.z);
K1 = -K1/cos_phi;
f1.x = K1*((r23.y*r23.y + r23.z*r23.z)*dsindC.x
- r23.x*r23.y*dsindC.y
- r23.x*r23.z*dsindC.z);
f1.y = K1*((r23.z*r23.z + r23.x*r23.x)*dsindC.y
- r23.y*r23.z*dsindC.z
- r23.y*r23.x*dsindC.x);
f1.z = K1*((r23.x*r23.x + r23.y*r23.y)*dsindC.z
- r23.z*r23.x*dsindC.x
- r23.z*r23.y*dsindC.y);
f3 = cross(K1,dsindB,r23);
f2.x = K1*(-(r23.y*r12.y + r23.z*r12.z)*dsindC.x
+(2.0*r23.x*r12.y - r12.x*r23.y)*dsindC.y
+(2.0*r23.x*r12.z - r12.x*r23.z)*dsindC.z
+dsindB.z*r34.y - dsindB.y*r34.z);
f2.y = K1*(-(r23.z*r12.z + r23.x*r12.x)*dsindC.y
+(2.0*r23.y*r12.z - r12.y*r23.z)*dsindC.z
+(2.0*r23.y*r12.x - r12.y*r23.x)*dsindC.x
+dsindB.x*r34.z - dsindB.z*r34.x);
f2.z = K1*(-(r23.x*r12.x + r23.y*r12.y)*dsindC.z
+(2.0*r23.z*r12.x - r12.z*r23.x)*dsindC.x
+(2.0*r23.z*r12.y - r12.z*r23.y)*dsindC.y
+dsindB.y*r34.x - dsindB.x*r34.y);
}
/* store the forces */
// p[0]->f[localIndex[0]] += f1;
// p[1]->f[localIndex[1]] += f2 - f1;
// p[2]->f[localIndex[2]] += f3 - f2;
// p[3]->f[localIndex[3]] += -f3;
p[0]->f[localIndex[0]].x += f1.x;
p[0]->f[localIndex[0]].y += f1.y;
p[0]->f[localIndex[0]].z += f1.z;
p[1]->f[localIndex[1]].x += f2.x - f1.x;
p[1]->f[localIndex[1]].y += f2.y - f1.y;
p[1]->f[localIndex[1]].z += f2.z - f1.z;
p[2]->f[localIndex[2]].x += f3.x - f2.x;
p[2]->f[localIndex[2]].y += f3.y - f2.y;
p[2]->f[localIndex[2]].z += f3.z - f2.z;
p[3]->f[localIndex[3]].x += -f3.x;
p[3]->f[localIndex[3]].y += -f3.y;
p[3]->f[localIndex[3]].z += -f3.z;
/* store the force for dihedral-only accelMD */
if ( p[0]->af ) {
p[0]->af[localIndex[0]].x += f1.x;
p[0]->af[localIndex[0]].y += f1.y;
p[0]->af[localIndex[0]].z += f1.z;
p[1]->af[localIndex[1]].x += f2.x - f1.x;
p[1]->af[localIndex[1]].y += f2.y - f1.y;
p[1]->af[localIndex[1]].z += f2.z - f1.z;
p[2]->af[localIndex[2]].x += f3.x - f2.x;
p[2]->af[localIndex[2]].y += f3.y - f2.y;
p[2]->af[localIndex[2]].z += f3.z - f2.z;
p[3]->af[localIndex[3]].x += -f3.x;
p[3]->af[localIndex[3]].y += -f3.y;
p[3]->af[localIndex[3]].z += -f3.z;
}
DebugM(3, "::computeForce() -- ending with delta energy " << K << std::endl);
reduction[dihedralEnergyIndex] += K;
reduction[virialIndex_XX] += ( f1.x * r12.x + f2.x * r23.x + f3.x * r34.x );
reduction[virialIndex_XY] += ( f1.x * r12.y + f2.x * r23.y + f3.x * r34.y );
reduction[virialIndex_XZ] += ( f1.x * r12.z + f2.x * r23.z + f3.x * r34.z );
reduction[virialIndex_YX] += ( f1.y * r12.x + f2.y * r23.x + f3.y * r34.x );
reduction[virialIndex_YY] += ( f1.y * r12.y + f2.y * r23.y + f3.y * r34.y );
reduction[virialIndex_YZ] += ( f1.y * r12.z + f2.y * r23.z + f3.y * r34.z );
reduction[virialIndex_ZX] += ( f1.z * r12.x + f2.z * r23.x + f3.z * r34.x );
reduction[virialIndex_ZY] += ( f1.z * r12.y + f2.z * r23.y + f3.z * r34.y );
reduction[virialIndex_ZZ] += ( f1.z * r12.z + f2.z * r23.z + f3.z * r34.z );
if (pressureProfileData) {
BigReal z1 = p[0]->x[localIndex[0]].position.z;
BigReal z2 = p[1]->x[localIndex[1]].position.z;
BigReal z3 = p[2]->x[localIndex[2]].position.z;
BigReal z4 = p[3]->x[localIndex[3]].position.z;
int n1 = (int)floor((z1-pressureProfileMin)/pressureProfileThickness);
int n2 = (int)floor((z2-pressureProfileMin)/pressureProfileThickness);
int n3 = (int)floor((z3-pressureProfileMin)/pressureProfileThickness);
int n4 = (int)floor((z4-pressureProfileMin)/pressureProfileThickness);
pp_clamp(n1, pressureProfileSlabs);
pp_clamp(n2, pressureProfileSlabs);
pp_clamp(n3, pressureProfileSlabs);
pp_clamp(n4, pressureProfileSlabs);
int p1 = p[0]->x[localIndex[0]].partition;
int p2 = p[1]->x[localIndex[1]].partition;
int p3 = p[2]->x[localIndex[2]].partition;
int p4 = p[3]->x[localIndex[3]].partition;
int pn = pressureProfileAtomTypes;
pp_reduction(pressureProfileSlabs, n1, n2,
p1, p2, pn,
f1.x * r12.x, f1.y * r12.y, f1.z * r12.z,
pressureProfileData);
pp_reduction(pressureProfileSlabs, n2, n3,
p2, p3, pn,
f2.x * r23.x, f2.y * r23.y, f2.z * r23.z,
pressureProfileData);
pp_reduction(pressureProfileSlabs, n3, n4,
p3, p4, pn,
f3.x * r34.x, f3.y * r34.y, f3.z * r34.z,
pressureProfileData);
}
}
void DihedralElem::submitReductionData(BigReal *data, SubmitReduction *reduction)
{
reduction->item(REDUCTION_DIHEDRAL_ENERGY) += data[dihedralEnergyIndex];
ADD_TENSOR(reduction,REDUCTION_VIRIAL_NORMAL,data,virialIndex);
ADD_TENSOR(reduction,REDUCTION_VIRIAL_AMD_DIHE,data,virialIndex);
}