ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/notifyCutoffs.F90
Revision: 1948
Committed: Fri Jan 14 20:31:16 2005 UTC (19 years, 5 months ago) by gezelter
File size: 4187 byte(s)
Log Message:
separating modules and C/Fortran interface subroutines

File Contents

# Content
1 !!
2 !! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3 !!
4 !! The University of Notre Dame grants you ("Licensee") a
5 !! non-exclusive, royalty free, license to use, modify and
6 !! redistribute this software in source and binary code form, provided
7 !! that the following conditions are met:
8 !!
9 !! 1. Acknowledgement of the program authors must be made in any
10 !! publication of scientific results based in part on use of the
11 !! program. An acceptable form of acknowledgement is citation of
12 !! the article in which the program was described (Matthew
13 !! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 !! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 !! Parallel Simulation Engine for Molecular Dynamics,"
16 !! J. Comput. Chem. 26, pp. 252-271 (2005))
17 !!
18 !! 2. Redistributions of source code must retain the above copyright
19 !! notice, this list of conditions and the following disclaimer.
20 !!
21 !! 3. Redistributions in binary form must reproduce the above copyright
22 !! notice, this list of conditions and the following disclaimer in the
23 !! documentation and/or other materials provided with the
24 !! distribution.
25 !!
26 !! This software is provided "AS IS," without a warranty of any
27 !! kind. All express or implied conditions, representations and
28 !! warranties, including any implied warranty of merchantability,
29 !! fitness for a particular purpose or non-infringement, are hereby
30 !! excluded. The University of Notre Dame and its licensors shall not
31 !! be liable for any damages suffered by licensee as a result of
32 !! using, modifying or distributing the software or its
33 !! derivatives. In no event will the University of Notre Dame or its
34 !! licensors be liable for any lost revenue, profit or data, or for
35 !! direct, indirect, special, consequential, incidental or punitive
36 !! damages, however caused and regardless of the theory of liability,
37 !! arising out of the use of or inability to use software, even if the
38 !! University of Notre Dame has been advised of the possibility of
39 !! such damages.
40 !!
41
42 module notifyCutoffs
43
44 use definitions
45 use doForces, only: setRlistDF
46 use reaction_field, only: setCutoffsRF
47 use lj, only: setCutoffLJ
48 use eam, only: setCutoffEAM
49 use switcheroo, only: set_switch
50 use status
51 implicit none
52
53 PRIVATE
54
55 character(len = statusMsgSize) :: errMsg
56
57 #define __FORTRAN90
58 #include "UseTheForce/fSwitchingFunction.h"
59
60 public::cutoffNotify
61
62 contains
63
64 subroutine cutoffNotify( this_rcut, this_rsw, this_rlist )
65
66 real(kind=dp), intent(in) :: this_rcut, this_rsw, this_rlist
67
68 real(kind=dp) :: rsw, rcut, rlist
69 integer :: localError
70 logical :: do_shift
71
72 rcut = this_rcut
73 rsw = this_rsw
74 rlist = this_rlist
75
76 if (rcut .lt. rsw) then
77
78 write(errMsg, *) 'cutoffRadius is ', rcut, newline // tab, &
79 'but switchingRadius is set larger at ', rsw , newline // tab, &
80 'That is probably not what you wanted to do!'
81
82 call handleWarning("cutoffNotify", errMsg)
83
84 endif
85
86 if (rlist .lt. rcut) then
87
88 write(errMsg, *) 'neighborListRadius is ', rlist, newline &
89 // tab, 'but cutoffRadius is set larger at ', rcut , newline &
90 // tab, 'That is probably a programming error!'
91
92 call handleWarning("cutoffNotify", errMsg)
93
94 endif
95
96 do_shift = .false.
97 if (abs(rcut-rsw) .lt. 0.0001) then
98
99 write(errMsg, *) &
100 'cutoffRadius and switchingRadius are set to the same', newline &
101 // tab, 'value. OOPSE will use shifted Lennard-Jones', newline &
102 // tab, 'potentials instead of switching functions.'
103
104 call handleInfo("cutoffNotify", errMsg)
105
106 do_shift = .true.
107
108 endif
109
110 call setRlistDF( rlist )
111 call setCutoffsRF( rcut, rsw )
112 call setCutoffLJ( rcut, do_shift, localError )
113 call setCutoffEAM( rcut, localError)
114 call set_switch(GROUP_SWITCH, rsw, rcut)
115
116 end subroutine cutoffNotify
117
118 end module notifyCutoffs