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 switcheroo |
43 |
|
44 |
use definitions |
45 |
use interpolation |
46 |
|
47 |
implicit none |
48 |
PRIVATE |
49 |
|
50 |
#define __FORTRAN90 |
51 |
#include "UseTheForce/fSwitchingFunction.h" |
52 |
#include "UseTheForce/DarkSide/fSwitchingFunctionType.h" |
53 |
|
54 |
real ( kind = dp ), dimension(NSWITCHTYPES) :: rin |
55 |
real ( kind = dp ), dimension(NSWITCHTYPES) :: rout |
56 |
real ( kind = dp ), dimension(NSWITCHTYPES) :: rin2 |
57 |
real ( kind = dp ), dimension(NSWITCHTYPES) :: rout2 |
58 |
real ( kind = dp ), save :: c0, c1, c2, c3, c4, c5 |
59 |
|
60 |
logical, dimension(NSWITCHTYPES) :: isOK |
61 |
logical, save :: haveFunctionType = .false. |
62 |
logical, save :: haveSqrtSpline = .false. |
63 |
logical, save :: useSpline = .true. |
64 |
integer, save :: functionType = CUBIC |
65 |
|
66 |
! spline structure |
67 |
type(cubicSpline), save :: splineSqrt |
68 |
|
69 |
public::set_switch |
70 |
public::set_function_type |
71 |
public::get_switch |
72 |
|
73 |
contains |
74 |
|
75 |
subroutine setupSqrtSpline(rmin, rmax, np) |
76 |
real( kind = dp ), intent(in) :: rmin, rmax |
77 |
integer, intent(in) :: np |
78 |
real( kind = dp ) :: rvals(np), yvals(np) |
79 |
real( kind = dp ) :: dr, r |
80 |
real( kind = dp ) :: yprime1, yprimen |
81 |
integer :: i |
82 |
|
83 |
dr = (rmax-rmin) / float(np-1) |
84 |
|
85 |
do i = 1, np |
86 |
r = rmin + float(i-1)*dr |
87 |
rvals(i) = r |
88 |
yvals(i) = dsqrt(r) |
89 |
enddo |
90 |
|
91 |
yprime1 = 1.0d0 / ( 2.0d0 * dsqrt( rmin ) ) |
92 |
yprimeN = 1.0d0 / ( 2.0d0 * dsqrt( rmax ) ) |
93 |
|
94 |
call newSpline(splineSqrt, rvals, yvals, yprime1, yprimen, .true.) |
95 |
|
96 |
return |
97 |
end subroutine setupSqrtSpline |
98 |
|
99 |
subroutine set_switch(SwitchType, rinner, router) |
100 |
|
101 |
real ( kind = dp ), intent(in):: rinner, router |
102 |
integer, intent(in) :: SwitchType |
103 |
integer :: i |
104 |
|
105 |
if (SwitchType .gt. NSWITCHTYPES) then |
106 |
write(default_error, *) & |
107 |
'set_switch: not that many switch types! ' |
108 |
return |
109 |
endif |
110 |
|
111 |
isOK(SwitchType) = .false. |
112 |
|
113 |
if (router .lt. rinner) then |
114 |
write(default_error, *) & |
115 |
'set_switch: router is less than rinner ' |
116 |
return |
117 |
endif |
118 |
|
119 |
if ((router .lt. 0.0d0) .or. (rinner .lt. 0.0d0)) then |
120 |
write(default_error, *) & |
121 |
'set_switch: one of the switches is negative!' |
122 |
return |
123 |
endif |
124 |
|
125 |
rin(SwitchType) = rinner |
126 |
rout(SwitchType) = router |
127 |
rin2(SwitchType) = rinner * rinner |
128 |
rout2(SwitchType) = router * router |
129 |
isOK(SwitchType) = .true. |
130 |
|
131 |
if (.not.haveSqrtSpline) then |
132 |
! fill arrays for building the spline |
133 |
call setupSqrtSpline(1.0d0, rout2(switchType), SPLINE_SEGMENTS) |
134 |
haveSqrtSpline = .true. |
135 |
endif |
136 |
|
137 |
end subroutine set_switch |
138 |
|
139 |
subroutine set_function_type(functionForm) |
140 |
integer, intent(in) :: functionForm |
141 |
functionType = functionForm |
142 |
|
143 |
if (functionType .eq. FIFTH_ORDER_POLY) then |
144 |
c0 = 1.0d0 |
145 |
c1 = 0.0d0 |
146 |
c2 = 0.0d0 |
147 |
c3 = -10.0d0 |
148 |
c4 = 15.0d0 |
149 |
c5 = -6.0d0 |
150 |
endif |
151 |
end subroutine set_function_type |
152 |
|
153 |
subroutine get_switch(r2, sw, dswdr, r, SwitchType, in_switching_region) |
154 |
|
155 |
real( kind = dp ), intent(in) :: r2 |
156 |
real( kind = dp ), intent(inout) :: sw, dswdr, r |
157 |
real( kind = dp ) :: ron, roff, a, b, c, d, dx |
158 |
real( kind = dp ) :: rval, rval2, rval3, rval4, rval5 |
159 |
real( kind = dp ) :: rvaldi, rvaldi2, rvaldi3, rvaldi4, rvaldi5 |
160 |
integer, intent(in) :: SwitchType |
161 |
logical, intent(inout) :: in_switching_region |
162 |
integer :: j |
163 |
|
164 |
sw = 0.0d0 |
165 |
dswdr = 0.0d0 |
166 |
in_switching_region = .false. |
167 |
|
168 |
if (.not.isOK(SwitchType)) then |
169 |
write(default_error, *) & |
170 |
'get_switch: this switching function has not been set up!' |
171 |
return |
172 |
endif |
173 |
|
174 |
if (r2.lt.rout2(SwitchType)) then |
175 |
if (r2.lt.rin2(SwitchType)) then |
176 |
|
177 |
sw = 1.0d0 |
178 |
dswdr = 0.0d0 |
179 |
return |
180 |
|
181 |
else |
182 |
if (useSpline) then |
183 |
j = MAX(1, MIN(splineSqrt%np, idint((r2-splineSqrt%x(1)) * splineSqrt%dx_i) + 1)) |
184 |
|
185 |
dx = r2 - splineSqrt%x(j) |
186 |
|
187 |
a = splineSqrt%c(1,j) |
188 |
b = splineSqrt%c(2,j) |
189 |
c = splineSqrt%c(3,j) |
190 |
d = splineSqrt%c(4,j) |
191 |
|
192 |
r = c + dx * d |
193 |
r = b + dx * r |
194 |
r = a + dx * r |
195 |
else |
196 |
r = dsqrt(r2) |
197 |
endif |
198 |
|
199 |
ron = rin(SwitchType) |
200 |
roff = rout(SwitchType) |
201 |
|
202 |
if (functionType .eq. FIFTH_ORDER_POLY) then |
203 |
rval = ( r - ron ) |
204 |
rval2 = rval*rval |
205 |
rval3 = rval2*rval |
206 |
rval4 = rval2*rval2 |
207 |
rval5 = rval3*rval2 |
208 |
rvaldi = 1.0d0/( roff - ron ) |
209 |
rvaldi2 = rvaldi*rvaldi |
210 |
rvaldi3 = rvaldi2*rvaldi |
211 |
rvaldi4 = rvaldi2*rvaldi2 |
212 |
rvaldi5 = rvaldi3*rvaldi2 |
213 |
sw = c0 + c1*rval*rvaldi + c2*rval2*rvaldi2 + c3*rval3*rvaldi3 & |
214 |
+ c4*rval4*rvaldi4 + c5*rval5*rvaldi5 |
215 |
dswdr = c1*rvaldi + 2.0d0*c2*rval*rvaldi2 & |
216 |
+ 3.0d0*c3*rval2*rvaldi3 + 4.0d0*c4*rval3*rvaldi4 & |
217 |
+ 5.0d0*c5*rval4*rvaldi5 |
218 |
|
219 |
else |
220 |
sw = (roff + 2.0d0*r - 3.0d0*ron)*(roff-r)**2/ ((roff-ron)**3) |
221 |
dswdr = 6.0d0*(r*r - r*ron - r*roff +roff*ron)/((roff-ron)**3) |
222 |
|
223 |
endif |
224 |
in_switching_region = .true. |
225 |
return |
226 |
endif |
227 |
else |
228 |
return |
229 |
endif |
230 |
|
231 |
end subroutine get_switch |
232 |
|
233 |
end module switcheroo |