1aaa7dc30SBarry Smith #include <../src/tao/bound/impls/tron/tron.h> 2aaa7dc30SBarry Smith #include <petsc-private/kspimpl.h> 3aaa7dc30SBarry 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); 28a7e14dcfSSatish Balay ierr = PetscFree(tao->data); 29a7e14dcfSSatish Balay PetscFunctionReturn(0); 30a7e14dcfSSatish Balay } 31a7e14dcfSSatish Balay 32a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 33a7e14dcfSSatish Balay #undef __FUNCT__ 34a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_TRON" 35441846f8SBarry Smith static PetscErrorCode TaoSetFromOptions_TRON(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; 42a7e14dcfSSatish Balay ierr = PetscOptionsHead("Newton Trust Region Method for bound constrained optimization");CHKERRQ(ierr); 4347a47007SBarry Smith ierr = PetscOptionsInt("-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; 107a7e14dcfSSatish Balay PetscInt iter=0,its; 108*d8ec8e84SJason Sarich PetscBool is_stcg,is_nash,is_gltr; 109e4cb33bbSBarry Smith TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 110e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING; 111a7e14dcfSSatish Balay PetscReal prered,actred,delta,f,f_new,rhok,gdx,xdiff,stepsize; 11247a47007SBarry Smith 113a7e14dcfSSatish Balay PetscFunctionBegin; 114a7e14dcfSSatish Balay tron->pgstepsize=1.0; 115a7e14dcfSSatish Balay tao->trust = tao->trust0; 116a7e14dcfSSatish Balay /* Project the current point onto the feasible set */ 117a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 118a7e14dcfSSatish Balay ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 119a7e14dcfSSatish Balay ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 120a7e14dcfSSatish Balay 121a7e14dcfSSatish Balay ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&tron->f,tao->gradient);CHKERRQ(ierr); 122a7e14dcfSSatish Balay ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 12353506e15SBarry Smith 124a7e14dcfSSatish Balay ierr = VecWhichBetween(tao->XL,tao->solution,tao->XU,&tron->Free_Local);CHKERRQ(ierr); 125a7e14dcfSSatish Balay 126a7e14dcfSSatish Balay /* Project the gradient and calculate the norm */ 127a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tao->gradient,tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr); 128a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 129a7e14dcfSSatish Balay 13053506e15SBarry Smith if (PetscIsInfOrNanReal(tron->f) || PetscIsInfOrNanReal(tron->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf pr NaN"); 131a7e14dcfSSatish Balay if (tao->trust <= 0) { 132a7e14dcfSSatish Balay tao->trust=PetscMax(tron->gnorm*tron->gnorm,1.0); 133a7e14dcfSSatish Balay } 134a7e14dcfSSatish Balay 135a7e14dcfSSatish Balay tron->stepsize=tao->trust; 136a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, tron->f, tron->gnorm, 0.0, tron->stepsize, &reason);CHKERRQ(ierr); 137a7e14dcfSSatish Balay while (reason==TAO_CONTINUE_ITERATING){ 138a7e14dcfSSatish Balay 139a7e14dcfSSatish Balay ierr = TronGradientProjections(tao,tron);CHKERRQ(ierr); 140a7e14dcfSSatish Balay f=tron->f; delta=tao->trust; 141a7e14dcfSSatish Balay tron->n_free_last = tron->n_free; 142ffad9901SBarry Smith ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr); 143a7e14dcfSSatish Balay 144a7e14dcfSSatish Balay ierr = ISGetSize(tron->Free_Local, &tron->n_free);CHKERRQ(ierr); 145a7e14dcfSSatish Balay 146a7e14dcfSSatish Balay /* If no free variables */ 147a7e14dcfSSatish Balay if (tron->n_free == 0) { 148a7e14dcfSSatish Balay actred=0; 149a7e14dcfSSatish Balay PetscInfo(tao,"No free variables in tron iteration."); 150a7e14dcfSSatish Balay break; 151a7e14dcfSSatish Balay 152a7e14dcfSSatish Balay } 153a7e14dcfSSatish Balay /* use free_local to mask/submat gradient, hessian, stepdirection */ 154b98f30f2SJason Sarich ierr = TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->R);CHKERRQ(ierr); 155b98f30f2SJason Sarich ierr = TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->DXFree);CHKERRQ(ierr); 156a7e14dcfSSatish Balay ierr = VecSet(tron->DXFree,0.0);CHKERRQ(ierr); 157a7e14dcfSSatish Balay ierr = VecScale(tron->R, -1.0);CHKERRQ(ierr); 158b98f30f2SJason Sarich ierr = TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub);CHKERRQ(ierr); 159a7e14dcfSSatish Balay if (tao->hessian == tao->hessian_pre) { 160a7e14dcfSSatish Balay ierr = MatDestroy(&tron->Hpre_sub);CHKERRQ(ierr); 161a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)(tron->H_sub));CHKERRQ(ierr); 162a7e14dcfSSatish Balay tron->Hpre_sub = tron->H_sub; 163a7e14dcfSSatish Balay } else { 164b98f30f2SJason Sarich ierr = TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type,&tron->Hpre_sub);CHKERRQ(ierr); 165a7e14dcfSSatish Balay } 166a7e14dcfSSatish Balay ierr = KSPReset(tao->ksp);CHKERRQ(ierr); 16723ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub);CHKERRQ(ierr); 168*d8ec8e84SJason Sarich ierr = PetscObjectTypeCompare((PetscObject)(tao->ksp),KSPSTCG,&is_stcg);CHKERRQ(ierr); 169*d8ec8e84SJason Sarich ierr = PetscObjectTypeCompare((PetscObject)(tao->ksp),KSPNASH,&is_nash);CHKERRQ(ierr); 170*d8ec8e84SJason Sarich ierr = PetscObjectTypeCompare((PetscObject)(tao->ksp),KSPGLTR,&is_gltr);CHKERRQ(ierr); 171a7e14dcfSSatish Balay while (1) { 172a7e14dcfSSatish Balay 173a7e14dcfSSatish Balay /* Approximately solve the reduced linear system */ 174*d8ec8e84SJason Sarich if (is_stcg) { 175a7e14dcfSSatish Balay ierr = KSPSTCGSetRadius(tao->ksp,delta);CHKERRQ(ierr); 176*d8ec8e84SJason Sarich } else if (is_nash) { 177*d8ec8e84SJason Sarich ierr = KSPNASHSetRadius(tao->ksp,delta);CHKERRQ(ierr); 178*d8ec8e84SJason Sarich } else if (is_gltr) { 179*d8ec8e84SJason Sarich ierr = KSPGLTRSetRadius(tao->ksp,delta);CHKERRQ(ierr); 180*d8ec8e84SJason Sarich } 181*d8ec8e84SJason Sarich 182a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, tron->R, tron->DXFree);CHKERRQ(ierr); 183a7e14dcfSSatish Balay ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr); 184a7e14dcfSSatish Balay tao->ksp_its+=its; 185a7e14dcfSSatish Balay ierr = VecSet(tao->stepdirection,0.0);CHKERRQ(ierr); 186a7e14dcfSSatish Balay 187a7e14dcfSSatish Balay /* Add dxfree matrix to compute step direction vector */ 1884473680cSBarry Smith ierr = VecISAXPY(tao->stepdirection,tron->Free_Local,1.0,tron->DXFree);CHKERRQ(ierr); 189a7e14dcfSSatish Balay if (0) { 190a7e14dcfSSatish Balay PetscReal rhs,stepnorm; 191a7e14dcfSSatish Balay ierr = VecNorm(tron->R,NORM_2,&rhs);CHKERRQ(ierr); 192a7e14dcfSSatish Balay ierr = VecNorm(tron->DXFree,NORM_2,&stepnorm);CHKERRQ(ierr); 193335036cbSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"|rhs|=%g\t|s|=%g\n",(double)rhs,(double)stepnorm);CHKERRQ(ierr); 194a7e14dcfSSatish Balay } 195a7e14dcfSSatish Balay 196a7e14dcfSSatish Balay 197a7e14dcfSSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr); 198335036cbSBarry Smith ierr = PetscInfo1(tao,"Expected decrease in function value: %14.12e\n",(double)gdx);CHKERRQ(ierr); 199a7e14dcfSSatish Balay 200a7e14dcfSSatish Balay ierr = VecCopy(tao->solution, tron->X_New);CHKERRQ(ierr); 201a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tron->G_New);CHKERRQ(ierr); 202a7e14dcfSSatish Balay 203a7e14dcfSSatish Balay stepsize=1.0;f_new=f; 204a7e14dcfSSatish Balay 205a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr); 20647a47007SBarry Smith ierr = TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection,&stepsize,&ls_reason);CHKERRQ(ierr);CHKERRQ(ierr); 207a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 208a7e14dcfSSatish Balay 209a7e14dcfSSatish Balay ierr = MatMult(tao->hessian, tao->stepdirection, tron->Work);CHKERRQ(ierr); 210a7e14dcfSSatish Balay ierr = VecAYPX(tron->Work, 0.5, tao->gradient);CHKERRQ(ierr); 211a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, tron->Work, &prered);CHKERRQ(ierr); 212a7e14dcfSSatish Balay actred = f_new - f; 213a7e14dcfSSatish Balay if (actred<0) { 214a7e14dcfSSatish Balay rhok=PetscAbs(-actred/prered); 215a7e14dcfSSatish Balay } else { 216a7e14dcfSSatish Balay rhok=0.0; 217a7e14dcfSSatish Balay } 218a7e14dcfSSatish Balay 219a7e14dcfSSatish Balay /* Compare actual improvement to the quadratic model */ 220a7e14dcfSSatish Balay if (rhok > tron->eta1) { /* Accept the point */ 221a7e14dcfSSatish Balay /* d = x_new - x */ 222a7e14dcfSSatish Balay ierr = VecCopy(tron->X_New, tao->stepdirection);CHKERRQ(ierr); 223a7e14dcfSSatish Balay ierr = VecAXPY(tao->stepdirection, -1.0, tao->solution);CHKERRQ(ierr); 224a7e14dcfSSatish Balay 225a7e14dcfSSatish Balay ierr = VecNorm(tao->stepdirection, NORM_2, &xdiff);CHKERRQ(ierr); 226a7e14dcfSSatish Balay xdiff *= stepsize; 227a7e14dcfSSatish Balay 228a7e14dcfSSatish Balay /* Adjust trust region size */ 229a7e14dcfSSatish Balay if (rhok < tron->eta2 ){ 230a7e14dcfSSatish Balay delta = PetscMin(xdiff,delta)*tron->sigma1; 231a7e14dcfSSatish Balay } else if (rhok > tron->eta4 ){ 232a7e14dcfSSatish Balay delta= PetscMin(xdiff,delta)*tron->sigma3; 233a7e14dcfSSatish Balay } else if (rhok > tron->eta3 ){ 234a7e14dcfSSatish Balay delta=PetscMin(xdiff,delta)*tron->sigma2; 235a7e14dcfSSatish Balay } 236a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tron->G_New,tron->X_New, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr); 237a7e14dcfSSatish Balay if (tron->Free_Local) { 238a7e14dcfSSatish Balay ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 239a7e14dcfSSatish Balay } 240a7e14dcfSSatish Balay ierr = VecWhichBetween(tao->XL, tron->X_New, tao->XU, &tron->Free_Local);CHKERRQ(ierr); 241a7e14dcfSSatish Balay f=f_new; 242a7e14dcfSSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr); 243a7e14dcfSSatish Balay ierr = VecCopy(tron->X_New, tao->solution);CHKERRQ(ierr); 244a7e14dcfSSatish Balay ierr = VecCopy(tron->G_New, tao->gradient);CHKERRQ(ierr); 245a7e14dcfSSatish Balay break; 246a7e14dcfSSatish Balay } 247a7e14dcfSSatish Balay else if (delta <= 1e-30) { 248a7e14dcfSSatish Balay break; 249a7e14dcfSSatish Balay } 250a7e14dcfSSatish Balay else { 251a7e14dcfSSatish Balay delta /= 4.0; 252a7e14dcfSSatish Balay } 253a7e14dcfSSatish Balay } /* end linear solve loop */ 254a7e14dcfSSatish Balay 255a7e14dcfSSatish Balay 256a7e14dcfSSatish Balay tron->f=f; tron->actred=actred; tao->trust=delta; 257a7e14dcfSSatish Balay iter++; 258a7e14dcfSSatish Balay ierr = TaoMonitor(tao, iter, tron->f, tron->gnorm, 0.0, delta, &reason);CHKERRQ(ierr); 259a7e14dcfSSatish Balay } /* END MAIN LOOP */ 260a7e14dcfSSatish Balay 261a7e14dcfSSatish Balay PetscFunctionReturn(0); 262a7e14dcfSSatish Balay } 263a7e14dcfSSatish Balay 264a7e14dcfSSatish Balay 265a7e14dcfSSatish Balay #undef __FUNCT__ 266a7e14dcfSSatish Balay #define __FUNCT__ "TronGradientProjections" 267441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao tao,TAO_TRON *tron) 268a7e14dcfSSatish Balay { 269a7e14dcfSSatish Balay PetscErrorCode ierr; 270a7e14dcfSSatish Balay PetscInt i; 271e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 272a7e14dcfSSatish Balay PetscReal actred=-1.0,actred_max=0.0; 273a7e14dcfSSatish Balay PetscReal f_new; 274a7e14dcfSSatish Balay /* 275a7e14dcfSSatish Balay The gradient and function value passed into and out of this 276a7e14dcfSSatish Balay routine should be current and correct. 277a7e14dcfSSatish Balay 278a7e14dcfSSatish Balay The free, active, and binding variables should be already identified 279a7e14dcfSSatish Balay */ 280a7e14dcfSSatish Balay PetscFunctionBegin; 281a7e14dcfSSatish Balay if (tron->Free_Local) { 282a7e14dcfSSatish Balay ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 283a7e14dcfSSatish Balay } 284a7e14dcfSSatish Balay ierr = VecWhichBetween(tao->XL,tao->solution,tao->XU,&tron->Free_Local);CHKERRQ(ierr); 285a7e14dcfSSatish Balay 286a7e14dcfSSatish Balay for (i=0;i<tron->maxgpits;i++){ 287a7e14dcfSSatish Balay 288a7e14dcfSSatish Balay if ( -actred <= (tron->pg_ftol)*actred_max) break; 289a7e14dcfSSatish Balay 290a7e14dcfSSatish Balay tron->gp_iterates++; tron->total_gp_its++; 291a7e14dcfSSatish Balay f_new=tron->f; 292a7e14dcfSSatish Balay 293a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 294a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 295a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,tron->pgstepsize);CHKERRQ(ierr); 296a7e14dcfSSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection, 297a7e14dcfSSatish Balay &tron->pgstepsize, &ls_reason);CHKERRQ(ierr); 298a7e14dcfSSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 299a7e14dcfSSatish Balay 300a7e14dcfSSatish Balay 301a7e14dcfSSatish Balay /* Update the iterate */ 302a7e14dcfSSatish Balay actred = f_new - tron->f; 303a7e14dcfSSatish Balay actred_max = PetscMax(actred_max,-(f_new - tron->f)); 304a7e14dcfSSatish Balay tron->f = f_new; 305a7e14dcfSSatish Balay if (tron->Free_Local) { 306a7e14dcfSSatish Balay ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr); 307a7e14dcfSSatish Balay } 308a7e14dcfSSatish Balay ierr = VecWhichBetween(tao->XL,tao->solution,tao->XU,&tron->Free_Local);CHKERRQ(ierr); 309a7e14dcfSSatish Balay } 310a7e14dcfSSatish Balay 311a7e14dcfSSatish Balay PetscFunctionReturn(0); 312a7e14dcfSSatish Balay } 313a7e14dcfSSatish Balay 314a7e14dcfSSatish Balay #undef __FUNCT__ 315a7e14dcfSSatish Balay #define __FUNCT__ "TaoComputeDual_TRON" 316441846f8SBarry Smith static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU) { 317a7e14dcfSSatish Balay 318a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data; 319a7e14dcfSSatish Balay PetscErrorCode ierr; 320a7e14dcfSSatish Balay 321a7e14dcfSSatish Balay PetscFunctionBegin; 322441846f8SBarry Smith PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 323a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXL,VEC_CLASSID,2); 324a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXU,VEC_CLASSID,3); 32547a47007SBarry Smith if (!tron->Work || !tao->gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n"); 326a7e14dcfSSatish Balay 327a7e14dcfSSatish Balay ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tron->Work);CHKERRQ(ierr); 328a7e14dcfSSatish Balay ierr = VecCopy(tron->Work,DXL);CHKERRQ(ierr); 329a7e14dcfSSatish Balay ierr = VecAXPY(DXL,-1.0,tao->gradient);CHKERRQ(ierr); 330a7e14dcfSSatish Balay ierr = VecSet(DXU,0.0);CHKERRQ(ierr); 331a7e14dcfSSatish Balay ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr); 332a7e14dcfSSatish Balay 333a7e14dcfSSatish Balay ierr = VecCopy(tao->gradient,DXU);CHKERRQ(ierr); 334a7e14dcfSSatish Balay ierr = VecAXPY(DXU,-1.0,tron->Work);CHKERRQ(ierr); 335a7e14dcfSSatish Balay ierr = VecSet(tron->Work,0.0);CHKERRQ(ierr); 336a7e14dcfSSatish Balay ierr = VecPointwiseMin(DXU,tron->Work,DXU);CHKERRQ(ierr); 337a7e14dcfSSatish Balay PetscFunctionReturn(0); 338a7e14dcfSSatish Balay } 339a7e14dcfSSatish Balay 340a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 341a7e14dcfSSatish Balay EXTERN_C_BEGIN 342a7e14dcfSSatish Balay #undef __FUNCT__ 343a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_TRON" 344441846f8SBarry Smith PetscErrorCode TaoCreate_TRON(Tao tao) 345a7e14dcfSSatish Balay { 346a7e14dcfSSatish Balay TAO_TRON *tron; 347a7e14dcfSSatish Balay PetscErrorCode ierr; 3488caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT; 349a7e14dcfSSatish Balay 35047a47007SBarry Smith PetscFunctionBegin; 351a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_TRON; 352a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_TRON; 353a7e14dcfSSatish Balay tao->ops->view = TaoView_TRON; 354a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_TRON; 355a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_TRON; 356a7e14dcfSSatish Balay tao->ops->computedual = TaoComputeDual_TRON; 357a7e14dcfSSatish Balay 3583c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&tron);CHKERRQ(ierr); 359a7e14dcfSSatish Balay 360a7e14dcfSSatish Balay tao->max_it = 50; 3616f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 3626f4723b1SBarry Smith tao->fatol = 1e-5; 3636f4723b1SBarry Smith tao->frtol = 1e-5; 3646f4723b1SBarry Smith tao->steptol = 1e-6; 3656f4723b1SBarry Smith #else 366a7e14dcfSSatish Balay tao->fatol = 1e-10; 367a7e14dcfSSatish Balay tao->frtol = 1e-10; 368a7e14dcfSSatish Balay tao->steptol = 1e-12; 3696f4723b1SBarry Smith #endif 3706f4723b1SBarry Smith tao->data = (void*)tron; 371a7e14dcfSSatish Balay tao->trust0 = 1.0; 372a7e14dcfSSatish Balay 373a7e14dcfSSatish Balay /* Initialize pointers and variables */ 374a7e14dcfSSatish Balay tron->n = 0; 375a7e14dcfSSatish Balay tron->maxgpits = 3; 376a7e14dcfSSatish Balay tron->pg_ftol = 0.001; 377a7e14dcfSSatish Balay 378a7e14dcfSSatish Balay tron->eta1 = 1.0e-4; 379a7e14dcfSSatish Balay tron->eta2 = 0.25; 380a7e14dcfSSatish Balay tron->eta3 = 0.50; 381a7e14dcfSSatish Balay tron->eta4 = 0.90; 382a7e14dcfSSatish Balay 383a7e14dcfSSatish Balay tron->sigma1 = 0.5; 384a7e14dcfSSatish Balay tron->sigma2 = 2.0; 385a7e14dcfSSatish Balay tron->sigma3 = 4.0; 386a7e14dcfSSatish Balay 387a7e14dcfSSatish Balay tron->gp_iterates = 0; /* Cumulative number */ 388a7e14dcfSSatish Balay tron->total_gp_its = 0; 389a7e14dcfSSatish Balay tron->n_free = 0; 390a7e14dcfSSatish Balay 3916c23d075SBarry Smith tron->DXFree=NULL; 3926c23d075SBarry Smith tron->R=NULL; 3936c23d075SBarry Smith tron->X_New=NULL; 3946c23d075SBarry Smith tron->G_New=NULL; 3956c23d075SBarry Smith tron->Work=NULL; 3966c23d075SBarry Smith tron->Free_Local=NULL; 3976c23d075SBarry Smith tron->H_sub=NULL; 3986c23d075SBarry Smith tron->Hpre_sub=NULL; 399a7e14dcfSSatish Balay tao->subset_type = TAO_SUBSET_SUBVEC; 400a7e14dcfSSatish Balay 401a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 402a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr); 403441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr); 404a7e14dcfSSatish Balay 405a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr); 406a7e14dcfSSatish Balay ierr = KSPSetType(tao->ksp,KSPSTCG);CHKERRQ(ierr); 407a7e14dcfSSatish Balay PetscFunctionReturn(0); 408a7e14dcfSSatish Balay } 409a7e14dcfSSatish Balay EXTERN_C_END 410