11bb2a437SAlp Dener #include <petsctaolinesearch.h> /*I "petsctaolinesearch.h" I*/ 2f5766c09SAlp Dener #include <../src/tao/unconstrained/impls/lmvm/lmvm.h> 3f5766c09SAlp Dener #include <../src/tao/bound/impls/blmvm/blmvm.h> 4a7e14dcfSSatish Balay 5d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_BLMVM(Tao tao) 6d71ae5a4SJacob Faibussowitsch { 7f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 8f5766c09SAlp Dener TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 9f5766c09SAlp Dener PetscReal f, fold, gdx, gnorm, gnorm2; 10f5766c09SAlp Dener PetscReal stepsize = 1.0, delta; 11a7e14dcfSSatish Balay 12f5766c09SAlp Dener PetscFunctionBegin; 13f5766c09SAlp Dener /* Project initial point onto bounds */ 149566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 159566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution)); 169566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU)); 17f5766c09SAlp Dener 18f5766c09SAlp Dener /* Check convergence criteria */ 199566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, blmP->unprojected_gradient)); 209566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient)); 21f5766c09SAlp Dener 229566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm)); 23*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN"); 24f5766c09SAlp Dener 25f5766c09SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 269566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its)); 279566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize)); 28dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 293ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS); 30f5766c09SAlp Dener 31f5766c09SAlp Dener /* Set counter for gradient/reset steps */ 32f5766c09SAlp Dener if (!blmP->recycle) { 33f5766c09SAlp Dener blmP->grad = 0; 34f5766c09SAlp Dener blmP->reset = 0; 359566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(blmP->M, PETSC_FALSE)); 36f5766c09SAlp Dener } 37f5766c09SAlp Dener 38f5766c09SAlp Dener /* Have not converged; continue with Newton method */ 39f5766c09SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 40e1e80dc8SAlp Dener /* Call general purpose update function */ 41e1e80dc8SAlp Dener if (tao->ops->update) { 42dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, update, tao->niter, tao->user_update); 43270bebe6SStefano Zampini PetscCall(TaoComputeObjective(tao, tao->solution, &f)); 44e1e80dc8SAlp Dener } 45f5766c09SAlp Dener /* Compute direction */ 46f5766c09SAlp Dener gnorm2 = gnorm * gnorm; 478cabe928SAlp Dener if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON; 488cabe928SAlp Dener if (f == 0.0) { 498cabe928SAlp Dener delta = 2.0 / gnorm2; 508cabe928SAlp Dener } else { 518cabe928SAlp Dener delta = 2.0 * PetscAbsScalar(f) / gnorm2; 528cabe928SAlp Dener } 539566063dSJacob Faibussowitsch PetscCall(MatLMVMSymBroydenSetDelta(blmP->M, delta)); 549566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(blmP->M, tao->solution, tao->gradient)); 559566063dSJacob Faibussowitsch PetscCall(MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection)); 569566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->stepdirection, tao->solution, tao->XL, tao->XU, tao->gradient)); 57f5766c09SAlp Dener 58f5766c09SAlp Dener /* Check for success (descent direction) */ 599566063dSJacob Faibussowitsch PetscCall(VecDot(blmP->unprojected_gradient, tao->gradient, &gdx)); 60f5766c09SAlp Dener if (gdx <= 0) { 61f5766c09SAlp Dener /* Step is not descent or solve was not successful 62f5766c09SAlp Dener Use steepest descent direction (scaled) */ 63f5766c09SAlp Dener ++blmP->grad; 64f5766c09SAlp Dener 659566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(blmP->M, PETSC_FALSE)); 669566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient)); 679566063dSJacob Faibussowitsch PetscCall(MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection)); 68f5766c09SAlp Dener } 699566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0)); 70f5766c09SAlp Dener 71f5766c09SAlp Dener /* Perform the linesearch */ 72f5766c09SAlp Dener fold = f; 739566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, blmP->Xold)); 749566063dSJacob Faibussowitsch PetscCall(VecCopy(blmP->unprojected_gradient, blmP->Gold)); 759566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0)); 769566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status)); 779566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 78f5766c09SAlp Dener 79f5766c09SAlp Dener if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 80f5766c09SAlp Dener /* Linesearch failed 81f5766c09SAlp Dener Reset factors and use scaled (projected) gradient step */ 82f5766c09SAlp Dener ++blmP->reset; 83f5766c09SAlp Dener 84f5766c09SAlp Dener f = fold; 859566063dSJacob Faibussowitsch PetscCall(VecCopy(blmP->Xold, tao->solution)); 869566063dSJacob Faibussowitsch PetscCall(VecCopy(blmP->Gold, blmP->unprojected_gradient)); 87f5766c09SAlp Dener 889566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(blmP->M, PETSC_FALSE)); 899566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient)); 909566063dSJacob Faibussowitsch PetscCall(MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection)); 919566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0)); 92f5766c09SAlp Dener 93f5766c09SAlp Dener /* This may be incorrect; linesearch has values for stepmax and stepmin 94f5766c09SAlp Dener that should be reset. */ 959566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0)); 969566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status)); 979566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 98f5766c09SAlp Dener 99f5766c09SAlp Dener if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 100f5766c09SAlp Dener tao->reason = TAO_DIVERGED_LS_FAILURE; 101f5766c09SAlp Dener break; 102f5766c09SAlp Dener } 103f5766c09SAlp Dener } 104f5766c09SAlp Dener 105f5766c09SAlp Dener /* Check for converged */ 1069566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient)); 1079566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm)); 1083c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Not-a-Number"); 109f5766c09SAlp Dener tao->niter++; 1109566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its)); 1119566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize)); 112dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 113f5766c09SAlp Dener } 1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115f5766c09SAlp Dener } 116f5766c09SAlp Dener 117d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetup_BLMVM(Tao tao) 118d71ae5a4SJacob Faibussowitsch { 119f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 120f5766c09SAlp Dener 121f5766c09SAlp Dener PetscFunctionBegin; 122f5766c09SAlp Dener /* Existence of tao->solution checked in TaoSetup() */ 1239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &blmP->Xold)); 1249566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &blmP->Gold)); 1259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &blmP->unprojected_gradient)); 12648a46eb9SPierre Jolivet if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 12748a46eb9SPierre Jolivet if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 128f5766c09SAlp Dener /* Allocate matrix for the limited memory approximation */ 1299566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(blmP->M, tao->solution, blmP->unprojected_gradient)); 130f5766c09SAlp Dener 131f5766c09SAlp Dener /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 1321baa6e33SBarry Smith if (blmP->H0) PetscCall(MatLMVMSetJ0(blmP->M, blmP->H0)); 1333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 134f5766c09SAlp Dener } 135f5766c09SAlp Dener 136d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_BLMVM(Tao tao) 137d71ae5a4SJacob Faibussowitsch { 138f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 139f5766c09SAlp Dener 140f5766c09SAlp Dener PetscFunctionBegin; 141f5766c09SAlp Dener if (tao->setupcalled) { 1429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->unprojected_gradient)); 1439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->Xold)); 1449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->Gold)); 145f5766c09SAlp Dener } 1469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&blmP->M)); 1473ba16761SJacob Faibussowitsch if (blmP->H0) PetscCall(PetscObjectDereference((PetscObject)blmP->H0)); 1489566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150f5766c09SAlp Dener } 151f5766c09SAlp Dener 152ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_BLMVM(Tao tao, PetscOptionItems PetscOptionsObject) 153d71ae5a4SJacob Faibussowitsch { 154f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 155b94d7dedSBarry Smith PetscBool is_spd, is_set; 156a7e14dcfSSatish Balay 157a7e14dcfSSatish Balay PetscFunctionBegin; 158d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Limited-memory variable-metric method for bound constrained optimization"); 1599566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_blmvm_recycle", "enable recycling of the BFGS matrix between subsequent TaoSolve() calls", "", blmP->recycle, &blmP->recycle, NULL)); 160d0609cedSBarry Smith PetscOptionsHeadEnd(); 1619566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(blmP->M, ((PetscObject)tao)->prefix)); 1629566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(blmP->M, "tao_blmvm_")); 1639566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(blmP->M)); 164b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(blmP->M, &is_set, &is_spd)); 165b94d7dedSBarry Smith PetscCheck(is_set && is_spd, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite"); 1663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 167a7e14dcfSSatish Balay } 168a7e14dcfSSatish Balay 169d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer) 170d71ae5a4SJacob Faibussowitsch { 171f5766c09SAlp Dener TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data; 172f5766c09SAlp Dener PetscBool isascii; 173a7e14dcfSSatish Balay 174a7e14dcfSSatish Balay PetscFunctionBegin; 1759566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 176f5766c09SAlp Dener if (isascii) { 17763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", lmP->grad)); 1789566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 1799566063dSJacob Faibussowitsch PetscCall(MatView(lmP->M, viewer)); 1809566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 181f5766c09SAlp Dener } 1823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 183f5766c09SAlp Dener } 184f5766c09SAlp Dener 185d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU) 186d71ae5a4SJacob Faibussowitsch { 187f5766c09SAlp Dener TAO_BLMVM *blm = (TAO_BLMVM *)tao->data; 188f5766c09SAlp Dener 189f5766c09SAlp Dener PetscFunctionBegin; 190f5766c09SAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 191f5766c09SAlp Dener PetscValidHeaderSpecific(DXL, VEC_CLASSID, 2); 192f5766c09SAlp Dener PetscValidHeaderSpecific(DXU, VEC_CLASSID, 3); 1933c859ba3SBarry Smith PetscCheck(tao->gradient && blm->unprojected_gradient, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Dual variables don't exist yet or no longer exist."); 194f5766c09SAlp Dener 1959566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, DXL)); 1969566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXL, -1.0, blm->unprojected_gradient)); 1979566063dSJacob Faibussowitsch PetscCall(VecSet(DXU, 0.0)); 1989566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(DXL, DXL, DXU)); 199f5766c09SAlp Dener 2009566063dSJacob Faibussowitsch PetscCall(VecCopy(blm->unprojected_gradient, DXU)); 2019566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU, -1.0, tao->gradient)); 2029566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU, 1.0, DXL)); 2033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 204f5766c09SAlp Dener } 205f5766c09SAlp Dener 206f5766c09SAlp Dener /*MC 207f5766c09SAlp Dener TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method 208f5766c09SAlp Dener for nonlinear minimization with bound constraints. It is an extension 20920f4b53cSBarry Smith of `TAOLMVM` 210f5766c09SAlp Dener 21120f4b53cSBarry Smith Options Database Key: 21220f4b53cSBarry Smith . -tao_lmm_recycle - enable recycling of LMVM information between subsequent `TaoSolve()` calls 213f5766c09SAlp Dener 214f5766c09SAlp Dener Level: beginner 21520f4b53cSBarry Smith 216fe8e7dddSPierre Jolivet .seealso: `Tao`, `TAOLMVM`, `TAOBLMVM`, `TaoLMVMGetH0()`, `TaoLMVMGetH0KSP()` 217f5766c09SAlp Dener M*/ 218d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao) 219d71ae5a4SJacob Faibussowitsch { 220f5766c09SAlp Dener TAO_BLMVM *blmP; 221f5766c09SAlp Dener const char *morethuente_type = TAOLINESEARCHMT; 222f5766c09SAlp Dener 223f5766c09SAlp Dener PetscFunctionBegin; 224f5766c09SAlp Dener tao->ops->setup = TaoSetup_BLMVM; 225f5766c09SAlp Dener tao->ops->solve = TaoSolve_BLMVM; 226f5766c09SAlp Dener tao->ops->view = TaoView_BLMVM; 227a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_BLMVM; 228f5766c09SAlp Dener tao->ops->destroy = TaoDestroy_BLMVM; 229f5766c09SAlp Dener tao->ops->computedual = TaoComputeDual_BLMVM; 230f5766c09SAlp Dener 2314dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&blmP)); 232f5766c09SAlp Dener blmP->H0 = NULL; 233f5766c09SAlp Dener blmP->recycle = PETSC_FALSE; 234f5766c09SAlp Dener tao->data = (void *)blmP; 235f5766c09SAlp Dener 236f5766c09SAlp Dener /* Override default settings (unless already changed) */ 237606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao)); 238606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_it, 2000); 239606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_funcs, 4000); 240f5766c09SAlp Dener 2419566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 2429566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 2439566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type)); 2449566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao)); 245f5766c09SAlp Dener 2469566063dSJacob Faibussowitsch PetscCall(KSPInitializePackage()); 2479566063dSJacob Faibussowitsch PetscCall(MatCreate(((PetscObject)tao)->comm, &blmP->M)); 2489566063dSJacob Faibussowitsch PetscCall(MatSetType(blmP->M, MATLMVMBFGS)); 2499566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1)); 2503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 251f5766c09SAlp Dener } 252f5766c09SAlp Dener 2531bb2a437SAlp Dener /*@ 25420f4b53cSBarry Smith TaoLMVMRecycle - Enable/disable recycling of the QN history between subsequent `TaoSolve()` calls. 2551bb2a437SAlp Dener 2561bb2a437SAlp Dener Input Parameters: 25720f4b53cSBarry Smith + tao - the `Tao` solver context 25820f4b53cSBarry Smith - flg - Boolean flag for recycling (`PETSC_TRUE` or `PETSC_FALSE`) 2591bb2a437SAlp Dener 2601bb2a437SAlp Dener Level: intermediate 26120f4b53cSBarry Smith 26276fbde31SPierre Jolivet .seealso: `Tao`, `TAOLMVM`, `TAOBLMVM` 2631bb2a437SAlp Dener @*/ 264d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLMVMRecycle(Tao tao, PetscBool flg) 265d71ae5a4SJacob Faibussowitsch { 266b39c12a9SAlp Dener TAO_LMVM *lmP; 267b39c12a9SAlp Dener TAO_BLMVM *blmP; 268b39c12a9SAlp Dener PetscBool is_lmvm, is_blmvm; 269b39c12a9SAlp Dener 270b39c12a9SAlp Dener PetscFunctionBegin; 2719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOLMVM, &is_lmvm)); 2729566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBLMVM, &is_blmvm)); 273b39c12a9SAlp Dener if (is_lmvm) { 274b39c12a9SAlp Dener lmP = (TAO_LMVM *)tao->data; 275b39c12a9SAlp Dener lmP->recycle = flg; 276b39c12a9SAlp Dener } else if (is_blmvm) { 277b39c12a9SAlp Dener blmP = (TAO_BLMVM *)tao->data; 278b39c12a9SAlp Dener blmP->recycle = flg; 2798854b543SStefano Zampini } 2803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 281b39c12a9SAlp Dener } 282b39c12a9SAlp Dener 2831bb2a437SAlp Dener /*@ 2841bb2a437SAlp Dener TaoLMVMSetH0 - Set the initial Hessian for the QN approximation 2851bb2a437SAlp Dener 2861bb2a437SAlp Dener Input Parameters: 28720f4b53cSBarry Smith + tao - the `Tao` solver context 28820f4b53cSBarry Smith - H0 - `Mat` object for the initial Hessian 2891bb2a437SAlp Dener 2901bb2a437SAlp Dener Level: advanced 2911bb2a437SAlp Dener 29220f4b53cSBarry Smith .seealso: `Tao`, `TAOLMVM`, `TAOBLMVM`, `TaoLMVMGetH0()`, `TaoLMVMGetH0KSP()` 2931bb2a437SAlp Dener @*/ 294d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0) 295d71ae5a4SJacob Faibussowitsch { 296f5766c09SAlp Dener TAO_LMVM *lmP; 297f5766c09SAlp Dener TAO_BLMVM *blmP; 298f5766c09SAlp Dener PetscBool is_lmvm, is_blmvm; 299f5766c09SAlp Dener 300b39c12a9SAlp Dener PetscFunctionBegin; 3019566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOLMVM, &is_lmvm)); 3029566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBLMVM, &is_blmvm)); 303f5766c09SAlp Dener if (is_lmvm) { 304f5766c09SAlp Dener lmP = (TAO_LMVM *)tao->data; 3059566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H0)); 306f5766c09SAlp Dener lmP->H0 = H0; 307f5766c09SAlp Dener } else if (is_blmvm) { 308f5766c09SAlp Dener blmP = (TAO_BLMVM *)tao->data; 3099566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H0)); 310f5766c09SAlp Dener blmP->H0 = H0; 3118854b543SStefano Zampini } 3123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 313f5766c09SAlp Dener } 314f5766c09SAlp Dener 3151bb2a437SAlp Dener /*@ 3161bb2a437SAlp Dener TaoLMVMGetH0 - Get the matrix object for the QN initial Hessian 3171bb2a437SAlp Dener 31820f4b53cSBarry Smith Input Parameter: 31920f4b53cSBarry Smith . tao - the `Tao` solver context 3201bb2a437SAlp Dener 32120f4b53cSBarry Smith Output Parameter: 32220f4b53cSBarry Smith . H0 - `Mat` object for the initial Hessian 3231bb2a437SAlp Dener 3241bb2a437SAlp Dener Level: advanced 3251bb2a437SAlp Dener 32620f4b53cSBarry Smith .seealso: `Tao`, `TAOLMVM`, `TAOBLMVM`, `TaoLMVMSetH0()`, `TaoLMVMGetH0KSP()` 3271bb2a437SAlp Dener @*/ 328d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0) 329d71ae5a4SJacob Faibussowitsch { 330f5766c09SAlp Dener TAO_LMVM *lmP; 331f5766c09SAlp Dener TAO_BLMVM *blmP; 332f5766c09SAlp Dener PetscBool is_lmvm, is_blmvm; 333f5766c09SAlp Dener Mat M; 334f5766c09SAlp Dener 335b39c12a9SAlp Dener PetscFunctionBegin; 3369566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOLMVM, &is_lmvm)); 3379566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBLMVM, &is_blmvm)); 338f5766c09SAlp Dener if (is_lmvm) { 339f5766c09SAlp Dener lmP = (TAO_LMVM *)tao->data; 340f5766c09SAlp Dener M = lmP->M; 341f5766c09SAlp Dener } else if (is_blmvm) { 342f5766c09SAlp Dener blmP = (TAO_BLMVM *)tao->data; 343f5766c09SAlp Dener M = blmP->M; 3448854b543SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM."); 3459566063dSJacob Faibussowitsch PetscCall(MatLMVMGetJ0(M, H0)); 3463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 347f5766c09SAlp Dener } 348f5766c09SAlp Dener 3491bb2a437SAlp Dener /*@ 3501bb2a437SAlp Dener TaoLMVMGetH0KSP - Get the iterative solver for applying the inverse of the QN initial Hessian 3511bb2a437SAlp Dener 3522fe279fdSBarry Smith Input Parameter: 35320f4b53cSBarry Smith . tao - the `Tao` solver context 3541bb2a437SAlp Dener 3552fe279fdSBarry Smith Output Parameter: 35620f4b53cSBarry Smith . ksp - `KSP` solver context for the initial Hessian 3571bb2a437SAlp Dener 3581bb2a437SAlp Dener Level: advanced 3591bb2a437SAlp Dener 360fe8e7dddSPierre Jolivet .seealso: `Tao`, `TAOLMVM`, `TAOBLMVM`, `TaoLMVMGetH0()` 3611bb2a437SAlp Dener @*/ 362d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp) 363d71ae5a4SJacob Faibussowitsch { 364f5766c09SAlp Dener TAO_LMVM *lmP; 365f5766c09SAlp Dener TAO_BLMVM *blmP; 366f5766c09SAlp Dener PetscBool is_lmvm, is_blmvm; 367f5766c09SAlp Dener Mat M; 368f5766c09SAlp Dener 3698854b543SStefano Zampini PetscFunctionBegin; 3709566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOLMVM, &is_lmvm)); 3719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBLMVM, &is_blmvm)); 372f5766c09SAlp Dener if (is_lmvm) { 373f5766c09SAlp Dener lmP = (TAO_LMVM *)tao->data; 374f5766c09SAlp Dener M = lmP->M; 375f5766c09SAlp Dener } else if (is_blmvm) { 376f5766c09SAlp Dener blmP = (TAO_BLMVM *)tao->data; 377f5766c09SAlp Dener M = blmP->M; 3788854b543SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM."); 3799566063dSJacob Faibussowitsch PetscCall(MatLMVMGetJ0KSP(M, ksp)); 3803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 381a9603a14SPatrick Farrell } 382