forked from TinkerTools/tinker
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathbaoab.f
236 lines (236 loc) · 6.29 KB
/
baoab.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
c
c
c ###############################################################
c ## COPYRIGHT (C) 2016 by Charles Matthews & Ben Leimkuhler ##
c ## All Rights Reserved ##
c ###############################################################
c
c ############################################################
c ## ##
c ## subroutine baoab -- BAOAB stochastic dynamics step ##
c ## ##
c ############################################################
c
c
c "baoab" implements a constrained stochastic dynamics time
c step using the geodesic BAOAB scheme
c
c literature reference:
c
c B. Leimkuhler and C. Matthews, "Efficient Molecular Dynamics
c Using Geodesic Integration and Solvent-Solute Splitting",
c Proceedings of the Royal Society A, 472, 20160138 (2016)
c
c
subroutine baoab (istep,dt)
use atomid
use atoms
use freeze
use mdstuf
use moldyn
use units
use usage
use virial
use limits
use potent
implicit none
integer i,j,k
integer istep
integer nrattle
real*8 dt,dtr
real*8 dt_2,dt_4
real*8 etot,epot
real*8 eksum
real*8 temp,pres
real*8 ekin(3,3)
real*8 stress(3,3)
real*8, allocatable :: xold(:)
real*8, allocatable :: yold(:)
real*8, allocatable :: zold(:)
real*8, allocatable :: vfric(:)
real*8, allocatable :: vrand(:,:)
real*8, allocatable :: derivs(:,:)
c
c
c set some time values for the dynamics integration
c
dt_2 = 0.5d0 * dt
dt_4 = 0.25d0 * dt
nrattle = 1
dtr = dt_2 / dble(nrattle)
c
c make half-step temperature and pressure corrections
c
call temper (dt,eksum,ekin,temp)
c
c perform dynamic allocation of some local arrays
c
allocate (xold(n))
allocate (yold(n))
allocate (zold(n))
allocate (derivs(3,n))
allocate (vfric(n))
allocate (vrand(3,n))
c
c find half-step velocities via the Verlet recursion
c
do i = 1, nuse
k = iuse(i)
do j = 1, 3
v(j,k) = v(j,k) + a(j,k)*dt_2
end do
end do
c
c find the constraint-corrected full-step velocities
c
if (use_rattle) call rattle2 (dt)
c
c take first A step according to the BAOAB sequence
c
do j = 1, nrattle
do i = 1, nuse
k = iuse(i)
xold(k) = x(k)
yold(k) = y(k)
zold(k) = z(k)
x(k) = x(k) + v(1,k)*dtr
y(k) = y(k) + v(2,k)*dtr
z(k) = z(k) + v(3,k)*dtr
end do
if (use_rattle) call rattle (dtr,xold,yold,zold)
if (use_rattle) call rattle2 (dtr)
end do
c
c find the constraint-corrected full-step velocities
c
if (use_rattle) call rattle2 (dt)
c
c update velocities with frictional and random components
c
call oprep (dt,vfric,vrand)
do i = 1, nuse
k = iuse(i)
do j = 1, 3
v(j,k) = v(j,k)*vfric(k) + vrand(j,k)
end do
end do
if (use_rattle) call rattle2 (dt)
c
c take second A step according to the BAOAB sequence
c
do j = 1, nrattle
do i = 1, nuse
k = iuse(i)
xold(k) = x(k)
yold(k) = y(k)
zold(k) = z(k)
x(k) = x(k) + v(1,k)*dtr
y(k) = y(k) + v(2,k)*dtr
z(k) = z(k) + v(3,k)*dtr
end do
if (use_rattle) call rattle (dtr,xold,yold,zold)
if (use_rattle) call rattle2 (dtr)
end do
c
c get the potential energy and atomic forces
c
call gradient (epot,derivs)
c
c use Newton's second law to get the next accelerations;
c find the full-step velocities using the Verlet recursion
c
do i = 1, nuse
k = iuse(i)
do j = 1, 3
a(j,k) = -ekcal * derivs(j,k) / mass(k)
v(j,k) = v(j,k) + a(j,k)*dt_2
end do
end do
c
c perform deallocation of some local arrays
c
deallocate (xold)
deallocate (yold)
deallocate (zold)
deallocate (derivs)
deallocate (vfric)
deallocate (vrand)
c
c find the constraint-corrected full-step velocities
c
if (use_rattle) call rattle2 (dt)
c
c compute and control the temperature and pressure
c
call kinetic (eksum,ekin,temp)
temp = 2.0d0 * eksum / (dble(nfree) * gasconst)
call pressure (dt,epot,ekin,temp,pres,stress)
c
c total energy is sum of kinetic and potential energies
c
etot = eksum + epot
c
c compute statistics and save trajectory for this step
c
call mdstat (istep,dt,etot,epot,eksum,temp,pres)
call mdsave (istep,dt,epot,eksum)
call mdrest (istep)
return
end
c
c
c #################################################################
c ## ##
c ## subroutine oprep -- frictional & random terms for BAOAB ##
c ## ##
c #################################################################
c
c
c "oprep" sets up the frictional and random terms needed to
c update positions and velocities for the BAOAB integrator
c
c
subroutine oprep (dt,vfric,vrand)
use atoms
use atomid
use bath
use stodyn
use units
use usage
implicit none
integer i,j,k
real*8 dt,ktm
real*8 egdt,normal
real*8 vsig,vnorm
real*8 vfric(*)
real*8 vrand(3,*)
logical first
save first
data first / .true. /
c
c
c set the atomic friction coefficients to the global value
c
if (first) then
first = .false.
if (.not. allocated(fgamma)) allocate (fgamma(n))
do i = 1, n
fgamma(i) = friction
end do
end if
c
c get the frictional and random terms for a BAOAB step
c
egdt = exp(-(friction * dt))
do i = 1, nuse
k = iuse(i)
vfric(k) = egdt
ktm = boltzmann * kelvin / mass(k)
vsig = sqrt(ktm*(1.0-egdt*egdt))
do j = 1, 3
vnorm = normal ()
vrand(j,k) = vsig * vnorm
end do
end do
return
end