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 8*198282dbSAlp Dener ------------------------------------------------------------ 9*198282dbSAlp Dener 10*198282dbSAlp Dener initialize trust radius (default: BNK_INIT_INTERPOLATION) 11*198282dbSAlp Dener x_0 = VecMedian(x_0) 12*198282dbSAlp Dener f_0, g_0 = TaoComputeObjectiveAndGradient(x_0) 13*198282dbSAlp Dener pg_0 = VecBoundGradientProjection(g_0) 14*198282dbSAlp Dener check convergence at pg_0 15*198282dbSAlp Dener niter = 0 16*198282dbSAlp Dener step_accepted = true 17*198282dbSAlp Dener 18*198282dbSAlp Dener while niter <= max_it 19*198282dbSAlp Dener if step_accepted 20*198282dbSAlp Dener niter += 1 21*198282dbSAlp Dener H_k = TaoComputeHessian(x_k) 22*198282dbSAlp Dener if pc_type == BNK_PC_BFGS 23*198282dbSAlp Dener add correction to BFGS approx 24*198282dbSAlp Dener if scale_type == BNK_SCALE_AHESS 25*198282dbSAlp Dener D = VecMedian(1e-6, abs(diag(H_k)), 1e6) 26*198282dbSAlp Dener scale BFGS with VecReciprocal(D) 27*198282dbSAlp Dener end 28*198282dbSAlp Dener end 29*198282dbSAlp Dener end 30*198282dbSAlp Dener 31*198282dbSAlp Dener if pc_type = BNK_PC_BFGS 32*198282dbSAlp Dener B_k = BFGS 33*198282dbSAlp Dener else 34*198282dbSAlp Dener B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6) 35*198282dbSAlp Dener B_k = VecReciprocal(B_k) 36*198282dbSAlp Dener end 37*198282dbSAlp Dener w = x_k - VecMedian(x_k - 0.001*B_k*g_k) 38*198282dbSAlp Dener eps = min(eps, norm2(w)) 39*198282dbSAlp Dener determine the active and inactive index sets such that 40*198282dbSAlp Dener L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0} 41*198282dbSAlp Dener U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0} 42*198282dbSAlp Dener F = {i : l_i = (x_k)_i = u_i} 43*198282dbSAlp Dener A = {L + U + F} 44*198282dbSAlp Dener I = {i : i not in A} 45*198282dbSAlp Dener 46*198282dbSAlp Dener generate the reduced system Hr_k dr_k = -gr_k for variables in I 47*198282dbSAlp Dener if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS 48*198282dbSAlp Dener D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6) 49*198282dbSAlp Dener scale BFGS with VecReciprocal(D) 50*198282dbSAlp Dener end 51*198282dbSAlp Dener solve Hr_k dr_k = -gr_k 52*198282dbSAlp Dener set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F 53*198282dbSAlp Dener 54*198282dbSAlp Dener x_{k+1} = VecMedian(x_k + d_k) 55*198282dbSAlp Dener s = x_{k+1} - x_k 56*198282dbSAlp Dener prered = dot(s, 0.5*gr_k - Hr_k*s) 57*198282dbSAlp Dener f_{k+1} = TaoComputeObjective(x_{k+1}) 58*198282dbSAlp Dener actred = f_k - f_{k+1} 59*198282dbSAlp Dener 60*198282dbSAlp Dener oldTrust = trust 61*198282dbSAlp Dener step_accepted, trust = TaoBNKUpdateTrustRadius(default: BNK_UPDATE_REDUCTION) 62*198282dbSAlp Dener if step_accepted 63*198282dbSAlp Dener g_{k+1} = TaoComputeGradient(x_{k+1}) 64*198282dbSAlp Dener pg_{k+1} = VecBoundGradientProjection(g_{k+1}) 65*198282dbSAlp Dener count the accepted Newton step 66*198282dbSAlp Dener else 67*198282dbSAlp Dener f_{k+1} = f_k 68*198282dbSAlp Dener x_{k+1} = x_k 69*198282dbSAlp Dener g_{k+1} = g_k 70*198282dbSAlp Dener pg_{k+1} = pg_k 71*198282dbSAlp Dener if trust == oldTrust 72*198282dbSAlp Dener terminate because we cannot shrink the radius any further 73*198282dbSAlp Dener end 74*198282dbSAlp Dener end 75*198282dbSAlp Dener 76*198282dbSAlp Dener check convergence at pg_{k+1} 77*198282dbSAlp 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 868d5ead36SAlp Dener PetscReal 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); 130*198282dbSAlp 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); 149fed79b8eSAlp Dener } else { 150fed79b8eSAlp Dener /* Step is bad, revert old solution and re-solve with new radius*/ 1518d5ead36SAlp Dener steplen = 0.0; 152fed79b8eSAlp Dener bnk->f = bnk->fold; 153fed79b8eSAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 154fed79b8eSAlp Dener ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 155fed79b8eSAlp Dener ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 15673e4db90SAlp Dener if (oldTrust == tao->trust) { 15773e4db90SAlp Dener /* Can't change the radius anymore so just terminate */ 158fed79b8eSAlp Dener tao->reason = TAO_DIVERGED_TR_REDUCTION; 159fed79b8eSAlp Dener } 160fed79b8eSAlp Dener } 161fed79b8eSAlp Dener 162fed79b8eSAlp Dener /* Check for termination */ 163fed79b8eSAlp Dener ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr); 16473e4db90SAlp Dener if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 165fed79b8eSAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 1668d5ead36SAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,steplen);CHKERRQ(ierr); 167fed79b8eSAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 168fed79b8eSAlp Dener } 169fed79b8eSAlp Dener PetscFunctionReturn(0); 170fed79b8eSAlp Dener } 171fed79b8eSAlp Dener 172df278d8fSAlp Dener /*------------------------------------------------------------*/ 173df278d8fSAlp Dener 174fed79b8eSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao) 175fed79b8eSAlp Dener { 176fed79b8eSAlp Dener TAO_BNK *bnk; 177fed79b8eSAlp Dener PetscErrorCode ierr; 178fed79b8eSAlp Dener 179fed79b8eSAlp Dener PetscFunctionBegin; 180fed79b8eSAlp Dener ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 181fed79b8eSAlp Dener tao->ops->solve=TaoSolve_BNTR; 182fed79b8eSAlp Dener 183fed79b8eSAlp Dener bnk = (TAO_BNK *)tao->data; 18466ed3702SAlp Dener bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */ 18566ed3702SAlp Dener bnk->sval = 0.0; /* disable Hessian shifting */ 186fed79b8eSAlp Dener PetscFunctionReturn(0); 187fed79b8eSAlp Dener }