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. This version includes a 7 line search fall-back in the event of a trust region failure. 8 9 The linear system solve has to be done with a conjugate gradient method. 10 */ 11 12 static PetscErrorCode TaoSolve_BNTL(Tao tao) 13 { 14 PetscErrorCode ierr; 15 TAO_BNK *bnk = (TAO_BNK *)tao->data; 16 KSPConvergedReason ksp_reason; 17 TaoLineSearchConvergedReason ls_reason; 18 19 PetscReal oldTrust, prered, actred, stepNorm, steplen; 20 PetscBool stepAccepted = PETSC_TRUE; 21 PetscInt stepType, updateType; 22 23 PetscFunctionBegin; 24 /* Project the current point onto the feasible set */ 25 ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 26 ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 27 28 /* Project the initial point onto the feasible region */ 29 ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 30 31 /* Check convergence criteria */ 32 ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr); 33 ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 34 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 35 if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 36 37 tao->reason = TAO_CONTINUE_ITERATING; 38 ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 39 ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,tao->trust);CHKERRQ(ierr); 40 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 41 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 42 43 /* Initialize the preconditioner and trust radius */ 44 ierr = TaoBNKInitialize(tao);CHKERRQ(ierr); 45 46 /* Have not converged; continue with Newton method */ 47 while (tao->reason == TAO_CONTINUE_ITERATING) { 48 tao->niter++; 49 tao->ksp_its=0; 50 /* Compute the Hessian */ 51 ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 52 /* Update the BFGS preconditioner */ 53 if (BNK_PC_BFGS == bnk->pc_type) { 54 if (BFGS_SCALE_PHESS == bnk->bfgs_scale_type) { 55 /* Obtain diagonal for the bfgs preconditioner */ 56 ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr); 57 ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 58 ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 59 ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 60 } 61 /* Update the limited memory preconditioner and get existing # of updates */ 62 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 63 } 64 65 /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */ 66 ierr = TaoBNKComputeStep(tao, &ksp_reason);CHKERRQ(ierr); 67 68 /* Store current solution before it changes */ 69 oldTrust = tao->trust; 70 bnk->fold = bnk->f; 71 ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 72 ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr); 73 ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr); 74 75 /* Temporarily accept the step and project it into the bounds */ 76 ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr); 77 ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 78 79 /* Check if the projection changed the step direction */ 80 ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr); 81 ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr); 82 ierr = VecNorm(tao->stepdirection, NORM_2, &stepNorm);CHKERRQ(ierr); 83 if (stepNorm != bnk->dnorm) { 84 /* Projection changed the step, so we have to recompute predicted reduction. 85 However, we deliberately do not change the step norm and the trust radius 86 in order for the safeguard to more closely mimic a piece-wise linesearch 87 along the bounds. */ 88 ierr = MatMult(tao->hessian, tao->stepdirection, bnk->Xwork);CHKERRQ(ierr); 89 ierr = VecAYPX(bnk->Xwork, -0.5, tao->gradient);CHKERRQ(ierr); 90 ierr = VecDot(bnk->Xwork, tao->stepdirection, &prered); 91 } else { 92 /* Step did not change, so we can just recover the pre-computed prediction */ 93 ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr); 94 } 95 prered = -prered; 96 97 /* Compute the actual reduction and update the trust radius */ 98 ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr); 99 actred = bnk->fold - bnk->f; 100 ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, stepType, &stepAccepted);CHKERRQ(ierr); 101 102 if (stepAccepted) { 103 /* Step is good, evaluate the gradient and the hessian */ 104 steplen = 1.0; 105 ++bnk->newt; 106 ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 107 ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 108 } else { 109 /* Trust-region rejected the step. Revert the solution. */ 110 bnk->f = bnk->fold; 111 ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 112 113 /* Trigger the line search */ 114 ierr = TaoBNKSafeguardStep(tao, ksp_reason, &stepType);CHKERRQ(ierr); 115 ierr = TaoBNKPerformLineSearch(tao, stepType, &steplen, &ls_reason);CHKERRQ(ierr); 116 if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 117 /* Line search failed, revert solution and terminate */ 118 bnk->f = bnk->fold; 119 ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 120 ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 121 ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 122 tao->trust = 0.0; 123 tao->reason = TAO_DIVERGED_LS_FAILURE; 124 } else { 125 /* Line search succeeded so we should update the trust radius based on the LS step length */ 126 updateType = bnk->update_type; 127 bnk->update_type = BNK_UPDATE_STEP; 128 tao->trust = oldTrust; 129 ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, stepType, &stepAccepted);CHKERRQ(ierr); 130 bnk->update_type = updateType; 131 /* count the accepted step types */ 132 switch (stepType) { 133 case BNK_NEWTON: 134 ++bnk->newt; 135 break; 136 case BNK_BFGS: 137 ++bnk->bfgs; 138 break; 139 case BNK_SCALED_GRADIENT: 140 ++bnk->sgrad; 141 break; 142 case BNK_GRADIENT: 143 ++bnk->grad; 144 break; 145 default: 146 break; 147 } 148 } 149 } 150 151 /* Check for termination */ 152 ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr); 153 if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 154 ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 155 ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,steplen);CHKERRQ(ierr); 156 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 157 } 158 PetscFunctionReturn(0); 159 } 160 161 /*------------------------------------------------------------*/ 162 163 PETSC_EXTERN PetscErrorCode TaoCreate_BNTL(Tao tao) 164 { 165 TAO_BNK *bnk; 166 PetscErrorCode ierr; 167 168 PetscFunctionBegin; 169 ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 170 tao->ops->solve=TaoSolve_BNTL; 171 172 bnk = (TAO_BNK *)tao->data; 173 bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */ 174 bnk->sval = 0.0; /* disable Hessian shifting */ 175 PetscFunctionReturn(0); 176 }