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 653 by chuckv, Fri Jul 25 20:00:17 2003 UTC vs.
Revision 1169 by gezelter, Wed May 12 19:44:54 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    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(*,*) 'cutoffNotify warning: cutoffRadius is ', rcut
34 +         write(*,*) '       but switchingRadius is set larger at ', rsw
35 +         write(*,*) '       That is probably not what you wanted to do!'
36 +      endif
37  
38 +      if (rlist .lt. rcut) then
39 +         write(*,*) 'cutoffNotify warning: neighborListRadius is ', rlist
40 +         write(*,*) '       but cutoffRadius is set larger at ', rcut
41 +         write(*,*) '       That is probably a programming error!'
42 +      endif
43 +
44 +      do_shift = .false.
45 +      if (abs(rcut-rsw) .lt. 0.0001) then
46 +         write(*,*) 'cutoffNotify info: cutoffRadius and switchingRadius'
47 +         write(*,*) '       are set to the same value.  OOPSE will use'
48 +         write(*,*) '       shifted Lennard-Jones potentials instead of'
49 +         write(*,*) '       switching functions.'
50 +         do_shift = .true.
51 +      endif
52 +      
53        call setRlistDF( rlist )
54 <      call setCutoffsDipole( ecr, rtaper )
55 <      call setCutoffsRF( ecr, rtaper )
56 <      call setCutoffLJ( rcut, localError )
54 >      call setCutoffsRF( rcut, rsw )
55 >      call setCutoffLJ( rcut, do_shift, localError )
56 >      call setCutoffEAM( rcut, localError)
57 >      call set_switch(GROUP_SWITCH, rsw, rcut)
58  
59      end subroutine cutoffNotify
60  

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines