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 1160 by gezelter, Tue May 11 21:31:15 2004 UTC

# Line 2 | Line 2 | module notifyCutoffs
2    
3    use definitions
4    use do_Forces, only:      setRlistDF
5  use dipole_dipole, only:  setCutoffsDipole
5    use reaction_field, only: setCutoffsRF
6    use lj, only:             setCutoffLJ
7 <
7 >  use eam, only:            setCutoffEAM
8 >  use switcheroo, only:     set_switch
9    implicit none
10  
11    PRIVATE
12 +
13 + #define __FORTRAN90
14 + #include "fSwitchingFunction.h"
15    
16    public::cutoffNotify
17  
18    contains
19      
20 <    subroutine cutoffNotify( this_rcut, this_rlist, this_ecr, this_est )
20 >    subroutine cutoffNotify( this_rcut, this_rsw, this_rlist )
21        
22 <      real(kind=dp), intent(in) :: this_rcut, this_rlist, this_ecr, this_est
22 >      real(kind=dp), intent(in) :: this_rcut, this_rsw, this_rlist
23  
24 <      real(kind=dp) :: rtaper, rcut, rlist, ecr
24 >      real(kind=dp) :: rsw, rcut, rlist
25        integer :: localError
26 <      
27 <
26 >      logical :: do_shift
27 >    
28        rcut   = this_rcut
29 +      rsw    = this_rsw
30        rlist  = this_rlist
27      ecr    = this_ecr
28      rtaper = this_ecr - this_est
31  
32 +      if (rcut .lt. rsw) then
33 +         write(*,*) 'warning, rcut = ', rcut, ' but rsw = ', rsw
34 +      endif
35  
36 +      if (rlist .lt. rcut) then
37 +         write(*,*) 'warning, rlist = ', rlist, ' but rcut = ', rcut
38 +      endif
39 +
40 +      do_shift = .false.
41 +      if (abs(rcut-rsw) .lt. 0.0001) then
42 +         do_shift = .true.
43 +      endif
44 +      
45        call setRlistDF( rlist )
46 <      call setCutoffsDipole( ecr, rtaper )
47 <      call setCutoffsRF( ecr, rtaper )
48 <      call setCutoffLJ( rcut, localError )
46 >      call setCutoffsRF( rcut, rsw )
47 >      call setCutoffLJ( rcut, do_shift, localError )
48 >      call setCutoffEAM( rcut, localError)
49 >      call set_switch(GROUP_SWITCH, rsw, rcut)
50  
51      end subroutine cutoffNotify
52  

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines