xref: /petsc/src/tao/bound/impls/bnk/bntr.c (revision 62675bee48cd89cd8757a6de6d7b87f2edd3afcc)
1 #include <../src/tao/bound/impls/bnk/bnk.h>
2 #include <petscksp.h>
3 
4 /*
5  Implements Newton's Method with a trust region approach for solving
6  bound constrained minimization problems.
7 
8  The linear system solve has to be done with a conjugate gradient method.
9 */
10 
11 static PetscErrorCode TaoSolve_BNTR(Tao tao)
12 {
13   PetscErrorCode               ierr;
14   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
15   KSPConvergedReason           ksp_reason;
16 
17   PetscReal                    oldTrust, prered, actred, stepNorm, steplen;
18   PetscBool                    stepAccepted = PETSC_TRUE, shift = PETSC_FALSE;
19   PetscInt                     stepType = BNK_NEWTON;
20 
21   PetscFunctionBegin;
22   /* Initialize the preconditioner, KSP solver and trust radius/line search */
23   tao->reason = TAO_CONTINUE_ITERATING;
24   ierr = TaoBNKInitialize(tao, bnk->init_type);CHKERRQ(ierr);
25   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
26 
27   /* Have not converged; continue with Newton method */
28   while (tao->reason == TAO_CONTINUE_ITERATING) {
29 
30     if (stepAccepted) {
31       tao->niter++;
32       tao->ksp_its=0;
33       /* Compute the hessian and update the BFGS preconditioner at the new iterate*/
34       ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
35     }
36 
37     /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
38     ierr = TaoBNKComputeStep(tao, shift, &ksp_reason);CHKERRQ(ierr);
39 
40     /* Store current solution before it changes */
41     oldTrust = tao->trust;
42     bnk->fold = bnk->f;
43     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
44     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
45     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
46 
47     /* Temporarily accept the step and project it into the bounds */
48     ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr);
49     ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
50 
51     /* Check if the projection changed the step direction */
52     ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
53     ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
54     ierr = VecNorm(tao->stepdirection, NORM_2, &stepNorm);CHKERRQ(ierr);
55     if (stepNorm != bnk->dnorm) {
56       /* Projection changed the step, so we have to recompute predicted reduction.
57          However, we deliberately do not change the step norm and the trust radius
58          in order for the safeguard to more closely mimic a piece-wise linesearch
59          along the bounds. */
60       ierr = MatMult(bnk->H_inactive, tao->stepdirection, bnk->Xwork);CHKERRQ(ierr);
61       ierr = VecAYPX(bnk->Xwork, -0.5, tao->gradient);CHKERRQ(ierr);
62       ierr = VecDot(bnk->Xwork, tao->stepdirection, &prered);
63     } else {
64       /* Step did not change, so we can just recover the pre-computed prediction */
65       ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr);
66     }
67     prered = -prered;
68 
69     /* Compute the actual reduction and update the trust radius */
70     ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr);
71     actred = bnk->fold - bnk->f;
72     ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr);
73 
74     if (stepAccepted) {
75       /* Step is good, evaluate the gradient and the hessian */
76       steplen = 1.0;
77       ++bnk->newt;
78       ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
79       ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
80     } else {
81       /* Step is bad, revert old solution and re-solve with new radius*/
82       steplen = 0.0;
83       bnk->f = bnk->fold;
84       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
85       ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
86       ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
87       if (oldTrust == tao->trust) {
88         /* Can't change the radius anymore so just terminate */
89         tao->reason = TAO_DIVERGED_TR_REDUCTION;
90       }
91     }
92 
93     /*  Check for termination */
94     ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
95     if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
96     ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
97     ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,steplen);CHKERRQ(ierr);
98     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
99   }
100   PetscFunctionReturn(0);
101 }
102 
103 /*------------------------------------------------------------*/
104 
105 PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao)
106 {
107   TAO_BNK        *bnk;
108   PetscErrorCode ierr;
109 
110   PetscFunctionBegin;
111   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
112   tao->ops->solve=TaoSolve_BNTR;
113 
114   bnk = (TAO_BNK *)tao->data;
115   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
116   bnk->sval = 0.0; /* disable Hessian shifting */
117   PetscFunctionReturn(0);
118 }