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 TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 92e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING; 93a7e14dcfSSatish Balay PetscReal prered,actred,delta,f,f_new,rhok,gdx,xdiff,stepsize; 9447a47007SBarry Smith 95a7e14dcfSSatish Balay PetscFunctionBegin; 96a7e14dcfSSatish Balay tron->pgstepsize = 1.0; 97a7e14dcfSSatish Balay tao->trust = tao->trust0; 98a7e14dcfSSatish Balay /* Project the current point onto the feasible set */ 99a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 100a7e14dcfSSatish Balay ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 101a7e14dcfSSatish Balay 102ce902467SBarry Smith /* Project the initial point onto the feasible region */ 103ce902467SBarry Smith ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 10453506e15SBarry Smith 105ce902467SBarry Smith /* Compute the objective function and gradient */ 106ce902467SBarry Smith ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&tron->f,tao->gradient);CHKERRQ(ierr); 107ce902467SBarry Smith ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 108ce902467SBarry Smith if (PetscIsInfOrNanReal(tron->f) || PetscIsInfOrNanReal(tron->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 109a7e14dcfSSatish Balay 110a7e14dcfSSatish Balay /* Project the gradient and calculate the norm */ 111a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 112a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 113a7e14dcfSSatish Balay 114ce902467SBarry Smith /* Initialize trust region radius */ 115ce902467SBarry Smith tao->trust=tao->trust0; 116a7e14dcfSSatish Balay if (tao->trust <= 0) { 117a7e14dcfSSatish Balay tao->trust=PetscMax(tron->gnorm*tron->gnorm,1.0); 118a7e14dcfSSatish Balay } 119a7e14dcfSSatish Balay 120ce902467SBarry Smith /* Initialize step sizes for the line searches */ 121ce902467SBarry Smith tron->pgstepsize=1.0; 122a7e14dcfSSatish Balay tron->stepsize=tao->trust; 123ce902467SBarry Smith 1248931d482SJason Sarich ierr = TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,tron->stepsize,&reason);CHKERRQ(ierr); 125a7e14dcfSSatish Balay while (reason==TAO_CONTINUE_ITERATING){ 126ce902467SBarry Smith 127ce902467SBarry Smith /* Perform projected gradient iterations */ 128a7e14dcfSSatish Balay ierr = TronGradientProjections(tao,tron);CHKERRQ(ierr); 129ce902467SBarry Smith 130ce902467SBarry Smith ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 131ce902467SBarry Smith ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 132ce902467SBarry Smith 133ce902467SBarry Smith tao->ksp_its=0; 134a7e14dcfSSatish Balay f=tron->f; delta=tao->trust; 135a7e14dcfSSatish Balay tron->n_free_last = tron->n_free; 136ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 137a7e14dcfSSatish Balay 138ce902467SBarry Smith ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 139ce902467SBarry Smith ierr = VecWhichInactive(tao->XL,tao->solution,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local);CHKERRQ(ierr); 140a7e14dcfSSatish Balay ierr = ISGetSize(tron->Free_Local, &tron->n_free);CHKERRQ(ierr); 141a7e14dcfSSatish Balay 142a7e14dcfSSatish Balay /* If no free variables */ 143a7e14dcfSSatish Balay if (tron->n_free == 0) { 14455723ba5SJason Sarich ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 1458931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, delta, &reason);CHKERRQ(ierr); 14655723ba5SJason Sarich if (!reason) { 14755723ba5SJason Sarich reason = TAO_CONVERGED_STEPTOL; 14855723ba5SJason Sarich ierr = TaoSetConvergedReason(tao,reason);CHKERRQ(ierr); 14955723ba5SJason Sarich } 150a7e14dcfSSatish Balay break; 151a7e14dcfSSatish Balay } 152a7e14dcfSSatish Balay /* use free_local to mask/submat gradient, hessian, stepdirection */ 153b98f30f2SJason Sarich ierr = TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->R);CHKERRQ(ierr); 154b98f30f2SJason Sarich ierr = TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->DXFree);CHKERRQ(ierr); 155a7e14dcfSSatish Balay ierr = VecSet(tron->DXFree,0.0);CHKERRQ(ierr); 156a7e14dcfSSatish Balay ierr = VecScale(tron->R, -1.0);CHKERRQ(ierr); 157b98f30f2SJason Sarich ierr = TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub);CHKERRQ(ierr); 158a7e14dcfSSatish Balay if (tao->hessian == tao->hessian_pre) { 159a7e14dcfSSatish Balay ierr = MatDestroy(&tron->Hpre_sub);CHKERRQ(ierr); 160a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)(tron->H_sub));CHKERRQ(ierr); 161a7e14dcfSSatish Balay tron->Hpre_sub = tron->H_sub; 162a7e14dcfSSatish Balay } else { 163b98f30f2SJason Sarich ierr = TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type,&tron->Hpre_sub);CHKERRQ(ierr); 164a7e14dcfSSatish Balay } 165a7e14dcfSSatish Balay ierr = KSPReset(tao->ksp);CHKERRQ(ierr); 16623ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub);CHKERRQ(ierr); 167a7e14dcfSSatish Balay while (1) { 168a7e14dcfSSatish Balay 169a7e14dcfSSatish Balay /* Approximately solve the reduced linear system */ 1708f1e51d3STodd Munson ierr = KSPCGSetRadius(tao->ksp,delta);CHKERRQ(ierr); 171d8ec8e84SJason Sarich 172a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tron->R, tron->DXFree);CHKERRQ(ierr); 173a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 174a7e14dcfSSatish Balay tao->ksp_its+=its; 175ae93cb3cSJason Sarich tao->ksp_tot_its+=its; 176a7e14dcfSSatish Balay ierr = VecSet(tao->stepdirection,0.0);CHKERRQ(ierr); 177a7e14dcfSSatish Balay 178a7e14dcfSSatish Balay /* Add dxfree matrix to compute step direction vector */ 1794473680cSBarry Smith ierr = VecISAXPY(tao->stepdirection,tron->Free_Local,1.0,tron->DXFree);CHKERRQ(ierr); 180a7e14dcfSSatish Balay 181a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 182a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tron->X_New);CHKERRQ(ierr); 183a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tron->G_New);CHKERRQ(ierr); 184a7e14dcfSSatish Balay 185a7e14dcfSSatish Balay stepsize=1.0;f_new=f; 186a7e14dcfSSatish Balay 187a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 18847a47007SBarry Smith ierr = TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection,&stepsize,&ls_reason);CHKERRQ(ierr);CHKERRQ(ierr); 189a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 190a7e14dcfSSatish Balay 191a7e14dcfSSatish Balay ierr = MatMult(tao->hessian, tao->stepdirection, tron->Work);CHKERRQ(ierr); 192a7e14dcfSSatish Balay ierr = VecAYPX(tron->Work, 0.5, tao->gradient);CHKERRQ(ierr); 193a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tron->Work, &prered);CHKERRQ(ierr); 194a7e14dcfSSatish Balay actred = f_new - f; 195ce902467SBarry Smith if ((PetscAbsScalar(actred) <= 1e-6) && (PetscAbsScalar(prered) <= 1e-6)) { 196ce902467SBarry Smith rhok = 1.0; 197ce902467SBarry Smith } else if (actred<0) { 198a7e14dcfSSatish Balay rhok=PetscAbs(-actred/prered); 199a7e14dcfSSatish Balay } else { 200a7e14dcfSSatish Balay rhok=0.0; 201a7e14dcfSSatish Balay } 202a7e14dcfSSatish Balay 203a7e14dcfSSatish Balay /* Compare actual improvement to the quadratic model */ 204a7e14dcfSSatish Balay if (rhok > tron->eta1) { /* Accept the point */ 205a7e14dcfSSatish Balay /* d = x_new - x */ 206a7e14dcfSSatish Balay ierr = VecCopy(tron->X_New, tao->stepdirection);CHKERRQ(ierr); 207a7e14dcfSSatish Balay ierr = VecAXPY(tao->stepdirection, -1.0, tao->solution);CHKERRQ(ierr); 208a7e14dcfSSatish Balay 209a7e14dcfSSatish Balay ierr = VecNorm(tao->stepdirection, NORM_2, &xdiff);CHKERRQ(ierr); 210a7e14dcfSSatish Balay xdiff *= stepsize; 211a7e14dcfSSatish Balay 212a7e14dcfSSatish Balay /* Adjust trust region size */ 213a7e14dcfSSatish Balay if (rhok < tron->eta2 ){ 214a7e14dcfSSatish Balay delta = PetscMin(xdiff,delta)*tron->sigma1; 215a7e14dcfSSatish Balay } else if (rhok > tron->eta4 ){ 216a7e14dcfSSatish Balay delta= PetscMin(xdiff,delta)*tron->sigma3; 217a7e14dcfSSatish Balay } else if (rhok > tron->eta3 ){ 218a7e14dcfSSatish Balay delta=PetscMin(xdiff,delta)*tron->sigma2; 219a7e14dcfSSatish Balay } 220a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tron->G_New,tron->X_New, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr); 221a7e14dcfSSatish Balay ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 222ce902467SBarry Smith ierr = VecWhichInactive(tao->XL,tron->X_New,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local);CHKERRQ(ierr); 223a7e14dcfSSatish Balay f=f_new; 224a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 225a7e14dcfSSatish Balay ierr = VecCopy(tron->X_New, tao->solution);CHKERRQ(ierr); 226a7e14dcfSSatish Balay ierr = VecCopy(tron->G_New, tao->gradient);CHKERRQ(ierr); 227a7e14dcfSSatish Balay break; 228a7e14dcfSSatish Balay } 229a7e14dcfSSatish Balay else if (delta <= 1e-30) { 230a7e14dcfSSatish Balay break; 231a7e14dcfSSatish Balay } 232a7e14dcfSSatish Balay else { 233a7e14dcfSSatish Balay delta /= 4.0; 234a7e14dcfSSatish Balay } 235a7e14dcfSSatish Balay } /* end linear solve loop */ 236a7e14dcfSSatish Balay 237a7e14dcfSSatish Balay tron->f=f; tron->actred=actred; tao->trust=delta; 2388931d482SJason Sarich tao->niter++; 2398931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, delta, &reason);CHKERRQ(ierr); 240a7e14dcfSSatish Balay } /* END MAIN LOOP */ 241a7e14dcfSSatish Balay PetscFunctionReturn(0); 242a7e14dcfSSatish Balay } 243a7e14dcfSSatish Balay 244441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao tao,TAO_TRON *tron) 245a7e14dcfSSatish Balay { 246a7e14dcfSSatish Balay PetscErrorCode ierr; 247a7e14dcfSSatish Balay PetscInt i; 248e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 249a7e14dcfSSatish Balay PetscReal actred=-1.0,actred_max=0.0; 250a7e14dcfSSatish Balay PetscReal f_new; 251a7e14dcfSSatish Balay /* 252a7e14dcfSSatish Balay The gradient and function value passed into and out of this 253a7e14dcfSSatish Balay routine should be current and correct. 254a7e14dcfSSatish Balay 255a7e14dcfSSatish Balay The free, active, and binding variables should be already identified 256a7e14dcfSSatish Balay */ 257a7e14dcfSSatish Balay PetscFunctionBegin; 258a7e14dcfSSatish Balay 259ce902467SBarry Smith for (i=0;i<tron->maxgpits;++i){ 260a7e14dcfSSatish Balay 261a7e14dcfSSatish Balay if (-actred <= (tron->pg_ftol)*actred_max) break; 262a7e14dcfSSatish Balay 263ce902467SBarry Smith ++tron->gp_iterates; 264ce902467SBarry Smith ++tron->total_gp_its; 265a7e14dcfSSatish Balay f_new=tron->f; 266a7e14dcfSSatish Balay 267a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient,tao->stepdirection);CHKERRQ(ierr); 268a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr); 269a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,tron->pgstepsize);CHKERRQ(ierr); 270a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection, 271a7e14dcfSSatish Balay &tron->pgstepsize, &ls_reason);CHKERRQ(ierr); 272a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 273a7e14dcfSSatish Balay 274ce902467SBarry Smith ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 275ce902467SBarry Smith ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 276ce902467SBarry Smith 277a7e14dcfSSatish Balay /* Update the iterate */ 278a7e14dcfSSatish Balay actred = f_new - tron->f; 279a7e14dcfSSatish Balay actred_max = PetscMax(actred_max,-(f_new - tron->f)); 280a7e14dcfSSatish Balay tron->f = f_new; 281a7e14dcfSSatish Balay } 282a7e14dcfSSatish Balay PetscFunctionReturn(0); 283a7e14dcfSSatish Balay } 284a7e14dcfSSatish Balay 285441846f8SBarry Smith static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU) { 286a7e14dcfSSatish Balay 287a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 288a7e14dcfSSatish Balay PetscErrorCode ierr; 289a7e14dcfSSatish Balay 290a7e14dcfSSatish Balay PetscFunctionBegin; 291441846f8SBarry Smith PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 292a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 293a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 29447a47007SBarry Smith if (!tron->Work || !tao->gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n"); 295a7e14dcfSSatish Balay 296a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tron->Work);CHKERRQ(ierr); 297a7e14dcfSSatish Balay ierr = VecCopy(tron->Work,DXL);CHKERRQ(ierr); 298a7e14dcfSSatish Balay ierr = VecAXPY(DXL,-1.0,tao->gradient);CHKERRQ(ierr); 299a7e14dcfSSatish Balay ierr = VecSet(DXU,0.0);CHKERRQ(ierr); 300a7e14dcfSSatish Balay ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr); 301a7e14dcfSSatish Balay 302a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient,DXU);CHKERRQ(ierr); 303a7e14dcfSSatish Balay ierr = VecAXPY(DXU,-1.0,tron->Work);CHKERRQ(ierr); 304a7e14dcfSSatish Balay ierr = VecSet(tron->Work,0.0);CHKERRQ(ierr); 305a7e14dcfSSatish Balay ierr = VecPointwiseMin(DXU,tron->Work,DXU);CHKERRQ(ierr); 306a7e14dcfSSatish Balay PetscFunctionReturn(0); 307a7e14dcfSSatish Balay } 308a7e14dcfSSatish Balay 309a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 3101522df2eSJason Sarich /*MC 3111522df2eSJason Sarich TAOTRON - The TRON algorithm is an active-set Newton trust region method 3121522df2eSJason Sarich for bound-constrained minimization. 3131522df2eSJason Sarich 3141522df2eSJason Sarich Options Database Keys: 3151522df2eSJason Sarich + -tao_tron_maxgpits - maximum number of gradient projections per TRON iterate 3161522df2eSJason Sarich - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets 3171522df2eSJason Sarich 3181eb8069cSJason Sarich Level: beginner 3191522df2eSJason Sarich M*/ 320728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_TRON(Tao tao) 321a7e14dcfSSatish Balay { 322a7e14dcfSSatish Balay TAO_TRON *tron; 323a7e14dcfSSatish Balay PetscErrorCode ierr; 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 3343c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&tron);CHKERRQ(ierr); 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 370a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 371*63b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 372a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 373441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 3745d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 375a7e14dcfSSatish Balay 376a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr); 377*63b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 3785d527766SPatrick Farrell ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr); 3798f1e51d3STodd Munson ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr); 380a7e14dcfSSatish Balay PetscFunctionReturn(0); 381a7e14dcfSSatish Balay } 382