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)); 22*a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 239566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 24a7e14dcfSSatish Balay PetscFunctionReturn(0); 25a7e14dcfSSatish Balay } 26a7e14dcfSSatish Balay 27a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 284416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_TRON(PetscOptionItems *PetscOptionsObject,Tao tao) 29a7e14dcfSSatish Balay { 30a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 31a7e14dcfSSatish Balay PetscBool flg; 32a7e14dcfSSatish Balay 33a7e14dcfSSatish Balay PetscFunctionBegin; 34d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Newton Trust Region Method for bound constrained optimization"); 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_tron_maxgpits","maximum number of gradient projections per TRON iterate","TaoSetMaxGPIts",tron->maxgpits,&tron->maxgpits,&flg)); 36d0609cedSBarry Smith PetscOptionsHeadEnd(); 379566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp)); 38a7e14dcfSSatish Balay PetscFunctionReturn(0); 39a7e14dcfSSatish Balay } 40a7e14dcfSSatish Balay 41a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 42441846f8SBarry Smith static PetscErrorCode TaoView_TRON(Tao tao, PetscViewer viewer) 43a7e14dcfSSatish Balay { 44a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 45a7e14dcfSSatish Balay PetscBool isascii; 46a7e14dcfSSatish Balay 47a7e14dcfSSatish Balay PetscFunctionBegin; 489566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); 49a7e14dcfSSatish Balay if (isascii) { 5063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Total PG its: %" PetscInt_FMT ",",tron->total_gp_its)); 519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"PG tolerance: %g \n",(double)tron->pg_ftol)); 52a7e14dcfSSatish Balay } 53a7e14dcfSSatish Balay PetscFunctionReturn(0); 54a7e14dcfSSatish Balay } 55a7e14dcfSSatish Balay 56a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 57441846f8SBarry Smith static PetscErrorCode TaoSetup_TRON(Tao tao) 58a7e14dcfSSatish Balay { 59a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 60a7e14dcfSSatish Balay 61a7e14dcfSSatish Balay PetscFunctionBegin; 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 PetscFunctionReturn(0); 70a7e14dcfSSatish Balay } 71a7e14dcfSSatish Balay 72441846f8SBarry Smith static PetscErrorCode TaoSolve_TRON(Tao tao) 7347a47007SBarry Smith { 74a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 758931d482SJason Sarich PetscInt its; 76e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING; 77a7e14dcfSSatish Balay PetscReal prered,actred,delta,f,f_new,rhok,gdx,xdiff,stepsize; 7847a47007SBarry Smith 79a7e14dcfSSatish Balay PetscFunctionBegin; 80a7e14dcfSSatish Balay tron->pgstepsize = 1.0; 81a7e14dcfSSatish Balay tao->trust = tao->trust0; 82a7e14dcfSSatish Balay /* Project the current point onto the feasible set */ 839566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 849566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU)); 85a7e14dcfSSatish Balay 86ce902467SBarry Smith /* Project the initial point onto the feasible region */ 879566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL,tao->solution,tao->XU,tao->solution)); 8853506e15SBarry Smith 89ce902467SBarry Smith /* Compute the objective function and gradient */ 909566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao,tao->solution,&tron->f,tao->gradient)); 919566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 923c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(tron->f) && !PetscIsInfOrNanReal(tron->gnorm),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN"); 93a7e14dcfSSatish Balay 94a7e14dcfSSatish Balay /* Project the gradient and calculate the norm */ 959566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient)); 969566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 97a7e14dcfSSatish Balay 98ce902467SBarry Smith /* Initialize trust region radius */ 99ce902467SBarry Smith tao->trust=tao->trust0; 100a7e14dcfSSatish Balay if (tao->trust <= 0) { 101a7e14dcfSSatish Balay tao->trust=PetscMax(tron->gnorm*tron->gnorm,1.0); 102a7e14dcfSSatish Balay } 103a7e14dcfSSatish Balay 104ce902467SBarry Smith /* Initialize step sizes for the line searches */ 105ce902467SBarry Smith tron->pgstepsize=1.0; 106a7e14dcfSSatish Balay tron->stepsize=tao->trust; 107ce902467SBarry Smith 1083ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 1099566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its)); 1109566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,tron->stepsize)); 1119566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 1123ecd9318SAlp Dener while (tao->reason==TAO_CONTINUE_ITERATING) { 113e1e80dc8SAlp Dener /* Call general purpose update function */ 114e1e80dc8SAlp Dener if (tao->ops->update) { 1159566063dSJacob Faibussowitsch PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update)); 116e1e80dc8SAlp Dener } 117ce902467SBarry Smith 118ce902467SBarry Smith /* Perform projected gradient iterations */ 1199566063dSJacob Faibussowitsch PetscCall(TronGradientProjections(tao,tron)); 120ce902467SBarry Smith 1219566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient)); 1229566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 123ce902467SBarry Smith 124ce902467SBarry Smith tao->ksp_its=0; 125a7e14dcfSSatish Balay f=tron->f; delta=tao->trust; 126a7e14dcfSSatish Balay tron->n_free_last = tron->n_free; 1279566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre)); 128a7e14dcfSSatish Balay 129baa3a814SAlp Dener /* Generate index set (IS) of which bound constraints are active */ 1309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tron->Free_Local)); 1319566063dSJacob Faibussowitsch PetscCall(VecWhichInactive(tao->XL,tao->solution,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local)); 1329566063dSJacob Faibussowitsch PetscCall(ISGetSize(tron->Free_Local, &tron->n_free)); 133a7e14dcfSSatish Balay 134a7e14dcfSSatish Balay /* If no free variables */ 135a7e14dcfSSatish Balay if (tron->n_free == 0) { 1369566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 1379566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its)); 1389566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,delta)); 1399566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 1403ecd9318SAlp Dener if (!tao->reason) { 1413ecd9318SAlp Dener tao->reason = TAO_CONVERGED_STEPTOL; 14255723ba5SJason Sarich } 143a7e14dcfSSatish Balay break; 144a7e14dcfSSatish Balay } 145a7e14dcfSSatish Balay /* use free_local to mask/submat gradient, hessian, stepdirection */ 1469566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->R)); 1479566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->DXFree)); 1489566063dSJacob Faibussowitsch PetscCall(VecSet(tron->DXFree,0.0)); 1499566063dSJacob Faibussowitsch PetscCall(VecScale(tron->R, -1.0)); 1509566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub)); 151a7e14dcfSSatish Balay if (tao->hessian == tao->hessian_pre) { 1529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tron->Hpre_sub)); 1539566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(tron->H_sub))); 154a7e14dcfSSatish Balay tron->Hpre_sub = tron->H_sub; 155a7e14dcfSSatish Balay } else { 1569566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type,&tron->Hpre_sub)); 157a7e14dcfSSatish Balay } 1589566063dSJacob Faibussowitsch PetscCall(KSPReset(tao->ksp)); 1599566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub)); 160a7e14dcfSSatish Balay while (1) { 161a7e14dcfSSatish Balay 162a7e14dcfSSatish Balay /* Approximately solve the reduced linear system */ 1639566063dSJacob Faibussowitsch PetscCall(KSPCGSetRadius(tao->ksp,delta)); 164d8ec8e84SJason Sarich 1659566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, tron->R, tron->DXFree)); 1669566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&its)); 167a7e14dcfSSatish Balay tao->ksp_its+=its; 168ae93cb3cSJason Sarich tao->ksp_tot_its+=its; 1699566063dSJacob Faibussowitsch PetscCall(VecSet(tao->stepdirection,0.0)); 170a7e14dcfSSatish Balay 171a7e14dcfSSatish Balay /* Add dxfree matrix to compute step direction vector */ 1729566063dSJacob Faibussowitsch PetscCall(VecISAXPY(tao->stepdirection,tron->Free_Local,1.0,tron->DXFree)); 173a7e14dcfSSatish Balay 1749566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx)); 1759566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, tron->X_New)); 1769566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, tron->G_New)); 177a7e14dcfSSatish Balay 178a7e14dcfSSatish Balay stepsize=1.0;f_new=f; 179a7e14dcfSSatish Balay 1809566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0)); 1819566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection,&stepsize,&ls_reason)); 1829566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 183a7e14dcfSSatish Balay 1849566063dSJacob Faibussowitsch PetscCall(MatMult(tao->hessian, tao->stepdirection, tron->Work)); 1859566063dSJacob Faibussowitsch PetscCall(VecAYPX(tron->Work, 0.5, tao->gradient)); 1869566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, tron->Work, &prered)); 187a7e14dcfSSatish Balay actred = f_new - f; 188ce902467SBarry Smith if ((PetscAbsScalar(actred) <= 1e-6) && (PetscAbsScalar(prered) <= 1e-6)) { 189ce902467SBarry Smith rhok = 1.0; 190ce902467SBarry Smith } else if (actred<0) { 191a7e14dcfSSatish Balay rhok=PetscAbs(-actred/prered); 192a7e14dcfSSatish Balay } else { 193a7e14dcfSSatish Balay rhok=0.0; 194a7e14dcfSSatish Balay } 195a7e14dcfSSatish Balay 196a7e14dcfSSatish Balay /* Compare actual improvement to the quadratic model */ 197a7e14dcfSSatish Balay if (rhok > tron->eta1) { /* Accept the point */ 198a7e14dcfSSatish Balay /* d = x_new - x */ 1999566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->X_New, tao->stepdirection)); 2009566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->stepdirection, -1.0, tao->solution)); 201a7e14dcfSSatish Balay 2029566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->stepdirection, NORM_2, &xdiff)); 203a7e14dcfSSatish Balay xdiff *= stepsize; 204a7e14dcfSSatish Balay 205a7e14dcfSSatish Balay /* Adjust trust region size */ 206a7e14dcfSSatish Balay if (rhok < tron->eta2) { 207a7e14dcfSSatish Balay delta = PetscMin(xdiff,delta)*tron->sigma1; 208a7e14dcfSSatish Balay } else if (rhok > tron->eta4) { 209a7e14dcfSSatish Balay delta= PetscMin(xdiff,delta)*tron->sigma3; 210a7e14dcfSSatish Balay } else if (rhok > tron->eta3) { 211a7e14dcfSSatish Balay delta=PetscMin(xdiff,delta)*tron->sigma2; 212a7e14dcfSSatish Balay } 2139566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tron->G_New,tron->X_New, tao->XL, tao->XU, tao->gradient)); 2149566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tron->Free_Local)); 2159566063dSJacob Faibussowitsch PetscCall(VecWhichInactive(tao->XL,tron->X_New,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local)); 216a7e14dcfSSatish Balay f=f_new; 2179566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 2189566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->X_New, tao->solution)); 2199566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->G_New, tao->gradient)); 220a7e14dcfSSatish Balay break; 221a7e14dcfSSatish Balay } 222a7e14dcfSSatish Balay else if (delta <= 1e-30) { 223a7e14dcfSSatish Balay break; 224a7e14dcfSSatish Balay } 225a7e14dcfSSatish Balay else { 226a7e14dcfSSatish Balay delta /= 4.0; 227a7e14dcfSSatish Balay } 228a7e14dcfSSatish Balay } /* end linear solve loop */ 229a7e14dcfSSatish Balay 230a7e14dcfSSatish Balay tron->f=f; tron->actred=actred; tao->trust=delta; 2318931d482SJason Sarich tao->niter++; 2329566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its)); 2339566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,stepsize)); 2349566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 235a7e14dcfSSatish Balay } /* END MAIN LOOP */ 236a7e14dcfSSatish Balay PetscFunctionReturn(0); 237a7e14dcfSSatish Balay } 238a7e14dcfSSatish Balay 239441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao tao,TAO_TRON *tron) 240a7e14dcfSSatish Balay { 241a7e14dcfSSatish Balay PetscInt i; 242e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 243a7e14dcfSSatish Balay PetscReal actred=-1.0,actred_max=0.0; 244a7e14dcfSSatish Balay PetscReal f_new; 245a7e14dcfSSatish Balay /* 246a7e14dcfSSatish Balay The gradient and function value passed into and out of this 247a7e14dcfSSatish Balay routine should be current and correct. 248a7e14dcfSSatish Balay 249a7e14dcfSSatish Balay The free, active, and binding variables should be already identified 250a7e14dcfSSatish Balay */ 251a7e14dcfSSatish Balay PetscFunctionBegin; 252a7e14dcfSSatish Balay 253ce902467SBarry Smith for (i=0;i<tron->maxgpits;++i) { 254a7e14dcfSSatish Balay 255a7e14dcfSSatish Balay if (-actred <= (tron->pg_ftol)*actred_max) break; 256a7e14dcfSSatish Balay 257ce902467SBarry Smith ++tron->gp_iterates; 258ce902467SBarry Smith ++tron->total_gp_its; 259a7e14dcfSSatish Balay f_new=tron->f; 260a7e14dcfSSatish Balay 2619566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient,tao->stepdirection)); 2629566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection,-1.0)); 2639566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,tron->pgstepsize)); 264d0609cedSBarry Smith PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection,&tron->pgstepsize, &ls_reason)); 2659566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 266a7e14dcfSSatish Balay 2679566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient)); 2689566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 269ce902467SBarry Smith 270a7e14dcfSSatish Balay /* Update the iterate */ 271a7e14dcfSSatish Balay actred = f_new - tron->f; 272a7e14dcfSSatish Balay actred_max = PetscMax(actred_max,-(f_new - tron->f)); 273a7e14dcfSSatish Balay tron->f = f_new; 274a7e14dcfSSatish Balay } 275a7e14dcfSSatish Balay PetscFunctionReturn(0); 276a7e14dcfSSatish Balay } 277a7e14dcfSSatish Balay 278bdb10af2SPierre Jolivet static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU) 279bdb10af2SPierre Jolivet { 280a7e14dcfSSatish Balay 281a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 282a7e14dcfSSatish Balay 283a7e14dcfSSatish Balay PetscFunctionBegin; 284441846f8SBarry Smith PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 285a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 286a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 2873c859ba3SBarry Smith PetscCheck(tron->Work && tao->gradient,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist."); 288a7e14dcfSSatish Balay 2899566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tron->Work)); 2909566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->Work,DXL)); 2919566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXL,-1.0,tao->gradient)); 2929566063dSJacob Faibussowitsch PetscCall(VecSet(DXU,0.0)); 2939566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(DXL,DXL,DXU)); 294a7e14dcfSSatish Balay 2959566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient,DXU)); 2969566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU,-1.0,tron->Work)); 2979566063dSJacob Faibussowitsch PetscCall(VecSet(tron->Work,0.0)); 2989566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(DXU,tron->Work,DXU)); 299a7e14dcfSSatish Balay PetscFunctionReturn(0); 300a7e14dcfSSatish Balay } 301a7e14dcfSSatish Balay 302a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 3031522df2eSJason Sarich /*MC 3041522df2eSJason Sarich TAOTRON - The TRON algorithm is an active-set Newton trust region method 3051522df2eSJason Sarich for bound-constrained minimization. 3061522df2eSJason Sarich 3071522df2eSJason Sarich Options Database Keys: 3081522df2eSJason Sarich + -tao_tron_maxgpits - maximum number of gradient projections per TRON iterate 3091522df2eSJason Sarich - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets 3101522df2eSJason Sarich 3111eb8069cSJason Sarich Level: beginner 3121522df2eSJason Sarich M*/ 313728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_TRON(Tao tao) 314a7e14dcfSSatish Balay { 315a7e14dcfSSatish Balay TAO_TRON *tron; 3168caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 317a7e14dcfSSatish Balay 31847a47007SBarry Smith PetscFunctionBegin; 319a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_TRON; 320a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_TRON; 321a7e14dcfSSatish Balay tao->ops->view = TaoView_TRON; 322a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_TRON; 323a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_TRON; 324a7e14dcfSSatish Balay tao->ops->computedual = TaoComputeDual_TRON; 325a7e14dcfSSatish Balay 3269566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&tron)); 3276f4723b1SBarry Smith tao->data = (void*)tron; 3286552cf8aSJason Sarich 3296552cf8aSJason Sarich /* Override default settings (unless already changed) */ 3306552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 50; 3316552cf8aSJason Sarich if (!tao->trust0_changed) tao->trust0 = 1.0; 3324f1535bcSTodd Munson if (!tao->steptol_changed) tao->steptol = 0.0; 333a7e14dcfSSatish Balay 334a7e14dcfSSatish Balay /* Initialize pointers and variables */ 335a7e14dcfSSatish Balay tron->n = 0; 336a7e14dcfSSatish Balay tron->maxgpits = 3; 337a7e14dcfSSatish Balay tron->pg_ftol = 0.001; 338a7e14dcfSSatish Balay 339a7e14dcfSSatish Balay tron->eta1 = 1.0e-4; 340a7e14dcfSSatish Balay tron->eta2 = 0.25; 341a7e14dcfSSatish Balay tron->eta3 = 0.50; 342a7e14dcfSSatish Balay tron->eta4 = 0.90; 343a7e14dcfSSatish Balay 344a7e14dcfSSatish Balay tron->sigma1 = 0.5; 345a7e14dcfSSatish Balay tron->sigma2 = 2.0; 346a7e14dcfSSatish Balay tron->sigma3 = 4.0; 347a7e14dcfSSatish Balay 348a7e14dcfSSatish Balay tron->gp_iterates = 0; /* Cumulative number */ 349a7e14dcfSSatish Balay tron->total_gp_its = 0; 350a7e14dcfSSatish Balay tron->n_free = 0; 351a7e14dcfSSatish Balay 3526c23d075SBarry Smith tron->DXFree=NULL; 3536c23d075SBarry Smith tron->R=NULL; 3546c23d075SBarry Smith tron->X_New=NULL; 3556c23d075SBarry Smith tron->G_New=NULL; 3566c23d075SBarry Smith tron->Work=NULL; 3576c23d075SBarry Smith tron->Free_Local=NULL; 3586c23d075SBarry Smith tron->H_sub=NULL; 3596c23d075SBarry Smith tron->Hpre_sub=NULL; 360a7e14dcfSSatish Balay tao->subset_type = TAO_SUBSET_SUBVEC; 361a7e14dcfSSatish Balay 3629566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 3639566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 3649566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch,morethuente_type)); 3659566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch,tao)); 3669566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix)); 367a7e14dcfSSatish Balay 3689566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp)); 3699566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 3709566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix)); 3719566063dSJacob Faibussowitsch PetscCall(KSPSetType(tao->ksp,KSPSTCG)); 372a7e14dcfSSatish Balay PetscFunctionReturn(0); 373a7e14dcfSSatish Balay } 374