CAPI trlib_leftmost¶
Functions¶
-
trlib_int_t
trlib_leftmost
(trlib_int_t nirblk, trlib_int_t *irblk, trlib_flt_t *diag, trlib_flt_t *offdiag, trlib_int_t warm, trlib_flt_t leftmost_minor, trlib_int_t itmax, trlib_flt_t tol_abs, trlib_int_t verbose, trlib_int_t unicode, char *prefix, FILE *fout, trlib_int_t *timing, trlib_int_t *ileftmost, trlib_flt_t *leftmost)¶ Computes smallest eigenvalue of symmetric tridiagonal matrix \(T \in \mathbb R^{n\times n}\), using a iteration based on last-pivot function of Parlett and Reid.
Let \(T = \begin{pmatrix} T_1 & & \\ & \ddots & \\ & & T_\ell \end{pmatrix}\) be composed into irreducible blocks \(T_i\).
Calls
trlib_leftmost_irreducible()
on every irreducible block in case of coldstart, in case of warmstart just updates information on \(T_\ell\).Parameters: - nirblk (
trlib_int_t
, input) – number of irreducible blocks \(\ell\), ensure \(\ell > 0\) - irblk (
trlib_int_t
, input) –pointer to indices of irreducible blocks, length
nirblk+1
:irblk[i]
is the start index of block \(i\) indiag
andoffdiag
irblk[i+1] - 1
is the stop index of block \(i\)irblk[i+1] - irred[i]
the dimension \(n_\ell\) of block \(i\)irblk[nirred]
the dimension of \(T\)
- diag (
trlib_flt_t
, input) – pointer to array holding diagonal of \(T\), lengthirblk[nirblk]
- offdiag (
trlib_flt_t
, input) – pointer to array holding offdiagonal of \(T\), lengthirblk[nirblk]
- warm (
trlib_int_t
, input) – set \(\ge 1\) if you want to update information about block \(\ell\), provide values inleftmost_minor
,ileftmost
,leftmost
; else0
- leftmost_minor (
trlib_flt_t
, input) – smallest eigenvalue of principal \((n_\ell-1)\times (n_\ell-1)\) submatrix of \(T_\ell\) - itmax (
trlib_int_t
, input) – maximum number of iterations - tol_abs (
trlib_flt_t
, input) – absolute stopping tolerance in Reid-Parlett zero finding, good default may be \(\sqrt{\texttt{macheps}}^{3/4}\) (TRLIB_EPS_POW_75
) - verbose (
trlib_int_t
, input) – determines the verbosity level of output that is written tofout
- unicode (
trlib_int_t
, input) – set to1
iffout
can handle unicode, otherwise to0
- prefix (char, input) – string that is printed before iteration output
- fout (FILE, input) – output stream
- timing (
trlib_int_t
, input/output) –gives timing details, provide allocated zero initialized memory of length
trlib_leftmost_timing_size()
block description 0 total duration - ileftmost (
trlib_int_t
, input/output) – index of block that corresponds to absolute smallest eigenvalue - leftmost (
trlib_flt_t
, input/output) –smallest eigenvalue of \(T\), length \(\ell\)
- on entry: allocated memory
- on exit:
leftmost[i]
smallest eigenvalue of \(T_i\)
Returns: status
TRLIB_LMR_CONV
successTRLIB_LMR_ITMAX
iteration limit exceeded
Return type: trlib_int_t
- nirblk (
-
trlib_int_t
trlib_leftmost_irreducible
(trlib_int_t n, trlib_flt_t *diag, trlib_flt_t *offdiag, trlib_int_t warm, trlib_flt_t leftmost_minor, trlib_int_t itmax, trlib_flt_t tol_abs, trlib_int_t verbose, trlib_int_t unicode, char *prefix, FILE *fout, trlib_int_t *timing, trlib_flt_t *leftmost, trlib_int_t *iter_pr)¶ Computes smallest eigenvalue of irreducible symmetric tridiagonal matrix \(T \in \mathbb R^{n\times n}\), using a iteration based on last-pivot function of Parlett and Reid.
Method is sketched on p. 516 in [Gould1999].
Note that this function most likely will fail in the case of a reducible matrix (
offdiag
contains 0).Convergence
Convergence is reported if \(\texttt{up}-\texttt{low} \le \texttt{tol}\_\texttt{abs} * \max\{1, \vert \texttt{low} \vert, \vert \texttt{up} \vert \}\) or \(\texttt{prlp} \le \texttt{tol}\_\texttt{abs}\), \(\texttt{low}\) and \(\texttt{up}\) denote bracket values enclosing the leftmost eigenvalue and \(\texttt{prlp}\) denotes the last-pivot function value used in root finding.
Parameters: - n (
trlib_int_t
, input) – dimension, ensure \(n > 0\) - diag (
trlib_flt_t
, input) – pointer to array holding diagonal of \(T\), lengthn
- offdiag (
trlib_flt_t
, input) – pointer to array holding offdiagonal of \(T\), lengthn-1
- warm (
trlib_int_t
, input) –set \(\ge 1\) if you provide a valid value in
leftmost_minor
, else0
. Exact value determines which model will be used for zero-findingwarm model 0 linear model for rational function 1 sensible heuristic model choice for lifted rational function 2 asymptotic quadratic model \(\theta^2 + b \theta + c\) for lifted rational function 3 taylor quadratic model \(a \theta^2 + b \theta + c\) for lifted rational function 4 linear model \(b \theta + c\) for lifted rational function - leftmost_minor (
trlib_flt_t
, input) – smallest eigenvalue of principal \((n-1)\times (n-1)\) submatrix of \(T\) - itmax (
trlib_int_t
, input) – maximum number of iterations - tol_abs (
trlib_flt_t
, input) – absolute stopping tolerance in Reid-Parlett zero finding, good default may be \(\sqrt{\texttt{macheps}}^{3/4}\) (TRLIB_EPS_POW_75
) - verbose (
trlib_int_t
, input) – determines the verbosity level of output that is written tofout
- unicode (
trlib_int_t
, input) – set to1
iffout
can handle unicode, otherwise to0
- prefix (char, input) – string that is printed before iteration output
- fout (FILE, input) – output stream
- timing (
trlib_int_t
, input/output) –gives timing details, provide allocated zero initialized memory of length
trlib_leftmost_timing_size()
block description 0 total duration - leftmost (
trlib_flt_t
, output) – smallest eigenvalue of \(T\) - iter_pr (
trlib_int_t
, output) – number of Parlett-Reid iterations
Returns: status
TRLIB_LMR_CONV
successTRLIB_LMR_ITMAX
iteration limit exceeded
Return type: trlib_int_t
- n (
-
trlib_int_t
trlib_leftmost_timing_size
(void)¶ size that has to be allocated for
timing
intrlib_leftmost_irreducible()
andtrlib_leftmost()