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; 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);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 */ 34 ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 35 /* Update the BFGS preconditioner */ 36 if (BNK_PC_BFGS == bnk->pc_type) { 37 if (BFGS_SCALE_PHESS == bnk->bfgs_scale_type) { 38 /* Obtain diagonal for the bfgs preconditioner */ 39 ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr); 40 ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 41 ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 42 ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 43 } 44 /* Update the limited memory preconditioner and get existing # of updates */ 45 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 46 } 47 } 48 49 /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */ 50 ierr = TaoBNKComputeStep(tao, &ksp_reason);CHKERRQ(ierr); 51 52 /* Store current solution before it changes */ 53 oldTrust = tao->trust; 54 bnk->fold = bnk->f; 55 ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 56 ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr); 57 ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr); 58 59 /* Temporarily accept the step and project it into the bounds */ 60 ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr); 61 ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 62 63 /* Check if the projection changed the step direction */ 64 ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr); 65 ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr); 66 ierr = VecNorm(tao->stepdirection, NORM_2, &stepNorm);CHKERRQ(ierr); 67 if (stepNorm != bnk->dnorm) { 68 /* Projection changed the step, so we have to recompute predicted reduction. 69 However, we deliberately do not change the step norm and the trust radius 70 in order for the safeguard to more closely mimic a piece-wise linesearch 71 along the bounds. */ 72 ierr = MatMult(bnk->H_inactive, tao->stepdirection, bnk->Xwork);CHKERRQ(ierr); 73 ierr = VecAYPX(bnk->Xwork, -0.5, tao->gradient);CHKERRQ(ierr); 74 ierr = VecDot(bnk->Xwork, tao->stepdirection, &prered); 75 } else { 76 /* Step did not change, so we can just recover the pre-computed prediction */ 77 ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr); 78 } 79 prered = -prered; 80 81 /* Compute the actual reduction and update the trust radius */ 82 ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr); 83 actred = bnk->fold - bnk->f; 84 ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr); 85 86 if (stepAccepted) { 87 /* Step is good, evaluate the gradient and the hessian */ 88 steplen = 1.0; 89 ++bnk->newt; 90 ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 91 ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 92 } else { 93 /* Step is bad, revert old solution and re-solve with new radius*/ 94 steplen = 0.0; 95 bnk->f = bnk->fold; 96 ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 97 ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 98 ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 99 if (oldTrust == tao->trust) { 100 /* Can't change the radius anymore so just terminate */ 101 tao->reason = TAO_DIVERGED_TR_REDUCTION; 102 } 103 } 104 105 /* Check for termination */ 106 ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr); 107 if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 108 ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 109 ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,steplen);CHKERRQ(ierr); 110 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 111 } 112 PetscFunctionReturn(0); 113 } 114 115 /*------------------------------------------------------------*/ 116 117 PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao) 118 { 119 TAO_BNK *bnk; 120 PetscErrorCode ierr; 121 122 PetscFunctionBegin; 123 ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 124 tao->ops->solve=TaoSolve_BNTR; 125 126 bnk = (TAO_BNK *)tao->data; 127 bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */ 128 bnk->sval = 0.0; /* disable Hessian shifting */ 129 PetscFunctionReturn(0); 130 }