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