1aaa7dc30SBarry Smith #include <../src/tao/bound/impls/tron/tron.h> 2aaa7dc30SBarry Smith #include <../src/tao/matrix/submatfree.h> 3a7e14dcfSSatish Balay 4a7e14dcfSSatish Balay /* TRON Routines */ 5441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao,TAO_TRON*); 6a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 7441846f8SBarry Smith static PetscErrorCode TaoDestroy_TRON(Tao tao) 8a7e14dcfSSatish Balay { 9a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 10a7e14dcfSSatish Balay 11a7e14dcfSSatish Balay PetscFunctionBegin; 129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tron->X_New)); 139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tron->G_New)); 149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tron->Work)); 159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tron->DXFree)); 169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tron->R)); 179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tron->diag)); 189566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&tron->scatter)); 199566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tron->Free_Local)); 209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tron->H_sub)); 219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tron->Hpre_sub)); 229566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 23a7e14dcfSSatish Balay PetscFunctionReturn(0); 24a7e14dcfSSatish Balay } 25a7e14dcfSSatish Balay 26a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 274416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_TRON(PetscOptionItems *PetscOptionsObject,Tao tao) 28a7e14dcfSSatish Balay { 29a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 30a7e14dcfSSatish Balay PetscBool flg; 31a7e14dcfSSatish Balay 32a7e14dcfSSatish Balay PetscFunctionBegin; 33d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Newton Trust Region Method for bound constrained optimization"); 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_tron_maxgpits","maximum number of gradient projections per TRON iterate","TaoSetMaxGPIts",tron->maxgpits,&tron->maxgpits,&flg)); 35d0609cedSBarry Smith PetscOptionsHeadEnd(); 369566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp)); 37a7e14dcfSSatish Balay PetscFunctionReturn(0); 38a7e14dcfSSatish Balay } 39a7e14dcfSSatish Balay 40a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 41441846f8SBarry Smith static PetscErrorCode TaoView_TRON(Tao tao, PetscViewer viewer) 42a7e14dcfSSatish Balay { 43a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 44a7e14dcfSSatish Balay PetscBool isascii; 45a7e14dcfSSatish Balay 46a7e14dcfSSatish Balay PetscFunctionBegin; 479566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); 48a7e14dcfSSatish Balay if (isascii) { 49*63a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Total PG its: %" PetscInt_FMT ",",tron->total_gp_its)); 509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"PG tolerance: %g \n",(double)tron->pg_ftol)); 51a7e14dcfSSatish Balay } 52a7e14dcfSSatish Balay PetscFunctionReturn(0); 53a7e14dcfSSatish Balay } 54a7e14dcfSSatish Balay 55a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 56441846f8SBarry Smith static PetscErrorCode TaoSetup_TRON(Tao tao) 57a7e14dcfSSatish Balay { 58a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 59a7e14dcfSSatish Balay 60a7e14dcfSSatish Balay PetscFunctionBegin; 61a7e14dcfSSatish Balay 62a7e14dcfSSatish Balay /* Allocate some arrays */ 639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tron->diag)); 649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tron->X_New)); 659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tron->G_New)); 669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tron->Work)); 679566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 69a7e14dcfSSatish Balay if (!tao->XL) { 709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->XL)); 719566063dSJacob Faibussowitsch PetscCall(VecSet(tao->XL, PETSC_NINFINITY)); 72a7e14dcfSSatish Balay } 73a7e14dcfSSatish Balay if (!tao->XU) { 749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->XU)); 759566063dSJacob Faibussowitsch PetscCall(VecSet(tao->XU, PETSC_INFINITY)); 76a7e14dcfSSatish Balay } 77a7e14dcfSSatish Balay PetscFunctionReturn(0); 78a7e14dcfSSatish Balay } 79a7e14dcfSSatish Balay 80441846f8SBarry Smith static PetscErrorCode TaoSolve_TRON(Tao tao) 8147a47007SBarry Smith { 82a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 838931d482SJason Sarich PetscInt its; 84e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING; 85a7e14dcfSSatish Balay PetscReal prered,actred,delta,f,f_new,rhok,gdx,xdiff,stepsize; 8647a47007SBarry Smith 87a7e14dcfSSatish Balay PetscFunctionBegin; 88a7e14dcfSSatish Balay tron->pgstepsize = 1.0; 89a7e14dcfSSatish Balay tao->trust = tao->trust0; 90a7e14dcfSSatish Balay /* Project the current point onto the feasible set */ 919566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 929566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU)); 93a7e14dcfSSatish Balay 94ce902467SBarry Smith /* Project the initial point onto the feasible region */ 959566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL,tao->solution,tao->XU,tao->solution)); 9653506e15SBarry Smith 97ce902467SBarry Smith /* Compute the objective function and gradient */ 989566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao,tao->solution,&tron->f,tao->gradient)); 999566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 1003c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(tron->f) && !PetscIsInfOrNanReal(tron->gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 101a7e14dcfSSatish Balay 102a7e14dcfSSatish Balay /* Project the gradient and calculate the norm */ 1039566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient)); 1049566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 105a7e14dcfSSatish Balay 106ce902467SBarry Smith /* Initialize trust region radius */ 107ce902467SBarry Smith tao->trust=tao->trust0; 108a7e14dcfSSatish Balay if (tao->trust <= 0) { 109a7e14dcfSSatish Balay tao->trust=PetscMax(tron->gnorm*tron->gnorm,1.0); 110a7e14dcfSSatish Balay } 111a7e14dcfSSatish Balay 112ce902467SBarry Smith /* Initialize step sizes for the line searches */ 113ce902467SBarry Smith tron->pgstepsize=1.0; 114a7e14dcfSSatish Balay tron->stepsize=tao->trust; 115ce902467SBarry Smith 1163ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 1179566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its)); 1189566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,tron->stepsize)); 1199566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 1203ecd9318SAlp Dener while (tao->reason==TAO_CONTINUE_ITERATING) { 121e1e80dc8SAlp Dener /* Call general purpose update function */ 122e1e80dc8SAlp Dener if (tao->ops->update) { 1239566063dSJacob Faibussowitsch PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update)); 124e1e80dc8SAlp Dener } 125ce902467SBarry Smith 126ce902467SBarry Smith /* Perform projected gradient iterations */ 1279566063dSJacob Faibussowitsch PetscCall(TronGradientProjections(tao,tron)); 128ce902467SBarry Smith 1299566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient)); 1309566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 131ce902467SBarry Smith 132ce902467SBarry Smith tao->ksp_its=0; 133a7e14dcfSSatish Balay f=tron->f; delta=tao->trust; 134a7e14dcfSSatish Balay tron->n_free_last = tron->n_free; 1359566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre)); 136a7e14dcfSSatish Balay 137baa3a814SAlp Dener /* Generate index set (IS) of which bound constraints are active */ 1389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tron->Free_Local)); 1399566063dSJacob Faibussowitsch PetscCall(VecWhichInactive(tao->XL,tao->solution,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local)); 1409566063dSJacob Faibussowitsch PetscCall(ISGetSize(tron->Free_Local, &tron->n_free)); 141a7e14dcfSSatish Balay 142a7e14dcfSSatish Balay /* If no free variables */ 143a7e14dcfSSatish Balay if (tron->n_free == 0) { 1449566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 1459566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its)); 1469566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,delta)); 1479566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 1483ecd9318SAlp Dener if (!tao->reason) { 1493ecd9318SAlp Dener tao->reason = TAO_CONVERGED_STEPTOL; 15055723ba5SJason Sarich } 151a7e14dcfSSatish Balay break; 152a7e14dcfSSatish Balay } 153a7e14dcfSSatish Balay /* use free_local to mask/submat gradient, hessian, stepdirection */ 1549566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->R)); 1559566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->DXFree)); 1569566063dSJacob Faibussowitsch PetscCall(VecSet(tron->DXFree,0.0)); 1579566063dSJacob Faibussowitsch PetscCall(VecScale(tron->R, -1.0)); 1589566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub)); 159a7e14dcfSSatish Balay if (tao->hessian == tao->hessian_pre) { 1609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tron->Hpre_sub)); 1619566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(tron->H_sub))); 162a7e14dcfSSatish Balay tron->Hpre_sub = tron->H_sub; 163a7e14dcfSSatish Balay } else { 1649566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type,&tron->Hpre_sub)); 165a7e14dcfSSatish Balay } 1669566063dSJacob Faibussowitsch PetscCall(KSPReset(tao->ksp)); 1679566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub)); 168a7e14dcfSSatish Balay while (1) { 169a7e14dcfSSatish Balay 170a7e14dcfSSatish Balay /* Approximately solve the reduced linear system */ 1719566063dSJacob Faibussowitsch PetscCall(KSPCGSetRadius(tao->ksp,delta)); 172d8ec8e84SJason Sarich 1739566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, tron->R, tron->DXFree)); 1749566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&its)); 175a7e14dcfSSatish Balay tao->ksp_its+=its; 176ae93cb3cSJason Sarich tao->ksp_tot_its+=its; 1779566063dSJacob Faibussowitsch PetscCall(VecSet(tao->stepdirection,0.0)); 178a7e14dcfSSatish Balay 179a7e14dcfSSatish Balay /* Add dxfree matrix to compute step direction vector */ 1809566063dSJacob Faibussowitsch PetscCall(VecISAXPY(tao->stepdirection,tron->Free_Local,1.0,tron->DXFree)); 181a7e14dcfSSatish Balay 1829566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx)); 1839566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, tron->X_New)); 1849566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, tron->G_New)); 185a7e14dcfSSatish Balay 186a7e14dcfSSatish Balay stepsize=1.0;f_new=f; 187a7e14dcfSSatish Balay 1889566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0)); 1899566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection,&stepsize,&ls_reason)); 1909566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 191a7e14dcfSSatish Balay 1929566063dSJacob Faibussowitsch PetscCall(MatMult(tao->hessian, tao->stepdirection, tron->Work)); 1939566063dSJacob Faibussowitsch PetscCall(VecAYPX(tron->Work, 0.5, tao->gradient)); 1949566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, tron->Work, &prered)); 195a7e14dcfSSatish Balay actred = f_new - f; 196ce902467SBarry Smith if ((PetscAbsScalar(actred) <= 1e-6) && (PetscAbsScalar(prered) <= 1e-6)) { 197ce902467SBarry Smith rhok = 1.0; 198ce902467SBarry Smith } else if (actred<0) { 199a7e14dcfSSatish Balay rhok=PetscAbs(-actred/prered); 200a7e14dcfSSatish Balay } else { 201a7e14dcfSSatish Balay rhok=0.0; 202a7e14dcfSSatish Balay } 203a7e14dcfSSatish Balay 204a7e14dcfSSatish Balay /* Compare actual improvement to the quadratic model */ 205a7e14dcfSSatish Balay if (rhok > tron->eta1) { /* Accept the point */ 206a7e14dcfSSatish Balay /* d = x_new - x */ 2079566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->X_New, tao->stepdirection)); 2089566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->stepdirection, -1.0, tao->solution)); 209a7e14dcfSSatish Balay 2109566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->stepdirection, NORM_2, &xdiff)); 211a7e14dcfSSatish Balay xdiff *= stepsize; 212a7e14dcfSSatish Balay 213a7e14dcfSSatish Balay /* Adjust trust region size */ 214a7e14dcfSSatish Balay if (rhok < tron->eta2) { 215a7e14dcfSSatish Balay delta = PetscMin(xdiff,delta)*tron->sigma1; 216a7e14dcfSSatish Balay } else if (rhok > tron->eta4) { 217a7e14dcfSSatish Balay delta= PetscMin(xdiff,delta)*tron->sigma3; 218a7e14dcfSSatish Balay } else if (rhok > tron->eta3) { 219a7e14dcfSSatish Balay delta=PetscMin(xdiff,delta)*tron->sigma2; 220a7e14dcfSSatish Balay } 2219566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tron->G_New,tron->X_New, tao->XL, tao->XU, tao->gradient)); 2229566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tron->Free_Local)); 2239566063dSJacob Faibussowitsch PetscCall(VecWhichInactive(tao->XL,tron->X_New,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local)); 224a7e14dcfSSatish Balay f=f_new; 2259566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 2269566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->X_New, tao->solution)); 2279566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->G_New, tao->gradient)); 228a7e14dcfSSatish Balay break; 229a7e14dcfSSatish Balay } 230a7e14dcfSSatish Balay else if (delta <= 1e-30) { 231a7e14dcfSSatish Balay break; 232a7e14dcfSSatish Balay } 233a7e14dcfSSatish Balay else { 234a7e14dcfSSatish Balay delta /= 4.0; 235a7e14dcfSSatish Balay } 236a7e14dcfSSatish Balay } /* end linear solve loop */ 237a7e14dcfSSatish Balay 238a7e14dcfSSatish Balay tron->f=f; tron->actred=actred; tao->trust=delta; 2398931d482SJason Sarich tao->niter++; 2409566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its)); 2419566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,stepsize)); 2429566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 243a7e14dcfSSatish Balay } /* END MAIN LOOP */ 244a7e14dcfSSatish Balay PetscFunctionReturn(0); 245a7e14dcfSSatish Balay } 246a7e14dcfSSatish Balay 247441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao tao,TAO_TRON *tron) 248a7e14dcfSSatish Balay { 249a7e14dcfSSatish Balay PetscInt i; 250e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 251a7e14dcfSSatish Balay PetscReal actred=-1.0,actred_max=0.0; 252a7e14dcfSSatish Balay PetscReal f_new; 253a7e14dcfSSatish Balay /* 254a7e14dcfSSatish Balay The gradient and function value passed into and out of this 255a7e14dcfSSatish Balay routine should be current and correct. 256a7e14dcfSSatish Balay 257a7e14dcfSSatish Balay The free, active, and binding variables should be already identified 258a7e14dcfSSatish Balay */ 259a7e14dcfSSatish Balay PetscFunctionBegin; 260a7e14dcfSSatish Balay 261ce902467SBarry Smith for (i=0;i<tron->maxgpits;++i) { 262a7e14dcfSSatish Balay 263a7e14dcfSSatish Balay if (-actred <= (tron->pg_ftol)*actred_max) break; 264a7e14dcfSSatish Balay 265ce902467SBarry Smith ++tron->gp_iterates; 266ce902467SBarry Smith ++tron->total_gp_its; 267a7e14dcfSSatish Balay f_new=tron->f; 268a7e14dcfSSatish Balay 2699566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient,tao->stepdirection)); 2709566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection,-1.0)); 2719566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,tron->pgstepsize)); 272d0609cedSBarry Smith PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection,&tron->pgstepsize, &ls_reason)); 2739566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 274a7e14dcfSSatish Balay 2759566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient)); 2769566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 277ce902467SBarry Smith 278a7e14dcfSSatish Balay /* Update the iterate */ 279a7e14dcfSSatish Balay actred = f_new - tron->f; 280a7e14dcfSSatish Balay actred_max = PetscMax(actred_max,-(f_new - tron->f)); 281a7e14dcfSSatish Balay tron->f = f_new; 282a7e14dcfSSatish Balay } 283a7e14dcfSSatish Balay PetscFunctionReturn(0); 284a7e14dcfSSatish Balay } 285a7e14dcfSSatish Balay 286bdb10af2SPierre Jolivet static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU) 287bdb10af2SPierre Jolivet { 288a7e14dcfSSatish Balay 289a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 290a7e14dcfSSatish Balay 291a7e14dcfSSatish Balay PetscFunctionBegin; 292441846f8SBarry Smith PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 293a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 294a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 2953c859ba3SBarry Smith PetscCheck(tron->Work && tao->gradient,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist."); 296a7e14dcfSSatish Balay 2979566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tron->Work)); 2989566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->Work,DXL)); 2999566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXL,-1.0,tao->gradient)); 3009566063dSJacob Faibussowitsch PetscCall(VecSet(DXU,0.0)); 3019566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(DXL,DXL,DXU)); 302a7e14dcfSSatish Balay 3039566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient,DXU)); 3049566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU,-1.0,tron->Work)); 3059566063dSJacob Faibussowitsch PetscCall(VecSet(tron->Work,0.0)); 3069566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(DXU,tron->Work,DXU)); 307a7e14dcfSSatish Balay PetscFunctionReturn(0); 308a7e14dcfSSatish Balay } 309a7e14dcfSSatish Balay 310a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 3111522df2eSJason Sarich /*MC 3121522df2eSJason Sarich TAOTRON - The TRON algorithm is an active-set Newton trust region method 3131522df2eSJason Sarich for bound-constrained minimization. 3141522df2eSJason Sarich 3151522df2eSJason Sarich Options Database Keys: 3161522df2eSJason Sarich + -tao_tron_maxgpits - maximum number of gradient projections per TRON iterate 3171522df2eSJason Sarich - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets 3181522df2eSJason Sarich 3191eb8069cSJason Sarich Level: beginner 3201522df2eSJason Sarich M*/ 321728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_TRON(Tao tao) 322a7e14dcfSSatish Balay { 323a7e14dcfSSatish Balay TAO_TRON *tron; 3248caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 325a7e14dcfSSatish Balay 32647a47007SBarry Smith PetscFunctionBegin; 327a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_TRON; 328a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_TRON; 329a7e14dcfSSatish Balay tao->ops->view = TaoView_TRON; 330a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_TRON; 331a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_TRON; 332a7e14dcfSSatish Balay tao->ops->computedual = TaoComputeDual_TRON; 333a7e14dcfSSatish Balay 3349566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&tron)); 3356f4723b1SBarry Smith tao->data = (void*)tron; 3366552cf8aSJason Sarich 3376552cf8aSJason Sarich /* Override default settings (unless already changed) */ 3386552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 50; 3396552cf8aSJason Sarich if (!tao->trust0_changed) tao->trust0 = 1.0; 3404f1535bcSTodd Munson if (!tao->steptol_changed) tao->steptol = 0.0; 341a7e14dcfSSatish Balay 342a7e14dcfSSatish Balay /* Initialize pointers and variables */ 343a7e14dcfSSatish Balay tron->n = 0; 344a7e14dcfSSatish Balay tron->maxgpits = 3; 345a7e14dcfSSatish Balay tron->pg_ftol = 0.001; 346a7e14dcfSSatish Balay 347a7e14dcfSSatish Balay tron->eta1 = 1.0e-4; 348a7e14dcfSSatish Balay tron->eta2 = 0.25; 349a7e14dcfSSatish Balay tron->eta3 = 0.50; 350a7e14dcfSSatish Balay tron->eta4 = 0.90; 351a7e14dcfSSatish Balay 352a7e14dcfSSatish Balay tron->sigma1 = 0.5; 353a7e14dcfSSatish Balay tron->sigma2 = 2.0; 354a7e14dcfSSatish Balay tron->sigma3 = 4.0; 355a7e14dcfSSatish Balay 356a7e14dcfSSatish Balay tron->gp_iterates = 0; /* Cumulative number */ 357a7e14dcfSSatish Balay tron->total_gp_its = 0; 358a7e14dcfSSatish Balay tron->n_free = 0; 359a7e14dcfSSatish Balay 3606c23d075SBarry Smith tron->DXFree=NULL; 3616c23d075SBarry Smith tron->R=NULL; 3626c23d075SBarry Smith tron->X_New=NULL; 3636c23d075SBarry Smith tron->G_New=NULL; 3646c23d075SBarry Smith tron->Work=NULL; 3656c23d075SBarry Smith tron->Free_Local=NULL; 3666c23d075SBarry Smith tron->H_sub=NULL; 3676c23d075SBarry Smith tron->Hpre_sub=NULL; 368a7e14dcfSSatish Balay tao->subset_type = TAO_SUBSET_SUBVEC; 369a7e14dcfSSatish Balay 3709566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 3719566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 3729566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch,morethuente_type)); 3739566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch,tao)); 3749566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix)); 375a7e14dcfSSatish Balay 3769566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp)); 3779566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 3789566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix)); 3799566063dSJacob Faibussowitsch PetscCall(KSPSetType(tao->ksp,KSPSTCG)); 380a7e14dcfSSatish Balay PetscFunctionReturn(0); 381a7e14dcfSSatish Balay } 382