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)); 22a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 239566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 24a7e14dcfSSatish Balay PetscFunctionReturn(0); 25a7e14dcfSSatish Balay } 26a7e14dcfSSatish Balay 27a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 28*dbbe0bcdSBarry Smith static PetscErrorCode TaoSetFromOptions_TRON(Tao tao,PetscOptionItems *PetscOptionsObject) 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)); 111*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest ,tao->cnvP); 1123ecd9318SAlp Dener while (tao->reason==TAO_CONTINUE_ITERATING) { 113e1e80dc8SAlp Dener /* Call general purpose update function */ 114*dbbe0bcdSBarry Smith PetscTryTypeMethod(tao,update, tao->niter, tao->user_update); 115ce902467SBarry Smith 116ce902467SBarry Smith /* Perform projected gradient iterations */ 1179566063dSJacob Faibussowitsch PetscCall(TronGradientProjections(tao,tron)); 118ce902467SBarry Smith 1199566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient)); 1209566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 121ce902467SBarry Smith 122ce902467SBarry Smith tao->ksp_its=0; 123a7e14dcfSSatish Balay f=tron->f; delta=tao->trust; 124a7e14dcfSSatish Balay tron->n_free_last = tron->n_free; 1259566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre)); 126a7e14dcfSSatish Balay 127baa3a814SAlp Dener /* Generate index set (IS) of which bound constraints are active */ 1289566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tron->Free_Local)); 1299566063dSJacob Faibussowitsch PetscCall(VecWhichInactive(tao->XL,tao->solution,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local)); 1309566063dSJacob Faibussowitsch PetscCall(ISGetSize(tron->Free_Local, &tron->n_free)); 131a7e14dcfSSatish Balay 132a7e14dcfSSatish Balay /* If no free variables */ 133a7e14dcfSSatish Balay if (tron->n_free == 0) { 1349566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 1359566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its)); 1369566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,delta)); 137*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest ,tao->cnvP); 1383ecd9318SAlp Dener if (!tao->reason) { 1393ecd9318SAlp Dener tao->reason = TAO_CONVERGED_STEPTOL; 14055723ba5SJason Sarich } 141a7e14dcfSSatish Balay break; 142a7e14dcfSSatish Balay } 143a7e14dcfSSatish Balay /* use free_local to mask/submat gradient, hessian, stepdirection */ 1449566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->R)); 1459566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->DXFree)); 1469566063dSJacob Faibussowitsch PetscCall(VecSet(tron->DXFree,0.0)); 1479566063dSJacob Faibussowitsch PetscCall(VecScale(tron->R, -1.0)); 1489566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub)); 149a7e14dcfSSatish Balay if (tao->hessian == tao->hessian_pre) { 1509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tron->Hpre_sub)); 1519566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(tron->H_sub))); 152a7e14dcfSSatish Balay tron->Hpre_sub = tron->H_sub; 153a7e14dcfSSatish Balay } else { 1549566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type,&tron->Hpre_sub)); 155a7e14dcfSSatish Balay } 1569566063dSJacob Faibussowitsch PetscCall(KSPReset(tao->ksp)); 1579566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub)); 158a7e14dcfSSatish Balay while (1) { 159a7e14dcfSSatish Balay 160a7e14dcfSSatish Balay /* Approximately solve the reduced linear system */ 1619566063dSJacob Faibussowitsch PetscCall(KSPCGSetRadius(tao->ksp,delta)); 162d8ec8e84SJason Sarich 1639566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, tron->R, tron->DXFree)); 1649566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&its)); 165a7e14dcfSSatish Balay tao->ksp_its+=its; 166ae93cb3cSJason Sarich tao->ksp_tot_its+=its; 1679566063dSJacob Faibussowitsch PetscCall(VecSet(tao->stepdirection,0.0)); 168a7e14dcfSSatish Balay 169a7e14dcfSSatish Balay /* Add dxfree matrix to compute step direction vector */ 1709566063dSJacob Faibussowitsch PetscCall(VecISAXPY(tao->stepdirection,tron->Free_Local,1.0,tron->DXFree)); 171a7e14dcfSSatish Balay 1729566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx)); 1739566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, tron->X_New)); 1749566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, tron->G_New)); 175a7e14dcfSSatish Balay 176a7e14dcfSSatish Balay stepsize=1.0;f_new=f; 177a7e14dcfSSatish Balay 1789566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0)); 1799566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection,&stepsize,&ls_reason)); 1809566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 181a7e14dcfSSatish Balay 1829566063dSJacob Faibussowitsch PetscCall(MatMult(tao->hessian, tao->stepdirection, tron->Work)); 1839566063dSJacob Faibussowitsch PetscCall(VecAYPX(tron->Work, 0.5, tao->gradient)); 1849566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, tron->Work, &prered)); 185a7e14dcfSSatish Balay actred = f_new - f; 186ce902467SBarry Smith if ((PetscAbsScalar(actred) <= 1e-6) && (PetscAbsScalar(prered) <= 1e-6)) { 187ce902467SBarry Smith rhok = 1.0; 188ce902467SBarry Smith } else if (actred<0) { 189a7e14dcfSSatish Balay rhok=PetscAbs(-actred/prered); 190a7e14dcfSSatish Balay } else { 191a7e14dcfSSatish Balay rhok=0.0; 192a7e14dcfSSatish Balay } 193a7e14dcfSSatish Balay 194a7e14dcfSSatish Balay /* Compare actual improvement to the quadratic model */ 195a7e14dcfSSatish Balay if (rhok > tron->eta1) { /* Accept the point */ 196a7e14dcfSSatish Balay /* d = x_new - x */ 1979566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->X_New, tao->stepdirection)); 1989566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->stepdirection, -1.0, tao->solution)); 199a7e14dcfSSatish Balay 2009566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->stepdirection, NORM_2, &xdiff)); 201a7e14dcfSSatish Balay xdiff *= stepsize; 202a7e14dcfSSatish Balay 203a7e14dcfSSatish Balay /* Adjust trust region size */ 204a7e14dcfSSatish Balay if (rhok < tron->eta2) { 205a7e14dcfSSatish Balay delta = PetscMin(xdiff,delta)*tron->sigma1; 206a7e14dcfSSatish Balay } else if (rhok > tron->eta4) { 207a7e14dcfSSatish Balay delta= PetscMin(xdiff,delta)*tron->sigma3; 208a7e14dcfSSatish Balay } else if (rhok > tron->eta3) { 209a7e14dcfSSatish Balay delta=PetscMin(xdiff,delta)*tron->sigma2; 210a7e14dcfSSatish Balay } 2119566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tron->G_New,tron->X_New, tao->XL, tao->XU, tao->gradient)); 2129566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tron->Free_Local)); 2139566063dSJacob Faibussowitsch PetscCall(VecWhichInactive(tao->XL,tron->X_New,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local)); 214a7e14dcfSSatish Balay f=f_new; 2159566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 2169566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->X_New, tao->solution)); 2179566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->G_New, tao->gradient)); 218a7e14dcfSSatish Balay break; 219a7e14dcfSSatish Balay } 220a7e14dcfSSatish Balay else if (delta <= 1e-30) { 221a7e14dcfSSatish Balay break; 222a7e14dcfSSatish Balay } 223a7e14dcfSSatish Balay else { 224a7e14dcfSSatish Balay delta /= 4.0; 225a7e14dcfSSatish Balay } 226a7e14dcfSSatish Balay } /* end linear solve loop */ 227a7e14dcfSSatish Balay 228a7e14dcfSSatish Balay tron->f=f; tron->actred=actred; tao->trust=delta; 2298931d482SJason Sarich tao->niter++; 2309566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its)); 2319566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,stepsize)); 232*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest ,tao->cnvP); 233a7e14dcfSSatish Balay } /* END MAIN LOOP */ 234a7e14dcfSSatish Balay PetscFunctionReturn(0); 235a7e14dcfSSatish Balay } 236a7e14dcfSSatish Balay 237441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao tao,TAO_TRON *tron) 238a7e14dcfSSatish Balay { 239a7e14dcfSSatish Balay PetscInt i; 240e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 241a7e14dcfSSatish Balay PetscReal actred=-1.0,actred_max=0.0; 242a7e14dcfSSatish Balay PetscReal f_new; 243a7e14dcfSSatish Balay /* 244a7e14dcfSSatish Balay The gradient and function value passed into and out of this 245a7e14dcfSSatish Balay routine should be current and correct. 246a7e14dcfSSatish Balay 247a7e14dcfSSatish Balay The free, active, and binding variables should be already identified 248a7e14dcfSSatish Balay */ 249a7e14dcfSSatish Balay PetscFunctionBegin; 250a7e14dcfSSatish Balay 251ce902467SBarry Smith for (i=0;i<tron->maxgpits;++i) { 252a7e14dcfSSatish Balay 253a7e14dcfSSatish Balay if (-actred <= (tron->pg_ftol)*actred_max) break; 254a7e14dcfSSatish Balay 255ce902467SBarry Smith ++tron->gp_iterates; 256ce902467SBarry Smith ++tron->total_gp_its; 257a7e14dcfSSatish Balay f_new=tron->f; 258a7e14dcfSSatish Balay 2599566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient,tao->stepdirection)); 2609566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection,-1.0)); 2619566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,tron->pgstepsize)); 262d0609cedSBarry Smith PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection,&tron->pgstepsize, &ls_reason)); 2639566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao)); 264a7e14dcfSSatish Balay 2659566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient)); 2669566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_2,&tron->gnorm)); 267ce902467SBarry Smith 268a7e14dcfSSatish Balay /* Update the iterate */ 269a7e14dcfSSatish Balay actred = f_new - tron->f; 270a7e14dcfSSatish Balay actred_max = PetscMax(actred_max,-(f_new - tron->f)); 271a7e14dcfSSatish Balay tron->f = f_new; 272a7e14dcfSSatish Balay } 273a7e14dcfSSatish Balay PetscFunctionReturn(0); 274a7e14dcfSSatish Balay } 275a7e14dcfSSatish Balay 276bdb10af2SPierre Jolivet static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU) 277bdb10af2SPierre Jolivet { 278a7e14dcfSSatish Balay 279a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 280a7e14dcfSSatish Balay 281a7e14dcfSSatish Balay PetscFunctionBegin; 282441846f8SBarry Smith PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 283a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 284a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 2853c859ba3SBarry Smith PetscCheck(tron->Work && tao->gradient,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist."); 286a7e14dcfSSatish Balay 2879566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tron->Work)); 2889566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->Work,DXL)); 2899566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXL,-1.0,tao->gradient)); 2909566063dSJacob Faibussowitsch PetscCall(VecSet(DXU,0.0)); 2919566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(DXL,DXL,DXU)); 292a7e14dcfSSatish Balay 2939566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient,DXU)); 2949566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU,-1.0,tron->Work)); 2959566063dSJacob Faibussowitsch PetscCall(VecSet(tron->Work,0.0)); 2969566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(DXU,tron->Work,DXU)); 297a7e14dcfSSatish Balay PetscFunctionReturn(0); 298a7e14dcfSSatish Balay } 299a7e14dcfSSatish Balay 300a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 3011522df2eSJason Sarich /*MC 3021522df2eSJason Sarich TAOTRON - The TRON algorithm is an active-set Newton trust region method 3031522df2eSJason Sarich for bound-constrained minimization. 3041522df2eSJason Sarich 3051522df2eSJason Sarich Options Database Keys: 3061522df2eSJason Sarich + -tao_tron_maxgpits - maximum number of gradient projections per TRON iterate 3071522df2eSJason Sarich - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets 3081522df2eSJason Sarich 3091eb8069cSJason Sarich Level: beginner 3101522df2eSJason Sarich M*/ 311728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_TRON(Tao tao) 312a7e14dcfSSatish Balay { 313a7e14dcfSSatish Balay TAO_TRON *tron; 3148caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 315a7e14dcfSSatish Balay 31647a47007SBarry Smith PetscFunctionBegin; 317a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_TRON; 318a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_TRON; 319a7e14dcfSSatish Balay tao->ops->view = TaoView_TRON; 320a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_TRON; 321a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_TRON; 322a7e14dcfSSatish Balay tao->ops->computedual = TaoComputeDual_TRON; 323a7e14dcfSSatish Balay 3249566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&tron)); 3256f4723b1SBarry Smith tao->data = (void*)tron; 3266552cf8aSJason Sarich 3276552cf8aSJason Sarich /* Override default settings (unless already changed) */ 3286552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 50; 3296552cf8aSJason Sarich if (!tao->trust0_changed) tao->trust0 = 1.0; 3304f1535bcSTodd Munson if (!tao->steptol_changed) tao->steptol = 0.0; 331a7e14dcfSSatish Balay 332a7e14dcfSSatish Balay /* Initialize pointers and variables */ 333a7e14dcfSSatish Balay tron->n = 0; 334a7e14dcfSSatish Balay tron->maxgpits = 3; 335a7e14dcfSSatish Balay tron->pg_ftol = 0.001; 336a7e14dcfSSatish Balay 337a7e14dcfSSatish Balay tron->eta1 = 1.0e-4; 338a7e14dcfSSatish Balay tron->eta2 = 0.25; 339a7e14dcfSSatish Balay tron->eta3 = 0.50; 340a7e14dcfSSatish Balay tron->eta4 = 0.90; 341a7e14dcfSSatish Balay 342a7e14dcfSSatish Balay tron->sigma1 = 0.5; 343a7e14dcfSSatish Balay tron->sigma2 = 2.0; 344a7e14dcfSSatish Balay tron->sigma3 = 4.0; 345a7e14dcfSSatish Balay 346a7e14dcfSSatish Balay tron->gp_iterates = 0; /* Cumulative number */ 347a7e14dcfSSatish Balay tron->total_gp_its = 0; 348a7e14dcfSSatish Balay tron->n_free = 0; 349a7e14dcfSSatish Balay 3506c23d075SBarry Smith tron->DXFree=NULL; 3516c23d075SBarry Smith tron->R=NULL; 3526c23d075SBarry Smith tron->X_New=NULL; 3536c23d075SBarry Smith tron->G_New=NULL; 3546c23d075SBarry Smith tron->Work=NULL; 3556c23d075SBarry Smith tron->Free_Local=NULL; 3566c23d075SBarry Smith tron->H_sub=NULL; 3576c23d075SBarry Smith tron->Hpre_sub=NULL; 358a7e14dcfSSatish Balay tao->subset_type = TAO_SUBSET_SUBVEC; 359a7e14dcfSSatish Balay 3609566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 3619566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 3629566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch,morethuente_type)); 3639566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch,tao)); 3649566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix)); 365a7e14dcfSSatish Balay 3669566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp)); 3679566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 3689566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix)); 3699566063dSJacob Faibussowitsch PetscCall(KSPSetType(tao->ksp,KSPSTCG)); 370a7e14dcfSSatish Balay PetscFunctionReturn(0); 371a7e14dcfSSatish Balay } 372