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 388 by chuckv, Fri Mar 21 22:11:50 2003 UTC vs.
Revision 483 by gezelter, Wed Apr 9 04:06:43 2003 UTC

# Line 4 | Line 4 | module dipole_dipole
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 <  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 24 | Line 25 | contains
25      rrfsq = rrf * rrf
26      ! pre converts from mu in units of debye to kcal/mol
27      pre = 14.38362_dp
28 <    
28 >
29      dipole_initialized = .true.
30      
31      return
# Line 38 | 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
46  
47      real( kind = dp ) :: pot
48      real( kind = dp ), dimension(3) :: d
49 <    real( kind = dp ), dimension(:,:) :: u_l
50 <    real( kind = dp ), dimension(:,:) :: f
51 <    real( kind = dp ), dimension(:,:) :: t
49 >    real( kind = dp ), dimension(3,getNlocal()) :: u_l
50 >    real( kind = dp ), dimension(3,getNlocal()) :: f
51 >    real( kind = dp ), dimension(3,getNlocal()) :: t
52      
53      real (kind = dp), dimension(3) :: ul1
54      real (kind = dp), dimension(3) :: ul2
# Line 79 | Line 80 | contains
80      ul2(3) = u_l(3,atom2)
81   #endif
82  
83 +    if( atom1 .eq. 2 )then
84 +       write (0,*) 'ul =', ul1(1), ul1(2), ul1(3)
85 +    endif
86 +
87 +    if( atom2 .eq. 2 )then
88 +       write (0,*) 'ul =', ul2(1), ul2(2), ul2(3)
89 +    endif
90 +
91 +
92      call getElementProperty(atypes, me1, "dipole_moment", mu1)
93      call getElementProperty(atypes, me2, "dipole_moment", mu2)
94  
# Line 91 | Line 101 | contains
101            taper = (rrf + 2.0d0*rij - 3.0d0*rt)*(rrf-rij)**2/ ((rrf-rt)**3)
102            dtdr = 6.0d0*(rij*rij - rij*rt - rij*rrf +rrf*rt)/((rrf-rt)**3)
103         endif
104 <      
104 >
105         r3 = r2*rij
106         r5 = r3*r2
107        
108         rdotu1 = d(1)*ul1(1) + d(2)*ul1(2) + d(3)*ul1(3)
109         rdotu2 = d(1)*ul2(1) + d(2)*ul2(2) + d(3)*ul2(3)
110         u1dotu2 = ul1(1)*ul2(1) + ul1(2)*ul2(2) + ul1(3)*ul2(3)
111 <      
111 >
112         dip2 = pre * mu1 * mu2
103      
113         dfact1 = 3.0d0*dip2 / r2
114         dfact2 = 3.0d0*dip2 / r5
115 <      
115 >
116         vterm = dip2*((u1dotu2/r3) - 3.0d0*(rdotu1*rdotu2/r5))
117 <      
117 >
118         if (do_pot) then
119   #ifdef IS_MPI
120            pot_row(atom1) = pot_row(atom1) + 0.5d0*vterm*taper
# Line 129 | Line 138 | contains
138              (5.0d0*(rdotu1*rdotu2)/r5)) -  &
139              dfact2*(ul1(3)*rdotu2 + ul2(3)*rdotu1))*taper +  &
140              vterm*dtdr*d(3)/rij
141 <      
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 <      
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 <      
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
# Line 174 | Line 183 | contains
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))
186 >
187 >          if (molMembershipList(atom1) .ne. molMembershipList(atom2)) then
188 >
189 >             tau_Temp(1) = tau_Temp(1) + dudx * d(1)
190 >             tau_Temp(2) = tau_Temp(2) + dudx * d(2)
191 >             tau_Temp(3) = tau_Temp(3) + dudx * d(3)
192 >             tau_Temp(4) = tau_Temp(4) + dudy * d(1)
193 >             tau_Temp(5) = tau_Temp(5) + dudy * d(2)
194 >             tau_Temp(6) = tau_Temp(6) + dudy * d(3)
195 >             tau_Temp(7) = tau_Temp(7) + dudz * d(1)
196 >             tau_Temp(8) = tau_Temp(8) + dudz * d(2)
197 >             tau_Temp(9) = tau_Temp(9) + dudz * d(3)
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