ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/f_verlet_constrained.F90
Revision: 388
Committed: Fri Mar 21 22:11:50 2003 UTC (21 years, 3 months ago) by chuckv
File size: 10199 byte(s)
Log Message:
various write statements for debugging

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