xref: /petsc/src/tao/bound/impls/bqnls/bqnls.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 #include <../src/tao/bound/impls/bqnk/bqnk.h>
2 
3 static const char *BNK_AS[64] = {"none", "bertsekas"};
4 
5 static PetscErrorCode TaoBQNLSComputeHessian(Tao tao) {
6   TAO_BNK  *bnk  = (TAO_BNK *)tao->data;
7   TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
8   PetscReal gnorm2, delta;
9 
10   PetscFunctionBegin;
11   /* Compute the initial scaling and update the approximation */
12   gnorm2 = bnk->gnorm * bnk->gnorm;
13   if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON;
14   if (bnk->f == 0.0) delta = 2.0 / gnorm2;
15   else delta = 2.0 * PetscAbsScalar(bnk->f) / gnorm2;
16   PetscCall(MatLMVMSymBroydenSetDelta(bqnk->B, delta));
17   PetscCall(MatLMVMUpdate(bqnk->B, tao->solution, bnk->unprojected_gradient));
18   PetscFunctionReturn(0);
19 }
20 
21 static PetscErrorCode TaoBQNLSComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type) {
22   TAO_BNK  *bnk  = (TAO_BNK *)tao->data;
23   TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
24   PetscInt  nupdates;
25 
26   PetscFunctionBegin;
27   PetscCall(MatSolve(bqnk->B, tao->gradient, tao->stepdirection));
28   PetscCall(VecScale(tao->stepdirection, -1.0));
29   PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
30   *ksp_reason = KSP_CONVERGED_ATOL;
31   PetscCall(MatLMVMGetUpdateCount(bqnk->B, &nupdates));
32   if (nupdates == 0) *step_type = BNK_SCALED_GRADIENT;
33   else *step_type = BNK_BFGS;
34   PetscFunctionReturn(0);
35 }
36 
37 static PetscErrorCode TaoSetFromOptions_BQNLS(Tao tao, PetscOptionItems *PetscOptionsObject) {
38   TAO_BNK  *bnk  = (TAO_BNK *)tao->data;
39   TAO_BQNK *bqnk = (TAO_BQNK *)bnk->ctx;
40   PetscBool is_set, is_spd;
41 
42   PetscFunctionBegin;
43   PetscOptionsHeadBegin(PetscOptionsObject, "Quasi-Newton-Krylov method for bound constrained optimization");
44   PetscCall(PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, NULL));
45   PetscCall(PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon, NULL));
46   PetscCall(PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol, NULL));
47   PetscCall(PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step, NULL));
48   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));
49   PetscOptionsHeadEnd();
50 
51   PetscCall(TaoSetOptionsPrefix(bnk->bncg, ((PetscObject)(tao))->prefix));
52   PetscCall(TaoAppendOptionsPrefix(bnk->bncg, "tao_bnk_"));
53   PetscCall(TaoSetFromOptions(bnk->bncg));
54 
55   PetscCall(MatSetOptionsPrefix(bqnk->B, ((PetscObject)tao)->prefix));
56   PetscCall(MatAppendOptionsPrefix(bqnk->B, "tao_bqnls_"));
57   PetscCall(MatSetFromOptions(bqnk->B));
58   PetscCall(MatIsSPDKnown(bqnk->B, &is_set, &is_spd));
59   PetscCheck(is_set && is_spd, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite");
60   PetscFunctionReturn(0);
61 }
62 
63 /*MC
64   TAOBQNLS - Bounded Quasi-Newton Line Search method for nonlinear minimization with bound
65              constraints. This method approximates the action of the inverse-Hessian with a
66              limited memory quasi-Newton formula. The quasi-Newton matrix and its options are
67              accessible via the prefix `-tao_bqnls_`
68 
69   Option Database Keys:
70 + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
71 . -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
72 . -tao_bnk_epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance
73 . -tao_bnk_as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
74 - -tao_bnk_as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
75 
76   Level: beginner
77 .seealso: `TAOBNK`
78 M*/
79 PETSC_EXTERN PetscErrorCode TaoCreate_BQNLS(Tao tao) {
80   TAO_BNK  *bnk;
81   TAO_BQNK *bqnk;
82 
83   PetscFunctionBegin;
84   PetscCall(TaoCreate_BQNK(tao));
85   tao->ops->setfromoptions = TaoSetFromOptions_BQNLS;
86 
87   bnk                 = (TAO_BNK *)tao->data;
88   bnk->update_type    = BNK_UPDATE_STEP;
89   bnk->computehessian = TaoBQNLSComputeHessian;
90   bnk->computestep    = TaoBQNLSComputeStep;
91 
92   bqnk        = (TAO_BQNK *)bnk->ctx;
93   bqnk->solve = TaoSolve_BNLS;
94   PetscCall(MatSetType(bqnk->B, MATLMVMBFGS));
95   PetscFunctionReturn(0);
96 }
97