1aaa7dc30SBarry Smith #include <../src/tao/bound/impls/tron/tron.h> 2af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 3af0996ceSBarry Smith #include <petsc/private/matimpl.h> 4aaa7dc30SBarry Smith #include <../src/tao/matrix/submatfree.h> 5a7e14dcfSSatish Balay 6a7e14dcfSSatish Balay 7a7e14dcfSSatish Balay /* TRON Routines */ 8441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao,TAO_TRON*); 9a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 10a7e14dcfSSatish Balay #undef __FUNCT__ 11a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_TRON" 12441846f8SBarry Smith static PetscErrorCode TaoDestroy_TRON(Tao tao) 13a7e14dcfSSatish Balay { 14a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 15a7e14dcfSSatish Balay PetscErrorCode ierr; 16a7e14dcfSSatish Balay 17a7e14dcfSSatish Balay PetscFunctionBegin; 18a7e14dcfSSatish Balay ierr = VecDestroy(&tron->X_New);CHKERRQ(ierr); 19a7e14dcfSSatish Balay ierr = VecDestroy(&tron->G_New);CHKERRQ(ierr); 20a7e14dcfSSatish Balay ierr = VecDestroy(&tron->Work);CHKERRQ(ierr); 21a7e14dcfSSatish Balay ierr = VecDestroy(&tron->DXFree);CHKERRQ(ierr); 22a7e14dcfSSatish Balay ierr = VecDestroy(&tron->R);CHKERRQ(ierr); 23a7e14dcfSSatish Balay ierr = VecDestroy(&tron->diag);CHKERRQ(ierr); 24a7e14dcfSSatish Balay ierr = VecScatterDestroy(&tron->scatter);CHKERRQ(ierr); 25a7e14dcfSSatish Balay ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 26a7e14dcfSSatish Balay ierr = MatDestroy(&tron->H_sub);CHKERRQ(ierr); 27a7e14dcfSSatish Balay ierr = MatDestroy(&tron->Hpre_sub);CHKERRQ(ierr); 28302440fdSBarry Smith ierr = PetscFree(tao->data);CHKERRQ(ierr); 29a7e14dcfSSatish Balay PetscFunctionReturn(0); 30a7e14dcfSSatish Balay } 31a7e14dcfSSatish Balay 32a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 33a7e14dcfSSatish Balay #undef __FUNCT__ 34a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_TRON" 358c34d3f5SBarry Smith static PetscErrorCode TaoSetFromOptions_TRON(PetscOptions *PetscOptionsObject,Tao tao) 36a7e14dcfSSatish Balay { 37a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 38a7e14dcfSSatish Balay PetscErrorCode ierr; 39a7e14dcfSSatish Balay PetscBool flg; 40a7e14dcfSSatish Balay 41a7e14dcfSSatish Balay PetscFunctionBegin; 421a1499c8SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Newton Trust Region Method for bound constrained optimization");CHKERRQ(ierr); 431522df2eSJason Sarich ierr = PetscOptionsInt("-tao_tron_maxgpits","maximum number of gradient projections per TRON iterate","TaoSetMaxGPIts",tron->maxgpits,&tron->maxgpits,&flg);CHKERRQ(ierr); 44a7e14dcfSSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 45a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 46a7e14dcfSSatish Balay ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 47a7e14dcfSSatish Balay PetscFunctionReturn(0); 48a7e14dcfSSatish Balay } 49a7e14dcfSSatish Balay 50a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 51a7e14dcfSSatish Balay #undef __FUNCT__ 52a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_TRON" 53441846f8SBarry Smith static PetscErrorCode TaoView_TRON(Tao tao, PetscViewer viewer) 54a7e14dcfSSatish Balay { 55a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 56a7e14dcfSSatish Balay PetscBool isascii; 57a7e14dcfSSatish Balay PetscErrorCode ierr; 58a7e14dcfSSatish Balay 59a7e14dcfSSatish Balay PetscFunctionBegin; 60a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 61a7e14dcfSSatish Balay if (isascii) { 62a7e14dcfSSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 63a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"Total PG its: %D,",tron->total_gp_its);CHKERRQ(ierr); 6447a47007SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"PG tolerance: %g \n",(double)tron->pg_ftol);CHKERRQ(ierr); 65a7e14dcfSSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 66a7e14dcfSSatish Balay } 67a7e14dcfSSatish Balay PetscFunctionReturn(0); 68a7e14dcfSSatish Balay } 69a7e14dcfSSatish Balay 70a7e14dcfSSatish Balay 71a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 72a7e14dcfSSatish Balay #undef __FUNCT__ 73a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetup_TRON" 74441846f8SBarry Smith static PetscErrorCode TaoSetup_TRON(Tao tao) 75a7e14dcfSSatish Balay { 76a7e14dcfSSatish Balay PetscErrorCode ierr; 77a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 78a7e14dcfSSatish Balay 79a7e14dcfSSatish Balay PetscFunctionBegin; 80a7e14dcfSSatish Balay 81a7e14dcfSSatish Balay /* Allocate some arrays */ 82a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tron->diag);CHKERRQ(ierr); 83a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tron->X_New);CHKERRQ(ierr); 84a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tron->G_New);CHKERRQ(ierr); 85a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tron->Work);CHKERRQ(ierr); 86a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); 87a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr); 88a7e14dcfSSatish Balay if (!tao->XL) { 89a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->XL);CHKERRQ(ierr); 90e270355aSBarry Smith ierr = VecSet(tao->XL, PETSC_NINFINITY);CHKERRQ(ierr); 91a7e14dcfSSatish Balay } 92a7e14dcfSSatish Balay if (!tao->XU) { 93a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &tao->XU);CHKERRQ(ierr); 94e270355aSBarry Smith ierr = VecSet(tao->XU, PETSC_INFINITY);CHKERRQ(ierr); 95a7e14dcfSSatish Balay } 96a7e14dcfSSatish Balay PetscFunctionReturn(0); 97a7e14dcfSSatish Balay } 98a7e14dcfSSatish Balay 99a7e14dcfSSatish Balay 100a7e14dcfSSatish Balay 101a7e14dcfSSatish Balay #undef __FUNCT__ 102a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_TRON" 103441846f8SBarry Smith static PetscErrorCode TaoSolve_TRON(Tao tao) 10447a47007SBarry Smith { 105a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 106a7e14dcfSSatish Balay PetscErrorCode ierr; 1078931d482SJason Sarich PetscInt its; 108e4cb33bbSBarry Smith TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 109e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING; 110a7e14dcfSSatish Balay PetscReal prered,actred,delta,f,f_new,rhok,gdx,xdiff,stepsize; 11147a47007SBarry Smith 112a7e14dcfSSatish Balay PetscFunctionBegin; 113a7e14dcfSSatish Balay tron->pgstepsize=1.0; 114a7e14dcfSSatish Balay tao->trust = tao->trust0; 115a7e14dcfSSatish Balay /* Project the current point onto the feasible set */ 116a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 117a7e14dcfSSatish Balay ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 118a7e14dcfSSatish Balay ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 119a7e14dcfSSatish Balay 120a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&tron->f,tao->gradient);CHKERRQ(ierr); 121a7e14dcfSSatish Balay ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 12253506e15SBarry Smith 123a7e14dcfSSatish Balay ierr = VecWhichBetween(tao->XL,tao->solution,tao->XU,&tron->Free_Local);CHKERRQ(ierr); 124a7e14dcfSSatish Balay 125a7e14dcfSSatish Balay /* Project the gradient and calculate the norm */ 126a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tao->gradient,tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr); 127a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 128a7e14dcfSSatish Balay 12953506e15SBarry Smith if (PetscIsInfOrNanReal(tron->f) || PetscIsInfOrNanReal(tron->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf pr NaN"); 130a7e14dcfSSatish Balay if (tao->trust <= 0) { 131a7e14dcfSSatish Balay tao->trust=PetscMax(tron->gnorm*tron->gnorm,1.0); 132a7e14dcfSSatish Balay } 133a7e14dcfSSatish Balay 134a7e14dcfSSatish Balay tron->stepsize=tao->trust; 1358931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, tron->stepsize, &reason);CHKERRQ(ierr); 136a7e14dcfSSatish Balay while (reason==TAO_CONTINUE_ITERATING){ 137ae93cb3cSJason Sarich tao->ksp_its=0; 138a7e14dcfSSatish Balay ierr = TronGradientProjections(tao,tron);CHKERRQ(ierr); 139a7e14dcfSSatish Balay f=tron->f; delta=tao->trust; 140a7e14dcfSSatish Balay tron->n_free_last = tron->n_free; 141ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 142a7e14dcfSSatish Balay 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) { 147a7e14dcfSSatish Balay actred=0; 148955c1f14SBarry Smith ierr = PetscInfo(tao,"No free variables in tron iteration.\n");CHKERRQ(ierr); 14955723ba5SJason Sarich ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 1508931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, delta, &reason);CHKERRQ(ierr); 15155723ba5SJason Sarich if (!reason) { 15255723ba5SJason Sarich reason = TAO_CONVERGED_STEPTOL; 15355723ba5SJason Sarich ierr = TaoSetConvergedReason(tao,reason);CHKERRQ(ierr); 15455723ba5SJason Sarich } 15555723ba5SJason Sarich 156a7e14dcfSSatish Balay break; 157a7e14dcfSSatish Balay 158a7e14dcfSSatish Balay } 159a7e14dcfSSatish Balay /* use free_local to mask/submat gradient, hessian, stepdirection */ 160b98f30f2SJason Sarich ierr = TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->R);CHKERRQ(ierr); 161b98f30f2SJason Sarich ierr = TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->DXFree);CHKERRQ(ierr); 162a7e14dcfSSatish Balay ierr = VecSet(tron->DXFree,0.0);CHKERRQ(ierr); 163a7e14dcfSSatish Balay ierr = VecScale(tron->R, -1.0);CHKERRQ(ierr); 164b98f30f2SJason Sarich ierr = TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub);CHKERRQ(ierr); 165a7e14dcfSSatish Balay if (tao->hessian == tao->hessian_pre) { 166a7e14dcfSSatish Balay ierr = MatDestroy(&tron->Hpre_sub);CHKERRQ(ierr); 167a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)(tron->H_sub));CHKERRQ(ierr); 168a7e14dcfSSatish Balay tron->Hpre_sub = tron->H_sub; 169a7e14dcfSSatish Balay } else { 170b98f30f2SJason Sarich ierr = TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type,&tron->Hpre_sub);CHKERRQ(ierr); 171a7e14dcfSSatish Balay } 172a7e14dcfSSatish Balay ierr = KSPReset(tao->ksp);CHKERRQ(ierr); 17323ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub);CHKERRQ(ierr); 174a7e14dcfSSatish Balay while (1) { 175a7e14dcfSSatish Balay 176a7e14dcfSSatish Balay /* Approximately solve the reduced linear system */ 177a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(tao->ksp,delta);CHKERRQ(ierr); 178d8ec8e84SJason Sarich 179a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tron->R, tron->DXFree);CHKERRQ(ierr); 180a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 181a7e14dcfSSatish Balay tao->ksp_its+=its; 182ae93cb3cSJason Sarich tao->ksp_tot_its+=its; 183a7e14dcfSSatish Balay ierr = VecSet(tao->stepdirection,0.0);CHKERRQ(ierr); 184a7e14dcfSSatish Balay 185a7e14dcfSSatish Balay /* Add dxfree matrix to compute step direction vector */ 1864473680cSBarry Smith ierr = VecISAXPY(tao->stepdirection,tron->Free_Local,1.0,tron->DXFree);CHKERRQ(ierr); 187a7e14dcfSSatish Balay if (0) { 188a7e14dcfSSatish Balay PetscReal rhs,stepnorm; 189a7e14dcfSSatish Balay ierr = VecNorm(tron->R,NORM_2,&rhs);CHKERRQ(ierr); 190a7e14dcfSSatish Balay ierr = VecNorm(tron->DXFree,NORM_2,&stepnorm);CHKERRQ(ierr); 191335036cbSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"|rhs|=%g\t|s|=%g\n",(double)rhs,(double)stepnorm);CHKERRQ(ierr); 192a7e14dcfSSatish Balay } 193a7e14dcfSSatish Balay 194a7e14dcfSSatish Balay 195a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 196335036cbSBarry Smith ierr = PetscInfo1(tao,"Expected decrease in function value: %14.12e\n",(double)gdx);CHKERRQ(ierr); 197a7e14dcfSSatish Balay 198a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tron->X_New);CHKERRQ(ierr); 199a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tron->G_New);CHKERRQ(ierr); 200a7e14dcfSSatish Balay 201a7e14dcfSSatish Balay stepsize=1.0;f_new=f; 202a7e14dcfSSatish Balay 203a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 20447a47007SBarry Smith ierr = TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection,&stepsize,&ls_reason);CHKERRQ(ierr);CHKERRQ(ierr); 205a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 206a7e14dcfSSatish Balay 207a7e14dcfSSatish Balay ierr = MatMult(tao->hessian, tao->stepdirection, tron->Work);CHKERRQ(ierr); 208a7e14dcfSSatish Balay ierr = VecAYPX(tron->Work, 0.5, tao->gradient);CHKERRQ(ierr); 209a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tron->Work, &prered);CHKERRQ(ierr); 210a7e14dcfSSatish Balay actred = f_new - f; 211a7e14dcfSSatish Balay if (actred<0) { 212a7e14dcfSSatish Balay rhok=PetscAbs(-actred/prered); 213a7e14dcfSSatish Balay } else { 214a7e14dcfSSatish Balay rhok=0.0; 215a7e14dcfSSatish Balay } 216a7e14dcfSSatish Balay 217a7e14dcfSSatish Balay /* Compare actual improvement to the quadratic model */ 218a7e14dcfSSatish Balay if (rhok > tron->eta1) { /* Accept the point */ 219a7e14dcfSSatish Balay /* d = x_new - x */ 220a7e14dcfSSatish Balay ierr = VecCopy(tron->X_New, tao->stepdirection);CHKERRQ(ierr); 221a7e14dcfSSatish Balay ierr = VecAXPY(tao->stepdirection, -1.0, tao->solution);CHKERRQ(ierr); 222a7e14dcfSSatish Balay 223a7e14dcfSSatish Balay ierr = VecNorm(tao->stepdirection, NORM_2, &xdiff);CHKERRQ(ierr); 224a7e14dcfSSatish Balay xdiff *= stepsize; 225a7e14dcfSSatish Balay 226a7e14dcfSSatish Balay /* Adjust trust region size */ 227a7e14dcfSSatish Balay if (rhok < tron->eta2 ){ 228a7e14dcfSSatish Balay delta = PetscMin(xdiff,delta)*tron->sigma1; 229a7e14dcfSSatish Balay } else if (rhok > tron->eta4 ){ 230a7e14dcfSSatish Balay delta= PetscMin(xdiff,delta)*tron->sigma3; 231a7e14dcfSSatish Balay } else if (rhok > tron->eta3 ){ 232a7e14dcfSSatish Balay delta=PetscMin(xdiff,delta)*tron->sigma2; 233a7e14dcfSSatish Balay } 234a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tron->G_New,tron->X_New, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr); 235a7e14dcfSSatish Balay if (tron->Free_Local) { 236a7e14dcfSSatish Balay ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 237a7e14dcfSSatish Balay } 238a7e14dcfSSatish Balay ierr = VecWhichBetween(tao->XL, tron->X_New, tao->XU, &tron->Free_Local);CHKERRQ(ierr); 239a7e14dcfSSatish Balay f=f_new; 240a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 241a7e14dcfSSatish Balay ierr = VecCopy(tron->X_New, tao->solution);CHKERRQ(ierr); 242a7e14dcfSSatish Balay ierr = VecCopy(tron->G_New, tao->gradient);CHKERRQ(ierr); 243a7e14dcfSSatish Balay break; 244a7e14dcfSSatish Balay } 245a7e14dcfSSatish Balay else if (delta <= 1e-30) { 246a7e14dcfSSatish Balay break; 247a7e14dcfSSatish Balay } 248a7e14dcfSSatish Balay else { 249a7e14dcfSSatish Balay delta /= 4.0; 250a7e14dcfSSatish Balay } 251a7e14dcfSSatish Balay } /* end linear solve loop */ 252a7e14dcfSSatish Balay 253a7e14dcfSSatish Balay 254a7e14dcfSSatish Balay tron->f=f; tron->actred=actred; tao->trust=delta; 2558931d482SJason Sarich tao->niter++; 2568931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, delta, &reason);CHKERRQ(ierr); 257a7e14dcfSSatish Balay } /* END MAIN LOOP */ 258a7e14dcfSSatish Balay 259a7e14dcfSSatish Balay PetscFunctionReturn(0); 260a7e14dcfSSatish Balay } 261a7e14dcfSSatish Balay 262a7e14dcfSSatish Balay 263a7e14dcfSSatish Balay #undef __FUNCT__ 264a7e14dcfSSatish Balay #define __FUNCT__ "TronGradientProjections" 265441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao tao,TAO_TRON *tron) 266a7e14dcfSSatish Balay { 267a7e14dcfSSatish Balay PetscErrorCode ierr; 268a7e14dcfSSatish Balay PetscInt i; 269e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 270a7e14dcfSSatish Balay PetscReal actred=-1.0,actred_max=0.0; 271a7e14dcfSSatish Balay PetscReal f_new; 272a7e14dcfSSatish Balay /* 273a7e14dcfSSatish Balay The gradient and function value passed into and out of this 274a7e14dcfSSatish Balay routine should be current and correct. 275a7e14dcfSSatish Balay 276a7e14dcfSSatish Balay The free, active, and binding variables should be already identified 277a7e14dcfSSatish Balay */ 278a7e14dcfSSatish Balay PetscFunctionBegin; 279a7e14dcfSSatish Balay if (tron->Free_Local) { 280a7e14dcfSSatish Balay ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 281a7e14dcfSSatish Balay } 282a7e14dcfSSatish Balay ierr = VecWhichBetween(tao->XL,tao->solution,tao->XU,&tron->Free_Local);CHKERRQ(ierr); 283a7e14dcfSSatish Balay 284a7e14dcfSSatish Balay for (i=0;i<tron->maxgpits;i++){ 285a7e14dcfSSatish Balay 286a7e14dcfSSatish Balay if ( -actred <= (tron->pg_ftol)*actred_max) break; 287a7e14dcfSSatish Balay 288a7e14dcfSSatish Balay tron->gp_iterates++; tron->total_gp_its++; 289a7e14dcfSSatish Balay f_new=tron->f; 290a7e14dcfSSatish Balay 291a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 292a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 293a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,tron->pgstepsize);CHKERRQ(ierr); 294a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection, 295a7e14dcfSSatish Balay &tron->pgstepsize, &ls_reason);CHKERRQ(ierr); 296a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 297a7e14dcfSSatish Balay 298a7e14dcfSSatish Balay 299a7e14dcfSSatish Balay /* Update the iterate */ 300a7e14dcfSSatish Balay actred = f_new - tron->f; 301a7e14dcfSSatish Balay actred_max = PetscMax(actred_max,-(f_new - tron->f)); 302a7e14dcfSSatish Balay tron->f = f_new; 303a7e14dcfSSatish Balay if (tron->Free_Local) { 304a7e14dcfSSatish Balay ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 305a7e14dcfSSatish Balay } 306a7e14dcfSSatish Balay ierr = VecWhichBetween(tao->XL,tao->solution,tao->XU,&tron->Free_Local);CHKERRQ(ierr); 307a7e14dcfSSatish Balay } 308a7e14dcfSSatish Balay 309a7e14dcfSSatish Balay PetscFunctionReturn(0); 310a7e14dcfSSatish Balay } 311a7e14dcfSSatish Balay 312a7e14dcfSSatish Balay #undef __FUNCT__ 313a7e14dcfSSatish Balay #define __FUNCT__ "TaoComputeDual_TRON" 314441846f8SBarry Smith static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU) { 315a7e14dcfSSatish Balay 316a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 317a7e14dcfSSatish Balay PetscErrorCode ierr; 318a7e14dcfSSatish Balay 319a7e14dcfSSatish Balay PetscFunctionBegin; 320441846f8SBarry Smith PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 321a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 322a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 32347a47007SBarry Smith if (!tron->Work || !tao->gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n"); 324a7e14dcfSSatish Balay 325a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tron->Work);CHKERRQ(ierr); 326a7e14dcfSSatish Balay ierr = VecCopy(tron->Work,DXL);CHKERRQ(ierr); 327a7e14dcfSSatish Balay ierr = VecAXPY(DXL,-1.0,tao->gradient);CHKERRQ(ierr); 328a7e14dcfSSatish Balay ierr = VecSet(DXU,0.0);CHKERRQ(ierr); 329a7e14dcfSSatish Balay ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr); 330a7e14dcfSSatish Balay 331a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient,DXU);CHKERRQ(ierr); 332a7e14dcfSSatish Balay ierr = VecAXPY(DXU,-1.0,tron->Work);CHKERRQ(ierr); 333a7e14dcfSSatish Balay ierr = VecSet(tron->Work,0.0);CHKERRQ(ierr); 334a7e14dcfSSatish Balay ierr = VecPointwiseMin(DXU,tron->Work,DXU);CHKERRQ(ierr); 335a7e14dcfSSatish Balay PetscFunctionReturn(0); 336a7e14dcfSSatish Balay } 337a7e14dcfSSatish Balay 338a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 3391522df2eSJason Sarich /*MC 3401522df2eSJason Sarich TAOTRON - The TRON algorithm is an active-set Newton trust region method 3411522df2eSJason Sarich for bound-constrained minimization. 3421522df2eSJason Sarich 3431522df2eSJason Sarich Options Database Keys: 3441522df2eSJason Sarich + -tao_tron_maxgpits - maximum number of gradient projections per TRON iterate 3451522df2eSJason Sarich - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets 3461522df2eSJason Sarich 3471eb8069cSJason Sarich Level: beginner 3481522df2eSJason Sarich M*/ 349a7e14dcfSSatish Balay #undef __FUNCT__ 350a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_TRON" 351728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_TRON(Tao tao) 352a7e14dcfSSatish Balay { 353a7e14dcfSSatish Balay TAO_TRON *tron; 354a7e14dcfSSatish Balay PetscErrorCode ierr; 3558caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 356a7e14dcfSSatish Balay 35747a47007SBarry Smith PetscFunctionBegin; 358a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_TRON; 359a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_TRON; 360a7e14dcfSSatish Balay tao->ops->view = TaoView_TRON; 361a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_TRON; 362a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_TRON; 363a7e14dcfSSatish Balay tao->ops->computedual = TaoComputeDual_TRON; 364a7e14dcfSSatish Balay 3653c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&tron);CHKERRQ(ierr); 366a7e14dcfSSatish Balay 367a7e14dcfSSatish Balay tao->max_it = 50; 3686f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 3696f4723b1SBarry Smith tao->fatol = 1e-5; 3706f4723b1SBarry Smith tao->frtol = 1e-5; 3716f4723b1SBarry Smith tao->steptol = 1e-6; 3726f4723b1SBarry Smith #else 373a7e14dcfSSatish Balay tao->fatol = 1e-10; 374a7e14dcfSSatish Balay tao->frtol = 1e-10; 375a7e14dcfSSatish Balay tao->steptol = 1e-12; 3766f4723b1SBarry Smith #endif 3776f4723b1SBarry Smith tao->data = (void*)tron; 378a7e14dcfSSatish Balay tao->trust0 = 1.0; 379a7e14dcfSSatish Balay 380a7e14dcfSSatish Balay /* Initialize pointers and variables */ 381a7e14dcfSSatish Balay tron->n = 0; 382a7e14dcfSSatish Balay tron->maxgpits = 3; 383a7e14dcfSSatish Balay tron->pg_ftol = 0.001; 384a7e14dcfSSatish Balay 385a7e14dcfSSatish Balay tron->eta1 = 1.0e-4; 386a7e14dcfSSatish Balay tron->eta2 = 0.25; 387a7e14dcfSSatish Balay tron->eta3 = 0.50; 388a7e14dcfSSatish Balay tron->eta4 = 0.90; 389a7e14dcfSSatish Balay 390a7e14dcfSSatish Balay tron->sigma1 = 0.5; 391a7e14dcfSSatish Balay tron->sigma2 = 2.0; 392a7e14dcfSSatish Balay tron->sigma3 = 4.0; 393a7e14dcfSSatish Balay 394a7e14dcfSSatish Balay tron->gp_iterates = 0; /* Cumulative number */ 395a7e14dcfSSatish Balay tron->total_gp_its = 0; 396a7e14dcfSSatish Balay tron->n_free = 0; 397a7e14dcfSSatish Balay 3986c23d075SBarry Smith tron->DXFree=NULL; 3996c23d075SBarry Smith tron->R=NULL; 4006c23d075SBarry Smith tron->X_New=NULL; 4016c23d075SBarry Smith tron->G_New=NULL; 4026c23d075SBarry Smith tron->Work=NULL; 4036c23d075SBarry Smith tron->Free_Local=NULL; 4046c23d075SBarry Smith tron->H_sub=NULL; 4056c23d075SBarry Smith tron->Hpre_sub=NULL; 406a7e14dcfSSatish Balay tao->subset_type = TAO_SUBSET_SUBVEC; 407a7e14dcfSSatish Balay 408a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 409a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 410441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 411*5d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 412a7e14dcfSSatish Balay 413a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr); 414*5d527766SPatrick Farrell ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr); 415a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp,KSPSTCG);CHKERRQ(ierr); 416a7e14dcfSSatish Balay PetscFunctionReturn(0); 417a7e14dcfSSatish Balay } 418728e0ed0SBarry Smith 419