ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/notifyCutoffs.F90
Revision: 1169
Committed: Wed May 12 19:44:54 2004 UTC (20 years, 1 month ago) by gezelter
File size: 1844 byte(s)
Log Message:
MPI fixes and removal of extraneous write statements

File Contents

# User Rev Content
1 mmeineke 626 module notifyCutoffs
2    
3     use definitions
4     use do_Forces, only: setRlistDF
5     use reaction_field, only: setCutoffsRF
6     use lj, only: setCutoffLJ
7 chuckv 653 use eam, only: setCutoffEAM
8 gezelter 1150 use switcheroo, only: set_switch
9 mmeineke 626 implicit none
10    
11     PRIVATE
12 gezelter 1150
13     #define __FORTRAN90
14     #include "fSwitchingFunction.h"
15 mmeineke 626
16     public::cutoffNotify
17    
18     contains
19    
20 gezelter 1154 subroutine cutoffNotify( this_rcut, this_rsw, this_rlist )
21 mmeineke 626
22 gezelter 1154 real(kind=dp), intent(in) :: this_rcut, this_rsw, this_rlist
23 mmeineke 626
24 gezelter 1154 real(kind=dp) :: rsw, rcut, rlist
25 mmeineke 626 integer :: localError
26 gezelter 1160 logical :: do_shift
27 gezelter 845
28 mmeineke 626 rcut = this_rcut
29 gezelter 1154 rsw = this_rsw
30 mmeineke 626 rlist = this_rlist
31    
32 gezelter 1154 if (rcut .lt. rsw) then
33 gezelter 1169 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 gezelter 845 endif
37 mmeineke 626
38 gezelter 1154 if (rlist .lt. rcut) then
39 gezelter 1169 write(*,*) 'cutoffNotify warning: neighborListRadius is ', rlist
40     write(*,*) ' but cutoffRadius is set larger at ', rcut
41     write(*,*) ' That is probably a programming error!'
42 gezelter 1154 endif
43 gezelter 1160
44     do_shift = .false.
45     if (abs(rcut-rsw) .lt. 0.0001) then
46 gezelter 1169 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 gezelter 1160 do_shift = .true.
51     endif
52 gezelter 1154
53 mmeineke 626 call setRlistDF( rlist )
54 gezelter 1154 call setCutoffsRF( rcut, rsw )
55 gezelter 1160 call setCutoffLJ( rcut, do_shift, localError )
56 gezelter 1154 call setCutoffEAM( rcut, localError)
57     call set_switch(GROUP_SWITCH, rsw, rcut)
58 mmeineke 626
59     end subroutine cutoffNotify
60    
61     end module notifyCutoffs