1 #include <../src/tao/bound/impls/bnk/bnk.h> 2 #include <petscksp.h> 3 4 /* 5 Implements Newton's Method with a line search approach for solving 6 bound constrained minimization problems. A projected More'-Thuente line 7 search is used to guarantee that the bfgs preconditioner remains positive 8 definite. 9 10 The method can shift the Hessian matrix. The shifting procedure is 11 adapted from the PATH algorithm for solving complementarity 12 problems. 13 14 The linear system solve should be done with a conjugate gradient 15 method, although any method can be used. 16 */ 17 18 static PetscErrorCode TaoSolve_BNLS(Tao tao) 19 { 20 PetscErrorCode ierr; 21 TAO_BNK *bnk = (TAO_BNK *)tao->data; 22 TaoLineSearchConvergedReason ls_reason; 23 24 PetscReal f_full, gdx; 25 PetscReal step = 1.0; 26 PetscReal delta; 27 PetscReal e_min; 28 29 PetscBool trustAccept; 30 PetscInt stepType; 31 PetscInt bfgsUpdates = 0; 32 PetscInt needH = 1; 33 34 PetscFunctionBegin; 35 /* Project the current point onto the feasible set */ 36 ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 37 ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 38 39 /* Project the initial point onto the feasible region */ 40 ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 41 42 /* Check convergence criteria */ 43 ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr); 44 ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 45 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 46 if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 47 48 tao->reason = TAO_CONTINUE_ITERATING; 49 ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 50 ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr); 51 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 52 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 53 54 /* Initialize the preconditioner and trust radius */ 55 ierr = TaoBNKInitialize(tao);CHKERRQ(ierr); 56 57 /* Have not converged; continue with Newton method */ 58 while (tao->reason == TAO_CONTINUE_ITERATING) { 59 ++tao->niter; 60 tao->ksp_its=0; 61 62 /* Compute the Hessian */ 63 ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 64 65 /* Update the BFGS preconditioner */ 66 if (BNK_PC_BFGS == bnk->pc_type) { 67 if (BFGS_SCALE_PHESS == bnk->bfgs_scale_type) { 68 /* Obtain diagonal for the bfgs preconditioner */ 69 ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr); 70 ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 71 ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 72 ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 73 } 74 /* Update the limited memory preconditioner and get existing # of updates */ 75 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 76 ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 77 } 78 79 /* Use the common BNK kernel to compute the safeguarded Newton step */ 80 ierr = TaoBNKComputeStep(tao, PETSC_TRUE, &stepType);CHKERRQ(ierr); 81 82 /* Store current solution before it changes */ 83 bnk->fold = bnk->f; 84 ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 85 ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr); 86 ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr); 87 88 /* Perform the linesearch */ 89 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 90 ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 91 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 92 93 while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != BNK_GRADIENT) { 94 /* Linesearch failed, revert solution */ 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 100 switch(stepType) { 101 case BNK_NEWTON: 102 /* Failed to obtain acceptable iterate with Newton 1step 103 Update the perturbation for next time */ 104 if (bnk->pert <= 0.0) { 105 /* Initialize the perturbation */ 106 bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 107 if (bnk->is_gltr) { 108 ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 109 bnk->pert = PetscMax(bnk->pert, -e_min); 110 } 111 } else { 112 /* Increase the perturbation */ 113 bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 114 } 115 116 if (BNK_PC_BFGS != bnk->pc_type) { 117 /* We don't have the bfgs matrix around and being updated 118 Must use gradient direction in this case */ 119 ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 120 ++bnk->grad; 121 stepType = BNK_GRADIENT; 122 } else { 123 /* Attempt to use the BFGS direction */ 124 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 125 ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 126 /* Check for success (descent direction) */ 127 ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 128 if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 129 /* BFGS direction is not descent or direction produced not a number 130 We can assert bfgsUpdates > 1 in this case 131 Use steepest descent direction (scaled) */ 132 133 if (bnk->f != 0.0) { 134 delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 135 } else { 136 delta = 2.0 / (bnk->gnorm*bnk->gnorm); 137 } 138 ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 139 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 140 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 141 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 142 ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 143 144 bfgsUpdates = 1; 145 ++bnk->sgrad; 146 stepType = BNK_SCALED_GRADIENT; 147 } else { 148 ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 149 if (1 == bfgsUpdates) { 150 /* The first BFGS direction is always the scaled gradient */ 151 ++bnk->sgrad; 152 stepType = BNK_SCALED_GRADIENT; 153 } else { 154 ++bnk->bfgs; 155 stepType = BNK_BFGS; 156 } 157 } 158 } 159 break; 160 161 case BNK_BFGS: 162 /* Can only enter if pc_type == BNK_PC_BFGS 163 Failed to obtain acceptable iterate with BFGS step 164 Attempt to use the scaled gradient direction */ 165 166 if (bnk->f != 0.0) { 167 delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 168 } else { 169 delta = 2.0 / (bnk->gnorm*bnk->gnorm); 170 } 171 ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 172 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 173 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 174 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 175 ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 176 177 bfgsUpdates = 1; 178 ++bnk->sgrad; 179 stepType = BNK_SCALED_GRADIENT; 180 break; 181 182 case BNK_SCALED_GRADIENT: 183 /* Can only enter if pc_type == BNK_PC_BFGS 184 The scaled gradient step did not produce a new iterate; 185 attemp to use the gradient direction. 186 Need to make sure we are not using a different diagonal scaling */ 187 188 ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr); 189 ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr); 190 ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 191 ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 192 ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 193 ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 194 195 bfgsUpdates = 1; 196 ++bnk->grad; 197 stepType = BNK_GRADIENT; 198 break; 199 } 200 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 201 202 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 203 ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 204 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 205 } 206 207 if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 208 /* Failed to find an improving point */ 209 bnk->f = bnk->fold; 210 ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 211 ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 212 ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 213 step = 0.0; 214 tao->reason = TAO_DIVERGED_LS_FAILURE; 215 break; 216 } 217 218 /* update trust radius */ 219 ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr); 220 ierr = TaoBNKUpdateTrustRadius(tao, bnk->fold, f_full, stepType, &trustAccept);CHKERRQ(ierr); 221 222 /* Check for termination */ 223 ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 224 if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 225 needH = 1; 226 ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 227 ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr); 228 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 229 } 230 PetscFunctionReturn(0); 231 } 232 233 PETSC_EXTERN PetscErrorCode TaoCreate_BNLS(Tao tao) 234 { 235 TAO_BNK *bnk; 236 PetscErrorCode ierr; 237 238 PetscFunctionBegin; 239 ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 240 tao->ops->solve = TaoSolve_BNLS; 241 242 bnk = (TAO_BNK *)tao->data; 243 bnk->update_type = BNK_UPDATE_STEP; /* trust region updates based on line search step length */ 244 PetscFunctionReturn(0); 245 }