1c14b763aSAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h> 2c14b763aSAlp Dener #include <petscksp.h> 3c14b763aSAlp Dener 4c14b763aSAlp Dener /* 5c14b763aSAlp Dener Implements Newton's Method with a trust region approach for solving 6c14b763aSAlp Dener bound constrained minimization problems. This version includes a 7c14b763aSAlp Dener line search fall-back in the event of a trust region failure. 8c14b763aSAlp Dener 9*df278d8fSAlp Dener The linear system solve has to be done with a conjugate gradient method. 10c14b763aSAlp Dener */ 11c14b763aSAlp Dener 12c14b763aSAlp Dener static PetscErrorCode TaoSolve_BNTL(Tao tao) 13c14b763aSAlp Dener { 14c14b763aSAlp Dener PetscErrorCode ierr; 15c14b763aSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 16c14b763aSAlp Dener TaoLineSearchConvergedReason ls_reason; 17c14b763aSAlp Dener 18c14b763aSAlp Dener PetscReal oldTrust, prered, actred, stepNorm, gdx, delta, steplen; 19c14b763aSAlp Dener PetscBool stepAccepted = PETSC_TRUE; 20c14b763aSAlp Dener PetscInt stepType, bfgsUpdates, updateType; 21c14b763aSAlp Dener 22c14b763aSAlp Dener PetscFunctionBegin; 23c14b763aSAlp Dener /* Project the current point onto the feasible set */ 24c14b763aSAlp Dener ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 25c14b763aSAlp Dener ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 26c14b763aSAlp Dener 27c14b763aSAlp Dener /* Project the initial point onto the feasible region */ 28c14b763aSAlp Dener ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 29c14b763aSAlp Dener 30c14b763aSAlp Dener /* Check convergence criteria */ 31c14b763aSAlp Dener ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr); 32c14b763aSAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 33c14b763aSAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 34c14b763aSAlp Dener if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 35c14b763aSAlp Dener 36c14b763aSAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 37c14b763aSAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 38c14b763aSAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,tao->trust);CHKERRQ(ierr); 39c14b763aSAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 40c14b763aSAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 41c14b763aSAlp Dener 42c14b763aSAlp Dener /* Initialize the preconditioner and trust radius */ 43c14b763aSAlp Dener ierr = TaoBNKInitialize(tao);CHKERRQ(ierr); 44c14b763aSAlp Dener 45c14b763aSAlp Dener /* Have not converged; continue with Newton method */ 46c14b763aSAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 47c14b763aSAlp Dener 48c14b763aSAlp Dener if (stepAccepted) { 49c14b763aSAlp Dener tao->niter++; 50c14b763aSAlp Dener tao->ksp_its=0; 51c14b763aSAlp Dener /* Compute the Hessian */ 52c14b763aSAlp Dener ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 53c14b763aSAlp Dener /* Update the BFGS preconditioner */ 54c14b763aSAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { 55c14b763aSAlp Dener if (BFGS_SCALE_PHESS == bnk->bfgs_scale_type) { 56c14b763aSAlp Dener /* Obtain diagonal for the bfgs preconditioner */ 57c14b763aSAlp Dener ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr); 58c14b763aSAlp Dener ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 59c14b763aSAlp Dener ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 60c14b763aSAlp Dener ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 61c14b763aSAlp Dener } 62c14b763aSAlp Dener /* Update the limited memory preconditioner and get existing # of updates */ 63c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 64c14b763aSAlp Dener } 65c14b763aSAlp Dener } 66c14b763aSAlp Dener 67c14b763aSAlp Dener /* Use the common BNK kernel to compute the raw Newton step */ 68c14b763aSAlp Dener ierr = TaoBNKComputeStep(tao, PETSC_FALSE, &stepType);CHKERRQ(ierr); 69c14b763aSAlp Dener 70c14b763aSAlp Dener /* Store current solution before it changes */ 71c14b763aSAlp Dener oldTrust = tao->trust; 72c14b763aSAlp Dener bnk->fold = bnk->f; 73c14b763aSAlp Dener ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 74c14b763aSAlp Dener ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr); 75c14b763aSAlp Dener ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr); 76c14b763aSAlp Dener 77c14b763aSAlp Dener /* Temporarily accept the step and project it into the bounds */ 78c14b763aSAlp Dener ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr); 79c14b763aSAlp Dener ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 80c14b763aSAlp Dener 81c14b763aSAlp Dener /* Check if the projection changed the step direction */ 82c14b763aSAlp Dener ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr); 83c14b763aSAlp Dener ierr = VecAXPBY(tao->stepdirection, -1.0, 1.0, bnk->Xold);CHKERRQ(ierr); 84c14b763aSAlp Dener ierr = VecNorm(tao->stepdirection, NORM_2, &stepNorm);CHKERRQ(ierr); 85c14b763aSAlp Dener if (stepNorm != bnk->dnorm) { 86c14b763aSAlp Dener /* Projection changed the step, so we have to adjust trust radius and recompute predicted reduction */ 87c14b763aSAlp Dener bnk->dnorm = stepNorm; 88c14b763aSAlp Dener tao->trust = bnk->dnorm; 89c14b763aSAlp Dener ierr = MatMult(tao->hessian, tao->stepdirection, bnk->Xwork);CHKERRQ(ierr); 90c14b763aSAlp Dener ierr = VecAYPX(bnk->Xwork, -0.5, tao->gradient);CHKERRQ(ierr); 91c14b763aSAlp Dener ierr = VecDot(bnk->Xwork, tao->stepdirection, &prered); 92c14b763aSAlp Dener } else { 93c14b763aSAlp Dener /* Step did not change, so we can just recover the pre-computed prediction */ 94c14b763aSAlp Dener ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr); 95c14b763aSAlp Dener } 96c14b763aSAlp Dener prered = -prered; 97c14b763aSAlp Dener 98c14b763aSAlp Dener /* Compute the actual reduction and update the trust radius */ 99c14b763aSAlp Dener ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr); 100c14b763aSAlp Dener actred = bnk->fold - bnk->f; 101c14b763aSAlp Dener ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, stepType, &stepAccepted);CHKERRQ(ierr); 102c14b763aSAlp Dener 103c14b763aSAlp Dener if (stepAccepted) { 104c14b763aSAlp Dener /* Step is good, evaluate the gradient and the hessian */ 105c14b763aSAlp Dener ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 106c14b763aSAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 107c14b763aSAlp Dener } else { 108c14b763aSAlp Dener /* Trust-region rejected the step. Revert the solution. */ 109c14b763aSAlp Dener bnk->f = bnk->fold; 110c14b763aSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 111c14b763aSAlp Dener 112c14b763aSAlp Dener /* Now check to make sure the Newton step is a descent direction... */ 113c14b763aSAlp Dener ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 114c14b763aSAlp Dener if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 115c14b763aSAlp Dener /* Newton step is not descent or direction produced Inf or NaN */ 116c14b763aSAlp Dener --bnk->newt; 117c14b763aSAlp Dener if (BNK_PC_BFGS != bnk->pc_type) { 118c14b763aSAlp Dener /* We don't have the BFGS matrix around and updated 119c14b763aSAlp Dener Must use gradient direction in this case */ 120c14b763aSAlp Dener ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 121c14b763aSAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 122c14b763aSAlp Dener ++bnk->grad; 123c14b763aSAlp Dener stepType = BNK_GRADIENT; 124c14b763aSAlp Dener } else { 125c14b763aSAlp Dener /* We have the BFGS matrix, so attempt to use the BFGS direction */ 126c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 127c14b763aSAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 128c14b763aSAlp Dener 129c14b763aSAlp Dener /* Check for success (descent direction) */ 130c14b763aSAlp Dener ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 131c14b763aSAlp Dener if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) { 132c14b763aSAlp Dener /* BFGS direction is not descent or direction produced not a number 133c14b763aSAlp Dener We can assert bfgsUpdates > 1 in this case because 134c14b763aSAlp Dener the first solve produces the scaled gradient direction, 135c14b763aSAlp Dener which is guaranteed to be descent */ 136c14b763aSAlp Dener 137c14b763aSAlp Dener /* Use steepest descent direction (scaled) */ 138c14b763aSAlp Dener if (bnk->f != 0.0) { 139c14b763aSAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 140c14b763aSAlp Dener } else { 141c14b763aSAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 142c14b763aSAlp Dener } 143c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 144c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 145c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, tao->gradient);CHKERRQ(ierr); 146c14b763aSAlp Dener ierr = MatLMVMSolve(bnk->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 147c14b763aSAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 148c14b763aSAlp Dener 149c14b763aSAlp Dener ++bnk->sgrad; 150c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 151c14b763aSAlp Dener } else { 152c14b763aSAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 153c14b763aSAlp Dener if (1 == bfgsUpdates) { 154c14b763aSAlp Dener /* The first BFGS direction is always the scaled gradient */ 155c14b763aSAlp Dener ++bnk->sgrad; 156c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 157c14b763aSAlp Dener } else { 158c14b763aSAlp Dener ++bnk->bfgs; 159c14b763aSAlp Dener stepType = BNK_BFGS; 160c14b763aSAlp Dener } 161c14b763aSAlp Dener } 162c14b763aSAlp Dener } 163c14b763aSAlp Dener } 164c14b763aSAlp Dener 165c14b763aSAlp Dener /* Trigger the line search */ 166c14b763aSAlp Dener ierr = TaoBNKPerformLineSearch(tao, stepType, &steplen, &ls_reason);CHKERRQ(ierr); 167c14b763aSAlp Dener if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 168c14b763aSAlp Dener /* Line search failed, revert solution and terminate */ 169c14b763aSAlp Dener bnk->f = bnk->fold; 170c14b763aSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 171c14b763aSAlp Dener ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 172c14b763aSAlp Dener ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 173c14b763aSAlp Dener tao->trust = 0.0; 174c14b763aSAlp Dener tao->reason = TAO_DIVERGED_LS_FAILURE; 175c14b763aSAlp Dener } else { 176c14b763aSAlp Dener /* Line search succeeded so we should update the trust radius based on the LS step length */ 177c14b763aSAlp Dener updateType = bnk->update_type; 178c14b763aSAlp Dener bnk->update_type = BNK_UPDATE_STEP; 179c14b763aSAlp Dener ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, stepType, &stepAccepted);CHKERRQ(ierr); 180c14b763aSAlp Dener bnk->update_type = updateType; 181c14b763aSAlp Dener } 182c14b763aSAlp Dener } 183c14b763aSAlp Dener 184c14b763aSAlp Dener /* Check for termination */ 185c14b763aSAlp Dener ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr); 186c14b763aSAlp Dener if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 187c14b763aSAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 188c14b763aSAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,tao->trust);CHKERRQ(ierr); 189c14b763aSAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 190c14b763aSAlp Dener } 191c14b763aSAlp Dener PetscFunctionReturn(0); 192c14b763aSAlp Dener } 193c14b763aSAlp Dener 194*df278d8fSAlp Dener /*------------------------------------------------------------*/ 195*df278d8fSAlp Dener 196c14b763aSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNTL(Tao tao) 197c14b763aSAlp Dener { 198c14b763aSAlp Dener TAO_BNK *bnk; 199c14b763aSAlp Dener PetscErrorCode ierr; 200c14b763aSAlp Dener 201c14b763aSAlp Dener PetscFunctionBegin; 202c14b763aSAlp Dener ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 203c14b763aSAlp Dener tao->ops->solve=TaoSolve_BNTL; 204c14b763aSAlp Dener 205c14b763aSAlp Dener bnk = (TAO_BNK *)tao->data; 206c14b763aSAlp Dener bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */ 207c14b763aSAlp Dener bnk->sval = 0.0; /* disable Hessian shifting */ 208c14b763aSAlp Dener PetscFunctionReturn(0); 209c14b763aSAlp Dener }