1aaa7dc30SBarry Smith #include <petsc-private/taolinesearchimpl.h> 2aaa7dc30SBarry Smith #include <../src/tao/linesearch/impls/gpcglinesearch/gpcglinesearch.h> 3a7e14dcfSSatish Balay 4a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 5a7e14dcfSSatish Balay 6a7e14dcfSSatish Balay #undef __FUNCT__ 7a7e14dcfSSatish Balay #define __FUNCT__ "TaoLineSearchDestroy_GPCG" 8a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchDestroy_GPCG(TaoLineSearch ls) 9a7e14dcfSSatish Balay { 10a7e14dcfSSatish Balay PetscErrorCode ierr; 118caf6e8cSBarry Smith TaoLineSearch_GPCG *ctx = (TaoLineSearch_GPCG *)ls->data; 12a7e14dcfSSatish Balay 13a7e14dcfSSatish Balay PetscFunctionBegin; 14a7e14dcfSSatish Balay ierr = VecDestroy(&ctx->W1);CHKERRQ(ierr); 15a7e14dcfSSatish Balay ierr = VecDestroy(&ctx->W2);CHKERRQ(ierr); 16a7e14dcfSSatish Balay ierr = VecDestroy(&ctx->Gold);CHKERRQ(ierr); 17f06e3bfaSBarry Smith ierr = VecDestroy(&ctx->x);CHKERRQ(ierr); 18f06e3bfaSBarry Smith ierr = PetscFree(ls->data);CHKERRQ(ierr); 19a7e14dcfSSatish Balay PetscFunctionReturn(0); 20a7e14dcfSSatish Balay } 21a7e14dcfSSatish Balay 22a7e14dcfSSatish Balay 23a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 24a7e14dcfSSatish Balay #undef __FUNCT__ 25a7e14dcfSSatish Balay #define __FUNCT__ "TaoLineSearchView_GPCG" 26a7e14dcfSSatish Balay static PetscErrorCode TaoLineSearchView_GPCG(TaoLineSearch ls, PetscViewer viewer) 27a7e14dcfSSatish Balay { 28a7e14dcfSSatish Balay PetscBool isascii; 29a7e14dcfSSatish Balay PetscErrorCode ierr; 30f06e3bfaSBarry Smith 31a7e14dcfSSatish Balay PetscFunctionBegin; 32a7e14dcfSSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 33a7e14dcfSSatish Balay if (isascii) { 34a7e14dcfSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," GPCG Line search");CHKERRQ(ierr); 35a7e14dcfSSatish Balay } 36a7e14dcfSSatish Balay PetscFunctionReturn(0); 37a7e14dcfSSatish Balay } 38a7e14dcfSSatish Balay 39a7e14dcfSSatish Balay /*------------------------------------------------------------*/ 40a7e14dcfSSatish Balay #undef __FUNCT__ 41a7e14dcfSSatish Balay #define __FUNCT__ "TaoLineSearchApply_GPCG" 42f06e3bfaSBarry Smith static PetscErrorCode TaoLineSearchApply_GPCG(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, Vec s) 43a7e14dcfSSatish Balay { 448caf6e8cSBarry Smith TaoLineSearch_GPCG *neP = (TaoLineSearch_GPCG *)ls->data; 45a7e14dcfSSatish Balay PetscErrorCode ierr; 46a7e14dcfSSatish Balay PetscInt i; 47a7e14dcfSSatish Balay PetscBool g_computed=PETSC_FALSE; /* to prevent extra gradient computation */ 48a7e14dcfSSatish Balay PetscReal d1,finit,actred,prered,rho, gdx; 49a7e14dcfSSatish Balay 50a7e14dcfSSatish Balay PetscFunctionBegin; 51a7e14dcfSSatish Balay /* ls->stepmin - lower bound for step */ 52a7e14dcfSSatish Balay /* ls->stepmax - upper bound for step */ 53a7e14dcfSSatish Balay /* ls->rtol - relative tolerance for an acceptable step */ 54a7e14dcfSSatish Balay /* ls->ftol - tolerance for sufficient decrease condition */ 55a7e14dcfSSatish Balay /* ls->gtol - tolerance for curvature condition */ 56a7e14dcfSSatish Balay /* ls->nfeval - number of function evaluations */ 57a7e14dcfSSatish Balay /* ls->nfeval - number of function/gradient evaluations */ 58a7e14dcfSSatish Balay /* ls->max_funcs - maximum number of function evaluations */ 59a7e14dcfSSatish Balay 60a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_CONTINUE_ITERATING; 61a7e14dcfSSatish Balay ls->step = ls->initstep; 62a7e14dcfSSatish Balay if (!neP->W2) { 63a7e14dcfSSatish Balay ierr = VecDuplicate(x,&neP->W2);CHKERRQ(ierr); 64a7e14dcfSSatish Balay ierr = VecDuplicate(x,&neP->W1);CHKERRQ(ierr); 65a7e14dcfSSatish Balay ierr = VecDuplicate(x,&neP->Gold);CHKERRQ(ierr); 66a7e14dcfSSatish Balay neP->x = x; 67a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)neP->x);CHKERRQ(ierr); 68f06e3bfaSBarry Smith } else if (x != neP->x) { 69a7e14dcfSSatish Balay ierr = VecDestroy(&neP->x);CHKERRQ(ierr); 70a7e14dcfSSatish Balay ierr = VecDestroy(&neP->W1);CHKERRQ(ierr); 71a7e14dcfSSatish Balay ierr = VecDestroy(&neP->W2);CHKERRQ(ierr); 72a7e14dcfSSatish Balay ierr = VecDestroy(&neP->Gold);CHKERRQ(ierr); 73a7e14dcfSSatish Balay ierr = VecDuplicate(x,&neP->W1);CHKERRQ(ierr); 74a7e14dcfSSatish Balay ierr = VecDuplicate(x,&neP->W2);CHKERRQ(ierr); 75a7e14dcfSSatish Balay ierr = VecDuplicate(x,&neP->Gold);CHKERRQ(ierr); 76a7e14dcfSSatish Balay ierr = PetscObjectDereference((PetscObject)neP->x);CHKERRQ(ierr); 77a7e14dcfSSatish Balay neP->x = x; 78a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)neP->x);CHKERRQ(ierr); 79a7e14dcfSSatish Balay } 80a7e14dcfSSatish Balay 81a7e14dcfSSatish Balay ierr = VecDot(g,s,&gdx);CHKERRQ(ierr); 82a7e14dcfSSatish Balay if (gdx > 0) { 83*c783eb9eSBarry Smith ierr = PetscInfo1(ls,"Line search error: search direction is not descent direction. dot(g,s) = %g\n",(double)gdx);CHKERRQ(ierr); 84a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_FAILED_ASCENT; 85a7e14dcfSSatish Balay PetscFunctionReturn(0); 86a7e14dcfSSatish Balay } 87a7e14dcfSSatish Balay ierr = VecCopy(x,neP->W2);CHKERRQ(ierr); 88a7e14dcfSSatish Balay ierr = VecCopy(g,neP->Gold);CHKERRQ(ierr); 89a7e14dcfSSatish Balay if (ls->bounded) { 90f06e3bfaSBarry Smith /* Compute the smallest steplength that will make one nonbinding variable equal the bound */ 91e270355aSBarry Smith ierr = VecStepBoundInfo(x,s,ls->lower,ls->upper,&rho,&actred,&d1);CHKERRQ(ierr); 92a7e14dcfSSatish Balay ls->step = PetscMin(ls->step,d1); 93a7e14dcfSSatish Balay } 94a7e14dcfSSatish Balay rho=0; actred=0; 95a7e14dcfSSatish Balay 96a7e14dcfSSatish Balay if (ls->step < 0) { 97f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Line search error: initial step parameter %g< 0\n",(double)ls->step);CHKERRQ(ierr); 98a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_OTHER; 99a7e14dcfSSatish Balay PetscFunctionReturn(0); 100a7e14dcfSSatish Balay } 101a7e14dcfSSatish Balay 102a7e14dcfSSatish Balay /* Initialization */ 103a7e14dcfSSatish Balay finit = *f; 104a7e14dcfSSatish Balay for (i=0; i< ls->max_funcs; i++) { 105a7e14dcfSSatish Balay /* Force the step to be within the bounds */ 106a7e14dcfSSatish Balay ls->step = PetscMax(ls->step,ls->stepmin); 107a7e14dcfSSatish Balay ls->step = PetscMin(ls->step,ls->stepmax); 108a7e14dcfSSatish Balay 109a7e14dcfSSatish Balay ierr = VecCopy(x,neP->W2);CHKERRQ(ierr); 110a7e14dcfSSatish Balay ierr = VecAXPY(neP->W2,ls->step,s);CHKERRQ(ierr); 111a7e14dcfSSatish Balay if (ls->bounded) { 112a7e14dcfSSatish Balay /* Make sure new vector is numerically within bounds */ 113a7e14dcfSSatish Balay ierr = VecMedian(neP->W2,ls->lower,ls->upper,neP->W2);CHKERRQ(ierr); 114a7e14dcfSSatish Balay } 115a7e14dcfSSatish Balay 116a7e14dcfSSatish Balay /* Gradient is not needed here. Unless there is a separate 117a7e14dcfSSatish Balay gradient routine, compute it here anyway to prevent recomputing at 118a7e14dcfSSatish Balay the end of the line search */ 119a7e14dcfSSatish Balay if (ls->hasobjective) { 120a7e14dcfSSatish Balay ierr = TaoLineSearchComputeObjective(ls,neP->W2,f);CHKERRQ(ierr); 121a7e14dcfSSatish Balay g_computed=PETSC_FALSE; 122a7e14dcfSSatish Balay } else if (ls->usegts){ 123a7e14dcfSSatish Balay ierr = TaoLineSearchComputeObjectiveAndGTS(ls,neP->W2,f,&gdx);CHKERRQ(ierr); 124a7e14dcfSSatish Balay g_computed=PETSC_FALSE; 125a7e14dcfSSatish Balay } else { 126a7e14dcfSSatish Balay ierr = TaoLineSearchComputeObjectiveAndGradient(ls,neP->W2,f,g);CHKERRQ(ierr); 127a7e14dcfSSatish Balay g_computed=PETSC_TRUE; 128a7e14dcfSSatish Balay } 129a7e14dcfSSatish Balay 130a7e14dcfSSatish Balay if (0 == i) { 131a7e14dcfSSatish Balay ls->f_fullstep = *f; 132a7e14dcfSSatish Balay } 133a7e14dcfSSatish Balay 134a7e14dcfSSatish Balay actred = *f - finit; 135a7e14dcfSSatish Balay ierr = VecCopy(neP->W2,neP->W1);CHKERRQ(ierr); 136a7e14dcfSSatish Balay ierr = VecAXPY(neP->W1,-1.0,x);CHKERRQ(ierr); /* W1 = W2 - X */ 137a7e14dcfSSatish Balay ierr = VecDot(neP->W1,neP->Gold,&prered);CHKERRQ(ierr); 138a7e14dcfSSatish Balay 139a7e14dcfSSatish Balay if (fabs(prered)<1.0e-100) prered=1.0e-12; 140a7e14dcfSSatish Balay rho = actred/prered; 141a7e14dcfSSatish Balay 142a7e14dcfSSatish Balay /* 143a7e14dcfSSatish Balay If sufficient progress has been obtained, accept the 144a7e14dcfSSatish Balay point. Otherwise, backtrack. 145a7e14dcfSSatish Balay */ 146a7e14dcfSSatish Balay 147a7e14dcfSSatish Balay if (actred > 0) { 148a7e14dcfSSatish Balay ierr = PetscInfo(ls,"Step resulted in ascent, rejecting.\n");CHKERRQ(ierr); 149a7e14dcfSSatish Balay ls->step = (ls->step)/2; 150a7e14dcfSSatish Balay } else if (rho > ls->ftol){ 151a7e14dcfSSatish Balay break; 152a7e14dcfSSatish Balay } else{ 153a7e14dcfSSatish Balay ls->step = (ls->step)/2; 154a7e14dcfSSatish Balay } 155a7e14dcfSSatish Balay 156a7e14dcfSSatish Balay /* Convergence testing */ 157a7e14dcfSSatish Balay 158a7e14dcfSSatish Balay if (ls->step <= ls->stepmin || ls->step >= ls->stepmax) { 159a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_OTHER; 160a7e14dcfSSatish Balay ierr = PetscInfo(ls,"Rounding errors may prevent further progress. May not be a step satisfying\n");CHKERRQ(ierr); 161a7e14dcfSSatish Balay ierr = PetscInfo(ls,"sufficient decrease and curvature conditions. Tolerances may be too small.\n");CHKERRQ(ierr); 162a7e14dcfSSatish Balay break; 163a7e14dcfSSatish Balay } 164a7e14dcfSSatish Balay if (ls->step == ls->stepmax) { 165f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Step is at the upper bound, stepmax (%g)\n",(double)ls->stepmax);CHKERRQ(ierr); 166a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_UPPERBOUND; 167a7e14dcfSSatish Balay break; 168a7e14dcfSSatish Balay } 169a7e14dcfSSatish Balay if (ls->step == ls->stepmin) { 170f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Step is at the lower bound, stepmin (%g)\n",(double)ls->stepmin);CHKERRQ(ierr); 171a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_LOWERBOUND; 172a7e14dcfSSatish Balay break; 173a7e14dcfSSatish Balay } 174a7e14dcfSSatish Balay if ((ls->nfeval+ls->nfgeval) >= ls->max_funcs) { 175a7e14dcfSSatish Balay ierr = PetscInfo2(ls,"Number of line search function evals (%D) > maximum (%D)\n",ls->nfeval+ls->nfgeval,ls->max_funcs);CHKERRQ(ierr); 176a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_MAXFCN; 177a7e14dcfSSatish Balay break; 178a7e14dcfSSatish Balay } 179a7e14dcfSSatish Balay if ((neP->bracket) && (ls->stepmax - ls->stepmin <= ls->rtol*ls->stepmax)){ 180f06e3bfaSBarry Smith ierr = PetscInfo1(ls,"Relative width of interval of uncertainty is at most rtol (%g)\n",(double)ls->rtol);CHKERRQ(ierr); 181a7e14dcfSSatish Balay ls->reason = TAOLINESEARCH_HALTED_RTOL; 182a7e14dcfSSatish Balay break; 183a7e14dcfSSatish Balay } 184a7e14dcfSSatish Balay } 185f06e3bfaSBarry Smith ierr = PetscInfo2(ls,"%D function evals in line search, step = %g\n",ls->nfeval+ls->nfgeval,(double)ls->step);CHKERRQ(ierr); 186a7e14dcfSSatish Balay /* set new solution vector and compute gradient if necessary */ 187a7e14dcfSSatish Balay ierr = VecCopy(neP->W2, x);CHKERRQ(ierr); 188a7e14dcfSSatish Balay if (!g_computed) { 189a7e14dcfSSatish Balay ierr = TaoLineSearchComputeGradient(ls,x,g);CHKERRQ(ierr); 190a7e14dcfSSatish Balay } 191a7e14dcfSSatish Balay PetscFunctionReturn(0); 192a7e14dcfSSatish Balay } 193a7e14dcfSSatish Balay 194a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 195a7e14dcfSSatish Balay #undef __FUNCT__ 196a7e14dcfSSatish Balay #define __FUNCT__ "TaoLineSearchCreate_GPCG" 197728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoLineSearchCreate_GPCG(TaoLineSearch ls) 198a7e14dcfSSatish Balay { 199a7e14dcfSSatish Balay PetscErrorCode ierr; 2008caf6e8cSBarry Smith TaoLineSearch_GPCG *neP; 201a7e14dcfSSatish Balay 202a7e14dcfSSatish Balay PetscFunctionBegin; 203a7e14dcfSSatish Balay ls->ftol = 0.05; 204a7e14dcfSSatish Balay ls->rtol = 0.0; 205a7e14dcfSSatish Balay ls->gtol = 0.0; 206a7e14dcfSSatish Balay ls->stepmin = 1.0e-20; 207a7e14dcfSSatish Balay ls->stepmax = 1.0e+20; 208a7e14dcfSSatish Balay ls->nfeval = 0; 209a7e14dcfSSatish Balay ls->max_funcs = 30; 210a7e14dcfSSatish Balay ls->step = 1.0; 211a7e14dcfSSatish Balay 2123c9e27cfSGeoffrey Irving ierr = PetscNewLog(ls,&neP);CHKERRQ(ierr); 213a7e14dcfSSatish Balay neP->bracket = 0; 214a7e14dcfSSatish Balay neP->infoc = 1; 215a7e14dcfSSatish Balay ls->data = (void*)neP; 216a7e14dcfSSatish Balay 217a7e14dcfSSatish Balay ls->ops->setup = 0; 21824517cd3SJason Sarich ls->ops->reset = 0; 219a7e14dcfSSatish Balay ls->ops->apply=TaoLineSearchApply_GPCG; 220a7e14dcfSSatish Balay ls->ops->view =TaoLineSearchView_GPCG; 221a7e14dcfSSatish Balay ls->ops->destroy=TaoLineSearchDestroy_GPCG; 222a7e14dcfSSatish Balay ls->ops->setfromoptions=0; 223a7e14dcfSSatish Balay PetscFunctionReturn(0); 224a7e14dcfSSatish Balay } 225728e0ed0SBarry Smith 226