ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_dipole_dipole.F90
Revision: 443
Committed: Wed Apr 2 22:19:03 2003 UTC (21 years, 3 months ago) by mmeineke
File size: 6303 byte(s)
Log Message:
dipoles mostly work, but there is a memory leak somewhere.

File Contents

# Content
1 module dipole_dipole
2
3 use force_globals
4 use definitions
5 use atype_module
6 use vector_class
7 #ifdef IS_MPI
8 use mpiSimulation
9 #endif
10 implicit none
11
12 real(kind=dp), save :: rrf
13 real(kind=dp), save :: rt
14 real(kind=dp), save :: rrfsq
15 real(kind=dp), save :: pre
16 logical, save :: dipole_initialized = .false.
17
18 contains
19
20 subroutine initialize_dipole(this_rrf, this_rt)
21 real(kind=dp), intent(in) :: this_rrf, this_rt
22 rrf = this_rrf
23 rt = this_rt
24 rrfsq = rrf * rrf
25 ! pre converts from mu in units of debye to kcal/mol
26 pre = 14.38362_dp
27
28 dipole_initialized = .true.
29
30 return
31 end subroutine initialize_dipole
32
33
34 subroutine do_dipole_pair(atom1, atom2, d, rij, r2, pot, u_l, f, t, &
35 do_pot, do_stress)
36
37 logical :: do_pot, do_stress
38
39 integer atom1, atom2, me1, me2
40 real(kind=dp) :: rij, mu1, mu2
41 real(kind=dp) :: dfact1, dfact2, dip2, r2, r3, r5, pre
42 real(kind=dp) :: dudx, dudy, dudz, dudu1x, dudu1y, dudu1z
43 real(kind=dp) :: dudu2x, dudu2y, dudu2z, rdotu1, rdotu2, u1dotu2
44 real(kind=dp) :: taper, dtdr, vterm
45
46 real( kind = dp ) :: pot
47 real( kind = dp ), dimension(3) :: d
48 real( kind = dp ), dimension(:,:) :: u_l
49 real( kind = dp ), dimension(:,:) :: f
50 real( kind = dp ), dimension(:,:) :: t
51
52 real (kind = dp), dimension(3) :: ul1
53 real (kind = dp), dimension(3) :: ul2
54
55 if (.not.dipole_initialized) then
56 write(default_error,*) 'Dipole-dipole not initialized!'
57 return
58 endif
59
60 #ifdef IS_MPI
61 me1 = atid_Row(atom1)
62 ul1(1) = u_l_Row(1,atom1)
63 ul1(2) = u_l_Row(2,atom1)
64 ul1(3) = u_l_Row(3,atom1)
65
66 me2 = atid_Col(atom2)
67 ul2(1) = u_l_Col(1,atom2)
68 ul2(2) = u_l_Col(2,atom2)
69 ul2(3) = u_l_Col(3,atom2)
70 #else
71 me1 = atid(atom1)
72 ul1(1) = u_l(1,atom1)
73 ul1(2) = u_l(2,atom1)
74 ul1(3) = u_l(3,atom1)
75
76 me2 = atid(atom2)
77 ul2(1) = u_l(1,atom2)
78 ul2(2) = u_l(2,atom2)
79 ul2(3) = u_l(3,atom2)
80 #endif
81
82 if( atom1 .eq. 2 )then
83 write (*,*) 'ul =', ul1(1), ul1(2), ul1(3)
84 endif
85
86 if( atom2 .eq. 2 )then
87 write (*,*) 'ul =', ul2(1), ul2(2), ul2(3)
88 endif
89
90
91 call getElementProperty(atypes, me1, "dipole_moment", mu1)
92 call getElementProperty(atypes, me2, "dipole_moment", mu2)
93
94 if (rij.le.rrf) then
95
96 if (rij.lt.rt) then
97 taper = 1.0d0
98 dtdr = 0.0d0
99 else
100 taper = (rrf + 2.0d0*rij - 3.0d0*rt)*(rrf-rij)**2/ ((rrf-rt)**3)
101 dtdr = 6.0d0*(rij*rij - rij*rt - rij*rrf +rrf*rt)/((rrf-rt)**3)
102 endif
103
104 r3 = r2*rij
105 r5 = r3*r2
106
107 rdotu1 = d(1)*ul1(1) + d(2)*ul1(2) + d(3)*ul1(3)
108 rdotu2 = d(1)*ul2(1) + d(2)*ul2(2) + d(3)*ul2(3)
109 u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
110
111 dip2 = pre * mu1 * mu2
112
113 dfact1 = 3.0d0*dip2 / r2
114 dfact2 = 3.0d0*dip2 / r5
115
116 vterm = dip2*((u1dotu2/r3) - 3.0d0*(rdotu1*rdotu2/r5))
117
118 if (do_pot) then
119 #ifdef IS_MPI
120 pot_row(atom1) = pot_row(atom1) + 0.5d0*vterm*taper
121 pot_col(atom2) = pot_col(atom2) + 0.5d0*vterm*taper
122 #else
123 pot = pot + vterm*taper
124 #endif
125 endif
126
127 dudx = (-dfact1 * d(1) * ((u1dotu2/r3) - &
128 (5.0d0*(rdotu1*rdotu2)/r5)) - &
129 dfact2*(ul1(1)*rdotu2 + ul2(1)*rdotu1))*taper + &
130 vterm*dtdr*d(1)/rij
131
132 dudy = (-dfact1 * d(2) * ((u1dotu2/r3) - &
133 (5.0d0*(rdotu1*rdotu2)/r5)) - &
134 dfact2*(ul1(2)*rdotu2 + ul2(2)*rdotu1))*taper + &
135 vterm*dtdr*d(2)/rij
136
137 dudz = (-dfact1 * d(3) * ((u1dotu2/r3) - &
138 (5.0d0*(rdotu1*rdotu2)/r5)) - &
139 dfact2*(ul1(3)*rdotu2 + ul2(3)*rdotu1))*taper + &
140 vterm*dtdr*d(3)/rij
141
142 dudu1x = (dip2*((ul2(1)/r3) - (3.0d0*d(1)*rdotu2/r5)))*taper
143 dudu1y = (dip2*((ul2(2)/r3) - (3.0d0*d(2)*rdotu2/r5)))*taper
144 dudu1z = (dip2*((ul2(3)/r3) - (3.0d0*d(3)*rdotu2/r5)))*taper
145
146 dudu2x = (dip2*((ul1(1)/r3) - (3.0d0*d(1)*rdotu1/r5)))*taper
147 dudu2y = (dip2*((ul1(2)/r3) - (3.0d0*d(2)*rdotu1/r5)))*taper
148 dudu2z = (dip2*((ul1(3)/r3) - (3.0d0*d(3)*rdotu1/r5)))*taper
149
150
151 #ifdef IS_MPI
152 f_Row(1,atom1) = f_Row(1,atom1) - dudx
153 f_Row(2,atom1) = f_Row(2,atom1) - dudy
154 f_Row(3,atom1) = f_Row(3,atom1) - dudz
155
156 f_Col(1,atom2) = f_Col(1,atom2) + dudx
157 f_Col(2,atom2) = f_Col(2,atom2) + dudy
158 f_Col(3,atom2) = f_Col(3,atom2) + dudz
159
160 t_Row(1,atom1) = t_Row(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
161 t_Row(2,atom1) = t_Row(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
162 t_Row(3,atom1) = t_Row(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
163
164 t_Col(1,atom2) = t_Col(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
165 t_Col(2,atom2) = t_Col(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
166 t_Col(3,atom2) = t_Col(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
167 #else
168 f(1,atom1) = f(1,atom1) - dudx
169 f(2,atom1) = f(2,atom1) - dudy
170 f(3,atom1) = f(3,atom1) - dudz
171
172 f(1,atom2) = f(1,atom2) + dudx
173 f(2,atom2) = f(2,atom2) + dudy
174 f(3,atom2) = f(3,atom2) + dudz
175
176 t(1,atom1) = t(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
177 t(2,atom1) = t(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
178 t(3,atom1) = t(3,atom1) - ul1(1)*dudu1y + ul1(2)*dudu1x
179
180 t(1,atom2) = t(1,atom2) - ul2(2)*dudu2z + ul2(3)*dudu2y
181 t(2,atom2) = t(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
182 t(3,atom2) = t(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
183 #endif
184
185 if (do_stress) then
186 tau_Temp(1) = tau_Temp(1) + dudx * d(1)
187 tau_Temp(2) = tau_Temp(2) + dudx * d(2)
188 tau_Temp(3) = tau_Temp(3) + dudx * d(3)
189 tau_Temp(4) = tau_Temp(4) + dudy * d(1)
190 tau_Temp(5) = tau_Temp(5) + dudy * d(2)
191 tau_Temp(6) = tau_Temp(6) + dudy * d(3)
192 tau_Temp(7) = tau_Temp(7) + dudz * d(1)
193 tau_Temp(8) = tau_Temp(8) + dudz * d(2)
194 tau_Temp(9) = tau_Temp(9) + dudz * d(3)
195 virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
196 endif
197
198 endif
199
200 return
201 end subroutine do_dipole_pair
202
203 end module dipole_dipole