1aaa7dc30SBarry Smith #include <../src/tao/bound/impls/tron/tron.h> 2aaa7dc30SBarry Smith #include <../src/tao/matrix/submatfree.h> 3a7e14dcfSSatish Balay 4a7e14dcfSSatish Balay 5a7e14dcfSSatish Balay /* TRON Routines */ 6441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao,TAO_TRON*); 7a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 8441846f8SBarry Smith static PetscErrorCode TaoDestroy_TRON(Tao tao) 9a7e14dcfSSatish Balay { 10a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 11a7e14dcfSSatish Balay PetscErrorCode ierr; 12a7e14dcfSSatish Balay 13a7e14dcfSSatish Balay PetscFunctionBegin; 14a7e14dcfSSatish Balay ierr = VecDestroy(&tron->X_New);CHKERRQ(ierr); 15a7e14dcfSSatish Balay ierr = VecDestroy(&tron->G_New);CHKERRQ(ierr); 16a7e14dcfSSatish Balay ierr = VecDestroy(&tron->Work);CHKERRQ(ierr); 17a7e14dcfSSatish Balay ierr = VecDestroy(&tron->DXFree);CHKERRQ(ierr); 18a7e14dcfSSatish Balay ierr = VecDestroy(&tron->R);CHKERRQ(ierr); 19a7e14dcfSSatish Balay ierr = VecDestroy(&tron->diag);CHKERRQ(ierr); 20a7e14dcfSSatish Balay ierr = VecScatterDestroy(&tron->scatter);CHKERRQ(ierr); 21a7e14dcfSSatish Balay ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 22a7e14dcfSSatish Balay ierr = MatDestroy(&tron->H_sub);CHKERRQ(ierr); 23a7e14dcfSSatish Balay ierr = MatDestroy(&tron->Hpre_sub);CHKERRQ(ierr); 24302440fdSBarry Smith ierr = PetscFree(tao->data);CHKERRQ(ierr); 25a7e14dcfSSatish Balay PetscFunctionReturn(0); 26a7e14dcfSSatish Balay } 27a7e14dcfSSatish Balay 28a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 294416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_TRON(PetscOptionItems *PetscOptionsObject,Tao tao) 30a7e14dcfSSatish Balay { 31a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 32a7e14dcfSSatish Balay PetscErrorCode ierr; 33a7e14dcfSSatish Balay PetscBool flg; 34a7e14dcfSSatish Balay 35a7e14dcfSSatish Balay PetscFunctionBegin; 361a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Newton Trust Region Method for bound constrained optimization");CHKERRQ(ierr); 371522df2eSJason Sarich ierr = PetscOptionsInt("-tao_tron_maxgpits","maximum number of gradient projections per TRON iterate","TaoSetMaxGPIts",tron->maxgpits,&tron->maxgpits,&flg);CHKERRQ(ierr); 38a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 39a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 40a7e14dcfSSatish Balay ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 41a7e14dcfSSatish Balay PetscFunctionReturn(0); 42a7e14dcfSSatish Balay } 43a7e14dcfSSatish Balay 44a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 45441846f8SBarry Smith static PetscErrorCode TaoView_TRON(Tao tao, PetscViewer viewer) 46a7e14dcfSSatish Balay { 47a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 48a7e14dcfSSatish Balay PetscBool isascii; 49a7e14dcfSSatish Balay PetscErrorCode ierr; 50a7e14dcfSSatish Balay 51a7e14dcfSSatish Balay PetscFunctionBegin; 52a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 53a7e14dcfSSatish Balay if (isascii) { 54a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"Total PG its: %D,",tron->total_gp_its);CHKERRQ(ierr); 5547a47007SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"PG tolerance: %g \n",(double)tron->pg_ftol);CHKERRQ(ierr); 56a7e14dcfSSatish Balay } 57a7e14dcfSSatish Balay PetscFunctionReturn(0); 58a7e14dcfSSatish Balay } 59a7e14dcfSSatish Balay 60a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 61441846f8SBarry Smith static PetscErrorCode TaoSetup_TRON(Tao tao) 62a7e14dcfSSatish Balay { 63a7e14dcfSSatish Balay PetscErrorCode ierr; 64a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 65a7e14dcfSSatish Balay 66a7e14dcfSSatish Balay PetscFunctionBegin; 67a7e14dcfSSatish Balay 68a7e14dcfSSatish Balay /* Allocate some arrays */ 69a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tron->diag);CHKERRQ(ierr); 70a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tron->X_New);CHKERRQ(ierr); 71a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tron->G_New);CHKERRQ(ierr); 72a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tron->Work);CHKERRQ(ierr); 73a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); 74a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr); 75a7e14dcfSSatish Balay if (!tao->XL) { 76a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->XL);CHKERRQ(ierr); 77e270355aSBarry Smith ierr = VecSet(tao->XL, PETSC_NINFINITY);CHKERRQ(ierr); 78a7e14dcfSSatish Balay } 79a7e14dcfSSatish Balay if (!tao->XU) { 80a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->XU);CHKERRQ(ierr); 81e270355aSBarry Smith ierr = VecSet(tao->XU, PETSC_INFINITY);CHKERRQ(ierr); 82a7e14dcfSSatish Balay } 83a7e14dcfSSatish Balay PetscFunctionReturn(0); 84a7e14dcfSSatish Balay } 85a7e14dcfSSatish Balay 86441846f8SBarry Smith static PetscErrorCode TaoSolve_TRON(Tao tao) 8747a47007SBarry Smith { 88a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 89a7e14dcfSSatish Balay PetscErrorCode ierr; 908931d482SJason Sarich PetscInt its; 91e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING; 92a7e14dcfSSatish Balay PetscReal prered,actred,delta,f,f_new,rhok,gdx,xdiff,stepsize; 9347a47007SBarry Smith 94a7e14dcfSSatish Balay PetscFunctionBegin; 95a7e14dcfSSatish Balay tron->pgstepsize = 1.0; 96a7e14dcfSSatish Balay tao->trust = tao->trust0; 97a7e14dcfSSatish Balay /* Project the current point onto the feasible set */ 98a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 99a7e14dcfSSatish Balay ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 100a7e14dcfSSatish Balay 101ce902467SBarry Smith /* Project the initial point onto the feasible region */ 102ce902467SBarry Smith ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 10353506e15SBarry Smith 104ce902467SBarry Smith /* Compute the objective function and gradient */ 105ce902467SBarry Smith ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&tron->f,tao->gradient);CHKERRQ(ierr); 106ce902467SBarry Smith ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 107ce902467SBarry Smith if (PetscIsInfOrNanReal(tron->f) || PetscIsInfOrNanReal(tron->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 108a7e14dcfSSatish Balay 109a7e14dcfSSatish Balay /* Project the gradient and calculate the norm */ 110a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 111a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 112a7e14dcfSSatish Balay 113ce902467SBarry Smith /* Initialize trust region radius */ 114ce902467SBarry Smith tao->trust=tao->trust0; 115a7e14dcfSSatish Balay if (tao->trust <= 0) { 116a7e14dcfSSatish Balay tao->trust=PetscMax(tron->gnorm*tron->gnorm,1.0); 117a7e14dcfSSatish Balay } 118a7e14dcfSSatish Balay 119ce902467SBarry Smith /* Initialize step sizes for the line searches */ 120ce902467SBarry Smith tron->pgstepsize=1.0; 121a7e14dcfSSatish Balay tron->stepsize=tao->trust; 122ce902467SBarry Smith 1233ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 1243ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 1253ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,tron->stepsize);CHKERRQ(ierr); 1263ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 1273ecd9318SAlp Dener while (tao->reason==TAO_CONTINUE_ITERATING){ 128ce902467SBarry Smith 129ce902467SBarry Smith /* Perform projected gradient iterations */ 130a7e14dcfSSatish Balay ierr = TronGradientProjections(tao,tron);CHKERRQ(ierr); 131ce902467SBarry Smith 132ce902467SBarry Smith ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 133ce902467SBarry Smith ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 134ce902467SBarry Smith 135ce902467SBarry Smith tao->ksp_its=0; 136a7e14dcfSSatish Balay f=tron->f; delta=tao->trust; 137a7e14dcfSSatish Balay tron->n_free_last = tron->n_free; 138ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 139a7e14dcfSSatish Balay 140baa3a814SAlp Dener /* Generate index set (IS) of which bound constraints are active */ 141ce902467SBarry Smith ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 142ce902467SBarry Smith ierr = VecWhichInactive(tao->XL,tao->solution,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local);CHKERRQ(ierr); 143a7e14dcfSSatish Balay ierr = ISGetSize(tron->Free_Local, &tron->n_free);CHKERRQ(ierr); 144a7e14dcfSSatish Balay 145a7e14dcfSSatish Balay /* If no free variables */ 146a7e14dcfSSatish Balay if (tron->n_free == 0) { 14755723ba5SJason Sarich ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 1483ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 1493ecd9318SAlp Dener ierr = TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,delta);CHKERRQ(ierr); 1503ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 1513ecd9318SAlp Dener if (!tao->reason) { 1523ecd9318SAlp Dener tao->reason = TAO_CONVERGED_STEPTOL; 15355723ba5SJason Sarich } 154a7e14dcfSSatish Balay break; 155a7e14dcfSSatish Balay } 156a7e14dcfSSatish Balay /* use free_local to mask/submat gradient, hessian, stepdirection */ 157b98f30f2SJason Sarich ierr = TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->R);CHKERRQ(ierr); 158b98f30f2SJason Sarich ierr = TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->DXFree);CHKERRQ(ierr); 159a7e14dcfSSatish Balay ierr = VecSet(tron->DXFree,0.0);CHKERRQ(ierr); 160a7e14dcfSSatish Balay ierr = VecScale(tron->R, -1.0);CHKERRQ(ierr); 161b98f30f2SJason Sarich ierr = TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub);CHKERRQ(ierr); 162a7e14dcfSSatish Balay if (tao->hessian == tao->hessian_pre) { 163a7e14dcfSSatish Balay ierr = MatDestroy(&tron->Hpre_sub);CHKERRQ(ierr); 164a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)(tron->H_sub));CHKERRQ(ierr); 165a7e14dcfSSatish Balay tron->Hpre_sub = tron->H_sub; 166a7e14dcfSSatish Balay } else { 167b98f30f2SJason Sarich ierr = TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type,&tron->Hpre_sub);CHKERRQ(ierr); 168a7e14dcfSSatish Balay } 169a7e14dcfSSatish Balay ierr = KSPReset(tao->ksp);CHKERRQ(ierr); 17023ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub);CHKERRQ(ierr); 171a7e14dcfSSatish Balay while (1) { 172a7e14dcfSSatish Balay 173a7e14dcfSSatish Balay /* Approximately solve the reduced linear system */ 1748f1e51d3STodd Munson ierr = KSPCGSetRadius(tao->ksp,delta);CHKERRQ(ierr); 175d8ec8e84SJason Sarich 176a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tron->R, tron->DXFree);CHKERRQ(ierr); 177a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 178a7e14dcfSSatish Balay tao->ksp_its+=its; 179ae93cb3cSJason Sarich tao->ksp_tot_its+=its; 180a7e14dcfSSatish Balay ierr = VecSet(tao->stepdirection,0.0);CHKERRQ(ierr); 181a7e14dcfSSatish Balay 182a7e14dcfSSatish Balay /* Add dxfree matrix to compute step direction vector */ 1834473680cSBarry Smith ierr = VecISAXPY(tao->stepdirection,tron->Free_Local,1.0,tron->DXFree);CHKERRQ(ierr); 184a7e14dcfSSatish Balay 185a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 186a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tron->X_New);CHKERRQ(ierr); 187a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tron->G_New);CHKERRQ(ierr); 188a7e14dcfSSatish Balay 189a7e14dcfSSatish Balay stepsize=1.0;f_new=f; 190a7e14dcfSSatish Balay 191a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 19247a47007SBarry Smith ierr = TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection,&stepsize,&ls_reason);CHKERRQ(ierr);CHKERRQ(ierr); 193a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 194a7e14dcfSSatish Balay 195a7e14dcfSSatish Balay ierr = MatMult(tao->hessian, tao->stepdirection, tron->Work);CHKERRQ(ierr); 196a7e14dcfSSatish Balay ierr = VecAYPX(tron->Work, 0.5, tao->gradient);CHKERRQ(ierr); 197a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tron->Work, &prered);CHKERRQ(ierr); 198a7e14dcfSSatish Balay actred = f_new - f; 199ce902467SBarry Smith if ((PetscAbsScalar(actred) <= 1e-6) && (PetscAbsScalar(prered) <= 1e-6)) { 200ce902467SBarry Smith rhok = 1.0; 201ce902467SBarry Smith } else if (actred<0) { 202a7e14dcfSSatish Balay rhok=PetscAbs(-actred/prered); 203a7e14dcfSSatish Balay } else { 204a7e14dcfSSatish Balay rhok=0.0; 205a7e14dcfSSatish Balay } 206a7e14dcfSSatish Balay 207a7e14dcfSSatish Balay /* Compare actual improvement to the quadratic model */ 208a7e14dcfSSatish Balay if (rhok > tron->eta1) { /* Accept the point */ 209a7e14dcfSSatish Balay /* d = x_new - x */ 210a7e14dcfSSatish Balay ierr = VecCopy(tron->X_New, tao->stepdirection);CHKERRQ(ierr); 211a7e14dcfSSatish Balay ierr = VecAXPY(tao->stepdirection, -1.0, tao->solution);CHKERRQ(ierr); 212a7e14dcfSSatish Balay 213a7e14dcfSSatish Balay ierr = VecNorm(tao->stepdirection, NORM_2, &xdiff);CHKERRQ(ierr); 214a7e14dcfSSatish Balay xdiff *= stepsize; 215a7e14dcfSSatish Balay 216a7e14dcfSSatish Balay /* Adjust trust region size */ 217a7e14dcfSSatish Balay if (rhok < tron->eta2 ){ 218a7e14dcfSSatish Balay delta = PetscMin(xdiff,delta)*tron->sigma1; 219a7e14dcfSSatish Balay } else if (rhok > tron->eta4 ){ 220a7e14dcfSSatish Balay delta= PetscMin(xdiff,delta)*tron->sigma3; 221a7e14dcfSSatish Balay } else if (rhok > tron->eta3 ){ 222a7e14dcfSSatish Balay delta=PetscMin(xdiff,delta)*tron->sigma2; 223a7e14dcfSSatish Balay } 224a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tron->G_New,tron->X_New, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr); 225a7e14dcfSSatish Balay ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 226ce902467SBarry Smith ierr = VecWhichInactive(tao->XL,tron->X_New,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local);CHKERRQ(ierr); 227a7e14dcfSSatish Balay f=f_new; 228a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 229a7e14dcfSSatish Balay ierr = VecCopy(tron->X_New, tao->solution);CHKERRQ(ierr); 230a7e14dcfSSatish Balay ierr = VecCopy(tron->G_New, tao->gradient);CHKERRQ(ierr); 231a7e14dcfSSatish Balay break; 232a7e14dcfSSatish Balay } 233a7e14dcfSSatish Balay else if (delta <= 1e-30) { 234a7e14dcfSSatish Balay break; 235a7e14dcfSSatish Balay } 236a7e14dcfSSatish Balay else { 237a7e14dcfSSatish Balay delta /= 4.0; 238a7e14dcfSSatish Balay } 239a7e14dcfSSatish Balay } /* end linear solve loop */ 240a7e14dcfSSatish Balay 241a7e14dcfSSatish Balay tron->f=f; tron->actred=actred; tao->trust=delta; 2428931d482SJason Sarich tao->niter++; 2433ecd9318SAlp Dener ierr = TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr); 244*770b7498SAlp Dener ierr = TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,stepsize);CHKERRQ(ierr); 2453ecd9318SAlp Dener ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 246a7e14dcfSSatish Balay } /* END MAIN LOOP */ 247a7e14dcfSSatish Balay PetscFunctionReturn(0); 248a7e14dcfSSatish Balay } 249a7e14dcfSSatish Balay 250441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao tao,TAO_TRON *tron) 251a7e14dcfSSatish Balay { 252a7e14dcfSSatish Balay PetscErrorCode ierr; 253a7e14dcfSSatish Balay PetscInt i; 254e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 255a7e14dcfSSatish Balay PetscReal actred=-1.0,actred_max=0.0; 256a7e14dcfSSatish Balay PetscReal f_new; 257a7e14dcfSSatish Balay /* 258a7e14dcfSSatish Balay The gradient and function value passed into and out of this 259a7e14dcfSSatish Balay routine should be current and correct. 260a7e14dcfSSatish Balay 261a7e14dcfSSatish Balay The free, active, and binding variables should be already identified 262a7e14dcfSSatish Balay */ 263a7e14dcfSSatish Balay PetscFunctionBegin; 264a7e14dcfSSatish Balay 265ce902467SBarry Smith for (i=0;i<tron->maxgpits;++i){ 266a7e14dcfSSatish Balay 267a7e14dcfSSatish Balay if (-actred <= (tron->pg_ftol)*actred_max) break; 268a7e14dcfSSatish Balay 269ce902467SBarry Smith ++tron->gp_iterates; 270ce902467SBarry Smith ++tron->total_gp_its; 271a7e14dcfSSatish Balay f_new=tron->f; 272a7e14dcfSSatish Balay 273a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient,tao->stepdirection);CHKERRQ(ierr); 274a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr); 275a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,tron->pgstepsize);CHKERRQ(ierr); 276a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection, 277a7e14dcfSSatish Balay &tron->pgstepsize, &ls_reason);CHKERRQ(ierr); 278a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 279a7e14dcfSSatish Balay 280ce902467SBarry Smith ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 281ce902467SBarry Smith ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 282ce902467SBarry Smith 283a7e14dcfSSatish Balay /* Update the iterate */ 284a7e14dcfSSatish Balay actred = f_new - tron->f; 285a7e14dcfSSatish Balay actred_max = PetscMax(actred_max,-(f_new - tron->f)); 286a7e14dcfSSatish Balay tron->f = f_new; 287a7e14dcfSSatish Balay } 288a7e14dcfSSatish Balay PetscFunctionReturn(0); 289a7e14dcfSSatish Balay } 290a7e14dcfSSatish Balay 291441846f8SBarry Smith static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU) { 292a7e14dcfSSatish Balay 293a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 294a7e14dcfSSatish Balay PetscErrorCode ierr; 295a7e14dcfSSatish Balay 296a7e14dcfSSatish Balay PetscFunctionBegin; 297441846f8SBarry Smith PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 298a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 299a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 30047a47007SBarry Smith if (!tron->Work || !tao->gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n"); 301a7e14dcfSSatish Balay 302a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tron->Work);CHKERRQ(ierr); 303a7e14dcfSSatish Balay ierr = VecCopy(tron->Work,DXL);CHKERRQ(ierr); 304a7e14dcfSSatish Balay ierr = VecAXPY(DXL,-1.0,tao->gradient);CHKERRQ(ierr); 305a7e14dcfSSatish Balay ierr = VecSet(DXU,0.0);CHKERRQ(ierr); 306a7e14dcfSSatish Balay ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr); 307a7e14dcfSSatish Balay 308a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient,DXU);CHKERRQ(ierr); 309a7e14dcfSSatish Balay ierr = VecAXPY(DXU,-1.0,tron->Work);CHKERRQ(ierr); 310a7e14dcfSSatish Balay ierr = VecSet(tron->Work,0.0);CHKERRQ(ierr); 311a7e14dcfSSatish Balay ierr = VecPointwiseMin(DXU,tron->Work,DXU);CHKERRQ(ierr); 312a7e14dcfSSatish Balay PetscFunctionReturn(0); 313a7e14dcfSSatish Balay } 314a7e14dcfSSatish Balay 315a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 3161522df2eSJason Sarich /*MC 3171522df2eSJason Sarich TAOTRON - The TRON algorithm is an active-set Newton trust region method 3181522df2eSJason Sarich for bound-constrained minimization. 3191522df2eSJason Sarich 3201522df2eSJason Sarich Options Database Keys: 3211522df2eSJason Sarich + -tao_tron_maxgpits - maximum number of gradient projections per TRON iterate 3221522df2eSJason Sarich - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets 3231522df2eSJason Sarich 3241eb8069cSJason Sarich Level: beginner 3251522df2eSJason Sarich M*/ 326728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_TRON(Tao tao) 327a7e14dcfSSatish Balay { 328a7e14dcfSSatish Balay TAO_TRON *tron; 329a7e14dcfSSatish Balay PetscErrorCode ierr; 3308caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 331a7e14dcfSSatish Balay 33247a47007SBarry Smith PetscFunctionBegin; 333a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_TRON; 334a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_TRON; 335a7e14dcfSSatish Balay tao->ops->view = TaoView_TRON; 336a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_TRON; 337a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_TRON; 338a7e14dcfSSatish Balay tao->ops->computedual = TaoComputeDual_TRON; 339a7e14dcfSSatish Balay 3403c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&tron);CHKERRQ(ierr); 3416f4723b1SBarry Smith tao->data = (void*)tron; 3426552cf8aSJason Sarich 3436552cf8aSJason Sarich /* Override default settings (unless already changed) */ 3446552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 50; 3456552cf8aSJason Sarich if (!tao->trust0_changed) tao->trust0 = 1.0; 3464f1535bcSTodd Munson if (!tao->steptol_changed) tao->steptol = 0.0; 347a7e14dcfSSatish Balay 348a7e14dcfSSatish Balay /* Initialize pointers and variables */ 349a7e14dcfSSatish Balay tron->n = 0; 350a7e14dcfSSatish Balay tron->maxgpits = 3; 351a7e14dcfSSatish Balay tron->pg_ftol = 0.001; 352a7e14dcfSSatish Balay 353a7e14dcfSSatish Balay tron->eta1 = 1.0e-4; 354a7e14dcfSSatish Balay tron->eta2 = 0.25; 355a7e14dcfSSatish Balay tron->eta3 = 0.50; 356a7e14dcfSSatish Balay tron->eta4 = 0.90; 357a7e14dcfSSatish Balay 358a7e14dcfSSatish Balay tron->sigma1 = 0.5; 359a7e14dcfSSatish Balay tron->sigma2 = 2.0; 360a7e14dcfSSatish Balay tron->sigma3 = 4.0; 361a7e14dcfSSatish Balay 362a7e14dcfSSatish Balay tron->gp_iterates = 0; /* Cumulative number */ 363a7e14dcfSSatish Balay tron->total_gp_its = 0; 364a7e14dcfSSatish Balay tron->n_free = 0; 365a7e14dcfSSatish Balay 3666c23d075SBarry Smith tron->DXFree=NULL; 3676c23d075SBarry Smith tron->R=NULL; 3686c23d075SBarry Smith tron->X_New=NULL; 3696c23d075SBarry Smith tron->G_New=NULL; 3706c23d075SBarry Smith tron->Work=NULL; 3716c23d075SBarry Smith tron->Free_Local=NULL; 3726c23d075SBarry Smith tron->H_sub=NULL; 3736c23d075SBarry Smith tron->Hpre_sub=NULL; 374a7e14dcfSSatish Balay tao->subset_type = TAO_SUBSET_SUBVEC; 375a7e14dcfSSatish Balay 376a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 37763b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 378a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 379441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 3805d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 381a7e14dcfSSatish Balay 382a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr); 38363b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 3845d527766SPatrick Farrell ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr); 3858f1e51d3STodd Munson ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr); 386a7e14dcfSSatish Balay PetscFunctionReturn(0); 387a7e14dcfSSatish Balay } 388