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 6eb910715SAlp Dener unconstrained minimization problems. A More'-Thuente line search 7eb910715SAlp Dener 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*080d2917SAlp Dener PetscReal f_full, gdx; 25eb910715SAlp Dener PetscReal step = 1.0; 26eb910715SAlp Dener PetscReal delta; 27*080d2917SAlp Dener PetscReal e_min; 28eb910715SAlp Dener 29*080d2917SAlp Dener PetscBool trustAccept; 30eb910715SAlp Dener PetscInt stepType; 31eb910715SAlp Dener PetscInt bfgsUpdates = 0; 32eb910715SAlp Dener PetscInt needH = 1; 33eb910715SAlp Dener 34eb910715SAlp Dener PetscFunctionBegin; 35eb910715SAlp Dener if (tao->XL || tao->XU || tao->ops->computebounds) { 36eb910715SAlp Dener ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by nls algorithm\n");CHKERRQ(ierr); 37eb910715SAlp Dener } 38eb910715SAlp Dener 39eb910715SAlp Dener /* Check convergence criteria */ 40eb910715SAlp Dener ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, tao->gradient);CHKERRQ(ierr); 41eb910715SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 42eb910715SAlp Dener if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 43eb910715SAlp Dener 44eb910715SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 45eb910715SAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 46eb910715SAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr); 47eb910715SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 48eb910715SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 49eb910715SAlp Dener 50eb910715SAlp Dener /* Initialize the preconditioner and trust radius */ 51eb910715SAlp Dener ierr = TaoBNKInitialize(tao);CHKERRQ(ierr); 52eb910715SAlp Dener 53eb910715SAlp Dener /* Have not converged; continue with Newton method */ 54eb910715SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 55eb910715SAlp Dener ++tao->niter; 56eb910715SAlp Dener tao->ksp_its=0; 57eb910715SAlp Dener 58eb910715SAlp Dener /* Use the common BNK kernel to compute the step */ 59eb910715SAlp Dener ierr = TaoBNKComputeStep(tao, &stepType);CHKERRQ(ierr); 60eb910715SAlp Dener 61*080d2917SAlp Dener /* Store current solution before it changes */ 62*080d2917SAlp Dener bnk->fold = bnk->f; 63eb910715SAlp Dener ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 64eb910715SAlp Dener ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr); 65eb910715SAlp Dener 66*080d2917SAlp Dener /* Perform the linesearch */ 67*080d2917SAlp Dener ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, tao->gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 68eb910715SAlp Dener ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 69eb910715SAlp Dener 70eb910715SAlp Dener while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != BNK_GRADIENT) { 71eb910715SAlp Dener /* Linesearch failed, revert solution */ 72*080d2917SAlp Dener bnk->f = bnk->fold; 73eb910715SAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 74eb910715SAlp Dener ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 75eb910715SAlp Dener 76eb910715SAlp Dener switch(stepType) { 77eb910715SAlp Dener case BNK_NEWTON: 78eb910715SAlp Dener /* Failed to obtain acceptable iterate with Newton 1step 79eb910715SAlp Dener Update the perturbation for next time */ 80eb910715SAlp Dener if (bnk->pert <= 0.0) { 81eb910715SAlp Dener /* Initialize the perturbation */ 82eb910715SAlp Dener bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm)); 83eb910715SAlp Dener if (bnk->is_gltr) { 84eb910715SAlp Dener ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr); 85eb910715SAlp Dener bnk->pert = PetscMax(bnk->pert, -e_min); 86eb910715SAlp Dener } 87eb910715SAlp Dener } else { 88eb910715SAlp Dener /* Increase the perturbation */ 89eb910715SAlp Dener bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm)); 90eb910715SAlp Dener } 91eb910715SAlp Dener 92eb910715SAlp Dener if (BNK_PC_BFGS != bnk->pc_type) { 93eb910715SAlp Dener /* We don't have the bfgs matrix around and being updated 94eb910715SAlp Dener Must use gradient direction in this case */ 95*080d2917SAlp Dener ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 96eb910715SAlp Dener ++bnk->grad; 97eb910715SAlp Dener stepType = BNK_GRADIENT; 98eb910715SAlp Dener } else { 99eb910715SAlp Dener /* Attempt to use the BFGS direction */ 100*080d2917SAlp Dener ierr = MatLMVMSolve(bnk->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 101eb910715SAlp Dener /* Check for success (descent direction) */ 102*080d2917SAlp Dener ierr = VecDot(tao->solution, tao->stepdirection, &gdx);CHKERRQ(ierr); 103eb910715SAlp Dener if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) { 104eb910715SAlp Dener /* BFGS direction is not descent or direction produced not a number 105eb910715SAlp Dener We can assert bfgsUpdates > 1 in this case 106eb910715SAlp Dener Use steepest descent direction (scaled) */ 107eb910715SAlp Dener 108eb910715SAlp Dener if (bnk->f != 0.0) { 109eb910715SAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 110eb910715SAlp Dener } else { 111eb910715SAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 112eb910715SAlp Dener } 113eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 114eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 115eb910715SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, tao->gradient);CHKERRQ(ierr); 116*080d2917SAlp Dener ierr = MatLMVMSolve(bnk->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 117eb910715SAlp Dener 118eb910715SAlp Dener bfgsUpdates = 1; 119eb910715SAlp Dener ++bnk->sgrad; 120eb910715SAlp Dener stepType = BNK_SCALED_GRADIENT; 121eb910715SAlp Dener } else { 122eb910715SAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 123eb910715SAlp Dener if (1 == bfgsUpdates) { 124eb910715SAlp Dener /* The first BFGS direction is always the scaled gradient */ 125eb910715SAlp Dener ++bnk->sgrad; 126eb910715SAlp Dener stepType = BNK_SCALED_GRADIENT; 127eb910715SAlp Dener } else { 128eb910715SAlp Dener ++bnk->bfgs; 129eb910715SAlp Dener stepType = BNK_BFGS; 130eb910715SAlp Dener } 131eb910715SAlp Dener } 132eb910715SAlp Dener } 133eb910715SAlp Dener break; 134eb910715SAlp Dener 135eb910715SAlp Dener case BNK_BFGS: 136eb910715SAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS 137eb910715SAlp Dener Failed to obtain acceptable iterate with BFGS step 138eb910715SAlp Dener Attempt to use the scaled gradient direction */ 139eb910715SAlp Dener 140eb910715SAlp Dener if (bnk->f != 0.0) { 141eb910715SAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm); 142eb910715SAlp Dener } else { 143eb910715SAlp Dener delta = 2.0 / (bnk->gnorm*bnk->gnorm); 144eb910715SAlp Dener } 145eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr); 146eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 147eb910715SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, tao->gradient);CHKERRQ(ierr); 148*080d2917SAlp Dener ierr = MatLMVMSolve(bnk->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 149eb910715SAlp Dener 150eb910715SAlp Dener bfgsUpdates = 1; 151eb910715SAlp Dener ++bnk->sgrad; 152eb910715SAlp Dener stepType = BNK_SCALED_GRADIENT; 153eb910715SAlp Dener break; 154eb910715SAlp Dener 155eb910715SAlp Dener case BNK_SCALED_GRADIENT: 156eb910715SAlp Dener /* Can only enter if pc_type == BNK_PC_BFGS 157eb910715SAlp Dener The scaled gradient step did not produce a new iterate; 158eb910715SAlp Dener attemp to use the gradient direction. 159eb910715SAlp Dener Need to make sure we are not using a different diagonal scaling */ 160eb910715SAlp Dener 161eb910715SAlp Dener ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr); 162eb910715SAlp Dener ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr); 163eb910715SAlp Dener ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr); 164eb910715SAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, tao->gradient);CHKERRQ(ierr); 165*080d2917SAlp Dener ierr = MatLMVMSolve(bnk->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 166eb910715SAlp Dener 167eb910715SAlp Dener bfgsUpdates = 1; 168eb910715SAlp Dener ++bnk->grad; 169eb910715SAlp Dener stepType = BNK_GRADIENT; 170eb910715SAlp Dener break; 171eb910715SAlp Dener } 172*080d2917SAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 173eb910715SAlp Dener 174*080d2917SAlp Dener ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, tao->gradient, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr); 175eb910715SAlp Dener ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 176eb910715SAlp Dener } 177eb910715SAlp Dener 178eb910715SAlp Dener if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 179eb910715SAlp Dener /* Failed to find an improving point */ 180*080d2917SAlp Dener bnk->f = bnk->fold; 181eb910715SAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 182eb910715SAlp Dener ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 183eb910715SAlp Dener step = 0.0; 184eb910715SAlp Dener tao->reason = TAO_DIVERGED_LS_FAILURE; 185eb910715SAlp Dener break; 186eb910715SAlp Dener } 187eb910715SAlp Dener 188*080d2917SAlp Dener /* update trust radius */ 189*080d2917SAlp Dener ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr); 190*080d2917SAlp Dener ierr = TaoBNKUpdateTrustRadius(tao, bnk->fold, f_full, stepType, &trustAccept);CHKERRQ(ierr); 191eb910715SAlp Dener 192eb910715SAlp Dener /* Check for termination */ 193eb910715SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 194eb910715SAlp Dener if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 195eb910715SAlp Dener needH = 1; 196eb910715SAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 197eb910715SAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr); 198eb910715SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 199eb910715SAlp Dener } 200eb910715SAlp Dener PetscFunctionReturn(0); 201eb910715SAlp Dener } 202eb910715SAlp Dener 203eb910715SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNLS(Tao tao) 204eb910715SAlp Dener { 205eb910715SAlp Dener PetscErrorCode ierr; 206eb910715SAlp Dener 207eb910715SAlp Dener PetscFunctionBegin; 208eb910715SAlp Dener ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 209eb910715SAlp Dener tao->ops->solve=TaoSolve_BNLS; 210eb910715SAlp Dener PetscFunctionReturn(0); 211eb910715SAlp Dener }