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 1233 by gezelter, Fri Jun 4 02:38:23 2004 UTC

# Line 12 | Line 12 | module notifyCutoffs
12  
13   #define __FORTRAN90
14   #include "fSwitchingFunction.h"
15 + #include "simError.h"
16    
17    public::cutoffNotify
18  
# Line 30 | Line 31 | module notifyCutoffs
31        rlist  = this_rlist
32  
33        if (rcut .lt. rsw) then
34 <         write(*,*) 'cutoffNotify warning: cutoffRadius is ', rcut
35 <         write(*,*) '       but switchingRadius is set larger at ', rsw
36 <         write(*,*) '       That is probably not what you wanted to do!'
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        if (rlist .lt. rcut) then
49 <         write(*,*) 'cutoffNotify warning: neighborListRadius is ', rlist
50 <         write(*,*) '       but cutoffRadius is set larger at ', rcut
51 <         write(*,*) '       That is probably a programming error!'
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 >
60        endif
61  
62        do_shift = .false.
63        if (abs(rcut-rsw) .lt. 0.0001) then
64 <         write(*,*) 'cutoffNotify info: cutoffRadius and switchingRadius'
65 <         write(*,*) '       are set to the same value.  OOPSE will use'
66 <         write(*,*) '       shifted Lennard-Jones potentials instead of'
67 <         write(*,*) '       switching functions.'
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 )

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines