ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/f_verlet_constrained.F90
Revision: 542
Committed: Fri May 30 21:31:48 2003 UTC (21 years, 1 month ago) by mmeineke
File size: 10382 byte(s)
Log Message:
currently modifiying Symplectic to become the basic integrator.

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