ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/calc_reaction_field.F90
Revision: 307
Committed: Mon Mar 10 19:37:48 2003 UTC (21 years, 5 months ago) by chuckv
File size: 3599 byte(s)
Log Message:
Changed names of calculation modules for do_Forces.

File Contents

# User Rev Content
1 chuckv 307 subroutine accumulate_rf(atom1, atom2, rij)
2    
3     include 'sizes.inc'
4     include 'simulation.inc'
5    
6     integer atom1, atom2
7     double precision taper, rij
8    
9     if (rij.le.rrf) then
10    
11     if (rij.lt.rt) then
12     taper = 1.0d0
13     else
14     taper = (rrf + 2.0d0*rij - 3.0d0*rt)*(rrf-rij)**2/ ((rrf-rt)**3)
15     endif
16    
17     rflx(atom1) = rflx(atom1) + ulx(atom2)*mu(atom2)*taper
18     rfly(atom1) = rfly(atom1) + uly(atom2)*mu(atom2)*taper
19     rflz(atom1) = rflz(atom1) + ulz(atom2)*mu(atom2)*taper
20    
21     rflx(atom2) = rflx(atom2) + ulx(atom1)*mu(atom1)*taper
22     rfly(atom2) = rfly(atom2) + uly(atom1)*mu(atom1)*taper
23     rflz(atom2) = rflz(atom2) + ulz(atom1)*mu(atom1)*taper
24    
25     endif
26     return
27    
28     end subroutine accumulate_rf
29    
30     subroutine accumulate_self_rf()
31    
32     include 'sizes.inc'
33     include 'simulation.inc'
34    
35     integer i, ia, a1, atype1
36     logical is_dipole_atype
37     external is_dipole_atype
38    
39     do i = 1, nmol
40     do ia = 1, na(i)
41     a1 = atom_index(i,ia)
42    
43     atype1 = atype(a1)
44    
45     if (is_dipole_atype(atype1)) then
46    
47     rflx(a1) = rflx(a1) + ulx(a1)*mu(a1)
48     rfly(a1) = rfly(a1) + uly(a1)*mu(a1)
49     rflz(a1) = rflz(a1) + ulz(a1)*mu(a1)
50    
51     endif
52     enddo
53     enddo
54    
55     return
56     end subroutine accumulate_self_rf
57    
58     subroutine reaction_field(pot)
59    
60     include 'sizes.inc'
61     include 'simulation.inc'
62    
63     double precision rrfsq, pre
64     integer i, ia, a1, atype1
65     logical is_dipole_atype
66     external is_dipole_atype
67    
68     ! do single loop to compute torques on dipoles:
69     ! pre converts from mu in units of debye to kcal/mol
70    
71     rrfsq = rrf * rrf
72     pre = 14.38362d0*2.0d0*(dielect-1.0d0)/((2.0d0*dielect+1.0d0)*rrfsq*rrf)
73    
74     do i = 1, nmol
75     do ia = 1, na(i)
76     a1 = atom_index(i,ia)
77    
78     atype1 = atype(a1)
79    
80     if (is_dipole_atype(atype1)) then
81    
82     ! The torque contribution is dipole cross reaction_field
83    
84     tlx(a1) = tlx(a1) + pre*mu(a1)*(uly(a1)*rflz(a1) - ulz(a1)*rfly(a1))
85     tly(a1) = tly(a1) + pre*mu(a1)*(ulz(a1)*rflx(a1) - ulx(a1)*rflz(a1))
86     tlz(a1) = tlz(a1) + pre*mu(a1)*(ulx(a1)*rfly(a1) - uly(a1)*rflx(a1))
87    
88     ! the potential contribution is -1/2 dipole dot reaction_field
89    
90     pot = pot - 0.5d0 * pre * mu(a1) * &
91     (rflx(a1)*ulx(a1) + rfly(a1)*uly(a1) + rflz(a1)*ulz(a1))
92    
93     endif
94    
95     enddo
96     enddo
97    
98     end subroutine reaction_field
99    
100     subroutine rf_correct_forces(atom1, atom2, dx, dy, dz, rij)
101     include 'sizes.inc'
102     include 'simulation.inc'
103    
104     integer atom1, atom2
105     double precision dtdr, rrfsq, prerf, rij
106     double precision dudx, dudy, dudz, u1dotu2, dx, dy, dz
107    
108     rrfsq = rrf * rrf
109     prerf = 14.38362d0*2.0d0*(dielect-1.0d0)/((2.0d0*dielect+1.0d0)*rrfsq*rrf)
110    
111     if (rij.le.rrf) then
112    
113     ! cubic taper
114     if (rij.lt.rt) then
115     dtdr = 0.0d0
116     else
117     dtdr = 6.0d0*(rij*rij - rij*rt - rij*rrf +rrf*rt)/((rrf-rt)**3)
118     endif
119    
120     u1dotu2 = ulx(atom1)*ulx(atom2) + uly(atom1)*uly(atom2) + &
121     ulz(atom1)*ulz(atom2)
122    
123     dudx = - prerf*mu(atom1)*mu(atom2)*u1dotu2*dtdr*dx/rij
124     dudy = - prerf*mu(atom1)*mu(atom2)*u1dotu2*dtdr*dy/rij
125     dudz = - prerf*mu(atom1)*mu(atom2)*u1dotu2*dtdr*dz/rij
126    
127     flx(atom1) = flx(atom1) + dudx
128     fly(atom1) = fly(atom1) + dudy
129     flz(atom1) = flz(atom1) + dudz
130    
131     flx(atom2) = flx(atom2) - dudx
132     fly(atom2) = fly(atom2) - dudy
133     flz(atom2) = flz(atom2) - dudz
134    
135     ! add contribution to the virial as well
136     virial = virial + ( dx*dudx + dy*dudy + dz*dudz )
137    
138     endif
139    
140     return
141     end subroutine rf_correct_forces