16b591159SAlp Dener #include <../src/tao/bound/impls/bqnk/bqnk.h>
26b591159SAlp Dener
370a3f44bSAlp Dener static const char *BNK_AS[64] = {"none", "bertsekas"};
470a3f44bSAlp Dener
TaoBQNLSComputeHessian(Tao tao)5d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoBQNLSComputeHessian(Tao tao)
6d71ae5a4SJacob Faibussowitsch {
76b591159SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
86b591159SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
9d5ae2380SAlp Dener PetscReal gnorm2, delta;
106b591159SAlp Dener
116b591159SAlp Dener PetscFunctionBegin;
12f5766c09SAlp Dener /* Compute the initial scaling and update the approximation */
13d5ae2380SAlp Dener gnorm2 = bnk->gnorm * bnk->gnorm;
148cabe928SAlp Dener if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON;
158ebe3e4eSStefano Zampini if (bnk->f == 0.0) delta = 2.0 / gnorm2;
168ebe3e4eSStefano Zampini else delta = 2.0 * PetscAbsScalar(bnk->f) / gnorm2;
179566063dSJacob Faibussowitsch PetscCall(MatLMVMSymBroydenSetDelta(bqnk->B, delta));
189566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(bqnk->B, tao->solution, bnk->unprojected_gradient));
193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
206b591159SAlp Dener }
216b591159SAlp Dener
TaoBQNLSComputeStep(Tao tao,PetscBool shift,KSPConvergedReason * ksp_reason,PetscInt * step_type)22d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoBQNLSComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
23d71ae5a4SJacob Faibussowitsch {
246b591159SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
256b591159SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
2665f5217aSAlp Dener PetscInt nupdates;
276b591159SAlp Dener
286b591159SAlp Dener PetscFunctionBegin;
299566063dSJacob Faibussowitsch PetscCall(MatSolve(bqnk->B, tao->gradient, tao->stepdirection));
309566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0));
319566063dSJacob Faibussowitsch PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
326b591159SAlp Dener *ksp_reason = KSP_CONVERGED_ATOL;
339566063dSJacob Faibussowitsch PetscCall(MatLMVMGetUpdateCount(bqnk->B, &nupdates));
348ebe3e4eSStefano Zampini if (nupdates == 0) *step_type = BNK_SCALED_GRADIENT;
358ebe3e4eSStefano Zampini else *step_type = BNK_BFGS;
363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
376b591159SAlp Dener }
386b591159SAlp Dener
TaoSetFromOptions_BQNLS(Tao tao,PetscOptionItems PetscOptionsObject)39*ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_BQNLS(Tao tao, PetscOptionItems PetscOptionsObject)
40d71ae5a4SJacob Faibussowitsch {
416b591159SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data;
426b591159SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
43b94d7dedSBarry Smith PetscBool is_set, is_spd;
446b591159SAlp Dener
456b591159SAlp Dener PetscFunctionBegin;
46d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Quasi-Newton-Krylov method for bound constrained optimization");
479566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, NULL));
489566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon, NULL));
499566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol, NULL));
509566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step, NULL));
519566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_bnk_max_cg_its", "number of BNCG iterations to take for each Newton step", "", bnk->max_cg_its, &bnk->max_cg_its, NULL));
52d0609cedSBarry Smith PetscOptionsHeadEnd();
538ebe3e4eSStefano Zampini
54f4f49eeaSPierre Jolivet PetscCall(TaoSetOptionsPrefix(bnk->bncg, ((PetscObject)tao)->prefix));
559566063dSJacob Faibussowitsch PetscCall(TaoAppendOptionsPrefix(bnk->bncg, "tao_bnk_"));
569566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(bnk->bncg));
578ebe3e4eSStefano Zampini
589566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(bqnk->B, ((PetscObject)tao)->prefix));
599566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(bqnk->B, "tao_bqnls_"));
609566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(bqnk->B));
61b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(bqnk->B, &is_set, &is_spd));
62b94d7dedSBarry Smith PetscCheck(is_set && is_spd, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite");
633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
646b591159SAlp Dener }
656b591159SAlp Dener
663850be85SAlp Dener /*MC
673850be85SAlp Dener TAOBQNLS - Bounded Quasi-Newton Line Search method for nonlinear minimization with bound
683850be85SAlp Dener constraints. This method approximates the action of the inverse-Hessian with a
693850be85SAlp Dener limited memory quasi-Newton formula. The quasi-Newton matrix and its options are
703850be85SAlp Dener accessible via the prefix `-tao_bqnls_`
713850be85SAlp Dener
722fe279fdSBarry Smith Options Database Keys:
739fa2b5dcSStefano Zampini + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
749fa2b5dcSStefano Zampini . -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
759fa2b5dcSStefano Zampini . -tao_bnk_epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance
769fa2b5dcSStefano Zampini . -tao_bnk_as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
77f1a722f8SMatthew G. Knepley - -tao_bnk_as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
783850be85SAlp Dener
793850be85SAlp Dener Level: beginner
802fe279fdSBarry Smith
81db781477SPatrick Sanan .seealso: `TAOBNK`
823850be85SAlp Dener M*/
TaoCreate_BQNLS(Tao tao)83d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BQNLS(Tao tao)
84d71ae5a4SJacob Faibussowitsch {
856b591159SAlp Dener TAO_BNK *bnk;
866b591159SAlp Dener TAO_BQNK *bqnk;
876b591159SAlp Dener
886b591159SAlp Dener PetscFunctionBegin;
899566063dSJacob Faibussowitsch PetscCall(TaoCreate_BQNK(tao));
906b591159SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BQNLS;
916b591159SAlp Dener
926b591159SAlp Dener bnk = (TAO_BNK *)tao->data;
936b591159SAlp Dener bnk->update_type = BNK_UPDATE_STEP;
946b591159SAlp Dener bnk->computehessian = TaoBQNLSComputeHessian;
956b591159SAlp Dener bnk->computestep = TaoBQNLSComputeStep;
966b591159SAlp Dener
976b591159SAlp Dener bqnk = (TAO_BQNK *)bnk->ctx;
98414d97d3SAlp Dener bqnk->solve = TaoSolve_BNLS;
999566063dSJacob Faibussowitsch PetscCall(MatSetType(bqnk->B, MATLMVMBFGS));
1003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1016b591159SAlp Dener }
102