xref: /petsc/src/tao/bound/impls/bqnls/bqnls.c (revision f5766c097db85df4ee62e2acb750cb76ef9efd11)
1 #include <../src/tao/bound/impls/bqnk/bqnk.h>
2 
3 static PetscErrorCode TaoBQNLSComputeHessian(Tao tao)
4 {
5   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
6   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
7   PetscErrorCode ierr;
8   PetscReal      gnorm2, delta;
9 
10   PetscFunctionBegin;
11   /* Compute the initial scaling and update the approximation */
12   gnorm2 = bnk->gnorm*bnk->gnorm;
13   delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / PetscMax(gnorm2, PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0));
14   ierr = MatSymBrdnSetDelta(bqnk->B, delta);CHKERRQ(ierr);
15   ierr = MatLMVMUpdate(bqnk->B, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
16   PetscFunctionReturn(0);
17 }
18 
19 static PetscErrorCode TaoBQNLSComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
20 {
21   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
22   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
23   PetscErrorCode ierr;
24   PetscInt       nupdates;
25 
26   PetscFunctionBegin;
27   ierr = MatSolve(bqnk->B, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
28   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
29   ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
30   *ksp_reason = KSP_CONVERGED_ATOL;
31   ierr = MatLMVMGetUpdateCount(bqnk->B, &nupdates);CHKERRQ(ierr);
32   if (nupdates == 0) {
33     *step_type = BNK_SCALED_GRADIENT;
34   } else {
35     *step_type = BNK_BFGS;
36   }
37   PetscFunctionReturn(0);
38 }
39 
40 static PetscErrorCode TaoSetFromOptions_BQNLS(PetscOptionItems *PetscOptionsObject,Tao tao)
41 {
42   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
43   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
44   PetscErrorCode ierr;
45   KSPType        ksp_type;
46   PetscBool      is_spd;
47 
48   PetscFunctionBegin;
49   ierr = PetscOptionsHead(PetscOptionsObject,"Quasi-Newton-Krylov method for bound constrained optimization");CHKERRQ(ierr);
50   ierr = PetscOptionsEList("-tao_bqnls_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, 0);CHKERRQ(ierr);
51   ierr = PetscOptionsReal("-tao_bqnls_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
52   ierr = PetscOptionsReal("-tao_bqnls_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr);
53   ierr = PetscOptionsReal("-tao_bqnls_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr);
54   ierr = PetscOptionsInt("-tao_bqnls_max_cg_its", "number of BNCG iterations to take for each Newton step", "", bnk->max_cg_its, &bnk->max_cg_its,NULL);CHKERRQ(ierr);
55   ierr = PetscOptionsTail();CHKERRQ(ierr);
56   ierr = TaoSetFromOptions(bnk->bncg);CHKERRQ(ierr);
57   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
58   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
59   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
60   bnk->is_nash = bnk->is_gltr = bnk->is_stcg = PETSC_FALSE;
61   ierr = MatSetFromOptions(bqnk->B);CHKERRQ(ierr);
62   ierr = MatGetOption(bqnk->B, MAT_SPD, &is_spd);CHKERRQ(ierr);
63   if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite");
64   PetscFunctionReturn(0);
65 }
66 
67 PETSC_EXTERN PetscErrorCode TaoCreate_BQNLS(Tao tao)
68 {
69   TAO_BNK        *bnk;
70   TAO_BQNK       *bqnk;
71   PetscErrorCode ierr;
72 
73   PetscFunctionBegin;
74   ierr = TaoCreate_BQNK(tao);CHKERRQ(ierr);
75   ierr = KSPSetOptionsPrefix(tao->ksp, "unused");CHKERRQ(ierr);
76   tao->ops->solve = TaoSolve_BNLS;
77   tao->ops->setfromoptions = TaoSetFromOptions_BQNLS;
78 
79   bnk = (TAO_BNK*)tao->data;
80   bnk->update_type = BNK_UPDATE_STEP;
81   bnk->computehessian = TaoBQNLSComputeHessian;
82   bnk->computestep = TaoBQNLSComputeStep;
83 
84   bqnk = (TAO_BQNK*)bnk->ctx;
85   ierr = MatSetOptionsPrefix(bqnk->B, "tao_bqnls_");CHKERRQ(ierr);
86   ierr = MatSetType(bqnk->B, MATLMVMBFGS);CHKERRQ(ierr);
87   PetscFunctionReturn(0);
88 }