Actual source code: bqnktr.c
1: #include <../src/tao/bound/impls/bqnk/bqnk.h>
2: #include <petscksp.h>
4: static PetscErrorCode TaoSetUp_BQNKTR(Tao tao)
5: {
6: KSP ksp;
7: PetscBool valid;
9: PetscFunctionBegin;
10: PetscCall(TaoSetUp_BQNK(tao));
11: PetscCall(TaoGetKSP(tao, &ksp));
12: PetscCall(PetscObjectHasFunction((PetscObject)ksp, "KSPCGSetRadius_C", &valid));
13: 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);
14: PetscFunctionReturn(PETSC_SUCCESS);
15: }
17: /*MC
18: TAOBQNKTR - Bounded Quasi-Newton-Krylov Trust Region method for nonlinear minimization with
19: bound constraints. This method approximates the Hessian-vector product using a
20: limited-memory quasi-Newton formula, and iteratively inverts the Hessian with a
21: Krylov solver. The quasi-Newton matrix and its settings can be accessed via the
22: prefix `-tao_bqnk_`. For options database, see TAOBNK
24: Level: beginner
26: Notes:
27: The base class for this method is `TAOBNK`
29: The various algorithmic factors can only be supplied via the options database
31: .seealso: `Tao`, `TAOBNK`, `TAONLS`, `TAONTL`, `TAONM`, `TaoType`, `TaoCreate()`, `TAOBQNKTR`, `TAOBQNKLS`
32: M*/
33: PETSC_EXTERN PetscErrorCode TaoCreate_BQNKTR(Tao tao)
34: {
35: TAO_BNK *bnk;
36: TAO_BQNK *bqnk;
38: PetscFunctionBegin;
39: PetscCall(TaoCreate_BQNK(tao));
40: tao->ops->setup = TaoSetUp_BQNKTR;
41: bnk = (TAO_BNK *)tao->data;
42: bqnk = (TAO_BQNK *)bnk->ctx;
43: bqnk->solve = TaoSolve_BNTR;
44: PetscFunctionReturn(PETSC_SUCCESS);
45: }