1414d97d3SAlp Dener #include <../src/tao/bound/impls/bqnk/bqnk.h> /*I "petsctao.h" I*/ /*I "petscmat.h" I*/ 2e0ed867bSAlp Dener #include <petscksp.h> 3e0ed867bSAlp Dener 470a3f44bSAlp Dener static const char *BQNK_INIT[64] = {"constant", "direction"}; 570a3f44bSAlp Dener static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"}; 670a3f44bSAlp Dener static const char *BNK_AS[64] = {"none", "bertsekas"}; 770a3f44bSAlp Dener 8e0ed867bSAlp Dener static PetscErrorCode TaoBQNKComputeHessian(Tao tao) 9e0ed867bSAlp Dener { 10e0ed867bSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 11e0ed867bSAlp Dener TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 12e0ed867bSAlp Dener PetscErrorCode ierr; 1365f5217aSAlp Dener PetscReal gnorm2, delta; 14e0ed867bSAlp Dener 15e0ed867bSAlp Dener PetscFunctionBegin; 16e0ed867bSAlp Dener /* Alias the LMVM matrix into the TAO hessian */ 17e0ed867bSAlp Dener if (tao->hessian) { 18e0ed867bSAlp Dener ierr = MatDestroy(&tao->hessian);CHKERRQ(ierr); 19e0ed867bSAlp Dener } 20e0ed867bSAlp Dener if (tao->hessian_pre) { 21e0ed867bSAlp Dener ierr = MatDestroy(&tao->hessian_pre);CHKERRQ(ierr); 22e0ed867bSAlp Dener } 23e0ed867bSAlp Dener ierr = PetscObjectReference((PetscObject)bqnk->B);CHKERRQ(ierr); 24e0ed867bSAlp Dener tao->hessian = bqnk->B; 25e0ed867bSAlp Dener ierr = PetscObjectReference((PetscObject)bqnk->B);CHKERRQ(ierr); 26e0ed867bSAlp Dener tao->hessian_pre = bqnk->B; 27e0ed867bSAlp Dener /* Update the Hessian with the latest solution */ 28f5766c09SAlp Dener if (bqnk->is_spd) { 2965f5217aSAlp Dener gnorm2 = bnk->gnorm*bnk->gnorm; 308cabe928SAlp Dener if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON; 318cabe928SAlp Dener if (bnk->f == 0.0) { 328cabe928SAlp Dener delta = 2.0 / gnorm2; 338cabe928SAlp Dener } else { 348cabe928SAlp Dener delta = 2.0 * PetscAbsScalar(bnk->f) / gnorm2; 358cabe928SAlp Dener } 36864588a7SAlp Dener ierr = MatLMVMSymBroydenSetDelta(bqnk->B, delta);CHKERRQ(ierr); 37f5766c09SAlp Dener } 38e0ed867bSAlp Dener ierr = MatLMVMUpdate(tao->hessian, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 39e0ed867bSAlp Dener ierr = MatLMVMResetShift(tao->hessian);CHKERRQ(ierr); 40e0ed867bSAlp Dener /* Prepare the reduced sub-matrices for the inactive set */ 414f4fdda4SAlp Dener ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr); 42e0ed867bSAlp Dener if (bnk->active_idx) { 43e0ed867bSAlp Dener ierr = MatCreateSubMatrixVirtual(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, &bnk->H_inactive);CHKERRQ(ierr); 444f4fdda4SAlp Dener ierr = PCLMVMSetIS(bqnk->pc, bnk->inactive_idx);CHKERRQ(ierr); 45e0ed867bSAlp Dener } else { 46e0ed867bSAlp Dener ierr = PetscObjectReference((PetscObject)tao->hessian);CHKERRQ(ierr); 47e0ed867bSAlp Dener bnk->H_inactive = tao->hessian; 484f4fdda4SAlp Dener ierr = PCLMVMClearIS(bqnk->pc);CHKERRQ(ierr); 49e0ed867bSAlp Dener } 504f4fdda4SAlp Dener ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr); 514f4fdda4SAlp Dener ierr = PetscObjectReference((PetscObject)bnk->H_inactive);CHKERRQ(ierr); 52e0ed867bSAlp Dener bnk->Hpre_inactive = bnk->H_inactive; 53e0ed867bSAlp Dener PetscFunctionReturn(0); 54e0ed867bSAlp Dener } 55e0ed867bSAlp Dener 566b591159SAlp Dener static PetscErrorCode TaoBQNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type) 57e0ed867bSAlp Dener { 58e0ed867bSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 59e0ed867bSAlp Dener TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 60e0ed867bSAlp Dener PetscErrorCode ierr; 61e0ed867bSAlp Dener 62e0ed867bSAlp Dener PetscFunctionBegin; 636b591159SAlp Dener ierr = TaoBNKComputeStep(tao, shift, ksp_reason, step_type);CHKERRQ(ierr); 64e0ed867bSAlp Dener if (*ksp_reason < 0) { 65e0ed867bSAlp Dener /* Krylov solver failed to converge so reset the LMVM matrix */ 66e0ed867bSAlp Dener ierr = MatLMVMReset(bqnk->B, PETSC_FALSE);CHKERRQ(ierr); 676b591159SAlp Dener ierr = MatLMVMUpdate(bqnk->B, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr); 68e0ed867bSAlp Dener } 69e0ed867bSAlp Dener PetscFunctionReturn(0); 70e0ed867bSAlp Dener } 71e0ed867bSAlp Dener 72414d97d3SAlp Dener PetscErrorCode TaoSolve_BQNK(Tao tao) 73414d97d3SAlp Dener { 74414d97d3SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 75414d97d3SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 76414d97d3SAlp Dener Mat_LMVM *lmvm = (Mat_LMVM*)bqnk->B->data; 77414d97d3SAlp Dener Mat_LMVM *J0; 78414d97d3SAlp Dener Mat_SymBrdn *diag_ctx; 79414d97d3SAlp Dener PetscBool flg = PETSC_FALSE; 80414d97d3SAlp Dener PetscErrorCode ierr; 81414d97d3SAlp Dener 82414d97d3SAlp Dener PetscFunctionBegin; 83414d97d3SAlp Dener if (!tao->recycle) { 84414d97d3SAlp Dener ierr = MatLMVMReset(bqnk->B, PETSC_FALSE);CHKERRQ(ierr); 85414d97d3SAlp Dener lmvm->nresets = 0; 86414d97d3SAlp Dener if (lmvm->J0) { 87414d97d3SAlp Dener ierr = PetscObjectBaseTypeCompare((PetscObject)lmvm->J0, MATLMVM, &flg);CHKERRQ(ierr); 88414d97d3SAlp Dener if (flg) { 89414d97d3SAlp Dener J0 = (Mat_LMVM*)lmvm->J0->data; 90414d97d3SAlp Dener J0->nresets = 0; 91414d97d3SAlp Dener } 92414d97d3SAlp Dener } 93414d97d3SAlp Dener flg = PETSC_FALSE; 94*1e1ea65dSPierre Jolivet ierr = PetscObjectTypeCompareAny((PetscObject)bqnk->B, &flg, MATLMVMSYMBROYDEN, MATLMVMSYMBADBROYDEN, MATLMVMBFGS, MATLMVMDFP, "");CHKERRQ(ierr); 95414d97d3SAlp Dener if (flg) { 96414d97d3SAlp Dener diag_ctx = (Mat_SymBrdn*)lmvm->ctx; 97414d97d3SAlp Dener J0 = (Mat_LMVM*)diag_ctx->D->data; 98414d97d3SAlp Dener J0->nresets = 0; 99414d97d3SAlp Dener } 100414d97d3SAlp Dener } 101414d97d3SAlp Dener ierr = (*bqnk->solve)(tao);CHKERRQ(ierr); 102414d97d3SAlp Dener PetscFunctionReturn(0); 103414d97d3SAlp Dener } 104414d97d3SAlp Dener 1056b591159SAlp Dener PetscErrorCode TaoSetUp_BQNK(Tao tao) 1064f4fdda4SAlp Dener { 1074f4fdda4SAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 1084f4fdda4SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 1094f4fdda4SAlp Dener PetscErrorCode ierr; 1104f4fdda4SAlp Dener PetscInt n, N; 1116b591159SAlp Dener PetscBool is_lmvm, is_sym, is_spd; 1124f4fdda4SAlp Dener 1134f4fdda4SAlp Dener PetscFunctionBegin; 1144f4fdda4SAlp Dener ierr = TaoSetUp_BNK(tao);CHKERRQ(ierr); 1154f4fdda4SAlp Dener ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 1164f4fdda4SAlp Dener ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 1174f4fdda4SAlp Dener ierr = MatSetSizes(bqnk->B, n, n, N, N);CHKERRQ(ierr); 1184f4fdda4SAlp Dener ierr = MatLMVMAllocate(bqnk->B,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr); 1194f4fdda4SAlp Dener ierr = PetscObjectBaseTypeCompare((PetscObject)bqnk->B, MATLMVM, &is_lmvm);CHKERRQ(ierr); 1204f4fdda4SAlp Dener if (!is_lmvm) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Matrix must be an LMVM-type"); 1214f4fdda4SAlp Dener ierr = MatGetOption(bqnk->B, MAT_SYMMETRIC, &is_sym);CHKERRQ(ierr); 1224f4fdda4SAlp Dener if (!is_sym) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric"); 1236b591159SAlp Dener ierr = MatGetOption(bqnk->B, MAT_SPD, &is_spd);CHKERRQ(ierr); 1244f4fdda4SAlp Dener ierr = KSPGetPC(tao->ksp, &bqnk->pc);CHKERRQ(ierr); 1254f4fdda4SAlp Dener ierr = PCSetType(bqnk->pc, PCLMVM);CHKERRQ(ierr); 1264f4fdda4SAlp Dener ierr = PCLMVMSetMatLMVM(bqnk->pc, bqnk->B);CHKERRQ(ierr); 1274f4fdda4SAlp Dener PetscFunctionReturn(0); 1284f4fdda4SAlp Dener } 1294f4fdda4SAlp Dener 130e0ed867bSAlp Dener static PetscErrorCode TaoSetFromOptions_BQNK(PetscOptionItems *PetscOptionsObject,Tao tao) 131e0ed867bSAlp Dener { 132e0ed867bSAlp Dener TAO_BNK *bnk = (TAO_BNK *)tao->data; 133e0ed867bSAlp Dener TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 134e0ed867bSAlp Dener PetscErrorCode ierr; 135e0ed867bSAlp Dener 136e0ed867bSAlp Dener PetscFunctionBegin; 1374f4fdda4SAlp Dener ierr = PetscOptionsHead(PetscOptionsObject,"Quasi-Newton-Krylov method for bound constrained optimization");CHKERRQ(ierr); 13883c8fe1dSLisandro Dalcin ierr = PetscOptionsEList("-tao_bqnk_init_type", "radius initialization type", "", BQNK_INIT, BQNK_INIT_TYPES, BQNK_INIT[bnk->init_type], &bnk->init_type, NULL);CHKERRQ(ierr); 13983c8fe1dSLisandro Dalcin ierr = PetscOptionsEList("-tao_bqnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, NULL);CHKERRQ(ierr); 14083c8fe1dSLisandro Dalcin ierr = PetscOptionsEList("-tao_bqnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, NULL);CHKERRQ(ierr); 1416b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr); 1426b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr); 1436b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr); 1446b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr); 1456b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr); 1466b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr); 1476b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr); 1486b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr); 1496b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr); 1506b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr); 1516b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_eta1", "(developer) threshold for rejecting step (-tao_bqnk_update_type reduction)", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr); 1526b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_eta2", "(developer) threshold for accepting marginal step (-tao_bqnk_update_type reduction)", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr); 1536b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bqnk_update_type reduction)", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr); 1546b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_eta4", "(developer) threshold for accepting good step (-tao_bqnk_update_type reduction)", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr); 1556b591159SAlp 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); 1566b591159SAlp 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); 1576b591159SAlp 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); 1586b591159SAlp 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); 1596b591159SAlp 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); 1606b591159SAlp 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); 1616b591159SAlp 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); 1626b591159SAlp 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); 1636b591159SAlp 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); 1646b591159SAlp 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); 1656b591159SAlp 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); 1666b591159SAlp 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); 1676b591159SAlp 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); 1686b591159SAlp 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); 1696b591159SAlp 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); 1706b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_mu2", "(developer) threshold for accepting good step (-tao_bqnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr); 1716b591159SAlp 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); 1726b591159SAlp 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); 1736b591159SAlp 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); 1746b591159SAlp 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); 1756b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_theta", "(developer) trust region interpolation factor (-tao_bqnk_update_type interpolation)", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr); 1766b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr); 1776b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr); 1786b591159SAlp Dener ierr = PetscOptionsReal("-tao_bqnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr); 1796b591159SAlp 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); 1806b591159SAlp 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); 1816b591159SAlp 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); 1824f4fdda4SAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 1836b591159SAlp Dener ierr = TaoSetFromOptions(bnk->bncg);CHKERRQ(ierr); 1846b591159SAlp Dener ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 1856b591159SAlp Dener ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 186e0ed867bSAlp Dener if (bnk->init_type == BNK_INIT_INTERPOLATION) bnk->init_type = BNK_INIT_DIRECTION; 187e0ed867bSAlp Dener ierr = MatSetFromOptions(bqnk->B);CHKERRQ(ierr); 188f5766c09SAlp Dener ierr = MatGetOption(bqnk->B, MAT_SPD, &bqnk->is_spd);CHKERRQ(ierr); 189e0ed867bSAlp Dener PetscFunctionReturn(0); 190e0ed867bSAlp Dener } 191e0ed867bSAlp Dener 192e0ed867bSAlp Dener static PetscErrorCode TaoView_BQNK(Tao tao, PetscViewer viewer) 193e0ed867bSAlp Dener { 194e0ed867bSAlp Dener TAO_BNK *bnk = (TAO_BNK*)tao->data; 195e0ed867bSAlp Dener TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 196e0ed867bSAlp Dener PetscErrorCode ierr; 197e0ed867bSAlp Dener PetscBool isascii; 198e0ed867bSAlp Dener 199e0ed867bSAlp Dener PetscFunctionBegin; 200e0ed867bSAlp Dener ierr = TaoView_BNK(tao, viewer);CHKERRQ(ierr); 201e0ed867bSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 202e0ed867bSAlp Dener if (isascii) { 203e0ed867bSAlp Dener ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 204e0ed867bSAlp Dener ierr = MatView(bqnk->B, viewer);CHKERRQ(ierr); 205e0ed867bSAlp Dener ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 206e0ed867bSAlp Dener } 207e0ed867bSAlp Dener PetscFunctionReturn(0); 208e0ed867bSAlp Dener } 209e0ed867bSAlp Dener 210e0ed867bSAlp Dener static PetscErrorCode TaoDestroy_BQNK(Tao tao) 211e0ed867bSAlp Dener { 212e0ed867bSAlp Dener TAO_BNK *bnk = (TAO_BNK*)tao->data; 213e0ed867bSAlp Dener TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 214e0ed867bSAlp Dener PetscErrorCode ierr; 215e0ed867bSAlp Dener 216e0ed867bSAlp Dener PetscFunctionBegin; 2174f4fdda4SAlp Dener ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr); 218cb384e1eSAlp Dener ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr); 219e0ed867bSAlp Dener ierr = MatDestroy(&bqnk->B);CHKERRQ(ierr); 220e0ed867bSAlp Dener ierr = PetscFree(bnk->ctx);CHKERRQ(ierr); 221e0ed867bSAlp Dener ierr = TaoDestroy_BNK(tao);CHKERRQ(ierr); 222e0ed867bSAlp Dener PetscFunctionReturn(0); 223e0ed867bSAlp Dener } 224e0ed867bSAlp Dener 225e0ed867bSAlp Dener PETSC_INTERN PetscErrorCode TaoCreate_BQNK(Tao tao) 226e0ed867bSAlp Dener { 227e0ed867bSAlp Dener TAO_BNK *bnk; 228e0ed867bSAlp Dener TAO_BQNK *bqnk; 229e0ed867bSAlp Dener PetscErrorCode ierr; 230e0ed867bSAlp Dener 231e0ed867bSAlp Dener PetscFunctionBegin; 232e0ed867bSAlp Dener ierr = TaoCreate_BNK(tao);CHKERRQ(ierr); 233f4a59f9fSAlp Dener ierr = KSPSetOptionsPrefix(tao->ksp, "tao_bqnk_");CHKERRQ(ierr); 234414d97d3SAlp Dener tao->ops->solve = TaoSolve_BQNK; 235e0ed867bSAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BQNK; 236e0ed867bSAlp Dener tao->ops->destroy = TaoDestroy_BQNK; 237e0ed867bSAlp Dener tao->ops->view = TaoView_BQNK; 2384f4fdda4SAlp Dener tao->ops->setup = TaoSetUp_BQNK; 239e0ed867bSAlp Dener 240e0ed867bSAlp Dener bnk = (TAO_BNK *)tao->data; 241e0ed867bSAlp Dener bnk->computehessian = TaoBQNKComputeHessian; 242e0ed867bSAlp Dener bnk->computestep = TaoBQNKComputeStep; 243e0ed867bSAlp Dener bnk->init_type = BNK_INIT_DIRECTION; 244e0ed867bSAlp Dener 245e0ed867bSAlp Dener ierr = PetscNewLog(tao,&bqnk);CHKERRQ(ierr); 246e0ed867bSAlp Dener bnk->ctx = (void*)bqnk; 247f5766c09SAlp Dener bqnk->is_spd = PETSC_TRUE; 248e0ed867bSAlp Dener 249e0ed867bSAlp Dener ierr = MatCreate(PetscObjectComm((PetscObject)tao), &bqnk->B);CHKERRQ(ierr); 250e0ed867bSAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)bqnk->B, (PetscObject)tao, 1);CHKERRQ(ierr); 251e0ed867bSAlp Dener ierr = MatSetOptionsPrefix(bqnk->B, "tao_bqnk_");CHKERRQ(ierr); 252e0ed867bSAlp Dener ierr = MatSetType(bqnk->B, MATLMVMSR1);CHKERRQ(ierr); 253e0ed867bSAlp Dener PetscFunctionReturn(0); 254e0ed867bSAlp Dener } 255f5766c09SAlp Dener 256414d97d3SAlp Dener /*@ 257414d97d3SAlp Dener TaoGetLMVMMatrix - Returns a pointer to the internal LMVM matrix. Valid 258414d97d3SAlp Dener only for quasi-Newton family of methods. 259414d97d3SAlp Dener 260414d97d3SAlp Dener Input Parameters: 261414d97d3SAlp Dener . tao - Tao solver context 262414d97d3SAlp Dener 263414d97d3SAlp Dener Output Parameters: 264414d97d3SAlp Dener . B - LMVM matrix 265414d97d3SAlp Dener 266414d97d3SAlp Dener Level: advanced 267414d97d3SAlp Dener 268414d97d3SAlp Dener .seealso: TAOBQNLS, TAOBQNKLS, TAOBQNKTL, TAOBQNKTR, MATLMVM, TaoSetLMVMMatrix() 269414d97d3SAlp Dener @*/ 270f5766c09SAlp Dener PetscErrorCode TaoGetLMVMMatrix(Tao tao, Mat *B) 271f5766c09SAlp Dener { 272f5766c09SAlp Dener TAO_BNK *bnk = (TAO_BNK*)tao->data; 273f5766c09SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 274f5766c09SAlp Dener PetscErrorCode ierr; 275414d97d3SAlp Dener PetscBool flg = PETSC_FALSE; 276f5766c09SAlp Dener 277f5766c09SAlp Dener PetscFunctionBegin; 278414d97d3SAlp Dener ierr = PetscObjectTypeCompareAny((PetscObject)tao, &flg, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, "");CHKERRQ(ierr); 279414d97d3SAlp Dener if (!flg) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM Matrix only exists for quasi-Newton algorithms"); 280f5766c09SAlp Dener *B = bqnk->B; 281f5766c09SAlp Dener PetscFunctionReturn(0); 282f5766c09SAlp Dener } 283414d97d3SAlp Dener 284414d97d3SAlp Dener /*@ 285414d97d3SAlp Dener TaoSetLMVMMatrix - Sets an external LMVM matrix into the Tao solver. Valid 286414d97d3SAlp Dener only for quasi-Newton family of methods. 287414d97d3SAlp Dener 288414d97d3SAlp Dener QN family of methods create their own LMVM matrices and users who wish to 289414d97d3SAlp Dener manipulate this matrix should use TaoGetLMVMMatrix() instead. 290414d97d3SAlp Dener 291414d97d3SAlp Dener Input Parameters: 292414d97d3SAlp Dener + tao - Tao solver context 293414d97d3SAlp Dener - B - LMVM matrix 294414d97d3SAlp Dener 295414d97d3SAlp Dener Level: advanced 296414d97d3SAlp Dener 297414d97d3SAlp Dener .seealso: TAOBQNLS, TAOBQNKLS, TAOBQNKTL, TAOBQNKTR, MATLMVM, TaoGetLMVMMatrix() 298414d97d3SAlp Dener @*/ 299414d97d3SAlp Dener PetscErrorCode TaoSetLMVMMatrix(Tao tao, Mat B) 300414d97d3SAlp Dener { 301414d97d3SAlp Dener TAO_BNK *bnk = (TAO_BNK*)tao->data; 302414d97d3SAlp Dener TAO_BQNK *bqnk = (TAO_BQNK*)bnk->ctx; 303414d97d3SAlp Dener PetscErrorCode ierr; 304414d97d3SAlp Dener PetscBool flg = PETSC_FALSE; 305414d97d3SAlp Dener 306414d97d3SAlp Dener PetscFunctionBegin; 307414d97d3SAlp Dener ierr = PetscObjectTypeCompareAny((PetscObject)tao, &flg, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, "");CHKERRQ(ierr); 308414d97d3SAlp Dener if (!flg) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM Matrix only exists for quasi-Newton algorithms"); 309414d97d3SAlp Dener ierr = PetscObjectBaseTypeCompare((PetscObject)B, MATLMVM, &flg);CHKERRQ(ierr); 310414d97d3SAlp Dener if (!flg) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given matrix is not an LMVM matrix"); 311414d97d3SAlp Dener if (bqnk->B) { 312414d97d3SAlp Dener ierr = MatDestroy(&bqnk->B);CHKERRQ(ierr); 313414d97d3SAlp Dener } 314414d97d3SAlp Dener ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 315414d97d3SAlp Dener bqnk->B = B; 316414d97d3SAlp Dener PetscFunctionReturn(0); 317414d97d3SAlp Dener } 318