xref: /petsc/src/tao/bound/impls/bnk/bntr.c (revision 66ed3702dc8f7c8e3a7cbd4d5287a373dbf44efa)
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 
17fed79b8eSAlp Dener   PetscReal                    oldTrust;
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) {
46*66ed3702SAlp 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 
76fed79b8eSAlp Dener     /* Test the new step for acceptance */
77fed79b8eSAlp Dener     ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr);
78fed79b8eSAlp Dener     ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr);
79fed79b8eSAlp Dener     ierr = TaoBNKUpdateTrustRadius(tao, bnk->fold, bnk->f, stepType, &stepAccepted);CHKERRQ(ierr);
80fed79b8eSAlp Dener 
81fed79b8eSAlp Dener     if (stepAccepted) {
82*66ed3702SAlp Dener       /* Step is good, evaluate the gradient and the hessian */
83fed79b8eSAlp Dener       ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
84fed79b8eSAlp Dener       ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
85fed79b8eSAlp Dener     } else {
86fed79b8eSAlp Dener       /* Step is bad, revert old solution and re-solve with new radius*/
87fed79b8eSAlp Dener       bnk->f = bnk->fold;
88fed79b8eSAlp Dener       ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
89fed79b8eSAlp Dener       ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
90fed79b8eSAlp Dener       ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
91*66ed3702SAlp Dener       if (oldTrust == tao->trust == bnk->min_radius) {
92fed79b8eSAlp Dener         /* Can't shrink trust radius any further, so we have to terminate */
93fed79b8eSAlp Dener         tao->reason = TAO_DIVERGED_TR_REDUCTION;
94fed79b8eSAlp Dener       }
95fed79b8eSAlp Dener     }
96fed79b8eSAlp Dener 
97fed79b8eSAlp Dener     /*  Check for termination */
98fed79b8eSAlp Dener     ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
99fed79b8eSAlp Dener     if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number");
100fed79b8eSAlp Dener     ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
101fed79b8eSAlp Dener     ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,tao->trust);CHKERRQ(ierr);
102fed79b8eSAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
103fed79b8eSAlp Dener   }
104fed79b8eSAlp Dener   PetscFunctionReturn(0);
105fed79b8eSAlp Dener }
106fed79b8eSAlp Dener 
107fed79b8eSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao)
108fed79b8eSAlp Dener {
109fed79b8eSAlp Dener   TAO_BNK        *bnk;
110fed79b8eSAlp Dener   PetscErrorCode ierr;
111fed79b8eSAlp Dener 
112fed79b8eSAlp Dener   PetscFunctionBegin;
113fed79b8eSAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
114fed79b8eSAlp Dener   tao->ops->solve=TaoSolve_BNTR;
115fed79b8eSAlp Dener 
116fed79b8eSAlp Dener   bnk = (TAO_BNK *)tao->data;
117*66ed3702SAlp Dener   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
118*66ed3702SAlp Dener   bnk->sval = 0.0; /* disable Hessian shifting */
119fed79b8eSAlp Dener   PetscFunctionReturn(0);
120fed79b8eSAlp Dener }