1 #include <../src/tao/bound/impls/bqnk/bqnk.h> 2 3 static PetscErrorCode TaoBQNLSComputeHessian(Tao tao) 4 { 5 TAO_BNK *bnk = (TAO_BNK *)tao->data; 6 TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 7 PetscErrorCode ierr; 8 PetscReal gnorm2, delta; 9 10 PetscFunctionBegin; 11 gnorm2 = bnk->gnorm*bnk->gnorm; 12 if (gnorm2 == 0.0) gnorm2 = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); 13 if (bnk->f != 0.0) { 14 delta = 2.0*PetscAbsScalar(bnk->f) / gnorm2; 15 } else { 16 delta = 2.0 / gnorm2; 17 } 18 ierr = MatSymBrdnSetDelta(bqnk->B, delta);CHKERRQ(ierr); 19 ierr = MatLMVMUpdate(bqnk->B, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 20 PetscFunctionReturn(0); 21 } 22 23 static PetscErrorCode TaoBQNLSComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type) 24 { 25 TAO_BNK *bnk = (TAO_BNK *)tao->data; 26 TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 27 PetscErrorCode ierr; 28 29 PetscFunctionBegin; 30 ierr = MatSolve(bqnk->B, tao->gradient, tao->stepdirection);CHKERRQ(ierr); 31 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 32 ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr); 33 *ksp_reason = KSP_CONVERGED_ATOL; 34 *step_type = BNK_BFGS; 35 PetscFunctionReturn(0); 36 } 37 38 static PetscErrorCode TaoSetFromOptions_BQNLS(PetscOptionItems *PetscOptionsObject,Tao tao) 39 { 40 TAO_BNK *bnk = (TAO_BNK *)tao->data; 41 TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 42 PetscErrorCode ierr; 43 KSPType ksp_type; 44 PetscBool is_spd; 45 46 PetscFunctionBegin; 47 ierr = PetscOptionsHead(PetscOptionsObject,"Quasi-Newton-Krylov method for bound constrained optimization");CHKERRQ(ierr); 48 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); 49 ierr = PetscOptionsReal("-tao_bqnls_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr); 50 ierr = PetscOptionsReal("-tao_bqnls_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr); 51 ierr = PetscOptionsReal("-tao_bqnls_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr); 52 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); 53 ierr = PetscOptionsTail();CHKERRQ(ierr); 54 ierr = TaoSetFromOptions(bnk->bncg);CHKERRQ(ierr); 55 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 56 ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 57 ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr); 58 bnk->is_nash = bnk->is_gltr = bnk->is_stcg = PETSC_FALSE; 59 ierr = MatSetFromOptions(bqnk->B);CHKERRQ(ierr); 60 ierr = MatGetOption(bqnk->B, MAT_SPD, &is_spd);CHKERRQ(ierr); 61 if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite"); 62 PetscFunctionReturn(0); 63 } 64 65 PETSC_EXTERN PetscErrorCode TaoCreate_BQNLS(Tao tao) 66 { 67 TAO_BNK *bnk; 68 TAO_BQNK *bqnk; 69 PetscErrorCode ierr; 70 71 PetscFunctionBegin; 72 ierr = TaoCreate_BQNK(tao);CHKERRQ(ierr); 73 ierr = KSPSetOptionsPrefix(tao->ksp, "unused");CHKERRQ(ierr); 74 tao->ops->solve = TaoSolve_BNLS; 75 tao->ops->setfromoptions = TaoSetFromOptions_BQNLS; 76 77 bnk = (TAO_BNK*)tao->data; 78 bnk->update_type = BNK_UPDATE_STEP; 79 bnk->computehessian = TaoBQNLSComputeHessian; 80 bnk->computestep = TaoBQNLSComputeStep; 81 82 bqnk = (TAO_BQNK*)bnk->ctx; 83 ierr = MatSetOptionsPrefix(bqnk->B, "tao_bqnls_");CHKERRQ(ierr); 84 ierr = MatSetType(bqnk->B, MATLMVMBFGS);CHKERRQ(ierr); 85 PetscFunctionReturn(0); 86 }