1ba92ff59SBarry Smith #include <petsctaolinesearch.h>
2aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/cg/taocg.h>
35e1affc2SSatish Balay
45e1affc2SSatish Balay #define CG_FletcherReeves 0
55e1affc2SSatish Balay #define CG_PolakRibiere 1
65e1affc2SSatish Balay #define CG_PolakRibierePlus 2
75e1affc2SSatish Balay #define CG_HestenesStiefel 3
85e1affc2SSatish Balay #define CG_DaiYuan 4
95e1affc2SSatish Balay #define CG_Types 5
105e1affc2SSatish Balay
1187f595a5SBarry Smith static const char *CG_Table[64] = {"fr", "pr", "prp", "hs", "dy"};
125e1affc2SSatish Balay
TaoSolve_CG(Tao tao)13d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_CG(Tao tao)
14d71ae5a4SJacob Faibussowitsch {
155e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG *)tao->data;
16e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
175e1affc2SSatish Balay PetscReal step = 1.0, f, gnorm, gnorm2, delta, gd, ginner, beta;
185e1affc2SSatish Balay PetscReal gd_old, gnorm2_old, f_old;
195e1affc2SSatish Balay
205e1affc2SSatish Balay PetscFunctionBegin;
2148a46eb9SPierre Jolivet if (tao->XL || tao->XU || tao->ops->computebounds) PetscCall(PetscInfo(tao, "WARNING: Variable bounds have been set but will be ignored by cg algorithm\n"));
225e1affc2SSatish Balay
235e1affc2SSatish Balay /* Check convergence criteria */
249566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
259566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm));
2676c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
275e1affc2SSatish Balay
283ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING;
299566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
309566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
31dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
323ba16761SJacob Faibussowitsch if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
335e1affc2SSatish Balay
345e1affc2SSatish Balay /* Set initial direction to -gradient */
359566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, tao->stepdirection));
369566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0));
375e1affc2SSatish Balay gnorm2 = gnorm * gnorm;
385e1affc2SSatish Balay
395e1affc2SSatish Balay /* Set initial scaling for the function */
405e1affc2SSatish Balay if (f != 0.0) {
415e1affc2SSatish Balay delta = 2.0 * PetscAbsScalar(f) / gnorm2;
425e1affc2SSatish Balay delta = PetscMax(delta, cgP->delta_min);
435e1affc2SSatish Balay delta = PetscMin(delta, cgP->delta_max);
445e1affc2SSatish Balay } else {
455e1affc2SSatish Balay delta = 2.0 / gnorm2;
465e1affc2SSatish Balay delta = PetscMax(delta, cgP->delta_min);
475e1affc2SSatish Balay delta = PetscMin(delta, cgP->delta_max);
485e1affc2SSatish Balay }
495e1affc2SSatish Balay /* Set counter for gradient and reset steps */
505e1affc2SSatish Balay cgP->ngradsteps = 0;
515e1affc2SSatish Balay cgP->nresetsteps = 0;
525e1affc2SSatish Balay
535e1affc2SSatish Balay while (1) {
54e1e80dc8SAlp Dener /* Call general purpose update function */
55270bebe6SStefano Zampini if (tao->ops->update) {
56270bebe6SStefano Zampini PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
57270bebe6SStefano Zampini PetscCall(TaoComputeObjective(tao, tao->solution, &f));
58270bebe6SStefano Zampini }
595e1affc2SSatish Balay /* Save the current gradient information */
605e1affc2SSatish Balay f_old = f;
615e1affc2SSatish Balay gnorm2_old = gnorm2;
629566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, cgP->X_old));
639566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, cgP->G_old));
649566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
655e1affc2SSatish Balay if ((gd >= 0) || PetscIsInfOrNanReal(gd)) {
665e1affc2SSatish Balay ++cgP->ngradsteps;
675e1affc2SSatish Balay if (f != 0.0) {
685e1affc2SSatish Balay delta = 2.0 * PetscAbsScalar(f) / gnorm2;
695e1affc2SSatish Balay delta = PetscMax(delta, cgP->delta_min);
705e1affc2SSatish Balay delta = PetscMin(delta, cgP->delta_max);
715e1affc2SSatish Balay } else {
725e1affc2SSatish Balay delta = 2.0 / gnorm2;
735e1affc2SSatish Balay delta = PetscMax(delta, cgP->delta_min);
745e1affc2SSatish Balay delta = PetscMin(delta, cgP->delta_max);
755e1affc2SSatish Balay }
765e1affc2SSatish Balay
779566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, tao->stepdirection));
789566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0));
795e1affc2SSatish Balay }
805e1affc2SSatish Balay
815e1affc2SSatish Balay /* Search direction for improving point */
829566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, delta));
839566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status));
849566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao));
8587f595a5SBarry Smith if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
865e1affc2SSatish Balay /* Linesearch failed */
875e1affc2SSatish Balay /* Reset factors and use scaled gradient step */
885e1affc2SSatish Balay ++cgP->nresetsteps;
895e1affc2SSatish Balay f = f_old;
905e1affc2SSatish Balay gnorm2 = gnorm2_old;
919566063dSJacob Faibussowitsch PetscCall(VecCopy(cgP->X_old, tao->solution));
929566063dSJacob Faibussowitsch PetscCall(VecCopy(cgP->G_old, tao->gradient));
935e1affc2SSatish Balay
945e1affc2SSatish Balay if (f != 0.0) {
955e1affc2SSatish Balay delta = 2.0 * PetscAbsScalar(f) / gnorm2;
965e1affc2SSatish Balay delta = PetscMax(delta, cgP->delta_min);
975e1affc2SSatish Balay delta = PetscMin(delta, cgP->delta_max);
985e1affc2SSatish Balay } else {
995e1affc2SSatish Balay delta = 2.0 / gnorm2;
1005e1affc2SSatish Balay delta = PetscMax(delta, cgP->delta_min);
1015e1affc2SSatish Balay delta = PetscMin(delta, cgP->delta_max);
1025e1affc2SSatish Balay }
1035e1affc2SSatish Balay
1049566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, tao->stepdirection));
1059566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0));
1065e1affc2SSatish Balay
1079566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, delta));
1089566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status));
1099566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao));
1105e1affc2SSatish Balay
11187f595a5SBarry Smith if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
1125e1affc2SSatish Balay /* Linesearch failed again */
1135e1affc2SSatish Balay /* switch to unscaled gradient */
1145e1affc2SSatish Balay f = f_old;
1159566063dSJacob Faibussowitsch PetscCall(VecCopy(cgP->X_old, tao->solution));
1169566063dSJacob Faibussowitsch PetscCall(VecCopy(cgP->G_old, tao->gradient));
1175e1affc2SSatish Balay delta = 1.0;
118*79813928SHan Liu PetscCall(VecCopy(tao->gradient, tao->stepdirection));
1199566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0));
1205e1affc2SSatish Balay
1219566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, delta));
1229566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status));
1239566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao));
1248caf6e8cSBarry Smith if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
1255e1affc2SSatish Balay /* Line search failed for last time -- give up */
1265e1affc2SSatish Balay f = f_old;
1279566063dSJacob Faibussowitsch PetscCall(VecCopy(cgP->X_old, tao->solution));
1289566063dSJacob Faibussowitsch PetscCall(VecCopy(cgP->G_old, tao->gradient));
1295e1affc2SSatish Balay step = 0.0;
1305e1affc2SSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE;
1315e1affc2SSatish Balay }
1325e1affc2SSatish Balay }
1335e1affc2SSatish Balay }
1345e1affc2SSatish Balay
1355e1affc2SSatish Balay /* Check for bad value */
1369566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient, NORM_2, &gnorm));
13776c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User-provided compute function generated infinity or NaN");
1385e1affc2SSatish Balay
1395e1affc2SSatish Balay /* Check for termination */
1405e1affc2SSatish Balay gnorm2 = gnorm * gnorm;
1418931d482SJason Sarich tao->niter++;
1429566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
1439566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
144dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
145ad540459SPierre Jolivet if (tao->reason != TAO_CONTINUE_ITERATING) break;
1465e1affc2SSatish Balay
1475e1affc2SSatish Balay /* Check for restart condition */
1489566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, cgP->G_old, &ginner));
1495e1affc2SSatish Balay if (PetscAbsScalar(ginner) >= cgP->eta * gnorm2) {
150a5b23f4aSJose E. Roman /* Gradients far from orthogonal; use steepest descent direction */
1515e1affc2SSatish Balay beta = 0.0;
1525e1affc2SSatish Balay } else {
1535e1affc2SSatish Balay /* Gradients close to orthogonal; use conjugate gradient formula */
1545e1affc2SSatish Balay switch (cgP->cg_type) {
155d71ae5a4SJacob Faibussowitsch case CG_FletcherReeves:
156d71ae5a4SJacob Faibussowitsch beta = gnorm2 / gnorm2_old;
157d71ae5a4SJacob Faibussowitsch break;
1585e1affc2SSatish Balay
159d71ae5a4SJacob Faibussowitsch case CG_PolakRibiere:
160d71ae5a4SJacob Faibussowitsch beta = (gnorm2 - ginner) / gnorm2_old;
161d71ae5a4SJacob Faibussowitsch break;
1625e1affc2SSatish Balay
163d71ae5a4SJacob Faibussowitsch case CG_PolakRibierePlus:
164d71ae5a4SJacob Faibussowitsch beta = PetscMax((gnorm2 - ginner) / gnorm2_old, 0.0);
165d71ae5a4SJacob Faibussowitsch break;
1665e1affc2SSatish Balay
1675e1affc2SSatish Balay case CG_HestenesStiefel:
1689566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
1699566063dSJacob Faibussowitsch PetscCall(VecDot(cgP->G_old, tao->stepdirection, &gd_old));
1705e1affc2SSatish Balay beta = (gnorm2 - ginner) / (gd - gd_old);
1715e1affc2SSatish Balay break;
1725e1affc2SSatish Balay
1735e1affc2SSatish Balay case CG_DaiYuan:
1749566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gd));
1759566063dSJacob Faibussowitsch PetscCall(VecDot(cgP->G_old, tao->stepdirection, &gd_old));
1765e1affc2SSatish Balay beta = gnorm2 / (gd - gd_old);
1775e1affc2SSatish Balay break;
1785e1affc2SSatish Balay
179d71ae5a4SJacob Faibussowitsch default:
180d71ae5a4SJacob Faibussowitsch beta = 0.0;
181d71ae5a4SJacob Faibussowitsch break;
1825e1affc2SSatish Balay }
1835e1affc2SSatish Balay }
1845e1affc2SSatish Balay
1855e1affc2SSatish Balay /* Compute the direction d=-g + beta*d */
1869566063dSJacob Faibussowitsch PetscCall(VecAXPBY(tao->stepdirection, -1.0, beta, tao->gradient));
1875e1affc2SSatish Balay
1885e1affc2SSatish Balay /* update initial steplength choice */
1895e1affc2SSatish Balay delta = 1.0;
1905e1affc2SSatish Balay delta = PetscMax(delta, cgP->delta_min);
1915e1affc2SSatish Balay delta = PetscMin(delta, cgP->delta_max);
1925e1affc2SSatish Balay }
1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1945e1affc2SSatish Balay }
1955e1affc2SSatish Balay
TaoSetUp_CG(Tao tao)196d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_CG(Tao tao)
197d71ae5a4SJacob Faibussowitsch {
1985e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG *)tao->data;
1995e1affc2SSatish Balay
2005e1affc2SSatish Balay PetscFunctionBegin;
2019566063dSJacob Faibussowitsch if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
2029566063dSJacob Faibussowitsch if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
2039566063dSJacob Faibussowitsch if (!cgP->X_old) PetscCall(VecDuplicate(tao->solution, &cgP->X_old));
2049566063dSJacob Faibussowitsch if (!cgP->G_old) PetscCall(VecDuplicate(tao->gradient, &cgP->G_old));
2053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2065e1affc2SSatish Balay }
2075e1affc2SSatish Balay
TaoDestroy_CG(Tao tao)208d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_CG(Tao tao)
209d71ae5a4SJacob Faibussowitsch {
2105e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG *)tao->data;
2115e1affc2SSatish Balay
2125e1affc2SSatish Balay PetscFunctionBegin;
2135e1affc2SSatish Balay if (tao->setupcalled) {
2149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cgP->X_old));
2159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cgP->G_old));
2165e1affc2SSatish Balay }
2179566063dSJacob Faibussowitsch PetscCall(TaoLineSearchDestroy(&tao->linesearch));
2189566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data));
2193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2205e1affc2SSatish Balay }
2215e1affc2SSatish Balay
TaoSetFromOptions_CG(Tao tao,PetscOptionItems PetscOptionsObject)222ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_CG(Tao tao, PetscOptionItems PetscOptionsObject)
223d71ae5a4SJacob Faibussowitsch {
2245e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG *)tao->data;
2255e1affc2SSatish Balay
2265e1affc2SSatish Balay PetscFunctionBegin;
2279566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
228d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Nonlinear Conjugate Gradient method for unconstrained optimization");
2299566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_cg_eta", "restart tolerance", "", cgP->eta, &cgP->eta, NULL));
2309566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-tao_cg_type", "cg formula", "", CG_Table, CG_Types, CG_Table[cgP->cg_type], &cgP->cg_type, NULL));
2319566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_cg_delta_min", "minimum delta value", "", cgP->delta_min, &cgP->delta_min, NULL));
2329566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_cg_delta_max", "maximum delta value", "", cgP->delta_max, &cgP->delta_max, NULL));
233d0609cedSBarry Smith PetscOptionsHeadEnd();
2343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2355e1affc2SSatish Balay }
2365e1affc2SSatish Balay
TaoView_CG(Tao tao,PetscViewer viewer)237d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_CG(Tao tao, PetscViewer viewer)
238d71ae5a4SJacob Faibussowitsch {
2395e1affc2SSatish Balay PetscBool isascii;
2405e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG *)tao->data;
2415e1affc2SSatish Balay
2425e1affc2SSatish Balay PetscFunctionBegin;
2439566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2445e1affc2SSatish Balay if (isascii) {
2459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer));
2469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", CG_Table[cgP->cg_type]));
24763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", cgP->ngradsteps));
24863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Reset steps: %" PetscInt_FMT "\n", cgP->nresetsteps));
2499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
2505e1affc2SSatish Balay }
2513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2525e1affc2SSatish Balay }
2534aa34175SJason Sarich
2544aa34175SJason Sarich /*MC
2554aa34175SJason Sarich TAOCG - Nonlinear conjugate gradient method is an extension of the
2564aa34175SJason Sarich nonlinear conjugate gradient solver for nonlinear optimization.
2574aa34175SJason Sarich
2584aa34175SJason Sarich Options Database Keys:
2594aa34175SJason Sarich + -tao_cg_eta <r> - restart tolerance
2604aa34175SJason Sarich . -tao_cg_type <taocg_type> - cg formula
2614aa34175SJason Sarich . -tao_cg_delta_min <r> - minimum delta value
2624aa34175SJason Sarich - -tao_cg_delta_max <r> - maximum delta value
2634aa34175SJason Sarich
2644aa34175SJason Sarich Notes:
2654aa34175SJason Sarich CG formulas are:
2664aa34175SJason Sarich "fr" - Fletcher-Reeves
2674aa34175SJason Sarich "pr" - Polak-Ribiere
2684aa34175SJason Sarich "prp" - Polak-Ribiere-Plus
2694aa34175SJason Sarich "hs" - Hestenes-Steifel
2704aa34175SJason Sarich "dy" - Dai-Yuan
2711eb8069cSJason Sarich Level: beginner
2724aa34175SJason Sarich M*/
2735e1affc2SSatish Balay
TaoCreate_CG(Tao tao)274d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_CG(Tao tao)
275d71ae5a4SJacob Faibussowitsch {
2765e1affc2SSatish Balay TAO_CG *cgP;
2778caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT;
27887f595a5SBarry Smith
2795e1affc2SSatish Balay PetscFunctionBegin;
2805e1affc2SSatish Balay tao->ops->setup = TaoSetUp_CG;
2815e1affc2SSatish Balay tao->ops->solve = TaoSolve_CG;
2825e1affc2SSatish Balay tao->ops->view = TaoView_CG;
2835e1affc2SSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_CG;
2845e1affc2SSatish Balay tao->ops->destroy = TaoDestroy_CG;
2855e1affc2SSatish Balay
2866552cf8aSJason Sarich /* Override default settings (unless already changed) */
287606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao));
288606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_it, 2000);
289606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_funcs, 4000);
2905e1affc2SSatish Balay
2915e1affc2SSatish Balay /* Note: nondefault values should be used for nonlinear conjugate gradient */
2925e1affc2SSatish Balay /* method. In particular, gtol should be less that 0.5; the value used in */
2935e1affc2SSatish Balay /* Nocedal and Wright is 0.10. We use the default values for the */
2945e1affc2SSatish Balay /* linesearch because it seems to work better. */
2959566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
2969566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
2979566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
2989566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
2999566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
3005e1affc2SSatish Balay
3014dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&cgP));
3025e1affc2SSatish Balay tao->data = (void *)cgP;
3035e1affc2SSatish Balay cgP->eta = 0.1;
3045e1affc2SSatish Balay cgP->delta_min = 1e-7;
3055e1affc2SSatish Balay cgP->delta_max = 100;
3065e1affc2SSatish Balay cgP->cg_type = CG_PolakRibierePlus;
3073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3085e1affc2SSatish Balay }
309