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 626 by mmeineke, Wed Jul 16 21:30:56 2003 UTC vs.
Revision 1154 by gezelter, Tue May 11 16:00:22 2004 UTC

# Line 2 | Line 2 | module notifyCutoffs
2    
3    use definitions
4    use do_Forces, only:      setRlistDF
5 +  use charge_charge, only:  setCutoffsCharge
6    use dipole_dipole, only:  setCutoffsDipole
7    use reaction_field, only: setCutoffsRF
8    use lj, only:             setCutoffLJ
9 <
9 >  use eam, only:            setCutoffEAM
10 >  use switcheroo, only:     set_switch
11    implicit none
12  
13    PRIVATE
14 +
15 + #define __FORTRAN90
16 + #include "fSwitchingFunction.h"
17    
18    public::cutoffNotify
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 <      
24 <
28 >    
29        rcut   = this_rcut
30 +      rsw    = this_rsw
31        rlist  = this_rlist
27      ecr    = this_ecr
28      rtaper = this_ecr - this_est
32  
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 setCutoffsDipole( ecr, rtaper )
43 <      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)
47 +      call set_switch(GROUP_SWITCH, rsw, rcut)
48  
49      end subroutine cutoffNotify
50  

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines