xref: /petsc/src/tao/bound/impls/bnk/bntr.c (revision 73e4db900df47a59ed4aa6fd03eeb9faefcb6f8a)
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 
8fed79b8eSAlp Dener  The linear system solve should be done with a conjugate gradient
9fed79b8eSAlp Dener  method, although any method can be used.
10fed79b8eSAlp Dener */
11fed79b8eSAlp Dener 
12fed79b8eSAlp Dener static PetscErrorCode TaoSolve_BNTR(Tao tao)
13fed79b8eSAlp Dener {
14fed79b8eSAlp Dener   PetscErrorCode               ierr;
15fed79b8eSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
16fed79b8eSAlp Dener 
17b1c2d0e3SAlp Dener   PetscReal                    oldTrust, prered, actred, stepNorm;
18fed79b8eSAlp Dener   PetscBool                    stepAccepted = PETSC_TRUE;
19fed79b8eSAlp Dener   PetscInt                     stepType;
20fed79b8eSAlp Dener 
21fed79b8eSAlp Dener   PetscFunctionBegin;
22fed79b8eSAlp Dener   /*   Project the current point onto the feasible set */
23fed79b8eSAlp Dener   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
24fed79b8eSAlp Dener   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
25fed79b8eSAlp Dener 
26fed79b8eSAlp Dener   /* Project the initial point onto the feasible region */
27fed79b8eSAlp Dener   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
28fed79b8eSAlp Dener 
29fed79b8eSAlp Dener   /* Check convergence criteria */
30fed79b8eSAlp Dener   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr);
31fed79b8eSAlp Dener   ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
32fed79b8eSAlp Dener   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
33fed79b8eSAlp Dener   if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
34fed79b8eSAlp Dener 
35fed79b8eSAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
36fed79b8eSAlp Dener   ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
37fed79b8eSAlp Dener   ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,tao->trust);CHKERRQ(ierr);
38fed79b8eSAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
39fed79b8eSAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
40fed79b8eSAlp Dener 
41fed79b8eSAlp Dener   /* Initialize the preconditioner and trust radius */
42fed79b8eSAlp Dener   ierr = TaoBNKInitialize(tao);CHKERRQ(ierr);
43fed79b8eSAlp Dener 
44fed79b8eSAlp Dener   /* Have not converged; continue with Newton method */
45fed79b8eSAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
4666ed3702SAlp Dener 
47fed79b8eSAlp Dener     if (stepAccepted) {
48fed79b8eSAlp Dener       tao->niter++;
49fed79b8eSAlp Dener       tao->ksp_its=0;
50fed79b8eSAlp Dener       /* Compute the Hessian */
51fed79b8eSAlp Dener       ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
52fed79b8eSAlp Dener       /* Update the BFGS preconditioner */
53fed79b8eSAlp Dener       if (BNK_PC_BFGS == bnk->pc_type) {
54fed79b8eSAlp Dener         if (BFGS_SCALE_PHESS == bnk->bfgs_scale_type) {
55fed79b8eSAlp Dener           /* Obtain diagonal for the bfgs preconditioner  */
56fed79b8eSAlp Dener           ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr);
57fed79b8eSAlp Dener           ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
58fed79b8eSAlp Dener           ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
59fed79b8eSAlp Dener           ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
60fed79b8eSAlp Dener         }
61fed79b8eSAlp Dener         /* Update the limited memory preconditioner and get existing # of updates */
62fed79b8eSAlp Dener         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
63fed79b8eSAlp Dener       }
64fed79b8eSAlp Dener     }
65fed79b8eSAlp Dener 
66fed79b8eSAlp Dener     /* Use the common BNK kernel to compute the raw Newton step */
67fed79b8eSAlp Dener     ierr = TaoBNKComputeStep(tao, PETSC_FALSE, &stepType);CHKERRQ(ierr);
68fed79b8eSAlp Dener 
69fed79b8eSAlp Dener     /* Store current solution before it changes */
70fed79b8eSAlp Dener     oldTrust = tao->trust;
71fed79b8eSAlp Dener     bnk->fold = bnk->f;
72fed79b8eSAlp Dener     ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
73fed79b8eSAlp Dener     ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr);
74fed79b8eSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr);
75fed79b8eSAlp Dener 
76b1c2d0e3SAlp Dener     /* Temporarily accept the step and project it into the bounds */
77fed79b8eSAlp Dener     ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr);
78b1c2d0e3SAlp Dener     ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
79b1c2d0e3SAlp Dener 
80b1c2d0e3SAlp Dener     /* Check if the projection changed the step direction */
81b1c2d0e3SAlp Dener     ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
82b1c2d0e3SAlp Dener     ierr = VecAXPBY(tao->stepdirection, -1.0, 1.0, bnk->Xold);CHKERRQ(ierr);
83b1c2d0e3SAlp Dener     ierr = VecNorm(tao->stepdirection, NORM_2, &stepNorm);CHKERRQ(ierr);
84b1c2d0e3SAlp Dener     if (stepNorm != bnk->dnorm) {
85b1c2d0e3SAlp Dener       /* Projection changed the step, so we have to adjust trust radius and recompute predicted reduction */
86b1c2d0e3SAlp Dener       bnk->dnorm = stepNorm;
87b1c2d0e3SAlp Dener       tao->trust = bnk->dnorm;
88b1c2d0e3SAlp Dener       ierr = MatMult(tao->hessian, tao->stepdirection, bnk->Xwork);CHKERRQ(ierr);
89b1c2d0e3SAlp Dener       ierr = VecAYPX(bnk->Xwork, -0.5, tao->gradient);CHKERRQ(ierr);
90b1c2d0e3SAlp Dener       ierr = VecDot(bnk->Xwork, tao->stepdirection, &prered);
91b1c2d0e3SAlp Dener     } else {
92b1c2d0e3SAlp Dener       /* Step did not change, so we can just recover the pre-computed prediction */
93b1c2d0e3SAlp Dener       ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr);
94b1c2d0e3SAlp Dener     }
95b1c2d0e3SAlp Dener     prered = -prered;
96b1c2d0e3SAlp Dener 
97b1c2d0e3SAlp Dener     /* Compute the actual reduction and update the trust radius */
98fed79b8eSAlp Dener     ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr);
99b1c2d0e3SAlp Dener     actred = bnk->fold - bnk->f;
100b1c2d0e3SAlp Dener     ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, stepType, &stepAccepted);CHKERRQ(ierr);
101fed79b8eSAlp Dener 
102fed79b8eSAlp Dener     if (stepAccepted) {
10366ed3702SAlp Dener       /* Step is good, evaluate the gradient and the hessian */
104fed79b8eSAlp Dener       ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
105fed79b8eSAlp Dener       ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
106fed79b8eSAlp Dener     } else {
107fed79b8eSAlp Dener       /* Step is bad, revert old solution and re-solve with new radius*/
108fed79b8eSAlp Dener       bnk->f = bnk->fold;
109fed79b8eSAlp Dener       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
110fed79b8eSAlp Dener       ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
111fed79b8eSAlp Dener       ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
112*73e4db90SAlp Dener       if (oldTrust == tao->trust) {
113*73e4db90SAlp Dener         /* Can't change the radius anymore so just terminate */
114fed79b8eSAlp Dener         tao->reason = TAO_DIVERGED_TR_REDUCTION;
115fed79b8eSAlp Dener       }
116fed79b8eSAlp Dener     }
117fed79b8eSAlp Dener 
118fed79b8eSAlp Dener     /*  Check for termination */
119fed79b8eSAlp Dener     ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
120*73e4db90SAlp Dener     if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
121fed79b8eSAlp Dener     ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
122fed79b8eSAlp Dener     ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,tao->trust);CHKERRQ(ierr);
123fed79b8eSAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
124fed79b8eSAlp Dener   }
125fed79b8eSAlp Dener   PetscFunctionReturn(0);
126fed79b8eSAlp Dener }
127fed79b8eSAlp Dener 
128fed79b8eSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao)
129fed79b8eSAlp Dener {
130fed79b8eSAlp Dener   TAO_BNK        *bnk;
131fed79b8eSAlp Dener   PetscErrorCode ierr;
132fed79b8eSAlp Dener 
133fed79b8eSAlp Dener   PetscFunctionBegin;
134fed79b8eSAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
135fed79b8eSAlp Dener   tao->ops->solve=TaoSolve_BNTR;
136fed79b8eSAlp Dener 
137fed79b8eSAlp Dener   bnk = (TAO_BNK *)tao->data;
13866ed3702SAlp Dener   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
13966ed3702SAlp Dener   bnk->sval = 0.0; /* disable Hessian shifting */
140fed79b8eSAlp Dener   PetscFunctionReturn(0);
141fed79b8eSAlp Dener }