-
Notifications
You must be signed in to change notification settings - Fork 70
/
functions.jl
393 lines (300 loc) · 10.5 KB
/
functions.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
# This file is part of the IntervalArithmetic.jl package; MIT licensed
## Powers
# CRlibm does not contain a correctly-rounded ^ function for Float64
# Use the BigFloat version from MPFR instead, which is correctly-rounded:
# Write explicitly like this to avoid ambiguity warnings:
for T in (:Integer, :Float64, :BigFloat, :Interval)
@eval ^(a::Interval{Float64}, x::$T) = atomic(Interval{Float64}, bigequiv(a)^x)
end
^(a::Interval{Float32}, x::Interval) = atomic(Interval{Float32}, bigequiv(a)^x)
# Integer power:
# overwrite new behaviour for small integer powers from
# https://github.com/JuliaLang/julia/pull/24240:
Base.literal_pow(::typeof(^), x::Interval{T}, ::Val{p}) where {T,p} = x^p
Base.eltype(x::Interval{T}) where {T<:Real} = Interval{T}
"""
numtype(::Interval{T}) where {T<:Real} = T
Returns the type of the bounds of the interval.
### Example
```julia
julia> numtype(1..2)
Float64
```
"""
numtype(::Interval{T}) where {T<:Real} = T
function ^(a::Interval{BigFloat}, n::Integer)
isempty(a) && return a
n == 0 && return one(a)
n == 1 && return a
# n == 2 && return sqr(a)
n < 0 && a == zero(a) && return emptyinterval(a)
T = BigFloat
if isodd(n) # odd power
isentire(a) && return a
if n > 0
a.lo == 0 && return @round(0, a.hi^n)
a.hi == 0 && return @round(a.lo^n, 0)
return @round(a.lo^n, a.hi^n)
else
if a.lo ≥ 0
a.lo == 0 && return @round(a.hi^n, Inf)
return @round(a.hi^n, a.lo^n)
elseif a.hi ≤ 0
a.hi == 0 && return @round(-Inf, a.lo^n)
return @round(a.hi^n, a.lo^n)
else
return entireinterval(a)
end
end
else # even power
if n > 0
if a.lo ≥ 0
return @round(a.lo^n, a.hi^n)
elseif a.hi ≤ 0
return @round(a.hi^n, a.lo^n)
else
return @round(mig(a)^n, mag(a)^n)
end
else
if a.lo ≥ 0
return @round(a.hi^n, a.lo^n)
elseif a.hi ≤ 0
return @round(a.lo^n, a.hi^n)
else
return @round(mag(a)^n, mig(a)^n)
end
end
end
end
function sqr(a::Interval{T}) where T<:Real
return a^2
# isempty(a) && return a
# if a.lo ≥ zero(T)
# return @round(a.lo^2, a.hi^2)
#
# elseif a.hi ≤ zero(T)
# return @round(a.hi^2, a.lo^2)
# end
#
# return @round(mig(a)^2, mag(a)^2)
end
^(a::Interval{BigFloat}, x::AbstractFloat) = a^big(x)
# Floating-point power of a BigFloat interval:
function ^(a::Interval{BigFloat}, x::BigFloat)
domain = Interval{BigFloat}(0, Inf)
if a == zero(a)
a = a ∩ domain
x > zero(x) && return zero(a)
return emptyinterval(a)
end
isinteger(x) && return a^(round(Int, x))
x == 0.5 && return sqrt(a)
a = a ∩ domain
(isempty(x) || isempty(a)) && return emptyinterval(a)
xx = atomic(Interval{BigFloat}, x)
# @round() can't be used directly, because both arguments may
# Inf or -Inf, which throws an error
# lo = @round(a.lo^xx.lo, a.lo^xx.lo)
lolod = @round_down(a.lo^xx.lo)
lolou = @round_up(a.lo^xx.lo)
lo = (lolod == Inf || lolou == -Inf) ?
wideinterval(lolod) : Interval(lolod, lolou)
# lo1 = @round(a.lo^xx.hi, a.lo^xx.hi)
lohid = @round_down(a.lo^xx.hi)
lohiu = @round_up(a.lo^xx.hi)
lo1 = (lohid == Inf || lohiu == -Inf) ?
wideinterval(lohid) : Interval(lohid, lohiu)
# hi = @round(a.hi^xx.lo, a.hi^xx.lo)
hilod = @round_down(a.hi^xx.lo)
hilou = @round_up(a.hi^xx.lo)
hi = (hilod == Inf || hilou == -Inf) ?
wideinterval(hilod) : Interval(hilod, hilou)
# hi1 = @round(a.hi^xx.hi, a.hi^xx.hi)
hihid = @round_down(a.hi^xx.hi)
hihiu = @round_up(a.hi^xx.hi)
hi1 = (hihid == Inf || hihiu == -Inf) ?
wideinterval(hihid) : Interval(hihid, hihiu)
lo = hull(lo, lo1)
hi = hull(hi, hi1)
return hull(lo, hi)
end
function ^(a::Interval{Rational{T}}, x::AbstractFloat) where T<:Integer
a = Interval(a.lo.num/a.lo.den, a.hi.num/a.hi.den)
a = a^x
atomic(Interval{Rational{T}}, a)
end
# Rational power
function ^(a::Interval{T}, x::Rational) where T
domain = Interval{T}(0, Inf)
p = x.num
q = x.den
isempty(a) && return emptyinterval(a)
x == 0 && return one(a)
if a == zero(a)
x > zero(x) && return zero(a)
return emptyinterval(a)
end
x == (1//2) && return sqrt(a)
if x >= 0
if a.lo ≥ 0
isinteger(x) && return a ^ Int64(x)
a = @biginterval(a)
ui = convert(Culong, q)
low = BigFloat()
high = BigFloat()
ccall((:mpfr_rootn_ui, :libmpfr) , Int32 , (Ref{BigFloat}, Ref{BigFloat}, Culong, MPFRRoundingMode) , low , a.lo , ui, MPFRRoundDown)
ccall((:mpfr_rootn_ui, :libmpfr) , Int32 , (Ref{BigFloat}, Ref{BigFloat}, Culong, MPFRRoundingMode) , high , a.hi , ui, MPFRRoundUp)
b = interval(low, high)
b = convert(Interval{T}, b)
return b^p
end
if a.lo < 0 && a.hi ≥ 0
isinteger(x) && return a ^ Int64(x)
a = a ∩ Interval{T}(0, Inf)
a = @biginterval(a)
ui = convert(Culong, q)
low = BigFloat()
high = BigFloat()
ccall((:mpfr_rootn_ui, :libmpfr) , Int32 , (Ref{BigFloat}, Ref{BigFloat}, Culong, MPFRRoundingMode) , low , a.lo , ui, MPFRRoundDown)
ccall((:mpfr_rootn_ui, :libmpfr) , Int32 , (Ref{BigFloat}, Ref{BigFloat}, Culong, MPFRRoundingMode) , high , a.hi , ui, MPFRRoundUp)
b = interval(low, high)
b = convert(Interval{T}, b)
return b^p
end
if a.hi < 0
isinteger(x) && return a ^ Int64(x)
return emptyinterval(a)
end
end
if x < 0
return inv(a^(-x))
end
end
# Interval power of an interval:
function ^(a::Interval{BigFloat}, x::Interval)
T = BigFloat
domain = Interval{T}(0, Inf)
a = a ∩ domain
(isempty(x) || isempty(a)) && return emptyinterval(a)
lo1 = a^x.lo
lo2 = a^x.hi
lo1 = hull(lo1, lo2)
hi1 = a^x.lo
hi2 = a^x.hi
hi1 = hull(hi1, hi2)
hull(lo1, hi1)
end
function sqrt(a::Interval{T}) where T
domain = Interval{T}(0, Inf)
a = a ∩ domain
isempty(a) && return a
@round(sqrt(a.lo), sqrt(a.hi)) # `sqrt` is correctly-rounded
end
"""
pow(x::Interval, n::Integer)
A faster implementation of `x^n`, currently using `power_by_squaring`.
`pow(x, n)` will usually return an interval that is slightly larger than that calculated by `x^n`, but is guaranteed to be a correct
enclosure when using multiplication with correct rounding.
"""
function pow(x::Interval, n::Integer) # fast integer power
if n < 0
return 1/pow(x, -n)
end
isempty(x) && return x
if iseven(n) && 0 ∈ x
return hull(zero(x),
hull(Base.power_by_squaring(Interval(mig(x)), n),
Base.power_by_squaring(Interval(mag(x)), n))
)
else
return hull( Base.power_by_squaring(Interval(x.lo), n),
Base.power_by_squaring(Interval(x.hi), n) )
end
end
function pow(x::Interval, y::Real) # fast real power, including for y an Interval
isempty(x) && return x
isinteger(y) && return pow(x, Int(y.lo))
return exp(y * log(x))
end
for f in (:exp, :expm1)
@eval begin
function ($f)(a::Interval{T}) where T
isempty(a) && return a
@round( ($f)(a.lo), ($f)(a.hi) )
end
end
end
for f in (:exp2, :exp10, :cbrt)
@eval function ($f)(x::BigFloat, r::RoundingMode) # add BigFloat functions with rounding:
setrounding(BigFloat, r) do
($f)(x)
end
end
@eval ($f)(a::Interval{T}) where T = atomic(Interval{T}, $f(bigequiv(a))) # no CRlibm version
@eval function ($f)(a::Interval{BigFloat})
isempty(a) && return a
@round( ($f)(a.lo), ($f)(a.hi) )
end
end
for f in (:log, :log2, :log10)
@eval function ($f)(a::Interval{T}) where T
domain = Interval{T}(0, Inf)
a = a ∩ domain
(isempty(a) || a.hi ≤ zero(T)) && return emptyinterval(a)
@round( ($f)(a.lo), ($f)(a.hi) )
end
end
function log1p(a::Interval{T}) where T
domain = Interval{T}(-1, Inf)
a = a ∩ domain
(isempty(a) || a.hi ≤ -one(T)) && return emptyinterval(a)
@round( log1p(a.lo), log1p(a.hi) )
end
hypot(a::Interval{BigFloat}, b::Interval{BigFloat}) = sqrt(a^2 + b^2)
hypot(a::Interval{T}, b::Interval{T}) where T= atomic(Interval{T}, hypot(bigequiv(a), bigequiv(b)))
"""
nthroot(a::Interval{BigFloat}, n::Integer)
Compute the real n-th root of Interval.
"""
function nthroot(a::Interval{BigFloat}, n::Integer)
n == 1 && return a
n < 0 && a == zero(a) && return emptyinterval(a)
isempty(a) && return a
if n > 0
a.hi < 0 && iseven(n) && return emptyinterval(BigFloat)
if a.lo < 0 && a.hi >= 0 && iseven(n)
a = a ∩ Interval{BigFloat}(0, Inf)
end
ui = convert(Culong, n)
low = BigFloat()
high = BigFloat()
ccall((:mpfr_rootn_ui, :libmpfr), Int32 , (Ref{BigFloat}, Ref{BigFloat}, Culong, MPFRRoundingMode) , low , a.lo , ui, MPFRRoundDown)
ccall((:mpfr_rootn_ui, :libmpfr), Int32 , (Ref{BigFloat}, Ref{BigFloat}, Culong, MPFRRoundingMode) , high , a.hi , ui, MPFRRoundUp)
b = interval(low , high)
return b
elseif n < 0
return inv(nthroot(a, -n))
elseif n == 0
return emptyinterval(a)
end
end
function nthroot(a::Interval{T}, n::Integer) where T
b = nthroot(bigequiv(a), n)
return convert(Interval{T}, b)
end
"""
Calculate `x::Interval mod y::Real` where y != zero(y) and y is not inteval`.
"""
function mod(x::Interval, y::Real)
@assert y != zero(y) """mod(x::Interval, y::Real)
is currently implemented only for a strictly positive or negative divisor y."""
division = x / y
fl = floor(division)
if !isthin(fl)
return y > zero(y) ? Interval(zero(y), y) : Interval(y, zero(y))
else
return y * (division - fl)
end
end
mod(x::Interval, y::Interval) = throw(ArgumentError("mod not defined for interval as divisor `y`"))
mod(x::Real, y::Interval) = throw(ArgumentError("mod not defined for interval as divisor `y`"))