1*5e1affc2SSatish Balay #include "taolinesearch.h" 2*5e1affc2SSatish Balay #include "taocg.h" 3*5e1affc2SSatish Balay 4*5e1affc2SSatish Balay #define CG_FletcherReeves 0 5*5e1affc2SSatish Balay #define CG_PolakRibiere 1 6*5e1affc2SSatish Balay #define CG_PolakRibierePlus 2 7*5e1affc2SSatish Balay #define CG_HestenesStiefel 3 8*5e1affc2SSatish Balay #define CG_DaiYuan 4 9*5e1affc2SSatish Balay #define CG_Types 5 10*5e1affc2SSatish Balay 11*5e1affc2SSatish Balay static const char *CG_Table[64] = { 12*5e1affc2SSatish Balay "fr", "pr", "prp", "hs", "dy" 13*5e1affc2SSatish Balay }; 14*5e1affc2SSatish Balay 15*5e1affc2SSatish Balay 16*5e1affc2SSatish Balay #undef __FUNCT__ 17*5e1affc2SSatish Balay #define __FUNCT__ "TaoSolve_CG" 18*5e1affc2SSatish Balay static PetscErrorCode TaoSolve_CG(TaoSolver tao) 19*5e1affc2SSatish Balay { 20*5e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG*)tao->data; 21*5e1affc2SSatish Balay PetscErrorCode ierr; 22*5e1affc2SSatish Balay TaoSolverTerminationReason reason = TAO_CONTINUE_ITERATING; 23*5e1affc2SSatish Balay TaoLineSearchTerminationReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 24*5e1affc2SSatish Balay PetscReal step=1.0,f,gnorm,gnorm2,delta,gd,ginner,beta; 25*5e1affc2SSatish Balay PetscReal gd_old,gnorm2_old,f_old; 26*5e1affc2SSatish Balay PetscInt iter=0; 27*5e1affc2SSatish Balay 28*5e1affc2SSatish Balay 29*5e1affc2SSatish Balay PetscFunctionBegin; 30*5e1affc2SSatish Balay 31*5e1affc2SSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 32*5e1affc2SSatish Balay ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by cg algorithm\n"); CHKERRQ(ierr); 33*5e1affc2SSatish Balay } 34*5e1affc2SSatish Balay 35*5e1affc2SSatish Balay /* Check convergence criteria */ 36*5e1affc2SSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient); CHKERRQ(ierr); 37*5e1affc2SSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&gnorm); CHKERRQ(ierr); 38*5e1affc2SSatish Balay if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) { 39*5e1affc2SSatish Balay SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 40*5e1affc2SSatish Balay } 41*5e1affc2SSatish Balay 42*5e1affc2SSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason); CHKERRQ(ierr); 43*5e1affc2SSatish Balay if (reason != TAO_CONTINUE_ITERATING) { 44*5e1affc2SSatish Balay PetscFunctionReturn(0); 45*5e1affc2SSatish Balay } 46*5e1affc2SSatish Balay 47*5e1affc2SSatish Balay 48*5e1affc2SSatish Balay /* Set initial direction to -gradient */ 49*5e1affc2SSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection); CHKERRQ(ierr); 50*5e1affc2SSatish Balay ierr = VecScale(tao->stepdirection, -1.0); CHKERRQ(ierr); 51*5e1affc2SSatish Balay gnorm2 = gnorm*gnorm; 52*5e1affc2SSatish Balay 53*5e1affc2SSatish Balay /* Set initial scaling for the function */ 54*5e1affc2SSatish Balay if (f != 0.0) { 55*5e1affc2SSatish Balay delta = 2.0*PetscAbsScalar(f) / gnorm2; 56*5e1affc2SSatish Balay delta = PetscMax(delta,cgP->delta_min); 57*5e1affc2SSatish Balay delta = PetscMin(delta,cgP->delta_max); 58*5e1affc2SSatish Balay } else { 59*5e1affc2SSatish Balay delta = 2.0 / gnorm2; 60*5e1affc2SSatish Balay delta = PetscMax(delta,cgP->delta_min); 61*5e1affc2SSatish Balay delta = PetscMin(delta,cgP->delta_max); 62*5e1affc2SSatish Balay } 63*5e1affc2SSatish Balay /* Set counter for gradient and reset steps */ 64*5e1affc2SSatish Balay cgP->ngradsteps = 0; 65*5e1affc2SSatish Balay cgP->nresetsteps = 0; 66*5e1affc2SSatish Balay 67*5e1affc2SSatish Balay while (1) { 68*5e1affc2SSatish Balay /* Save the current gradient information */ 69*5e1affc2SSatish Balay f_old = f; 70*5e1affc2SSatish Balay gnorm2_old = gnorm2; 71*5e1affc2SSatish Balay ierr = VecCopy(tao->solution, cgP->X_old); CHKERRQ(ierr); 72*5e1affc2SSatish Balay ierr = VecCopy(tao->gradient, cgP->G_old); CHKERRQ(ierr); 73*5e1affc2SSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gd); CHKERRQ(ierr); 74*5e1affc2SSatish Balay if ((gd >= 0) || PetscIsInfOrNanReal(gd)) { 75*5e1affc2SSatish Balay ++cgP->ngradsteps; 76*5e1affc2SSatish Balay if (f != 0.0) { 77*5e1affc2SSatish Balay delta = 2.0*PetscAbsScalar(f) / gnorm2; 78*5e1affc2SSatish Balay delta = PetscMax(delta,cgP->delta_min); 79*5e1affc2SSatish Balay delta = PetscMin(delta,cgP->delta_max); 80*5e1affc2SSatish Balay } else { 81*5e1affc2SSatish Balay delta = 2.0 / gnorm2; 82*5e1affc2SSatish Balay delta = PetscMax(delta,cgP->delta_min); 83*5e1affc2SSatish Balay delta = PetscMin(delta,cgP->delta_max); 84*5e1affc2SSatish Balay } 85*5e1affc2SSatish Balay 86*5e1affc2SSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection); CHKERRQ(ierr); 87*5e1affc2SSatish Balay ierr = VecScale(tao->stepdirection, -1.0); CHKERRQ(ierr); 88*5e1affc2SSatish Balay } 89*5e1affc2SSatish Balay 90*5e1affc2SSatish Balay /* Search direction for improving point */ 91*5e1affc2SSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,delta); 92*5e1affc2SSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status); CHKERRQ(ierr); 93*5e1affc2SSatish Balay ierr = TaoAddLineSearchCounts(tao); CHKERRQ(ierr); 94*5e1affc2SSatish Balay if (ls_status != TAOLINESEARCH_SUCCESS && 95*5e1affc2SSatish Balay ls_status != TAOLINESEARCH_SUCCESS_USER) { 96*5e1affc2SSatish Balay /* Linesearch failed */ 97*5e1affc2SSatish Balay /* Reset factors and use scaled gradient step */ 98*5e1affc2SSatish Balay ++cgP->nresetsteps; 99*5e1affc2SSatish Balay f = f_old; 100*5e1affc2SSatish Balay gnorm2 = gnorm2_old; 101*5e1affc2SSatish Balay ierr = VecCopy(cgP->X_old, tao->solution); CHKERRQ(ierr); 102*5e1affc2SSatish Balay ierr = VecCopy(cgP->G_old, tao->gradient); CHKERRQ(ierr); 103*5e1affc2SSatish Balay 104*5e1affc2SSatish Balay if (f != 0.0) { 105*5e1affc2SSatish Balay delta = 2.0*PetscAbsScalar(f) / gnorm2; 106*5e1affc2SSatish Balay delta = PetscMax(delta,cgP->delta_min); 107*5e1affc2SSatish Balay delta = PetscMin(delta,cgP->delta_max); 108*5e1affc2SSatish Balay } else { 109*5e1affc2SSatish Balay delta = 2.0 / gnorm2; 110*5e1affc2SSatish Balay delta = PetscMax(delta,cgP->delta_min); 111*5e1affc2SSatish Balay delta = PetscMin(delta,cgP->delta_max); 112*5e1affc2SSatish Balay } 113*5e1affc2SSatish Balay 114*5e1affc2SSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection); CHKERRQ(ierr); 115*5e1affc2SSatish Balay ierr = VecScale(tao->stepdirection, -1.0); CHKERRQ(ierr); 116*5e1affc2SSatish Balay 117*5e1affc2SSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,delta); 118*5e1affc2SSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status); CHKERRQ(ierr); 119*5e1affc2SSatish Balay ierr = TaoAddLineSearchCounts(tao); CHKERRQ(ierr); 120*5e1affc2SSatish Balay 121*5e1affc2SSatish Balay if (ls_status != TAOLINESEARCH_SUCCESS && 122*5e1affc2SSatish Balay ls_status != TAOLINESEARCH_SUCCESS_USER) { 123*5e1affc2SSatish Balay /* Linesearch failed again */ 124*5e1affc2SSatish Balay /* switch to unscaled gradient */ 125*5e1affc2SSatish Balay f = f_old; 126*5e1affc2SSatish Balay gnorm2 = gnorm2_old; 127*5e1affc2SSatish Balay ierr = VecCopy(cgP->X_old, tao->solution); CHKERRQ(ierr); 128*5e1affc2SSatish Balay ierr = VecCopy(cgP->G_old, tao->gradient); CHKERRQ(ierr); 129*5e1affc2SSatish Balay delta = 1.0; 130*5e1affc2SSatish Balay ierr = VecCopy(tao->solution, tao->stepdirection); CHKERRQ(ierr); 131*5e1affc2SSatish Balay ierr = VecScale(tao->stepdirection, -1.0); CHKERRQ(ierr); 132*5e1affc2SSatish Balay 133*5e1affc2SSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,delta); 134*5e1affc2SSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status); CHKERRQ(ierr); 135*5e1affc2SSatish Balay ierr = TaoAddLineSearchCounts(tao); CHKERRQ(ierr); 136*5e1affc2SSatish Balay if (ls_status != TAOLINESEARCH_SUCCESS && 137*5e1affc2SSatish Balay ls_status != TAOLINESEARCH_SUCCESS_USER) { 138*5e1affc2SSatish Balay 139*5e1affc2SSatish Balay /* Line search failed for last time -- give up */ 140*5e1affc2SSatish Balay f = f_old; 141*5e1affc2SSatish Balay gnorm2 = gnorm2_old; 142*5e1affc2SSatish Balay ierr = VecCopy(cgP->X_old, tao->solution); CHKERRQ(ierr); 143*5e1affc2SSatish Balay ierr = VecCopy(cgP->G_old, tao->gradient); CHKERRQ(ierr); 144*5e1affc2SSatish Balay step = 0.0; 145*5e1affc2SSatish Balay reason = TAO_DIVERGED_LS_FAILURE; 146*5e1affc2SSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 147*5e1affc2SSatish Balay } 148*5e1affc2SSatish Balay } 149*5e1affc2SSatish Balay } 150*5e1affc2SSatish Balay 151*5e1affc2SSatish Balay /* Check for bad value */ 152*5e1affc2SSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&gnorm); CHKERRQ(ierr); 153*5e1affc2SSatish Balay if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) { 154*5e1affc2SSatish Balay SETERRQ(PETSC_COMM_SELF,1,"User-provided compute function generated Inf or NaN"); 155*5e1affc2SSatish Balay } 156*5e1affc2SSatish Balay 157*5e1affc2SSatish Balay /* Check for termination */ 158*5e1affc2SSatish Balay gnorm2 =gnorm * gnorm; 159*5e1affc2SSatish Balay iter++; 160*5e1affc2SSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason); CHKERRQ(ierr); 161*5e1affc2SSatish Balay if (reason != TAO_CONTINUE_ITERATING) { 162*5e1affc2SSatish Balay break; 163*5e1affc2SSatish Balay } 164*5e1affc2SSatish Balay 165*5e1affc2SSatish Balay /* Check for restart condition */ 166*5e1affc2SSatish Balay ierr = VecDot(tao->gradient, cgP->G_old, &ginner); CHKERRQ(ierr); 167*5e1affc2SSatish Balay if (PetscAbsScalar(ginner) >= cgP->eta * gnorm2) { 168*5e1affc2SSatish Balay /* Gradients far from orthognal; use steepest descent direction */ 169*5e1affc2SSatish Balay beta = 0.0; 170*5e1affc2SSatish Balay } else { 171*5e1affc2SSatish Balay /* Gradients close to orthogonal; use conjugate gradient formula */ 172*5e1affc2SSatish Balay switch (cgP->cg_type) { 173*5e1affc2SSatish Balay case CG_FletcherReeves: 174*5e1affc2SSatish Balay beta = gnorm2 / gnorm2_old; 175*5e1affc2SSatish Balay break; 176*5e1affc2SSatish Balay 177*5e1affc2SSatish Balay case CG_PolakRibiere: 178*5e1affc2SSatish Balay beta = (gnorm2 - ginner) / gnorm2_old; 179*5e1affc2SSatish Balay break; 180*5e1affc2SSatish Balay 181*5e1affc2SSatish Balay case CG_PolakRibierePlus: 182*5e1affc2SSatish Balay beta = PetscMax((gnorm2-ginner)/gnorm2_old, 0.0); 183*5e1affc2SSatish Balay break; 184*5e1affc2SSatish Balay 185*5e1affc2SSatish Balay case CG_HestenesStiefel: 186*5e1affc2SSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gd); CHKERRQ(ierr); 187*5e1affc2SSatish Balay ierr = VecDot(cgP->G_old, tao->stepdirection, &gd_old); CHKERRQ(ierr); 188*5e1affc2SSatish Balay beta = (gnorm2 - ginner) / (gd - gd_old); 189*5e1affc2SSatish Balay break; 190*5e1affc2SSatish Balay 191*5e1affc2SSatish Balay case CG_DaiYuan: 192*5e1affc2SSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gd); CHKERRQ(ierr); 193*5e1affc2SSatish Balay ierr = VecDot(cgP->G_old, tao->stepdirection, &gd_old); CHKERRQ(ierr); 194*5e1affc2SSatish Balay beta = gnorm2 / (gd - gd_old); 195*5e1affc2SSatish Balay break; 196*5e1affc2SSatish Balay 197*5e1affc2SSatish Balay default: 198*5e1affc2SSatish Balay beta = 0.0; 199*5e1affc2SSatish Balay break; 200*5e1affc2SSatish Balay } 201*5e1affc2SSatish Balay } 202*5e1affc2SSatish Balay 203*5e1affc2SSatish Balay /* Compute the direction d=-g + beta*d */ 204*5e1affc2SSatish Balay ierr = VecAXPBY(tao->stepdirection, -1.0, beta, tao->gradient); CHKERRQ(ierr); 205*5e1affc2SSatish Balay 206*5e1affc2SSatish Balay /* update initial steplength choice */ 207*5e1affc2SSatish Balay delta = 1.0; 208*5e1affc2SSatish Balay delta = PetscMax(delta, cgP->delta_min); 209*5e1affc2SSatish Balay delta = PetscMin(delta, cgP->delta_max); 210*5e1affc2SSatish Balay } 211*5e1affc2SSatish Balay PetscFunctionReturn(0); 212*5e1affc2SSatish Balay } 213*5e1affc2SSatish Balay 214*5e1affc2SSatish Balay #undef __FUNCT__ 215*5e1affc2SSatish Balay #define __FUNCT__ "TaoSetUp_CG" 216*5e1affc2SSatish Balay static PetscErrorCode TaoSetUp_CG(TaoSolver tao) 217*5e1affc2SSatish Balay { 218*5e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG*)tao->data; 219*5e1affc2SSatish Balay PetscErrorCode ierr; 220*5e1affc2SSatish Balay 221*5e1affc2SSatish Balay PetscFunctionBegin; 222*5e1affc2SSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient); CHKERRQ(ierr);} 223*5e1affc2SSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection); CHKERRQ(ierr); } 224*5e1affc2SSatish Balay if (!cgP->X_old) {ierr = VecDuplicate(tao->solution,&cgP->X_old); CHKERRQ(ierr);} 225*5e1affc2SSatish Balay if (!cgP->G_old) {ierr = VecDuplicate(tao->gradient,&cgP->G_old); CHKERRQ(ierr); } 226*5e1affc2SSatish Balay 227*5e1affc2SSatish Balay PetscFunctionReturn(0); 228*5e1affc2SSatish Balay } 229*5e1affc2SSatish Balay 230*5e1affc2SSatish Balay 231*5e1affc2SSatish Balay #undef __FUNCT__ 232*5e1affc2SSatish Balay #define __FUNCT__ "TaoDestroy_CG" 233*5e1affc2SSatish Balay static PetscErrorCode TaoDestroy_CG(TaoSolver tao) 234*5e1affc2SSatish Balay { 235*5e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG*) tao->data; 236*5e1affc2SSatish Balay PetscErrorCode ierr; 237*5e1affc2SSatish Balay 238*5e1affc2SSatish Balay PetscFunctionBegin; 239*5e1affc2SSatish Balay 240*5e1affc2SSatish Balay if (tao->setupcalled) { 241*5e1affc2SSatish Balay ierr = VecDestroy(&cgP->X_old); CHKERRQ(ierr); 242*5e1affc2SSatish Balay ierr = VecDestroy(&cgP->G_old); CHKERRQ(ierr); 243*5e1affc2SSatish Balay } 244*5e1affc2SSatish Balay ierr = TaoLineSearchDestroy(&tao->linesearch); CHKERRQ(ierr); 245*5e1affc2SSatish Balay 246*5e1affc2SSatish Balay ierr = PetscFree(tao->data); CHKERRQ(ierr); 247*5e1affc2SSatish Balay tao->data = PETSC_NULL; 248*5e1affc2SSatish Balay 249*5e1affc2SSatish Balay PetscFunctionReturn(0); 250*5e1affc2SSatish Balay } 251*5e1affc2SSatish Balay 252*5e1affc2SSatish Balay #undef __FUNCT__ 253*5e1affc2SSatish Balay #define __FUNCT__ "TaoSetFromOptions_CG" 254*5e1affc2SSatish Balay static PetscErrorCode TaoSetFromOptions_CG(TaoSolver tao) 255*5e1affc2SSatish Balay { 256*5e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG*)tao->data; 257*5e1affc2SSatish Balay PetscErrorCode ierr; 258*5e1affc2SSatish Balay 259*5e1affc2SSatish Balay PetscFunctionBegin; 260*5e1affc2SSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch); CHKERRQ(ierr); 261*5e1affc2SSatish Balay ierr = PetscOptionsHead("Nonlinear Conjugate Gradient method for unconstrained optimization"); CHKERRQ(ierr); 262*5e1affc2SSatish Balay ierr = PetscOptionsReal("-tao_cg_eta","restart tolerance", "", cgP->eta,&cgP->eta,0);CHKERRQ(ierr); 263*5e1affc2SSatish Balay ierr = PetscOptionsEList("-tao_cg_type","cg formula", "", CG_Table, CG_Types, CG_Table[cgP->cg_type], &cgP->cg_type, 0); CHKERRQ(ierr); 264*5e1affc2SSatish Balay ierr = PetscOptionsReal("-tao_cg_delta_min","minimum delta value", "", cgP-> 265*5e1affc2SSatish Balay delta_min,&cgP->delta_min,0); CHKERRQ(ierr); 266*5e1affc2SSatish Balay ierr = PetscOptionsReal("-tao_cg_delta_max","maximum delta value", "", cgP-> 267*5e1affc2SSatish Balay delta_max,&cgP->delta_max,0); CHKERRQ(ierr); 268*5e1affc2SSatish Balay ierr = PetscOptionsTail(); CHKERRQ(ierr); 269*5e1affc2SSatish Balay PetscFunctionReturn(0); 270*5e1affc2SSatish Balay } 271*5e1affc2SSatish Balay 272*5e1affc2SSatish Balay #undef __FUNCT__ 273*5e1affc2SSatish Balay #define __FUNCT__ "TaoView_CG" 274*5e1affc2SSatish Balay static PetscErrorCode TaoView_CG(TaoSolver tao, PetscViewer viewer) 275*5e1affc2SSatish Balay { 276*5e1affc2SSatish Balay PetscBool isascii; 277*5e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG*)tao->data; 278*5e1affc2SSatish Balay PetscErrorCode ierr; 279*5e1affc2SSatish Balay 280*5e1affc2SSatish Balay PetscFunctionBegin; 281*5e1affc2SSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii); CHKERRQ(ierr); 282*5e1affc2SSatish Balay if (isascii) { 283*5e1affc2SSatish Balay ierr = PetscViewerASCIIPushTab(viewer); CHKERRQ(ierr); 284*5e1affc2SSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", CG_Table[cgP->cg_type]); CHKERRQ(ierr); 285*5e1affc2SSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", cgP->ngradsteps); CHKERRQ(ierr); 286*5e1affc2SSatish Balay ierr= PetscViewerASCIIPrintf(viewer, "Reset steps: %D\n", cgP->nresetsteps); CHKERRQ(ierr); 287*5e1affc2SSatish Balay ierr = PetscViewerASCIIPopTab(viewer); CHKERRQ(ierr); 288*5e1affc2SSatish Balay } else { 289*5e1affc2SSatish Balay SETERRQ1(((PetscObject)tao)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for TAO CG",((PetscObject)viewer)->type_name); 290*5e1affc2SSatish Balay } 291*5e1affc2SSatish Balay PetscFunctionReturn(0); 292*5e1affc2SSatish Balay } 293*5e1affc2SSatish Balay 294*5e1affc2SSatish Balay 295*5e1affc2SSatish Balay EXTERN_C_BEGIN 296*5e1affc2SSatish Balay #undef __FUNCT__ 297*5e1affc2SSatish Balay #define __FUNCT__ "TaoCreate_CG" 298*5e1affc2SSatish Balay PetscErrorCode TaoCreate_CG(TaoSolver tao) 299*5e1affc2SSatish Balay { 300*5e1affc2SSatish Balay TAO_CG *cgP; 301*5e1affc2SSatish Balay const char *morethuente_type = TAOLINESEARCH_MT; 302*5e1affc2SSatish Balay PetscErrorCode ierr; 303*5e1affc2SSatish Balay PetscFunctionBegin; 304*5e1affc2SSatish Balay tao->ops->setup = TaoSetUp_CG; 305*5e1affc2SSatish Balay tao->ops->solve = TaoSolve_CG; 306*5e1affc2SSatish Balay tao->ops->view = TaoView_CG; 307*5e1affc2SSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_CG; 308*5e1affc2SSatish Balay tao->ops->destroy = TaoDestroy_CG; 309*5e1affc2SSatish Balay 310*5e1affc2SSatish Balay tao->max_it = 2000; 311*5e1affc2SSatish Balay tao->max_funcs = 4000; 312*5e1affc2SSatish Balay tao->fatol = 1e-4; 313*5e1affc2SSatish Balay tao->frtol = 1e-4; 314*5e1affc2SSatish Balay 315*5e1affc2SSatish Balay /* Note: nondefault values should be used for nonlinear conjugate gradient */ 316*5e1affc2SSatish Balay /* method. In particular, gtol should be less that 0.5; the value used in */ 317*5e1affc2SSatish Balay /* Nocedal and Wright is 0.10. We use the default values for the */ 318*5e1affc2SSatish Balay /* linesearch because it seems to work better. */ 319*5e1affc2SSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch); CHKERRQ(ierr); 320*5e1affc2SSatish Balay ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type); CHKERRQ(ierr); 321*5e1affc2SSatish Balay ierr = TaoLineSearchUseTaoSolverRoutines(tao->linesearch, tao); CHKERRQ(ierr); 322*5e1affc2SSatish Balay 323*5e1affc2SSatish Balay 324*5e1affc2SSatish Balay ierr = PetscNewLog(tao,&cgP); CHKERRQ(ierr); 325*5e1affc2SSatish Balay tao->data = (void*)cgP; 326*5e1affc2SSatish Balay cgP->eta = 0.1; 327*5e1affc2SSatish Balay cgP->delta_min = 1e-7; 328*5e1affc2SSatish Balay cgP->delta_max = 100; 329*5e1affc2SSatish Balay cgP->cg_type = CG_PolakRibierePlus; 330*5e1affc2SSatish Balay 331*5e1affc2SSatish Balay PetscFunctionReturn(0); 332*5e1affc2SSatish Balay } 333*5e1affc2SSatish Balay 334*5e1affc2SSatish Balay EXTERN_C_END 335*5e1affc2SSatish Balay 336*5e1affc2SSatish Balay 337*5e1affc2SSatish Balay 338*5e1affc2SSatish Balay 339*5e1affc2SSatish Balay 340