90 lines
26 KiB
HTML
90 lines
26 KiB
HTML
|
<!DOCTYPE html>
|
|||
|
<html lang="en"><head><meta charset="UTF-8"/><meta name="viewport" content="width=device-width, initial-scale=1.0"/><title>Linear algebra · The Julia Language</title><script>(function(i,s,o,g,r,a,m){i['GoogleAnalyticsObject']=r;i[r]=i[r]||function(){
|
|||
|
(i[r].q=i[r].q||[]).push(arguments)},i[r].l=1*new Date();a=s.createElement(o),
|
|||
|
m=s.getElementsByTagName(o)[0];a.async=1;a.src=g;m.parentNode.insertBefore(a,m)
|
|||
|
})(window,document,'script','https://www.google-analytics.com/analytics.js','ga');
|
|||
|
|
|||
|
ga('create', 'UA-28835595-6', 'auto');
|
|||
|
ga('send', 'pageview');
|
|||
|
</script><link href="https://cdnjs.cloudflare.com/ajax/libs/normalize/4.2.0/normalize.min.css" rel="stylesheet" type="text/css"/><link href="https://fonts.googleapis.com/css?family=Lato|Roboto+Mono" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/4.6.3/css/font-awesome.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/highlight.js/9.12.0/styles/default.min.css" rel="stylesheet" type="text/css"/><script>documenterBaseURL=".."</script><script src="https://cdnjs.cloudflare.com/ajax/libs/require.js/2.2.0/require.min.js" data-main="../assets/documenter.js"></script><script src="../siteinfo.js"></script><script src="../../versions.js"></script><link href="../assets/documenter.css" rel="stylesheet" type="text/css"/><link href="../assets/julia-manual.css" rel="stylesheet" type="text/css"/></head><body><nav class="toc"><a href="../index.html"><img class="logo" src="../assets/logo.png" alt="The Julia Language logo"/></a><h1>The Julia Language</h1><select id="version-selector" onChange="window.location.href=this.value" style="visibility: hidden"></select><form class="search" id="search-form" action="../search.html"><input id="search-query" name="q" type="text" placeholder="Search docs"/></form><ul><li><a class="toctext" href="../index.html">Home</a></li><li><span class="toctext">Manual</span><ul><li><a class="toctext" href="introduction.html">Introduction</a></li><li><a class="toctext" href="getting-started.html">Getting Started</a></li><li><a class="toctext" href="variables.html">Variables</a></li><li><a class="toctext" href="integers-and-floating-point-numbers.html">Integers and Floating-Point Numbers</a></li><li><a class="toctext" href="mathematical-operations.html">Mathematical Operations and Elementary Functions</a></li><li><a class="toctext" href="complex-and-rational-numbers.html">Complex and Rational Numbers</a></li><li><a class="toctext" href="strings.html">Strings</a></li><li><a class="toctext" href="functions.html">Functions</a></li><li><a class="toctext" href="control-flow.html">Control Flow</a></li><li><a class="toctext" href="variables-and-scoping.html">Scope of Variables</a></li><li><a class="toctext" href="types.html">Types</a></li><li><a class="toctext" href="methods.html">Methods</a></li><li><a class="toctext" href="constructors.html">Constructors</a></li><li><a class="toctext" href="conversion-and-promotion.html">Conversion and Promotion</a></li><li><a class="toctext" href="interfaces.html">Interfaces</a></li><li><a class="toctext" href="modules.html">Modules</a></li><li><a class="toctext" href="documentation.html">Documentation</a></li><li><a class="toctext" href="metaprogramming.html">Metaprogramming</a></li><li><a class="toctext" href="arrays.html">Multi-dimensional Arrays</a></li><li class="current"><a class="toctext" href="linear-algebra.html">Linear algebra</a><ul class="internal"><li><a class="toctext" href="#Special-matrices-1">Special matrices</a></li><li><a class="toctext" href="#man-linalg-factorizations-1">Matrix factorizations</a></li></ul></li><li><a class="toctext" href="networking-and-streams.html">Networking and Streams</a></li><li><a class="toctext" href="parallel-computing.html">Parallel Computing</a></li><li><a class="toctext" href="dates.html">Date and DateTime</a></li><li><a class="toctext" href="interacting-with-julia.html">Interacting With Julia</a></li><li><a class="toctext" href="running-external-programs.html">Running External Programs</a></li><li><a class="toctext" href="calling-c-and-fortran-code.html">Calling C and Fortran Code</a></li><li><a class="toctext" href="handling-operating-system-variation.html">Handling Operating System Variation</a></li><li><a class="toctext" href="environment-variables.html">Environment Variables</a></li><li><a class="toctext" href="embedding.html">Embedding Julia</a></li><li><a class="toctext" href="packages.html">Packages</a></li><li><a class="toctext" href="profile.html">Profiling</a></li><li><a class="toctext" href="stacktraces.html">Stack Traces</a></li><l
|
|||
|
3×3 Array{Int64,2}:
|
|||
|
1 2 3
|
|||
|
4 1 6
|
|||
|
7 8 1
|
|||
|
|
|||
|
julia> trace(A)
|
|||
|
3
|
|||
|
|
|||
|
julia> det(A)
|
|||
|
104.0
|
|||
|
|
|||
|
julia> inv(A)
|
|||
|
3×3 Array{Float64,2}:
|
|||
|
-0.451923 0.211538 0.0865385
|
|||
|
0.365385 -0.192308 0.0576923
|
|||
|
0.240385 0.0576923 -0.0673077</code></pre><p>As well as other useful operations, such as finding eigenvalues or eigenvectors:</p><pre><code class="language-julia-repl">julia> A = [1.5 2 -4; 3 -1 -6; -10 2.3 4]
|
|||
|
3×3 Array{Float64,2}:
|
|||
|
1.5 2.0 -4.0
|
|||
|
3.0 -1.0 -6.0
|
|||
|
-10.0 2.3 4.0
|
|||
|
|
|||
|
julia> eigvals(A)
|
|||
|
3-element Array{Complex{Float64},1}:
|
|||
|
9.31908+0.0im
|
|||
|
-2.40954+2.72095im
|
|||
|
-2.40954-2.72095im
|
|||
|
|
|||
|
julia> eigvecs(A)
|
|||
|
3×3 Array{Complex{Float64},2}:
|
|||
|
-0.488645+0.0im 0.182546-0.39813im 0.182546+0.39813im
|
|||
|
-0.540358+0.0im 0.692926+0.0im 0.692926-0.0im
|
|||
|
0.68501+0.0im 0.254058-0.513301im 0.254058+0.513301im</code></pre><p>In addition, Julia provides many <a href="linear-algebra.html#man-linalg-factorizations-1">factorizations</a> which can be used to speed up problems such as linear solve or matrix exponentiation by pre-factorizing a matrix into a form more amenable (for performance or memory reasons) to the problem. See the documentation on <a href="../stdlib/linalg.html#Base.LinAlg.factorize"><code>factorize</code></a> for more information. As an example:</p><pre><code class="language-julia-repl">julia> A = [1.5 2 -4; 3 -1 -6; -10 2.3 4]
|
|||
|
3×3 Array{Float64,2}:
|
|||
|
1.5 2.0 -4.0
|
|||
|
3.0 -1.0 -6.0
|
|||
|
-10.0 2.3 4.0
|
|||
|
|
|||
|
julia> factorize(A)
|
|||
|
Base.LinAlg.LU{Float64,Array{Float64,2}} with factors L and U:
|
|||
|
[1.0 0.0 0.0; -0.15 1.0 0.0; -0.3 -0.132196 1.0]
|
|||
|
[-10.0 2.3 4.0; 0.0 2.345 -3.4; 0.0 0.0 -5.24947]</code></pre><p>Since <code>A</code> is not Hermitian, symmetric, triangular, tridiagonal, or bidiagonal, an LU factorization may be the best we can do. Compare with:</p><pre><code class="language-julia-repl">julia> B = [1.5 2 -4; 2 -1 -3; -4 -3 5]
|
|||
|
3×3 Array{Float64,2}:
|
|||
|
1.5 2.0 -4.0
|
|||
|
2.0 -1.0 -3.0
|
|||
|
-4.0 -3.0 5.0
|
|||
|
|
|||
|
julia> factorize(B)
|
|||
|
Base.LinAlg.BunchKaufman{Float64,Array{Float64,2}}([-1.64286 0.142857 -0.8; 2.0 -2.8 -0.6; -4.0 -3.0 5.0], [1, 2, 3], 'U', true, false, 0)</code></pre><p>Here, Julia was able to detect that <code>B</code> is in fact symmetric, and used a more appropriate factorization. Often it's possible to write more efficient code for a matrix that is known to have certain properties e.g. it is symmetric, or tridiagonal. Julia provides some special types so that you can "tag" matrices as having these properties. For instance:</p><pre><code class="language-julia-repl">julia> B = [1.5 2 -4; 2 -1 -3; -4 -3 5]
|
|||
|
3×3 Array{Float64,2}:
|
|||
|
1.5 2.0 -4.0
|
|||
|
2.0 -1.0 -3.0
|
|||
|
-4.0 -3.0 5.0
|
|||
|
|
|||
|
julia> sB = Symmetric(B)
|
|||
|
3×3 Symmetric{Float64,Array{Float64,2}}:
|
|||
|
1.5 2.0 -4.0
|
|||
|
2.0 -1.0 -3.0
|
|||
|
-4.0 -3.0 5.0</code></pre><p><code>sB</code> has been tagged as a matrix that's (real) symmetric, so for later operations we might perform on it, such as eigenfactorization or computing matrix-vector products, efficiencies can be found by only referencing half of it. For example:</p><pre><code class="language-julia-repl">julia> B = [1.5 2 -4; 2 -1 -3; -4 -3 5]
|
|||
|
3×3 Array{Float64,2}:
|
|||
|
1.5 2.0 -4.0
|
|||
|
2.0 -1.0 -3.0
|
|||
|
-4.0 -3.0 5.0
|
|||
|
|
|||
|
julia> sB = Symmetric(B)
|
|||
|
3×3 Symmetric{Float64,Array{Float64,2}}:
|
|||
|
1.5 2.0 -4.0
|
|||
|
2.0 -1.0 -3.0
|
|||
|
-4.0 -3.0 5.0
|
|||
|
|
|||
|
julia> x = [1; 2; 3]
|
|||
|
3-element Array{Int64,1}:
|
|||
|
1
|
|||
|
2
|
|||
|
3
|
|||
|
|
|||
|
julia> sB\x
|
|||
|
3-element Array{Float64,1}:
|
|||
|
-1.73913
|
|||
|
-1.1087
|
|||
|
-1.45652</code></pre><p>The <code>\</code> operation here performs the linear solution. Julia's parser provides convenient dispatch to specialized methods for the <em>transpose</em> of a matrix left-divided by a vector, or for the various combinations of transpose operations in matrix-matrix solutions. Many of these are further specialized for certain special matrix types. For example, <code>A\B</code> will end up calling <a href="../stdlib/linalg.html#Base.LinAlg.A_ldiv_B!"><code>Base.LinAlg.A_ldiv_B!</code></a> while <code>A'\B</code> will end up calling <a href="../stdlib/linalg.html#Base.Ac_ldiv_B"><code>Base.LinAlg.Ac_ldiv_B</code></a>, even though we used the same left-division operator. This works for matrices too: <code>A.'\B.'</code> would call <a href="../stdlib/linalg.html#Base.At_ldiv_Bt"><code>Base.LinAlg.At_ldiv_Bt</code></a>. The left-division operator is pretty powerful and it's easy to write compact, readable code that is flexible enough to solve all sorts of systems of linear equations.</p><h2><a class="nav-anchor" id="Special-matrices-1" href="#Special-matrices-1">Special matrices</a></h2><p><a href="http://www2.imm.dtu.dk/pubdb/views/publication_details.php?id=3274">Matrices with special symmetries and structures</a> arise often in linear algebra and are frequently associated with various matrix factorizations. Julia features a rich collection of special matrix types, which allow for fast computation with specialized routines that are specially developed for particular matrix types.</p><p>The following tables summarize the types of special matrices that have been implemented in Julia, as well as whether hooks to various optimized methods for them in LAPACK are available.</p><table><tr><th>Type</th><th>Description</th></tr><tr><td><a href="../stdlib/linalg.html#Base.LinAlg.Hermitian"><code>Hermitian</code></a></td><td><a href="https://en.wikipedia.org/wiki/Hermitian_matrix">Hermitian matrix</a></td></tr><tr><td><a href="../stdlib/linalg.html#Base.LinAlg.UpperTriangular"><code>UpperTriangular</code></a></td><td>Upper <a href="https://en.wikipedia.org/wiki/Triangular_matrix">triangular matrix</a></td></tr><tr><td><a href="../stdlib/linalg.html#Base.LinAlg.LowerTriangular"><code>LowerTriangular</code></a></td><td>Lower <a href="https://en.wikipedia.org/wiki/Triangular_matrix">triangular matrix</a></td></tr><tr><td><a href="../stdlib/linalg.html#Base.LinAlg.Tridiagonal"><code>Tridiagonal</code></a></td><td><a href="https://en.wikipedia.org/wiki/Tridiagonal_matrix">Tridiagonal matrix</a></td></tr><tr><td><a href="../stdlib/linalg.html#Base.LinAlg.SymTridiagonal"><code>SymTridiagonal</code></a></td><td>Symmetric tridiagonal matrix</td></tr><tr><td><a href="../stdlib/linalg.html#Base.LinAlg.Bidiagonal"><code>Bidiagonal</code></a></td><td>Upper/lower <a href="https://en.wikipedia.org/wiki/Bidiagonal_matrix">bidiagonal matrix</a></td></tr><tr><td><a href="../stdlib/linalg.html#Base.LinAlg.Diagonal"><code>Diagonal</code></a></td><td><a href="https://en.wikipedia.org/wiki/Diagonal_matrix">Diagonal matrix</a></td></tr><tr><td><code>UniformScaling</code></td><td><a href="https://en.wikipedia.org/wiki/Uniform_scaling">Uniform scaling operator</a></td></tr></table><h3><a class="nav-anchor" id="Elementary-operations-1" href="#Elementary-operations-1">Elementary operations</a></h3><table><tr><th>Matrix type</th><th><code>+</code></th><th><code>-</code></th><th><code>*</code></th><th><code>\</code></th><th>Other functions with optimized methods</th></tr><tr><td><a href="../stdlib/linalg.html#Base.LinAlg.Hermitian"><code>Hermitian</code></a></td><td> </td><td> </td><td> </td><td>MV</td><td><a href="../stdlib/linalg.html#Base.inv"><code>inv()</code></a>, <a href="../stdlib/linalg.html#Base.LinAlg.sqrtm"><code>sqrtm()</code></a>, <a href="../stdlib/linalg.html#Base.LinAlg.expm"><code>expm()</code></a></td></tr><tr><td><a href="../stdlib/linalg.html#Base.LinAlg.UpperTriangular"><code>UpperTriangular</code></a></td><td> </td><td> </td><td>MV</td><td>MV</td><td><a href="../stdlib/linalg.html#Base.in
|