1 #include <../src/tao/bound/impls/bqnk/bqnk.h> 2 #include <petscksp.h> 3 4 static PetscErrorCode TaoBQNKComputeHessian(Tao tao) 5 { 6 TAO_BNK *bnk = (TAO_BNK *)tao->data; 7 TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 8 PetscErrorCode ierr; 9 PetscReal gnorm2, delta; 10 11 PetscFunctionBegin; 12 /* Alias the LMVM matrix into the TAO hessian */ 13 if (tao->hessian) { 14 ierr = MatDestroy(&tao->hessian);CHKERRQ(ierr); 15 } 16 if (tao->hessian_pre) { 17 ierr = MatDestroy(&tao->hessian_pre);CHKERRQ(ierr); 18 } 19 ierr = PetscObjectReference((PetscObject)bqnk->B);CHKERRQ(ierr); 20 tao->hessian = bqnk->B; 21 ierr = PetscObjectReference((PetscObject)bqnk->B);CHKERRQ(ierr); 22 tao->hessian_pre = bqnk->B; 23 /* Update the Hessian with the latest solution */ 24 gnorm2 = bnk->gnorm*bnk->gnorm; 25 delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / PetscMax(gnorm2, PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0)); 26 ierr = MatSymBrdnSetDelta(bqnk->B, delta);CHKERRQ(ierr); 27 ierr = MatLMVMUpdate(tao->hessian, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 28 ierr = MatLMVMResetShift(tao->hessian);CHKERRQ(ierr); 29 /* Prepare the reduced sub-matrices for the inactive set */ 30 ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr); 31 if (bnk->active_idx) { 32 ierr = MatCreateSubMatrixVirtual(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, &bnk->H_inactive);CHKERRQ(ierr); 33 ierr = PCLMVMSetIS(bqnk->pc, bnk->inactive_idx);CHKERRQ(ierr); 34 } else { 35 ierr = PetscObjectReference((PetscObject)tao->hessian);CHKERRQ(ierr); 36 bnk->H_inactive = tao->hessian; 37 ierr = PCLMVMClearIS(bqnk->pc);CHKERRQ(ierr); 38 } 39 ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr); 40 ierr = PetscObjectReference((PetscObject)bnk->H_inactive);CHKERRQ(ierr); 41 bnk->Hpre_inactive = bnk->H_inactive; 42 PetscFunctionReturn(0); 43 } 44 45 static PetscErrorCode TaoBQNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type) 46 { 47 TAO_BNK *bnk = (TAO_BNK *)tao->data; 48 TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 49 PetscErrorCode ierr; 50 51 PetscFunctionBegin; 52 ierr = TaoBNKComputeStep(tao, shift, ksp_reason, step_type);CHKERRQ(ierr); 53 if (*ksp_reason < 0) { 54 /* Krylov solver failed to converge so reset the LMVM matrix */ 55 ierr = MatLMVMReset(bqnk->B, PETSC_FALSE);CHKERRQ(ierr); 56 ierr = MatLMVMUpdate(bqnk->B, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 57 } 58 PetscFunctionReturn(0); 59 } 60 61 PetscErrorCode TaoSetUp_BQNK(Tao tao) 62 { 63 TAO_BNK *bnk = (TAO_BNK *)tao->data; 64 TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 65 PetscErrorCode ierr; 66 PetscInt n, N; 67 PetscBool is_lmvm, is_sym, is_spd; 68 69 PetscFunctionBegin; 70 ierr = TaoSetUp_BNK(tao);CHKERRQ(ierr); 71 ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 72 ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 73 ierr = MatSetSizes(bqnk->B, n, n, N, N);CHKERRQ(ierr); 74 ierr = MatLMVMAllocate(bqnk->B,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr); 75 ierr = PetscObjectBaseTypeCompare((PetscObject)bqnk->B, MATLMVM, &is_lmvm);CHKERRQ(ierr); 76 if (!is_lmvm) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Matrix must be an LMVM-type"); 77 ierr = MatGetOption(bqnk->B, MAT_SYMMETRIC, &is_sym);CHKERRQ(ierr); 78 if (!is_sym) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric"); 79 ierr = MatGetOption(bqnk->B, MAT_SPD, &is_spd);CHKERRQ(ierr); 80 ierr = KSPGetPC(tao->ksp, &bqnk->pc);CHKERRQ(ierr); 81 ierr = PCSetType(bqnk->pc, PCLMVM);CHKERRQ(ierr); 82 ierr = PCLMVMSetMatLMVM(bqnk->pc, bqnk->B);CHKERRQ(ierr); 83 PetscFunctionReturn(0); 84 } 85 86 static PetscErrorCode TaoSetFromOptions_BQNK(PetscOptionItems *PetscOptionsObject,Tao tao) 87 { 88 TAO_BNK *bnk = (TAO_BNK *)tao->data; 89 TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 90 PetscErrorCode ierr; 91 KSPType ksp_type; 92 93 PetscFunctionBegin; 94 ierr = PetscOptionsHead(PetscOptionsObject,"Quasi-Newton-Krylov method for bound constrained optimization");CHKERRQ(ierr); 95 ierr = PetscOptionsEList("-tao_bqnk_init_type", "radius initialization type", "", BQNK_INIT, BQNK_INIT_TYPES, BQNK_INIT[bnk->init_type], &bnk->init_type, 0);CHKERRQ(ierr); 96 ierr = PetscOptionsEList("-tao_bqnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, 0);CHKERRQ(ierr); 97 ierr = PetscOptionsEList("-tao_bqnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, 0);CHKERRQ(ierr); 98 ierr = PetscOptionsReal("-tao_bqnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr); 99 ierr = PetscOptionsReal("-tao_bqnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr); 100 ierr = PetscOptionsReal("-tao_bqnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr); 101 ierr = PetscOptionsReal("-tao_bqnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr); 102 ierr = PetscOptionsReal("-tao_bqnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr); 103 ierr = PetscOptionsReal("-tao_bqnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr); 104 ierr = PetscOptionsReal("-tao_bqnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr); 105 ierr = PetscOptionsReal("-tao_bqnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr); 106 ierr = PetscOptionsReal("-tao_bqnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr); 107 ierr = PetscOptionsReal("-tao_bqnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr); 108 ierr = PetscOptionsReal("-tao_bqnk_eta1", "(developer) threshold for rejecting step (-tao_bqnk_update_type reduction)", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr); 109 ierr = PetscOptionsReal("-tao_bqnk_eta2", "(developer) threshold for accepting marginal step (-tao_bqnk_update_type reduction)", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr); 110 ierr = PetscOptionsReal("-tao_bqnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bqnk_update_type reduction)", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr); 111 ierr = PetscOptionsReal("-tao_bqnk_eta4", "(developer) threshold for accepting good step (-tao_bqnk_update_type reduction)", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr); 112 ierr = PetscOptionsReal("-tao_bqnk_alpha1", "(developer) radius reduction factor for rejected step (-tao_bqnk_update_type reduction)", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr); 113 ierr = PetscOptionsReal("-tao_bqnk_alpha2", "(developer) radius reduction factor for marginally accepted bad step (-tao_bqnk_update_type reduction)", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr); 114 ierr = PetscOptionsReal("-tao_bqnk_alpha3", "(developer) radius increase factor for reasonable accepted step (-tao_bqnk_update_type reduction)", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr); 115 ierr = PetscOptionsReal("-tao_bqnk_alpha4", "(developer) radius increase factor for good accepted step (-tao_bqnk_update_type reduction)", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr); 116 ierr = PetscOptionsReal("-tao_bqnk_alpha5", "(developer) radius increase factor for very good accepted step (-tao_bqnk_update_type reduction)", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr); 117 ierr = PetscOptionsReal("-tao_bqnk_nu1", "(developer) threshold for small line-search step length (-tao_bqnk_update_type step)", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr); 118 ierr = PetscOptionsReal("-tao_bqnk_nu2", "(developer) threshold for reasonable line-search step length (-tao_bqnk_update_type step)", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr); 119 ierr = PetscOptionsReal("-tao_bqnk_nu3", "(developer) threshold for large line-search step length (-tao_bqnk_update_type step)", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr); 120 ierr = PetscOptionsReal("-tao_bqnk_nu4", "(developer) threshold for very large line-search step length (-tao_bqnk_update_type step)", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr); 121 ierr = PetscOptionsReal("-tao_bqnk_omega1", "(developer) radius reduction factor for very small line-search step length (-tao_bqnk_update_type step)", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr); 122 ierr = PetscOptionsReal("-tao_bqnk_omega2", "(developer) radius reduction factor for small line-search step length (-tao_bqnk_update_type step)", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr); 123 ierr = PetscOptionsReal("-tao_bqnk_omega3", "(developer) radius factor for decent line-search step length (-tao_bqnk_update_type step)", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr); 124 ierr = PetscOptionsReal("-tao_bqnk_omega4", "(developer) radius increase factor for large line-search step length (-tao_bqnk_update_type step)", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr); 125 ierr = PetscOptionsReal("-tao_bqnk_omega5", "(developer) radius increase factor for very large line-search step length (-tao_bqnk_update_type step)", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr); 126 ierr = PetscOptionsReal("-tao_bqnk_mu1", "(developer) threshold for accepting very good step (-tao_bqnk_update_type interpolation)", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr); 127 ierr = PetscOptionsReal("-tao_bqnk_mu2", "(developer) threshold for accepting good step (-tao_bqnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr); 128 ierr = PetscOptionsReal("-tao_bqnk_gamma1", "(developer) radius reduction factor for rejected very bad step (-tao_bqnk_update_type interpolation)", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr); 129 ierr = PetscOptionsReal("-tao_bqnk_gamma2", "(developer) radius reduction factor for rejected bad step (-tao_bqnk_update_type interpolation)", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr); 130 ierr = PetscOptionsReal("-tao_bqnk_gamma3", "(developer) radius increase factor for accepted good step (-tao_bqnk_update_type interpolation)", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr); 131 ierr = PetscOptionsReal("-tao_bqnk_gamma4", "(developer) radius increase factor for accepted very good step (-tao_bqnk_update_type interpolation)", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr); 132 ierr = PetscOptionsReal("-tao_bqnk_theta", "(developer) trust region interpolation factor (-tao_bqnk_update_type interpolation)", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr); 133 ierr = PetscOptionsReal("-tao_bqnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr); 134 ierr = PetscOptionsReal("-tao_bqnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr); 135 ierr = PetscOptionsReal("-tao_bqnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr); 136 ierr = PetscOptionsReal("-tao_bqnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr); 137 ierr = PetscOptionsReal("-tao_bqnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr); 138 ierr = PetscOptionsInt("-tao_bqnk_max_cg_its", "number of BNCG iterations to take for each Newton step", "", bnk->max_cg_its, &bnk->max_cg_its,NULL);CHKERRQ(ierr); 139 ierr = PetscOptionsTail();CHKERRQ(ierr); 140 ierr = TaoSetFromOptions(bnk->bncg);CHKERRQ(ierr); 141 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 142 ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 143 ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr); 144 ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr); 145 ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr); 146 ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr); 147 if (bnk->init_type == BNK_INIT_INTERPOLATION) bnk->init_type = BNK_INIT_DIRECTION; 148 ierr = MatSetFromOptions(bqnk->B);CHKERRQ(ierr); 149 PetscFunctionReturn(0); 150 } 151 152 static PetscErrorCode TaoView_BQNK(Tao tao, PetscViewer viewer) 153 { 154 TAO_BNK *bnk = (TAO_BNK*)tao->data; 155 TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 156 PetscErrorCode ierr; 157 PetscBool isascii; 158 159 PetscFunctionBegin; 160 ierr = TaoView_BNK(tao, viewer);CHKERRQ(ierr); 161 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 162 if (isascii) { 163 ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 164 ierr = MatView(bqnk->B, viewer);CHKERRQ(ierr); 165 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 166 } 167 PetscFunctionReturn(0); 168 } 169 170 static PetscErrorCode TaoDestroy_BQNK(Tao tao) 171 { 172 TAO_BNK *bnk = (TAO_BNK*)tao->data; 173 TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 174 PetscErrorCode ierr; 175 176 PetscFunctionBegin; 177 ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr); 178 ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr); 179 ierr = MatDestroy(&bqnk->B);CHKERRQ(ierr); 180 ierr = PetscFree(bnk->ctx);CHKERRQ(ierr); 181 ierr = TaoDestroy_BNK(tao);CHKERRQ(ierr); 182 PetscFunctionReturn(0); 183 } 184 185 PETSC_INTERN PetscErrorCode TaoCreate_BQNK(Tao tao) 186 { 187 TAO_BNK *bnk; 188 TAO_BQNK *bqnk; 189 PetscErrorCode ierr; 190 191 PetscFunctionBegin; 192 ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 193 ierr = KSPSetOptionsPrefix(tao->ksp, "tao_bqnk_");CHKERRQ(ierr); 194 tao->ops->setfromoptions = TaoSetFromOptions_BQNK; 195 tao->ops->destroy = TaoDestroy_BQNK; 196 tao->ops->view = TaoView_BQNK; 197 tao->ops->setup = TaoSetUp_BQNK; 198 199 bnk = (TAO_BNK *)tao->data; 200 bnk->computehessian = TaoBQNKComputeHessian; 201 bnk->computestep = TaoBQNKComputeStep; 202 bnk->init_type = BNK_INIT_DIRECTION; 203 204 ierr = PetscNewLog(tao,&bqnk);CHKERRQ(ierr); 205 bnk->ctx = (void*)bqnk; 206 207 ierr = MatCreate(PetscObjectComm((PetscObject)tao), &bqnk->B);CHKERRQ(ierr); 208 ierr = PetscObjectIncrementTabLevel((PetscObject)bqnk->B, (PetscObject)tao, 1);CHKERRQ(ierr); 209 ierr = MatSetOptionsPrefix(bqnk->B, "tao_bqnk_");CHKERRQ(ierr); 210 ierr = MatSetType(bqnk->B, MATLMVMSR1);CHKERRQ(ierr); 211 PetscFunctionReturn(0); 212 }