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 /*------------------------------------------------------------*/ 6f5766c09SAlp Dener static PetscErrorCode TaoSolve_BLMVM(Tao tao) 7f5766c09SAlp Dener { 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)); 299566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 30f5766c09SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 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) { 439566063dSJacob Faibussowitsch PetscCall((*tao->ops->update)(tao, 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)); 1139566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 114f5766c09SAlp Dener } 115f5766c09SAlp Dener PetscFunctionReturn(0); 116f5766c09SAlp Dener } 117f5766c09SAlp Dener 118f5766c09SAlp Dener static PetscErrorCode TaoSetup_BLMVM(Tao tao) 119f5766c09SAlp Dener { 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)); 127f5766c09SAlp Dener if (!tao->stepdirection) { 1289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 129f5766c09SAlp Dener } 130f5766c09SAlp Dener if (!tao->gradient) { 1319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 132f5766c09SAlp Dener } 133f5766c09SAlp Dener /* Allocate matrix for the limited memory approximation */ 1349566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(blmP->M,tao->solution,blmP->unprojected_gradient)); 135f5766c09SAlp Dener 136f5766c09SAlp Dener /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 1371baa6e33SBarry Smith if (blmP->H0) PetscCall(MatLMVMSetJ0(blmP->M, blmP->H0)); 138f5766c09SAlp Dener PetscFunctionReturn(0); 139f5766c09SAlp Dener } 140f5766c09SAlp Dener 141f5766c09SAlp Dener /* ---------------------------------------------------------- */ 142f5766c09SAlp Dener static PetscErrorCode TaoDestroy_BLMVM(Tao tao) 143f5766c09SAlp Dener { 144f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 145f5766c09SAlp Dener 146f5766c09SAlp Dener PetscFunctionBegin; 147f5766c09SAlp Dener if (tao->setupcalled) { 1489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->unprojected_gradient)); 1499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->Xold)); 1509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->Gold)); 151f5766c09SAlp Dener } 1529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&blmP->M)); 153f5766c09SAlp Dener if (blmP->H0) { 154f5766c09SAlp Dener PetscObjectDereference((PetscObject)blmP->H0); 155f5766c09SAlp Dener } 1569566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 157f5766c09SAlp Dener PetscFunctionReturn(0); 158f5766c09SAlp Dener } 159f5766c09SAlp Dener 160f5766c09SAlp Dener /*------------------------------------------------------------*/ 1612ec5c1acSAlp Dener static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao) 162a7e14dcfSSatish Balay { 163f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 164*b94d7dedSBarry Smith PetscBool is_spd,is_set; 165a7e14dcfSSatish Balay 166a7e14dcfSSatish Balay PetscFunctionBegin; 167d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization"); 1689566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL)); 169d0609cedSBarry Smith PetscOptionsHeadEnd(); 1709566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(blmP->M, ((PetscObject)tao)->prefix)); 1719566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(blmP->M, "tao_blmvm_")); 1729566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(blmP->M)); 173*b94d7dedSBarry Smith PetscCall(MatIsSPDKnown(blmP->M, &is_set, &is_spd)); 174*b94d7dedSBarry Smith PetscCheck(is_set && is_spd,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite"); 175a7e14dcfSSatish Balay PetscFunctionReturn(0); 176a7e14dcfSSatish Balay } 177a7e14dcfSSatish Balay 178f5766c09SAlp Dener /*------------------------------------------------------------*/ 1795bd1e576SStefano Zampini static PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer) 180a7e14dcfSSatish Balay { 181f5766c09SAlp Dener TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data; 182f5766c09SAlp Dener PetscBool isascii; 183a7e14dcfSSatish Balay 184a7e14dcfSSatish Balay PetscFunctionBegin; 1859566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 186f5766c09SAlp Dener if (isascii) { 18763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", lmP->grad)); 1889566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 1899566063dSJacob Faibussowitsch PetscCall(MatView(lmP->M, viewer)); 1909566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 191f5766c09SAlp Dener } 192f5766c09SAlp Dener PetscFunctionReturn(0); 193f5766c09SAlp Dener } 194f5766c09SAlp Dener 195f5766c09SAlp Dener static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU) 196f5766c09SAlp Dener { 197f5766c09SAlp Dener TAO_BLMVM *blm = (TAO_BLMVM *) tao->data; 198f5766c09SAlp Dener 199f5766c09SAlp Dener PetscFunctionBegin; 200f5766c09SAlp Dener PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 201f5766c09SAlp Dener PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 202f5766c09SAlp Dener PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 2033c859ba3SBarry Smith PetscCheck(tao->gradient && blm->unprojected_gradient,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist."); 204f5766c09SAlp Dener 2059566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient,DXL)); 2069566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXL,-1.0,blm->unprojected_gradient)); 2079566063dSJacob Faibussowitsch PetscCall(VecSet(DXU,0.0)); 2089566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(DXL,DXL,DXU)); 209f5766c09SAlp Dener 2109566063dSJacob Faibussowitsch PetscCall(VecCopy(blm->unprojected_gradient,DXU)); 2119566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU,-1.0,tao->gradient)); 2129566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU,1.0,DXL)); 213f5766c09SAlp Dener PetscFunctionReturn(0); 214f5766c09SAlp Dener } 215f5766c09SAlp Dener 216f5766c09SAlp Dener /* ---------------------------------------------------------- */ 217f5766c09SAlp Dener /*MC 218f5766c09SAlp Dener TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method 219f5766c09SAlp Dener for nonlinear minimization with bound constraints. It is an extension 220f5766c09SAlp Dener of TAOLMVM 221f5766c09SAlp Dener 222f5766c09SAlp Dener Options Database Keys: 223f5766c09SAlp Dener . -tao_lmm_recycle - enable recycling of LMVM information between subsequent TaoSolve calls 224f5766c09SAlp Dener 225f5766c09SAlp Dener Level: beginner 226f5766c09SAlp Dener M*/ 227f5766c09SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao) 228f5766c09SAlp Dener { 229f5766c09SAlp Dener TAO_BLMVM *blmP; 230f5766c09SAlp Dener const char *morethuente_type = TAOLINESEARCHMT; 231f5766c09SAlp Dener 232f5766c09SAlp Dener PetscFunctionBegin; 233f5766c09SAlp Dener tao->ops->setup = TaoSetup_BLMVM; 234f5766c09SAlp Dener tao->ops->solve = TaoSolve_BLMVM; 235f5766c09SAlp Dener tao->ops->view = TaoView_BLMVM; 236a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_BLMVM; 237f5766c09SAlp Dener tao->ops->destroy = TaoDestroy_BLMVM; 238f5766c09SAlp Dener tao->ops->computedual = TaoComputeDual_BLMVM; 239f5766c09SAlp Dener 2409566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&blmP)); 241f5766c09SAlp Dener blmP->H0 = NULL; 242f5766c09SAlp Dener blmP->recycle = PETSC_FALSE; 243f5766c09SAlp Dener tao->data = (void*)blmP; 244f5766c09SAlp Dener 245f5766c09SAlp Dener /* Override default settings (unless already changed) */ 246f5766c09SAlp Dener if (!tao->max_it_changed) tao->max_it = 2000; 247f5766c09SAlp Dener if (!tao->max_funcs_changed) tao->max_funcs = 4000; 248f5766c09SAlp Dener 2499566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 2509566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 2519566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type)); 2529566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch,tao)); 253f5766c09SAlp Dener 2549566063dSJacob Faibussowitsch PetscCall(KSPInitializePackage()); 2559566063dSJacob Faibussowitsch PetscCall(MatCreate(((PetscObject)tao)->comm, &blmP->M)); 2569566063dSJacob Faibussowitsch PetscCall(MatSetType(blmP->M, MATLMVMBFGS)); 2579566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1)); 258f5766c09SAlp Dener PetscFunctionReturn(0); 259f5766c09SAlp Dener } 260f5766c09SAlp Dener 2611bb2a437SAlp Dener /*@ 2621bb2a437SAlp Dener TaoLMVMRecycle - Enable/disable recycling of the QN history between subsequent TaoSolve calls. 2631bb2a437SAlp Dener 2641bb2a437SAlp Dener Input Parameters: 2651bb2a437SAlp Dener + tao - the Tao solver context 2661bb2a437SAlp Dener - flg - Boolean flag for recycling (PETSC_TRUE or PETSC_FALSE) 2671bb2a437SAlp Dener 2681bb2a437SAlp Dener Level: intermediate 2691bb2a437SAlp Dener @*/ 270b39c12a9SAlp Dener PetscErrorCode TaoLMVMRecycle(Tao tao, PetscBool flg) 271b39c12a9SAlp Dener { 272b39c12a9SAlp Dener TAO_LMVM *lmP; 273b39c12a9SAlp Dener TAO_BLMVM *blmP; 274b39c12a9SAlp Dener PetscBool is_lmvm, is_blmvm; 275b39c12a9SAlp Dener 276b39c12a9SAlp Dener PetscFunctionBegin; 2779566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm)); 2789566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm)); 279b39c12a9SAlp Dener if (is_lmvm) { 280b39c12a9SAlp Dener lmP = (TAO_LMVM *)tao->data; 281b39c12a9SAlp Dener lmP->recycle = flg; 282b39c12a9SAlp Dener } else if (is_blmvm) { 283b39c12a9SAlp Dener blmP = (TAO_BLMVM *)tao->data; 284b39c12a9SAlp Dener blmP->recycle = flg; 2858854b543SStefano Zampini } 286b39c12a9SAlp Dener PetscFunctionReturn(0); 287b39c12a9SAlp Dener } 288b39c12a9SAlp Dener 2891bb2a437SAlp Dener /*@ 2901bb2a437SAlp Dener TaoLMVMSetH0 - Set the initial Hessian for the QN approximation 2911bb2a437SAlp Dener 2921bb2a437SAlp Dener Input Parameters: 2931bb2a437SAlp Dener + tao - the Tao solver context 2941bb2a437SAlp Dener - H0 - Mat object for the initial Hessian 2951bb2a437SAlp Dener 2961bb2a437SAlp Dener Level: advanced 2971bb2a437SAlp Dener 298db781477SPatrick Sanan .seealso: `TaoLMVMGetH0()`, `TaoLMVMGetH0KSP()` 2991bb2a437SAlp Dener @*/ 300f5766c09SAlp Dener PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0) 301f5766c09SAlp Dener { 302f5766c09SAlp Dener TAO_LMVM *lmP; 303f5766c09SAlp Dener TAO_BLMVM *blmP; 304f5766c09SAlp Dener PetscBool is_lmvm, is_blmvm; 305f5766c09SAlp Dener 306b39c12a9SAlp Dener PetscFunctionBegin; 3079566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm)); 3089566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm)); 309f5766c09SAlp Dener if (is_lmvm) { 310f5766c09SAlp Dener lmP = (TAO_LMVM *)tao->data; 3119566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H0)); 312f5766c09SAlp Dener lmP->H0 = H0; 313f5766c09SAlp Dener } else if (is_blmvm) { 314f5766c09SAlp Dener blmP = (TAO_BLMVM *)tao->data; 3159566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H0)); 316f5766c09SAlp Dener blmP->H0 = H0; 3178854b543SStefano Zampini } 318f5766c09SAlp Dener PetscFunctionReturn(0); 319f5766c09SAlp Dener } 320f5766c09SAlp Dener 3211bb2a437SAlp Dener /*@ 3221bb2a437SAlp Dener TaoLMVMGetH0 - Get the matrix object for the QN initial Hessian 3231bb2a437SAlp Dener 3241bb2a437SAlp Dener Input Parameters: 3251bb2a437SAlp Dener . tao - the Tao solver context 3261bb2a437SAlp Dener 3271bb2a437SAlp Dener Output Parameters: 3281bb2a437SAlp Dener . H0 - Mat object for the initial Hessian 3291bb2a437SAlp Dener 3301bb2a437SAlp Dener Level: advanced 3311bb2a437SAlp Dener 332db781477SPatrick Sanan .seealso: `TaoLMVMSetH0()`, `TaoLMVMGetH0KSP()` 3331bb2a437SAlp Dener @*/ 334f5766c09SAlp Dener PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0) 335f5766c09SAlp Dener { 336f5766c09SAlp Dener TAO_LMVM *lmP; 337f5766c09SAlp Dener TAO_BLMVM *blmP; 338f5766c09SAlp Dener PetscBool is_lmvm, is_blmvm; 339f5766c09SAlp Dener Mat M; 340f5766c09SAlp Dener 341b39c12a9SAlp Dener PetscFunctionBegin; 3429566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm)); 3439566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm)); 344f5766c09SAlp Dener if (is_lmvm) { 345f5766c09SAlp Dener lmP = (TAO_LMVM *)tao->data; 346f5766c09SAlp Dener M = lmP->M; 347f5766c09SAlp Dener } else if (is_blmvm) { 348f5766c09SAlp Dener blmP = (TAO_BLMVM *)tao->data; 349f5766c09SAlp Dener M = blmP->M; 3508854b543SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM."); 3519566063dSJacob Faibussowitsch PetscCall(MatLMVMGetJ0(M, H0)); 352f5766c09SAlp Dener PetscFunctionReturn(0); 353f5766c09SAlp Dener } 354f5766c09SAlp Dener 3551bb2a437SAlp Dener /*@ 3561bb2a437SAlp Dener TaoLMVMGetH0KSP - Get the iterative solver for applying the inverse of the QN initial Hessian 3571bb2a437SAlp Dener 3581bb2a437SAlp Dener Input Parameters: 3591bb2a437SAlp Dener . tao - the Tao solver context 3601bb2a437SAlp Dener 3611bb2a437SAlp Dener Output Parameters: 3621bb2a437SAlp Dener . ksp - KSP solver context for the initial Hessian 3631bb2a437SAlp Dener 3641bb2a437SAlp Dener Level: advanced 3651bb2a437SAlp Dener 366db781477SPatrick Sanan .seealso: `TaoLMVMGetH0()`, `TaoLMVMGetH0KSP()` 3671bb2a437SAlp Dener @*/ 368f5766c09SAlp Dener PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp) 369f5766c09SAlp Dener { 370f5766c09SAlp Dener TAO_LMVM *lmP; 371f5766c09SAlp Dener TAO_BLMVM *blmP; 372f5766c09SAlp Dener PetscBool is_lmvm, is_blmvm; 373f5766c09SAlp Dener Mat M; 374f5766c09SAlp Dener 3758854b543SStefano Zampini PetscFunctionBegin; 3769566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm)); 3779566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm)); 378f5766c09SAlp Dener if (is_lmvm) { 379f5766c09SAlp Dener lmP = (TAO_LMVM *)tao->data; 380f5766c09SAlp Dener M = lmP->M; 381f5766c09SAlp Dener } else if (is_blmvm) { 382f5766c09SAlp Dener blmP = (TAO_BLMVM *)tao->data; 383f5766c09SAlp Dener M = blmP->M; 3848854b543SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM."); 3859566063dSJacob Faibussowitsch PetscCall(MatLMVMGetJ0KSP(M, ksp)); 386a9603a14SPatrick Farrell PetscFunctionReturn(0); 387a9603a14SPatrick Farrell } 388