xref: /petsc/src/tao/bound/impls/bqnls/bqnls.c (revision 8ebe3e4e9e00d86ece2e9fcd0cc84910b0ad437c)
16b591159SAlp Dener #include <../src/tao/bound/impls/bqnk/bqnk.h>
26b591159SAlp Dener 
370a3f44bSAlp Dener static const char *BNK_AS[64] = {"none", "bertsekas"};
470a3f44bSAlp Dener 
56b591159SAlp Dener static PetscErrorCode TaoBQNLSComputeHessian(Tao tao)
66b591159SAlp Dener {
76b591159SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
86b591159SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
96b591159SAlp Dener   PetscErrorCode ierr;
10d5ae2380SAlp Dener   PetscReal      gnorm2, delta;
116b591159SAlp Dener 
126b591159SAlp Dener   PetscFunctionBegin;
13f5766c09SAlp Dener   /* Compute the initial scaling and update the approximation */
14d5ae2380SAlp Dener   gnorm2 = bnk->gnorm*bnk->gnorm;
158cabe928SAlp Dener   if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON;
16*8ebe3e4eSStefano Zampini   if (bnk->f == 0.0) delta = 2.0 / gnorm2;
17*8ebe3e4eSStefano Zampini   else delta = 2.0 * PetscAbsScalar(bnk->f) / gnorm2;
18864588a7SAlp Dener   ierr = MatLMVMSymBroydenSetDelta(bqnk->B, delta);CHKERRQ(ierr);
196b591159SAlp Dener   ierr = MatLMVMUpdate(bqnk->B, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
206b591159SAlp Dener   PetscFunctionReturn(0);
216b591159SAlp Dener }
226b591159SAlp Dener 
236b591159SAlp Dener static PetscErrorCode TaoBQNLSComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
246b591159SAlp Dener {
256b591159SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
266b591159SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
276b591159SAlp Dener   PetscErrorCode ierr;
2865f5217aSAlp Dener   PetscInt       nupdates;
296b591159SAlp Dener 
306b591159SAlp Dener   PetscFunctionBegin;
319515a401SAlp Dener   ierr = MatSolve(bqnk->B, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
326b591159SAlp Dener   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
336b591159SAlp Dener   ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
346b591159SAlp Dener   *ksp_reason = KSP_CONVERGED_ATOL;
3565f5217aSAlp Dener   ierr = MatLMVMGetUpdateCount(bqnk->B, &nupdates);CHKERRQ(ierr);
36*8ebe3e4eSStefano Zampini   if (nupdates == 0) *step_type = BNK_SCALED_GRADIENT;
37*8ebe3e4eSStefano Zampini   else *step_type = BNK_BFGS;
386b591159SAlp Dener   PetscFunctionReturn(0);
396b591159SAlp Dener }
406b591159SAlp Dener 
416b591159SAlp Dener static PetscErrorCode TaoSetFromOptions_BQNLS(PetscOptionItems *PetscOptionsObject,Tao tao)
426b591159SAlp Dener {
436b591159SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
446b591159SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
456b591159SAlp Dener   PetscErrorCode ierr;
466b591159SAlp Dener   PetscBool      is_spd;
476b591159SAlp Dener 
486b591159SAlp Dener   PetscFunctionBegin;
496b591159SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Quasi-Newton-Krylov method for bound constrained optimization");CHKERRQ(ierr);
509fa2b5dcSStefano Zampini   ierr = PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, NULL);CHKERRQ(ierr);
519fa2b5dcSStefano Zampini   ierr = PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
529fa2b5dcSStefano Zampini   ierr = PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr);
539fa2b5dcSStefano Zampini   ierr = PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr);
549fa2b5dcSStefano Zampini   ierr = 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);CHKERRQ(ierr);
556b591159SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
56*8ebe3e4eSStefano Zampini 
57*8ebe3e4eSStefano Zampini   ierr = TaoSetOptionsPrefix(bnk->bncg,((PetscObject)(tao))->prefix);CHKERRQ(ierr);
58*8ebe3e4eSStefano Zampini   ierr = TaoAppendOptionsPrefix(bnk->bncg,"tao_bnk_");CHKERRQ(ierr);
596b591159SAlp Dener   ierr = TaoSetFromOptions(bnk->bncg);CHKERRQ(ierr);
60*8ebe3e4eSStefano Zampini 
61*8ebe3e4eSStefano Zampini   ierr = MatSetOptionsPrefix(bqnk->B, ((PetscObject)tao)->prefix);CHKERRQ(ierr);
62*8ebe3e4eSStefano Zampini   ierr = MatAppendOptionsPrefix(bqnk->B, "tao_bqnls_");CHKERRQ(ierr);
636b591159SAlp Dener   ierr = MatSetFromOptions(bqnk->B);CHKERRQ(ierr);
646b591159SAlp Dener   ierr = MatGetOption(bqnk->B, MAT_SPD, &is_spd);CHKERRQ(ierr);
656b591159SAlp Dener   if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite");
666b591159SAlp Dener   PetscFunctionReturn(0);
676b591159SAlp Dener }
686b591159SAlp Dener 
693850be85SAlp Dener /*MC
703850be85SAlp Dener   TAOBQNLS - Bounded Quasi-Newton Line Search method for nonlinear minimization with bound
713850be85SAlp Dener              constraints. This method approximates the action of the inverse-Hessian with a
723850be85SAlp Dener              limited memory quasi-Newton formula. The quasi-Newton matrix and its options are
733850be85SAlp Dener              accessible via the prefix `-tao_bqnls_`
743850be85SAlp Dener 
759fa2b5dcSStefano Zampini   Option Database Keys:
769fa2b5dcSStefano Zampini + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
779fa2b5dcSStefano Zampini . -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
789fa2b5dcSStefano Zampini . -tao_bnk_epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance
799fa2b5dcSStefano Zampini . -tao_bnk_as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
809fa2b5dcSStefano Zampini . -tao_bnk_as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
813850be85SAlp Dener 
823850be85SAlp Dener   Level: beginner
839fa2b5dcSStefano Zampini .seealso: TAOBNK
843850be85SAlp Dener M*/
856b591159SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BQNLS(Tao tao)
866b591159SAlp Dener {
876b591159SAlp Dener   TAO_BNK        *bnk;
886b591159SAlp Dener   TAO_BQNK       *bqnk;
896b591159SAlp Dener   PetscErrorCode ierr;
906b591159SAlp Dener 
916b591159SAlp Dener   PetscFunctionBegin;
926b591159SAlp Dener   ierr = TaoCreate_BQNK(tao);CHKERRQ(ierr);
936b591159SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BQNLS;
946b591159SAlp Dener 
956b591159SAlp Dener   bnk = (TAO_BNK*)tao->data;
966b591159SAlp Dener   bnk->update_type = BNK_UPDATE_STEP;
976b591159SAlp Dener   bnk->computehessian = TaoBQNLSComputeHessian;
986b591159SAlp Dener   bnk->computestep = TaoBQNLSComputeStep;
996b591159SAlp Dener 
1006b591159SAlp Dener   bqnk = (TAO_BQNK*)bnk->ctx;
101414d97d3SAlp Dener   bqnk->solve = TaoSolve_BNLS;
1026b591159SAlp Dener   ierr = MatSetType(bqnk->B, MATLMVMBFGS);CHKERRQ(ierr);
1036b591159SAlp Dener   PetscFunctionReturn(0);
1046b591159SAlp Dener }
105