Actual source code: bqnktl.c
1: #include <../src/tao/bound/impls/bqnk/bqnk.h>
3: static PetscErrorCode TaoSetUp_BQNKTL(Tao tao)
4: {
5: KSP ksp;
6: PetscVoidFn *valid;
8: PetscFunctionBegin;
9: PetscCall(TaoSetUp_BQNK(tao));
10: PetscCall(TaoGetKSP(tao, &ksp));
11: PetscCall(PetscObjectQueryFunction((PetscObject)ksp, "KSPCGSetRadius_C", &valid));
12: PetscCheck(valid, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Not for KSP type %s. Must use a trust-region CG method for KSP (e.g. KSPNASH, KSPSTCG, KSPGLTR)", ((PetscObject)ksp)->type_name);
13: PetscFunctionReturn(PETSC_SUCCESS);
14: }
16: /*MC
17: TAOBQNKTL - Bounded Quasi-Newton-Krylov Trust-region with Line-search fallback, for nonlinear
18: minimization with bound constraints. This method approximates the Hessian-vector
19: product using a limited-memory quasi-Newton formula, and iteratively inverts the
20: Hessian with a Krylov solver. The quasi-Newton matrix and its settings can be
21: accessed via the prefix `-tao_bqnk_`. For options database, see `TAOBNK`
23: Level: beginner
25: .seealso: `Tao`, `TaoType`, `TAOBNK`, `TAOBQNKTR`, `TAOBQNKLS`
26: M*/
27: PETSC_EXTERN PetscErrorCode TaoCreate_BQNKTL(Tao tao)
28: {
29: TAO_BNK *bnk;
30: TAO_BQNK *bqnk;
32: PetscFunctionBegin;
33: PetscCall(TaoCreate_BQNK(tao));
34: tao->ops->setup = TaoSetUp_BQNKTL;
35: bnk = (TAO_BNK *)tao->data;
36: bqnk = (TAO_BQNK *)bnk->ctx;
37: bqnk->solve = TaoSolve_BNTL;
38: PetscFunctionReturn(PETSC_SUCCESS);
39: }