ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_dipole_dipole.F90
Revision: 394
Committed: Mon Mar 24 21:55:34 2003 UTC (21 years, 3 months ago) by gezelter
File size: 6217 byte(s)
Log Message:
electrostatic changes for dipole / RF separation

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