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 944 by gezelter, Tue Jan 13 22:34:55 2004 UTC vs.
Revision 945 by gezelter, Thu Jan 15 01:14:55 2004 UTC

# Line 12 | Line 12 | module charge_charge
12    implicit none
13  
14    PRIVATE
15 <  real(kind=dp), save :: rcut = 0.0_DP
15 >  real(kind=dp), save :: ecr = 0.0_DP
16    real(kind=dp), save :: rt  = 0.0_DP
17    real(kind=dp), save :: pre = 0.0_DP
18    logical, save :: haveCutoffs = .false.
19    logical, save :: haveChargeMap = .false.
20  
21 <  public::setCutoffs
21 >  public::setCutoffsCharge
22    public::do_charge_pair
23  
24    type :: ChargeList
# Line 29 | Line 29 | contains
29  
30   contains
31      
32 <  subroutine setCutoffs(this_rcut, this_rt)
33 <    real(kind=dp), intent(in) :: this_rcut, this_rt
34 <    rcut = this_rcut
32 >  subroutine setCutoffsCharge(this_ecr, this_rt)
33 >    real(kind=dp), intent(in) :: this_ecr, this_rt
34 >    ecr = this_ecr
35      rt = this_rt
36  
37      ! pre converts from fundamental charge to kcal/mol
# Line 40 | Line 40 | contains
40      haveCutoffs = .true.
41      
42      return
43 <  end subroutine setCutoffs
43 >  end subroutine setCutoffsCharge
44  
45    subroutine createChargeMap(status)
46      integer :: nAtypes
# Line 118 | Line 118 | contains
118      q1 = ChargeMap(me1)%charge
119      q2 = ChargeMap(me2)%charge
120  
121 <    if (rij.le.rcut) then
121 >    if (rij.le.ecr) then
122        
123         if (rij.lt.rt) then
124            taper = 1.0d0
125            dtdr = 0.0d0
126         else
127 <          taper = (rcut + 2.0d0*rij - 3.0d0*rt)*(rcut-rij)**2/ ((rcut-rt)**3)
128 <          dtdr = 6.0d0*(rij*rij - rij*rt - rij*rcut +rcut*rt)/((rcut-rt)**3)
127 >          taper = (ecr + 2.0d0*rij - 3.0d0*rt)*(ecr-rij)**2/ ((ecr-rt)**3)
128 >          dtdr = 6.0d0*(rij*rij - rij*rt - rij*ecr +ecr*rt)/((ecr-rt)**3)
129         endif
130  
131         vterm = pre * q1 * q2 / rij

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines