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 135e1affc2SSatish Balay #undef __FUNCT__ 145e1affc2SSatish Balay #define __FUNCT__ "TaoSolve_CG" 15441846f8SBarry Smith static PetscErrorCode TaoSolve_CG(Tao tao) 165e1affc2SSatish Balay { 175e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG*)tao->data; 185e1affc2SSatish Balay PetscErrorCode ierr; 19*e4cb33bbSBarry Smith TaoConvergedReason reason = TAO_CONTINUE_ITERATING; 20*e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 215e1affc2SSatish Balay PetscReal step=1.0,f,gnorm,gnorm2,delta,gd,ginner,beta; 225e1affc2SSatish Balay PetscReal gd_old,gnorm2_old,f_old; 235e1affc2SSatish Balay PetscInt iter=0; 245e1affc2SSatish Balay 255e1affc2SSatish Balay PetscFunctionBegin; 265e1affc2SSatish Balay if (tao->XL || tao->XU || tao->ops->computebounds) { 275e1affc2SSatish Balay ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by cg algorithm\n");CHKERRQ(ierr); 285e1affc2SSatish Balay } 295e1affc2SSatish Balay 305e1affc2SSatish Balay /* Check convergence criteria */ 315e1affc2SSatish Balay ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr); 325e1affc2SSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 3387f595a5SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 345e1affc2SSatish Balay 355e1affc2SSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr); 3687f595a5SBarry Smith if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 375e1affc2SSatish Balay 385e1affc2SSatish Balay /* Set initial direction to -gradient */ 395e1affc2SSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 405e1affc2SSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 415e1affc2SSatish Balay gnorm2 = gnorm*gnorm; 425e1affc2SSatish Balay 435e1affc2SSatish Balay /* Set initial scaling for the function */ 445e1affc2SSatish Balay if (f != 0.0) { 455e1affc2SSatish Balay delta = 2.0*PetscAbsScalar(f) / gnorm2; 465e1affc2SSatish Balay delta = PetscMax(delta,cgP->delta_min); 475e1affc2SSatish Balay delta = PetscMin(delta,cgP->delta_max); 485e1affc2SSatish Balay } else { 495e1affc2SSatish Balay delta = 2.0 / gnorm2; 505e1affc2SSatish Balay delta = PetscMax(delta,cgP->delta_min); 515e1affc2SSatish Balay delta = PetscMin(delta,cgP->delta_max); 525e1affc2SSatish Balay } 535e1affc2SSatish Balay /* Set counter for gradient and reset steps */ 545e1affc2SSatish Balay cgP->ngradsteps = 0; 555e1affc2SSatish Balay cgP->nresetsteps = 0; 565e1affc2SSatish Balay 575e1affc2SSatish Balay while (1) { 585e1affc2SSatish Balay /* Save the current gradient information */ 595e1affc2SSatish Balay f_old = f; 605e1affc2SSatish Balay gnorm2_old = gnorm2; 615e1affc2SSatish Balay ierr = VecCopy(tao->solution, cgP->X_old);CHKERRQ(ierr); 625e1affc2SSatish Balay ierr = VecCopy(tao->gradient, cgP->G_old);CHKERRQ(ierr); 635e1affc2SSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 645e1affc2SSatish Balay if ((gd >= 0) || PetscIsInfOrNanReal(gd)) { 655e1affc2SSatish Balay ++cgP->ngradsteps; 665e1affc2SSatish Balay if (f != 0.0) { 675e1affc2SSatish Balay delta = 2.0*PetscAbsScalar(f) / gnorm2; 685e1affc2SSatish Balay delta = PetscMax(delta,cgP->delta_min); 695e1affc2SSatish Balay delta = PetscMin(delta,cgP->delta_max); 705e1affc2SSatish Balay } else { 715e1affc2SSatish Balay delta = 2.0 / gnorm2; 725e1affc2SSatish Balay delta = PetscMax(delta,cgP->delta_min); 735e1affc2SSatish Balay delta = PetscMin(delta,cgP->delta_max); 745e1affc2SSatish Balay } 755e1affc2SSatish Balay 765e1affc2SSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 775e1affc2SSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 785e1affc2SSatish Balay } 795e1affc2SSatish Balay 805e1affc2SSatish Balay /* Search direction for improving point */ 815e1affc2SSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,delta); 825e1affc2SSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 835e1affc2SSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 8487f595a5SBarry Smith if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 855e1affc2SSatish Balay /* Linesearch failed */ 865e1affc2SSatish Balay /* Reset factors and use scaled gradient step */ 875e1affc2SSatish Balay ++cgP->nresetsteps; 885e1affc2SSatish Balay f = f_old; 895e1affc2SSatish Balay gnorm2 = gnorm2_old; 905e1affc2SSatish Balay ierr = VecCopy(cgP->X_old, tao->solution);CHKERRQ(ierr); 915e1affc2SSatish Balay ierr = VecCopy(cgP->G_old, tao->gradient);CHKERRQ(ierr); 925e1affc2SSatish Balay 935e1affc2SSatish Balay if (f != 0.0) { 945e1affc2SSatish Balay delta = 2.0*PetscAbsScalar(f) / gnorm2; 955e1affc2SSatish Balay delta = PetscMax(delta,cgP->delta_min); 965e1affc2SSatish Balay delta = PetscMin(delta,cgP->delta_max); 975e1affc2SSatish Balay } else { 985e1affc2SSatish Balay delta = 2.0 / gnorm2; 995e1affc2SSatish Balay delta = PetscMax(delta,cgP->delta_min); 1005e1affc2SSatish Balay delta = PetscMin(delta,cgP->delta_max); 1015e1affc2SSatish Balay } 1025e1affc2SSatish Balay 1035e1affc2SSatish Balay ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 1045e1affc2SSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 1055e1affc2SSatish Balay 1065e1affc2SSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,delta); 1075e1affc2SSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 1085e1affc2SSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 1095e1affc2SSatish Balay 11087f595a5SBarry Smith if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 1115e1affc2SSatish Balay /* Linesearch failed again */ 1125e1affc2SSatish Balay /* switch to unscaled gradient */ 1135e1affc2SSatish Balay f = f_old; 1145e1affc2SSatish Balay gnorm2 = gnorm2_old; 1155e1affc2SSatish Balay ierr = VecCopy(cgP->X_old, tao->solution);CHKERRQ(ierr); 1165e1affc2SSatish Balay ierr = VecCopy(cgP->G_old, tao->gradient);CHKERRQ(ierr); 1175e1affc2SSatish Balay delta = 1.0; 1185e1affc2SSatish Balay ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr); 1195e1affc2SSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 1205e1affc2SSatish Balay 1215e1affc2SSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,delta); 1225e1affc2SSatish Balay ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 1235e1affc2SSatish Balay ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 1245e1affc2SSatish Balay if (ls_status != TAOLINESEARCH_SUCCESS && 1255e1affc2SSatish Balay ls_status != TAOLINESEARCH_SUCCESS_USER) { 1265e1affc2SSatish Balay 1275e1affc2SSatish Balay /* Line search failed for last time -- give up */ 1285e1affc2SSatish Balay f = f_old; 1295e1affc2SSatish Balay gnorm2 = gnorm2_old; 1305e1affc2SSatish Balay ierr = VecCopy(cgP->X_old, tao->solution);CHKERRQ(ierr); 1315e1affc2SSatish Balay ierr = VecCopy(cgP->G_old, tao->gradient);CHKERRQ(ierr); 1325e1affc2SSatish Balay step = 0.0; 1335e1affc2SSatish Balay reason = TAO_DIVERGED_LS_FAILURE; 1345e1affc2SSatish Balay tao->reason = TAO_DIVERGED_LS_FAILURE; 1355e1affc2SSatish Balay } 1365e1affc2SSatish Balay } 1375e1affc2SSatish Balay } 1385e1affc2SSatish Balay 1395e1affc2SSatish Balay /* Check for bad value */ 1405e1affc2SSatish Balay ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 14187f595a5SBarry Smith if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1,"User-provided compute function generated Inf or NaN"); 1425e1affc2SSatish Balay 1435e1affc2SSatish Balay /* Check for termination */ 1445e1affc2SSatish Balay gnorm2 =gnorm * gnorm; 1455e1affc2SSatish Balay iter++; 1465e1affc2SSatish Balay ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr); 1475e1affc2SSatish Balay if (reason != TAO_CONTINUE_ITERATING) { 1485e1affc2SSatish Balay break; 1495e1affc2SSatish Balay } 1505e1affc2SSatish Balay 1515e1affc2SSatish Balay /* Check for restart condition */ 1525e1affc2SSatish Balay ierr = VecDot(tao->gradient, cgP->G_old, &ginner);CHKERRQ(ierr); 1535e1affc2SSatish Balay if (PetscAbsScalar(ginner) >= cgP->eta * gnorm2) { 1545e1affc2SSatish Balay /* Gradients far from orthognal; use steepest descent direction */ 1555e1affc2SSatish Balay beta = 0.0; 1565e1affc2SSatish Balay } else { 1575e1affc2SSatish Balay /* Gradients close to orthogonal; use conjugate gradient formula */ 1585e1affc2SSatish Balay switch (cgP->cg_type) { 1595e1affc2SSatish Balay case CG_FletcherReeves: 1605e1affc2SSatish Balay beta = gnorm2 / gnorm2_old; 1615e1affc2SSatish Balay break; 1625e1affc2SSatish Balay 1635e1affc2SSatish Balay case CG_PolakRibiere: 1645e1affc2SSatish Balay beta = (gnorm2 - ginner) / gnorm2_old; 1655e1affc2SSatish Balay break; 1665e1affc2SSatish Balay 1675e1affc2SSatish Balay case CG_PolakRibierePlus: 1685e1affc2SSatish Balay beta = PetscMax((gnorm2-ginner)/gnorm2_old, 0.0); 1695e1affc2SSatish Balay break; 1705e1affc2SSatish Balay 1715e1affc2SSatish Balay case CG_HestenesStiefel: 1725e1affc2SSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 1735e1affc2SSatish Balay ierr = VecDot(cgP->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 1745e1affc2SSatish Balay beta = (gnorm2 - ginner) / (gd - gd_old); 1755e1affc2SSatish Balay break; 1765e1affc2SSatish Balay 1775e1affc2SSatish Balay case CG_DaiYuan: 1785e1affc2SSatish Balay ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 1795e1affc2SSatish Balay ierr = VecDot(cgP->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 1805e1affc2SSatish Balay beta = gnorm2 / (gd - gd_old); 1815e1affc2SSatish Balay break; 1825e1affc2SSatish Balay 1835e1affc2SSatish Balay default: 1845e1affc2SSatish Balay beta = 0.0; 1855e1affc2SSatish Balay break; 1865e1affc2SSatish Balay } 1875e1affc2SSatish Balay } 1885e1affc2SSatish Balay 1895e1affc2SSatish Balay /* Compute the direction d=-g + beta*d */ 1905e1affc2SSatish Balay ierr = VecAXPBY(tao->stepdirection, -1.0, beta, tao->gradient);CHKERRQ(ierr); 1915e1affc2SSatish Balay 1925e1affc2SSatish Balay /* update initial steplength choice */ 1935e1affc2SSatish Balay delta = 1.0; 1945e1affc2SSatish Balay delta = PetscMax(delta, cgP->delta_min); 1955e1affc2SSatish Balay delta = PetscMin(delta, cgP->delta_max); 1965e1affc2SSatish Balay } 1975e1affc2SSatish Balay PetscFunctionReturn(0); 1985e1affc2SSatish Balay } 1995e1affc2SSatish Balay 2005e1affc2SSatish Balay #undef __FUNCT__ 2015e1affc2SSatish Balay #define __FUNCT__ "TaoSetUp_CG" 202441846f8SBarry Smith static PetscErrorCode TaoSetUp_CG(Tao tao) 2035e1affc2SSatish Balay { 2045e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG*)tao->data; 2055e1affc2SSatish Balay PetscErrorCode ierr; 2065e1affc2SSatish Balay 2075e1affc2SSatish Balay PetscFunctionBegin; 2085e1affc2SSatish Balay if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);} 2095e1affc2SSatish Balay if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); } 2105e1affc2SSatish Balay if (!cgP->X_old) {ierr = VecDuplicate(tao->solution,&cgP->X_old);CHKERRQ(ierr);} 2115e1affc2SSatish Balay if (!cgP->G_old) {ierr = VecDuplicate(tao->gradient,&cgP->G_old);CHKERRQ(ierr); } 2125e1affc2SSatish Balay PetscFunctionReturn(0); 2135e1affc2SSatish Balay } 2145e1affc2SSatish Balay 2155e1affc2SSatish Balay #undef __FUNCT__ 2165e1affc2SSatish Balay #define __FUNCT__ "TaoDestroy_CG" 217441846f8SBarry Smith static PetscErrorCode TaoDestroy_CG(Tao tao) 2185e1affc2SSatish Balay { 2195e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG*) tao->data; 2205e1affc2SSatish Balay PetscErrorCode ierr; 2215e1affc2SSatish Balay 2225e1affc2SSatish Balay PetscFunctionBegin; 2235e1affc2SSatish Balay if (tao->setupcalled) { 2245e1affc2SSatish Balay ierr = VecDestroy(&cgP->X_old);CHKERRQ(ierr); 2255e1affc2SSatish Balay ierr = VecDestroy(&cgP->G_old);CHKERRQ(ierr); 2265e1affc2SSatish Balay } 2275e1affc2SSatish Balay ierr = TaoLineSearchDestroy(&tao->linesearch);CHKERRQ(ierr); 2285e1affc2SSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 2295e1affc2SSatish Balay PetscFunctionReturn(0); 2305e1affc2SSatish Balay } 2315e1affc2SSatish Balay 2325e1affc2SSatish Balay #undef __FUNCT__ 2335e1affc2SSatish Balay #define __FUNCT__ "TaoSetFromOptions_CG" 234441846f8SBarry Smith static PetscErrorCode TaoSetFromOptions_CG(Tao tao) 2355e1affc2SSatish Balay { 2365e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG*)tao->data; 2375e1affc2SSatish Balay PetscErrorCode ierr; 2385e1affc2SSatish Balay 2395e1affc2SSatish Balay PetscFunctionBegin; 2405e1affc2SSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 2415e1affc2SSatish Balay ierr = PetscOptionsHead("Nonlinear Conjugate Gradient method for unconstrained optimization");CHKERRQ(ierr); 2425e1affc2SSatish Balay ierr = PetscOptionsReal("-tao_cg_eta","restart tolerance", "", cgP->eta,&cgP->eta,0);CHKERRQ(ierr); 2435e1affc2SSatish Balay ierr = PetscOptionsEList("-tao_cg_type","cg formula", "", CG_Table, CG_Types, CG_Table[cgP->cg_type], &cgP->cg_type, 0);CHKERRQ(ierr); 2446c23d075SBarry Smith ierr = PetscOptionsReal("-tao_cg_delta_min","minimum delta value", "", cgP->delta_min,&cgP->delta_min,0);CHKERRQ(ierr); 2456c23d075SBarry Smith ierr = PetscOptionsReal("-tao_cg_delta_max","maximum delta value", "", cgP->delta_max,&cgP->delta_max,0);CHKERRQ(ierr); 2465e1affc2SSatish Balay ierr = PetscOptionsTail();CHKERRQ(ierr); 2475e1affc2SSatish Balay PetscFunctionReturn(0); 2485e1affc2SSatish Balay } 2495e1affc2SSatish Balay 2505e1affc2SSatish Balay #undef __FUNCT__ 2515e1affc2SSatish Balay #define __FUNCT__ "TaoView_CG" 252441846f8SBarry Smith static PetscErrorCode TaoView_CG(Tao tao, PetscViewer viewer) 2535e1affc2SSatish Balay { 2545e1affc2SSatish Balay PetscBool isascii; 2555e1affc2SSatish Balay TAO_CG *cgP = (TAO_CG*)tao->data; 2565e1affc2SSatish Balay PetscErrorCode ierr; 2575e1affc2SSatish Balay 2585e1affc2SSatish Balay PetscFunctionBegin; 2595e1affc2SSatish Balay ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 2605e1affc2SSatish Balay if (isascii) { 2615e1affc2SSatish Balay ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 2625e1affc2SSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", CG_Table[cgP->cg_type]);CHKERRQ(ierr); 2635e1affc2SSatish Balay ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", cgP->ngradsteps);CHKERRQ(ierr); 2645e1affc2SSatish Balay ierr= PetscViewerASCIIPrintf(viewer, "Reset steps: %D\n", cgP->nresetsteps);CHKERRQ(ierr); 2655e1affc2SSatish Balay ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2665e1affc2SSatish Balay } 2675e1affc2SSatish Balay PetscFunctionReturn(0); 2685e1affc2SSatish Balay } 2695e1affc2SSatish Balay 2705e1affc2SSatish Balay 2715e1affc2SSatish Balay EXTERN_C_BEGIN 2725e1affc2SSatish Balay #undef __FUNCT__ 2735e1affc2SSatish Balay #define __FUNCT__ "TaoCreate_CG" 274441846f8SBarry Smith PetscErrorCode TaoCreate_CG(Tao tao) 2755e1affc2SSatish Balay { 2765e1affc2SSatish Balay TAO_CG *cgP; 2775e1affc2SSatish Balay const char *morethuente_type = TAOLINESEARCH_MT; 2785e1affc2SSatish Balay PetscErrorCode ierr; 27987f595a5SBarry Smith 2805e1affc2SSatish Balay PetscFunctionBegin; 2815e1affc2SSatish Balay tao->ops->setup = TaoSetUp_CG; 2825e1affc2SSatish Balay tao->ops->solve = TaoSolve_CG; 2835e1affc2SSatish Balay tao->ops->view = TaoView_CG; 2845e1affc2SSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_CG; 2855e1affc2SSatish Balay tao->ops->destroy = TaoDestroy_CG; 2865e1affc2SSatish Balay 2875e1affc2SSatish Balay tao->max_it = 2000; 2885e1affc2SSatish Balay tao->max_funcs = 4000; 2895e1affc2SSatish Balay tao->fatol = 1e-4; 2905e1affc2SSatish Balay tao->frtol = 1e-4; 2915e1affc2SSatish Balay 2925e1affc2SSatish Balay /* Note: nondefault values should be used for nonlinear conjugate gradient */ 2935e1affc2SSatish Balay /* method. In particular, gtol should be less that 0.5; the value used in */ 2945e1affc2SSatish Balay /* Nocedal and Wright is 0.10. We use the default values for the */ 2955e1affc2SSatish Balay /* linesearch because it seems to work better. */ 2965e1affc2SSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 2975e1affc2SSatish Balay ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 298441846f8SBarry Smith ierr = TaoLineSearchUseTaoRoutines(tao->linesearch, tao);CHKERRQ(ierr); 2995e1affc2SSatish Balay 3005e1affc2SSatish Balay ierr = PetscNewLog(tao,&cgP);CHKERRQ(ierr); 3015e1affc2SSatish Balay tao->data = (void*)cgP; 3025e1affc2SSatish Balay cgP->eta = 0.1; 3035e1affc2SSatish Balay cgP->delta_min = 1e-7; 3045e1affc2SSatish Balay cgP->delta_max = 100; 3055e1affc2SSatish Balay cgP->cg_type = CG_PolakRibierePlus; 3065e1affc2SSatish Balay PetscFunctionReturn(0); 3075e1affc2SSatish Balay } 3085e1affc2SSatish Balay EXTERN_C_END 309