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 24c14b763aSAlp Dener PetscReal f_full, prered, actred; 25c14b763aSAlp Dener PetscReal steplen = 1.0; 26eb910715SAlp Dener 27080d2917SAlp Dener PetscBool trustAccept; 28eb910715SAlp Dener PetscInt stepType; 29eb910715SAlp Dener PetscInt bfgsUpdates = 0; 30eb910715SAlp Dener 31eb910715SAlp Dener PetscFunctionBegin; 3209164190SAlp Dener /* Project the current point onto the feasible set */ 3309164190SAlp Dener ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 3409164190SAlp Dener ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 3509164190SAlp Dener 3609164190SAlp Dener /* Project the initial point onto the feasible region */ 3709164190SAlp Dener ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 38eb910715SAlp Dener 39eb910715SAlp Dener /* Check convergence criteria */ 4009164190SAlp Dener ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr); 4109164190SAlp Dener ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 42eb910715SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 43eb910715SAlp Dener if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 44eb910715SAlp Dener 45eb910715SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 46eb910715SAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 47c14b763aSAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,steplen);CHKERRQ(ierr); 48eb910715SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 49eb910715SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 50eb910715SAlp Dener 51eb910715SAlp Dener /* Initialize the preconditioner and trust radius */ 52eb910715SAlp Dener ierr = TaoBNKInitialize(tao);CHKERRQ(ierr); 53eb910715SAlp Dener 54eb910715SAlp Dener /* Have not converged; continue with Newton method */ 55eb910715SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 56eb910715SAlp Dener ++tao->niter; 57eb910715SAlp Dener tao->ksp_its=0; 58eb910715SAlp Dener 5909164190SAlp Dener /* Compute the Hessian */ 6009164190SAlp Dener ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 6109164190SAlp Dener 62fed79b8eSAlp Dener /* Update the BFGS preconditioner */ 63fed79b8eSAlp Dener if (BNK_PC_BFGS == bnk->pc_type) { 64fed79b8eSAlp Dener if (BFGS_SCALE_PHESS == bnk->bfgs_scale_type) { 65fed79b8eSAlp Dener /* Obtain diagonal for the bfgs preconditioner */ 66fed79b8eSAlp Dener ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr); 67fed79b8eSAlp Dener ierr = VecAbs(bnk->Diag);CHKERRQ(ierr); 68fed79b8eSAlp Dener ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr); 69fed79b8eSAlp Dener ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr); 70fed79b8eSAlp Dener } 71fed79b8eSAlp Dener /* Update the limited memory preconditioner and get existing # of updates */ 72fed79b8eSAlp Dener ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 73fed79b8eSAlp Dener ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr); 74fed79b8eSAlp Dener } 75fed79b8eSAlp Dener 76fed79b8eSAlp Dener /* Use the common BNK kernel to compute the safeguarded Newton step */ 77fed79b8eSAlp Dener ierr = TaoBNKComputeStep(tao, PETSC_TRUE, &stepType);CHKERRQ(ierr); 78eb910715SAlp Dener 79080d2917SAlp Dener /* Store current solution before it changes */ 80080d2917SAlp Dener bnk->fold = bnk->f; 81eb910715SAlp Dener ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr); 82eb910715SAlp Dener ierr = VecCopy(tao->gradient, bnk->Gold);CHKERRQ(ierr); 8309164190SAlp Dener ierr = VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);CHKERRQ(ierr); 84eb910715SAlp Dener 85c14b763aSAlp Dener /* Trigger the line search */ 86c14b763aSAlp Dener ierr = TaoBNKPerformLineSearch(tao, stepType, &steplen, &ls_reason);CHKERRQ(ierr); 87eb910715SAlp Dener 88eb910715SAlp Dener if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) { 89eb910715SAlp Dener /* Failed to find an improving point */ 90080d2917SAlp Dener bnk->f = bnk->fold; 91eb910715SAlp Dener ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr); 92eb910715SAlp Dener ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr); 9309164190SAlp Dener ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr); 94c14b763aSAlp Dener steplen = 0.0; 95eb910715SAlp Dener tao->reason = TAO_DIVERGED_LS_FAILURE; 96eb910715SAlp Dener break; 97eb910715SAlp Dener } 98eb910715SAlp Dener 99080d2917SAlp Dener /* update trust radius */ 100080d2917SAlp Dener ierr = TaoLineSearchGetFullStepObjective(tao->linesearch, &f_full);CHKERRQ(ierr); 101b1c2d0e3SAlp Dener ierr = KSPCGGetObjFcn(tao->ksp, &prered);CHKERRQ(ierr); 102b1c2d0e3SAlp Dener prered = -prered; 103b1c2d0e3SAlp Dener actred = bnk->fold - f_full; 104b1c2d0e3SAlp Dener ierr = TaoBNKUpdateTrustRadius(tao, prered, actred, stepType, &trustAccept);CHKERRQ(ierr); 105eb910715SAlp Dener 106eb910715SAlp Dener /* Check for termination */ 107eb910715SAlp Dener ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr); 108eb910715SAlp Dener if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User provided compute function generated Not-a-Number"); 109eb910715SAlp Dener ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 110c14b763aSAlp Dener ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,steplen);CHKERRQ(ierr); 111eb910715SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 112eb910715SAlp Dener } 113eb910715SAlp Dener PetscFunctionReturn(0); 114eb910715SAlp Dener } 115eb910715SAlp Dener 116*df278d8fSAlp Dener /*------------------------------------------------------------*/ 117*df278d8fSAlp Dener 118eb910715SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BNLS(Tao tao) 119eb910715SAlp Dener { 120fed79b8eSAlp Dener TAO_BNK *bnk; 121eb910715SAlp Dener PetscErrorCode ierr; 122eb910715SAlp Dener 123eb910715SAlp Dener PetscFunctionBegin; 124eb910715SAlp Dener ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 125eb910715SAlp Dener tao->ops->solve = TaoSolve_BNLS; 126fed79b8eSAlp Dener 127fed79b8eSAlp Dener bnk = (TAO_BNK *)tao->data; 12866ed3702SAlp Dener bnk->update_type = BNK_UPDATE_STEP; /* trust region updates based on line search step length */ 129eb910715SAlp Dener PetscFunctionReturn(0); 130eb910715SAlp Dener }