1e0ed867bSAlp Dener #include <../src/tao/bound/impls/bqnk/bqnk.h>
2e0ed867bSAlp Dener
TaoSetUp_BQNKTL(Tao tao)3d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_BQNKTL(Tao tao)
4d71ae5a4SJacob Faibussowitsch {
52e6e4ca1SStefano Zampini KSP ksp;
6*0cd8b6e2SStefano Zampini PetscBool valid;
75eb5f4d6SAlp Dener
85eb5f4d6SAlp Dener PetscFunctionBegin;
99566063dSJacob Faibussowitsch PetscCall(TaoSetUp_BQNK(tao));
109566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp));
11*0cd8b6e2SStefano Zampini PetscCall(PetscObjectHasFunction((PetscObject)ksp, "KSPCGSetRadius_C", &valid));
123c859ba3SBarry Smith 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);
133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
145eb5f4d6SAlp Dener }
155eb5f4d6SAlp Dener
163850be85SAlp Dener /*MC
173850be85SAlp Dener TAOBQNKTL - Bounded Quasi-Newton-Krylov Trust-region with Line-search fallback, for nonlinear
183850be85SAlp Dener minimization with bound constraints. This method approximates the Hessian-vector
193850be85SAlp Dener product using a limited-memory quasi-Newton formula, and iteratively inverts the
203850be85SAlp Dener Hessian with a Krylov solver. The quasi-Newton matrix and its settings can be
21a1cb98faSBarry Smith accessed via the prefix `-tao_bqnk_`. For options database, see `TAOBNK`
223850be85SAlp Dener
233850be85SAlp Dener Level: beginner
24a1cb98faSBarry Smith
25a1cb98faSBarry Smith .seealso: `Tao`, `TaoType`, `TAOBNK`, `TAOBQNKTR`, `TAOBQNKLS`
263850be85SAlp Dener M*/
TaoCreate_BQNKTL(Tao tao)27d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BQNKTL(Tao tao)
28d71ae5a4SJacob Faibussowitsch {
29414d97d3SAlp Dener TAO_BNK *bnk;
30414d97d3SAlp Dener TAO_BQNK *bqnk;
31e0ed867bSAlp Dener
32e0ed867bSAlp Dener PetscFunctionBegin;
339566063dSJacob Faibussowitsch PetscCall(TaoCreate_BQNK(tao));
345eb5f4d6SAlp Dener tao->ops->setup = TaoSetUp_BQNKTL;
35414d97d3SAlp Dener bnk = (TAO_BNK *)tao->data;
36414d97d3SAlp Dener bqnk = (TAO_BQNK *)bnk->ctx;
37414d97d3SAlp Dener bqnk->solve = TaoSolve_BNTL;
383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
39e0ed867bSAlp Dener }
40