forked from TinkerTools/tinker
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathetors.f
442 lines (442 loc) · 12.9 KB
/
etors.f
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
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
c
c
c ###################################################
c ## COPYRIGHT (C) 1990 by Jay William Ponder ##
c ## All Rights Reserved ##
c ###################################################
c
c ##############################################################
c ## ##
c ## subroutine etors -- torsional angle potential energy ##
c ## ##
c ##############################################################
c
c
c "etors" calculates the torsional potential energy
c
c
subroutine etors
use warp
implicit none
c
c
c choose standard or potential energy smoothing version
c
if (use_smooth) then
call etors0b
else
call etors0a
end if
return
end
c
c
c ###############################################################
c ## ##
c ## subroutine etors0a -- standard torsional angle energy ##
c ## ##
c ###############################################################
c
c
c "etors0a" calculates the torsional potential energy
c using a standard sum of Fourier terms
c
c
subroutine etors0a
use atoms
use bound
use energi
use group
use torpot
use tors
use usage
implicit none
integer i,ia,ib,ic,id
real*8 e,eps,rcb,fgrp
real*8 xt,yt,zt,rt2
real*8 xu,yu,zu,ru2
real*8 xtu,ytu,ztu,rtru
real*8 v1,v2,v3,v4,v5,v6
real*8 c1,c2,c3,c4,c5,c6
real*8 s1,s2,s3,s4,s5,s6
real*8 cosine,cosine2
real*8 cosine3,cosine4
real*8 cosine5,cosine6
real*8 sine,sine2,sine3
real*8 sine4,sine5,sine6
real*8 phi1,phi2,phi3
real*8 phi4,phi5,phi6
real*8 xia,yia,zia
real*8 xib,yib,zib
real*8 xic,yic,zic
real*8 xid,yid,zid
real*8 xba,yba,zba
real*8 xdc,ydc,zdc
real*8 xcb,ycb,zcb
logical proceed
c
c
c zero out the torsional potential energy
c
et = 0.0d0
if (ntors .eq. 0) return
c
c set tolerance for minimum distance and angle values
c
eps = 0.0001d0
c
c OpenMP directives for the major loop structure
c
!$OMP PARALLEL default(private) shared(ntors,itors,tors1,tors2,tors3,
!$OMP& tors4,tors5,tors6,use,x,y,z,torsunit,eps,use_group,use_polymer)
!$OMP& shared(et)
!$OMP DO reduction(+:et) schedule(guided)
c
c calculate the torsional angle energy term
c
do i = 1, ntors
ia = itors(1,i)
ib = itors(2,i)
ic = itors(3,i)
id = itors(4,i)
c
c decide whether to compute the current interaction
c
proceed = .true.
if (use_group) call groups (proceed,fgrp,ia,ib,ic,id,0,0)
if (proceed) proceed = (use(ia) .or. use(ib) .or.
& use(ic) .or. use(id))
c
c compute the value of the torsional angle
c
if (proceed) then
xia = x(ia)
yia = y(ia)
zia = z(ia)
xib = x(ib)
yib = y(ib)
zib = z(ib)
xic = x(ic)
yic = y(ic)
zic = z(ic)
xid = x(id)
yid = y(id)
zid = z(id)
xba = xib - xia
yba = yib - yia
zba = zib - zia
xcb = xic - xib
ycb = yic - yib
zcb = zic - zib
xdc = xid - xic
ydc = yid - yic
zdc = zid - zic
if (use_polymer) then
call image (xba,yba,zba)
call image (xcb,ycb,zcb)
call image (xdc,ydc,zdc)
end if
rcb = sqrt(max(xcb*xcb+ycb*ycb+zcb*zcb,eps))
xt = yba*zcb - ycb*zba
yt = zba*xcb - zcb*xba
zt = xba*ycb - xcb*yba
xu = ycb*zdc - ydc*zcb
yu = zcb*xdc - zdc*xcb
zu = xcb*ydc - xdc*ycb
xtu = yt*zu - yu*zt
ytu = zt*xu - zu*xt
ztu = xt*yu - xu*yt
rt2 = max(xt*xt+yt*yt+zt*zt,eps)
ru2 = max(xu*xu+yu*yu+zu*zu,eps)
rtru = sqrt(rt2 * ru2)
cosine = (xt*xu + yt*yu + zt*zu) / rtru
sine = (xcb*xtu + ycb*ytu + zcb*ztu) / (rcb*rtru)
c
c set the torsional parameters for this angle
c
v1 = tors1(1,i)
c1 = tors1(3,i)
s1 = tors1(4,i)
v2 = tors2(1,i)
c2 = tors2(3,i)
s2 = tors2(4,i)
v3 = tors3(1,i)
c3 = tors3(3,i)
s3 = tors3(4,i)
v4 = tors4(1,i)
c4 = tors4(3,i)
s4 = tors4(4,i)
v5 = tors5(1,i)
c5 = tors5(3,i)
s5 = tors5(4,i)
v6 = tors6(1,i)
c6 = tors6(3,i)
s6 = tors6(4,i)
c
c compute the multiple angle trigonometry and the phase terms
c
cosine2 = cosine*cosine - sine*sine
sine2 = 2.0d0 * cosine * sine
cosine3 = cosine*cosine2 - sine*sine2
sine3 = cosine*sine2 + sine*cosine2
cosine4 = cosine*cosine3 - sine*sine3
sine4 = cosine*sine3 + sine*cosine3
cosine5 = cosine*cosine4 - sine*sine4
sine5 = cosine*sine4 + sine*cosine4
cosine6 = cosine*cosine5 - sine*sine5
sine6 = cosine*sine5 + sine*cosine5
phi1 = 1.0d0 + (cosine*c1 + sine*s1)
phi2 = 1.0d0 + (cosine2*c2 + sine2*s2)
phi3 = 1.0d0 + (cosine3*c3 + sine3*s3)
phi4 = 1.0d0 + (cosine4*c4 + sine4*s4)
phi5 = 1.0d0 + (cosine5*c5 + sine5*s5)
phi6 = 1.0d0 + (cosine6*c6 + sine6*s6)
c
c calculate the torsional energy for this angle
c
e = torsunit * (v1*phi1 + v2*phi2 + v3*phi3
& + v4*phi4 + v5*phi5 + v6*phi6)
c
c scale the interaction based on its group membership
c
if (use_group) e = e * fgrp
c
c increment the total torsional angle energy
c
et = et + e
end if
end do
c
c OpenMP directives for the major loop structure
c
!$OMP END DO
!$OMP END PARALLEL
return
end
c
c ##############################################################
c ## ##
c ## subroutine etors0b -- torsional energy for smoothing ##
c ## ##
c ##############################################################
c
c
c "etors0b" calculates the torsional potential energy
c for use with potential energy smoothing methods
c
c
subroutine etors0b
use atoms
use energi
use group
use math
use torpot
use tors
use usage
use warp
implicit none
integer i,ia,ib,ic,id
real*8 e,eps,rcb,fgrp
real*8 width,wterm
real*8 xt,yt,zt,rt2
real*8 xu,yu,zu,ru2
real*8 xtu,ytu,ztu,rtru
real*8 v1,v2,v3,v4,v5,v6
real*8 c1,c2,c3,c4,c5,c6
real*8 s1,s2,s3,s4,s5,s6
real*8 cosine,cosine2
real*8 cosine3,cosine4
real*8 cosine5,cosine6
real*8 sine,sine2,sine3
real*8 sine4,sine5,sine6
real*8 damp1,damp2,damp3
real*8 damp4,damp5,damp6
real*8 phi1,phi2,phi3
real*8 phi4,phi5,phi6
real*8 xia,yia,zia
real*8 xib,yib,zib
real*8 xic,yic,zic
real*8 xid,yid,zid
real*8 xba,yba,zba
real*8 xdc,ydc,zdc
real*8 xcb,ycb,zcb
logical proceed
c
c
c zero out the torsional potential energy
c
et = 0.0d0
if (ntors .eq. 0) return
c
c set tolerance for minimum distance and angle values
c
eps = 0.0001d0
c
c set the extent of smoothing to be performed
c
width = difft * deform
if (width .le. 0.0d0) then
damp1 = 1.0d0
damp2 = 1.0d0
damp3 = 1.0d0
damp4 = 1.0d0
damp5 = 1.0d0
damp6 = 1.0d0
else if (use_dem) then
damp1 = exp(-width)
damp2 = exp(-4.0d0*width)
damp3 = exp(-9.0d0*width)
damp4 = exp(-16.0d0*width)
damp5 = exp(-25.0d0*width)
damp6 = exp(-36.0d0*width)
else if (use_gda) then
wterm = difft / 12.0d0
else if (use_tophat .or. use_stophat) then
damp1 = 0.0d0
damp2 = 0.0d0
damp3 = 0.0d0
damp4 = 0.0d0
damp5 = 0.0d0
damp6 = 0.0d0
if (width .lt. pi) damp1 = sin(width) / width
wterm = 2.0d0 * width
if (wterm .lt. pi) damp2 = sin(wterm) / wterm
wterm = 3.0d0 * width
if (wterm .lt. pi) damp3 = sin(wterm) / wterm
wterm = 4.0d0 * width
if (wterm .lt. pi) damp4 = sin(wterm) / wterm
wterm = 5.0d0 * width
if (wterm .lt. pi) damp5 = sin(wterm) / wterm
wterm = 6.0d0 * width
if (wterm .lt. pi) damp6 = sin(wterm) / wterm
end if
c
c calculate the torsional angle energy term
c
do i = 1, ntors
ia = itors(1,i)
ib = itors(2,i)
ic = itors(3,i)
id = itors(4,i)
c
c decide whether to compute the current interaction
c
proceed = .true.
if (use_group) call groups (proceed,fgrp,ia,ib,ic,id,0,0)
if (proceed) proceed = (use(ia) .or. use(ib) .or.
& use(ic) .or. use(id))
c
c compute the value of the torsional angle
c
if (proceed) then
xia = x(ia)
yia = y(ia)
zia = z(ia)
xib = x(ib)
yib = y(ib)
zib = z(ib)
xic = x(ic)
yic = y(ic)
zic = z(ic)
xid = x(id)
yid = y(id)
zid = z(id)
xba = xib - xia
yba = yib - yia
zba = zib - zia
xcb = xic - xib
ycb = yic - yib
zcb = zic - zib
xdc = xid - xic
ydc = yid - yic
zdc = zid - zic
rcb = sqrt(max(xcb*xcb+ycb*ycb+zcb*zcb,eps))
xt = yba*zcb - ycb*zba
yt = zba*xcb - zcb*xba
zt = xba*ycb - xcb*yba
xu = ycb*zdc - ydc*zcb
yu = zcb*xdc - zdc*xcb
zu = xcb*ydc - xdc*ycb
xtu = yt*zu - yu*zt
ytu = zt*xu - zu*xt
ztu = xt*yu - xu*yt
rt2 = max(xt*xt+yt*yt+zt*zt,eps)
ru2 = max(xu*xu+yu*yu+zu*zu,eps)
rtru = sqrt(rt2 * ru2)
cosine = (xt*xu + yt*yu + zt*zu) / rtru
sine = (xcb*xtu + ycb*ytu + zcb*ztu) / (rcb*rtru)
c
c set the torsional parameters for this angle
c
v1 = tors1(1,i)
c1 = tors1(3,i)
s1 = tors1(4,i)
v2 = tors2(1,i)
c2 = tors2(3,i)
s2 = tors2(4,i)
v3 = tors3(1,i)
c3 = tors3(3,i)
s3 = tors3(4,i)
v4 = tors4(1,i)
c4 = tors4(3,i)
s4 = tors4(4,i)
v5 = tors5(1,i)
c5 = tors5(3,i)
s5 = tors5(4,i)
v6 = tors6(1,i)
c6 = tors6(3,i)
s6 = tors6(4,i)
c
c compute the multiple angle trigonometry and the phase terms
c
cosine2 = cosine*cosine - sine*sine
sine2 = 2.0d0 * cosine * sine
cosine3 = cosine*cosine2 - sine*sine2
sine3 = cosine*sine2 + sine*cosine2
cosine4 = cosine*cosine3 - sine*sine3
sine4 = cosine*sine3 + sine*cosine3
cosine5 = cosine*cosine4 - sine*sine4
sine5 = cosine*sine4 + sine*cosine4
cosine6 = cosine*cosine5 - sine*sine5
sine6 = cosine*sine5 + sine*cosine5
phi1 = 1.0d0 + (cosine*c1 + sine*s1)
phi2 = 1.0d0 + (cosine2*c2 + sine2*s2)
phi3 = 1.0d0 + (cosine3*c3 + sine3*s3)
phi4 = 1.0d0 + (cosine4*c4 + sine4*s4)
phi5 = 1.0d0 + (cosine5*c5 + sine5*s5)
phi6 = 1.0d0 + (cosine6*c6 + sine6*s6)
c
c transform the potential function via smoothing
c
if (use_gda) then
width = wterm * (m2(ia)+m2(ib)+m2(ic)+m2(id))
damp1 = exp(-width)
damp2 = exp(-4.0d0*width)
damp3 = exp(-9.0d0*width)
damp4 = exp(-16.0d0*width)
damp5 = exp(-25.0d0*width)
damp6 = exp(-36.0d0*width)
end if
phi1 = phi1 * damp1
phi2 = phi2 * damp2
phi3 = phi3 * damp3
phi4 = phi4 * damp4
phi5 = phi5 * damp5
phi6 = phi6 * damp6
c
c calculate the torsional energy for this angle
c
e = torsunit * (v1*phi1 + v2*phi2 + v3*phi3
& + v4*phi4 + v5*phi5 + v6*phi6)
c
c scale the interaction based on its group membership
c
if (use_group) e = e * fgrp
c
c increment the total torsional angle energy
c
et = et + e
end if
end do
return
end