xref: /petsc/src/tao/bound/impls/bqnk/bqnk.c (revision f5766c097db85df4ee62e2acb750cb76ef9efd11)
1e0ed867bSAlp Dener #include <../src/tao/bound/impls/bqnk/bqnk.h>
2e0ed867bSAlp Dener #include <petscksp.h>
3e0ed867bSAlp Dener 
4e0ed867bSAlp Dener static PetscErrorCode TaoBQNKComputeHessian(Tao tao)
5e0ed867bSAlp Dener {
6e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
7e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
8e0ed867bSAlp Dener   PetscErrorCode ierr;
965f5217aSAlp Dener   PetscReal      gnorm2, delta;
10e0ed867bSAlp Dener 
11e0ed867bSAlp Dener   PetscFunctionBegin;
12e0ed867bSAlp Dener   /* Alias the LMVM matrix into the TAO hessian */
13e0ed867bSAlp Dener   if (tao->hessian) {
14e0ed867bSAlp Dener     ierr = MatDestroy(&tao->hessian);CHKERRQ(ierr);
15e0ed867bSAlp Dener   }
16e0ed867bSAlp Dener   if (tao->hessian_pre) {
17e0ed867bSAlp Dener     ierr = MatDestroy(&tao->hessian_pre);CHKERRQ(ierr);
18e0ed867bSAlp Dener   }
19e0ed867bSAlp Dener   ierr = PetscObjectReference((PetscObject)bqnk->B);CHKERRQ(ierr);
20e0ed867bSAlp Dener   tao->hessian = bqnk->B;
21e0ed867bSAlp Dener   ierr = PetscObjectReference((PetscObject)bqnk->B);CHKERRQ(ierr);
22e0ed867bSAlp Dener   tao->hessian_pre = bqnk->B;
23e0ed867bSAlp Dener   /* Update the Hessian with the latest solution */
24*f5766c09SAlp Dener   if (bqnk->is_spd) {
2565f5217aSAlp Dener     gnorm2 = bnk->gnorm*bnk->gnorm;
2665f5217aSAlp Dener     delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / PetscMax(gnorm2, PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0));
2765f5217aSAlp Dener     ierr = MatSymBrdnSetDelta(bqnk->B, delta);CHKERRQ(ierr);
28*f5766c09SAlp Dener   }
29e0ed867bSAlp Dener   ierr = MatLMVMUpdate(tao->hessian, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
30e0ed867bSAlp Dener   ierr = MatLMVMResetShift(tao->hessian);CHKERRQ(ierr);
31e0ed867bSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
324f4fdda4SAlp Dener   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
33e0ed867bSAlp Dener   if (bnk->active_idx) {
34e0ed867bSAlp Dener     ierr = MatCreateSubMatrixVirtual(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, &bnk->H_inactive);CHKERRQ(ierr);
354f4fdda4SAlp Dener     ierr = PCLMVMSetIS(bqnk->pc, bnk->inactive_idx);CHKERRQ(ierr);
36e0ed867bSAlp Dener   } else {
37e0ed867bSAlp Dener     ierr = PetscObjectReference((PetscObject)tao->hessian);CHKERRQ(ierr);
38e0ed867bSAlp Dener     bnk->H_inactive = tao->hessian;
394f4fdda4SAlp Dener     ierr = PCLMVMClearIS(bqnk->pc);CHKERRQ(ierr);
40e0ed867bSAlp Dener   }
414f4fdda4SAlp Dener   ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
424f4fdda4SAlp Dener   ierr = PetscObjectReference((PetscObject)bnk->H_inactive);CHKERRQ(ierr);
43e0ed867bSAlp Dener   bnk->Hpre_inactive = bnk->H_inactive;
44e0ed867bSAlp Dener   PetscFunctionReturn(0);
45e0ed867bSAlp Dener }
46e0ed867bSAlp Dener 
476b591159SAlp Dener static PetscErrorCode TaoBQNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
48e0ed867bSAlp Dener {
49e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
50e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
51e0ed867bSAlp Dener   PetscErrorCode ierr;
52e0ed867bSAlp Dener 
53e0ed867bSAlp Dener   PetscFunctionBegin;
546b591159SAlp Dener   ierr = TaoBNKComputeStep(tao, shift, ksp_reason, step_type);CHKERRQ(ierr);
55e0ed867bSAlp Dener   if (*ksp_reason < 0) {
56e0ed867bSAlp Dener     /* Krylov solver failed to converge so reset the LMVM matrix */
57e0ed867bSAlp Dener     ierr = MatLMVMReset(bqnk->B, PETSC_FALSE);CHKERRQ(ierr);
586b591159SAlp Dener     ierr = MatLMVMUpdate(bqnk->B, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
59e0ed867bSAlp Dener   }
60e0ed867bSAlp Dener   PetscFunctionReturn(0);
61e0ed867bSAlp Dener }
62e0ed867bSAlp Dener 
636b591159SAlp Dener PetscErrorCode TaoSetUp_BQNK(Tao tao)
644f4fdda4SAlp Dener {
654f4fdda4SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
664f4fdda4SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
674f4fdda4SAlp Dener   PetscErrorCode ierr;
684f4fdda4SAlp Dener   PetscInt       n, N;
696b591159SAlp Dener   PetscBool      is_lmvm, is_sym, is_spd;
704f4fdda4SAlp Dener 
714f4fdda4SAlp Dener   PetscFunctionBegin;
724f4fdda4SAlp Dener   ierr = TaoSetUp_BNK(tao);CHKERRQ(ierr);
734f4fdda4SAlp Dener   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
744f4fdda4SAlp Dener   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
754f4fdda4SAlp Dener   ierr = MatSetSizes(bqnk->B, n, n, N, N);CHKERRQ(ierr);
764f4fdda4SAlp Dener   ierr = MatLMVMAllocate(bqnk->B,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
774f4fdda4SAlp Dener   ierr = PetscObjectBaseTypeCompare((PetscObject)bqnk->B, MATLMVM, &is_lmvm);CHKERRQ(ierr);
784f4fdda4SAlp Dener   if (!is_lmvm) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Matrix must be an LMVM-type");
794f4fdda4SAlp Dener   ierr = MatGetOption(bqnk->B, MAT_SYMMETRIC, &is_sym);CHKERRQ(ierr);
804f4fdda4SAlp Dener   if (!is_sym) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric");
816b591159SAlp Dener   ierr = MatGetOption(bqnk->B, MAT_SPD, &is_spd);CHKERRQ(ierr);
824f4fdda4SAlp Dener   ierr = KSPGetPC(tao->ksp, &bqnk->pc);CHKERRQ(ierr);
834f4fdda4SAlp Dener   ierr = PCSetType(bqnk->pc, PCLMVM);CHKERRQ(ierr);
844f4fdda4SAlp Dener   ierr = PCLMVMSetMatLMVM(bqnk->pc, bqnk->B);CHKERRQ(ierr);
854f4fdda4SAlp Dener   PetscFunctionReturn(0);
864f4fdda4SAlp Dener }
874f4fdda4SAlp Dener 
88e0ed867bSAlp Dener static PetscErrorCode TaoSetFromOptions_BQNK(PetscOptionItems *PetscOptionsObject,Tao tao)
89e0ed867bSAlp Dener {
90e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
91e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
92e0ed867bSAlp Dener   PetscErrorCode ierr;
936b591159SAlp Dener   KSPType        ksp_type;
94e0ed867bSAlp Dener 
95e0ed867bSAlp Dener   PetscFunctionBegin;
964f4fdda4SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Quasi-Newton-Krylov method for bound constrained optimization");CHKERRQ(ierr);
976b591159SAlp Dener   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);
986b591159SAlp Dener   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);
996b591159SAlp Dener   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);
1006b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
1016b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
1026b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
1036b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
1046b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
1056b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
1066b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
1076b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
1086b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
1096b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
1106b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_eta1", "(developer) threshold for rejecting step (-tao_bqnk_update_type reduction)", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
1116b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_eta2", "(developer) threshold for accepting marginal step (-tao_bqnk_update_type reduction)", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
1126b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bqnk_update_type reduction)", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
1136b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_eta4", "(developer) threshold for accepting good step (-tao_bqnk_update_type reduction)", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
1146b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_alpha1", "(developer) radius reduction factor for rejected step (-tao_bqnk_update_type reduction)", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
1156b591159SAlp Dener   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);
1166b591159SAlp Dener   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);
1176b591159SAlp Dener   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);
1186b591159SAlp Dener   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);
1196b591159SAlp Dener   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);
1206b591159SAlp Dener   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);
1216b591159SAlp Dener   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);
1226b591159SAlp Dener   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);
1236b591159SAlp Dener   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);
1246b591159SAlp Dener   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);
1256b591159SAlp Dener   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);
1266b591159SAlp Dener   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);
1276b591159SAlp Dener   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);
1286b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_mu1", "(developer) threshold for accepting very good step (-tao_bqnk_update_type interpolation)", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
1296b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_mu2", "(developer) threshold for accepting good step (-tao_bqnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
1306b591159SAlp Dener   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);
1316b591159SAlp Dener   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);
1326b591159SAlp Dener   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);
1336b591159SAlp Dener   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);
1346b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_theta", "(developer) trust region interpolation factor (-tao_bqnk_update_type interpolation)", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
1356b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
1366b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
1376b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
1386b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr);
1396b591159SAlp Dener   ierr = PetscOptionsReal("-tao_bqnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr);
1406b591159SAlp Dener   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);
1414f4fdda4SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
1426b591159SAlp Dener   ierr = TaoSetFromOptions(bnk->bncg);CHKERRQ(ierr);
1436b591159SAlp Dener   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1446b591159SAlp Dener   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1456b591159SAlp Dener   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
1466b591159SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
1476b591159SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
1486b591159SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
149e0ed867bSAlp Dener   if (bnk->init_type == BNK_INIT_INTERPOLATION) bnk->init_type = BNK_INIT_DIRECTION;
150e0ed867bSAlp Dener   ierr = MatSetFromOptions(bqnk->B);CHKERRQ(ierr);
151*f5766c09SAlp Dener   ierr = MatGetOption(bqnk->B, MAT_SPD, &bqnk->is_spd);CHKERRQ(ierr);
152e0ed867bSAlp Dener   PetscFunctionReturn(0);
153e0ed867bSAlp Dener }
154e0ed867bSAlp Dener 
155e0ed867bSAlp Dener static PetscErrorCode TaoView_BQNK(Tao tao, PetscViewer viewer)
156e0ed867bSAlp Dener {
157e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
158e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
159e0ed867bSAlp Dener   PetscErrorCode ierr;
160e0ed867bSAlp Dener   PetscBool      isascii;
161e0ed867bSAlp Dener 
162e0ed867bSAlp Dener   PetscFunctionBegin;
163e0ed867bSAlp Dener   ierr = TaoView_BNK(tao, viewer);CHKERRQ(ierr);
164e0ed867bSAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
165e0ed867bSAlp Dener   if (isascii) {
166e0ed867bSAlp Dener     ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
167e0ed867bSAlp Dener     ierr = MatView(bqnk->B, viewer);CHKERRQ(ierr);
168e0ed867bSAlp Dener     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
169e0ed867bSAlp Dener   }
170e0ed867bSAlp Dener   PetscFunctionReturn(0);
171e0ed867bSAlp Dener }
172e0ed867bSAlp Dener 
173e0ed867bSAlp Dener static PetscErrorCode TaoDestroy_BQNK(Tao tao)
174e0ed867bSAlp Dener {
175e0ed867bSAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
176e0ed867bSAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
177e0ed867bSAlp Dener   PetscErrorCode ierr;
178e0ed867bSAlp Dener 
179e0ed867bSAlp Dener   PetscFunctionBegin;
1804f4fdda4SAlp Dener   ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
181cb384e1eSAlp Dener   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
182e0ed867bSAlp Dener   ierr = MatDestroy(&bqnk->B);CHKERRQ(ierr);
183e0ed867bSAlp Dener   ierr = PetscFree(bnk->ctx);CHKERRQ(ierr);
184e0ed867bSAlp Dener   ierr = TaoDestroy_BNK(tao);CHKERRQ(ierr);
185e0ed867bSAlp Dener   PetscFunctionReturn(0);
186e0ed867bSAlp Dener }
187e0ed867bSAlp Dener 
188e0ed867bSAlp Dener PETSC_INTERN PetscErrorCode TaoCreate_BQNK(Tao tao)
189e0ed867bSAlp Dener {
190e0ed867bSAlp Dener   TAO_BNK        *bnk;
191e0ed867bSAlp Dener   TAO_BQNK       *bqnk;
192e0ed867bSAlp Dener   PetscErrorCode ierr;
193e0ed867bSAlp Dener 
194e0ed867bSAlp Dener   PetscFunctionBegin;
195e0ed867bSAlp Dener   ierr = TaoCreate_BNK(tao);CHKERRQ(ierr);
196f4a59f9fSAlp Dener   ierr = KSPSetOptionsPrefix(tao->ksp, "tao_bqnk_");CHKERRQ(ierr);
197e0ed867bSAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BQNK;
198e0ed867bSAlp Dener   tao->ops->destroy = TaoDestroy_BQNK;
199e0ed867bSAlp Dener   tao->ops->view = TaoView_BQNK;
2004f4fdda4SAlp Dener   tao->ops->setup = TaoSetUp_BQNK;
201e0ed867bSAlp Dener 
202e0ed867bSAlp Dener   bnk = (TAO_BNK *)tao->data;
203e0ed867bSAlp Dener   bnk->computehessian = TaoBQNKComputeHessian;
204e0ed867bSAlp Dener   bnk->computestep = TaoBQNKComputeStep;
205e0ed867bSAlp Dener   bnk->init_type = BNK_INIT_DIRECTION;
206e0ed867bSAlp Dener 
207e0ed867bSAlp Dener   ierr = PetscNewLog(tao,&bqnk);CHKERRQ(ierr);
208e0ed867bSAlp Dener   bnk->ctx = (void*)bqnk;
209*f5766c09SAlp Dener   bqnk->is_spd = PETSC_TRUE;
210e0ed867bSAlp Dener 
211e0ed867bSAlp Dener   ierr = MatCreate(PetscObjectComm((PetscObject)tao), &bqnk->B);CHKERRQ(ierr);
212e0ed867bSAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)bqnk->B, (PetscObject)tao, 1);CHKERRQ(ierr);
213e0ed867bSAlp Dener   ierr = MatSetOptionsPrefix(bqnk->B, "tao_bqnk_");CHKERRQ(ierr);
214e0ed867bSAlp Dener   ierr = MatSetType(bqnk->B, MATLMVMSR1);CHKERRQ(ierr);
215e0ed867bSAlp Dener   PetscFunctionReturn(0);
216e0ed867bSAlp Dener }
217*f5766c09SAlp Dener 
218*f5766c09SAlp Dener PetscErrorCode TaoGetLMVMMatrix(Tao tao, Mat *B)
219*f5766c09SAlp Dener {
220*f5766c09SAlp Dener   TAO_BNK        *bnk = (TAO_BNK*)tao->data;
221*f5766c09SAlp Dener   TAO_BQNK       *bqnk = (TAO_BQNK*)bnk->ctx;
222*f5766c09SAlp Dener   PetscErrorCode ierr;
223*f5766c09SAlp Dener   PetscBool      is_bqnls, is_bqnkls, is_bqnktr, is_bqnktl;
224*f5766c09SAlp Dener 
225*f5766c09SAlp Dener   PetscFunctionBegin;
226*f5766c09SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)tao, TAOBQNLS, &is_bqnls);CHKERRQ(ierr);
227*f5766c09SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)tao, TAOBQNKLS, &is_bqnkls);CHKERRQ(ierr);
228*f5766c09SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)tao, TAOBQNKTR, &is_bqnktr);CHKERRQ(ierr);
229*f5766c09SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)tao, TAOBQNKTL, &is_bqnktl);CHKERRQ(ierr);
230*f5766c09SAlp Dener   if (!is_bqnls && !is_bqnkls && !is_bqnktr && is_bqnktl) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM Matrix only exists for quasi-Newton algorithms");
231*f5766c09SAlp Dener   *B = bqnk->B;
232*f5766c09SAlp Dener   PetscFunctionReturn(0);
233*f5766c09SAlp Dener }