ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/mdtools/md_code/f_verlet_constrained.f
Revision: 11
Committed: Tue Jul 9 18:40:59 2002 UTC (21 years, 11 months ago) by mmeineke
File size: 11048 byte(s)
Log Message:
This commit was generated by cvs2svn to compensate for changes in r10, which
included commits to RCS files with non-trunk default branches.

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