Skip to content

Commit 0b49126

Browse files
committed
Added logabspfaffian(::AbstractMatrix{<:Complex})
I realised that the algorithm used to compute `pfaffian(::AbstractMatrix{<:Complex})` can be used to compute the `logabspfaffian` with little modifications and decided to add that too for completeness sake
1 parent 3c93770 commit 0b49126

File tree

3 files changed

+60
-6
lines changed

3 files changed

+60
-6
lines changed

docs/src/pfaffian.md

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -49,7 +49,7 @@ julia> pfaffian!(BigInt[0 2 -7; -2 0 -8; 7 8 0])
4949
```
5050
(Note that the Pfaffian is *always zero* for any *odd* size skew-symmetric matrix.)
5151

52-
Since the computation of the pfaffian can easily overflow/underflow the maximum/minimum representable floating-point value, we also provide a function `logabspfaffian` for real skew-symmetric matrices (along with an in-place variant `logabspfaffian!`) that returns a tuple `(logpf, sign)` such
52+
Since the computation of the pfaffian can easily overflow/underflow the maximum/minimum representable floating-point value, we also provide a function `logabspfaffian` (along with an in-place variant `logabspfaffian!`) that returns a tuple `(logpf, sign)` such
5353
that the Pfaffian is `sign * exp(logpf)`. (This is similar to the [`logabsdet`](https://docs.julialang.org/en/v1/stdlib/LinearAlgebra/#LinearAlgebra.logabsdet) function in Julia's `LinearAlgebra` library to compute the log of the determinant.)
5454

5555
```jl

src/pfaffian.jl

Lines changed: 56 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -171,6 +171,7 @@ function _logabspfaffian!(A::SkewHermitian{<:Real})
171171
end
172172
return logpf, sgn*sign(det(H.Q))
173173
end
174+
174175
function _logabspfaffian!(A::SkewHermTridiagonal{<:Real})
175176
n = size(A, 1)
176177
isodd(n) && return convert(eltype(A.ev), -Inf), zero(eltype(A.ev))
@@ -182,18 +183,62 @@ function _logabspfaffian!(A::SkewHermTridiagonal{<:Real})
182183
end
183184
return logpf, sgn
184185
end
186+
187+
function _logabspfaffian!(A::AbstractMatrix{<:Complex})
188+
n = size(A, 1)
189+
isodd(n) && return convert(real(eltype(A)), -Inf), zero(eltype(A))
190+
tau = Array{eltype(A)}(undef, n - 2)
191+
logpf = zero(real(eltype(A)))
192+
phase = zero(real(eltype(A)))
193+
@inbounds for k in 1:2:n-1
194+
tauk = @view tau[k:end]
195+
196+
# Pivot if neccessary
197+
@views kp = k + argmax(Iterators.map(abs, A[k+1:end, k]))
198+
if kp != k + 1
199+
@inbounds @simd for l in k:n
200+
A[k+1,l], A[kp,l] = A[kp,l], A[k+1,l]
201+
end
202+
@inbounds @simd for l in k:n
203+
A[l,k+1], A[l,kp] = A[l,kp], A[l,k+1]
204+
end
205+
phase += π
206+
end
207+
208+
# Apply Gauss transformation and update `pf`
209+
@inbounds if A[k+1,k] != zero(eltype(A))
210+
@inbounds @views tauk .= A[k,k+2:end] ./ A[k,k+1]
211+
logpf += @inbounds log(abs(A[k,k+1]))
212+
phase += @inbounds angle(A[k,k+1])
213+
if k + 2 <= n
214+
@inbounds for l1 in eachindex(tauk)
215+
@simd for l2 in eachindex(tauk)
216+
@fastmath A[k+1+l2,k+1+l1] +=
217+
tauk[l2] * A[k+1+l1,k+1] -
218+
tauk[l1] * A[k+1+l2,k+1]
219+
end
220+
end
221+
end
222+
else
223+
return convert(real(eltype(A)), -Inf), zero(eltype(A))
224+
end
225+
end
226+
227+
return logpf, cis(phase)
228+
end
229+
185230
logabspfaffian!(A::Union{SkewHermitian{<:Real},SkewHermTridiagonal{<:Real}})= _logabspfaffian!(A)
186231
logabspfaffian(A::Union{SkewHermitian{<:Real},SkewHermTridiagonal{<:Real}})= logabspfaffian!(copyeigtype(A))
187232

188233
"""
189234
logabspfaffian(A)
190235
191-
Returns a tuple `(log|Pf A|, sign)`, with the log of the absolute value of the pfaffian of `A` as first output
192-
and the sign (±1) of the pfaffian as second output. A must be a real skew-Hermitian matrix.
193-
If `A` is not of type `SkewHermitian{<:Real}`, then `isskewhermitian(A)`
194-
is checked to ensure that `A == -A'`
236+
Returns a tuple `(log|pf(A)|, sign)`, with the log of the absolute value of the pfaffian of
237+
`A` as first output and the sign of the pfaffian as second output such that
238+
`pfaffian(A) ≈ sign * exp(log|pf(A)|)`. `A` must be a skew-symmetric matrix. If `A` is not of type
239+
`SkewHermitian{<:Real}`, then `isskewsymmetric(A)` is checked to ensure that `A == -transpose(A)`
195240
"""
196-
logabspfaffian(A::AbstractMatrix{<:Real}) = logabspfaffian!(copyeigtype(A))
241+
logabspfaffian(A::AbstractMatrix) = logabspfaffian!(copyeigtype(A))
197242

198243

199244
"""
@@ -205,3 +250,9 @@ function logabspfaffian!(A::AbstractMatrix{<:Real})
205250
isskewhermitian(A) || throw(ArgumentError("Pfaffian requires a skew-Hermitian matrix"))
206251
return _logabspfaffian!(SkewHermitian(A))
207252
end
253+
254+
function logabspfaffian!(A::AbstractMatrix{<:Complex})
255+
LinearAlgebra.require_one_based_indexing(A)
256+
isskewsymmetric(A) || throw(ArgumentError("Pfaffian requires a skew-symmetric matrix"))
257+
return _logabspfaffian!(A)
258+
end

test/runtests.jl

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -414,6 +414,9 @@ end
414414
if VERSION v"1.7" # for exact det of BigInt matrices
415415
@test SkewLinearAlgebra.exactpfaffian(Abig)^2 == det(Abig)
416416
end
417+
418+
logpf, sign = logabspfaffian(A)
419+
@test pfaffian(A) sign * exp(logpf)
417420
end
418421

419422
# issue #49

0 commit comments

Comments
 (0)