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 /*------------------------------------------------------------*/ 6d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_BLMVM(Tao tao) 7d71ae5a4SJacob Faibussowitsch { 8f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 9f5766c09SAlp Dener TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 10f5766c09SAlp Dener PetscReal f, fold, gdx, gnorm, gnorm2; 11f5766c09SAlp Dener PetscReal stepsize = 1.0, delta; 12a7e14dcfSSatish Balay 13f5766c09SAlp Dener PetscFunctionBegin; 14f5766c09SAlp Dener /* Project initial point onto bounds */ 159566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 169566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution)); 179566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU)); 18f5766c09SAlp Dener 19f5766c09SAlp Dener /* Check convergence criteria */ 209566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, blmP->unprojected_gradient)); 219566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient)); 22f5766c09SAlp Dener 239566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm)); 243c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 25f5766c09SAlp Dener 26f5766c09SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 279566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its)); 289566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize)); 29dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 30*3ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS); 31f5766c09SAlp Dener 32f5766c09SAlp Dener /* Set counter for gradient/reset steps */ 33f5766c09SAlp Dener if (!blmP->recycle) { 34f5766c09SAlp Dener blmP->grad = 0; 35f5766c09SAlp Dener blmP->reset = 0; 369566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(blmP->M, PETSC_FALSE)); 37f5766c09SAlp Dener } 38f5766c09SAlp Dener 39f5766c09SAlp Dener /* Have not converged; continue with Newton method */ 40f5766c09SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 41e1e80dc8SAlp Dener /* Call general purpose update function */ 42e1e80dc8SAlp Dener if (tao->ops->update) { 43dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, update, tao->niter, tao->user_update); 447494f0b1SStefano Zampini PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient)); 45e1e80dc8SAlp Dener } 46f5766c09SAlp Dener /* Compute direction */ 47f5766c09SAlp Dener gnorm2 = gnorm * gnorm; 488cabe928SAlp Dener if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON; 498cabe928SAlp Dener if (f == 0.0) { 508cabe928SAlp Dener delta = 2.0 / gnorm2; 518cabe928SAlp Dener } else { 528cabe928SAlp Dener delta = 2.0 * PetscAbsScalar(f) / gnorm2; 538cabe928SAlp Dener } 549566063dSJacob Faibussowitsch PetscCall(MatLMVMSymBroydenSetDelta(blmP->M, delta)); 559566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(blmP->M, tao->solution, tao->gradient)); 569566063dSJacob Faibussowitsch PetscCall(MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection)); 579566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->stepdirection, tao->solution, tao->XL, tao->XU, tao->gradient)); 58f5766c09SAlp Dener 59f5766c09SAlp Dener /* Check for success (descent direction) */ 609566063dSJacob Faibussowitsch PetscCall(VecDot(blmP->unprojected_gradient, tao->gradient, &gdx)); 61f5766c09SAlp Dener if (gdx <= 0) { 62f5766c09SAlp Dener /* Step is not descent or solve was not successful 63f5766c09SAlp Dener Use steepest descent direction (scaled) */ 64f5766c09SAlp Dener ++blmP->grad; 65f5766c09SAlp Dener 669566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(blmP->M, PETSC_FALSE)); 679566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient)); 689566063dSJacob Faibussowitsch PetscCall(MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection)); 69f5766c09SAlp Dener } 709566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0)); 71f5766c09SAlp Dener 72f5766c09SAlp Dener /* Perform the linesearch */ 73f5766c09SAlp Dener fold = f; 749566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, blmP->Xold)); 759566063dSJacob Faibussowitsch PetscCall(VecCopy(blmP->unprojected_gradient, blmP->Gold)); 769566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0)); 779566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status)); 789566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 79f5766c09SAlp Dener 80f5766c09SAlp Dener if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 81f5766c09SAlp Dener /* Linesearch failed 82f5766c09SAlp Dener Reset factors and use scaled (projected) gradient step */ 83f5766c09SAlp Dener ++blmP->reset; 84f5766c09SAlp Dener 85f5766c09SAlp Dener f = fold; 869566063dSJacob Faibussowitsch PetscCall(VecCopy(blmP->Xold, tao->solution)); 879566063dSJacob Faibussowitsch PetscCall(VecCopy(blmP->Gold, blmP->unprojected_gradient)); 88f5766c09SAlp Dener 899566063dSJacob Faibussowitsch PetscCall(MatLMVMReset(blmP->M, PETSC_FALSE)); 909566063dSJacob Faibussowitsch PetscCall(MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient)); 919566063dSJacob Faibussowitsch PetscCall(MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection)); 929566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0)); 93f5766c09SAlp Dener 94f5766c09SAlp Dener /* This may be incorrect; linesearch has values for stepmax and stepmin 95f5766c09SAlp Dener that should be reset. */ 969566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0)); 979566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status)); 989566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 99f5766c09SAlp Dener 100f5766c09SAlp Dener if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 101f5766c09SAlp Dener tao->reason = TAO_DIVERGED_LS_FAILURE; 102f5766c09SAlp Dener break; 103f5766c09SAlp Dener } 104f5766c09SAlp Dener } 105f5766c09SAlp Dener 106f5766c09SAlp Dener /* Check for converged */ 1079566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient)); 1089566063dSJacob Faibussowitsch PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm)); 1093c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Not-a-Number"); 110f5766c09SAlp Dener tao->niter++; 1119566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its)); 1129566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize)); 113dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 114f5766c09SAlp Dener } 115*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 116f5766c09SAlp Dener } 117f5766c09SAlp Dener 118d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetup_BLMVM(Tao tao) 119d71ae5a4SJacob Faibussowitsch { 120f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 121f5766c09SAlp Dener 122f5766c09SAlp Dener PetscFunctionBegin; 123f5766c09SAlp Dener /* Existence of tao->solution checked in TaoSetup() */ 1249566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &blmP->Xold)); 1259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &blmP->Gold)); 1269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &blmP->unprojected_gradient)); 12748a46eb9SPierre Jolivet if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 12848a46eb9SPierre Jolivet if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 129f5766c09SAlp Dener /* Allocate matrix for the limited memory approximation */ 1309566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(blmP->M, tao->solution, blmP->unprojected_gradient)); 131f5766c09SAlp Dener 132f5766c09SAlp Dener /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 1331baa6e33SBarry Smith if (blmP->H0) PetscCall(MatLMVMSetJ0(blmP->M, blmP->H0)); 134*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 135f5766c09SAlp Dener } 136f5766c09SAlp Dener 137f5766c09SAlp Dener /* ---------------------------------------------------------- */ 138d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_BLMVM(Tao tao) 139d71ae5a4SJacob Faibussowitsch { 140f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 141f5766c09SAlp Dener 142f5766c09SAlp Dener PetscFunctionBegin; 143f5766c09SAlp Dener if (tao->setupcalled) { 1449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->unprojected_gradient)); 1459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->Xold)); 1469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->Gold)); 147f5766c09SAlp Dener } 1489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&blmP->M)); 149*3ba16761SJacob Faibussowitsch if (blmP->H0) PetscCall(PetscObjectDereference((PetscObject)blmP->H0)); 1509566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 151*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 152f5766c09SAlp Dener } 153f5766c09SAlp Dener 154f5766c09SAlp Dener /*------------------------------------------------------------*/ 155d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetFromOptions_BLMVM(Tao tao, PetscOptionItems *PetscOptionsObject) 156d71ae5a4SJacob Faibussowitsch { 157f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 158b94d7dedSBarry Smith PetscBool is_spd, is_set; 159a7e14dcfSSatish Balay 160a7e14dcfSSatish Balay PetscFunctionBegin; 161d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Limited-memory variable-metric method for bound constrained optimization"); 1629566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_blmvm_recycle", "enable recycling of the BFGS matrix between subsequent TaoSolve() calls", "", blmP->recycle, &blmP->recycle, NULL)); 163d0609cedSBarry Smith PetscOptionsHeadEnd(); 1649566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(blmP->M, ((PetscObject)tao)->prefix)); 1659566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(blmP->M, "tao_blmvm_")); 1669566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(blmP->M)); 167b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(blmP->M, &is_set, &is_spd)); 168b94d7dedSBarry Smith PetscCheck(is_set && is_spd, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite"); 169*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 170a7e14dcfSSatish Balay } 171a7e14dcfSSatish Balay 172f5766c09SAlp Dener /*------------------------------------------------------------*/ 173d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer) 174d71ae5a4SJacob Faibussowitsch { 175f5766c09SAlp Dener TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data; 176f5766c09SAlp Dener PetscBool isascii; 177a7e14dcfSSatish Balay 178a7e14dcfSSatish Balay PetscFunctionBegin; 1799566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 180f5766c09SAlp Dener if (isascii) { 18163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", lmP->grad)); 1829566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 1839566063dSJacob Faibussowitsch PetscCall(MatView(lmP->M, viewer)); 1849566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 185f5766c09SAlp Dener } 186*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 187f5766c09SAlp Dener } 188f5766c09SAlp Dener 189d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU) 190d71ae5a4SJacob Faibussowitsch { 191f5766c09SAlp Dener TAO_BLMVM *blm = (TAO_BLMVM *)tao->data; 192f5766c09SAlp Dener 193f5766c09SAlp Dener PetscFunctionBegin; 194f5766c09SAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 195f5766c09SAlp Dener PetscValidHeaderSpecific(DXL, VEC_CLASSID, 2); 196f5766c09SAlp Dener PetscValidHeaderSpecific(DXU, VEC_CLASSID, 3); 1973c859ba3SBarry Smith PetscCheck(tao->gradient && blm->unprojected_gradient, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Dual variables don't exist yet or no longer exist."); 198f5766c09SAlp Dener 1999566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, DXL)); 2009566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXL, -1.0, blm->unprojected_gradient)); 2019566063dSJacob Faibussowitsch PetscCall(VecSet(DXU, 0.0)); 2029566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(DXL, DXL, DXU)); 203f5766c09SAlp Dener 2049566063dSJacob Faibussowitsch PetscCall(VecCopy(blm->unprojected_gradient, DXU)); 2059566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU, -1.0, tao->gradient)); 2069566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU, 1.0, DXL)); 207*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 208f5766c09SAlp Dener } 209f5766c09SAlp Dener 210f5766c09SAlp Dener /* ---------------------------------------------------------- */ 211f5766c09SAlp Dener /*MC 212f5766c09SAlp Dener TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method 213f5766c09SAlp Dener for nonlinear minimization with bound constraints. It is an extension 214f5766c09SAlp Dener of TAOLMVM 215f5766c09SAlp Dener 216f5766c09SAlp Dener Options Database Keys: 217f5766c09SAlp Dener . -tao_lmm_recycle - enable recycling of LMVM information between subsequent TaoSolve calls 218f5766c09SAlp Dener 219f5766c09SAlp Dener Level: beginner 220f5766c09SAlp Dener M*/ 221d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao) 222d71ae5a4SJacob Faibussowitsch { 223f5766c09SAlp Dener TAO_BLMVM *blmP; 224f5766c09SAlp Dener const char *morethuente_type = TAOLINESEARCHMT; 225f5766c09SAlp Dener 226f5766c09SAlp Dener PetscFunctionBegin; 227f5766c09SAlp Dener tao->ops->setup = TaoSetup_BLMVM; 228f5766c09SAlp Dener tao->ops->solve = TaoSolve_BLMVM; 229f5766c09SAlp Dener tao->ops->view = TaoView_BLMVM; 230a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_BLMVM; 231f5766c09SAlp Dener tao->ops->destroy = TaoDestroy_BLMVM; 232f5766c09SAlp Dener tao->ops->computedual = TaoComputeDual_BLMVM; 233f5766c09SAlp Dener 2344dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&blmP)); 235f5766c09SAlp Dener blmP->H0 = NULL; 236f5766c09SAlp Dener blmP->recycle = PETSC_FALSE; 237f5766c09SAlp Dener tao->data = (void *)blmP; 238f5766c09SAlp Dener 239f5766c09SAlp Dener /* Override default settings (unless already changed) */ 240f5766c09SAlp Dener if (!tao->max_it_changed) tao->max_it = 2000; 241f5766c09SAlp Dener if (!tao->max_funcs_changed) tao->max_funcs = 4000; 242f5766c09SAlp Dener 2439566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 2459566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type)); 2469566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao)); 247f5766c09SAlp Dener 2489566063dSJacob Faibussowitsch PetscCall(KSPInitializePackage()); 2499566063dSJacob Faibussowitsch PetscCall(MatCreate(((PetscObject)tao)->comm, &blmP->M)); 2509566063dSJacob Faibussowitsch PetscCall(MatSetType(blmP->M, MATLMVMBFGS)); 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1)); 252*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 253f5766c09SAlp Dener } 254f5766c09SAlp Dener 2551bb2a437SAlp Dener /*@ 2561bb2a437SAlp Dener TaoLMVMRecycle - Enable/disable recycling of the QN history between subsequent TaoSolve calls. 2571bb2a437SAlp Dener 2581bb2a437SAlp Dener Input Parameters: 2591bb2a437SAlp Dener + tao - the Tao solver context 2601bb2a437SAlp Dener - flg - Boolean flag for recycling (PETSC_TRUE or PETSC_FALSE) 2611bb2a437SAlp Dener 2621bb2a437SAlp Dener Level: intermediate 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 } 280*3ba16761SJacob 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: 2871bb2a437SAlp Dener + tao - the Tao solver context 2881bb2a437SAlp Dener - H0 - Mat object for the initial Hessian 2891bb2a437SAlp Dener 2901bb2a437SAlp Dener Level: advanced 2911bb2a437SAlp Dener 292db781477SPatrick Sanan .seealso: `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 } 312*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 313f5766c09SAlp Dener } 314f5766c09SAlp Dener 3151bb2a437SAlp Dener /*@ 3161bb2a437SAlp Dener TaoLMVMGetH0 - Get the matrix object for the QN initial Hessian 3171bb2a437SAlp Dener 3181bb2a437SAlp Dener Input Parameters: 3191bb2a437SAlp Dener . tao - the Tao solver context 3201bb2a437SAlp Dener 3211bb2a437SAlp Dener Output Parameters: 3221bb2a437SAlp Dener . H0 - Mat object for the initial Hessian 3231bb2a437SAlp Dener 3241bb2a437SAlp Dener Level: advanced 3251bb2a437SAlp Dener 326db781477SPatrick Sanan .seealso: `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)); 346*3ba16761SJacob 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 3521bb2a437SAlp Dener Input Parameters: 3531bb2a437SAlp Dener . tao - the Tao solver context 3541bb2a437SAlp Dener 3551bb2a437SAlp Dener Output Parameters: 3561bb2a437SAlp Dener . ksp - KSP solver context for the initial Hessian 3571bb2a437SAlp Dener 3581bb2a437SAlp Dener Level: advanced 3591bb2a437SAlp Dener 360db781477SPatrick Sanan .seealso: `TaoLMVMGetH0()`, `TaoLMVMGetH0KSP()` 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)); 380*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 381a9603a14SPatrick Farrell } 382