xref: /petsc/src/tao/bound/impls/bnk/bntr.c (revision 62675bee48cd89cd8757a6de6d7b87f2edd3afcc)
1fed79b8eSAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
2fed79b8eSAlp Dener #include <petscksp.h>
3fed79b8eSAlp Dener 
4fed79b8eSAlp Dener /*
5fed79b8eSAlp Dener  Implements Newton's Method with a trust region approach for solving
6fed79b8eSAlp Dener  bound constrained minimization problems.
7fed79b8eSAlp Dener 
8df278d8fSAlp Dener  The linear system solve has to be done with a conjugate gradient method.
9fed79b8eSAlp Dener */
10fed79b8eSAlp Dener 
11fed79b8eSAlp Dener static PetscErrorCode TaoSolve_BNTR(Tao tao)
12fed79b8eSAlp Dener {
13fed79b8eSAlp Dener   PetscErrorCode               ierr;
14fed79b8eSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
15e465cd6fSAlp Dener   KSPConvergedReason           ksp_reason;
16fed79b8eSAlp Dener 
178d5ead36SAlp Dener   PetscReal                    oldTrust, prered, actred, stepNorm, steplen;
18*62675beeSAlp Dener   PetscBool                    stepAccepted = PETSC_TRUE, shift = PETSC_FALSE;
19e465cd6fSAlp Dener   PetscInt                     stepType = BNK_NEWTON;
20fed79b8eSAlp Dener 
21fed79b8eSAlp Dener   PetscFunctionBegin;
2228017e9fSAlp Dener   /* Initialize the preconditioner, KSP solver and trust radius/line search */
23fed79b8eSAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
24*62675beeSAlp Dener   ierr = TaoBNKInitialize(tao, bnk->init_type);CHKERRQ(ierr);
2528017e9fSAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
26fed79b8eSAlp Dener 
27fed79b8eSAlp Dener   /* Have not converged; continue with Newton method */
28fed79b8eSAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
2966ed3702SAlp Dener 
30fed79b8eSAlp Dener     if (stepAccepted) {
31fed79b8eSAlp Dener       tao->niter++;
32fed79b8eSAlp Dener       tao->ksp_its=0;
33*62675beeSAlp Dener       /* Compute the hessian and update the BFGS preconditioner at the new iterate*/
34*62675beeSAlp Dener       ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
35fed79b8eSAlp Dener     }
36fed79b8eSAlp 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);
39fed79b8eSAlp Dener 
40fed79b8eSAlp Dener     /* Store current solution before it changes */
41fed79b8eSAlp Dener     oldTrust = tao->trust;
42fed79b8eSAlp Dener     bnk->fold = bnk->f;
43fed79b8eSAlp Dener     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
44fed79b8eSAlp Dener     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
45fed79b8eSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
46fed79b8eSAlp Dener 
47b1c2d0e3SAlp Dener     /* Temporarily accept the step and project it into the bounds */
48fed79b8eSAlp Dener     ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr);
49b1c2d0e3SAlp Dener     ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
50b1c2d0e3SAlp Dener 
51b1c2d0e3SAlp Dener     /* Check if the projection changed the step direction */
52b1c2d0e3SAlp Dener     ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
538d5ead36SAlp Dener     ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
54b1c2d0e3SAlp Dener     ierr = VecNorm(tao->stepdirection, NORM_2, &stepNorm);CHKERRQ(ierr);
55b1c2d0e3SAlp 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);
61b1c2d0e3SAlp Dener       ierr = VecAYPX(bnk->Xwork, -0.5, tao->gradient);CHKERRQ(ierr);
62b1c2d0e3SAlp Dener       ierr = VecDot(bnk->Xwork, tao->stepdirection, &prered);
63b1c2d0e3SAlp Dener     } else {
64b1c2d0e3SAlp Dener       /* Step did not change, so we can just recover the pre-computed prediction */
65b1c2d0e3SAlp Dener       ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr);
66b1c2d0e3SAlp Dener     }
67b1c2d0e3SAlp Dener     prered = -prered;
68b1c2d0e3SAlp Dener 
69b1c2d0e3SAlp Dener     /* Compute the actual reduction and update the trust radius */
70fed79b8eSAlp Dener     ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr);
71b1c2d0e3SAlp Dener     actred = bnk->fold - bnk->f;
7228017e9fSAlp Dener     ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr);
73fed79b8eSAlp Dener 
74fed79b8eSAlp Dener     if (stepAccepted) {
7566ed3702SAlp Dener       /* Step is good, evaluate the gradient and the hessian */
768d5ead36SAlp Dener       steplen = 1.0;
77e465cd6fSAlp Dener       ++bnk->newt;
78fed79b8eSAlp Dener       ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
79fed79b8eSAlp Dener       ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
80fed79b8eSAlp Dener     } else {
81fed79b8eSAlp Dener       /* Step is bad, revert old solution and re-solve with new radius*/
828d5ead36SAlp Dener       steplen = 0.0;
83fed79b8eSAlp Dener       bnk->f = bnk->fold;
84fed79b8eSAlp Dener       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
85fed79b8eSAlp Dener       ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
86fed79b8eSAlp Dener       ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
8773e4db90SAlp Dener       if (oldTrust == tao->trust) {
8873e4db90SAlp Dener         /* Can't change the radius anymore so just terminate */
89fed79b8eSAlp Dener         tao->reason = TAO_DIVERGED_TR_REDUCTION;
90fed79b8eSAlp Dener       }
91fed79b8eSAlp Dener     }
92fed79b8eSAlp Dener 
93fed79b8eSAlp Dener     /*  Check for termination */
94fed79b8eSAlp Dener     ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
9573e4db90SAlp Dener     if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
96fed79b8eSAlp Dener     ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
978d5ead36SAlp Dener     ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,steplen);CHKERRQ(ierr);
98fed79b8eSAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
99fed79b8eSAlp Dener   }
100fed79b8eSAlp Dener   PetscFunctionReturn(0);
101fed79b8eSAlp Dener }
102fed79b8eSAlp Dener 
103df278d8fSAlp Dener /*------------------------------------------------------------*/
104df278d8fSAlp Dener 
105fed79b8eSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao)
106fed79b8eSAlp Dener {
107fed79b8eSAlp Dener   TAO_BNK        *bnk;
108fed79b8eSAlp Dener   PetscErrorCode ierr;
109fed79b8eSAlp Dener 
110fed79b8eSAlp Dener   PetscFunctionBegin;
111fed79b8eSAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
112fed79b8eSAlp Dener   tao->ops->solve=TaoSolve_BNTR;
113fed79b8eSAlp Dener 
114fed79b8eSAlp Dener   bnk = (TAO_BNK *)tao->data;
11566ed3702SAlp Dener   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
11666ed3702SAlp Dener   bnk->sval = 0.0; /* disable Hessian shifting */
117fed79b8eSAlp Dener   PetscFunctionReturn(0);
118fed79b8eSAlp Dener }