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