@@ -164,6 +164,28 @@ Same as [`eigvals`](@ref), but saves space by overwriting the input `A`, instead
164
164
The option `permute=true` permutes the matrix to become
165
165
closer to upper triangular, and `scale=true` scales the matrix by its diagonal elements to
166
166
make rows and columns more equal in norm.
167
+
168
+ !!! note
169
+ The input matrix `A` will not contain its eigenvalues after `eigvals!` is
170
+ called on it - `A` is used as a workspace.
171
+
172
+ # Examples
173
+ ```jldoctest
174
+ julia> A = [1. 2.; 3. 4.]
175
+ 2×2 Array{Float64,2}:
176
+ 1.0 2.0
177
+ 3.0 4.0
178
+
179
+ julia> eigvals!(A)
180
+ 2-element Array{Float64,1}:
181
+ -0.3722813232690143
182
+ 5.372281323269014
183
+
184
+ julia> A
185
+ 2×2 Array{Float64,2}:
186
+ -0.372281 -1.0
187
+ 0.0 5.37228
188
+ ```
167
189
"""
168
190
function eigvals! (A:: StridedMatrix{<:BlasReal} ; permute:: Bool = true , scale:: Bool = true )
169
191
issymmetric (A) && return eigvals! (Symmetric (A))
@@ -185,6 +207,19 @@ before the eigenvalue calculation. The option `permute=true` permutes the matrix
185
207
become closer to upper triangular, and `scale=true` scales the matrix by its diagonal
186
208
elements to make rows and columns more equal in norm. The default is `true` for both
187
209
options.
210
+
211
+ # Examples
212
+ ```jldoctest
213
+ julia> diag_matrix = [1 0; 0 4]
214
+ 2×2 Array{Int64,2}:
215
+ 1 0
216
+ 0 4
217
+
218
+ julia> eigvals(diag_matrix)
219
+ 2-element Array{Float64,1}:
220
+ 1.0
221
+ 4.0
222
+ ```
188
223
"""
189
224
function eigvals (A:: StridedMatrix{T} ; permute:: Bool = true , scale:: Bool = true ) where T
190
225
S = promote_type (Float32, typeof (one (T)/ norm (one (T))))
@@ -317,6 +352,31 @@ Computes the generalized eigenvalue decomposition of `A` and `B`, returning a
317
352
`GeneralizedEigen` factorization object `F` which contains the generalized eigenvalues in
318
353
`F[:values]` and the generalized eigenvectors in the columns of the matrix `F[:vectors]`.
319
354
(The `k`th generalized eigenvector can be obtained from the slice `F[:vectors][:, k]`.)
355
+
356
+ # Examples
357
+ ```jldoctest
358
+ julia> A = [1 0; 0 -1]
359
+ 2×2 Array{Int64,2}:
360
+ 1 0
361
+ 0 -1
362
+
363
+ julia> B = [0 1; 1 0]
364
+ 2×2 Array{Int64,2}:
365
+ 0 1
366
+ 1 0
367
+
368
+ julia> F = eigfact(A, B);
369
+
370
+ julia> F[:values]
371
+ 2-element Array{Complex{Float64},1}:
372
+ 0.0 + 1.0im
373
+ 0.0 - 1.0im
374
+
375
+ julia> F[:vectors]
376
+ 2×2 Array{Complex{Float64},2}:
377
+ 0.0-1.0im 0.0+1.0im
378
+ -1.0-0.0im -1.0+0.0im
379
+ ```
320
380
"""
321
381
function eigfact (A:: AbstractMatrix{TA} , B:: AbstractMatrix{TB} ) where {TA,TB}
322
382
S = promote_type (Float32, typeof (one (TA)/ norm (one (TA))),TB)
361
421
"""
362
422
eigvals!(A, B) -> values
363
423
364
- Same as [`eigvals`](@ref), but saves space by overwriting the input `A` (and `B`), instead of creating copies.
424
+ Same as [`eigvals`](@ref), but saves space by overwriting the input `A` (and `B`),
425
+ instead of creating copies.
426
+
427
+ !!! note
428
+ The input matrices `A` and `B` will not contain their eigenvalues after
429
+ `eigvals!` is called. They are used as workspaces.
430
+
431
+ # Examples
432
+ ```jldoctest
433
+ julia> A = [1. 0.; 0. -1.]
434
+ 2×2 Array{Float64,2}:
435
+ 1.0 0.0
436
+ 0.0 -1.0
437
+
438
+ julia> B = [0. 1.; 1. 0.]
439
+ 2×2 Array{Float64,2}:
440
+ 0.0 1.0
441
+ 1.0 0.0
442
+
443
+ julia> eigvals!(A, B)
444
+ 2-element Array{Complex{Float64},1}:
445
+ 0.0 + 1.0im
446
+ 0.0 - 1.0im
447
+
448
+ julia> A
449
+ 2×2 Array{Float64,2}:
450
+ -0.0 -1.0
451
+ 1.0 -0.0
452
+
453
+ julia> B
454
+ 2×2 Array{Float64,2}:
455
+ 1.0 0.0
456
+ 0.0 1.0
457
+ ```
365
458
"""
366
459
function eigvals! (A:: StridedMatrix{T} , B:: StridedMatrix{T} ) where T<: BlasReal
367
460
issymmetric (A) && isposdef (B) && return eigvals! (Symmetric (A), Symmetric (B))
0 commit comments