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 9df278d8fSAlp 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 tao->niter++; 48c14b763aSAlp Dener tao->ksp_its=0; 49c14b763aSAlp Dener /* Compute the Hessian */ 50c14b763aSAlp Dener ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 51c14b763aSAlp Dener /* Update the BFGS preconditioner */ 52c14b763aSAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { 53c14b763aSAlp Dener if (BFGS_SCALE_PHESS == bnk->bfgs_scale_type) { 54c14b763aSAlp Dener /* Obtain diagonal for the bfgs preconditioner */ 55c14b763aSAlp Dener ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr); 56c14b763aSAlp Dener ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 57c14b763aSAlp Dener ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 58c14b763aSAlp Dener ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 59c14b763aSAlp Dener } 60c14b763aSAlp Dener /* Update the limited memory preconditioner and get existing # of updates */ 61c14b763aSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 62c14b763aSAlp Dener } 63c14b763aSAlp Dener 648d5ead36SAlp Dener /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */ 65c14b763aSAlp Dener ierr = TaoBNKComputeStep(tao, PETSC_FALSE, &stepType);CHKERRQ(ierr); 66c14b763aSAlp Dener 67c14b763aSAlp Dener /* Store current solution before it changes */ 68c14b763aSAlp Dener oldTrust = tao->trust; 69c14b763aSAlp Dener bnk->fold = bnk->f; 70c14b763aSAlp Dener ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 71c14b763aSAlp Dener ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr); 72c14b763aSAlp Dener ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr); 73c14b763aSAlp Dener 74c14b763aSAlp Dener /* Temporarily accept the step and project it into the bounds */ 75c14b763aSAlp Dener ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr); 76c14b763aSAlp Dener ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 77c14b763aSAlp Dener 78c14b763aSAlp Dener /* Check if the projection changed the step direction */ 79c14b763aSAlp Dener ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr); 808d5ead36SAlp Dener ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr); 81c14b763aSAlp Dener ierr = VecNorm(tao->stepdirection, NORM_2, &stepNorm);CHKERRQ(ierr); 82c14b763aSAlp Dener if (stepNorm != bnk->dnorm) { 838d5ead36SAlp Dener /* Projection changed the step, so we have to recompute predicted reduction. 848d5ead36SAlp Dener However, we deliberately do not change the step norm and the trust radius 858d5ead36SAlp Dener in order for the safeguard to more closely mimic a piece-wise linesearch 868d5ead36SAlp Dener along the bounds. */ 87c14b763aSAlp Dener ierr = MatMult(tao->hessian, tao->stepdirection, bnk->Xwork);CHKERRQ(ierr); 88c14b763aSAlp Dener ierr = VecAYPX(bnk->Xwork, -0.5, tao->gradient);CHKERRQ(ierr); 89c14b763aSAlp Dener ierr = VecDot(bnk->Xwork, tao->stepdirection, &prered); 90c14b763aSAlp Dener } else { 91c14b763aSAlp Dener /* Step did not change, so we can just recover the pre-computed prediction */ 92c14b763aSAlp Dener ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr); 93c14b763aSAlp Dener } 94c14b763aSAlp Dener prered = -prered; 95c14b763aSAlp Dener 96c14b763aSAlp Dener /* Compute the actual reduction and update the trust radius */ 97c14b763aSAlp Dener ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr); 98c14b763aSAlp Dener actred = bnk->fold - bnk->f; 99c14b763aSAlp Dener ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, stepType, &stepAccepted);CHKERRQ(ierr); 100c14b763aSAlp Dener 101c14b763aSAlp Dener if (stepAccepted) { 102c14b763aSAlp Dener /* Step is good, evaluate the gradient and the hessian */ 1038d5ead36SAlp Dener steplen = 1.0; 104c14b763aSAlp Dener ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 105c14b763aSAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 106c14b763aSAlp Dener } else { 107c14b763aSAlp Dener /* Trust-region rejected the step. Revert the solution. */ 108c14b763aSAlp Dener bnk->f = bnk->fold; 109c14b763aSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 110c14b763aSAlp Dener 111c14b763aSAlp Dener /* Now check to make sure the Newton step is a descent direction... */ 112c14b763aSAlp Dener ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 113c14b763aSAlp Dener if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) { 114c14b763aSAlp Dener /* Newton step is not descent or direction produced Inf or NaN */ 115c14b763aSAlp Dener --bnk->newt; 116c14b763aSAlp Dener if (BNK_PC_BFGS != bnk->pc_type) { 117c14b763aSAlp Dener /* We don't have the BFGS matrix around and updated 118c14b763aSAlp Dener Must use gradient direction in this case */ 119c14b763aSAlp Dener ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 120c14b763aSAlp Dener ++bnk->grad; 121c14b763aSAlp Dener stepType = BNK_GRADIENT; 122c14b763aSAlp Dener } else { 123c14b763aSAlp Dener /* We have the BFGS matrix, so attempt to use the BFGS direction */ 124a41f356dSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 125c14b763aSAlp Dener 1268d5ead36SAlp Dener /* Check for success (descent direction) 1278d5ead36SAlp Dener NOTE: Negative gdx here means not a descent direction because 1288d5ead36SAlp Dener the fall-back step is missing a negative sign. */ 129c14b763aSAlp Dener ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr); 1308d5ead36SAlp Dener if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 131c14b763aSAlp Dener /* BFGS direction is not descent or direction produced not a number 132c14b763aSAlp Dener We can assert bfgsUpdates > 1 in this case because 133c14b763aSAlp Dener the first solve produces the scaled gradient direction, 134c14b763aSAlp Dener which is guaranteed to be descent */ 135c14b763aSAlp Dener 136c14b763aSAlp Dener /* Use steepest descent direction (scaled) */ 137c14b763aSAlp Dener if (bnk->f != 0.0) { 138c14b763aSAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 139c14b763aSAlp Dener } else { 140c14b763aSAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 141c14b763aSAlp Dener } 142c14b763aSAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 143c14b763aSAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 144a41f356dSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 145a41f356dSAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 146c14b763aSAlp Dener 147c14b763aSAlp Dener ++bnk->sgrad; 148c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 149c14b763aSAlp Dener } else { 150c14b763aSAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 151c14b763aSAlp Dener if (1 == bfgsUpdates) { 152c14b763aSAlp Dener /* The first BFGS direction is always the scaled gradient */ 153c14b763aSAlp Dener ++bnk->sgrad; 154c14b763aSAlp Dener stepType = BNK_SCALED_GRADIENT; 155c14b763aSAlp Dener } else { 156c14b763aSAlp Dener ++bnk->bfgs; 157c14b763aSAlp Dener stepType = BNK_BFGS; 158c14b763aSAlp Dener } 159c14b763aSAlp Dener } 160c14b763aSAlp Dener } 161c14b763aSAlp Dener } 162*770b7498SAlp Dener /* Make sure the safeguarded fall-back step is zero for actively bounded variables */ 163*770b7498SAlp Dener ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 1648d5ead36SAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 165c14b763aSAlp Dener 166c14b763aSAlp Dener /* Trigger the line search */ 167c14b763aSAlp Dener ierr = TaoBNKPerformLineSearch(tao, stepType, &steplen, &ls_reason);CHKERRQ(ierr); 168c14b763aSAlp Dener if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 169c14b763aSAlp Dener /* Line search failed, revert solution and terminate */ 170c14b763aSAlp Dener bnk->f = bnk->fold; 171c14b763aSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 172c14b763aSAlp Dener ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 173c14b763aSAlp Dener ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 174c14b763aSAlp Dener tao->trust = 0.0; 175c14b763aSAlp Dener tao->reason = TAO_DIVERGED_LS_FAILURE; 176c14b763aSAlp Dener } else { 177c14b763aSAlp Dener /* Line search succeeded so we should update the trust radius based on the LS step length */ 178c14b763aSAlp Dener updateType = bnk->update_type; 179c14b763aSAlp Dener bnk->update_type = BNK_UPDATE_STEP; 180*770b7498SAlp Dener tao->trust = oldTrust; 181c14b763aSAlp Dener ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, stepType, &stepAccepted);CHKERRQ(ierr); 182c14b763aSAlp Dener bnk->update_type = updateType; 183c14b763aSAlp Dener } 184c14b763aSAlp Dener } 185c14b763aSAlp Dener 186c14b763aSAlp Dener /* Check for termination */ 187c14b763aSAlp Dener ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr); 188c14b763aSAlp Dener if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 189c14b763aSAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 1908d5ead36SAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,steplen);CHKERRQ(ierr); 191c14b763aSAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 192c14b763aSAlp Dener } 193c14b763aSAlp Dener PetscFunctionReturn(0); 194c14b763aSAlp Dener } 195c14b763aSAlp Dener 196df278d8fSAlp Dener /*------------------------------------------------------------*/ 197df278d8fSAlp Dener 198c14b763aSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNTL(Tao tao) 199c14b763aSAlp Dener { 200c14b763aSAlp Dener TAO_BNK *bnk; 201c14b763aSAlp Dener PetscErrorCode ierr; 202c14b763aSAlp Dener 203c14b763aSAlp Dener PetscFunctionBegin; 204c14b763aSAlp Dener ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 205c14b763aSAlp Dener tao->ops->solve=TaoSolve_BNTL; 206c14b763aSAlp Dener 207c14b763aSAlp Dener bnk = (TAO_BNK *)tao->data; 208c14b763aSAlp Dener bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */ 209c14b763aSAlp Dener bnk->sval = 0.0; /* disable Hessian shifting */ 210c14b763aSAlp Dener PetscFunctionReturn(0); 211c14b763aSAlp Dener }