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 5f5766c09SAlp Dener /*------------------------------------------------------------*/ 69371c9d4SSatish Balay static PetscErrorCode TaoSolve_BLMVM(Tao tao) { 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)); 233c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf 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); 29f5766c09SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 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); 437494f0b1SStefano Zampini PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient)); 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 } 114f5766c09SAlp Dener PetscFunctionReturn(0); 115f5766c09SAlp Dener } 116f5766c09SAlp Dener 1179371c9d4SSatish Balay static PetscErrorCode TaoSetup_BLMVM(Tao tao) { 118f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 119f5766c09SAlp Dener 120f5766c09SAlp Dener PetscFunctionBegin; 121f5766c09SAlp Dener /* Existence of tao->solution checked in TaoSetup() */ 1229566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &blmP->Xold)); 1239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &blmP->Gold)); 1249566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &blmP->unprojected_gradient)); 125*48a46eb9SPierre Jolivet if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 126*48a46eb9SPierre Jolivet if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 127f5766c09SAlp Dener /* Allocate matrix for the limited memory approximation */ 1289566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(blmP->M, tao->solution, blmP->unprojected_gradient)); 129f5766c09SAlp Dener 130f5766c09SAlp Dener /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 1311baa6e33SBarry Smith if (blmP->H0) PetscCall(MatLMVMSetJ0(blmP->M, blmP->H0)); 132f5766c09SAlp Dener PetscFunctionReturn(0); 133f5766c09SAlp Dener } 134f5766c09SAlp Dener 135f5766c09SAlp Dener /* ---------------------------------------------------------- */ 1369371c9d4SSatish Balay static PetscErrorCode TaoDestroy_BLMVM(Tao tao) { 137f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 138f5766c09SAlp Dener 139f5766c09SAlp Dener PetscFunctionBegin; 140f5766c09SAlp Dener if (tao->setupcalled) { 1419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->unprojected_gradient)); 1429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->Xold)); 1439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->Gold)); 144f5766c09SAlp Dener } 1459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&blmP->M)); 1469371c9d4SSatish Balay if (blmP->H0) { PetscObjectDereference((PetscObject)blmP->H0); } 1479566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 148f5766c09SAlp Dener PetscFunctionReturn(0); 149f5766c09SAlp Dener } 150f5766c09SAlp Dener 151f5766c09SAlp Dener /*------------------------------------------------------------*/ 1529371c9d4SSatish Balay static PetscErrorCode TaoSetFromOptions_BLMVM(Tao tao, PetscOptionItems *PetscOptionsObject) { 153f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 154b94d7dedSBarry Smith PetscBool is_spd, is_set; 155a7e14dcfSSatish Balay 156a7e14dcfSSatish Balay PetscFunctionBegin; 157d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Limited-memory variable-metric method for bound constrained optimization"); 1589566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_blmvm_recycle", "enable recycling of the BFGS matrix between subsequent TaoSolve() calls", "", blmP->recycle, &blmP->recycle, NULL)); 159d0609cedSBarry Smith PetscOptionsHeadEnd(); 1609566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(blmP->M, ((PetscObject)tao)->prefix)); 1619566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(blmP->M, "tao_blmvm_")); 1629566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(blmP->M)); 163b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(blmP->M, &is_set, &is_spd)); 164b94d7dedSBarry Smith PetscCheck(is_set && is_spd, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite"); 165a7e14dcfSSatish Balay PetscFunctionReturn(0); 166a7e14dcfSSatish Balay } 167a7e14dcfSSatish Balay 168f5766c09SAlp Dener /*------------------------------------------------------------*/ 1699371c9d4SSatish Balay static PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer) { 170f5766c09SAlp Dener TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data; 171f5766c09SAlp Dener PetscBool isascii; 172a7e14dcfSSatish Balay 173a7e14dcfSSatish Balay PetscFunctionBegin; 1749566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 175f5766c09SAlp Dener if (isascii) { 17663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", lmP->grad)); 1779566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 1789566063dSJacob Faibussowitsch PetscCall(MatView(lmP->M, viewer)); 1799566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 180f5766c09SAlp Dener } 181f5766c09SAlp Dener PetscFunctionReturn(0); 182f5766c09SAlp Dener } 183f5766c09SAlp Dener 1849371c9d4SSatish Balay static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU) { 185f5766c09SAlp Dener TAO_BLMVM *blm = (TAO_BLMVM *)tao->data; 186f5766c09SAlp Dener 187f5766c09SAlp Dener PetscFunctionBegin; 188f5766c09SAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 189f5766c09SAlp Dener PetscValidHeaderSpecific(DXL, VEC_CLASSID, 2); 190f5766c09SAlp Dener PetscValidHeaderSpecific(DXU, VEC_CLASSID, 3); 1913c859ba3SBarry Smith PetscCheck(tao->gradient && blm->unprojected_gradient, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Dual variables don't exist yet or no longer exist."); 192f5766c09SAlp Dener 1939566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, DXL)); 1949566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXL, -1.0, blm->unprojected_gradient)); 1959566063dSJacob Faibussowitsch PetscCall(VecSet(DXU, 0.0)); 1969566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(DXL, DXL, DXU)); 197f5766c09SAlp Dener 1989566063dSJacob Faibussowitsch PetscCall(VecCopy(blm->unprojected_gradient, DXU)); 1999566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU, -1.0, tao->gradient)); 2009566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU, 1.0, DXL)); 201f5766c09SAlp Dener PetscFunctionReturn(0); 202f5766c09SAlp Dener } 203f5766c09SAlp Dener 204f5766c09SAlp Dener /* ---------------------------------------------------------- */ 205f5766c09SAlp Dener /*MC 206f5766c09SAlp Dener TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method 207f5766c09SAlp Dener for nonlinear minimization with bound constraints. It is an extension 208f5766c09SAlp Dener of TAOLMVM 209f5766c09SAlp Dener 210f5766c09SAlp Dener Options Database Keys: 211f5766c09SAlp Dener . -tao_lmm_recycle - enable recycling of LMVM information between subsequent TaoSolve calls 212f5766c09SAlp Dener 213f5766c09SAlp Dener Level: beginner 214f5766c09SAlp Dener M*/ 2159371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao) { 216f5766c09SAlp Dener TAO_BLMVM *blmP; 217f5766c09SAlp Dener const char *morethuente_type = TAOLINESEARCHMT; 218f5766c09SAlp Dener 219f5766c09SAlp Dener PetscFunctionBegin; 220f5766c09SAlp Dener tao->ops->setup = TaoSetup_BLMVM; 221f5766c09SAlp Dener tao->ops->solve = TaoSolve_BLMVM; 222f5766c09SAlp Dener tao->ops->view = TaoView_BLMVM; 223a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_BLMVM; 224f5766c09SAlp Dener tao->ops->destroy = TaoDestroy_BLMVM; 225f5766c09SAlp Dener tao->ops->computedual = TaoComputeDual_BLMVM; 226f5766c09SAlp Dener 2279566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao, &blmP)); 228f5766c09SAlp Dener blmP->H0 = NULL; 229f5766c09SAlp Dener blmP->recycle = PETSC_FALSE; 230f5766c09SAlp Dener tao->data = (void *)blmP; 231f5766c09SAlp Dener 232f5766c09SAlp Dener /* Override default settings (unless already changed) */ 233f5766c09SAlp Dener if (!tao->max_it_changed) tao->max_it = 2000; 234f5766c09SAlp Dener if (!tao->max_funcs_changed) tao->max_funcs = 4000; 235f5766c09SAlp Dener 2369566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 2379566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 2389566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type)); 2399566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao)); 240f5766c09SAlp Dener 2419566063dSJacob Faibussowitsch PetscCall(KSPInitializePackage()); 2429566063dSJacob Faibussowitsch PetscCall(MatCreate(((PetscObject)tao)->comm, &blmP->M)); 2439566063dSJacob Faibussowitsch PetscCall(MatSetType(blmP->M, MATLMVMBFGS)); 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1)); 245f5766c09SAlp Dener PetscFunctionReturn(0); 246f5766c09SAlp Dener } 247f5766c09SAlp Dener 2481bb2a437SAlp Dener /*@ 2491bb2a437SAlp Dener TaoLMVMRecycle - Enable/disable recycling of the QN history between subsequent TaoSolve calls. 2501bb2a437SAlp Dener 2511bb2a437SAlp Dener Input Parameters: 2521bb2a437SAlp Dener + tao - the Tao solver context 2531bb2a437SAlp Dener - flg - Boolean flag for recycling (PETSC_TRUE or PETSC_FALSE) 2541bb2a437SAlp Dener 2551bb2a437SAlp Dener Level: intermediate 2561bb2a437SAlp Dener @*/ 2579371c9d4SSatish Balay PetscErrorCode TaoLMVMRecycle(Tao tao, PetscBool flg) { 258b39c12a9SAlp Dener TAO_LMVM *lmP; 259b39c12a9SAlp Dener TAO_BLMVM *blmP; 260b39c12a9SAlp Dener PetscBool is_lmvm, is_blmvm; 261b39c12a9SAlp Dener 262b39c12a9SAlp Dener PetscFunctionBegin; 2639566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOLMVM, &is_lmvm)); 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBLMVM, &is_blmvm)); 265b39c12a9SAlp Dener if (is_lmvm) { 266b39c12a9SAlp Dener lmP = (TAO_LMVM *)tao->data; 267b39c12a9SAlp Dener lmP->recycle = flg; 268b39c12a9SAlp Dener } else if (is_blmvm) { 269b39c12a9SAlp Dener blmP = (TAO_BLMVM *)tao->data; 270b39c12a9SAlp Dener blmP->recycle = flg; 2718854b543SStefano Zampini } 272b39c12a9SAlp Dener PetscFunctionReturn(0); 273b39c12a9SAlp Dener } 274b39c12a9SAlp Dener 2751bb2a437SAlp Dener /*@ 2761bb2a437SAlp Dener TaoLMVMSetH0 - Set the initial Hessian for the QN approximation 2771bb2a437SAlp Dener 2781bb2a437SAlp Dener Input Parameters: 2791bb2a437SAlp Dener + tao - the Tao solver context 2801bb2a437SAlp Dener - H0 - Mat object for the initial Hessian 2811bb2a437SAlp Dener 2821bb2a437SAlp Dener Level: advanced 2831bb2a437SAlp Dener 284db781477SPatrick Sanan .seealso: `TaoLMVMGetH0()`, `TaoLMVMGetH0KSP()` 2851bb2a437SAlp Dener @*/ 2869371c9d4SSatish Balay PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0) { 287f5766c09SAlp Dener TAO_LMVM *lmP; 288f5766c09SAlp Dener TAO_BLMVM *blmP; 289f5766c09SAlp Dener PetscBool is_lmvm, is_blmvm; 290f5766c09SAlp Dener 291b39c12a9SAlp Dener PetscFunctionBegin; 2929566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOLMVM, &is_lmvm)); 2939566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBLMVM, &is_blmvm)); 294f5766c09SAlp Dener if (is_lmvm) { 295f5766c09SAlp Dener lmP = (TAO_LMVM *)tao->data; 2969566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H0)); 297f5766c09SAlp Dener lmP->H0 = H0; 298f5766c09SAlp Dener } else if (is_blmvm) { 299f5766c09SAlp Dener blmP = (TAO_BLMVM *)tao->data; 3009566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H0)); 301f5766c09SAlp Dener blmP->H0 = H0; 3028854b543SStefano Zampini } 303f5766c09SAlp Dener PetscFunctionReturn(0); 304f5766c09SAlp Dener } 305f5766c09SAlp Dener 3061bb2a437SAlp Dener /*@ 3071bb2a437SAlp Dener TaoLMVMGetH0 - Get the matrix object for the QN initial Hessian 3081bb2a437SAlp Dener 3091bb2a437SAlp Dener Input Parameters: 3101bb2a437SAlp Dener . tao - the Tao solver context 3111bb2a437SAlp Dener 3121bb2a437SAlp Dener Output Parameters: 3131bb2a437SAlp Dener . H0 - Mat object for the initial Hessian 3141bb2a437SAlp Dener 3151bb2a437SAlp Dener Level: advanced 3161bb2a437SAlp Dener 317db781477SPatrick Sanan .seealso: `TaoLMVMSetH0()`, `TaoLMVMGetH0KSP()` 3181bb2a437SAlp Dener @*/ 3199371c9d4SSatish Balay PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0) { 320f5766c09SAlp Dener TAO_LMVM *lmP; 321f5766c09SAlp Dener TAO_BLMVM *blmP; 322f5766c09SAlp Dener PetscBool is_lmvm, is_blmvm; 323f5766c09SAlp Dener Mat M; 324f5766c09SAlp Dener 325b39c12a9SAlp Dener PetscFunctionBegin; 3269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOLMVM, &is_lmvm)); 3279566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBLMVM, &is_blmvm)); 328f5766c09SAlp Dener if (is_lmvm) { 329f5766c09SAlp Dener lmP = (TAO_LMVM *)tao->data; 330f5766c09SAlp Dener M = lmP->M; 331f5766c09SAlp Dener } else if (is_blmvm) { 332f5766c09SAlp Dener blmP = (TAO_BLMVM *)tao->data; 333f5766c09SAlp Dener M = blmP->M; 3348854b543SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM."); 3359566063dSJacob Faibussowitsch PetscCall(MatLMVMGetJ0(M, H0)); 336f5766c09SAlp Dener PetscFunctionReturn(0); 337f5766c09SAlp Dener } 338f5766c09SAlp Dener 3391bb2a437SAlp Dener /*@ 3401bb2a437SAlp Dener TaoLMVMGetH0KSP - Get the iterative solver for applying the inverse of the QN initial Hessian 3411bb2a437SAlp Dener 3421bb2a437SAlp Dener Input Parameters: 3431bb2a437SAlp Dener . tao - the Tao solver context 3441bb2a437SAlp Dener 3451bb2a437SAlp Dener Output Parameters: 3461bb2a437SAlp Dener . ksp - KSP solver context for the initial Hessian 3471bb2a437SAlp Dener 3481bb2a437SAlp Dener Level: advanced 3491bb2a437SAlp Dener 350db781477SPatrick Sanan .seealso: `TaoLMVMGetH0()`, `TaoLMVMGetH0KSP()` 3511bb2a437SAlp Dener @*/ 3529371c9d4SSatish Balay PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp) { 353f5766c09SAlp Dener TAO_LMVM *lmP; 354f5766c09SAlp Dener TAO_BLMVM *blmP; 355f5766c09SAlp Dener PetscBool is_lmvm, is_blmvm; 356f5766c09SAlp Dener Mat M; 357f5766c09SAlp Dener 3588854b543SStefano Zampini PetscFunctionBegin; 3599566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOLMVM, &is_lmvm)); 3609566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOBLMVM, &is_blmvm)); 361f5766c09SAlp Dener if (is_lmvm) { 362f5766c09SAlp Dener lmP = (TAO_LMVM *)tao->data; 363f5766c09SAlp Dener M = lmP->M; 364f5766c09SAlp Dener } else if (is_blmvm) { 365f5766c09SAlp Dener blmP = (TAO_BLMVM *)tao->data; 366f5766c09SAlp Dener M = blmP->M; 3678854b543SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM."); 3689566063dSJacob Faibussowitsch PetscCall(MatLMVMGetJ0KSP(M, ksp)); 369a9603a14SPatrick Farrell PetscFunctionReturn(0); 370a9603a14SPatrick Farrell } 371