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