1e0ed867bSAlp Dener #include <../src/tao/bound/impls/bqnk/bqnk.h> 2e0ed867bSAlp Dener 35eb5f4d6SAlp Dener static PetscErrorCode TaoSetUp_BQNKTL(Tao tao) 45eb5f4d6SAlp Dener { 55eb5f4d6SAlp Dener PetscErrorCode ierr; 62e6e4ca1SStefano Zampini KSP ksp; 72e6e4ca1SStefano Zampini PetscVoidFunction valid; 85eb5f4d6SAlp Dener 95eb5f4d6SAlp Dener PetscFunctionBegin; 105eb5f4d6SAlp Dener ierr = TaoSetUp_BQNK(tao);CHKERRQ(ierr); 112e6e4ca1SStefano Zampini ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr); 122e6e4ca1SStefano Zampini ierr = PetscObjectQueryFunction((PetscObject)ksp,"KSPCGSetRadius_C",&valid);CHKERRQ(ierr); 13*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); 145eb5f4d6SAlp Dener PetscFunctionReturn(0); 155eb5f4d6SAlp Dener } 165eb5f4d6SAlp Dener 173850be85SAlp Dener /*MC 183850be85SAlp Dener TAOBQNKTL - Bounded Quasi-Newton-Krylov Trust-region with Line-search fallback, for nonlinear 193850be85SAlp Dener minimization with bound constraints. This method approximates the Hessian-vector 203850be85SAlp Dener product using a limited-memory quasi-Newton formula, and iteratively inverts the 213850be85SAlp Dener Hessian with a Krylov solver. The quasi-Newton matrix and its settings can be 229fa2b5dcSStefano Zampini accessed via the prefix `-tao_bqnk_`. For options database, see TAOBNK 233850be85SAlp Dener 243850be85SAlp Dener Level: beginner 259fa2b5dcSStefano Zampini .seealso TAOBNK, TAOBQNKTR, TAOBQNKLS 263850be85SAlp Dener M*/ 27e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BQNKTL(Tao tao) 28e0ed867bSAlp Dener { 29414d97d3SAlp Dener TAO_BNK *bnk; 30414d97d3SAlp Dener TAO_BQNK *bqnk; 31e0ed867bSAlp Dener PetscErrorCode ierr; 32e0ed867bSAlp Dener 33e0ed867bSAlp Dener PetscFunctionBegin; 34e0ed867bSAlp Dener ierr = TaoCreate_BQNK(tao);CHKERRQ(ierr); 355eb5f4d6SAlp Dener tao->ops->setup = TaoSetUp_BQNKTL; 36414d97d3SAlp Dener bnk = (TAO_BNK*)tao->data; 37414d97d3SAlp Dener bqnk = (TAO_BQNK*)bnk->ctx; 38414d97d3SAlp Dener bqnk->solve = TaoSolve_BNTL; 39e0ed867bSAlp Dener PetscFunctionReturn(0); 40e0ed867bSAlp Dener } 41