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 6*2b97c8d8SAlp Dener bound constrained minimization problems. A projected More'-Thuente line 7*2b97c8d8SAlp 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 24080d2917SAlp Dener PetscReal f_full, gdx; 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 PetscInt needH = 1; 33eb910715SAlp Dener 34eb910715SAlp Dener PetscFunctionBegin; 3509164190SAlp Dener /* Project the current point onto the feasible set */ 3609164190SAlp Dener ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 3709164190SAlp Dener ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 3809164190SAlp Dener 3909164190SAlp Dener /* Project the initial point onto the feasible region */ 4009164190SAlp Dener ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 41eb910715SAlp Dener 42eb910715SAlp Dener /* Check convergence criteria */ 4309164190SAlp Dener ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr); 4409164190SAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 45eb910715SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 46eb910715SAlp Dener if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 47eb910715SAlp Dener 48eb910715SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 49eb910715SAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 50eb910715SAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr); 51eb910715SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 52eb910715SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 53eb910715SAlp Dener 54eb910715SAlp Dener /* Initialize the preconditioner and trust radius */ 55eb910715SAlp Dener ierr = TaoBNKInitialize(tao);CHKERRQ(ierr); 56eb910715SAlp Dener 57eb910715SAlp Dener /* Have not converged; continue with Newton method */ 58eb910715SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 59eb910715SAlp Dener ++tao->niter; 60eb910715SAlp Dener tao->ksp_its=0; 61eb910715SAlp Dener 6209164190SAlp Dener /* Compute the Hessian */ 6309164190SAlp Dener ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 6409164190SAlp Dener 65eb910715SAlp Dener /* Use the common BNK kernel to compute the step */ 66eb910715SAlp Dener ierr = TaoBNKComputeStep(tao, &stepType);CHKERRQ(ierr); 67eb910715SAlp Dener 68080d2917SAlp Dener /* Store current solution before it changes */ 69080d2917SAlp Dener bnk->fold = bnk->f; 70eb910715SAlp Dener ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 71eb910715SAlp Dener ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr); 7209164190SAlp Dener ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr); 73eb910715SAlp Dener 74080d2917SAlp Dener /* Perform the linesearch */ 7509164190SAlp Dener ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 7609164190SAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 77eb910715SAlp Dener ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 78eb910715SAlp Dener 79eb910715SAlp Dener while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != BNK_GRADIENT) { 80eb910715SAlp Dener /* Linesearch failed, revert solution */ 81080d2917SAlp Dener bnk->f = bnk->fold; 82eb910715SAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 83eb910715SAlp Dener ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 8409164190SAlp Dener ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 85eb910715SAlp Dener 86eb910715SAlp Dener switch(stepType) { 87eb910715SAlp Dener case BNK_NEWTON: 88eb910715SAlp Dener /* Failed to obtain acceptable iterate with Newton 1step 89eb910715SAlp Dener Update the perturbation for next time */ 90eb910715SAlp Dener if (bnk->pert <= 0.0) { 91eb910715SAlp Dener /* Initialize the perturbation */ 92eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 93eb910715SAlp Dener if (bnk->is_gltr) { 94eb910715SAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 95eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 96eb910715SAlp Dener } 97eb910715SAlp Dener } else { 98eb910715SAlp Dener /* Increase the perturbation */ 99eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 100eb910715SAlp Dener } 101eb910715SAlp Dener 102eb910715SAlp Dener if (BNK_PC_BFGS != bnk->pc_type) { 103eb910715SAlp Dener /* We don't have the bfgs matrix around and being updated 104eb910715SAlp Dener Must use gradient direction in this case */ 105080d2917SAlp Dener ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 106eb910715SAlp Dener ++bnk->grad; 107eb910715SAlp Dener stepType = BNK_GRADIENT; 108eb910715SAlp Dener } else { 109eb910715SAlp Dener /* Attempt to use the BFGS direction */ 11009164190SAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 11109164190SAlp Dener ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 112eb910715SAlp Dener /* Check for success (descent direction) */ 11309164190SAlp Dener ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 114eb910715SAlp Dener if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 115eb910715SAlp Dener /* BFGS direction is not descent or direction produced not a number 116eb910715SAlp Dener We can assert bfgsUpdates > 1 in this case 117eb910715SAlp Dener Use steepest descent direction (scaled) */ 118eb910715SAlp Dener 119eb910715SAlp Dener if (bnk->f != 0.0) { 120eb910715SAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 121eb910715SAlp Dener } else { 122eb910715SAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 123eb910715SAlp Dener } 124eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 125eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 12609164190SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 12709164190SAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 12809164190SAlp Dener ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 129eb910715SAlp Dener 130eb910715SAlp Dener bfgsUpdates = 1; 131eb910715SAlp Dener ++bnk->sgrad; 132eb910715SAlp Dener stepType = BNK_SCALED_GRADIENT; 133eb910715SAlp Dener } else { 134eb910715SAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 135eb910715SAlp Dener if (1 == bfgsUpdates) { 136eb910715SAlp Dener /* The first BFGS direction is always the scaled gradient */ 137eb910715SAlp Dener ++bnk->sgrad; 138eb910715SAlp Dener stepType = BNK_SCALED_GRADIENT; 139eb910715SAlp Dener } else { 140eb910715SAlp Dener ++bnk->bfgs; 141eb910715SAlp Dener stepType = BNK_BFGS; 142eb910715SAlp Dener } 143eb910715SAlp Dener } 144eb910715SAlp Dener } 145eb910715SAlp Dener break; 146eb910715SAlp Dener 147eb910715SAlp Dener case BNK_BFGS: 148eb910715SAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS 149eb910715SAlp Dener Failed to obtain acceptable iterate with BFGS step 150eb910715SAlp Dener Attempt to use the scaled gradient direction */ 151eb910715SAlp Dener 152eb910715SAlp Dener if (bnk->f != 0.0) { 153eb910715SAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 154eb910715SAlp Dener } else { 155eb910715SAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 156eb910715SAlp Dener } 157eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 158eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 15909164190SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 16009164190SAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 16109164190SAlp Dener ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 162eb910715SAlp Dener 163eb910715SAlp Dener bfgsUpdates = 1; 164eb910715SAlp Dener ++bnk->sgrad; 165eb910715SAlp Dener stepType = BNK_SCALED_GRADIENT; 166eb910715SAlp Dener break; 167eb910715SAlp Dener 168eb910715SAlp Dener case BNK_SCALED_GRADIENT: 169eb910715SAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS 170eb910715SAlp Dener The scaled gradient step did not produce a new iterate; 171eb910715SAlp Dener attemp to use the gradient direction. 172eb910715SAlp Dener Need to make sure we are not using a different diagonal scaling */ 173eb910715SAlp Dener 174eb910715SAlp Dener ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr); 175eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr); 176eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 17709164190SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 17809164190SAlp Dener ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr); 17909164190SAlp Dener ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr); 180eb910715SAlp Dener 181eb910715SAlp Dener bfgsUpdates = 1; 182eb910715SAlp Dener ++bnk->grad; 183eb910715SAlp Dener stepType = BNK_GRADIENT; 184eb910715SAlp Dener break; 185eb910715SAlp Dener } 186080d2917SAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 187eb910715SAlp Dener 18809164190SAlp Dener ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 18909164190SAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 190eb910715SAlp Dener ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 191eb910715SAlp Dener } 192eb910715SAlp Dener 193eb910715SAlp Dener if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 194eb910715SAlp Dener /* Failed to find an improving point */ 195080d2917SAlp Dener bnk->f = bnk->fold; 196eb910715SAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 197eb910715SAlp Dener ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 19809164190SAlp Dener ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 199eb910715SAlp Dener step = 0.0; 200eb910715SAlp Dener tao->reason = TAO_DIVERGED_LS_FAILURE; 201eb910715SAlp Dener break; 202eb910715SAlp Dener } 203eb910715SAlp Dener 204080d2917SAlp Dener /* update trust radius */ 205080d2917SAlp Dener ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr); 206080d2917SAlp Dener ierr = TaoBNKUpdateTrustRadius(tao, bnk->fold, f_full, stepType, &trustAccept);CHKERRQ(ierr); 207eb910715SAlp Dener 208eb910715SAlp Dener /* Check for termination */ 209eb910715SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 210eb910715SAlp Dener if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 211eb910715SAlp Dener needH = 1; 212eb910715SAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 213eb910715SAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr); 214eb910715SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 215eb910715SAlp Dener } 216eb910715SAlp Dener PetscFunctionReturn(0); 217eb910715SAlp Dener } 218eb910715SAlp Dener 219eb910715SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNLS(Tao tao) 220eb910715SAlp Dener { 221eb910715SAlp Dener PetscErrorCode ierr; 222eb910715SAlp Dener 223eb910715SAlp Dener PetscFunctionBegin; 224eb910715SAlp Dener ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 225eb910715SAlp Dener tao->ops->solve=TaoSolve_BNLS; 226eb910715SAlp Dener PetscFunctionReturn(0); 227eb910715SAlp Dener }