cuBLAS
======

.. raw:: html

    <p>Provides basic linear algebra building blocks. See <a class="reference external" href="http://docs.nvidia.com/cuda/cublas/index.html">NVIDIA cuBLAS</a>.</p>
    <p>The cuBLAS binding provides an interface that accepts NumPy arrays and Numba&#8217;s
    CUDA device arrays. The binding automatically transfers NumPy array arguments to
    the device as required. This automatic transfer may generate some unnecessary
    transfers, so optimal performance is likely to be obtained by the manual
    transfer for NumPy arrays into device arrays and using the cuBLAS to manipulate
    device arrays where possible.</p>
    <p>No special naming convention is used to identify the data
    type, unlike in the BLAS C and Fortran APIs. Arguments for array storage
    information which are part of the cuBLAS C API are also not necessary since
    NumPy arrays and device arrays contain this information.</p>
    <p>All functions are accessed through the <a class="reference internal" href="#accelerate.cuda.blas.Blas" title="accelerate.cuda.blas.Blas"><code class="xref py py-class docutils literal"><span class="pre">accelerate.cuda.blas.Blas</span></code></a> class:</p>
    <dl class="class">
    <dt id="accelerate.cuda.blas.Blas">
    <em class="property">class </em><code class="descclassname">accelerate.cuda.blas.</code><code class="descname">Blas</code><span class="sig-paren">(</span><em>*args</em>, <em>**kws</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas" title="Permalink to this definition">¶</a></dt>
    <dd><p>All BLAS subprograms are available under the Blas object.</p>
    <table class="docutils field-list" frame="void" rules="none">
    <col class="field-name" />
    <col class="field-body" />
    <tbody valign="top">
    <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><strong>stream</strong> &#8211; Optional. A CUDA Stream.</td>
    </tr>
    </tbody>
    </table>
    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.amax">
    <code class="descname">amax</code><span class="sig-paren">(</span><em>x</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.amax" title="Permalink to this definition">¶</a></dt>
    <dd><p>Same as np.argmax(x)</p>
    </dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.amin">
    <code class="descname">amin</code><span class="sig-paren">(</span><em>x</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.amin" title="Permalink to this definition">¶</a></dt>
    <dd><p>Same as np.argmin(x)</p>
    </dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.asum">
    <code class="descname">asum</code><span class="sig-paren">(</span><em>x</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.asum" title="Permalink to this definition">¶</a></dt>
    <dd><p>Same as np.sum(x)</p>
    </dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.axpy">
    <code class="descname">axpy</code><span class="sig-paren">(</span><em>alpha</em>, <em>x</em>, <em>y</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.axpy" title="Permalink to this definition">¶</a></dt>
    <dd><p>Same as y = alpha * x + y</p>
    </dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.dgmm">
    <code class="descname">dgmm</code><span class="sig-paren">(</span><em>side</em>, <em>m</em>, <em>n</em>, <em>A:r</em>, <em>x:r</em>, <em>C:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.dgmm" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.dot">
    <code class="descname">dot</code><span class="sig-paren">(</span><em>x</em>, <em>y</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.dot" title="Permalink to this definition">¶</a></dt>
    <dd><p>Same as np.dot</p>
    </dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.dotc">
    <code class="descname">dotc</code><span class="sig-paren">(</span><em>x</em>, <em>y</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.dotc" title="Permalink to this definition">¶</a></dt>
    <dd><p>Same as np.vdot</p>
    </dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.gbmv">
    <code class="descname">gbmv</code><span class="sig-paren">(</span><em>trans</em>, <em>m</em>, <em>n</em>, <em>kl</em>, <em>ku</em>, <em>alpha</em>, <em>A:r</em>, <em>x:r</em>, <em>beta</em>, <em>y:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.gbmv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.geam">
    <code class="descname">geam</code><span class="sig-paren">(</span><em>transa</em>, <em>transb</em>, <em>m</em>, <em>n</em>, <em>alpha</em>, <em>A:r</em>, <em>beta</em>, <em>B:r</em>, <em>C:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.geam" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.gemm">
    <code class="descname">gemm</code><span class="sig-paren">(</span><em>transa</em>, <em>transb</em>, <em>m</em>, <em>n</em>, <em>k</em>, <em>alpha</em>, <em>A:r</em>, <em>B:r</em>, <em>beta</em>, <em>C:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.gemm" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.gemv">
    <code class="descname">gemv</code><span class="sig-paren">(</span><em>trans</em>, <em>m</em>, <em>n</em>, <em>alpha</em>, <em>A:r</em>, <em>x:r</em>, <em>beta</em>, <em>y:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.gemv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.ger">
    <code class="descname">ger</code><span class="sig-paren">(</span><em>m</em>, <em>n</em>, <em>alpha</em>, <em>x:r</em>, <em>y:r</em>, <em>A:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.ger" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.gerc">
    <code class="descname">gerc</code><span class="sig-paren">(</span><em>m</em>, <em>n</em>, <em>alpha</em>, <em>x:r</em>, <em>y:r</em>, <em>A:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.gerc" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.geru">
    <code class="descname">geru</code><span class="sig-paren">(</span><em>m</em>, <em>n</em>, <em>alpha</em>, <em>x:r</em>, <em>y:r</em>, <em>A:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.geru" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.hbmv">
    <code class="descname">hbmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>k</em>, <em>alpha</em>, <em>A:r</em>, <em>x:r</em>, <em>beta</em>, <em>y:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.hbmv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.hemm">
    <code class="descname">hemm</code><span class="sig-paren">(</span><em>side</em>, <em>uplo</em>, <em>m</em>, <em>n</em>, <em>alpha</em>, <em>A:r</em>, <em>B:r</em>, <em>beta</em>, <em>C:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.hemm" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.hemv">
    <code class="descname">hemv</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>A:r</em>, <em>x:r</em>, <em>beta</em>, <em>y:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.hemv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.her">
    <code class="descname">her</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x:r</em>, <em>A:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.her" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.her2">
    <code class="descname">her2</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x:r</em>, <em>y:r</em>, <em>A:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.her2" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.herk">
    <code class="descname">herk</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>n</em>, <em>k</em>, <em>alpha</em>, <em>A:r</em>, <em>beta</em>, <em>C:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.herk" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.hpmv">
    <code class="descname">hpmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>AP:r</em>, <em>x:r</em>, <em>beta</em>, <em>y:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.hpmv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.hpr">
    <code class="descname">hpr</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x:r</em>, <em>AP:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.hpr" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.hpr2">
    <code class="descname">hpr2</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x:r</em>, <em>y:r</em>, <em>A:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.hpr2" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.nrm2">
    <code class="descname">nrm2</code><span class="sig-paren">(</span><em>x</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.nrm2" title="Permalink to this definition">¶</a></dt>
    <dd><p>Same as np.linalg.norm</p>
    </dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.rot">
    <code class="descname">rot</code><span class="sig-paren">(</span><em>x</em>, <em>y</em>, <em>c</em>, <em>s</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.rot" title="Permalink to this definition">¶</a></dt>
    <dd><p>Same as x, y = c * x + s * y, -s * x + c * y</p>
    </dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.rotg">
    <code class="descname">rotg</code><span class="sig-paren">(</span><em>a</em>, <em>b</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.rotg" title="Permalink to this definition">¶</a></dt>
    <dd><p>Compute the given rotation matrix given a column vector (a, b).
    Returns r, z, c, s.</p>
    <p>r &#8212; r = a ** 2 + b ** 2
    z &#8212; Use to recover c and s.</p>
    <blockquote>
    <div><dl class="docutils">
    <dt>if abs(z) &lt; 1:</dt>
    <dd>c, s = 1 - z ** 2, z</dd>
    <dt>elif abs(z) == 1:</dt>
    <dd>c, s = 0, 1</dd>
    <dt>else:</dt>
    <dd>c, s = 1 / z, 1 - z ** 2</dd>
    </dl>
    </div></blockquote>
    <p>c &#8212; Cosine element of the rotation matrix.
    s &#8212; Sine element of the rotation matrix.</p>
    </dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.rotm">
    <code class="descname">rotm</code><span class="sig-paren">(</span><em>x</em>, <em>y</em>, <em>param</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.rotm" title="Permalink to this definition">¶</a></dt>
    <dd><p>Applies the modified Givens transformation.</p>
    <p>x, y = h11 * x + h12 * y, h21 * x + h22 * y</p>
    <p>param &#8212; [flag, h11, h21, h12, h22]</p>
    <p>Refer to cuBLAS documentation for detail.</p>
    </dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.rotmg">
    <code class="descname">rotmg</code><span class="sig-paren">(</span><em>d1</em>, <em>d2</em>, <em>x1</em>, <em>y1</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.rotmg" title="Permalink to this definition">¶</a></dt>
    <dd><p>Constructs the modified Givens transformation.</p>
    <p>Returns param that is usable in rotm.</p>
    <p>Refer to cuBLAS documentation for detail.</p>
    </dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.sbmv">
    <code class="descname">sbmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>k</em>, <em>alpha</em>, <em>A:r</em>, <em>x:r</em>, <em>beta</em>, <em>y:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.sbmv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.scal">
    <code class="descname">scal</code><span class="sig-paren">(</span><em>alpha</em>, <em>x</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.scal" title="Permalink to this definition">¶</a></dt>
    <dd><p>Same as x = alpha * x</p>
    </dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.spmv">
    <code class="descname">spmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>AP:r</em>, <em>x:r</em>, <em>beta</em>, <em>y:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.spmv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.spr">
    <code class="descname">spr</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x:r</em>, <em>AP:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.spr" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.spr2">
    <code class="descname">spr2</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x:r</em>, <em>y:r</em>, <em>A:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.spr2" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.symm">
    <code class="descname">symm</code><span class="sig-paren">(</span><em>side</em>, <em>uplo</em>, <em>m</em>, <em>n</em>, <em>alpha</em>, <em>A:r</em>, <em>B:r</em>, <em>beta</em>, <em>C:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.symm" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.symv">
    <code class="descname">symv</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>A:r</em>, <em>x:r</em>, <em>beta</em>, <em>y:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.symv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.syr">
    <code class="descname">syr</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x:r</em>, <em>A:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.syr" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.syr2">
    <code class="descname">syr2</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x:r</em>, <em>y:r</em>, <em>A:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.syr2" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.syrk">
    <code class="descname">syrk</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>n</em>, <em>k</em>, <em>alpha</em>, <em>A:r</em>, <em>beta</em>, <em>C:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.syrk" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.tbmv">
    <code class="descname">tbmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>n</em>, <em>k</em>, <em>A:r</em>, <em>x:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.tbmv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.tbsv">
    <code class="descname">tbsv</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>n</em>, <em>k</em>, <em>A:r</em>, <em>x:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.tbsv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.tpmv">
    <code class="descname">tpmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>n</em>, <em>AP:r</em>, <em>x:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.tpmv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.tpsv">
    <code class="descname">tpsv</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>n</em>, <em>AP:r</em>, <em>x:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.tpsv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.trmm">
    <code class="descname">trmm</code><span class="sig-paren">(</span><em>side</em>, <em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>m</em>, <em>n</em>, <em>alpha</em>, <em>A:r</em>, <em>B:r</em>, <em>C:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.trmm" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.trmv">
    <code class="descname">trmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>n</em>, <em>A:r</em>, <em>x:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.trmv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.trsm">
    <code class="descname">trsm</code><span class="sig-paren">(</span><em>side</em>, <em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>m</em>, <em>n</em>, <em>alpha</em>, <em>A:r</em>, <em>B:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.trsm" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    <dl class="method">
    <dt id="accelerate.cuda.blas.Blas.trsv">
    <code class="descname">trsv</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>n</em>, <em>A:r</em>, <em>x:w</em><span class="sig-paren">)</span><a class="headerlink" href="#accelerate.cuda.blas.Blas.trsv" title="Permalink to this definition">¶</a></dt>
    <dd></dd></dl>

    </dd></dl>

    <div class="section" id="blas-level-1">
    <h2>BLAS Level 1<a class="headerlink" href="#blas-level-1" title="Permalink to this headline">¶</a></h2>
    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">nrm2</code><span class="sig-paren">(</span><em>x</em><span class="sig-paren">)</span></dt>
    <dd><p>Computes the L2 norm for array <cite>x</cite>. Same as <cite>numpy.linalg.norm(x)</cite>.</p>
    <table class="docutils field-list" frame="void" rules="none">
    <col class="field-name" />
    <col class="field-body" />
    <tbody valign="top">
    <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><strong>x</strong> (<em>array</em>) &#8211; input vector</td>
    </tr>
    <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body">resulting norm.</td>
    </tr>
    </tbody>
    </table>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">dot</code><span class="sig-paren">(</span><em>x</em>, <em>y</em><span class="sig-paren">)</span></dt>
    <dd><p>Compute the dot product of array <cite>x</cite> and array <cite>y</cite>.  Same as <cite>np.dot(x, y)</cite>.</p>
    <table class="docutils field-list" frame="void" rules="none">
    <col class="field-name" />
    <col class="field-body" />
    <tbody valign="top">
    <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><ul class="first simple">
    <li><strong>x</strong> (<em>array</em>) &#8211; vector</li>
    <li><strong>y</strong> (<em>array</em>) &#8211; vector</li>
    </ul>
    </td>
    </tr>
    <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body"><p class="first last">dot product of <cite>x</cite> and <cite>y</cite></p>
    </td>
    </tr>
    </tbody>
    </table>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">dotc</code><span class="sig-paren">(</span><em>x</em>, <em>y</em><span class="sig-paren">)</span></dt>
    <dd><p>Uses the conjugate of the element of the vectors to compute the dot product
    of array <cite>x</cite> and array <cite>y</cite> for complex dtype only.  Same as <cite>np.vdot(x, y)</cite>.</p>
    <table class="docutils field-list" frame="void" rules="none">
    <col class="field-name" />
    <col class="field-body" />
    <tbody valign="top">
    <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><ul class="first simple">
    <li><strong>x</strong> (<em>array</em>) &#8211; vector</li>
    <li><strong>y</strong> (<em>array</em>) &#8211; vector</li>
    </ul>
    </td>
    </tr>
    <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body"><p class="first last">dot product of <cite>x</cite> and <cite>y</cite></p>
    </td>
    </tr>
    </tbody>
    </table>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">scal</code><span class="sig-paren">(</span><em>alpha</em>, <em>x</em><span class="sig-paren">)</span></dt>
    <dd><p>Scale <cite>x</cite> inplace by alpha.  Same as <cite>x = alpha * x</cite></p>
    <table class="docutils field-list" frame="void" rules="none">
    <col class="field-name" />
    <col class="field-body" />
    <tbody valign="top">
    <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><ul class="first last simple">
    <li><strong>alpha</strong> &#8211; scalar</li>
    <li><strong>x</strong> (<em>array</em>) &#8211; vector</li>
    </ul>
    </td>
    </tr>
    </tbody>
    </table>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">axpy</code><span class="sig-paren">(</span><em>alpha</em>, <em>x</em><span class="sig-paren">)</span></dt>
    <dd><p>Compute <cite>y = alpha * x + y</cite> inplace.</p>
    <table class="docutils field-list" frame="void" rules="none">
    <col class="field-name" />
    <col class="field-body" />
    <tbody valign="top">
    <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><ul class="first last simple">
    <li><strong>alpha</strong> &#8211; scalar</li>
    <li><strong>x</strong> (<em>array</em>) &#8211; vector</li>
    </ul>
    </td>
    </tr>
    </tbody>
    </table>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">amax</code><span class="sig-paren">(</span><em>x</em><span class="sig-paren">)</span></dt>
    <dd><p>Find the index of the first largest element in array <cite>x</cite>.
    Same as <cite>np.argmax(x)</cite></p>
    <table class="docutils field-list" frame="void" rules="none">
    <col class="field-name" />
    <col class="field-body" />
    <tbody valign="top">
    <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><strong>x</strong> (<em>array</em>) &#8211; vector</td>
    </tr>
    <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body">index (start from 0).</td>
    </tr>
    </tbody>
    </table>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">amin</code><span class="sig-paren">(</span><em>x</em><span class="sig-paren">)</span></dt>
    <dd><p>Find the index of the first largest element in array <cite>x</cite>.
    Same as <cite>np.argmin(x)</cite></p>
    <table class="docutils field-list" frame="void" rules="none">
    <col class="field-name" />
    <col class="field-body" />
    <tbody valign="top">
    <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><strong>x</strong> (<em>array</em>) &#8211; vector</td>
    </tr>
    <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body">index (start from 0).</td>
    </tr>
    </tbody>
    </table>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">asum</code><span class="sig-paren">(</span><em>x</em><span class="sig-paren">)</span></dt>
    <dd><p>Compute the sum of all element in array <cite>x</cite>.</p>
    <table class="docutils field-list" frame="void" rules="none">
    <col class="field-name" />
    <col class="field-body" />
    <tbody valign="top">
    <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><strong>x</strong> (<em>array</em>) &#8211; vector</td>
    </tr>
    <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body"><cite>x.sum()</cite></td>
    </tr>
    </tbody>
    </table>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">rot</code><span class="sig-paren">(</span><em>x</em>, <em>y</em>, <em>c</em>, <em>s</em><span class="sig-paren">)</span></dt>
    <dd><p>Apply the Givens rotation matrix specified by the cosine element <cite>c</cite> and the
    sine element <cite>s</cite> inplace on vector element <cite>x</cite> and <cite>y</cite>.</p>
    <p>Same as <cite>x, y = c * x + s * y, -s * x + c * y</cite></p>
    <table class="docutils field-list" frame="void" rules="none">
    <col class="field-name" />
    <col class="field-body" />
    <tbody valign="top">
    <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><ul class="first last simple">
    <li><strong>x</strong> (<em>array</em>) &#8211; vector</li>
    <li><strong>y</strong> (<em>array</em>) &#8211; vector</li>
    </ul>
    </td>
    </tr>
    </tbody>
    </table>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">rotg</code><span class="sig-paren">(</span><em>a</em>, <em>b</em><span class="sig-paren">)</span></dt>
    <dd><p>Constructs the Givens rotation matrix with the column vector (a, b).</p>
    <table class="docutils field-list" frame="void" rules="none">
    <col class="field-name" />
    <col class="field-body" />
    <tbody valign="top">
    <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><ul class="first simple">
    <li><strong>a</strong> &#8211; first element of the column vector</li>
    <li><strong>b</strong> &#8211; second element of the column vector</li>
    </ul>
    </td>
    </tr>
    <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body"><p class="first last"><p>a tuple (r, z, c, s)</p>
    <p>r &#8211; <cite>r = a**2 + b**2</cite></p>
    <dl class="docutils">
    <dt>z &#8211; Use to reconstruct <cite>c</cite> and <cite>s</cite>.</dt>
    <dd><p class="first last">Refer to cuBLAS documentation for detail.</p>
    </dd>
    </dl>
    <p>c &#8211; The consine element.</p>
    <p>s &#8211; The sine element.</p>
    </p>
    </td>
    </tr>
    </tbody>
    </table>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">rotm</code><span class="sig-paren">(</span><em>x</em>, <em>y</em>, <em>param</em><span class="sig-paren">)</span></dt>
    <dd><p>Applies the modified Givens transformation inplace.</p>
    <p>Same as:</p>
    <div class="highlight-python"><div class="highlight"><pre><span></span><span class="n">param</span> <span class="o">=</span> <span class="n">flag</span><span class="p">,</span> <span class="n">h11</span><span class="p">,</span> <span class="n">h21</span><span class="p">,</span> <span class="n">h12</span><span class="p">,</span> <span class="n">h22</span>
    <span class="n">x</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">=</span> <span class="n">h11</span> <span class="o">*</span> <span class="n">x</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">+</span> <span class="n">h12</span> <span class="o">*</span> <span class="n">y</span><span class="p">[</span><span class="n">i</span><span class="p">]</span>
    <span class="n">y</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">=</span> <span class="n">h21</span> <span class="o">*</span> <span class="n">x</span><span class="p">[</span><span class="n">i</span><span class="p">]</span> <span class="o">+</span> <span class="n">h22</span> <span class="o">*</span> <span class="n">y</span><span class="p">[</span><span class="n">i</span><span class="p">]</span>
    </pre></div>
    </div>
    <p>Refer to the cuBLAS documentation for the use of <cite>flag</cite>.</p>
    <table class="docutils field-list" frame="void" rules="none">
    <col class="field-name" />
    <col class="field-body" />
    <tbody valign="top">
    <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><ul class="first last simple">
    <li><strong>x</strong> (<em>array</em>) &#8211; vector</li>
    <li><strong>y</strong> (<em>array</em>) &#8211; vector</li>
    </ul>
    </td>
    </tr>
    </tbody>
    </table>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">rotmg</code><span class="sig-paren">(</span><em>d1</em>, <em>d2</em>, <em>x1</em>, <em>y1</em><span class="sig-paren">)</span></dt>
    <dd><p>Constructs the modified Givens transformation <cite>H</cite> that zeros out the second
    entry of a column vector <cite>(d1 * x1, d2 * y1)</cite>.</p>
    <table class="docutils field-list" frame="void" rules="none">
    <col class="field-name" />
    <col class="field-body" />
    <tbody valign="top">
    <tr class="field-odd field"><th class="field-name">Parameters:</th><td class="field-body"><ul class="first simple">
    <li><strong>d1</strong> &#8211; scaling factor for the x-coordinate of the input vector</li>
    <li><strong>d2</strong> &#8211; scaling factor for the y-coordinate of the input vector</li>
    <li><strong>x1</strong> &#8211; x-coordinate of the input vector</li>
    <li><strong>y1</strong> &#8211; y-coordinate of the input vector</li>
    </ul>
    </td>
    </tr>
    <tr class="field-even field"><th class="field-name">Returns:</th><td class="field-body"><p class="first last">A 1D array that is usable in <cite>rotm</cite>.
    The first element is the flag for <cite>rotm</cite>.
    The rest of the elements corresponds to the <cite>h11, h21, h12, h22</cite>
    elements of <cite>H</cite>.</p>
    </td>
    </tr>
    </tbody>
    </table>
    </dd></dl>

    </div>
    <div class="section" id="blas-level-2">
    <h2>BLAS Level 2<a class="headerlink" href="#blas-level-2" title="Permalink to this headline">¶</a></h2>
    <p>All level 2 routines follow the following naming convention for all arguments:</p>
    <ul>
    <li><dl class="first docutils">
    <dt>A, B, C, AP &#8211; (2D array) Matrix argument.</dt>
    <dd><p class="first last"><cite>AP</cite> implies packed storage for banded matrix.</p>
    </dd>
    </dl>
    </li>
    <li><p class="first">x, y, z &#8211; (1D arrays)  Vector argument.</p>
    </li>
    <li><p class="first">alpha, beta &#8211; (scalar) Can be floats or complex numbers depending.</p>
    </li>
    <li><p class="first">m &#8211; (scalar)  Number of rows of matrix <cite>A</cite>.</p>
    </li>
    <li><dl class="first docutils">
    <dt>n &#8211; (scalar)  Number of columns of matrix <cite>A</cite>.  If <cite>m</cite> is not needed,</dt>
    <dd><p class="first last"><cite>n</cite> also means the number of rows of the matrix <cite>A</cite>; thus,
    implying a square matrix.</p>
    </dd>
    </dl>
    </li>
    <li><dl class="first docutils">
    <dt>trans, transa, transb &#8211; (string)</dt>
    <dd><p class="first">Select the operation <cite>op</cite> to apply to a matrix:</p>
    <ul class="simple">
    <li>&#8216;N&#8217;: <cite>op(X) = X</cite>, the identity operation;</li>
    <li>&#8216;T&#8217;: <cite>op(X) = X**T</cite>, the transpose;</li>
    <li>&#8216;C&#8217;: <cite>op(X) = X**H</cite>, the conjugate transpose.</li>
    </ul>
    <p class="last"><cite>trans</cite> only applies to the only matrix argument.
    <cite>transa</cite> and <cite>transb</cite> apply to matrix <cite>A</cite> and matrix <cite>B</cite>,
    respectively.</p>
    </dd>
    </dl>
    </li>
    <li><dl class="first docutils">
    <dt>uplo &#8211; (string) Can be &#8216;U&#8217; for filling the upper trianglar matrix; or &#8216;L&#8217; for</dt>
    <dd><p class="first last">filling the lower trianglar matrix.</p>
    </dd>
    </dl>
    </li>
    <li><p class="first">diag &#8211; (boolean)  Whether the matrix diagonal has unit elements.</p>
    </li>
    <li><dl class="first docutils">
    <dt>mode &#8211; (string) &#8216;L&#8217; means the matrix is on the left side in the equation.</dt>
    <dd><p class="first last">&#8216;R&#8217; means the matrix is on the right side in the equation.</p>
    </dd>
    </dl>
    </li>
    </ul>
    <div class="admonition note">
    <p class="first admonition-title">Note</p>
    <p class="last">The last array argument is always overwritten with the result.</p>
    </div>
    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">gbmv</code><span class="sig-paren">(</span><em>trans</em>, <em>m</em>, <em>n</em>, <em>kl</em>, <em>ku</em>, <em>alpha</em>, <em>A</em>, <em>x</em>, <em>beta</em>, <em>y</em><span class="sig-paren">)</span></dt>
    <dd><p>banded matrix-vector multiplication <cite>y = alpha * op(A) * x + beta * y</cite> where
    <cite>A</cite> has <cite>kl</cite> sub-diagonals and <cite>ku</cite> super-diagonals.</p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">gemv</code><span class="sig-paren">(</span><em>trans</em>, <em>m</em>, <em>n</em>, <em>alpha</em>, <em>A</em>, <em>x</em>, <em>beta</em>, <em>y</em><span class="sig-paren">)</span></dt>
    <dd><p>matrix-vector multiplication <cite>y = alpha * op(A) * x + beta * y</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">trmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>n</em>, <em>A</em>, <em>x</em><span class="sig-paren">)</span></dt>
    <dd><p>triangular matrix-vector multiplication <cite>x = op(A) * x</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">tbmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>n</em>, <em>k</em>, <em>A</em>, <em>x</em><span class="sig-paren">)</span></dt>
    <dd><p>triangular banded matrix-vector <cite>x = op(A) * x</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">tpmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>n</em>, <em>AP</em>, <em>x</em><span class="sig-paren">)</span></dt>
    <dd><p>triangular packed matrix-vector multiplication <cite>x = op(A) * x</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">trsv</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>n</em>, <em>A</em>, <em>x</em><span class="sig-paren">)</span></dt>
    <dd><p>Solves the triangular linear system with a single right-hand-side.
    <cite>op(A) * x = b</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">tpsv</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>n</em>, <em>AP</em>, <em>x</em><span class="sig-paren">)</span></dt>
    <dd><p>Solves the packed triangular linear system with a single right-hand-side.
    <cite>op(A) * x = b</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">tbsv</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>n</em>, <em>k</em>, <em>A</em>, <em>x</em><span class="sig-paren">)</span></dt>
    <dd><p>Solves the triangular banded linear system with a single right-hand-side.
    <cite>op(A) * x = b</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">symv</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>A</em>, <em>x</em>, <em>beta</em>, <em>y</em><span class="sig-paren">)</span></dt>
    <dd><p>symmetric matrix-vector multiplication <cite>y = alpha * A * x + beta * y</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">hemv</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>A</em>, <em>x</em>, <em>beta</em>, <em>y</em><span class="sig-paren">)</span></dt>
    <dd><p>Hermitian matrix-vector multiplication <cite>y = alpha * A * x + beta * y</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">sbmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>k</em>, <em>alpha</em>, <em>A</em>, <em>x</em>, <em>beta</em>, <em>y</em><span class="sig-paren">)</span></dt>
    <dd><p>symmetric banded matrix-vector multiplication  <cite>y = alpha * A * x + beta * y</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">hbmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>k</em>, <em>alpha</em>, <em>A</em>, <em>x</em>, <em>beta</em>, <em>y</em><span class="sig-paren">)</span></dt>
    <dd><p>Hermitian banded matrix-vector multiplication  <cite>y = alpha * A * x + beta * y</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">spmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>AP</em>, <em>x</em>, <em>beta</em>, <em>y</em><span class="sig-paren">)</span></dt>
    <dd><p>symmetric packed matrix-vector multiplication <cite>y = alpha * A * x + beta * y</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">hpmv</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>AP</em>, <em>x</em>, <em>beta</em>, <em>y</em><span class="sig-paren">)</span></dt>
    <dd><p>Hermitian packed matrix-vector multiplication <cite>y = alpha * A * x + beta * y</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">ger</code><span class="sig-paren">(</span><em>m</em>, <em>n</em>, <em>alpha</em>, <em>x</em>, <em>y</em>, <em>A</em><span class="sig-paren">)</span></dt>
    <dd><p>the rank-1 update <cite>A := alpha * x * y ** T + A</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">geru</code><span class="sig-paren">(</span><em>m</em>, <em>n</em>, <em>alpha</em>, <em>x</em>, <em>y</em>, <em>A</em><span class="sig-paren">)</span></dt>
    <dd><p>the rank-1 update <cite>A := alpha * x * y ** T + A</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">gerc</code><span class="sig-paren">(</span><em>m</em>, <em>n</em>, <em>alpha</em>, <em>x</em>, <em>y</em>, <em>A</em><span class="sig-paren">)</span></dt>
    <dd><p>the rank-1 update <cite>A := alpha * x * y ** H + A</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">syr</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x</em>, <em>A</em><span class="sig-paren">)</span></dt>
    <dd><p>symmetric rank 1 operation <cite>A := alpha * x * x ** T + A</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">her</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x</em>, <em>A</em><span class="sig-paren">)</span></dt>
    <dd><p>hermitian rank 1 operation  <cite>A := alpha * x * x ** H + A</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">spr</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x</em>, <em>AP</em><span class="sig-paren">)</span></dt>
    <dd><p>the symmetric rank 1 operation <cite>A := alpha * x * x ** T + A</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">hpr</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x</em>, <em>AP</em><span class="sig-paren">)</span></dt>
    <dd><p>hermitian rank 1 operation <cite>A := alpha * x * x ** H + A</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">syr2</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x</em>, <em>y</em>, <em>A</em><span class="sig-paren">)</span></dt>
    <dd><p>symmetric rank-2 update <cite>A = alpha * x * y ** T + y * x ** T + A</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">her2</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x</em>, <em>y</em>, <em>A</em><span class="sig-paren">)</span></dt>
    <dd><p>Hermitian rank-2 update <cite>A = alpha * x * y ** H + alpha * y * x ** H + A</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">spr2</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x</em>, <em>y</em>, <em>A</em><span class="sig-paren">)</span></dt>
    <dd><p>packed symmetric rank-2 update <cite>A = alpha * x * y ** T + y * x ** T + A</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">hpr2</code><span class="sig-paren">(</span><em>uplo</em>, <em>n</em>, <em>alpha</em>, <em>x</em>, <em>y</em>, <em>A</em><span class="sig-paren">)</span></dt>
    <dd><p>packed Hermitian rank-2 update <cite>A = alpha * x * y ** H + alpha * y * x ** H + A</cite></p>
    </dd></dl>

    </div>
    <div class="section" id="blas-level-3">
    <h2>BLAS Level 3<a class="headerlink" href="#blas-level-3" title="Permalink to this headline">¶</a></h2>
    <p>All level 3 routines follow the same naming convention for arguments as in
    level 2 routines.</p>
    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">gemm</code><span class="sig-paren">(</span><em>transa</em>, <em>transb</em>, <em>m</em>, <em>n</em>, <em>k</em>, <em>alpha</em>, <em>A</em>, <em>B</em>, <em>beta</em>, <em>C</em><span class="sig-paren">)</span></dt>
    <dd><p>matrix-matrix multiplication <cite>C = alpha * op(A) * op(B) + beta * C</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">syrk</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>n</em>, <em>k</em>, <em>alpha</em>, <em>A</em>, <em>beta</em>, <em>C</em><span class="sig-paren">)</span></dt>
    <dd><p>symmetric rank- k update <cite>C = alpha * op(A) * op(A) ** T + beta * C</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">herk</code><span class="sig-paren">(</span><em>uplo</em>, <em>trans</em>, <em>n</em>, <em>k</em>, <em>alpha</em>, <em>A</em>, <em>beta</em>, <em>C</em><span class="sig-paren">)</span></dt>
    <dd><p>Hermitian rank- k update <cite>C = alpha * op(A) * op(A) ** H + beta * C</cite></p>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">symm</code><span class="sig-paren">(</span><em>side</em>, <em>uplo</em>, <em>m</em>, <em>n</em>, <em>alpha</em>, <em>A</em>, <em>B</em>, <em>beta</em>, <em>C</em><span class="sig-paren">)</span></dt>
    <dd><p>symmetric matrix-matrix multiplication:</p>
    <div class="highlight-python"><div class="highlight"><pre><span></span><span class="k">if</span>  <span class="n">side</span> <span class="o">==</span> <span class="s1">&#39;L&#39;</span><span class="p">:</span>
        <span class="n">C</span> <span class="o">=</span> <span class="n">alpha</span> <span class="o">*</span> <span class="n">A</span> <span class="o">*</span> <span class="n">B</span> <span class="o">+</span> <span class="n">beta</span> <span class="o">*</span> <span class="n">C</span>
    <span class="k">else</span><span class="p">:</span>  <span class="c1"># side == &#39;R&#39;</span>
        <span class="n">C</span> <span class="o">=</span> <span class="n">alpha</span> <span class="o">*</span> <span class="n">B</span> <span class="o">*</span> <span class="n">A</span> <span class="o">+</span> <span class="n">beta</span> <span class="o">*</span> <span class="n">C</span>
    </pre></div>
    </div>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">hemm</code><span class="sig-paren">(</span><em>side</em>, <em>uplo</em>, <em>m</em>, <em>n</em>, <em>alpha</em>, <em>A</em>, <em>B</em>, <em>beta</em>, <em>C</em><span class="sig-paren">)</span></dt>
    <dd><p>Hermitian matrix-matrix multiplication:</p>
    <div class="highlight-python"><div class="highlight"><pre><span></span><span class="k">if</span>  <span class="n">side</span> <span class="o">==</span> <span class="s1">&#39;L&#39;</span><span class="p">:</span>
        <span class="n">C</span> <span class="o">=</span> <span class="n">alpha</span> <span class="o">*</span> <span class="n">A</span> <span class="o">*</span> <span class="n">B</span> <span class="o">+</span> <span class="n">beta</span> <span class="o">*</span> <span class="n">C</span>
    <span class="k">else</span><span class="p">:</span>   <span class="c1">#  side == &#39;R&#39;:</span>
        <span class="n">C</span> <span class="o">=</span> <span class="n">alpha</span> <span class="o">*</span> <span class="n">B</span> <span class="o">*</span> <span class="n">A</span> <span class="o">+</span> <span class="n">beta</span> <span class="o">*</span> <span class="n">C</span>
    </pre></div>
    </div>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">trsm</code><span class="sig-paren">(</span><em>side</em>, <em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>m</em>, <em>n</em>, <em>alpha</em>, <em>A</em>, <em>B</em><span class="sig-paren">)</span></dt>
    <dd><p>Solves the triangular linear system with multiple right-hand-sides:</p>
    <div class="highlight-python"><div class="highlight"><pre><span></span><span class="k">if</span>  <span class="n">side</span> <span class="o">==</span> <span class="s1">&#39;L&#39;</span><span class="p">:</span>
        <span class="n">op</span><span class="p">(</span><span class="n">A</span><span class="p">)</span> <span class="o">*</span> <span class="n">X</span> <span class="o">=</span> <span class="n">alpha</span> <span class="o">*</span> <span class="n">B</span>
    <span class="k">else</span><span class="p">:</span>       <span class="c1"># side == &#39;R&#39;</span>
        <span class="n">X</span> <span class="o">*</span> <span class="n">op</span><span class="p">(</span><span class="n">A</span><span class="p">)</span> <span class="o">=</span> <span class="n">alpha</span> <span class="o">*</span> <span class="n">B</span>
    </pre></div>
    </div>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">trmm</code><span class="sig-paren">(</span><em>side</em>, <em>uplo</em>, <em>trans</em>, <em>diag</em>, <em>m</em>, <em>n</em>, <em>alpha</em>, <em>A</em>, <em>B</em>, <em>C</em><span class="sig-paren">)</span></dt>
    <dd><p>triangular matrix-matrix multiplication:</p>
    <div class="highlight-python"><div class="highlight"><pre><span></span>if  side == &#39;:&#39;
        C = alpha * op(A) * B
    else:   # side == &#39;R&#39;
        C = alpha * B * op(A)
    </pre></div>
    </div>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">dgmm</code><span class="sig-paren">(</span><em>side</em>, <em>m</em>, <em>n</em>, <em>A</em>, <em>x</em>, <em>C</em><span class="sig-paren">)</span></dt>
    <dd><p>matrix-matrix multiplication:</p>
    <div class="highlight-python"><div class="highlight"><pre><span></span><span class="k">if</span>  <span class="n">mode</span> <span class="o">==</span> <span class="s1">&#39;R&#39;</span><span class="p">:</span>
        <span class="n">C</span> <span class="o">=</span> <span class="n">A</span> <span class="o">*</span> <span class="n">x</span> <span class="o">*</span> <span class="n">diag</span><span class="p">(</span><span class="n">X</span><span class="p">)</span>
    <span class="k">else</span><span class="p">:</span>       <span class="c1"># mode == &#39;L&#39;</span>
        <span class="n">C</span> <span class="o">=</span> <span class="n">diag</span><span class="p">(</span><span class="n">X</span><span class="p">)</span> <span class="o">*</span> <span class="n">x</span> <span class="o">*</span> <span class="n">A</span>
    </pre></div>
    </div>
    </dd></dl>

    <dl class="method">
    <dt>
    <code class="descclassname">accelerate.cuda.blas.Blas.</code><code class="descname">geam</code><span class="sig-paren">(</span><em>transa</em>, <em>transb</em>, <em>m</em>, <em>n</em>, <em>alpha</em>, <em>A</em>, <em>beta</em>, <em>B</em>, <em>C</em><span class="sig-paren">)</span></dt>
    <dd><p>matrix-matrix addition/transposition <cite>C = alpha * op(A) + beta * op(B)</cite></p>
    </dd></dl>
