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 1233 by gezelter, Fri Jun 4 02:38:23 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
5    use reaction_field, only: setCutoffsRF
6    use lj, only:             setCutoffLJ
7    use eam, only:            setCutoffEAM
# Line 14 | Line 12 | module notifyCutoffs
12  
13   #define __FORTRAN90
14   #include "fSwitchingFunction.h"
15 + #include "simError.h"
16    
17    public::cutoffNotify
18  
19    contains
20      
21 <    subroutine cutoffNotify( this_rcut, this_rlist, this_ecr, this_est )
21 >    subroutine cutoffNotify( this_rcut, this_rsw, this_rlist )
22        
23 <      real(kind=dp), intent(in) :: this_rcut, this_rlist, this_ecr, this_est
23 >      real(kind=dp), intent(in) :: this_rcut, this_rsw, this_rlist
24  
25 <      real(kind=dp) :: rtaper, rcut, rlist, ecr
25 >      real(kind=dp) :: rsw, rcut, rlist
26        integer :: localError
27 +      logical :: do_shift
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 >
35 >         write(painCave%errMsg, *) 'cutoffRadius is ', rcut, &
36 >              achar(10) // achar(9), &
37 >              'but switchingRadius is set larger at ', rsw , &
38 >              achar(10) // achar(9) , &
39 >              'That is probably not what you wanted to do!', &
40 >              achar(10) // achar(0)
41 >        
42 >         painCave%severity = OOPSE_WARNING
43 >         painCave%isFatal = .false.
44 >         call c_simError(painCave)
45 >
46        endif
47  
48 <      call setRlistDF( rlist )
49 <      call setCutoffsCharge( ecr, rtaper )
50 <      call setCutoffsDipole( ecr, rtaper )
51 <      call setCutoffsRF( ecr, rtaper )
52 <      call setCutoffLJ( rcut, localError )
53 <      call setCutoffEAM(rcut, localError)
48 >      if (rlist .lt. rcut) then
49 >        
50 >         write(painCave%errMsg, *) 'neighborListRadius is ', rlist, &
51 >              achar(10) // achar(9), &
52 >              'but cutoffRadius is set larger at ', rcut , &
53 >              achar(10) // achar(9) , &
54 >              'That is probably a programming error!', &
55 >              achar(10) // achar(0)
56 >         painCave%severity = OOPSE_WARNING
57 >         painCave%isFatal = .false.
58 >         call c_simError(painCave)
59  
45      if (ecr.gt.rcut) then
46         call set_switch(GROUP_SWITCH, rtaper, ecr)
47      else
48         call set_switch(GROUP_SWITCH, rcut, rcut)
60        endif
61  
62 +      do_shift = .false.
63 +      if (abs(rcut-rsw) .lt. 0.0001) then
64 +
65 +         write(painCave%errMsg,*) 'cutoffRadius and switchingRadius', &
66 +              achar(10) // achar(9), &
67 +              'are set to the same value.  OOPSE will use', &
68 +              achar(10) // achar(9), &
69 +              'shifted Lennard-Jones potentials instead of', &
70 +              achar(10) // achar(9), &
71 +              'switching functions.', &
72 +              achar(10) // achar(0)
73 +         painCave%severity = OOPSE_INFO
74 +         painCave%isFatal = .false.
75 +         call c_simError(painCave)
76 +
77 +         do_shift = .true.
78 +
79 +      endif
80 +      
81 +      call setRlistDF( rlist )
82 +      call setCutoffsRF( rcut, rsw )
83 +      call setCutoffLJ( rcut, do_shift, localError )
84 +      call setCutoffEAM( rcut, localError)
85 +      call set_switch(GROUP_SWITCH, rsw, rcut)
86 +
87      end subroutine cutoffNotify
88  
89    end module notifyCutoffs

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines