ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_dipole_dipole.F90
Revision: 626
Committed: Wed Jul 16 21:30:56 2003 UTC (20 years, 11 months ago) by mmeineke
File size: 6500 byte(s)
Log Message:
Changed how cutoffs were handled from C. Now notifyCutoffs in Fortran notifies those who need the information of any changes to cutoffs.

File Contents

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