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 1169 by gezelter, Wed May 12 19:44:54 2004 UTC vs.
Revision 1239 by gezelter, Fri Jun 4 14:59:27 2004 UTC

# Line 6 | Line 6 | module notifyCutoffs
6    use lj, only:             setCutoffLJ
7    use eam, only:            setCutoffEAM
8    use switcheroo, only:     set_switch
9 +  use status
10    implicit none
11  
12    PRIVATE
13  
14 +  character(len = statusMsgSize) :: errMsg
15 +
16   #define __FORTRAN90
17   #include "fSwitchingFunction.h"
15  
16  public::cutoffNotify
18  
19 <  contains
19 >  public::cutoffNotify
20 >  
21 > contains
22 >  
23 >  subroutine cutoffNotify( this_rcut, this_rsw, this_rlist )
24      
25 <    subroutine cutoffNotify( this_rcut, this_rsw, this_rlist )
26 <      
27 <      real(kind=dp), intent(in) :: this_rcut, this_rsw, this_rlist
25 >    real(kind=dp), intent(in) :: this_rcut, this_rsw, this_rlist
26 >    
27 >    real(kind=dp) :: rsw, rcut, rlist
28 >    integer :: localError
29 >    logical :: do_shift
30 >    
31 >    rcut   = this_rcut
32 >    rsw    = this_rsw
33 >    rlist  = this_rlist
34 >    
35 >    if (rcut .lt. rsw) then
36 >      
37 >       write(errMsg, *) 'cutoffRadius is ', rcut, newline // tab, &            
38 >              'but switchingRadius is set larger at ', rsw , newline // tab, &
39 >              'That is probably not what you wanted to do!'
40 >        
41 >       call handleWarning("cutoffNotify", errMsg)
42  
24      real(kind=dp) :: rsw, rcut, rlist
25      integer :: localError
26      logical :: do_shift
27    
28      rcut   = this_rcut
29      rsw    = this_rsw
30      rlist  = this_rlist
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!'
43        endif
44 <
44 >      
45        if (rlist .lt. rcut) then
46 <         write(*,*) 'cutoffNotify warning: neighborListRadius is ', rlist
47 <         write(*,*) '       but cutoffRadius is set larger at ', rcut
48 <         write(*,*) '       That is probably a programming error!'
46 >        
47 >         write(errMsg, *) 'neighborListRadius is ', rlist, newline &
48 >              // tab,  'but cutoffRadius is set larger at ', rcut , newline &
49 >              // tab,  'That is probably a programming error!'
50 >        
51 >         call handleWarning("cutoffNotify", errMsg)
52 >        
53        endif
54 <
54 >      
55        do_shift = .false.
56        if (abs(rcut-rsw) .lt. 0.0001) then
57 <         write(*,*) 'cutoffNotify info: cutoffRadius and switchingRadius'
58 <         write(*,*) '       are set to the same value.  OOPSE will use'
59 <         write(*,*) '       shifted Lennard-Jones potentials instead of'
60 <         write(*,*) '       switching functions.'
57 >
58 >         write(errMsg, *) 'cutoffRadius and switchingRadius', newline &
59 >              // tab, 'are set to the same value.  OOPSE will use', newline &
60 >              // tab, 'shifted Lennard-Jones potentials instead of', newline &
61 >              // tab, 'switching functions.'
62 >
63 >         call handleInfo("cutoffNotify", errMsg)
64 >
65           do_shift = .true.
66 +
67        endif
68        
69        call setRlistDF( rlist )

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines