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

Comparing trunk/OOPSE/libmdtools/notifyCutoffs.F90 (file contents):
Revision 1150 by gezelter, Fri May 7 21:35:05 2004 UTC vs.
Revision 1154 by gezelter, Tue May 11 16:00:22 2004 UTC

# Line 19 | Line 19 | module notifyCutoffs
19  
20    contains
21      
22 <    subroutine cutoffNotify( this_rcut, this_rlist, this_ecr, this_est )
22 >    subroutine cutoffNotify( this_rcut, this_rsw, this_rlist )
23        
24 <      real(kind=dp), intent(in) :: this_rcut, this_rlist, this_ecr, this_est
24 >      real(kind=dp), intent(in) :: this_rcut, this_rsw, this_rlist
25  
26 <      real(kind=dp) :: rtaper, rcut, rlist, ecr
26 >      real(kind=dp) :: rsw, rcut, rlist
27        integer :: localError
28      
29        rcut   = this_rcut
30 +      rsw    = this_rsw
31        rlist  = this_rlist
31      ecr    = this_ecr
32      rtaper = this_ecr - this_est
32  
33 <      if ((rlist .lt. rcut) .or. (rlist .lt. ecr)) then
34 <         write(*,*) 'warning, rlist = ', rlist, ' but rcut, ecr = ', rcut, ecr
33 >      if (rcut .lt. rsw) then
34 >         write(*,*) 'warning, rcut = ', rcut, ' but rsw = ', rsw
35        endif
36  
37 +      if (rlist .lt. rcut) then
38 +         write(*,*) 'warning, rlist = ', rlist, ' but rcut = ', rcut
39 +      endif
40 +      
41        call setRlistDF( rlist )
42 <      call setCutoffsCharge( ecr, rtaper )
43 <      call setCutoffsDipole( ecr, rtaper )
44 <      call setCutoffsRF( ecr, rtaper )
42 >      call setCutoffsCharge( rcut, rsw )
43 >      call setCutoffsDipole( rcut, rsw )
44 >      call setCutoffsRF( rcut, rsw )
45        call setCutoffLJ( rcut, localError )
46 <      call setCutoffEAM(rcut, localError)
46 >      call setCutoffEAM( rcut, localError)
47 >      call set_switch(GROUP_SWITCH, rsw, rcut)
48  
45      if (ecr.gt.rcut) then
46         call set_switch(GROUP_SWITCH, rtaper, ecr)
47      else
48         call set_switch(GROUP_SWITCH, rcut, rcut)
49      endif
50
49      end subroutine cutoffNotify
50  
51    end module notifyCutoffs

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines