10a844d1aSPeter Brune #include <../src/snes/impls/ncg/snesncgimpl.h> /*I "petscsnes.h" I*/ 26a6fc655SJed Brown const char *const SNESNCGTypes[] = {"FR","PRP","HS","DY","CD","SNESNCGType","SNES_NCG_",0}; 3fef7b6d8SPeter Brune 4fef7b6d8SPeter Brune #undef __FUNCT__ 5fef7b6d8SPeter Brune #define __FUNCT__ "SNESReset_NCG" 6fef7b6d8SPeter Brune PetscErrorCode SNESReset_NCG(SNES snes) 7fef7b6d8SPeter Brune { 8fef7b6d8SPeter Brune PetscFunctionBegin; 9fef7b6d8SPeter Brune PetscFunctionReturn(0); 10fef7b6d8SPeter Brune } 11fef7b6d8SPeter Brune 121a4f838cSPeter Brune #define SNESLINESEARCHNCGLINEAR "linear" 13bbd5d0b3SPeter Brune 14fef7b6d8SPeter Brune /* 15fef7b6d8SPeter Brune SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG(). 16fef7b6d8SPeter Brune 17fef7b6d8SPeter Brune Input Parameter: 18fef7b6d8SPeter Brune . snes - the SNES context 19fef7b6d8SPeter Brune 20fef7b6d8SPeter Brune Application Interface Routine: SNESDestroy() 21fef7b6d8SPeter Brune */ 22fef7b6d8SPeter Brune #undef __FUNCT__ 23fef7b6d8SPeter Brune #define __FUNCT__ "SNESDestroy_NCG" 24fef7b6d8SPeter Brune PetscErrorCode SNESDestroy_NCG(SNES snes) 25fef7b6d8SPeter Brune { 26fef7b6d8SPeter Brune PetscErrorCode ierr; 27fef7b6d8SPeter Brune 28fef7b6d8SPeter Brune PetscFunctionBegin; 29fef7b6d8SPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 30fef7b6d8SPeter Brune PetscFunctionReturn(0); 31fef7b6d8SPeter Brune } 32fef7b6d8SPeter Brune 33fef7b6d8SPeter Brune /* 34fef7b6d8SPeter Brune SNESSetUp_NCG - Sets up the internal data structures for the later use 35fef7b6d8SPeter Brune of the SNESNCG nonlinear solver. 36fef7b6d8SPeter Brune 37fef7b6d8SPeter Brune Input Parameters: 38fef7b6d8SPeter Brune + snes - the SNES context 39fef7b6d8SPeter Brune - x - the solution vector 40fef7b6d8SPeter Brune 41fef7b6d8SPeter Brune Application Interface Routine: SNESSetUp() 42fef7b6d8SPeter Brune */ 43bbd5d0b3SPeter Brune 448cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch); 45bbd5d0b3SPeter Brune 46fef7b6d8SPeter Brune #undef __FUNCT__ 47fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetUp_NCG" 48fef7b6d8SPeter Brune PetscErrorCode SNESSetUp_NCG(SNES snes) 49fef7b6d8SPeter Brune { 50fef7b6d8SPeter Brune PetscErrorCode ierr; 51fef7b6d8SPeter Brune 52fef7b6d8SPeter Brune PetscFunctionBegin; 53fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes,2);CHKERRQ(ierr); 54bbd5d0b3SPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 55*a71552e2SPeter Brune if (snes->pcside == PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning"); 56bdf89e91SBarry Smith ierr = SNESLineSearchRegister(SNESLINESEARCHNCGLINEAR, SNESLineSearchCreate_NCGLinear);CHKERRQ(ierr); 57fef7b6d8SPeter Brune PetscFunctionReturn(0); 58fef7b6d8SPeter Brune } 59fef7b6d8SPeter Brune /* 6005b53524SPeter Brune SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method. 61fef7b6d8SPeter Brune 62fef7b6d8SPeter Brune Input Parameter: 63fef7b6d8SPeter Brune . snes - the SNES context 64fef7b6d8SPeter Brune 65fef7b6d8SPeter Brune Application Interface Routine: SNESSetFromOptions() 66fef7b6d8SPeter Brune */ 67fef7b6d8SPeter Brune #undef __FUNCT__ 68fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NCG" 69fef7b6d8SPeter Brune static PetscErrorCode SNESSetFromOptions_NCG(SNES snes) 70fef7b6d8SPeter Brune { 71dfb256c7SPeter Brune SNES_NCG *ncg = (SNES_NCG*)snes->data; 72fef7b6d8SPeter Brune PetscErrorCode ierr; 730a844d1aSPeter Brune PetscBool debug; 74f1c6b773SPeter Brune SNESLineSearch linesearch; 75a11d90f7SPeter Brune SNESNCGType ncgtype=ncg->type; 76fef7b6d8SPeter Brune 77fef7b6d8SPeter Brune PetscFunctionBegin; 78fef7b6d8SPeter Brune ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr); 790298fd71SBarry Smith ierr = PetscOptionsBool("-snes_ncg_monitor","Monitor NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);CHKERRQ(ierr); 800298fd71SBarry Smith ierr = PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);CHKERRQ(ierr); 810a844d1aSPeter Brune ierr = SNESNCGSetType(snes, ncgtype);CHKERRQ(ierr); 82dfb256c7SPeter Brune if (debug) { 83ce94432eSBarry Smith ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 84dfb256c7SPeter Brune } 85fef7b6d8SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 869e764e56SPeter Brune if (!snes->linesearch) { 877601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 889e764e56SPeter Brune if (!snes->pc) { 891a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 909e764e56SPeter Brune } else { 911a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 929e764e56SPeter Brune } 939e764e56SPeter Brune } 94fef7b6d8SPeter Brune PetscFunctionReturn(0); 95fef7b6d8SPeter Brune } 96fef7b6d8SPeter Brune 97fef7b6d8SPeter Brune /* 98fef7b6d8SPeter Brune SNESView_NCG - Prints info from the SNESNCG data structure. 99fef7b6d8SPeter Brune 100fef7b6d8SPeter Brune Input Parameters: 101fef7b6d8SPeter Brune + SNES - the SNES context 102fef7b6d8SPeter Brune - viewer - visualization context 103fef7b6d8SPeter Brune 104fef7b6d8SPeter Brune Application Interface Routine: SNESView() 105fef7b6d8SPeter Brune */ 106fef7b6d8SPeter Brune #undef __FUNCT__ 107fef7b6d8SPeter Brune #define __FUNCT__ "SNESView_NCG" 108fef7b6d8SPeter Brune static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer) 109fef7b6d8SPeter Brune { 110fef7b6d8SPeter Brune PetscBool iascii; 111fef7b6d8SPeter Brune PetscErrorCode ierr; 112fef7b6d8SPeter Brune 113fef7b6d8SPeter Brune PetscFunctionBegin; 114251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 115fef7b6d8SPeter Brune if (iascii) { 116fef7b6d8SPeter Brune } 117fef7b6d8SPeter Brune PetscFunctionReturn(0); 118fef7b6d8SPeter Brune } 119fef7b6d8SPeter Brune 120fef7b6d8SPeter Brune #undef __FUNCT__ 121f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply_NCGLinear" 122f1c6b773SPeter Brune PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch) 123ea630c6eSPeter Brune { 124ea630c6eSPeter Brune PetscScalar alpha, ptAp; 125bbd5d0b3SPeter Brune Vec X, Y, F, W; 126bbd5d0b3SPeter Brune SNES snes; 127ea630c6eSPeter Brune PetscErrorCode ierr; 128bbd5d0b3SPeter Brune PetscReal *fnorm, *xnorm, *ynorm; 129ea630c6eSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 130ea630c6eSPeter Brune 131ea630c6eSPeter Brune PetscFunctionBegin; 132f1c6b773SPeter Brune ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 133bbd5d0b3SPeter Brune X = linesearch->vec_sol; 134bbd5d0b3SPeter Brune W = linesearch->vec_sol_new; 135bbd5d0b3SPeter Brune F = linesearch->vec_func; 136bbd5d0b3SPeter Brune Y = linesearch->vec_update; 137bbd5d0b3SPeter Brune fnorm = &linesearch->fnorm; 138bbd5d0b3SPeter Brune xnorm = &linesearch->xnorm; 139bbd5d0b3SPeter Brune ynorm = &linesearch->ynorm; 140bbd5d0b3SPeter Brune 141ea630c6eSPeter Brune /* 142ea630c6eSPeter Brune 143d2140ca3SPeter Brune The exact step size for unpreconditioned linear CG is just: 144ea630c6eSPeter Brune alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) 145ea630c6eSPeter Brune */ 146a5892487SPeter Brune ierr = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr); 147ea630c6eSPeter Brune ierr = VecDot(F, F, &alpha);CHKERRQ(ierr); 148ea630c6eSPeter Brune ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr); 149ea630c6eSPeter Brune ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr); 150ea630c6eSPeter Brune alpha = alpha / ptAp; 151ce94432eSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)snes), "alpha: %G\n", PetscRealPart(alpha));CHKERRQ(ierr); 152bbd5d0b3SPeter Brune ierr = VecAXPY(X, alpha, Y);CHKERRQ(ierr); 153bbd5d0b3SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 154bbd5d0b3SPeter Brune 155bbd5d0b3SPeter Brune ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 156bbd5d0b3SPeter Brune ierr = VecNorm(X, NORM_2, xnorm);CHKERRQ(ierr); 157bbd5d0b3SPeter Brune ierr = VecNorm(Y, NORM_2, ynorm);CHKERRQ(ierr); 158ea630c6eSPeter Brune PetscFunctionReturn(0); 159ea630c6eSPeter Brune } 160ea630c6eSPeter Brune 161bbd5d0b3SPeter Brune #undef __FUNCT__ 162f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_NCGLinear" 1638cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch) 164bbd5d0b3SPeter Brune { 165bbd5d0b3SPeter Brune PetscFunctionBegin; 166f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_NCGLinear; 1670298fd71SBarry Smith linesearch->ops->destroy = NULL; 1680298fd71SBarry Smith linesearch->ops->setfromoptions = NULL; 1690298fd71SBarry Smith linesearch->ops->reset = NULL; 1700298fd71SBarry Smith linesearch->ops->view = NULL; 1710298fd71SBarry Smith linesearch->ops->setup = NULL; 172bbd5d0b3SPeter Brune PetscFunctionReturn(0); 173bbd5d0b3SPeter Brune } 174bbd5d0b3SPeter Brune 17588d374b2SPeter Brune #undef __FUNCT__ 1768c40d5fbSBarry Smith #define __FUNCT__ "SNESNCGComputeYtJtF_Private" 17788d374b2SPeter Brune /* 17888d374b2SPeter Brune 17988d374b2SPeter Brune Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian. 18088d374b2SPeter Brune 18188d374b2SPeter Brune */ 18267392de3SBarry Smith PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf) 18367392de3SBarry Smith { 18488d374b2SPeter Brune PetscErrorCode ierr; 18588d374b2SPeter Brune PetscScalar ftf, ftg, fty, h; 18667392de3SBarry Smith 18788d374b2SPeter Brune PetscFunctionBegin; 18888d374b2SPeter Brune ierr = VecDot(F, F, &ftf);CHKERRQ(ierr); 18988d374b2SPeter Brune ierr = VecDot(F, Y, &fty);CHKERRQ(ierr); 19088d374b2SPeter Brune h = 1e-5*fty / fty; 19188d374b2SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 19288d374b2SPeter Brune ierr = VecAXPY(W, -h, Y);CHKERRQ(ierr); /* this is arbitrary */ 19388d374b2SPeter Brune ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 19488d374b2SPeter Brune ierr = VecDot(G, F, &ftg);CHKERRQ(ierr); 19588d374b2SPeter Brune *ytJtf = (ftg - ftf) / h; 19688d374b2SPeter Brune PetscFunctionReturn(0); 19788d374b2SPeter Brune } 19888d374b2SPeter Brune 1990a844d1aSPeter Brune #undef __FUNCT__ 2000a844d1aSPeter Brune #define __FUNCT__ "SNESNCGSetType" 2010a844d1aSPeter Brune /*@ 2020a844d1aSPeter Brune SNESNCGSetType - Sets the conjugate update type for SNESNCG. 2030a844d1aSPeter Brune 2040a844d1aSPeter Brune Logically Collective on SNES 2050a844d1aSPeter Brune 2060a844d1aSPeter Brune Input Parameters: 2070a844d1aSPeter Brune + snes - the iterative context 2080a844d1aSPeter Brune - btype - update type 2090a844d1aSPeter Brune 2100a844d1aSPeter Brune Options Database: 2110a844d1aSPeter Brune . -snes_ncg_type<prp,fr,hs,dy,cd> 2120a844d1aSPeter Brune 2130a844d1aSPeter Brune Level: intermediate 2140a844d1aSPeter Brune 2150a844d1aSPeter Brune SNESNCGSelectTypes: 2160a844d1aSPeter Brune + SNES_NCG_FR - Fletcher-Reeves update 2170a844d1aSPeter Brune . SNES_NCG_PRP - Polak-Ribiere-Polyak update 2180a844d1aSPeter Brune . SNES_NCG_HS - Hestenes-Steifel update 2190a844d1aSPeter Brune . SNES_NCG_DY - Dai-Yuan update 2200a844d1aSPeter Brune - SNES_NCG_CD - Conjugate Descent update 2210a844d1aSPeter Brune 2220a844d1aSPeter Brune Notes: 2230a844d1aSPeter Brune PRP is the default, and the only one that tolerates generalized search directions. 2240a844d1aSPeter Brune 2250a844d1aSPeter Brune .keywords: SNES, SNESNCG, selection, type, set 2260a844d1aSPeter Brune @*/ 2270adebc6cSBarry Smith PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype) 2280adebc6cSBarry Smith { 2290a844d1aSPeter Brune PetscErrorCode ierr; 2306e111a19SKarl Rupp 2310a844d1aSPeter Brune PetscFunctionBegin; 2320a844d1aSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2330a844d1aSPeter Brune ierr = PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));CHKERRQ(ierr); 2340a844d1aSPeter Brune PetscFunctionReturn(0); 2350a844d1aSPeter Brune } 2360a844d1aSPeter Brune 2370a844d1aSPeter Brune #undef __FUNCT__ 2380a844d1aSPeter Brune #define __FUNCT__ "SNESNCGSetType_NCG" 2390adebc6cSBarry Smith PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype) 2400adebc6cSBarry Smith { 2410a844d1aSPeter Brune SNES_NCG *ncg = (SNES_NCG*)snes->data; 2426e111a19SKarl Rupp 2430a844d1aSPeter Brune PetscFunctionBegin; 2440a844d1aSPeter Brune ncg->type = btype; 2450a844d1aSPeter Brune PetscFunctionReturn(0); 2460a844d1aSPeter Brune } 2470a844d1aSPeter Brune 248fef7b6d8SPeter Brune /* 249fef7b6d8SPeter Brune SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method. 250fef7b6d8SPeter Brune 251fef7b6d8SPeter Brune Input Parameters: 252fef7b6d8SPeter Brune . snes - the SNES context 253fef7b6d8SPeter Brune 254fef7b6d8SPeter Brune Output Parameter: 255fef7b6d8SPeter Brune . outits - number of iterations until termination 256fef7b6d8SPeter Brune 257fef7b6d8SPeter Brune Application Interface Routine: SNESSolve() 258fef7b6d8SPeter Brune */ 259fef7b6d8SPeter Brune #undef __FUNCT__ 260fef7b6d8SPeter Brune #define __FUNCT__ "SNESSolve_NCG" 261fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes) 262fef7b6d8SPeter Brune { 263dfb256c7SPeter Brune SNES_NCG *ncg = (SNES_NCG*)snes->data; 264*a71552e2SPeter Brune Vec X,dX,lX,F,Fold; 265bbd5d0b3SPeter Brune PetscReal fnorm, ynorm, xnorm, beta = 0.0; 2665d115551SPeter Brune PetscScalar dXdotF, dXolddotFold, dXdotFold, lXdotF, lXdotFold; 267fef7b6d8SPeter Brune PetscInt maxits, i; 268fef7b6d8SPeter Brune PetscErrorCode ierr; 269fef7b6d8SPeter Brune PetscBool lsSuccess = PETSC_TRUE; 270f1c6b773SPeter Brune SNESLineSearch linesearch; 27188d374b2SPeter Brune 272fef7b6d8SPeter Brune PetscFunctionBegin; 273fef7b6d8SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 274fef7b6d8SPeter Brune 275fef7b6d8SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 276fef7b6d8SPeter Brune X = snes->vec_sol; /* X^n */ 2775d115551SPeter Brune Fold = snes->work[0]; /* The previous iterate of X */ 27888d374b2SPeter Brune dX = snes->work[1]; /* the preconditioned direction */ 279169e2be8SPeter Brune lX = snes->vec_sol_update; /* the conjugate direction */ 280fef7b6d8SPeter Brune F = snes->vec_func; /* residual vector */ 281fef7b6d8SPeter Brune 2827601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 2839e764e56SPeter Brune 284ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 285fef7b6d8SPeter Brune snes->iter = 0; 286fef7b6d8SPeter Brune snes->norm = 0.; 287ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 288fef7b6d8SPeter Brune 289bbd5d0b3SPeter Brune /* compute the initial function and preconditioned update dX */ 290*a71552e2SPeter Brune 291*a71552e2SPeter Brune if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 292*a71552e2SPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,dX);CHKERRQ(ierr); 293*a71552e2SPeter Brune ierr = VecCopy(dX,F);CHKERRQ(ierr); 294*a71552e2SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 295*a71552e2SPeter Brune } else { 296e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 297fef7b6d8SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 298fef7b6d8SPeter Brune if (snes->domainerror) { 299fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 300fef7b6d8SPeter Brune PetscFunctionReturn(0); 301fef7b6d8SPeter Brune } 3021aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 303e4ed7901SPeter Brune if (!snes->norm_init_set) { 304fef7b6d8SPeter Brune /* convergence test */ 305*a71552e2SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 306189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 307189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 308189a9710SBarry Smith PetscFunctionReturn(0); 309189a9710SBarry Smith } 310e4ed7901SPeter Brune } else { 311e4ed7901SPeter Brune fnorm = snes->norm_init; 312e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 313e4ed7901SPeter Brune } 314*a71552e2SPeter Brune ierr = VecCopy(F,dX);CHKERRQ(ierr); 315*a71552e2SPeter Brune } 316*a71552e2SPeter Brune if (snes->pc) { 317*a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 318*a71552e2SPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,dX);CHKERRQ(ierr); 319*a71552e2SPeter Brune } 320*a71552e2SPeter Brune } 321*a71552e2SPeter Brune ierr = VecCopy(dX,lX);CHKERRQ(ierr); 322*a71552e2SPeter Brune ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr); 323*a71552e2SPeter Brune 324*a71552e2SPeter Brune 325ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 326fef7b6d8SPeter Brune snes->norm = fnorm; 327ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 328a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 329fef7b6d8SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 330fef7b6d8SPeter Brune 331fef7b6d8SPeter Brune /* set parameter for default relative tolerance convergence test */ 332fef7b6d8SPeter Brune snes->ttol = fnorm*snes->rtol; 333fef7b6d8SPeter Brune /* test convergence */ 334fef7b6d8SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 335fef7b6d8SPeter Brune if (snes->reason) PetscFunctionReturn(0); 336fef7b6d8SPeter Brune 337fef7b6d8SPeter Brune /* Call general purpose update function */ 338fef7b6d8SPeter Brune if (snes->ops->update) { 339fef7b6d8SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 340fef7b6d8SPeter Brune } 341fef7b6d8SPeter Brune 342fef7b6d8SPeter Brune /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 343bbd5d0b3SPeter Brune 34409c08436SPeter Brune for (i = 1; i < maxits + 1; i++) { 345fef7b6d8SPeter Brune lsSuccess = PETSC_TRUE; 34629c94e34SPeter Brune /* some update types require the old update direction or conjugate direction */ 3470a844d1aSPeter Brune if (ncg->type != SNES_NCG_FR) { 3485d115551SPeter Brune ierr = VecCopy(F, Fold);CHKERRQ(ierr); 34988d374b2SPeter Brune } 350f1c6b773SPeter Brune ierr = SNESLineSearchApply(linesearch,X,F,&fnorm,lX);CHKERRQ(ierr); 351f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(linesearch, &lsSuccess);CHKERRQ(ierr); 352fef7b6d8SPeter Brune if (!lsSuccess) { 353fef7b6d8SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 354fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 3555d10c001SPeter Brune PetscFunctionReturn(0); 356fef7b6d8SPeter Brune } 357fef7b6d8SPeter Brune } 358fef7b6d8SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 359fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 3605d10c001SPeter Brune PetscFunctionReturn(0); 361fef7b6d8SPeter Brune } 362fef7b6d8SPeter Brune if (snes->domainerror) { 363fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 364fef7b6d8SPeter Brune PetscFunctionReturn(0); 365fef7b6d8SPeter Brune } 366f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 367fef7b6d8SPeter Brune /* Monitor convergence */ 368ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 369169e2be8SPeter Brune snes->iter = i; 370fef7b6d8SPeter Brune snes->norm = fnorm; 371ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 372a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 373fef7b6d8SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 374302dbbaeSPeter Brune 375fef7b6d8SPeter Brune /* Test for convergence */ 376d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3775d10c001SPeter Brune if (snes->reason) PetscFunctionReturn(0); 378302dbbaeSPeter Brune 379302dbbaeSPeter Brune /* Call general purpose update function */ 380302dbbaeSPeter Brune if (snes->ops->update) { 381302dbbaeSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 382302dbbaeSPeter Brune } 383*a71552e2SPeter Brune if (snes->pc) { 384*a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 385*a71552e2SPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,dX);CHKERRQ(ierr); 386*a71552e2SPeter Brune ierr = VecCopy(dX,F);CHKERRQ(ierr); 387*a71552e2SPeter Brune } else { 388*a71552e2SPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,dX);CHKERRQ(ierr); 389302dbbaeSPeter Brune } 390a3685310SPeter Brune } else { 391a3685310SPeter Brune ierr = VecCopy(F,dX);CHKERRQ(ierr); 392302dbbaeSPeter Brune } 393302dbbaeSPeter Brune 39429c94e34SPeter Brune /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 3950a844d1aSPeter Brune switch (ncg->type) { 3960a844d1aSPeter Brune case SNES_NCG_FR: /* Fletcher-Reeves */ 3975d115551SPeter Brune dXolddotFold = dXdotF; 398*a71552e2SPeter Brune ierr = VecDot(dX, F, &dXdotF);CHKERRQ(ierr); 3995d115551SPeter Brune beta = PetscRealPart(dXdotF / dXolddotFold); 400dfb256c7SPeter Brune break; 4010a844d1aSPeter Brune case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */ 4025d115551SPeter Brune dXolddotFold = dXdotF; 40315f0db4aSPeter Brune ierr = VecDotBegin(F, dX, &dXdotF);CHKERRQ(ierr); 40415f0db4aSPeter Brune ierr = VecDotBegin(Fold, dX, &dXdotFold);CHKERRQ(ierr); 40515f0db4aSPeter Brune ierr = VecDotEnd(F, dX, &dXdotF);CHKERRQ(ierr); 40615f0db4aSPeter Brune ierr = VecDotEnd(Fold, dX, &dXdotFold);CHKERRQ(ierr); 4075d115551SPeter Brune beta = PetscRealPart(((dXdotF - dXdotFold) / dXolddotFold)); 408dfb256c7SPeter Brune if (beta < 0.0) beta = 0.0; /* restart */ 409dfb256c7SPeter Brune break; 4100a844d1aSPeter Brune case SNES_NCG_HS: /* Hestenes-Stiefel */ 41115f0db4aSPeter Brune ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr); 41215f0db4aSPeter Brune ierr = VecDotBegin(dX, Fold, &dXdotFold);CHKERRQ(ierr); 41315f0db4aSPeter Brune ierr = VecDotBegin(lX, F, &lXdotF);CHKERRQ(ierr); 41415f0db4aSPeter Brune ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr); 41515f0db4aSPeter Brune ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr); 41615f0db4aSPeter Brune ierr = VecDotEnd(dX, Fold, &dXdotFold);CHKERRQ(ierr); 41715f0db4aSPeter Brune ierr = VecDotEnd(lX, F, &lXdotF);CHKERRQ(ierr); 41815f0db4aSPeter Brune ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr); 4195d115551SPeter Brune beta = PetscRealPart((dXdotF - dXdotFold) / (lXdotF - lXdotFold)); 420dfb256c7SPeter Brune break; 4210a844d1aSPeter Brune case SNES_NCG_DY: /* Dai-Yuan */ 42215f0db4aSPeter Brune ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr); 42315f0db4aSPeter Brune ierr = VecDotBegin(lX, F, &lXdotF);CHKERRQ(ierr); 42415f0db4aSPeter Brune ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr); 42515f0db4aSPeter Brune ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr); 42615f0db4aSPeter Brune ierr = VecDotEnd(lX, F, &lXdotF);CHKERRQ(ierr); 42715f0db4aSPeter Brune ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr); 4285d115551SPeter Brune beta = PetscRealPart(dXdotF / (lXdotFold - lXdotF));CHKERRQ(ierr); 429dfb256c7SPeter Brune break; 4300a844d1aSPeter Brune case SNES_NCG_CD: /* Conjugate Descent */ 43115f0db4aSPeter Brune ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr); 43215f0db4aSPeter Brune ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr); 43315f0db4aSPeter Brune ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr); 43415f0db4aSPeter Brune ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr); 4355d115551SPeter Brune beta = PetscRealPart(dXdotF / lXdotFold);CHKERRQ(ierr); 436dfb256c7SPeter Brune break; 437dfb256c7SPeter Brune } 438dfb256c7SPeter Brune if (ncg->monitor) { 439646217ecSPeter Brune ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);CHKERRQ(ierr); 440dfb256c7SPeter Brune } 441dfb256c7SPeter Brune ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr); 442fef7b6d8SPeter Brune } 443fef7b6d8SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 444fef7b6d8SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 445fef7b6d8SPeter Brune PetscFunctionReturn(0); 446fef7b6d8SPeter Brune } 447fef7b6d8SPeter Brune 4480a844d1aSPeter Brune 4490a844d1aSPeter Brune 450fef7b6d8SPeter Brune /*MC 451fef7b6d8SPeter Brune SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems. 452fef7b6d8SPeter Brune 453fef7b6d8SPeter Brune Level: beginner 454fef7b6d8SPeter Brune 455fef7b6d8SPeter Brune Options Database: 4561867fe5bSPeter Brune + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter. 45741a8cb50SPeter Brune . -snes_linesearch_type <cp,l2,basic> - Line search type. 4581867fe5bSPeter Brune - -snes_ncg_monitor - Print relevant information about the ncg iteration. 459fef7b6d8SPeter Brune 460fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate 461fef7b6d8SPeter Brune gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 462fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x. 463fef7b6d8SPeter Brune 464fef7b6d8SPeter Brune 46504d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN 466fef7b6d8SPeter Brune M*/ 467fef7b6d8SPeter Brune #undef __FUNCT__ 468fef7b6d8SPeter Brune #define __FUNCT__ "SNESCreate_NCG" 4698cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes) 470fef7b6d8SPeter Brune { 471fef7b6d8SPeter Brune PetscErrorCode ierr; 472ea630c6eSPeter Brune SNES_NCG * neP; 473fef7b6d8SPeter Brune 474fef7b6d8SPeter Brune PetscFunctionBegin; 475fef7b6d8SPeter Brune snes->ops->destroy = SNESDestroy_NCG; 476fef7b6d8SPeter Brune snes->ops->setup = SNESSetUp_NCG; 477fef7b6d8SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NCG; 478fef7b6d8SPeter Brune snes->ops->view = SNESView_NCG; 479fef7b6d8SPeter Brune snes->ops->solve = SNESSolve_NCG; 480fef7b6d8SPeter Brune snes->ops->reset = SNESReset_NCG; 481fef7b6d8SPeter Brune 482fef7b6d8SPeter Brune snes->usesksp = PETSC_FALSE; 483fef7b6d8SPeter Brune snes->usespc = PETSC_TRUE; 484*a71552e2SPeter Brune snes->pcside = PC_LEFT; 485*a71552e2SPeter Brune snes->functype = SNES_FUNCTION_PRECONDITIONED; 486fef7b6d8SPeter Brune 48788976e71SPeter Brune if (!snes->tolerancesset) { 4880e444f03SPeter Brune snes->max_funcs = 30000; 4890e444f03SPeter Brune snes->max_its = 10000; 490c522fa08SPeter Brune snes->stol = 1e-20; 49188976e71SPeter Brune } 4920e444f03SPeter Brune 493fef7b6d8SPeter Brune ierr = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr); 494fef7b6d8SPeter Brune snes->data = (void*) neP; 4950298fd71SBarry Smith neP->monitor = NULL; 4960a844d1aSPeter Brune neP->type = SNES_NCG_PRP; 497bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);CHKERRQ(ierr); 498fef7b6d8SPeter Brune PetscFunctionReturn(0); 499fef7b6d8SPeter Brune } 500