ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/switcheroo.F90
Revision: 2723
Committed: Thu Apr 20 21:02:00 2006 UTC (18 years, 4 months ago) by chrisfen
File size: 6765 byte(s)
Log Message:
splined up sticky

File Contents

# User Rev Content
1 gezelter 1930 !!
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 gezelter 2722 !!$
42 gezelter 1608 module switcheroo
43    
44     use definitions
45 chrisfen 2715 use interpolation
46 gezelter 2722 use status
47 gezelter 1608
48     implicit none
49     PRIVATE
50    
51     #define __FORTRAN90
52 chrisfen 2424 #include "UseTheForce/DarkSide/fSwitchingFunctionType.h"
53 gezelter 2722
54     !! number of points for the spline approximations
55     INTEGER, PARAMETER :: np = 150
56 gezelter 1608
57 gezelter 2722 real ( kind = dp ), save :: rin
58     real ( kind = dp ), save :: rout
59     real ( kind = dp ), save :: rin2
60     real ( kind = dp ), save :: rout2
61 gezelter 1608
62 gezelter 2722 logical, save :: haveSplines = .false.
63     logical, save :: switchIsCubic = .true.
64 chrisfen 2424 integer, save :: functionType = CUBIC
65 gezelter 1608
66 gezelter 2722 ! spline structures
67     type(cubicSpline), save :: r2Spline
68     type(cubicSpline), save :: switchSpline
69 gezelter 1608
70 gezelter 2722 public::set_switch_type
71 gezelter 1608 public::set_switch
72     public::get_switch
73 gezelter 2722 public::delete_switch
74 gezelter 1608
75     contains
76    
77 gezelter 2722 subroutine set_switch(rinner, router)
78 gezelter 2717
79 gezelter 1608 real ( kind = dp ), intent(in):: rinner, router
80 gezelter 2722 real ( kind = dp ), dimension(np) :: xvals, yvals
81     real ( kind = dp ), dimension(2) :: rCubVals, sCubVals
82     real ( kind = dp ) :: rval, rval2, rval3, rval4, rval5
83     real ( kind = dp ) :: rvaldi, rvaldi2, rvaldi3, rvaldi4, rvaldi5
84     real ( kind = dp ) :: c0, c3, c4, c5, dx, r, r2
85 chrisfen 2715 integer :: i
86 gezelter 1608
87     if (router .lt. rinner) then
88 gezelter 2722 call handleError("set_switch", "router is less than rinner")
89 gezelter 1608 return
90     endif
91    
92     if ((router .lt. 0.0d0) .or. (rinner .lt. 0.0d0)) then
93 gezelter 2722 call handleError("set_switch", "one of the switches is negative!")
94 gezelter 1608 return
95     endif
96 gezelter 2204
97 gezelter 2722 rin = rinner
98     rout = router
99     rin2 = rinner * rinner
100     rout2 = router * router
101 gezelter 1608
102 gezelter 2722 if ((router-rinner) .lt. 1e-8) then
103     ! no reason to set up splines if the switching region is tiny
104     return
105 chrisfen 2715 endif
106 gezelter 2722
107     dx = (rout2-rin2) / dble(np-1)
108 chrisfen 2715
109 gezelter 2722 do i = 1, np
110     r2 = rin2 + dble(i-1)*dx
111     xvals(i) = r2
112     yvals(i) = dsqrt(r2)
113     enddo
114 gezelter 1608
115 gezelter 2722 call newSpline(r2spline, xvals, yvals, .true.)
116 chrisfen 2424
117     if (functionType .eq. FIFTH_ORDER_POLY) then
118     c0 = 1.0d0
119     c3 = -10.0d0
120     c4 = 15.0d0
121     c5 = -6.0d0
122 gezelter 2722
123     dx = (rout-rin) / dble(np-1)
124    
125     do i = 1, np
126     r = rin + dble(i-1)*dx
127     xvals(i) = r
128    
129     rval = ( r - rin )
130     rval2 = rval*rval
131     rval3 = rval2*rval
132     rval4 = rval2*rval2
133     rval5 = rval3*rval2
134     rvaldi = 1.0d0/( rout - rin )
135     rvaldi2 = rvaldi*rvaldi
136     rvaldi3 = rvaldi2*rvaldi
137     rvaldi4 = rvaldi2*rvaldi2
138     rvaldi5 = rvaldi3*rvaldi2
139     yvals(i)= c0 + c3*rval3*rvaldi3 + c4*rval4*rvaldi4 + c5*rval5*rvaldi5
140     enddo
141    
142     call newSpline(switchSpline, xvals, yvals, .true.)
143    
144     switchIsCubic = .false.
145     else
146     rCubVals(1) = rin
147     rCubVals(2) = rout
148     sCubVals(1) = 1.0d0
149     sCubVals(2) = 0.0d0
150     call newSpline(switchSpline, rCubVals, sCubVals, .true.)
151 chrisfen 2424 endif
152 gezelter 2722
153     haveSplines = .true.
154     return
155     end subroutine set_switch
156 chrisfen 2424
157 gezelter 2722 subroutine set_switch_type(functionForm)
158     integer, intent(in) :: functionForm
159     functionType = functionForm
160 gezelter 1608
161 gezelter 2722 if ((functionType.eq.FIFTH_ORDER_POLY).or.(functionType.eq.CUBIC)) then
162     if (haveSplines) then
163     call delete_switch()
164     call set_switch(rin, rout)
165     endif
166     else
167     call handleError("set_switch_type", &
168     "Unknown type of switching function!")
169     return
170     endif
171     end subroutine set_switch_type
172    
173     subroutine delete_switch()
174     call deleteSpline(switchSpline)
175     call deleteSpline(r2spline)
176     return
177     end subroutine delete_switch
178    
179     subroutine get_switch(r2, sw, dswdr, r, in_switching_region)
180    
181 gezelter 1608 real( kind = dp ), intent(in) :: r2
182     real( kind = dp ), intent(inout) :: sw, dswdr, r
183     logical, intent(inout) :: in_switching_region
184 gezelter 2717 integer :: j
185 gezelter 2722 real ( kind = dp ) :: a, b, c, d, dx
186 gezelter 1608
187 chrisfen 2723 sw = 1.0d0
188 gezelter 1608 dswdr = 0.0d0
189     in_switching_region = .false.
190    
191 chrisfen 2723 if (r2.gt.rin2) then
192     if (r2.gt.rout2) then
193 gezelter 1608
194 chrisfen 2723 sw = 0.0d0
195 gezelter 1608 dswdr = 0.0d0
196     return
197 chrisfen 2723
198 gezelter 2722 else
199    
200     call lookupUniformSpline(r2Spline, r2, r)
201     if (switchIsCubic) then
202     ! super zippy automated use of precomputed spline coefficients
203     dx = r - rin
204     sw = switchSpline%y(1) + dx*(dx*(switchSpline%c(1) + &
205     dx*switchSpline%d(1)))
206     dswdr = dx*(2.0d0 * switchSpline%c(1) + &
207     3.0d0 * dx * switchSpline%d(1))
208 chrisfen 2715 else
209 gezelter 2722 call lookupUniformSpline1d(switchSpline, r, sw, dswdr)
210 chrisfen 2715 endif
211 gezelter 2722
212 gezelter 1608 in_switching_region = .true.
213 gezelter 2722
214 gezelter 1608 return
215     endif
216     else
217     return
218 gezelter 2204 endif
219 gezelter 2722
220 gezelter 1608 end subroutine get_switch
221 chrisfen 2424
222 gezelter 1608 end module switcheroo