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