16b591159SAlp Dener #include <../src/tao/bound/impls/bqnk/bqnk.h> 26b591159SAlp Dener 36b591159SAlp Dener static PetscErrorCode TaoBQNLSComputeHessian(Tao tao) 46b591159SAlp Dener { 56b591159SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 66b591159SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 76b591159SAlp Dener PetscErrorCode ierr; 8d5ae2380SAlp Dener PetscReal gnorm2, delta; 96b591159SAlp Dener 106b591159SAlp Dener PetscFunctionBegin; 11d5ae2380SAlp Dener gnorm2 = bnk->gnorm*bnk->gnorm; 12*65f5217aSAlp Dener delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / PetscMax(gnorm2, PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0)); 13d5ae2380SAlp Dener ierr = MatSymBrdnSetDelta(bqnk->B, delta);CHKERRQ(ierr); 146b591159SAlp Dener ierr = MatLMVMUpdate(bqnk->B, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 156b591159SAlp Dener PetscFunctionReturn(0); 166b591159SAlp Dener } 176b591159SAlp Dener 186b591159SAlp Dener static PetscErrorCode TaoBQNLSComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type) 196b591159SAlp Dener { 206b591159SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 216b591159SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 226b591159SAlp Dener PetscErrorCode ierr; 23*65f5217aSAlp Dener PetscInt nupdates; 246b591159SAlp Dener 256b591159SAlp Dener PetscFunctionBegin; 266b591159SAlp Dener ierr = MatSolve(bqnk->B, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 276b591159SAlp Dener ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 286b591159SAlp Dener ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr); 296b591159SAlp Dener *ksp_reason = KSP_CONVERGED_ATOL; 30*65f5217aSAlp Dener ierr = MatLMVMGetUpdateCount(bqnk->B, &nupdates);CHKERRQ(ierr); 31*65f5217aSAlp Dener if (nupdates == 0) { 32*65f5217aSAlp Dener *step_type = BNK_SCALED_GRADIENT; 33*65f5217aSAlp Dener } else { 346b591159SAlp Dener *step_type = BNK_BFGS; 35*65f5217aSAlp Dener } 366b591159SAlp Dener PetscFunctionReturn(0); 376b591159SAlp Dener } 386b591159SAlp Dener 396b591159SAlp Dener static PetscErrorCode TaoSetFromOptions_BQNLS(PetscOptionItems *PetscOptionsObject,Tao tao) 406b591159SAlp Dener { 416b591159SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 426b591159SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 436b591159SAlp Dener PetscErrorCode ierr; 446b591159SAlp Dener KSPType ksp_type; 456b591159SAlp Dener PetscBool is_spd; 466b591159SAlp Dener 476b591159SAlp Dener PetscFunctionBegin; 486b591159SAlp Dener ierr = PetscOptionsHead(PetscOptionsObject,"Quasi-Newton-Krylov method for bound constrained optimization");CHKERRQ(ierr); 496b591159SAlp Dener ierr = PetscOptionsEList("-tao_bqnls_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, 0);CHKERRQ(ierr); 506b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnls_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr); 516b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnls_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr); 526b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnls_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr); 536b591159SAlp Dener ierr = PetscOptionsInt("-tao_bqnls_max_cg_its", "number of BNCG iterations to take for each Newton step", "", bnk->max_cg_its, &bnk->max_cg_its,NULL);CHKERRQ(ierr); 546b591159SAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 556b591159SAlp Dener ierr = TaoSetFromOptions(bnk->bncg);CHKERRQ(ierr); 566b591159SAlp Dener ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 576b591159SAlp Dener ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 586b591159SAlp Dener ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr); 596b591159SAlp Dener bnk->is_nash = bnk->is_gltr = bnk->is_stcg = PETSC_FALSE; 606b591159SAlp Dener ierr = MatSetFromOptions(bqnk->B);CHKERRQ(ierr); 616b591159SAlp Dener ierr = MatGetOption(bqnk->B, MAT_SPD, &is_spd);CHKERRQ(ierr); 626b591159SAlp Dener if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite"); 636b591159SAlp Dener PetscFunctionReturn(0); 646b591159SAlp Dener } 656b591159SAlp Dener 666b591159SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BQNLS(Tao tao) 676b591159SAlp Dener { 686b591159SAlp Dener TAO_BNK *bnk; 696b591159SAlp Dener TAO_BQNK *bqnk; 706b591159SAlp Dener PetscErrorCode ierr; 716b591159SAlp Dener 726b591159SAlp Dener PetscFunctionBegin; 736b591159SAlp Dener ierr = TaoCreate_BQNK(tao);CHKERRQ(ierr); 746b591159SAlp Dener ierr = KSPSetOptionsPrefix(tao->ksp, "unused");CHKERRQ(ierr); 756b591159SAlp Dener tao->ops->solve = TaoSolve_BNLS; 766b591159SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BQNLS; 776b591159SAlp Dener 786b591159SAlp Dener bnk = (TAO_BNK*)tao->data; 796b591159SAlp Dener bnk->update_type = BNK_UPDATE_STEP; 806b591159SAlp Dener bnk->computehessian = TaoBQNLSComputeHessian; 816b591159SAlp Dener bnk->computestep = TaoBQNLSComputeStep; 826b591159SAlp Dener 836b591159SAlp Dener bqnk = (TAO_BQNK*)bnk->ctx; 846b591159SAlp Dener ierr = MatSetOptionsPrefix(bqnk->B, "tao_bqnls_");CHKERRQ(ierr); 856b591159SAlp Dener ierr = MatSetType(bqnk->B, MATLMVMBFGS);CHKERRQ(ierr); 866b591159SAlp Dener PetscFunctionReturn(0); 876b591159SAlp Dener }