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 16 PetscReal oldTrust, prered, actred, stepNorm, steplen; 17 PetscBool stepAccepted = PETSC_TRUE; 18 PetscInt stepType; 19 20 PetscFunctionBegin; 21 /* Project the current point onto the feasible set */ 22 ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 23 ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 24 25 /* Project the initial point onto the feasible region */ 26 ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 27 28 /* Check convergence criteria */ 29 ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr); 30 ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 31 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 32 if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 33 34 tao->reason = TAO_CONTINUE_ITERATING; 35 ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 36 ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,tao->trust);CHKERRQ(ierr); 37 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 38 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 39 40 /* Initialize the preconditioner and trust radius */ 41 ierr = TaoBNKInitialize(tao);CHKERRQ(ierr); 42 43 /* Have not converged; continue with Newton method */ 44 while (tao->reason == TAO_CONTINUE_ITERATING) { 45 46 if (stepAccepted) { 47 tao->niter++; 48 tao->ksp_its=0; 49 /* Compute the Hessian */ 50 ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 51 /* Update the BFGS preconditioner */ 52 if (BNK_PC_BFGS == bnk->pc_type) { 53 if (BFGS_SCALE_PHESS == bnk->bfgs_scale_type) { 54 /* Obtain diagonal for the bfgs preconditioner */ 55 ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr); 56 ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 57 ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 58 ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 59 } 60 /* Update the limited memory preconditioner and get existing # of updates */ 61 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 62 } 63 } 64 65 /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */ 66 ierr = TaoBNKComputeStep(tao, PETSC_FALSE, &stepType);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 ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 106 ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 107 } else { 108 /* Step is bad, revert old solution and re-solve with new radius*/ 109 steplen = 0.0; 110 bnk->f = bnk->fold; 111 ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 112 ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 113 ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 114 if (oldTrust == tao->trust) { 115 /* Can't change the radius anymore so just terminate */ 116 tao->reason = TAO_DIVERGED_TR_REDUCTION; 117 } 118 } 119 120 /* Check for termination */ 121 ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr); 122 if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 123 ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 124 ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,steplen);CHKERRQ(ierr); 125 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 126 } 127 PetscFunctionReturn(0); 128 } 129 130 /*------------------------------------------------------------*/ 131 132 PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao) 133 { 134 TAO_BNK *bnk; 135 PetscErrorCode ierr; 136 137 PetscFunctionBegin; 138 ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 139 tao->ops->solve=TaoSolve_BNTR; 140 141 bnk = (TAO_BNK *)tao->data; 142 bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */ 143 bnk->sval = 0.0; /* disable Hessian shifting */ 144 PetscFunctionReturn(0); 145 }