xref: /petsc/src/tao/bound/impls/bqnk/bqnktr.c (revision 5fedec97950c19de564efaecd0f125b1a6cb2b20)
1 #include <../src/tao/bound/impls/bqnk/bqnk.h>
2 #include <petscksp.h>
3 
4 static PetscErrorCode TaoSetUp_BQNKTR(Tao tao)
5 {
6   TAO_BNK         *bnk = (TAO_BNK*)tao->data;
7   PetscErrorCode ierr;
8 
9   PetscFunctionBegin;
10   ierr = TaoSetUp_BQNK(tao);CHKERRQ(ierr);
11   if (!bnk->is_nash && !bnk->is_stcg && !bnk->is_gltr) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Must use a trust-region CG method for KSP (KSPNASH, KSPSTCG, KSPGLTR)");
12   PetscFunctionReturn(0);
13 }
14 
15 /*MC
16   TAOBQNKTR - Bounded Quasi-Newton-Krylov Trust Region method for nonlinear minimization with
17               bound constraints. This method approximates the Hessian-vector product using a
18               limited-memory quasi-Newton formula, and iteratively inverts the Hessian with a
19               Krylov solver. The quasi-Newton matrix and its settings can be accessed via the
20               prefix `-tao_bqnk_`. For options database, see TAOBNK
21 
22   Level: beginner
23 .seealso TAOBNK, TAOBQNKTR, TAOBQNKLS
24 M*/
25 PETSC_EXTERN PetscErrorCode TaoCreate_BQNKTR(Tao tao)
26 {
27   TAO_BNK        *bnk;
28   TAO_BQNK       *bqnk;
29   PetscErrorCode ierr;
30 
31   PetscFunctionBegin;
32   ierr = TaoCreate_BQNK(tao);CHKERRQ(ierr);
33   tao->ops->setup = TaoSetUp_BQNKTR;
34   bnk = (TAO_BNK*)tao->data;
35   bqnk = (TAO_BQNK*)bnk->ctx;
36   bqnk->solve = TaoSolve_BNTR;
37   PetscFunctionReturn(0);
38 }
39