@@ -136,19 +136,15 @@ end
136
136
# cholesky!. Destructive methods for computing Cholesky factorization of real symmetric
137
137
# or Hermitian matrix
138
138
# # No pivoting (default)
139
- function cholesky! (A:: RealHermSymComplexHerm , :: Val{false} = Val (false ))
140
- if A. uplo == ' U'
141
- CU, info = _chol! (A. data, UpperTriangular)
142
- Cholesky (CU. data, ' U' , info)
143
- else
144
- CL, info = _chol! (A. data, LowerTriangular)
145
- Cholesky (CL. data, ' L' , info)
146
- end
139
+ function cholesky! (A:: RealHermSymComplexHerm , :: Val{false} = Val (false ); check:: Bool = true )
140
+ C, info = _chol! (A. data, A. uplo == ' U' ? UpperTriangular : LowerTriangular)
141
+ check && checkpositivedefinite (info)
142
+ return Cholesky (C. data, A. uplo, info)
147
143
end
148
144
149
145
# ## for StridedMatrices, check that matrix is symmetric/Hermitian
150
146
"""
151
- cholesky!(A, Val(false)) -> Cholesky
147
+ cholesky!(A, Val(false); check = true ) -> Cholesky
152
148
153
149
The same as [`cholesky`](@ref), but saves space by overwriting the input `A`,
154
150
instead of creating a copy. An [`InexactError`](@ref) exception is thrown if
@@ -168,53 +164,58 @@ Stacktrace:
168
164
[...]
169
165
```
170
166
"""
171
- function cholesky! (A:: StridedMatrix , :: Val{false} = Val (false ))
167
+ function cholesky! (A:: StridedMatrix , :: Val{false} = Val (false ); check :: Bool = true )
172
168
checksquare (A)
173
169
if ! ishermitian (A) # return with info = -1 if not Hermitian
170
+ check && checkpositivedefinite (- 1 )
174
171
return Cholesky (A, ' U' , convert (BlasInt, - 1 ))
175
172
else
176
- return cholesky! (Hermitian (A), Val (false ))
173
+ return cholesky! (Hermitian (A), Val (false ); check = check )
177
174
end
178
175
end
179
176
180
177
181
178
# # With pivoting
182
179
# ## BLAS/LAPACK element types
183
180
function cholesky! (A:: RealHermSymComplexHerm{<:BlasReal,<:StridedMatrix} ,
184
- :: Val{true} ; tol = 0.0 )
181
+ :: Val{true} ; tol = 0.0 , check :: Bool = true )
185
182
AA, piv, rank, info = LAPACK. pstrf! (A. uplo, A. data, tol)
186
- return CholeskyPivoted {eltype(AA),typeof(AA)} (AA, A. uplo, piv, rank, tol, info)
183
+ C = CholeskyPivoted {eltype(AA),typeof(AA)} (AA, A. uplo, piv, rank, tol, info)
184
+ check && chkfullrank (C)
185
+ return C
187
186
end
188
187
189
188
# ## Non BLAS/LAPACK element types (generic). Since generic fallback for pivoted Cholesky
190
189
# ## is not implemented yet we throw an error
191
- cholesky! (A:: RealHermSymComplexHerm{<:Real} , :: Val{true} ; tol = 0.0 ) =
190
+ cholesky! (A:: RealHermSymComplexHerm{<:Real} , :: Val{true} ; tol = 0.0 , check :: Bool = true ) =
192
191
throw (ArgumentError (" generic pivoted Cholesky factorization is not implemented yet" ))
193
192
194
193
# ## for StridedMatrices, check that matrix is symmetric/Hermitian
195
194
"""
196
- cholesky!(A, Val(true); tol = 0.0) -> CholeskyPivoted
195
+ cholesky!(A, Val(true); tol = 0.0, check = true ) -> CholeskyPivoted
197
196
198
197
The same as [`cholesky`](@ref), but saves space by overwriting the input `A`,
199
198
instead of creating a copy. An [`InexactError`](@ref) exception is thrown if the
200
199
factorization produces a number not representable by the element type of `A`,
201
200
e.g. for integer types.
202
201
"""
203
- function cholesky! (A:: StridedMatrix , :: Val{true} ; tol = 0.0 )
202
+ function cholesky! (A:: StridedMatrix , :: Val{true} ; tol = 0.0 , check :: Bool = true )
204
203
checksquare (A)
205
- if ! ishermitian (A) # return with info = -1 if not Hermitian
206
- return CholeskyPivoted (A, ' U' , Vector {BlasInt} (),convert (BlasInt, 1 ),
207
- tol, convert (BlasInt, - 1 ))
204
+ if ! ishermitian (A)
205
+ C = CholeskyPivoted (A, ' U' , Vector {BlasInt} (),convert (BlasInt, 1 ),
206
+ tol, convert (BlasInt, - 1 ))
207
+ check && chkfullrank (C)
208
+ return C
208
209
else
209
- return cholesky! (Hermitian (A), Val (true ); tol = tol)
210
+ return cholesky! (Hermitian (A), Val (true ); tol = tol, check = check )
210
211
end
211
212
end
212
213
213
214
# cholesky. Non-destructive methods for computing Cholesky factorization of real symmetric
214
215
# or Hermitian matrix
215
216
# # No pivoting (default)
216
217
"""
217
- cholesky(A, Val(false)) -> Cholesky
218
+ cholesky(A, Val(false); check = true ) -> Cholesky
218
219
219
220
Compute the Cholesky factorization of a dense symmetric positive definite matrix `A`
220
221
and return a `Cholesky` factorization. The matrix `A` can either be a [`Symmetric`](@ref) or [`Hermitian`](@ref)
@@ -223,6 +224,10 @@ The triangular Cholesky factor can be obtained from the factorization `F` with:
223
224
The following functions are available for `Cholesky` objects: [`size`](@ref), [`\\ `](@ref),
224
225
[`inv`](@ref), [`det`](@ref), [`logdet`](@ref) and [`isposdef`](@ref).
225
226
227
+ When `check = true`, an error is thrown if the decomposition fails.
228
+ When `check = false`, responsibility for checking the decomposition's
229
+ validity (via [`issuccess`](@ref)) lies with the user.
230
+
226
231
# Examples
227
232
```jldoctest
228
233
julia> A = [4. 12. -16.; 12. 37. -43.; -16. -43. 98.]
@@ -256,12 +261,12 @@ true
256
261
```
257
262
"""
258
263
cholesky (A:: Union{StridedMatrix,RealHermSymComplexHerm{<:Real,<:StridedMatrix}} ,
259
- :: Val{false} = Val (false )) = cholesky! (cholcopy (A))
264
+ :: Val{false} = Val (false ); check :: Bool = true ) = cholesky! (cholcopy (A); check = check )
260
265
261
266
262
267
# # With pivoting
263
268
"""
264
- cholesky(A, Val(true); tol = 0.0) -> CholeskyPivoted
269
+ cholesky(A, Val(true); tol = 0.0, check = true ) -> CholeskyPivoted
265
270
266
271
Compute the pivoted Cholesky factorization of a dense symmetric positive semi-definite matrix `A`
267
272
and return a `CholeskyPivoted` factorization. The matrix `A` can either be a [`Symmetric`](@ref)
@@ -271,9 +276,14 @@ The following functions are available for `PivotedCholesky` objects:
271
276
[`size`](@ref), [`\\ `](@ref), [`inv`](@ref), [`det`](@ref), and [`rank`](@ref).
272
277
The argument `tol` determines the tolerance for determining the rank.
273
278
For negative values, the tolerance is the machine precision.
279
+
280
+ When `check = true`, an error is thrown if the decomposition fails.
281
+ When `check = false`, responsibility for checking the decomposition's
282
+ validity (via [`issuccess`](@ref)) lies with the user.
274
283
"""
275
284
cholesky (A:: Union{StridedMatrix,RealHermSymComplexHerm{<:Real,<:StridedMatrix}} ,
276
- :: Val{true} ; tol = 0.0 ) = cholesky! (cholcopy (A), Val (true ); tol = tol)
285
+ :: Val{true} ; tol = 0.0 , check:: Bool = true ) =
286
+ cholesky! (cholcopy (A), Val (true ); tol = tol, check = check)
277
287
278
288
# # Number
279
289
function cholesky (x:: Number , uplo:: Symbol = :U )
@@ -319,11 +329,11 @@ function getproperty(C::Cholesky, d::Symbol)
319
329
Cuplo = getfield (C, :uplo )
320
330
info = getfield (C, :info )
321
331
if d == :U
322
- return @assertposdef UpperTriangular (Symbol (Cuplo) == d ? Cfactors : copy (Cfactors' )) info
332
+ return UpperTriangular (Symbol (Cuplo) == d ? Cfactors : copy (Cfactors' ))
323
333
elseif d == :L
324
- return @assertposdef LowerTriangular (Symbol (Cuplo) == d ? Cfactors : copy (Cfactors' )) info
334
+ return LowerTriangular (Symbol (Cuplo) == d ? Cfactors : copy (Cfactors' ))
325
335
elseif d == :UL
326
- return @assertposdef (Symbol (Cuplo) == :U ? UpperTriangular (Cfactors) : LowerTriangular (Cfactors)) info
336
+ return (Symbol (Cuplo) == :U ? UpperTriangular (Cfactors) : LowerTriangular (Cfactors))
327
337
else
328
338
return getfield (C, d)
329
339
end
@@ -335,16 +345,12 @@ function getproperty(C::CholeskyPivoted{T}, d::Symbol) where T<:BlasFloat
335
345
Cfactors = getfield (C, :factors )
336
346
Cuplo = getfield (C, :uplo )
337
347
if d == :U
338
- chkfullrank (C)
339
348
return UpperTriangular (Symbol (Cuplo) == d ? Cfactors : copy (Cfactors' ))
340
349
elseif d == :L
341
- chkfullrank (C)
342
350
return LowerTriangular (Symbol (Cuplo) == d ? Cfactors : copy (Cfactors' ))
343
351
elseif d == :p
344
- chkfullrank (C)
345
352
return getfield (C, :piv )
346
353
elseif d == :P
347
- chkfullrank (C)
348
354
n = size (C, 1 )
349
355
P = zeros (T, n, n)
350
356
for i = 1 : n
@@ -379,7 +385,7 @@ function show(io::IO, mime::MIME{Symbol("text/plain")}, C::CholeskyPivoted{<:Any
379
385
end
380
386
381
387
ldiv! (C:: Cholesky{T,<:AbstractMatrix} , B:: StridedVecOrMat{T} ) where {T<: BlasFloat } =
382
- @assertposdef LAPACK. potrs! (C. uplo, C. factors, B) C . info
388
+ LAPACK. potrs! (C. uplo, C. factors, B)
383
389
384
390
function ldiv! (C:: Cholesky{<:Any,<:AbstractMatrix} , B:: StridedVecOrMat )
385
391
if C. uplo == ' L'
@@ -390,11 +396,9 @@ function ldiv!(C::Cholesky{<:Any,<:AbstractMatrix}, B::StridedVecOrMat)
390
396
end
391
397
392
398
function ldiv! (C:: CholeskyPivoted{T} , B:: StridedVector{T} ) where T<: BlasFloat
393
- chkfullrank (C)
394
399
invpermute! (LAPACK. potrs! (C. uplo, C. factors, permute! (B, C. piv)), C. piv)
395
400
end
396
401
function ldiv! (C:: CholeskyPivoted{T} , B:: StridedMatrix{T} ) where T<: BlasFloat
397
- chkfullrank (C)
398
402
n = size (C, 1 )
399
403
for i= 1 : size (B, 2 )
400
404
permute! (view (B, 1 : n, i), C. piv)
429
433
isposdef (C:: Union{Cholesky,CholeskyPivoted} ) = C. info == 0
430
434
431
435
function det (C:: Cholesky )
432
- isposdef (C) || throw (PosDefException (C. info))
433
436
dd = one (real (eltype (C)))
434
437
@inbounds for i in 1 : size (C. factors,1 )
435
438
dd *= real (C. factors[i,i])^ 2
@@ -438,8 +441,6 @@ function det(C::Cholesky)
438
441
end
439
442
440
443
function logdet (C:: Cholesky )
441
- # need to check first, or log will throw DomainError
442
- isposdef (C) || throw (PosDefException (C. info))
443
444
dd = zero (real (eltype (C)))
444
445
@inbounds for i in 1 : size (C. factors,1 )
445
446
dd += log (real (C. factors[i,i]))
@@ -472,12 +473,11 @@ function logdet(C::CholeskyPivoted)
472
473
end
473
474
474
475
inv! (C:: Cholesky{<:BlasFloat,<:StridedMatrix} ) =
475
- @assertposdef copytri! (LAPACK. potri! (C. uplo, C. factors), C. uplo, true ) C . info
476
+ copytri! (LAPACK. potri! (C. uplo, C. factors), C. uplo, true )
476
477
477
478
inv (C:: Cholesky{<:BlasFloat,<:StridedMatrix} ) = inv! (copy (C))
478
479
479
480
function inv (C:: CholeskyPivoted )
480
- chkfullrank (C)
481
481
ipiv = invperm (C. piv)
482
482
copytri! (LAPACK. potri! (C. uplo, copy (C. factors)), C. uplo, true )[ipiv, ipiv]
483
483
end
0 commit comments