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

Comparing:
branches/mmeineke/OOPSE/libmdtools/calc_reaction_field.F90 (file contents), Revision 377 by mmeineke, Fri Mar 21 17:42:12 2003 UTC vs.
trunk/OOPSE/libmdtools/calc_reaction_field.F90 (file contents), Revision 621 by gezelter, Wed Jul 16 02:11:02 2003 UTC

# Line 3 | Line 3 | module reaction_field
3    use definitions
4    use atype_module
5    use vector_class
6 +  use simulation
7   #ifdef IS_MPI
8    use mpiSimulation
9   #endif
# Line 27 | Line 28 | contains
28      rrfsq = rrf * rrf
29      pre = 14.38362d0*2.0d0*(dielect-1.0d0)/((2.0d0*dielect+1.0d0)*rrfsq*rrf)
30      
31 +
32 +    write(*,*) 'rrf = ', rrf
33 +    write(*,*) 'rt = ', rt
34 +    write(*,*) 'dielect = ', dielect
35 +    write(*,*) 'pre = ', pre
36      rf_initialized = .true.
37      
38      return
# Line 36 | Line 42 | contains
42  
43      integer, intent(in) :: atom1, atom2
44      real (kind = dp), intent(in) :: rij
45 <    real (kind = dp), dimension(:,:) :: u_l    
45 >    real (kind = dp), dimension(3,getNlocal()) :: u_l    
46  
47      integer :: me1, me2
48      real (kind = dp) :: taper, mu1, mu2
# Line 53 | Line 59 | contains
59         if (rij.lt.rt) then
60            taper = 1.0d0
61         else
62 +          write(*,*) 'rf in taper region'
63            taper = (rrf + 2.0d0*rij - 3.0d0*rt)*(rrf-rij)**2/ ((rrf-rt)**3)
64         endif
65        
# Line 87 | Line 94 | contains
94         rf_Row(2,atom1) = rf_Row(2,atom1) + ul2(2)*mu2*taper
95         rf_Row(3,atom1) = rf_Row(3,atom1) + ul2(3)*mu2*taper
96        
97 <       rf_Col(1,atom2) = rf_Col(1,atom2) + ul1(1)*mu2*taper
98 <       rf_Col(2,atom2) = rf_Col(2,atom2) + ul1(2)*mu2*taper
99 <       rf_Col(3,atom2) = rf_Col(3,atom2) + ul1(3)*mu2*taper
97 >       rf_Col(1,atom2) = rf_Col(1,atom2) + ul1(1)*mu1*taper
98 >       rf_Col(2,atom2) = rf_Col(2,atom2) + ul1(2)*mu1*taper
99 >       rf_Col(3,atom2) = rf_Col(3,atom2) + ul1(3)*mu1*taper
100   #else
101         rf(1,atom1) = rf(1,atom1) + ul2(1)*mu2*taper
102         rf(2,atom1) = rf(2,atom1) + ul2(2)*mu2*taper
103         rf(3,atom1) = rf(3,atom1) + ul2(3)*mu2*taper
104        
105 <       rf(1,atom2) = rf(1,atom2) + ul1(1)*mu2*taper
106 <       rf(2,atom2) = rf(2,atom2) + ul1(2)*mu2*taper
107 <       rf(3,atom2) = rf(3,atom2) + ul1(3)*mu2*taper    
105 >       rf(1,atom2) = rf(1,atom2) + ul1(1)*mu1*taper
106 >       rf(2,atom2) = rf(2,atom2) + ul1(2)*mu1*taper
107 >       rf(3,atom2) = rf(3,atom2) + ul1(3)*mu1*taper    
108   #endif
109        
110      endif
# Line 108 | Line 115 | contains
115      
116      integer, intent(in) :: atom1
117      real(kind=dp), intent(in) :: mu1
118 <    real(kind=dp), dimension(:,:) :: u_l
118 >    real(kind=dp), dimension(3,getNlocal()) :: u_l
119      
120      !! should work for both MPI and non-MPI version since this is not pairwise.
121      rf(1,atom1) = rf(1,atom1) + u_l(1,atom1)*mu1
# Line 124 | Line 131 | contains
131      real (kind=dp), intent(in) :: mu1
132      real (kind=dp), intent(inout) :: rfpot
133      logical, intent(in) :: do_pot
134 <    real (kind = dp), dimension(:,:) :: u_l    
135 <    real (kind = dp), dimension(:,:) :: t
134 >    real (kind = dp), dimension(3,getNlocal()) :: u_l    
135 >    real (kind = dp), dimension(3,getNlocal()) :: t
136  
137      if (.not.rf_initialized) then
138         write(default_error,*) 'Reaction field not initialized!'
# Line 135 | Line 142 | contains
142      ! compute torques on dipoles:
143      ! pre converts from mu in units of debye to kcal/mol
144      
145 <    ! The torque contribution is dipole cross reaction_field
146 <    
145 >    ! The torque contribution is dipole cross reaction_field  
146 >
147      t(1,a1) = t(1,a1) + pre*mu1*(u_l(2,a1)*rf(3,a1) - u_l(3,a1)*rf(2,a1))
148      t(2,a1) = t(2,a1) + pre*mu1*(u_l(3,a1)*rf(1,a1) - u_l(1,a1)*rf(3,a1))
149      t(3,a1) = t(3,a1) + pre*mu1*(u_l(1,a1)*rf(2,a1) - u_l(2,a1)*rf(1,a1))
# Line 156 | Line 163 | contains
163      integer, intent(in) :: atom1, atom2
164      real(kind=dp), dimension(3), intent(in) :: d
165      real(kind=dp), intent(in) :: rij
166 <    real( kind = dp ), dimension(:,:) :: u_l
167 <    real( kind = dp ), dimension(:,:) :: f
166 >    real( kind = dp ), dimension(3,getNlocal()) :: u_l
167 >    real( kind = dp ), dimension(3,getNlocal()) :: f
168      logical, intent(in) :: do_stress
169      
170      real (kind = dp), dimension(3) :: ul1
# Line 177 | Line 184 | contains
184         if (rij.lt.rt) then
185            dtdr = 0.0d0
186         else
187 +          write(*,*) 'rf correct in taper region'
188            dtdr = 6.0d0*(rij*rij - rij*rt - rij*rrf +rrf*rt)/((rrf-rt)**3)
189         endif
190        
# Line 230 | Line 238 | contains
238   #endif
239        
240         if (do_stress) then          
241 <          tau_Temp(1) = tau_Temp(1) + dudx * d(1)
242 <          tau_Temp(2) = tau_Temp(2) + dudx * d(2)
243 <          tau_Temp(3) = tau_Temp(3) + dudx * d(3)
244 <          tau_Temp(4) = tau_Temp(4) + dudy * d(1)
245 <          tau_Temp(5) = tau_Temp(5) + dudy * d(2)
246 <          tau_Temp(6) = tau_Temp(6) + dudy * d(3)
247 <          tau_Temp(7) = tau_Temp(7) + dudz * d(1)
248 <          tau_Temp(8) = tau_Temp(8) + dudz * d(2)
249 <          tau_Temp(9) = tau_Temp(9) + dudz * d(3)
250 <          virial_Temp = virial_Temp + (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
241 >          if (molMembershipList(atom1) .ne. molMembershipList(atom2)) then
242 >
243 >             ! because the d vector is the rj - ri vector, and
244 >             ! because dudx, dudy, and dudz are the
245 >             ! (positive) force on atom i (negative on atom j) we need
246 >             ! a negative sign here:
247 >
248 >             tau_Temp(1) = tau_Temp(1) - d(1) * dudx
249 >             tau_Temp(2) = tau_Temp(2) - d(1) * dudy
250 >             tau_Temp(3) = tau_Temp(3) - d(1) * dudz
251 >             tau_Temp(4) = tau_Temp(4) - d(2) * dudx
252 >             tau_Temp(5) = tau_Temp(5) - d(2) * dudy
253 >             tau_Temp(6) = tau_Temp(6) - d(2) * dudz
254 >             tau_Temp(7) = tau_Temp(7) - d(3) * dudx
255 >             tau_Temp(8) = tau_Temp(8) - d(3) * dudy
256 >             tau_Temp(9) = tau_Temp(9) - d(3) * dudz
257 >             virial_Temp = virial_Temp + &
258 >                  (tau_Temp(1) + tau_Temp(5) + tau_Temp(9))
259 >          endif
260         endif      
261      endif
262      

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines