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)); 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)); 1129566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 113f5766c09SAlp Dener } 114f5766c09SAlp Dener PetscFunctionReturn(0); 115f5766c09SAlp Dener } 116f5766c09SAlp Dener 117f5766c09SAlp Dener static PetscErrorCode TaoSetup_BLMVM(Tao tao) 118f5766c09SAlp Dener { 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)); 126f5766c09SAlp Dener 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 if (!tao->XL) { 1349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->XL)); 1359566063dSJacob Faibussowitsch PetscCall(VecSet(tao->XL,PETSC_NINFINITY)); 136f5766c09SAlp Dener } 137f5766c09SAlp Dener if (!tao->XU) { 1389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->XU)); 1399566063dSJacob Faibussowitsch PetscCall(VecSet(tao->XU,PETSC_INFINITY)); 140f5766c09SAlp Dener } 141f5766c09SAlp Dener /* Allocate matrix for the limited memory approximation */ 1429566063dSJacob Faibussowitsch PetscCall(MatLMVMAllocate(blmP->M,tao->solution,blmP->unprojected_gradient)); 143f5766c09SAlp Dener 144f5766c09SAlp Dener /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */ 145f5766c09SAlp Dener if (blmP->H0) { 1469566063dSJacob Faibussowitsch PetscCall(MatLMVMSetJ0(blmP->M, blmP->H0)); 147f5766c09SAlp Dener } 148f5766c09SAlp Dener PetscFunctionReturn(0); 149f5766c09SAlp Dener } 150f5766c09SAlp Dener 151f5766c09SAlp Dener /* ---------------------------------------------------------- */ 152f5766c09SAlp Dener static PetscErrorCode TaoDestroy_BLMVM(Tao tao) 153f5766c09SAlp Dener { 154f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 155f5766c09SAlp Dener 156f5766c09SAlp Dener PetscFunctionBegin; 157f5766c09SAlp Dener if (tao->setupcalled) { 1589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->unprojected_gradient)); 1599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->Xold)); 1609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&blmP->Gold)); 161f5766c09SAlp Dener } 1629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&blmP->M)); 163f5766c09SAlp Dener if (blmP->H0) { 164f5766c09SAlp Dener PetscObjectDereference((PetscObject)blmP->H0); 165f5766c09SAlp Dener } 1669566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 167f5766c09SAlp Dener PetscFunctionReturn(0); 168f5766c09SAlp Dener } 169f5766c09SAlp Dener 170f5766c09SAlp Dener /*------------------------------------------------------------*/ 1712ec5c1acSAlp Dener static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao) 172a7e14dcfSSatish Balay { 173f5766c09SAlp Dener TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data; 174e5fecd4eSAlp Dener PetscBool is_spd; 175a7e14dcfSSatish Balay 176a7e14dcfSSatish Balay PetscFunctionBegin; 177*d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization"); 1789566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL)); 179*d0609cedSBarry Smith PetscOptionsHeadEnd(); 1809566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(blmP->M, ((PetscObject)tao)->prefix)); 1819566063dSJacob Faibussowitsch PetscCall(MatAppendOptionsPrefix(blmP->M, "tao_blmvm_")); 1829566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(blmP->M)); 1839566063dSJacob Faibussowitsch PetscCall(MatGetOption(blmP->M, MAT_SPD, &is_spd)); 1843c859ba3SBarry Smith PetscCheck(is_spd,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite"); 185a7e14dcfSSatish Balay PetscFunctionReturn(0); 186a7e14dcfSSatish Balay } 187a7e14dcfSSatish Balay 188f5766c09SAlp Dener /*------------------------------------------------------------*/ 1895bd1e576SStefano Zampini static PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer) 190a7e14dcfSSatish Balay { 191f5766c09SAlp Dener TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data; 192f5766c09SAlp Dener PetscBool isascii; 193a7e14dcfSSatish Balay 194a7e14dcfSSatish Balay PetscFunctionBegin; 1959566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 196f5766c09SAlp Dener if (isascii) { 1979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad)); 1989566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 1999566063dSJacob Faibussowitsch PetscCall(MatView(lmP->M, viewer)); 2009566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 201f5766c09SAlp Dener } 202f5766c09SAlp Dener PetscFunctionReturn(0); 203f5766c09SAlp Dener } 204f5766c09SAlp Dener 205f5766c09SAlp Dener static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU) 206f5766c09SAlp Dener { 207f5766c09SAlp Dener TAO_BLMVM *blm = (TAO_BLMVM *) tao->data; 208f5766c09SAlp Dener 209f5766c09SAlp Dener PetscFunctionBegin; 210f5766c09SAlp Dener PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 211f5766c09SAlp Dener PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 212f5766c09SAlp Dener PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 2133c859ba3SBarry Smith PetscCheck(tao->gradient && blm->unprojected_gradient,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist."); 214f5766c09SAlp Dener 2159566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient,DXL)); 2169566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXL,-1.0,blm->unprojected_gradient)); 2179566063dSJacob Faibussowitsch PetscCall(VecSet(DXU,0.0)); 2189566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(DXL,DXL,DXU)); 219f5766c09SAlp Dener 2209566063dSJacob Faibussowitsch PetscCall(VecCopy(blm->unprojected_gradient,DXU)); 2219566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU,-1.0,tao->gradient)); 2229566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU,1.0,DXL)); 223f5766c09SAlp Dener PetscFunctionReturn(0); 224f5766c09SAlp Dener } 225f5766c09SAlp Dener 226f5766c09SAlp Dener /* ---------------------------------------------------------- */ 227f5766c09SAlp Dener /*MC 228f5766c09SAlp Dener TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method 229f5766c09SAlp Dener for nonlinear minimization with bound constraints. It is an extension 230f5766c09SAlp Dener of TAOLMVM 231f5766c09SAlp Dener 232f5766c09SAlp Dener Options Database Keys: 233f5766c09SAlp Dener . -tao_lmm_recycle - enable recycling of LMVM information between subsequent TaoSolve calls 234f5766c09SAlp Dener 235f5766c09SAlp Dener Level: beginner 236f5766c09SAlp Dener M*/ 237f5766c09SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao) 238f5766c09SAlp Dener { 239f5766c09SAlp Dener TAO_BLMVM *blmP; 240f5766c09SAlp Dener const char *morethuente_type = TAOLINESEARCHMT; 241f5766c09SAlp Dener 242f5766c09SAlp Dener PetscFunctionBegin; 243f5766c09SAlp Dener tao->ops->setup = TaoSetup_BLMVM; 244f5766c09SAlp Dener tao->ops->solve = TaoSolve_BLMVM; 245f5766c09SAlp Dener tao->ops->view = TaoView_BLMVM; 246a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_BLMVM; 247f5766c09SAlp Dener tao->ops->destroy = TaoDestroy_BLMVM; 248f5766c09SAlp Dener tao->ops->computedual = TaoComputeDual_BLMVM; 249f5766c09SAlp Dener 2509566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&blmP)); 251f5766c09SAlp Dener blmP->H0 = NULL; 252f5766c09SAlp Dener blmP->recycle = PETSC_FALSE; 253f5766c09SAlp Dener tao->data = (void*)blmP; 254f5766c09SAlp Dener 255f5766c09SAlp Dener /* Override default settings (unless already changed) */ 256f5766c09SAlp Dener if (!tao->max_it_changed) tao->max_it = 2000; 257f5766c09SAlp Dener if (!tao->max_funcs_changed) tao->max_funcs = 4000; 258f5766c09SAlp Dener 2599566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 2619566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type)); 2629566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch,tao)); 263f5766c09SAlp Dener 2649566063dSJacob Faibussowitsch PetscCall(KSPInitializePackage()); 2659566063dSJacob Faibussowitsch PetscCall(MatCreate(((PetscObject)tao)->comm, &blmP->M)); 2669566063dSJacob Faibussowitsch PetscCall(MatSetType(blmP->M, MATLMVMBFGS)); 2679566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1)); 268f5766c09SAlp Dener PetscFunctionReturn(0); 269f5766c09SAlp Dener } 270f5766c09SAlp Dener 2711bb2a437SAlp Dener /*@ 2721bb2a437SAlp Dener TaoLMVMRecycle - Enable/disable recycling of the QN history between subsequent TaoSolve calls. 2731bb2a437SAlp Dener 2741bb2a437SAlp Dener Input Parameters: 2751bb2a437SAlp Dener + tao - the Tao solver context 2761bb2a437SAlp Dener - flg - Boolean flag for recycling (PETSC_TRUE or PETSC_FALSE) 2771bb2a437SAlp Dener 2781bb2a437SAlp Dener Level: intermediate 2791bb2a437SAlp Dener @*/ 280b39c12a9SAlp Dener PetscErrorCode TaoLMVMRecycle(Tao tao, PetscBool flg) 281b39c12a9SAlp Dener { 282b39c12a9SAlp Dener TAO_LMVM *lmP; 283b39c12a9SAlp Dener TAO_BLMVM *blmP; 284b39c12a9SAlp Dener PetscBool is_lmvm, is_blmvm; 285b39c12a9SAlp Dener 286b39c12a9SAlp Dener PetscFunctionBegin; 2879566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm)); 2889566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm)); 289b39c12a9SAlp Dener if (is_lmvm) { 290b39c12a9SAlp Dener lmP = (TAO_LMVM *)tao->data; 291b39c12a9SAlp Dener lmP->recycle = flg; 292b39c12a9SAlp Dener } else if (is_blmvm) { 293b39c12a9SAlp Dener blmP = (TAO_BLMVM *)tao->data; 294b39c12a9SAlp Dener blmP->recycle = flg; 2958854b543SStefano Zampini } 296b39c12a9SAlp Dener PetscFunctionReturn(0); 297b39c12a9SAlp Dener } 298b39c12a9SAlp Dener 2991bb2a437SAlp Dener /*@ 3001bb2a437SAlp Dener TaoLMVMSetH0 - Set the initial Hessian for the QN approximation 3011bb2a437SAlp Dener 3021bb2a437SAlp Dener Input Parameters: 3031bb2a437SAlp Dener + tao - the Tao solver context 3041bb2a437SAlp Dener - H0 - Mat object for the initial Hessian 3051bb2a437SAlp Dener 3061bb2a437SAlp Dener Level: advanced 3071bb2a437SAlp Dener 3081bb2a437SAlp Dener .seealso: TaoLMVMGetH0(), TaoLMVMGetH0KSP() 3091bb2a437SAlp Dener @*/ 310f5766c09SAlp Dener PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0) 311f5766c09SAlp Dener { 312f5766c09SAlp Dener TAO_LMVM *lmP; 313f5766c09SAlp Dener TAO_BLMVM *blmP; 314f5766c09SAlp Dener PetscBool is_lmvm, is_blmvm; 315f5766c09SAlp Dener 316b39c12a9SAlp Dener PetscFunctionBegin; 3179566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm)); 3189566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm)); 319f5766c09SAlp Dener if (is_lmvm) { 320f5766c09SAlp Dener lmP = (TAO_LMVM *)tao->data; 3219566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H0)); 322f5766c09SAlp Dener lmP->H0 = H0; 323f5766c09SAlp Dener } else if (is_blmvm) { 324f5766c09SAlp Dener blmP = (TAO_BLMVM *)tao->data; 3259566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)H0)); 326f5766c09SAlp Dener blmP->H0 = H0; 3278854b543SStefano Zampini } 328f5766c09SAlp Dener PetscFunctionReturn(0); 329f5766c09SAlp Dener } 330f5766c09SAlp Dener 3311bb2a437SAlp Dener /*@ 3321bb2a437SAlp Dener TaoLMVMGetH0 - Get the matrix object for the QN initial Hessian 3331bb2a437SAlp Dener 3341bb2a437SAlp Dener Input Parameters: 3351bb2a437SAlp Dener . tao - the Tao solver context 3361bb2a437SAlp Dener 3371bb2a437SAlp Dener Output Parameters: 3381bb2a437SAlp Dener . H0 - Mat object for the initial Hessian 3391bb2a437SAlp Dener 3401bb2a437SAlp Dener Level: advanced 3411bb2a437SAlp Dener 3421bb2a437SAlp Dener .seealso: TaoLMVMSetH0(), TaoLMVMGetH0KSP() 3431bb2a437SAlp Dener @*/ 344f5766c09SAlp Dener PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0) 345f5766c09SAlp Dener { 346f5766c09SAlp Dener TAO_LMVM *lmP; 347f5766c09SAlp Dener TAO_BLMVM *blmP; 348f5766c09SAlp Dener PetscBool is_lmvm, is_blmvm; 349f5766c09SAlp Dener Mat M; 350f5766c09SAlp Dener 351b39c12a9SAlp Dener PetscFunctionBegin; 3529566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm)); 3539566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm)); 354f5766c09SAlp Dener if (is_lmvm) { 355f5766c09SAlp Dener lmP = (TAO_LMVM *)tao->data; 356f5766c09SAlp Dener M = lmP->M; 357f5766c09SAlp Dener } else if (is_blmvm) { 358f5766c09SAlp Dener blmP = (TAO_BLMVM *)tao->data; 359f5766c09SAlp Dener M = blmP->M; 3608854b543SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM."); 3619566063dSJacob Faibussowitsch PetscCall(MatLMVMGetJ0(M, H0)); 362f5766c09SAlp Dener PetscFunctionReturn(0); 363f5766c09SAlp Dener } 364f5766c09SAlp Dener 3651bb2a437SAlp Dener /*@ 3661bb2a437SAlp Dener TaoLMVMGetH0KSP - Get the iterative solver for applying the inverse of the QN initial Hessian 3671bb2a437SAlp Dener 3681bb2a437SAlp Dener Input Parameters: 3691bb2a437SAlp Dener . tao - the Tao solver context 3701bb2a437SAlp Dener 3711bb2a437SAlp Dener Output Parameters: 3721bb2a437SAlp Dener . ksp - KSP solver context for the initial Hessian 3731bb2a437SAlp Dener 3741bb2a437SAlp Dener Level: advanced 3751bb2a437SAlp Dener 3761bb2a437SAlp Dener .seealso: TaoLMVMGetH0(), TaoLMVMGetH0KSP() 3771bb2a437SAlp Dener @*/ 378f5766c09SAlp Dener PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp) 379f5766c09SAlp Dener { 380f5766c09SAlp Dener TAO_LMVM *lmP; 381f5766c09SAlp Dener TAO_BLMVM *blmP; 382f5766c09SAlp Dener PetscBool is_lmvm, is_blmvm; 383f5766c09SAlp Dener Mat M; 384f5766c09SAlp Dener 3858854b543SStefano Zampini PetscFunctionBegin; 3869566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOLMVM,&is_lmvm)); 3879566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tao,TAOBLMVM,&is_blmvm)); 388f5766c09SAlp Dener if (is_lmvm) { 389f5766c09SAlp Dener lmP = (TAO_LMVM *)tao->data; 390f5766c09SAlp Dener M = lmP->M; 391f5766c09SAlp Dener } else if (is_blmvm) { 392f5766c09SAlp Dener blmP = (TAO_BLMVM *)tao->data; 393f5766c09SAlp Dener M = blmP->M; 3948854b543SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONG, "This routine applies to TAO_LMVM and TAO_BLMVM."); 3959566063dSJacob Faibussowitsch PetscCall(MatLMVMGetJ0KSP(M, ksp)); 396a9603a14SPatrick Farrell PetscFunctionReturn(0); 397a9603a14SPatrick Farrell } 398