ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/calc_dipole_dipole.F90
(Generate patch)

Comparing trunk/OOPSE/libmdtools/calc_dipole_dipole.F90 (file contents):
Revision 460 by chuckv, Fri Apr 4 22:22:30 2003 UTC vs.
Revision 619 by mmeineke, Tue Jul 15 22:22:41 2003 UTC

# Line 10 | Line 10 | module dipole_dipole
10   #endif
11    implicit none
12  
13 <  real(kind=dp), save :: rrf
14 <  real(kind=dp), save :: rt  
15 <  real(kind=dp), save :: rrfsq
16 <  real(kind=dp), save :: pre
13 >  real(kind=dp), save :: rrf = 0.0
14 >  real(kind=dp), save :: rt  = 0.0
15 >  real(kind=dp), save :: rrfsq = 0.0
16 >  real(kind=dp), save :: pre = 0.0
17    logical, save :: dipole_initialized = .false.
18  
19   contains
# Line 39 | Line 39 | contains
39  
40      integer atom1, atom2, me1, me2
41      real(kind=dp) :: rij, mu1, mu2
42 <    real(kind=dp) :: dfact1, dfact2, dip2, r2, r3, r5, pre
42 >    real(kind=dp) :: dfact1, dfact2, dip2, r2, r3, r5
43      real(kind=dp) :: dudx, dudy, dudz, dudu1x, dudu1y, dudu1z
44      real(kind=dp) :: dudu2x, dudu2y, dudu2z, rdotu1, rdotu2, u1dotu2
45      real(kind=dp) :: taper, dtdr, vterm
# Line 80 | Line 80 | contains
80      ul2(3) = u_l(3,atom2)
81   #endif
82  
83    if( atom1 .eq. 2 )then
84       write (*,*) 'ul =', ul1(1), ul1(2), ul1(3)
85    endif
83  
87    if( atom2 .eq. 2 )then
88       write (*,*) 'ul =', ul2(1), ul2(2), ul2(3)
89    endif
90
91
84      call getElementProperty(atypes, me1, "dipole_moment", mu1)
85      call getElementProperty(atypes, me2, "dipole_moment", mu2)
86  
87 +
88      if (rij.le.rrf) then
89 +
90        
91 +      
92         if (rij.lt.rt) then
93            taper = 1.0d0
94            dtdr = 0.0d0
# Line 101 | Line 96 | contains
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 <      
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 <      
106 >
107         dip2 = pre * mu1 * mu2
113      
108         dfact1 = 3.0d0*dip2 / r2
109         dfact2 = 3.0d0*dip2 / r5
110 <      
110 >
111         vterm = dip2*((u1dotu2/r3) - 3.0d0*(rdotu1*rdotu2/r5))
112 <      
112 >
113         if (do_pot) then
114   #ifdef IS_MPI
115            pot_row(atom1) = pot_row(atom1) + 0.5d0*vterm*taper
# Line 139 | Line 133 | contains
133              (5.0d0*(rdotu1*rdotu2)/r5)) -  &
134              dfact2*(ul1(3)*rdotu2 + ul2(3)*rdotu1))*taper +  &
135              vterm*dtdr*d(3)/rij
136 <      
136 >
137         dudu1x = (dip2*((ul2(1)/r3) - (3.0d0*d(1)*rdotu2/r5)))*taper
138         dudu1y = (dip2*((ul2(2)/r3) - (3.0d0*d(2)*rdotu2/r5)))*taper
139         dudu1z = (dip2*((ul2(3)/r3) - (3.0d0*d(3)*rdotu2/r5)))*taper
140 <      
140 >
141         dudu2x = (dip2*((ul1(1)/r3) - (3.0d0*d(1)*rdotu1/r5)))*taper
142         dudu2y = (dip2*((ul1(2)/r3) - (3.0d0*d(2)*rdotu1/r5)))*taper
143         dudu2z = (dip2*((ul1(3)/r3) - (3.0d0*d(3)*rdotu1/r5)))*taper
144        
145 <      
145 >
146   #ifdef IS_MPI
147 <       f_Row(1,atom1) = f_Row(1,atom1) - dudx
148 <       f_Row(2,atom1) = f_Row(2,atom1) - dudy
149 <       f_Row(3,atom1) = f_Row(3,atom1) - dudz
147 >       f_Row(1,atom1) = f_Row(1,atom1) + dudx
148 >       f_Row(2,atom1) = f_Row(2,atom1) + dudy
149 >       f_Row(3,atom1) = f_Row(3,atom1) + dudz
150  
151 <       f_Col(1,atom2) = f_Col(1,atom2) + dudx
152 <       f_Col(2,atom2) = f_Col(2,atom2) + dudy
153 <       f_Col(3,atom2) = f_Col(3,atom2) + dudz
151 >       f_Col(1,atom2) = f_Col(1,atom2) - dudx
152 >       f_Col(2,atom2) = f_Col(2,atom2) - dudy
153 >       f_Col(3,atom2) = f_Col(3,atom2) - dudz
154        
155         t_Row(1,atom1) = t_Row(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
156         t_Row(2,atom1) = t_Row(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
# Line 166 | Line 160 | contains
160         t_Col(2,atom2) = t_Col(2,atom2) - ul2(3)*dudu2x + ul2(1)*dudu2z
161         t_Col(3,atom2) = t_Col(3,atom2) - ul2(1)*dudu2y + ul2(2)*dudu2x
162   #else
163 <       f(1,atom1) = f(1,atom1) - dudx
164 <       f(2,atom1) = f(2,atom1) - dudy
165 <       f(3,atom1) = f(3,atom1) - dudz
163 >       f(1,atom1) = f(1,atom1) + dudx
164 >       f(2,atom1) = f(2,atom1) + dudy
165 >       f(3,atom1) = f(3,atom1) + dudz
166        
167 <       f(1,atom2) = f(1,atom2) + dudx
168 <       f(2,atom2) = f(2,atom2) + dudy
169 <       f(3,atom2) = f(3,atom2) + dudz
167 >       f(1,atom2) = f(1,atom2) - dudx
168 >       f(2,atom2) = f(2,atom2) - dudy
169 >       f(3,atom2) = f(3,atom2) - dudz
170        
171         t(1,atom1) = t(1,atom1) - ul1(2)*dudu1z + ul1(3)*dudu1y
172         t(2,atom1) = t(2,atom1) - ul1(3)*dudu1x + ul1(1)*dudu1z
# Line 184 | Line 178 | contains
178   #endif
179        
180         if (do_stress) then          
181 <          tau_Temp(1) = tau_Temp(1) + dudx * d(1)
182 <          tau_Temp(2) = tau_Temp(2) + dudx * d(2)
183 <          tau_Temp(3) = tau_Temp(3) + dudx * d(3)
184 <          tau_Temp(4) = tau_Temp(4) + dudy * d(1)
185 <          tau_Temp(5) = tau_Temp(5) + dudy * d(2)
186 <          tau_Temp(6) = tau_Temp(6) + dudy * d(3)
187 <          tau_Temp(7) = tau_Temp(7) + dudz * d(1)
188 <          tau_Temp(8) = tau_Temp(8) + dudz * d(2)
189 <          tau_Temp(9) = tau_Temp(9) + dudz * d(3)
190 <          virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
181 >
182 >          if (molMembershipList(atom1) .ne. molMembershipList(atom2)) then
183 >
184 >             ! because the d vector is the rj - ri vector, and
185 >             ! because dudx, dudy, dudz are the (positive) force on
186 >             ! atom i (negative on atom j) we need a negative sign here:
187 >
188 >             tau_Temp(1) = tau_Temp(1) - d(1) * dudx
189 >             tau_Temp(2) = tau_Temp(2) - d(1) * dudy
190 >             tau_Temp(3) = tau_Temp(3) - d(1) * dudz
191 >             tau_Temp(4) = tau_Temp(4) - d(2) * dudx
192 >             tau_Temp(5) = tau_Temp(5) - d(2) * dudy
193 >             tau_Temp(6) = tau_Temp(6) - d(2) * dudz
194 >             tau_Temp(7) = tau_Temp(7) - d(3) * dudx
195 >             tau_Temp(8) = tau_Temp(8) - d(3) * dudy
196 >             tau_Temp(9) = tau_Temp(9) - d(3) * dudz
197 >
198 >             virial_Temp = virial_Temp + &
199 >                  (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
200 >
201 >          endif          
202         endif
203        
204      endif

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines