CAPI trlib_krylov¶
Functions¶
-
trlib_int_t
trlib_krylov_min
(trlib_int_t init, trlib_flt_t radius, trlib_int_t equality, trlib_int_t itmax, trlib_int_t itmax_lanczos, trlib_flt_t tol_rel_i, trlib_flt_t tol_abs_i, trlib_flt_t tol_rel_b, trlib_flt_t tol_abs_b, trlib_flt_t zero, trlib_flt_t obj_lo, trlib_int_t ctl_invariant, trlib_int_t convexify, trlib_int_t earlyterm, trlib_flt_t g_dot_g, trlib_flt_t v_dot_g, trlib_flt_t p_dot_Hp, trlib_int_t *iwork, trlib_flt_t *fwork, trlib_int_t refine, trlib_int_t verbose, trlib_int_t unicode, char *prefix, FILE *fout, trlib_int_t *timing, trlib_int_t *action, trlib_int_t *iter, trlib_int_t *ityp, trlib_flt_t *flt1, trlib_flt_t *flt2, trlib_flt_t *flt3)¶ Solves trust region subproblem: Computes minimizer to
\(\min \frac 12 \langle s, H s \rangle + \langle g_0, s \rangle\) subject to the trust region constraint \(\Vert s \Vert_M \le r\), where
- \(H\) is available via matrix-vector products \(v \mapsto Hv\),
- \(\Vert s \Vert_M = \sqrt{\langle s, Ms \rangle}\) with \(M\) positive definite,
- \(M\) is available via matrix-vector products of its inverse, \(v \mapsto M^{-1}v\).
The minimizer is a global minimizer (modulo floating point), as long as the hard case does not occur. The hard case is characterized by the fact that the eigenspace corresponding to the smallest eigenvalue is degenerate.
In that case a global minimizer in the Krylov subspaces sampled so far is returned, sampling the complete search space can be forced by setting
ctl_invariant
toTRLIB_CLC_EXP_INV_GLO
, but if this is desired factorization based method will most likely be much more efficient.A preconditioned Conjugate Gradient / Lanczos method is used, that possibly employs a reduction to a tridiagnoal subproblem that is solved using Moré-Sorensens method by explicitly factorizing the tridiagonal matrix, calling
trlib_tri_factor_min()
.The method builds upon algorithm 5.1 in [Gould1999], details of the implementation can be found [Lenders2016]:
Convergence
The stopping criterion is based on the gradient of the Lagrangian. Convergence in iteration \(i\) is reported as soon as
- interior case: \(\Vert g_{i+1} \Vert_{M^{-1}} \le \max\{ \texttt{tol}\_\texttt{abs}\_\texttt{i}, \eta_i \, \Vert g_0 \Vert_{M^{-1}} \}\).
- boundary case: \(\Vert g_{i+1} \Vert_{M^{-1}} \le \max\{ \texttt{tol}\_\texttt{abs}\_\texttt{b}, \eta_b \, \Vert g_0 \Vert_{M^{-1}} \}\).
- hard case:
ctl_invariant
determines if this is the sole stopping criterion or if further invariant subspaces are going to be explored. See the documentation of this option.
Here \(\eta_i = \begin{cases} \texttt{tol}\_\texttt{rel}\_\texttt{i} & \texttt{tol}\_\texttt{rel}\_\texttt{i} > 0 \\ \min\{ 0.5, \sqrt{\Vert g_{i+1} \Vert_{M^{-1}}} \} & \texttt{tol}\_\texttt{rel}\_\texttt{i} = -1 \\ \min\{ 0.5, \Vert g_{i+1} \Vert_{M^{-1}} \} & \texttt{tol}\_\texttt{rel}\_\texttt{i} = -2 \end{cases},\) \(\eta_b = \begin{cases} \texttt{tol}\_\texttt{rel}\_\texttt{b} & \texttt{tol}\_\texttt{rel}\_\texttt{b} > 0 \\ \min\{ 0.5, \sqrt{\Vert g_{i+1} \Vert_{M^{-1}}} \} & \texttt{tol}\_\texttt{rel}\_\texttt{b} = -1 \\ \min\{ 0.5, \Vert g_{i+1} \Vert_{M^{-1}} \} & \texttt{tol}\_\texttt{rel}\_\texttt{b} = -2, \\ \max\{10^{-6}, \min\{ 0.5, \sqrt{\Vert g_{i+1} \Vert_{M^{-1}}} \}\} & \texttt{tol}\_\texttt{rel}\_\texttt{b} = -3 \\ \max\{10^{-6}, \min\{ 0.5, \Vert g_{i+1} \Vert_{M^{-1}}\} \} & \texttt{tol}\_\texttt{rel}\_\texttt{b} = -4 \end{cases}.\)
Choice of \(\eta_i\) and \(\eta_b\) depend on the application, the author has good overall experience in unconstrained optimization with \(\texttt{tol}\_\texttt{rel}\_\texttt{i} = -2\) and \(\texttt{tol}\_\texttt{rel}\_\texttt{b} = -3\). Some remarks to keep in mind when choosing the tolerances:
- Choosing fixed small values for \(\eta_i\) and \(\eta_b\), for example \(\texttt{tol}\_\texttt{rel}\_\texttt{i} = 10^{-8}\) and \(\texttt{tol}\_\texttt{rel}\_\texttt{b} = 10^{-6}\) lead to quick convergence in a NLP algorithm with lowest number of function evaluations but at a possible excessive amount of matrix-vector products as this also solves problems accurately far away from the solution.
- Choosing \(\eta \sim O(\sqrt{\Vert g \Vert_{M^{-1}}})\) leads to superlinear convergence in a nonlinear programming algorithm.
- Choosing \(\eta \sim O(\Vert g \Vert_{M^{-1}})\) leads to quadratic convergence in a nonlinear programming algorithm.
- It is questionable whether it makes sense to solve the boundary with very high accuracy, as the final solution of a nonlinear program should be interior.
Calling Scheme
This function employs a reverse communication paradigm. The functions exits whenever there is an action to be performed by the user as indicated by the
action
. The user should perform this action and continue with the other values unchanged as long as the return value is positive.User provided storage
The user has to manage 5/6 vectors refered by the names \(g\), \(g_-\), \(v\), \(s\), \(p\), \(\textit{Hp}\) and a matrix \(Q\) with
itmax
columns to store the Lanczos directions \(p_i\). The user must perform the actions on those as indicated by the return value.In the case of trivial preconditioner, \(M = I\), it always holds \(v = g\) and the vector \(v\) is not necessary.
\(s\) holds the solution candidate.
Note that the action \(s \leftarrow Q h\) will only sometimes be requested in the very final iteration before convergence. If memory and not computational load is an issue, you may very well save
iter
,ityp
,flt1
,flt2
flt2
instead of \(q_j\) and when \(s \leftarrow Q s\) is requested simultaniously recompute the directions \(q_j\) and update the direction \(s\) once \(q_j\) has been computed as \(s \leftarrow h_j q_j\) with initialization \(s \leftarrow 0\).Furthermore the user has to provide a integer and floating point workspace, details are described in
iwork
andfwork
.Resolves
Reentry with new radius
You can efficiently hotstart from old results if you have a new problem with decreased trust region radius. Just set status to
TRLIB_CLS_HOTSTART
. Furthermore hotstarting with increased trust region radius should be trivial as you should be able to just increase the radius and take off the computation from the previous stage. However as this is an untypical scenario, it has not been tested at all.Reentry to get minimizer of regularized problem
You can reenter to compute the solution of a convexified trust region problem. You may be interested to this once you see the effect of bad conditioning or also as a cheap retry in a trust region algorithm before continuing to the next point after a discarded step.
In this case, call the function with status set to
TRLIB_CLS_HOTSTART_P
.Reentry to get unconstrained minimizer of constant regularized problem
After a successful solve you can compute the norm of the unconstrained minimizer to the problem with regularized hessian \(H + \beta M\) (\(\beta\) large enough such that \(H + \beta M\) is positive definite) in the previously expanded Krylov subspace. You will certainly be interested in this if you want to implement the TRACE algorithm described in [Curtis2016].
In this case, call the function with status set to
TRLIB_CLS_HOTSTART_R
andradius
set to \(\beta\).On exit, the
obj
field offwork
yields the norm of the solution vector.If you not only want the norm but also the unconstrained minimizer, you can get it by one step of
TRLIB_CLA_RETRANSF()
. However since this is usually not the case, you are not asked to do this in the reverse communication scheme.Reentry to find suitable TRACE regularization parameter
For the aforementioned TRACE algorithm, there is also the need to compute \(\beta\) such that \(\sigma_{\text l} \le \frac{\beta}{\Vert s(\beta) \Vert} \le \sigma_{\text u}\), where \(\Vert s(\beta) \Vert\) denotes the norm of the regularized unconstrained minimizer.
To get this, set
status
toTRLIB_CLS_HOTSTART_T
andradius
to an initial guess for \(\beta\),tol_rel_i
to \(\sigma_{\text l}\) andtol_rel_b
to \(\sigma_{\text u}\). The fieldobj
offwork
contains the desired solution norm on exit andflt1
is the desired regularization parameter \(\beta\).If you not only want the norm but also the unconstrained minimizer, you can get it by one step of
TRLIB_CLA_RETRANSF
. However since this is usually not the case, you are not asked to do this in the reverse communication scheme.Reentry with new gradient
You can efficiently hotstart from old results if you have a new problem with changed \(g_0\) where the previous sampled Krylov subspace is large enough to contain the solution to the new problem, convergence will not be checked:
- get pointer
gt
to negative gradient of tridiagonal problem in fwork withtrlib_krylov_gt()
- store
gt[0]
and overwrite \(\texttt{gt} \leftarrow - Q_i^T \texttt{grad}\) - set
status
toTRLIB_CLS_HOTSTART_G
and start reverse communication process - to reset data for new problem make sure that you restore
gt[0]
and setgt[1:] = 0
for those elements previously overwritten
Hard case
If you want to investigate the problem also if the hard case, you can either sample through the invariant subspaces (
ctl_invariant
set toTRLIB_CLC_EXP_INV_GLO
) or solve the problem with a gradient for which the hard does not occur and then hotstart with the actual gradient.Parameters: - init (
trlib_int_t
, input) – set toTRLIB_CLS_INIT
on first call, toTRLIB_CLS_HOTSTART
on hotstart with smaller radius and otherwise to0
. - radius (
trlib_flt_t
, input) – trust region radius - equality (
trlib_int_t
, input) – set to1
if trust region constraint should be enforced as equality - itmax (
trlib_int_t
, input) – maximum number of iterations - itmax_lanczos (
trlib_int_t
, input) – maximum number of Lanczos type iterations. You are strongly advised to set this to a small number (say 25) unless you know better. Keep in mind that Lanczos type iteration are only performed when curvature information is flat and Lanczos may amplify rounding errors without reorthogonalization. If you allow a lot of Lanczos type iterations consider reorthogonalizing the new direction against the vector storage. - tol_rel_i (
trlib_flt_t
, input) – relative stopping tolerance for interior solution - tol_abs_i (
trlib_flt_t
, input) – absolute stopping tolerance for interior solution - tol_rel_b (
trlib_flt_t
, input) – relative stopping tolerance for boundary solution - tol_abs_b (
trlib_flt_t
, input) – absolute stopping tolerance for boundary solution - zero (
trlib_flt_t
, input) – threshold that determines when \(\langle p, Hp \rangle\) is considered to be zero and thus control eventual Lanczos switch - obj_lo (
trlib_flt_t
, input) – lower bound on objective, returns if a point is found with function value \(\le \texttt{obj}\_\texttt{lo}\). Set to a positive value to ignore this bound. - ctl_invariant (
trlib_int_t
, input) –- set to
TRLIB_CLC_NO_EXP_INV
if you want to search only in the first invariant Krylov subspace - set to
TRLIB_CLC_EXP_INV_LOC
if you want to continue to expand the Krylov subspaces but terminate if there is convergence indication in the subspaces sampled so far. - set to
TRLIB_CLC_EXP_INV_GLO
if you want to continue to expand the Krylov subspaces even in the case of convergence to a local minimizer within in the subspaces sampled so far. Reverse communication continues as long as ctl_invariant is set toTRLIB_CLC_EXP_INV_GLO
, so you should resetctl_invariant
to eitherTRLIB_CLC_EXP_INV_LOC
orTRLIB_CLC_EXP_INV_LOC
if there is no reason to continue, for example because you cannot find a new nonzero random vector orthogonal to the sampled directions ifaction
isTRLIB_CLA_NEW_KRYLOV
.
- set to
- convexify (
trlib_int_t
, input) – set to 1 if you like to monitor if the tridiagonal solution and the backtransformed solution match and if not resolve with a convexified problem, else set to 0 - earlyterm (
trlib_int_t
, input) – set to 1 if you like to terminate in the boundary case if it unlikely that much progress will be made fast but no convergence is reached, else set to 0 - g_dot_g (
trlib_flt_t
, input) – dot product \(\langle g, g \rangle\) - v_dot_g (
trlib_flt_t
, input) – dot product \(\langle v, g \rangle\) - p_dot_Hp (
trlib_flt_t
, input) – dot product \(\langle p, Hp \rangle\) - iwork (
trlib_int_t
, input/output) –integer workspace for problem, internal memory layout described in table below
- on first call, provide allocated memory for
iwork_size
entries provided bytrlib_krylov_memory_size()
- leave untouched between iterations
start end(exclusive) description 0 1 internal status flag 1 2 iteration counter 2 3 flag that indicates if \(H\) is positive definite in sampled Krylov subspace 3 4 flag that indicates if interior solution is expected 4 5 flag that indicates if warmstart information to leftmost
is available5 6 index to block with smallest leftmost 6 7 flag that indicates if warmstart information to \(\lambda_0\) is available 7 8 flag that indicates if warmstart information to \(\lambda\) is available 8 9 iteration in which switched to Lanczos iteration, -1
: no switch occured9 10 return code from trlib_tri_factor_min()
10 11 sub_fail
exit code from subrotines called intrlib_tri_factor_min()
11 12 number of newton iterations needed in trlib_tri_factor_min()
12 13 last iteration in which a headline has been printed 13 14 kind of iteration headline that has been printed last 14 15 status variable if convexified version is going to be resolved 15 16 number of irreducible blocks 16 17 + itmax
decomposition into irredecubile block, irblk
fortrlib_tri_factor_min()
- on first call, provide allocated memory for
- fwork (
trlib_flt_t
, input/output) –floating point workspace for problem, internal memory layout described in table below
- on very first call, provide allocated memory for at least
fwork_size
entries that has been initialized usingtrlib_prepare_memory()
,fwork_size
can be obtained by callingtrlib_krylov_memory_size()
- can be used for different problem instances with matching dimensions and itmax without reinitialization
- leave untouched between iterations
start end (exclusive) description 0 1 stopping tolerance in case of interior solution 1 2 stopping tolerance in case of boundary solution 2 3 dot product \(\langle v, g \rangle\) in current iteration 3 4 dot product \(\langle p, Hp \rangle\) in current iteration 4 5 ratio between projected CG gradient and Lanczos direction in current iteration 5 6 ratio between projected CG gradient and Lanczos direction in previous iteration 6 7 Lagrange multiplier \(\lambda_0\) for trust region constraint 7 8 Lagrange multiplier \(\lambda\) for trust region constraint 8 9 objective function value in current iterate 9 10 \(\langle s_i, Mp_i \rangle\) 10 11 \(\langle p_i, Mp_i \rangle\) 11 12 \(\langle s_i, Ms_i \rangle\) 12 13 \(\sigma_i\) 13 14 max Rayleigh quotient, \(\max_i \frac{\langle p, Hp \rangle}{\langle p, M p \rangle}\) 14 15 min Rayleigh quotient, \(\min_i \frac{\langle p, Hp \rangle}{\langle p, M p \rangle}\) 15 16 + itmax
\(\alpha_i, i \ge 0\), step length CG 16 + itmax
17 + 2 itmax
\(\beta_i, i \ge 0\), step update factor CG 17 + 2 itmax
18 + 3 itmax
neglin
fortrlib_tri_factor_min()
, just given by \(- \gamma_0 e_1\)18 + 3 itmax
19 + 4 itmax
solution \(h_0\) of tridiagonal subproblem provided as sol
bytrlib_tri_factor_min()
19 + 4 itmax
20 + 5 itmax
solution \(h\) of tridiagonal subproblem provided as sol
bytrlib_tri_factor_min()
20 + 5 itmax
21 + 6 itmax
\(\delta_i, i \ge 0\), curvature in Lanczos, diagonal of \(T\) in Lanczos tridiagonalization process 21 + 6 itmax
22 + 7 itmax
diagonal of Cholesky of \(T_0 + \lambda_0 I\) 22 + 7 itmax
23 + 8 itmax
diagonal of Cholesky of \(T + \lambda I\) 23 + 8 itmax
23 + 9 itmax
\(\gamma_i, i \ge 1\), norm of gradients in Lanczos; provides offdiagonal of \(T\) in Lanczos tridiagonalization process 23 + 9 itmax
23 + 10 itmax
offdiagonal of Cholesky factorization of \(T_0 + \lambda_0 I\) 23 + 10 itmax
23 + 11 itmax
offdiagonal of Cholesky factorization of \(T + \lambda I\) 23 + 11 itmax
24 + 12 itmax
ones
fortrlib_tri_factor_min()
andtrlib_eigen_inverse()
24 + 12 itmax
25 + 13 itmax
leftmost
fortrlib_tri_factor_min()
25 + 13 itmax
26 + 14 itmax
regdiag
fortrlib_tri_factor_regularize_posdef()
26 + 14 itmax
27 + 15 itmax
history of convergence criteria values 27 + 15 itmax
27 + 15 itmax
+fmz
fwork
fortrlib_tri_factor_min()
,fmz
is given bytrlib_tri_factor_memory_size()
with argumentitmax+1
- on very first call, provide allocated memory for at least
- refine (
trlib_int_t
, input) – set to1
if iterative refinement should be used on solving linear systems, otherwise to0
- 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, all values are multiples of nanoseconds, provide zero allocated memory of length
trlib_krylov_timing_size()
block description 0 total duration 1 timing
oftrlib_tri_factor_min()
- action (
trlib_int_t
, output) –The user should perform the following action depending on
action
andityp
on the vectors he manages, see the table below. The table makes use of the notation explained in the User provided storage above and the following:- \(i\):
iter
- \(q_j\): \(j\)-th column of \(Q\)
- \(Q_i\): matrix consisting of the first \(i+1\) columns of \(Q\), \(Q_i = (q_0, \ldots, q_i)\)
- \(h_i\): vector of length \(i+1\) stored in
fwork
at start positionh_pointer
provided bytrlib_krylov_memory_size()
- \(p \leftarrow \perp_M Q_j\): optionally \(M\)-reorthonormalize \(p\) against \(Q_j\)
- \(g \leftarrow \texttt{rand}\perp Q_j\) find nonzero random vector that is orthogonal to \(Q_j\)
- Note that
TRLIB_CLA_NEW_KRYLOV
is unlikely and only occurs on problems that employ the hard case and only ifctl_invariant
\(\neq\)TRLIB_CLC_NO_EXP_INV
. If you want to safe yourself from the trouble implementing this and are confident that you don’t need to expand several invariant Krylov subspaces, just ensure thatctl_invariant
is set toTRLIB_CLC_NO_EXP_INV
.
action ityp command TRLIB_CLA_TRIVIAL
TRLIB_CLT_CG
,TRLIB_CLT_L
do nothing TRLIB_CLA_RETRANSF
TRLIB_CLT_CG
,TRLIB_CLT_L
\(s \leftarrow Q_i h_i\) TRLIB_CLA_INIT
TRLIB_CLT_CG
,TRLIB_CLT_L
\(s \leftarrow 0\), \(g \leftarrow g_0\), \(g_- \leftarrow 0\), \(v \leftarrow M^{-1} g\), \(p \leftarrow -v\), \(\textit{Hp} \leftarrow Hp\), \(\texttt{g}\_\texttt{dot}\_\texttt{g} \leftarrow \langle g, g \rangle\), \(\texttt{v}\_\texttt{dot}\_\texttt{g} \leftarrow \langle v, g \rangle\), \(\texttt{p}\_\texttt{dot}\_\texttt{Hp} \leftarrow \langle p, \textit{Hp} \rangle\), \(q_0 \leftarrow \frac{1}{\sqrt{\texttt{v}\_\texttt{dot}\_\texttt{g}}} v\) TRLIB_CLA_UPDATE_STATIO
TRLIB_CLT_CG
\(s \leftarrow s + \texttt{flt1} \, p\) TRLIB_CLA_UPDATE_STATIO
TRLIB_CLT_L
do nothing TRLIB_CLA_UPDATE_GRAD
TRLIB_CLT_CG
\(q_i \leftarrow \texttt{flt2} \, v\), \(g_- \leftarrow g\), \(g \leftarrow g + \texttt{flt1} \, \textit{Hp}\), \(v \leftarrow M^{-1} g\), \(\texttt{g}\_\texttt{dot}\_\texttt{g} \leftarrow \langle g, g \rangle\), \(\texttt{v}\_\texttt{dot}\_\texttt{g} \leftarrow \langle v, g \rangle\) TRLIB_CLA_UPDATE_GRAD
TRLIB_CLT_L
\(s \leftarrow \textit{Hp} + \texttt{flt1}\, g + \texttt{flt2}\, g_-\), \(g_- \leftarrow \texttt{flt3}\, g\), \(g \leftarrow s\), \(v \leftarrow M^{-1} g\), \(\texttt{v}\_\texttt{dot}\_\texttt{g} \leftarrow \langle v, g \rangle\) TRLIB_CLA_UPDATE_DIR
TRLIB_CLT_CG
\(p \leftarrow \texttt{flt1} \, v + \texttt{flt2} \, p\) with \(\texttt{flt1} = -1\), \(\textit{Hp} \leftarrow Hp\), \(\texttt{p}\_\texttt{dot}\_\texttt{Hp} \leftarrow \langle p, \textit{Hp} \rangle\) TRLIB_CLA_UPDATE_DIR
TRLIB_CLT_L
\(p \leftarrow \texttt{flt1} \, v + \texttt{flt2} \, p\) with \(\texttt{flt2} = 0\), \(p \leftarrow \perp_M Q_{i-1}\), \(\textit{Hp} \leftarrow Hp\), \(\texttt{p}\_\texttt{dot}\_\texttt{Hp} \leftarrow \langle p, \textit{Hp} \rangle\), \(q_i \leftarrow p\) TRLIB_CLA_NEW_KRYLOV
TRLIB_CLT_CG
,TRLIB_CLT_L
\(g \leftarrow \texttt{rand}\perp Q_{i-1}\), \(g_- \leftarrow 0\), \(v \leftarrow M^{-1} g\), \(\texttt{v}\_\texttt{dot}\_\texttt{g} \leftarrow \langle v, g \rangle\), \(p \leftarrow \frac{1}{\sqrt{\texttt{v}\_\texttt{dot}\_\texttt{g}}} v\), \(\texttt{p}\_\texttt{dot}\_\texttt{Hp} \leftarrow \langle p, \textit{Hp} \rangle\), \(q_{i+1} \leftarrow p\) TRLIB_CLA_CONV_HARD
TRLIB_CLT_CG
,TRLIB_CLT_L
\(\texttt{v}\_\texttt{dot}\_\texttt{g} \leftarrow \langle Hs+g_0+\texttt{flt1}\, Ms, M^{-1}(Hs+g_0) + \texttt{flt1} s \rangle\) TRLIB_CLA_OBJVAL
TRLIB_CLT_CG
,TRLIB_CLT_L
\(\texttt{g}\_\texttt{dot}\_\texttt{g} \leftarrow \tfrac 12 \langle s, Hs \rangle + \langle g, s \rangle\) - \(i\):
- iter (
trlib_int_t
, output) – iteration counter to tell user position in vector storage - ityp (
trlib_int_t
, output) – iteration type, seeaction
- flt1 (
trlib_flt_t
, output) – floating point value that user needs for actions - flt2 (
trlib_flt_t
, output) – floating point value that user needs for actions - flt3 (
trlib_flt_t
, output) – floating point value that user needs for actions
Returns: status flag with following meaning:
TRLIB_CLR_CONTINUE
no convergence yet, continue in reverse communicationTRLIB_CLR_CONV_BOUND
successful exit with converged solution on boundary, end reverse communication processTRLIB_CLR_CONV_INTERIOR
successful exit with converged interior solution, end reverse communication processTRLIB_CLR_APPROX_HARD
succesful exit with approximate solution, hard case occured, end reverse communication processTRLIB_CLR_NEWTON_BREAK
exit with breakdown in Newton iteration intrlib_tri_factor_min()
, most likely converged to boundary solutionTRLIB_TTR_HARD_INIT_LAM
hard case encountered without being able to find suitable initial \(\lambda\) for Newton iteration, returned approximate stationary point that maybe suboptimalTRLIB_CLR_ITMAX
iteration limit exceeded, end reverse communication processTRLIB_CLR_UNBDBEL
problem seems to be unbounded from below, end reverse communication processTRLIB_CLR_FAIL_FACTOR
failure on factorization, end reverse communication processTRLIB_CLR_FAIL_LINSOLVE
failure on backsolve, end reverse communication processTRLIB_CLR_FAIL_NUMERIC
failure as one of the valuesv_dot_g
orp_dot_Hp
are not a floating point numberTRLIB_CLR_UNLIKE_CONV
exit early as convergence seems to be unlikelyTRLIB_CLR_PCINDEF
preconditioner apperas to be indefinite, end reverse communication processTRLIB_CLR_UNEXPECT_INT
unexpected interior solution found, expected boundary solution, end reverse communication processTRLIB_CLR_FAIL_TTR
failure occured intrlib_tri_factor_min()
, checkiwork[7]
andiwork[8]
for detailsTRLIB_CLR_FAIL_HARD
failure due to occurence of hard case: invariant subspace encountered without local convergence and request for early termination without exploring further invariant subspaces
Return type: trlib_int_t
-
trlib_int_t
trlib_krylov_prepare_memory
(trlib_int_t itmax, trlib_flt_t *fwork)¶ Prepares floating point workspace for :c:func::trlib_krylov_min
Initializes floating point workspace
fwork
fortrlib_krylov_min()
Parameters: - itmax (
trlib_int_t
, input) – maximum number of iterations - fwork (
trlib_flt_t
, input/output) – floating point workspace to be used bytrlib_krylov_min()
Returns: 0
Return type: trlib_int_t
- itmax (
-
trlib_int_t
trlib_krylov_memory_size
(trlib_int_t itmax, trlib_int_t *iwork_size, trlib_int_t *fwork_size, trlib_int_t *h_pointer)¶ Gives information on memory that has to be allocated for :c:func::trlib_krylov_min
Parameters: - itmax (
trlib_int_t
, input) – maximum number of iterations - iwork_size (
trlib_int_t
, output) – size of integer workspace iwork that has to be allocated fortrlib_krylov_min()
- fwork_size (
trlib_int_t
, output) – size of floating point workspace fwork that has to be allocated fortrlib_krylov_min()
- h_pointer (
trlib_int_t
, output) – start index of vector \(h\) that has to be used in reverse communication in actionTRLIB_CLA_RETRANSF
Returns: 0
Return type: trlib_int_t
- itmax (
-
trlib_int_t
trlib_krylov_timing_size
(void)¶ size that has to be allocated for
timing
intrlib_krylov_min()
Returns: 0
Return type: trlib_int_t
-
trlib_int_t
trlib_krylov_gt
(trlib_int_t itmax, trlib_int_t *gt_pointer)¶ Gives pointer to negative gradient of tridiagonal problem
Parameters: - itmax (
trlib_int_t
, input) – itmax maximum number of iterations - gt_pointer (
trlib_int_t
, output) – pointer to negative gradient of tridiagonal subproblem
Returns: 0
Return type: trlib_int_t
- itmax (
Definitions¶
-
TRLIB_CLR_CONV_BOUND
¶ 0
-
TRLIB_CLR_CONV_INTERIOR
¶ 1
-
TRLIB_CLR_APPROX_HARD
¶ 2
-
TRLIB_CLR_NEWTON_BREAK
¶ 3
-
TRLIB_CLR_HARD_INIT_LAM
¶ 4
-
TRLIB_CLR_FAIL_HARD
¶ 5
-
TRLIB_CLR_UNBDBEL
¶ 6
-
TRLIB_CLR_UNLIKE_CONV
¶ 7
-
TRLIB_CLR_CONTINUE
¶ 10
-
TRLIB_CLR_ITMAX
¶ -1
-
TRLIB_CLR_FAIL_FACTOR
¶ -3
-
TRLIB_CLR_FAIL_LINSOLVE
¶ -4
-
TRLIB_CLR_FAIL_NUMERIC
¶ -5
-
TRLIB_CLR_FAIL_TTR
¶ -7
-
TRLIB_CLR_PCINDEF
¶ -8
-
TRLIB_CLR_UNEXPECT_INT
¶ -9
-
TRLIB_CLT_CG
¶ 1
-
TRLIB_CLT_L
¶ 2
-
TRLIB_CLA_TRIVIAL
¶ 0
-
TRLIB_CLA_INIT
¶ 1
-
TRLIB_CLA_RETRANSF
¶ 2
-
TRLIB_CLA_UPDATE_STATIO
¶ 3
-
TRLIB_CLA_UPDATE_GRAD
¶ 4
-
TRLIB_CLA_UPDATE_DIR
¶ 5
-
TRLIB_CLA_NEW_KRYLOV
¶ 6
-
TRLIB_CLA_CONV_HARD
¶ 7
-
TRLIB_CLA_OBJVAL
¶ 8
-
TRLIB_CLS_INIT
¶ 1
-
TRLIB_CLS_HOTSTART
¶ 2
-
TRLIB_CLS_HOTSTART_G
¶ 3
-
TRLIB_CLS_HOTSTART_P
¶ 4
-
TRLIB_CLS_HOTSTART_R
¶ 5
-
TRLIB_CLS_HOTSTART_T
¶ 6
-
TRLIB_CLS_VEC_INIT
¶ 7
-
TRLIB_CLS_CG_NEW_ITER
¶ 8
-
TRLIB_CLS_CG_UPDATE_S
¶ 9
-
TRLIB_CLS_CG_UPDATE_GV
¶ 10
-
TRLIB_CLS_CG_UPDATE_P
¶ 11
-
TRLIB_CLS_LANCZOS_SWT
¶ 12
-
TRLIB_CLS_L_UPDATE_P
¶ 13
-
TRLIB_CLS_L_CMP_CONV
¶ 14
-
TRLIB_CLS_L_CMP_CONV_RT
¶ 15
-
TRLIB_CLS_L_CHK_CONV
¶ 16
-
TRLIB_CLS_L_NEW_ITER
¶ 17
-
TRLIB_CLS_CG_IF_IRBLK_P
¶ 18
-
TRLIB_CLS_CG_IF_IRBLK_C
¶ 19
-
TRLIB_CLS_CG_IF_IRBLK_N
¶ 20
-
TRLIB_CLC_NO_EXP_INV
¶ 0
-
TRLIB_CLC_EXP_INV_LOC
¶ 1
-
TRLIB_CLC_EXP_INV_GLO
¶ 2
-
TRLIB_CLT_CG_INT
¶ 0
-
TRLIB_CLT_CG_BOUND
¶ 1
-
TRLIB_CLT_LANCZOS
¶ 2
-
TRLIB_CLT_HOTSTART
¶ 3