ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/f_verlet_constrained.F90
Revision: 505
Committed: Fri Apr 25 02:00:11 2003 UTC (21 years, 2 months ago) by mmeineke
File size: 10249 byte(s)
Log Message:
added a new test for constraint failure

File Contents

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