ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/f_verlet_constrained.F90
Revision: 389
Committed: Mon Mar 24 15:26:05 2003 UTC (21 years, 3 months ago) by mmeineke
File size: 10149 byte(s)
Log Message:
fixed bug where short range interactions were not being calculated.

removed some debug print statements

File Contents

# Content
1 subroutine v_constrain_a(dt,natoms,mass,rx,ry,rz, &
2 vx,vy,vz,fx,fy,fz, n_constrained, constraints_sqr, c_i, c_j, &
3 box_x, box_y, box_z)
4 implicit none
5
6 ! *******************************************************************
7 ! ** FIRST PART OF VELOCITY VERLET ALGORITHM **
8 ! ** **
9 ! ** USAGE: **
10 ! ** **
11 ! ** THE FIRST PART OF THE ALGORITHM IS A TAYLOR SERIES WHICH **
12 ! ** ADVANCES POSITIONS FROM T TO T + DT AND VELOCITIES FROM **
13 ! ** T TO T + DT/2. AFTER THIS, THE FORCE ROUTINE IS CALLED. **
14 ! *******************************************************************
15
16
17 ! move part a calculate velocities
18
19 INTEGER I
20 double precision DT2, DTSQ2
21
22 double precision box_x, box_y, box_z
23
24 ! Global Parameters
25
26 INTEGER maxn,ndim,max_dim
27 parameter (maxn=2048,ndim=3,max_dim=3)
28 INTEGER natoms
29
30 ! Global arrays
31
32 double precision mass(natoms)
33 double precision RX(natoms), RY(natoms), RZ(natoms) !Position
34 double precision VX(natoms), VY(natoms), VZ(natoms) !Velocity
35 double precision FX(natoms), FY(natoms), FZ(natoms) !force
36
37 ! Paramaters for force subroutines
38
39 double precision dt
40
41 ! variables for the constraint methods
42
43 integer n_constrained
44 integer c_i(n_constrained), c_j(n_constrained)
45 double precision constraints_sqr(n_constrained)
46 logical done
47 logical moving(natoms), moved(natoms)
48 double precision rptol, tol, tol2
49 double precision pxab, pyab, pzab, pabsq, rabsq, diffsq
50 double precision rxab, ryab, rzab, rpab, gab
51 double precision dx, dy, dz, rma, rmb
52 integer a, b, it, maxit
53 parameter (rptol = 1.0d-6, tol = 1.0d-6 )
54 parameter ( maxit = 100 )
55
56 double precision px(natoms), py(natoms), pz(natoms)
57
58 double precision e_convert
59 parameter ( e_convert = 4.184d-4 )
60
61
62 DT2 = DT / 2.0d0
63 DTSQ2 = DT * DT2
64
65 DO I = 1, natoms
66
67 px(I) = RX(I) + DT * VX(I) + &
68 (DTSQ2 * FX(I) * e_convert / mass(i))
69 py(I) = RY(I) + DT * VY(I) + &
70 (DTSQ2 * FY(I) * e_convert / mass(i))
71 pz(I) = RZ(I) + DT * VZ(I) + &
72 (DTSQ2 * FZ(I) * e_convert / mass(i))
73 VX(I) = VX(I) + (DT2 * FX(I) * e_convert / mass(i))
74 VY(I) = VY(I) + (DT2 * FY(I) * e_convert / mass(i))
75 VZ(I) = VZ(I) + (DT2 * FZ(I) * e_convert / mass(i))
76
77 moving(i) = .false.
78 moved(i) = .true.
79
80 end do
81
82 it = 0
83 done = .false.
84
85 ! start iterative loop
86
87 do while( (.not. done) .and. (it .le. maxit))
88
89 done = .true.
90
91 do i = 1, n_constrained
92
93 a = c_i(i)
94 b = c_j(i)
95
96 if( moved(a) .or. moved(b) ) then
97
98 pxab = px(a) - px(b)
99 pyab = py(a) - py(b)
100 pzab = pz(a) - pz(b)
101
102 pxab = pxab - box_x * dsign( 1.0d0, pxab ) &
103 * int( dabs( pxab / box_x ) + 0.5d0 )
104 pyab = pyab - box_y * dsign( 1.0d0, pyab ) &
105 * int( dabs( pyab / box_y ) + 0.5d0 )
106 pzab = pzab - box_z * dsign( 1.0d0, pzab ) &
107 * int( dabs( pzab / box_z ) + 0.5d0 )
108
109
110 pabsq = pxab * pxab + pyab * pyab + pzab * pzab
111 rabsq = constraints_sqr(i)
112 diffsq = rabsq - pabsq
113 ! write ( *,* )diffsq
114
115 if( dabs(diffsq) .gt. tol ) then
116
117 rxab = rx(a) - rx(b)
118 ryab = ry(a) - ry(b)
119 rzab = rz(a) - rz(b)
120
121 rxab = rxab - box_x * dsign( 1.0d0, rxab ) &
122 * int( dabs( rxab / box_x ) + 0.5d0 )
123 ryab = ryab - box_y * dsign( 1.0d0, ryab ) &
124 * int( dabs( ryab / box_y ) + 0.5d0 )
125 rzab = rzab - box_z * dsign( 1.0d0, rzab ) &
126 * int( dabs( rzab / box_z ) + 0.5d0 )
127
128 rpab = rxab * pxab + ryab * pyab + rzab * pzab
129
130 if( dabs(rpab) .lt. ( rabsq * rptol ) ) then
131 write (*, '('' Constraint Failure '')' )
132 write (*,*) a-1, b-1,rpab, rabsq * rptol
133 stop
134 end if
135
136 rma = 1.0d0 / mass(a)
137 rmb = 1.0d0 / mass(b)
138
139 gab = diffsq / ( 2.0d0 * ( rma + rmb ) * rpab )
140 dx = rxab * gab
141 dy = ryab * gab
142 dz = rzab * gab
143
144 px(a) = px(a) + rma * dx
145 py(a) = py(a) + rma * dy
146 pz(a) = pz(a) + rma * dz
147
148 px(b) = px(b) - rmb * dx
149 py(b) = py(b) - rmb * dy
150 pz(b) = pz(b) - rmb * dz
151
152 dx = dx / dt
153 dy = dy / dt
154 dz = dz / dt
155
156 vx(a) = vx(a) + rma * dx
157 vy(a) = vy(a) + rma * dy
158 vz(a) = vz(a) + rma * dz
159
160 vx(b) = vx(b) - rmb * dx
161 vy(b) = vy(b) - rmb * dy
162 vz(b) = vz(b) - rmb * dz
163
164 moving(a) = .true.
165 moving(b) = .true.
166 done = .false.
167
168 endif
169 endif
170
171 enddo
172
173 do i = 1, natoms
174 moved(i) = moving(i)
175 moving(i) = .false.
176 enddo
177
178 it = it + 1
179 enddo
180
181 ! end of iterative loop
182
183 if( .not. done) then
184 write (*, '('' too many constraint iterations in move_a '')' )
185 stop
186 endif
187
188 ! store new values
189
190 do i = 1, natoms
191 rx(i) = px(i)
192 ry(i) = py(i)
193 rz(i) = pz(i)
194 enddo
195
196
197
198 RETURN
199 end subroutine v_constrain_a
200
201 Subroutine v_constrain_b(dt,natoms,mass,rx,ry,rz, &
202 vx,vy,vz,fx,fy,fz,k, n_constrained, constraints_sqr, &
203 c_i, c_j, box_x, box_y, box_z)
204 implicit none
205
206 ! *******************************************************************
207 ! ** SECOND PART OF VELOCITY VERLET ALGORITHM **
208 ! ** **
209 ! ** USAGE: **
210 ! ** **
211 ! ** THE SECOND PART OF THE ALGORITHM ADVANCES VELOCITIES FROM **
212 ! ** T + DT/2 TO T + DT. THIS ASSUMES THAT FORCES HAVE BEEN **
213 ! ** COMPUTED IN THE FORCE ROUTINE AND STORED IN FX, FY, FZ. **
214 ! *******************************************************************
215
216 ! declarations
217
218 integer i
219 double precision accvel2, dt2
220 double precision box_x, box_y, box_z
221
222 ! Global Parameters
223
224
225 INTEGER natoms
226
227 ! Global arrays
228
229 double precision mass(natoms)
230 double precision RX(natoms), RY(natoms), RZ(natoms) !position
231 double precision VX(natoms), VY(natoms), VZ(natoms) !velocity
232 double precision FX(natoms), FY(natoms), FZ(natoms) !force
233
234
235 ! Declaration for the kinetic energy
236
237 double precision k
238
239 ! Paramaters for force subroutines
240
241 double precision dt
242
243 ! constraint parameters and variables
244
245 double precision tol, rxab, ryab, rzab, gab
246
247 double precision vxab, vyab, vzab, dx, dy, dz, rma, rmb, rvab
248 integer a, b, it, maxit, n_constrained
249 parameter ( maxit = 100 )
250 logical done
251 logical moving(natoms), moved(natoms)
252 integer c_i(n_constrained), c_j(n_constrained)
253 double precision constraints_sqr(n_constrained)
254
255 double precision e_convert
256 parameter ( e_convert = 4.184d-4 )
257
258
259 ! *******************************************************************
260
261 tol = 1.0d-6 / dt
262 dt2 = dt / 2.0d0
263 accvel2 = 0.0d0
264 k = 0.0d0
265 DO i = 1, natoms
266
267 VX(I) = VX(I) + (DT2 * FX(I) * e_convert / mass(i))
268 VY(I) = VY(I) + (DT2 * FY(I) * e_convert / mass(i))
269 VZ(I) = VZ(I) + (DT2 * FZ(I) * e_convert / mass(i))
270
271 moving(i) = .false.
272 moved(i) = .true.
273 enddo
274
275 it = 0
276 done = .false.
277
278 do while( (.not. done) .and. ( it .le. maxit) )
279
280 done = .true.
281
282 do i=1, n_constrained
283
284 a = c_i(i)
285 b = c_j(i)
286
287 if( moved(a) .or. moved(b) ) then
288
289 vxab = vx(a) - vx(b)
290 vyab = vy(a) - vy(b)
291 vzab = vz(a) - vz(b)
292
293 rxab = rx(a) - rx(b)
294 ryab = ry(a) - ry(b)
295 rzab = rz(a) - rz(b)
296
297 rxab = rxab - box_x * dsign( 1.0d0, rxab ) &
298 * int( dabs( rxab / box_x ) + 0.5d0 )
299 ryab = ryab - box_y * dsign( 1.0d0, ryab ) &
300 * int( dabs( ryab / box_y ) + 0.5d0 )
301 rzab = rzab - box_z * dsign( 1.0d0, rzab ) &
302 * int( dabs( rzab / box_z ) + 0.5d0 )
303
304
305 rma = 1.0d0 / mass(a)
306 rmb = 1.0d0 / mass(b)
307
308 rvab = rxab * vxab + ryab * vyab + rzab * vzab
309
310 gab = -rvab / ( ( rma + rmb ) * constraints_sqr(i) )
311
312 if ( dabs(gab) .gt. tol) then
313
314 dx = rxab * gab
315 dy = ryab * gab
316 dz = rzab * gab
317
318 vx(a) = vx(a) + rma * dx
319 vy(a) = vy(a) + rma * dy
320 vz(a) = vz(a) + rma * dz
321
322 vx(b) = vx(b) - rmb * dx
323 vy(b) = vy(b) - rmb * dy
324 vz(b) = vz(b) - rmb * dz
325
326 moving(a) = .true.
327 moving(b) = .true.
328 done = .false.
329
330 endif
331
332 endif
333
334 enddo
335
336 do i = 1, natoms
337
338 moved(i) = moving(i)
339 moving(i) = .false.
340
341 enddo
342
343 it = it + 1
344
345 enddo
346
347 ! End of iterative loop
348
349 if (.not. done) then
350
351 write(*, '('' Too many constraint iterations in moveb '')')
352 stop
353
354 endif
355
356 do i = 1, natoms
357
358 accvel2 = VX(I) ** 2 + VY(I) ** 2 + VZ(I) ** 2
359 k = k + 0.5d0 * mass(i) * accvel2 / e_convert
360
361 end do
362
363
364 end subroutine v_constrain_b