xref: /petsc/src/tao/bound/impls/bnk/bntl.c (revision 62675bee48cd89cd8757a6de6d7b87f2edd3afcc)
1c14b763aSAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
2c14b763aSAlp Dener #include <petscksp.h>
3c14b763aSAlp Dener 
4c14b763aSAlp Dener /*
5c14b763aSAlp Dener  Implements Newton's Method with a trust region approach for solving
6c14b763aSAlp Dener  bound constrained minimization problems. This version includes a
7c14b763aSAlp Dener  line search fall-back in the event of a trust region failure.
8c14b763aSAlp Dener 
9df278d8fSAlp Dener  The linear system solve has to be done with a conjugate gradient method.
10c14b763aSAlp Dener */
11c14b763aSAlp Dener 
12c14b763aSAlp Dener static PetscErrorCode TaoSolve_BNTL(Tao tao)
13c14b763aSAlp Dener {
14c14b763aSAlp Dener   PetscErrorCode               ierr;
15c14b763aSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
16e465cd6fSAlp Dener   KSPConvergedReason           ksp_reason;
17c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
18c14b763aSAlp Dener 
19e465cd6fSAlp Dener   PetscReal                    oldTrust, prered, actred, stepNorm, steplen;
20*62675beeSAlp Dener   PetscBool                    stepAccepted = PETSC_TRUE, shift = PETSC_FALSE;
2128017e9fSAlp Dener   PetscInt                     stepType;
22c14b763aSAlp Dener 
23c14b763aSAlp Dener   PetscFunctionBegin;
2428017e9fSAlp Dener   /* Initialize the preconditioner, KSP solver and trust radius/line search */
25c14b763aSAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
26*62675beeSAlp Dener   ierr = TaoBNKInitialize(tao, bnk->init_type);CHKERRQ(ierr);
2728017e9fSAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
28c14b763aSAlp Dener 
29c14b763aSAlp Dener   /* Have not converged; continue with Newton method */
30c14b763aSAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
31c14b763aSAlp Dener     tao->niter++;
32c14b763aSAlp Dener     tao->ksp_its=0;
33*62675beeSAlp Dener 
34*62675beeSAlp Dener     /* Compute the hessian and update the BFGS preconditioner at the new iterate*/
35*62675beeSAlp Dener     ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
36c14b763aSAlp Dener 
378d5ead36SAlp Dener     /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
38*62675beeSAlp Dener     ierr = TaoBNKComputeStep(tao, shift, &ksp_reason);CHKERRQ(ierr);
39c14b763aSAlp Dener 
40c14b763aSAlp Dener     /* Store current solution before it changes */
41c14b763aSAlp Dener     oldTrust = tao->trust;
42c14b763aSAlp Dener     bnk->fold = bnk->f;
43c14b763aSAlp Dener     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
44c14b763aSAlp Dener     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
45c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
46c14b763aSAlp Dener 
47c14b763aSAlp Dener     /* Temporarily accept the step and project it into the bounds */
48c14b763aSAlp Dener     ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr);
49c14b763aSAlp Dener     ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
50c14b763aSAlp Dener 
51c14b763aSAlp Dener     /* Check if the projection changed the step direction */
52c14b763aSAlp Dener     ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
538d5ead36SAlp Dener     ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
54c14b763aSAlp Dener     ierr = VecNorm(tao->stepdirection, NORM_2, &stepNorm);CHKERRQ(ierr);
55c14b763aSAlp Dener     if (stepNorm != bnk->dnorm) {
568d5ead36SAlp Dener       /* Projection changed the step, so we have to recompute predicted reduction.
578d5ead36SAlp Dener          However, we deliberately do not change the step norm and the trust radius
588d5ead36SAlp Dener          in order for the safeguard to more closely mimic a piece-wise linesearch
598d5ead36SAlp Dener          along the bounds. */
6028017e9fSAlp Dener       ierr = MatMult(bnk->H_inactive, tao->stepdirection, bnk->Xwork);CHKERRQ(ierr);
61c14b763aSAlp Dener       ierr = VecAYPX(bnk->Xwork, -0.5, tao->gradient);CHKERRQ(ierr);
62c14b763aSAlp Dener       ierr = VecDot(bnk->Xwork, tao->stepdirection, &prered);
63c14b763aSAlp Dener     } else {
64c14b763aSAlp Dener       /* Step did not change, so we can just recover the pre-computed prediction */
65c14b763aSAlp Dener       ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr);
66c14b763aSAlp Dener     }
67c14b763aSAlp Dener     prered = -prered;
68c14b763aSAlp Dener 
69c14b763aSAlp Dener     /* Compute the actual reduction and update the trust radius */
70c14b763aSAlp Dener     ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr);
71c14b763aSAlp Dener     actred = bnk->fold - bnk->f;
7228017e9fSAlp Dener     ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr);
73c14b763aSAlp Dener 
74c14b763aSAlp Dener     if (stepAccepted) {
75c14b763aSAlp Dener       /* Step is good, evaluate the gradient and the hessian */
768d5ead36SAlp Dener       steplen = 1.0;
77e465cd6fSAlp Dener       ++bnk->newt;
78c14b763aSAlp Dener       ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
79c14b763aSAlp Dener       ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
80c14b763aSAlp Dener     } else {
81c14b763aSAlp Dener       /* Trust-region rejected the step. Revert the solution. */
82c14b763aSAlp Dener       bnk->f = bnk->fold;
83c14b763aSAlp Dener       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
84c14b763aSAlp Dener 
85c14b763aSAlp Dener       /* Trigger the line search */
86e465cd6fSAlp Dener       ierr = TaoBNKSafeguardStep(tao, ksp_reason, &stepType);CHKERRQ(ierr);
87c14b763aSAlp Dener       ierr = TaoBNKPerformLineSearch(tao, stepType, &steplen, &ls_reason);CHKERRQ(ierr);
88c14b763aSAlp Dener       if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
89c14b763aSAlp Dener         /* Line search failed, revert solution and terminate */
90c14b763aSAlp Dener         bnk->f = bnk->fold;
91c14b763aSAlp Dener         ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
92c14b763aSAlp Dener         ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
93c14b763aSAlp Dener         ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
94c14b763aSAlp Dener         tao->trust = 0.0;
95c14b763aSAlp Dener         tao->reason = TAO_DIVERGED_LS_FAILURE;
96c14b763aSAlp Dener       } else {
97c14b763aSAlp Dener         /* Line search succeeded so we should update the trust radius based on the LS step length */
98770b7498SAlp Dener         tao->trust = oldTrust;
9928017e9fSAlp Dener         ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, BNK_UPDATE_STEP, stepType, &stepAccepted);CHKERRQ(ierr);
100*62675beeSAlp Dener         /* count the accepted step type */
101*62675beeSAlp Dener         ierr = TaoBNKAddStepCounts(tao, stepType);CHKERRQ(ierr);
102c14b763aSAlp Dener       }
103c14b763aSAlp Dener     }
104c14b763aSAlp Dener 
105c14b763aSAlp Dener     /*  Check for termination */
106c14b763aSAlp Dener     ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
107c14b763aSAlp Dener     if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
108c14b763aSAlp Dener     ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
1098d5ead36SAlp Dener     ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,steplen);CHKERRQ(ierr);
110c14b763aSAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
111c14b763aSAlp Dener   }
112c14b763aSAlp Dener   PetscFunctionReturn(0);
113c14b763aSAlp Dener }
114c14b763aSAlp Dener 
115df278d8fSAlp Dener /*------------------------------------------------------------*/
116df278d8fSAlp Dener 
117c14b763aSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNTL(Tao tao)
118c14b763aSAlp Dener {
119c14b763aSAlp Dener   TAO_BNK        *bnk;
120c14b763aSAlp Dener   PetscErrorCode ierr;
121c14b763aSAlp Dener 
122c14b763aSAlp Dener   PetscFunctionBegin;
123c14b763aSAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
124c14b763aSAlp Dener   tao->ops->solve=TaoSolve_BNTL;
125c14b763aSAlp Dener 
126c14b763aSAlp Dener   bnk = (TAO_BNK *)tao->data;
127c14b763aSAlp Dener   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
128c14b763aSAlp Dener   bnk->sval = 0.0; /* disable Hessian shifting */
129c14b763aSAlp Dener   PetscFunctionReturn(0);
130c14b763aSAlp Dener }