1fed79b8eSAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h> 2fed79b8eSAlp Dener #include <petscksp.h> 3fed79b8eSAlp Dener 4fed79b8eSAlp Dener /* 5fed79b8eSAlp Dener Implements Newton's Method with a trust region approach for solving 6fed79b8eSAlp Dener bound constrained minimization problems. 7fed79b8eSAlp Dener 8198282dbSAlp Dener ------------------------------------------------------------ 9198282dbSAlp Dener 10198282dbSAlp Dener initialize trust radius (default: BNK_INIT_INTERPOLATION) 11198282dbSAlp Dener x_0 = VecMedian(x_0) 12198282dbSAlp Dener f_0, g_0 = TaoComputeObjectiveAndGradient(x_0) 13198282dbSAlp Dener pg_0 = VecBoundGradientProjection(g_0) 14198282dbSAlp Dener check convergence at pg_0 15198282dbSAlp Dener niter = 0 16198282dbSAlp Dener step_accepted = true 17198282dbSAlp Dener 18198282dbSAlp Dener while niter <= max_it 19198282dbSAlp Dener if step_accepted 20198282dbSAlp Dener niter += 1 21198282dbSAlp Dener H_k = TaoComputeHessian(x_k) 22198282dbSAlp Dener if pc_type == BNK_PC_BFGS 23198282dbSAlp Dener add correction to BFGS approx 24198282dbSAlp Dener if scale_type == BNK_SCALE_AHESS 25198282dbSAlp Dener D = VecMedian(1e-6, abs(diag(H_k)), 1e6) 26198282dbSAlp Dener scale BFGS with VecReciprocal(D) 27198282dbSAlp Dener end 28198282dbSAlp Dener end 29198282dbSAlp Dener end 30198282dbSAlp Dener 31198282dbSAlp Dener if pc_type = BNK_PC_BFGS 32198282dbSAlp Dener B_k = BFGS 33198282dbSAlp Dener else 34198282dbSAlp Dener B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6) 35198282dbSAlp Dener B_k = VecReciprocal(B_k) 36198282dbSAlp Dener end 37198282dbSAlp Dener w = x_k - VecMedian(x_k - 0.001*B_k*g_k) 38198282dbSAlp Dener eps = min(eps, norm2(w)) 39198282dbSAlp Dener determine the active and inactive index sets such that 40198282dbSAlp Dener L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0} 41198282dbSAlp Dener U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0} 42198282dbSAlp Dener F = {i : l_i = (x_k)_i = u_i} 43198282dbSAlp Dener A = {L + U + F} 44198282dbSAlp Dener I = {i : i not in A} 45198282dbSAlp Dener 46198282dbSAlp Dener generate the reduced system Hr_k dr_k = -gr_k for variables in I 47198282dbSAlp Dener if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS 48198282dbSAlp Dener D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6) 49198282dbSAlp Dener scale BFGS with VecReciprocal(D) 50198282dbSAlp Dener end 51198282dbSAlp Dener solve Hr_k dr_k = -gr_k 52198282dbSAlp Dener set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F 53198282dbSAlp Dener 54198282dbSAlp Dener x_{k+1} = VecMedian(x_k + d_k) 55198282dbSAlp Dener s = x_{k+1} - x_k 56198282dbSAlp Dener prered = dot(s, 0.5*gr_k - Hr_k*s) 57198282dbSAlp Dener f_{k+1} = TaoComputeObjective(x_{k+1}) 58198282dbSAlp Dener actred = f_k - f_{k+1} 59198282dbSAlp Dener 60198282dbSAlp Dener oldTrust = trust 61198282dbSAlp Dener step_accepted, trust = TaoBNKUpdateTrustRadius(default: BNK_UPDATE_REDUCTION) 62198282dbSAlp Dener if step_accepted 63198282dbSAlp Dener g_{k+1} = TaoComputeGradient(x_{k+1}) 64198282dbSAlp Dener pg_{k+1} = VecBoundGradientProjection(g_{k+1}) 65198282dbSAlp Dener count the accepted Newton step 66198282dbSAlp Dener else 67198282dbSAlp Dener f_{k+1} = f_k 68198282dbSAlp Dener x_{k+1} = x_k 69198282dbSAlp Dener g_{k+1} = g_k 70198282dbSAlp Dener pg_{k+1} = pg_k 71198282dbSAlp Dener if trust == oldTrust 72198282dbSAlp Dener terminate because we cannot shrink the radius any further 73198282dbSAlp Dener end 74198282dbSAlp Dener end 75198282dbSAlp Dener 76198282dbSAlp Dener check convergence at pg_{k+1} 77198282dbSAlp Dener end 78fed79b8eSAlp Dener */ 79fed79b8eSAlp Dener 80fed79b8eSAlp Dener static PetscErrorCode TaoSolve_BNTR(Tao tao) 81fed79b8eSAlp Dener { 82fed79b8eSAlp Dener PetscErrorCode ierr; 83fed79b8eSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 84e465cd6fSAlp Dener KSPConvergedReason ksp_reason; 85fed79b8eSAlp Dener 86*9b6ef848SAlp Dener PetscReal resnorm, oldTrust, prered, actred, stepNorm, steplen; 8762675beeSAlp Dener PetscBool stepAccepted = PETSC_TRUE, shift = PETSC_FALSE; 88e465cd6fSAlp Dener PetscInt stepType = BNK_NEWTON; 89fed79b8eSAlp Dener 90fed79b8eSAlp Dener PetscFunctionBegin; 9128017e9fSAlp Dener /* Initialize the preconditioner, KSP solver and trust radius/line search */ 92fed79b8eSAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 9362675beeSAlp Dener ierr = TaoBNKInitialize(tao, bnk->init_type);CHKERRQ(ierr); 9428017e9fSAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 95fed79b8eSAlp Dener 96fed79b8eSAlp Dener /* Have not converged; continue with Newton method */ 97fed79b8eSAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 9866ed3702SAlp Dener 99fed79b8eSAlp Dener if (stepAccepted) { 100fed79b8eSAlp Dener tao->niter++; 101fed79b8eSAlp Dener tao->ksp_its=0; 10262675beeSAlp Dener /* Compute the hessian and update the BFGS preconditioner at the new iterate*/ 10362675beeSAlp Dener ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr); 104fed79b8eSAlp Dener } 105fed79b8eSAlp Dener 1068d5ead36SAlp Dener /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */ 10762675beeSAlp Dener ierr = TaoBNKComputeStep(tao, shift, &ksp_reason);CHKERRQ(ierr); 108fed79b8eSAlp Dener 109fed79b8eSAlp Dener /* Store current solution before it changes */ 110fed79b8eSAlp Dener oldTrust = tao->trust; 111fed79b8eSAlp Dener bnk->fold = bnk->f; 112fed79b8eSAlp Dener ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 113fed79b8eSAlp Dener ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr); 114fed79b8eSAlp Dener ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr); 115fed79b8eSAlp Dener 116b1c2d0e3SAlp Dener /* Temporarily accept the step and project it into the bounds */ 117fed79b8eSAlp Dener ierr = VecAXPY(tao->solution, 1.0, tao->stepdirection);CHKERRQ(ierr); 118b1c2d0e3SAlp Dener ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr); 119b1c2d0e3SAlp Dener 120b1c2d0e3SAlp Dener /* Check if the projection changed the step direction */ 121b1c2d0e3SAlp Dener ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr); 1228d5ead36SAlp Dener ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr); 123b1c2d0e3SAlp Dener ierr = VecNorm(tao->stepdirection, NORM_2, &stepNorm);CHKERRQ(ierr); 124b1c2d0e3SAlp Dener if (stepNorm != bnk->dnorm) { 1258d5ead36SAlp Dener /* Projection changed the step, so we have to recompute predicted reduction. 1268d5ead36SAlp Dener However, we deliberately do not change the step norm and the trust radius 1278d5ead36SAlp Dener in order for the safeguard to more closely mimic a piece-wise linesearch 1288d5ead36SAlp Dener along the bounds. */ 12928017e9fSAlp Dener ierr = MatMult(bnk->H_inactive, tao->stepdirection, bnk->Xwork);CHKERRQ(ierr); 130198282dbSAlp Dener ierr = VecAYPX(bnk->Xwork, -0.5, bnk->G_inactive);CHKERRQ(ierr); 131b1c2d0e3SAlp Dener ierr = VecDot(bnk->Xwork, tao->stepdirection, &prered); 132b1c2d0e3SAlp Dener } else { 133b1c2d0e3SAlp Dener /* Step did not change, so we can just recover the pre-computed prediction */ 134b1c2d0e3SAlp Dener ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr); 135b1c2d0e3SAlp Dener } 136b1c2d0e3SAlp Dener prered = -prered; 137b1c2d0e3SAlp Dener 138b1c2d0e3SAlp Dener /* Compute the actual reduction and update the trust radius */ 139fed79b8eSAlp Dener ierr = TaoComputeObjective(tao, tao->solution, &bnk->f);CHKERRQ(ierr); 140b1c2d0e3SAlp Dener actred = bnk->fold - bnk->f; 14128017e9fSAlp Dener ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);CHKERRQ(ierr); 142fed79b8eSAlp Dener 143fed79b8eSAlp Dener if (stepAccepted) { 14466ed3702SAlp Dener /* Step is good, evaluate the gradient and the hessian */ 1458d5ead36SAlp Dener steplen = 1.0; 146e465cd6fSAlp Dener ++bnk->newt; 147fed79b8eSAlp Dener ierr = TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 148fed79b8eSAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 149*9b6ef848SAlp Dener ierr = VecNorm(tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr); 150*9b6ef848SAlp Dener if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 151fed79b8eSAlp Dener } else { 152fed79b8eSAlp Dener /* Step is bad, revert old solution and re-solve with new radius*/ 1538d5ead36SAlp Dener steplen = 0.0; 154fed79b8eSAlp Dener bnk->f = bnk->fold; 155fed79b8eSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 156fed79b8eSAlp Dener ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 157fed79b8eSAlp Dener ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 15873e4db90SAlp Dener if (oldTrust == tao->trust) { 15973e4db90SAlp Dener /* Can't change the radius anymore so just terminate */ 160fed79b8eSAlp Dener tao->reason = TAO_DIVERGED_TR_REDUCTION; 161fed79b8eSAlp Dener } 162fed79b8eSAlp Dener } 163fed79b8eSAlp Dener 164fed79b8eSAlp Dener /* Check for termination */ 165*9b6ef848SAlp Dener ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->Gwork);CHKERRQ(ierr); 166*9b6ef848SAlp Dener ierr = VecNorm(bnk->Gwork, NORM_2, &resnorm);CHKERRQ(ierr); 167*9b6ef848SAlp Dener ierr = TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 168*9b6ef848SAlp Dener ierr = TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);CHKERRQ(ierr); 169fed79b8eSAlp Dener ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr); 170fed79b8eSAlp Dener } 171fed79b8eSAlp Dener PetscFunctionReturn(0); 172fed79b8eSAlp Dener } 173fed79b8eSAlp Dener 174df278d8fSAlp Dener /*------------------------------------------------------------*/ 175df278d8fSAlp Dener 176*9b6ef848SAlp Dener PETSC_INTERN PetscErrorCode TaoSetUp_BNTR(Tao tao) 177*9b6ef848SAlp Dener { 178*9b6ef848SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 179*9b6ef848SAlp Dener PetscErrorCode ierr; 180*9b6ef848SAlp Dener 181*9b6ef848SAlp Dener PetscFunctionBegin; 182*9b6ef848SAlp Dener ierr = TaoSetUp_BNK(tao);CHKERRQ(ierr); 183*9b6ef848SAlp Dener if (!bnk->is_nash && !bnk->is_stcg && !bnk->is_gltr) SETERRQ(PETSC_COMM_SELF,1,"Must use a trust-region CG method for KSP (KSPNASH, KSPSTCG, KSPGLTR)"); 184*9b6ef848SAlp Dener PetscFunctionReturn(0); 185*9b6ef848SAlp Dener } 186*9b6ef848SAlp Dener 187*9b6ef848SAlp Dener /*------------------------------------------------------------*/ 188*9b6ef848SAlp Dener 189*9b6ef848SAlp Dener PETSC_INTERN PetscErrorCode TaoCreate_BNTR(Tao tao) 190fed79b8eSAlp Dener { 191fed79b8eSAlp Dener TAO_BNK *bnk; 192fed79b8eSAlp Dener PetscErrorCode ierr; 193fed79b8eSAlp Dener 194fed79b8eSAlp Dener PetscFunctionBegin; 195fed79b8eSAlp Dener ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 196fed79b8eSAlp Dener tao->ops->solve=TaoSolve_BNTR; 197*9b6ef848SAlp Dener tao->ops->setup=TaoSetUp_BNTR; 198fed79b8eSAlp Dener 199fed79b8eSAlp Dener bnk = (TAO_BNK *)tao->data; 20066ed3702SAlp Dener bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */ 20166ed3702SAlp Dener bnk->sval = 0.0; /* disable Hessian shifting */ 202fed79b8eSAlp Dener PetscFunctionReturn(0); 203fed79b8eSAlp Dener }