xref: /petsc/src/tao/bound/impls/bqnk/bqnktr.c (revision 3c859ba3a04a72e8efcb87bc7ffc046a6cbab413)
1e0ed867bSAlp Dener #include <../src/tao/bound/impls/bqnk/bqnk.h>
2e0ed867bSAlp Dener #include <petscksp.h>
35eb5f4d6SAlp Dener 
45eb5f4d6SAlp Dener static PetscErrorCode TaoSetUp_BQNKTR(Tao tao)
55eb5f4d6SAlp Dener {
65eb5f4d6SAlp Dener   PetscErrorCode    ierr;
72e6e4ca1SStefano Zampini   KSP               ksp;
82e6e4ca1SStefano Zampini   PetscVoidFunction valid;
95eb5f4d6SAlp Dener 
105eb5f4d6SAlp Dener   PetscFunctionBegin;
115eb5f4d6SAlp Dener   ierr = TaoSetUp_BQNK(tao);CHKERRQ(ierr);
122e6e4ca1SStefano Zampini   ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr);
132e6e4ca1SStefano Zampini   ierr = PetscObjectQueryFunction((PetscObject)ksp,"KSPCGSetRadius_C",&valid);CHKERRQ(ierr);
14*3c859ba3SBarry 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);
155eb5f4d6SAlp Dener   PetscFunctionReturn(0);
165eb5f4d6SAlp Dener }
175eb5f4d6SAlp Dener 
183850be85SAlp Dener /*MC
193850be85SAlp Dener   TAOBQNKTR - Bounded Quasi-Newton-Krylov Trust Region method for nonlinear minimization with
203850be85SAlp Dener               bound constraints. This method approximates the Hessian-vector product using a
213850be85SAlp Dener               limited-memory quasi-Newton formula, and iteratively inverts the Hessian with a
223850be85SAlp Dener               Krylov solver. The quasi-Newton matrix and its settings can be accessed via the
239fa2b5dcSStefano Zampini               prefix `-tao_bqnk_`. For options database, see TAOBNK
243850be85SAlp Dener 
253850be85SAlp Dener   Level: beginner
269fa2b5dcSStefano Zampini .seealso TAOBNK, TAOBQNKTR, TAOBQNKLS
273850be85SAlp Dener M*/
28e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BQNKTR(Tao tao)
29e0ed867bSAlp Dener {
30414d97d3SAlp Dener   TAO_BNK        *bnk;
31414d97d3SAlp Dener   TAO_BQNK       *bqnk;
32e0ed867bSAlp Dener   PetscErrorCode ierr;
33e0ed867bSAlp Dener 
34e0ed867bSAlp Dener   PetscFunctionBegin;
35e0ed867bSAlp Dener   ierr = TaoCreate_BQNK(tao);CHKERRQ(ierr);
365eb5f4d6SAlp Dener   tao->ops->setup = TaoSetUp_BQNKTR;
37414d97d3SAlp Dener   bnk = (TAO_BNK*)tao->data;
38414d97d3SAlp Dener   bqnk = (TAO_BQNK*)bnk->ctx;
39414d97d3SAlp Dener   bqnk->solve = TaoSolve_BNTR;
40e0ed867bSAlp Dener   PetscFunctionReturn(0);
41e0ed867bSAlp Dener }
42