forked from m3g/packmol
-
Notifications
You must be signed in to change notification settings - Fork 0
/
computeg.f90
309 lines (227 loc) · 7.83 KB
/
computeg.f90
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
!
! Written by Leandro Martínez, 2009-2011.
! Copyright (c) 2009-2018, Leandro Martínez, Jose Mario Martinez,
! Ernesto G. Birgin.
!
! Subroutine that computes the analytical derivatives
!
subroutine computeg(n,x,g)
use sizes
use compute_data
use input, only : fix
implicit none
integer :: n
integer :: idatom, iatom, irest
integer :: i, j, k, ilubar, ilugan, icart, itype, imol
integer :: ibox, iboxx, iboxy, iboxz
integer :: k1, k2
integer :: iratcount
double precision :: x(n), g(n)
double precision :: dv1beta(3), dv1gama(3), dv1teta(3),&
dv2beta(3), dv2gama(3), dv2teta(3),&
dv3beta(3), dv3gama(3), dv3teta(3)
double precision :: v1(3), v2(3), v3(3)
double precision :: xbar, ybar, zbar
double precision :: xtemp, ytemp, ztemp
double precision :: beta, gama, teta, cb, sb, cg, sg, ct, st
! Reset gradients
do i = 1, ntotat
do j = 1, 3
gxcar(i,j) = 0.d0
end do
end do
! Reset boxes
if(.not.init1) call resetboxes()
! Transform baricenter and angles into cartesian coordinates
! Computes cartesian coordinates from vector x and coor
ilubar = 0
ilugan = ntotmol*3
icart = 0
do itype = 1, ntype
if(.not.comptype(itype)) then
icart = icart + nmols(itype)*natoms(itype)
else
do imol = 1, nmols(itype)
xbar = x(ilubar + 1)
ybar = x(ilubar + 2)
zbar = x(ilubar + 3)
! Compute the rotation matrix
beta = x(ilugan + 1)
gama = x(ilugan + 2)
teta = x(ilugan + 3)
call eulerrmat(beta,gama,teta,v1,v2,v3)
idatom = idfirst(itype) - 1
do iatom = 1, natoms(itype)
icart = icart + 1
idatom = idatom + 1
call compcart(icart,xbar,ybar,zbar, &
coor(idatom,1),coor(idatom,2),coor(idatom,3), &
v1,v2,v3)
! Gradient relative to the wall distace
do iratcount = 1, nratom(icart)
irest = iratom(icart,iratcount)
call gwalls(icart,irest)
end do
if(.not.init1) then
xtemp = xcart(icart,1) - sizemin(1)
ytemp = xcart(icart,2) - sizemin(2)
ztemp = xcart(icart,3) - sizemin(3)
iboxx = int(xtemp/boxl(1)) + 1
iboxy = int(ytemp/boxl(2)) + 1
iboxz = int(ztemp/boxl(3)) + 1
if(xtemp.le.0) iboxx = 1
if(ytemp.le.0) iboxy = 1
if(ztemp.le.0) iboxz = 1
if(iboxx.gt.nboxes(1)) iboxx = nboxes(1)
if(iboxy.gt.nboxes(2)) iboxy = nboxes(2)
if(iboxz.gt.nboxes(3)) iboxz = nboxes(3)
! Atom linked list
latomnext(icart) = latomfirst(iboxx,iboxy,iboxz)
latomfirst(iboxx,iboxy,iboxz) = icart
! Box with atoms linked list
if ( .not. hasfree(iboxx,iboxy,iboxz) ) then
hasfree(iboxx,iboxy,iboxz) = .true.
call ijk_to_ibox(iboxx,iboxy,iboxz,ibox)
lboxnext(ibox) = lboxfirst
lboxfirst = ibox
! Add boxes with fixed atoms which are vicinal to this box, and
! are behind
if ( fix ) then
call add_box_behind(iboxx-1,iboxy,iboxz)
call add_box_behind(iboxx,iboxy-1,iboxz)
call add_box_behind(iboxx,iboxy,iboxz-1)
call add_box_behind(iboxx,iboxy-1,iboxz+1)
call add_box_behind(iboxx,iboxy-1,iboxz-1)
call add_box_behind(iboxx-1,iboxy+1,iboxz)
call add_box_behind(iboxx-1,iboxy,iboxz+1)
call add_box_behind(iboxx-1,iboxy-1,iboxz)
call add_box_behind(iboxx-1,iboxy,iboxz-1)
call add_box_behind(iboxx-1,iboxy+1,iboxz+1)
call add_box_behind(iboxx-1,iboxy+1,iboxz-1)
call add_box_behind(iboxx-1,iboxy-1,iboxz+1)
call add_box_behind(iboxx-1,iboxy-1,iboxz-1)
end if
end if
ibtype(icart) = itype
ibmol(icart) = imol
end if
end do
ilugan = ilugan + 3
ilubar = ilubar + 3
end do
end if
end do
if( .not. init1 ) then
!
! Gradient relative to minimum distance
!
ibox = lboxfirst
do while( ibox > 0 )
call ibox_to_ijk(ibox,i,j,k)
icart = latomfirst(i,j,k)
do while ( icart .ne. 0 )
if(comptype(ibtype(icart))) then
! Interactions inside box
call gparc(icart,latomnext(icart))
! Interactions of boxes that share faces
call gparc(icart,latomfirst(i+1,j,k))
call gparc(icart,latomfirst(i,j+1,k))
call gparc(icart,latomfirst(i,j,k+1))
! Interactions of boxes that share axes
call gparc(icart,latomfirst(i+1,j+1,k))
call gparc(icart,latomfirst(i+1,j,k+1))
call gparc(icart,latomfirst(i+1,j-1,k))
call gparc(icart,latomfirst(i+1,j,k-1))
call gparc(icart,latomfirst(i,j+1,k+1))
call gparc(icart,latomfirst(i,j+1,k-1))
! Interactions of boxes that share vertices
call gparc(icart,latomfirst(i+1,j+1,k+1))
call gparc(icart,latomfirst(i+1,j+1,k-1))
call gparc(icart,latomfirst(i+1,j-1,k+1))
call gparc(icart,latomfirst(i+1,j-1,k-1))
end if
icart = latomnext(icart)
end do
ibox = lboxnext(ibox)
end do
end if
! Computing the gradient using chain rule
do i = 1, n
g(i) = 0.d0
end do
k1 = 0
k2 = ntotmol * 3
icart = 0
do itype = 1, ntype
if(.not.comptype(itype)) then
icart = icart + nmols(itype)*natoms(itype)
else
do imol = 1, nmols(itype)
beta = x(k2 + 1)
gama = x(k2 + 2)
teta = x(k2 + 3)
cb = dcos(beta)
sb = dsin(beta)
cg = dcos(gama)
sg = dsin(gama)
ct = dcos(teta)
st = dsin(teta)
dv1beta(1) = - cb * sg * ct - sb * cg
dv2beta(1) = - sb * sg * ct + cb * cg
dv3beta(1) = 0.d0
dv1gama(1) = - sb * cg * ct - cb * sg
dv2gama(1) = cb * cg * ct - sb * sg
dv3gama(1) = cg * st
dv1teta(1) = sb * sg * st
dv2teta(1) = - cb * sg * st
dv3teta(1) = sg * ct
dv1beta(2) = - cb * cg * ct + sb * sg
dv2beta(2) = - sb * cg * ct - cb * sg
dv3beta(2) = 0.d0
dv1gama(2) = sb * sg * ct - cb * cg
dv2gama(2) = - sg * cb * ct - cg * sb
dv3gama(2) = - sg * st
dv1teta(2) = sb * cg * st
dv2teta(2) = - cb * cg * st
dv3teta(2) = cg * ct
dv1beta(3) = cb * st
dv2beta(3) = sb * st
dv3beta(3) = 0.d0
dv1gama(3) = 0.d0
dv2gama(3) = 0.d0
dv3gama(3) = 0.d0
dv1teta(3) = sb * ct
dv2teta(3) = - cb * ct
dv3teta(3) = - st
idatom = idfirst(itype) - 1
do iatom = 1, natoms(itype)
icart = icart + 1
idatom = idatom + 1
do k = 1, 3
g(k1+k) = g(k1+k) + gxcar(icart, k)
end do
do k = 1, 3
g(k2 + 1) = g(k2 + 1) &
+ (coor(idatom,1) * dv1beta(k) &
+ coor(idatom, 2) * dv2beta(k) &
+ coor(idatom, 3) * dv3beta(k)) &
* gxcar(icart, k)
g(k2 + 2) = g(k2 + 2) &
+ (coor(idatom,1) * dv1gama(k) &
+ coor(idatom, 2) * dv2gama(k) &
+ coor(idatom, 3) * dv3gama(k)) &
* gxcar(icart, k)
g(k2 + 3) = g(k2 + 3) &
+ (coor(idatom,1) * dv1teta(k) &
+ coor(idatom, 2) * dv2teta(k) &
+ coor(idatom, 3) * dv3teta(k)) &
* gxcar(icart, k)
end do
end do
k2 = k2 + 3
k1 = k1 + 3
end do
end if
end do
return
end subroutine computeg