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 = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 55a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"Total PG its: %D,",tron->total_gp_its);CHKERRQ(ierr); 5647a47007SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"PG tolerance: %g \n",(double)tron->pg_ftol);CHKERRQ(ierr); 57a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 58a7e14dcfSSatish Balay } 59a7e14dcfSSatish Balay PetscFunctionReturn(0); 60a7e14dcfSSatish Balay } 61a7e14dcfSSatish Balay 62a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 63441846f8SBarry Smith static PetscErrorCode TaoSetup_TRON(Tao tao) 64a7e14dcfSSatish Balay { 65a7e14dcfSSatish Balay PetscErrorCode ierr; 66a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 67a7e14dcfSSatish Balay 68a7e14dcfSSatish Balay PetscFunctionBegin; 69a7e14dcfSSatish Balay 70a7e14dcfSSatish Balay /* Allocate some arrays */ 71a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tron->diag);CHKERRQ(ierr); 72a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tron->X_New);CHKERRQ(ierr); 73a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tron->G_New);CHKERRQ(ierr); 74a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tron->Work);CHKERRQ(ierr); 75a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); 76a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr); 77a7e14dcfSSatish Balay if (!tao->XL) { 78a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->XL);CHKERRQ(ierr); 79e270355aSBarry Smith ierr = VecSet(tao->XL, PETSC_NINFINITY);CHKERRQ(ierr); 80a7e14dcfSSatish Balay } 81a7e14dcfSSatish Balay if (!tao->XU) { 82a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->XU);CHKERRQ(ierr); 83e270355aSBarry Smith ierr = VecSet(tao->XU, PETSC_INFINITY);CHKERRQ(ierr); 84a7e14dcfSSatish Balay } 85a7e14dcfSSatish Balay PetscFunctionReturn(0); 86a7e14dcfSSatish Balay } 87a7e14dcfSSatish Balay 88441846f8SBarry Smith static PetscErrorCode TaoSolve_TRON(Tao tao) 8947a47007SBarry Smith { 90a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 91a7e14dcfSSatish Balay PetscErrorCode ierr; 928931d482SJason Sarich PetscInt its; 93e4cb33bbSBarry Smith TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 94e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING; 95a7e14dcfSSatish Balay PetscReal prered,actred,delta,f,f_new,rhok,gdx,xdiff,stepsize; 9647a47007SBarry Smith 97a7e14dcfSSatish Balay PetscFunctionBegin; 98a7e14dcfSSatish Balay tron->pgstepsize = 1.0; 99a7e14dcfSSatish Balay tao->trust = tao->trust0; 100a7e14dcfSSatish Balay /* Project the current point onto the feasible set */ 101a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 102a7e14dcfSSatish Balay ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 103a7e14dcfSSatish Balay 104*ce902467SBarry Smith /* Project the initial point onto the feasible region */ 105*ce902467SBarry Smith ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 10653506e15SBarry Smith 107*ce902467SBarry Smith /* Compute the objective function and gradient */ 108*ce902467SBarry Smith ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&tron->f,tao->gradient);CHKERRQ(ierr); 109*ce902467SBarry Smith ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 110*ce902467SBarry Smith if (PetscIsInfOrNanReal(tron->f) || PetscIsInfOrNanReal(tron->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 111a7e14dcfSSatish Balay 112a7e14dcfSSatish Balay /* Project the gradient and calculate the norm */ 113a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 114a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 115a7e14dcfSSatish Balay 116*ce902467SBarry Smith /* Initialize trust region radius */ 117*ce902467SBarry Smith tao->trust=tao->trust0; 118a7e14dcfSSatish Balay if (tao->trust <= 0) { 119a7e14dcfSSatish Balay tao->trust=PetscMax(tron->gnorm*tron->gnorm,1.0); 120a7e14dcfSSatish Balay } 121a7e14dcfSSatish Balay 122*ce902467SBarry Smith /* Initialize step sizes for the line searches */ 123*ce902467SBarry Smith tron->pgstepsize=1.0; 124a7e14dcfSSatish Balay tron->stepsize=tao->trust; 125*ce902467SBarry Smith 1268931d482SJason Sarich ierr = TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,tron->stepsize,&reason);CHKERRQ(ierr); 127a7e14dcfSSatish Balay while (reason==TAO_CONTINUE_ITERATING){ 128*ce902467SBarry Smith 129*ce902467SBarry Smith /* Perform projected gradient iterations */ 130a7e14dcfSSatish Balay ierr = TronGradientProjections(tao,tron);CHKERRQ(ierr); 131*ce902467SBarry Smith 132*ce902467SBarry Smith ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 133*ce902467SBarry Smith ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 134*ce902467SBarry Smith 135*ce902467SBarry 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 140*ce902467SBarry Smith ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 141*ce902467SBarry Smith ierr = VecWhichInactive(tao->XL,tao->solution,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local);CHKERRQ(ierr); 142a7e14dcfSSatish Balay ierr = ISGetSize(tron->Free_Local, &tron->n_free);CHKERRQ(ierr); 143a7e14dcfSSatish Balay 144a7e14dcfSSatish Balay /* If no free variables */ 145a7e14dcfSSatish Balay if (tron->n_free == 0) { 14655723ba5SJason Sarich ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 1478931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, delta, &reason);CHKERRQ(ierr); 14855723ba5SJason Sarich if (!reason) { 14955723ba5SJason Sarich reason = TAO_CONVERGED_STEPTOL; 15055723ba5SJason Sarich ierr = TaoSetConvergedReason(tao,reason);CHKERRQ(ierr); 15155723ba5SJason Sarich } 152a7e14dcfSSatish Balay break; 153a7e14dcfSSatish Balay } 154a7e14dcfSSatish Balay /* use free_local to mask/submat gradient, hessian, stepdirection */ 155b98f30f2SJason Sarich ierr = TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->R);CHKERRQ(ierr); 156b98f30f2SJason Sarich ierr = TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->DXFree);CHKERRQ(ierr); 157a7e14dcfSSatish Balay ierr = VecSet(tron->DXFree,0.0);CHKERRQ(ierr); 158a7e14dcfSSatish Balay ierr = VecScale(tron->R, -1.0);CHKERRQ(ierr); 159b98f30f2SJason Sarich ierr = TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub);CHKERRQ(ierr); 160a7e14dcfSSatish Balay if (tao->hessian == tao->hessian_pre) { 161a7e14dcfSSatish Balay ierr = MatDestroy(&tron->Hpre_sub);CHKERRQ(ierr); 162a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)(tron->H_sub));CHKERRQ(ierr); 163a7e14dcfSSatish Balay tron->Hpre_sub = tron->H_sub; 164a7e14dcfSSatish Balay } else { 165b98f30f2SJason Sarich ierr = TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type,&tron->Hpre_sub);CHKERRQ(ierr); 166a7e14dcfSSatish Balay } 167a7e14dcfSSatish Balay ierr = KSPReset(tao->ksp);CHKERRQ(ierr); 16823ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub);CHKERRQ(ierr); 169a7e14dcfSSatish Balay while (1) { 170a7e14dcfSSatish Balay 171a7e14dcfSSatish Balay /* Approximately solve the reduced linear system */ 1728f1e51d3STodd Munson ierr = KSPCGSetRadius(tao->ksp,delta);CHKERRQ(ierr); 173d8ec8e84SJason Sarich 174a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tron->R, tron->DXFree);CHKERRQ(ierr); 175a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 176a7e14dcfSSatish Balay tao->ksp_its+=its; 177ae93cb3cSJason Sarich tao->ksp_tot_its+=its; 178a7e14dcfSSatish Balay ierr = VecSet(tao->stepdirection,0.0);CHKERRQ(ierr); 179a7e14dcfSSatish Balay 180a7e14dcfSSatish Balay /* Add dxfree matrix to compute step direction vector */ 1814473680cSBarry Smith ierr = VecISAXPY(tao->stepdirection,tron->Free_Local,1.0,tron->DXFree);CHKERRQ(ierr); 182a7e14dcfSSatish Balay 183a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 184a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tron->X_New);CHKERRQ(ierr); 185a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tron->G_New);CHKERRQ(ierr); 186a7e14dcfSSatish Balay 187a7e14dcfSSatish Balay stepsize=1.0;f_new=f; 188a7e14dcfSSatish Balay 189a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 19047a47007SBarry Smith ierr = TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection,&stepsize,&ls_reason);CHKERRQ(ierr);CHKERRQ(ierr); 191a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 192a7e14dcfSSatish Balay 193a7e14dcfSSatish Balay ierr = MatMult(tao->hessian, tao->stepdirection, tron->Work);CHKERRQ(ierr); 194a7e14dcfSSatish Balay ierr = VecAYPX(tron->Work, 0.5, tao->gradient);CHKERRQ(ierr); 195a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tron->Work, &prered);CHKERRQ(ierr); 196a7e14dcfSSatish Balay actred = f_new - f; 197*ce902467SBarry Smith if ((PetscAbsScalar(actred) <= 1e-6) && (PetscAbsScalar(prered) <= 1e-6)) { 198*ce902467SBarry Smith rhok = 1.0; 199*ce902467SBarry Smith } else if (actred<0) { 200a7e14dcfSSatish Balay rhok=PetscAbs(-actred/prered); 201a7e14dcfSSatish Balay } else { 202a7e14dcfSSatish Balay rhok=0.0; 203a7e14dcfSSatish Balay } 204a7e14dcfSSatish Balay 205a7e14dcfSSatish Balay /* Compare actual improvement to the quadratic model */ 206a7e14dcfSSatish Balay if (rhok > tron->eta1) { /* Accept the point */ 207a7e14dcfSSatish Balay /* d = x_new - x */ 208a7e14dcfSSatish Balay ierr = VecCopy(tron->X_New, tao->stepdirection);CHKERRQ(ierr); 209a7e14dcfSSatish Balay ierr = VecAXPY(tao->stepdirection, -1.0, tao->solution);CHKERRQ(ierr); 210a7e14dcfSSatish Balay 211a7e14dcfSSatish Balay ierr = VecNorm(tao->stepdirection, NORM_2, &xdiff);CHKERRQ(ierr); 212a7e14dcfSSatish Balay xdiff *= stepsize; 213a7e14dcfSSatish Balay 214a7e14dcfSSatish Balay /* Adjust trust region size */ 215a7e14dcfSSatish Balay if (rhok < tron->eta2 ){ 216a7e14dcfSSatish Balay delta = PetscMin(xdiff,delta)*tron->sigma1; 217a7e14dcfSSatish Balay } else if (rhok > tron->eta4 ){ 218a7e14dcfSSatish Balay delta= PetscMin(xdiff,delta)*tron->sigma3; 219a7e14dcfSSatish Balay } else if (rhok > tron->eta3 ){ 220a7e14dcfSSatish Balay delta=PetscMin(xdiff,delta)*tron->sigma2; 221a7e14dcfSSatish Balay } 222a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tron->G_New,tron->X_New, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr); 223a7e14dcfSSatish Balay ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 224*ce902467SBarry Smith ierr = VecWhichInactive(tao->XL,tron->X_New,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local);CHKERRQ(ierr); 225a7e14dcfSSatish Balay f=f_new; 226a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 227a7e14dcfSSatish Balay ierr = VecCopy(tron->X_New, tao->solution);CHKERRQ(ierr); 228a7e14dcfSSatish Balay ierr = VecCopy(tron->G_New, tao->gradient);CHKERRQ(ierr); 229a7e14dcfSSatish Balay break; 230a7e14dcfSSatish Balay } 231a7e14dcfSSatish Balay else if (delta <= 1e-30) { 232a7e14dcfSSatish Balay break; 233a7e14dcfSSatish Balay } 234a7e14dcfSSatish Balay else { 235a7e14dcfSSatish Balay delta /= 4.0; 236a7e14dcfSSatish Balay } 237a7e14dcfSSatish Balay } /* end linear solve loop */ 238a7e14dcfSSatish Balay 239a7e14dcfSSatish Balay tron->f=f; tron->actred=actred; tao->trust=delta; 2408931d482SJason Sarich tao->niter++; 2418931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, delta, &reason);CHKERRQ(ierr); 242a7e14dcfSSatish Balay } /* END MAIN LOOP */ 243a7e14dcfSSatish Balay PetscFunctionReturn(0); 244a7e14dcfSSatish Balay } 245a7e14dcfSSatish Balay 246441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao tao,TAO_TRON *tron) 247a7e14dcfSSatish Balay { 248a7e14dcfSSatish Balay PetscErrorCode ierr; 249a7e14dcfSSatish Balay PetscInt i; 250e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 251a7e14dcfSSatish Balay PetscReal actred=-1.0,actred_max=0.0; 252a7e14dcfSSatish Balay PetscReal f_new; 253a7e14dcfSSatish Balay /* 254a7e14dcfSSatish Balay The gradient and function value passed into and out of this 255a7e14dcfSSatish Balay routine should be current and correct. 256a7e14dcfSSatish Balay 257a7e14dcfSSatish Balay The free, active, and binding variables should be already identified 258a7e14dcfSSatish Balay */ 259a7e14dcfSSatish Balay PetscFunctionBegin; 260a7e14dcfSSatish Balay 261*ce902467SBarry Smith for (i=0;i<tron->maxgpits;++i){ 262a7e14dcfSSatish Balay 263a7e14dcfSSatish Balay if (-actred <= (tron->pg_ftol)*actred_max) break; 264a7e14dcfSSatish Balay 265*ce902467SBarry Smith ++tron->gp_iterates; 266*ce902467SBarry Smith ++tron->total_gp_its; 267a7e14dcfSSatish Balay f_new=tron->f; 268a7e14dcfSSatish Balay 269a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient,tao->stepdirection);CHKERRQ(ierr); 270a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr); 271a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,tron->pgstepsize);CHKERRQ(ierr); 272a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection, 273a7e14dcfSSatish Balay &tron->pgstepsize, &ls_reason);CHKERRQ(ierr); 274a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 275a7e14dcfSSatish Balay 276*ce902467SBarry Smith ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr); 277*ce902467SBarry Smith ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 278*ce902467SBarry Smith 279a7e14dcfSSatish Balay /* Update the iterate */ 280a7e14dcfSSatish Balay actred = f_new - tron->f; 281a7e14dcfSSatish Balay actred_max = PetscMax(actred_max,-(f_new - tron->f)); 282a7e14dcfSSatish Balay tron->f = f_new; 283a7e14dcfSSatish Balay } 284a7e14dcfSSatish Balay PetscFunctionReturn(0); 285a7e14dcfSSatish Balay } 286a7e14dcfSSatish Balay 287441846f8SBarry Smith static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU) { 288a7e14dcfSSatish Balay 289a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 290a7e14dcfSSatish Balay PetscErrorCode ierr; 291a7e14dcfSSatish Balay 292a7e14dcfSSatish Balay PetscFunctionBegin; 293441846f8SBarry Smith PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 294a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 295a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 29647a47007SBarry Smith if (!tron->Work || !tao->gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n"); 297a7e14dcfSSatish Balay 298a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tron->Work);CHKERRQ(ierr); 299a7e14dcfSSatish Balay ierr = VecCopy(tron->Work,DXL);CHKERRQ(ierr); 300a7e14dcfSSatish Balay ierr = VecAXPY(DXL,-1.0,tao->gradient);CHKERRQ(ierr); 301a7e14dcfSSatish Balay ierr = VecSet(DXU,0.0);CHKERRQ(ierr); 302a7e14dcfSSatish Balay ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr); 303a7e14dcfSSatish Balay 304a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient,DXU);CHKERRQ(ierr); 305a7e14dcfSSatish Balay ierr = VecAXPY(DXU,-1.0,tron->Work);CHKERRQ(ierr); 306a7e14dcfSSatish Balay ierr = VecSet(tron->Work,0.0);CHKERRQ(ierr); 307a7e14dcfSSatish Balay ierr = VecPointwiseMin(DXU,tron->Work,DXU);CHKERRQ(ierr); 308a7e14dcfSSatish Balay PetscFunctionReturn(0); 309a7e14dcfSSatish Balay } 310a7e14dcfSSatish Balay 311a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 3121522df2eSJason Sarich /*MC 3131522df2eSJason Sarich TAOTRON - The TRON algorithm is an active-set Newton trust region method 3141522df2eSJason Sarich for bound-constrained minimization. 3151522df2eSJason Sarich 3161522df2eSJason Sarich Options Database Keys: 3171522df2eSJason Sarich + -tao_tron_maxgpits - maximum number of gradient projections per TRON iterate 3181522df2eSJason Sarich - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets 3191522df2eSJason Sarich 3201eb8069cSJason Sarich Level: beginner 3211522df2eSJason Sarich M*/ 322728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_TRON(Tao tao) 323a7e14dcfSSatish Balay { 324a7e14dcfSSatish Balay TAO_TRON *tron; 325a7e14dcfSSatish Balay PetscErrorCode ierr; 3268caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 327a7e14dcfSSatish Balay 32847a47007SBarry Smith PetscFunctionBegin; 329a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_TRON; 330a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_TRON; 331a7e14dcfSSatish Balay tao->ops->view = TaoView_TRON; 332a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_TRON; 333a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_TRON; 334a7e14dcfSSatish Balay tao->ops->computedual = TaoComputeDual_TRON; 335a7e14dcfSSatish Balay 3363c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&tron);CHKERRQ(ierr); 3376f4723b1SBarry Smith tao->data = (void*)tron; 3386552cf8aSJason Sarich 3396552cf8aSJason Sarich /* Override default settings (unless already changed) */ 3406552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 50; 3416552cf8aSJason Sarich if (!tao->trust0_changed) tao->trust0 = 1.0; 3424f1535bcSTodd Munson if (!tao->steptol_changed) tao->steptol = 0.0; 343a7e14dcfSSatish Balay 344a7e14dcfSSatish Balay /* Initialize pointers and variables */ 345a7e14dcfSSatish Balay tron->n = 0; 346a7e14dcfSSatish Balay tron->maxgpits = 3; 347a7e14dcfSSatish Balay tron->pg_ftol = 0.001; 348a7e14dcfSSatish Balay 349a7e14dcfSSatish Balay tron->eta1 = 1.0e-4; 350a7e14dcfSSatish Balay tron->eta2 = 0.25; 351a7e14dcfSSatish Balay tron->eta3 = 0.50; 352a7e14dcfSSatish Balay tron->eta4 = 0.90; 353a7e14dcfSSatish Balay 354a7e14dcfSSatish Balay tron->sigma1 = 0.5; 355a7e14dcfSSatish Balay tron->sigma2 = 2.0; 356a7e14dcfSSatish Balay tron->sigma3 = 4.0; 357a7e14dcfSSatish Balay 358a7e14dcfSSatish Balay tron->gp_iterates = 0; /* Cumulative number */ 359a7e14dcfSSatish Balay tron->total_gp_its = 0; 360a7e14dcfSSatish Balay tron->n_free = 0; 361a7e14dcfSSatish Balay 3626c23d075SBarry Smith tron->DXFree=NULL; 3636c23d075SBarry Smith tron->R=NULL; 3646c23d075SBarry Smith tron->X_New=NULL; 3656c23d075SBarry Smith tron->G_New=NULL; 3666c23d075SBarry Smith tron->Work=NULL; 3676c23d075SBarry Smith tron->Free_Local=NULL; 3686c23d075SBarry Smith tron->H_sub=NULL; 3696c23d075SBarry Smith tron->Hpre_sub=NULL; 370a7e14dcfSSatish Balay tao->subset_type = TAO_SUBSET_SUBVEC; 371a7e14dcfSSatish Balay 372a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 373a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 374441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 3755d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 376a7e14dcfSSatish Balay 377a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);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