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

Comparing trunk/OOPSE/libmdtools/calc_charge_charge.F90 (file contents):
Revision 1131 by tim, Thu Apr 22 21:33:55 2004 UTC vs.
Revision 1132 by tim, Sat Apr 24 04:31:36 2004 UTC

# Line 14 | Line 14 | module charge_charge
14    PRIVATE
15    real(kind=dp), save :: ecr = 0.0_DP
16    real(kind=dp), save :: rt  = 0.0_DP
17 +  real(kind=dp), save :: dielect  = 0.0_DP
18    real(kind=dp), save :: pre = 0.0_DP
19 +  real(kind=dp), save :: pre_rf = 0.0_DP
20 +  real(kind=dp), save :: ecr3 = 0.0_DP
21 +  real(kind=dp), save :: b0 = 0.0_DP
22 +  real(kind=dp), save :: pre_correction = 0.0_DP
23 +  
24 +  
25    logical, save :: haveCutoffs = .false.
26    logical, save :: haveChargeMap = .false.
27 +  logical, save :: haveDielectric = .false.
28  
29    public::setCutoffsCharge
30    public::do_charge_pair
31 +  public::initialize_charge
32  
33    type :: ChargeList
34       real(kind=DP) :: charge = 0.0_DP
# Line 28 | Line 37 | contains
37    type(ChargeList), dimension(:), allocatable :: ChargeMap
38  
39   contains
40 +
41 +  subroutine initialize_charge(this_dielect)
42 +    real(kind=dp), intent(in) :: this_dielect
43      
44 +    dielect = this_dielect
45 +    haveDielectric = .true.
46 +    
47 +    ! because setCutoffsCharge is called before initialize_charge
48 +    ! we need to call it agin to make sure all of the precalculation
49 +    ! value is correct. for the time being, just a quick hack
50 +    call setCutoffsCharge(ecr, rt)
51 +    return
52 +  end subroutine initialize_charge
53 +  
54 +    
55    subroutine setCutoffsCharge(this_ecr, this_rt)
56      real(kind=dp), intent(in) :: this_ecr, this_rt
57      ecr = this_ecr
# Line 36 | Line 59 | contains
59  
60      ! pre converts from fundamental charge to kcal/mol
61      pre = 332.06508_DP
62 <
62 >    
63 >    !if (.not.haveDielectric)then
64 >    !   write(default_error,*) 'Dielect constant in charge module is not set!'
65 >    !endif
66 >    
67 >    ecr3 = ecr * ecr * ecr
68 >    b0 = 2.0d0*(dielect-1.0d0)/(2.0d0*dielect+1.0d0)
69 >    pre_rf = pre * b0 /(2.0d0* ecr3)
70 >    pre_correction = pre* (1 + 0.5d0 * b0) /ecr
71      haveCutoffs = .true.
72      
73      return
# Line 92 | Line 123 | contains
123      real( kind = dp ) :: pot
124      real( kind = dp ), dimension(3) :: d
125      real( kind = dp ), dimension(3,nLocal) :: f
126 +    real( kind = dp) :: ecr3, b0
127 +    real( kind = dp) :: vcolomb, vrf, vcorrection
128      
129      if (.not. haveCutoffs) then
130         write(default_error,*) 'charge-charge does not have cutoffs set!'
# Line 118 | Line 151 | contains
151  
152      q1 = ChargeMap(me1)%charge
153      q2 = ChargeMap(me2)%charge
154 <
154 >    
155      if (rij.le.ecr) then
156        
157         if (rij.lt.rt) then
# Line 129 | Line 162 | contains
162            dtdr = 6.0d0*(rij*rij - rij*rt - rij*ecr +ecr*rt)/((ecr-rt)**3)
163         endif
164  
165 <       vterm = pre * q1 * q2 / rij
165 >       vcolomb = pre * q1 * q2 /rij
166 >       vrf = pre_rf * q1 * q2 * r2
167 >       vcorrection = pre_correction * q1 * q2
168 >      
169 >       vterm = vcolomb + vrf - vcorrection
170 >       !dudr = vterm * ( dtdr  - taper / rij )
171 >      
172 >       dudr =  vterm * dtdr - taper*(vcolomb/rij + 2.0d0 * vrf * taper/rij)
173  
134       !!if(rij .le. 1.5) then
135       !!  write(*,*) atom1, atom2, q1, q2, rij, vterm
136       !!endif
137
138       dudr = vterm * ( dtdr  - taper / rij )
139
174         drdx = d(1) / rij
175         drdy = d(2) / rij
176         drdz = d(3) / rij
# Line 197 | Line 231 | contains
231               tau_Temp(4) = tau_Temp(4) - d(2) * fx
232               tau_Temp(5) = tau_Temp(5) - d(2) * fy
233               tau_Temp(6) = tau_Temp(6) - d(2) * fz
200             tau_Temp(7) = tau_Temp(7) - d(3) * fx
234               tau_Temp(8) = tau_Temp(8) - d(3) * fy
235               tau_Temp(9) = tau_Temp(9) - d(3) * fz
236              

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines