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 129105210eSPeter Brune #define SNESLINESEARCHNCGLINEAR "ncglinear" 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); 54a71552e2SPeter Brune if (snes->pcside == PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning"); 556c67d002SPeter Brune if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED; 56fef7b6d8SPeter Brune PetscFunctionReturn(0); 57fef7b6d8SPeter Brune } 58fef7b6d8SPeter Brune /* 5905b53524SPeter Brune SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method. 60fef7b6d8SPeter Brune 61fef7b6d8SPeter Brune Input Parameter: 62fef7b6d8SPeter Brune . snes - the SNES context 63fef7b6d8SPeter Brune 64fef7b6d8SPeter Brune Application Interface Routine: SNESSetFromOptions() 65fef7b6d8SPeter Brune */ 66fef7b6d8SPeter Brune #undef __FUNCT__ 67fef7b6d8SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NCG" 68fef7b6d8SPeter Brune static PetscErrorCode SNESSetFromOptions_NCG(SNES snes) 69fef7b6d8SPeter Brune { 70dfb256c7SPeter Brune SNES_NCG *ncg = (SNES_NCG*)snes->data; 71fef7b6d8SPeter Brune PetscErrorCode ierr; 720a844d1aSPeter Brune PetscBool debug; 73f1c6b773SPeter Brune SNESLineSearch linesearch; 74a11d90f7SPeter Brune SNESNCGType ncgtype=ncg->type; 75fef7b6d8SPeter Brune 76fef7b6d8SPeter Brune PetscFunctionBegin; 77fef7b6d8SPeter Brune ierr = PetscOptionsHead("SNES NCG options");CHKERRQ(ierr); 780298fd71SBarry Smith ierr = PetscOptionsBool("-snes_ncg_monitor","Monitor NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);CHKERRQ(ierr); 790298fd71SBarry Smith ierr = PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);CHKERRQ(ierr); 800a844d1aSPeter Brune ierr = SNESNCGSetType(snes, ncgtype);CHKERRQ(ierr); 81dfb256c7SPeter Brune if (debug) { 82ce94432eSBarry Smith ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr); 83dfb256c7SPeter Brune } 84fef7b6d8SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 859e764e56SPeter Brune if (!snes->linesearch) { 867601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 879e764e56SPeter Brune if (!snes->pc) { 881a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);CHKERRQ(ierr); 899e764e56SPeter Brune } else { 901a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 919e764e56SPeter Brune } 929e764e56SPeter Brune } 939105210eSPeter Brune ierr = SNESLineSearchRegister(SNESLINESEARCHNCGLINEAR, SNESLineSearchCreate_NCGLinear);CHKERRQ(ierr); 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 1419105210eSPeter Brune if (!snes->jacobian) { 1429105210eSPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 1439105210eSPeter Brune } 1449105210eSPeter Brune 145ea630c6eSPeter Brune /* 146ea630c6eSPeter Brune 147d2140ca3SPeter Brune The exact step size for unpreconditioned linear CG is just: 148ea630c6eSPeter Brune alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) 149ea630c6eSPeter Brune */ 150a5892487SPeter Brune ierr = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr); 151ea630c6eSPeter Brune ierr = VecDot(F, F, &alpha);CHKERRQ(ierr); 152ea630c6eSPeter Brune ierr = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr); 153ea630c6eSPeter Brune ierr = VecDot(Y, W, &ptAp);CHKERRQ(ierr); 154ea630c6eSPeter Brune alpha = alpha / ptAp; 1559105210eSPeter Brune ierr = VecAXPY(X,-alpha,Y);CHKERRQ(ierr); 156bbd5d0b3SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 157bbd5d0b3SPeter Brune 158bbd5d0b3SPeter Brune ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr); 159bbd5d0b3SPeter Brune ierr = VecNorm(X,NORM_2,xnorm);CHKERRQ(ierr); 160bbd5d0b3SPeter Brune ierr = VecNorm(Y,NORM_2,ynorm);CHKERRQ(ierr); 161ea630c6eSPeter Brune PetscFunctionReturn(0); 162ea630c6eSPeter Brune } 163ea630c6eSPeter Brune 164bbd5d0b3SPeter Brune #undef __FUNCT__ 165f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_NCGLinear" 1668cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch) 167bbd5d0b3SPeter Brune { 168bbd5d0b3SPeter Brune PetscFunctionBegin; 169f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_NCGLinear; 1700298fd71SBarry Smith linesearch->ops->destroy = NULL; 1710298fd71SBarry Smith linesearch->ops->setfromoptions = NULL; 1720298fd71SBarry Smith linesearch->ops->reset = NULL; 1730298fd71SBarry Smith linesearch->ops->view = NULL; 1740298fd71SBarry Smith linesearch->ops->setup = NULL; 175bbd5d0b3SPeter Brune PetscFunctionReturn(0); 176bbd5d0b3SPeter Brune } 177bbd5d0b3SPeter Brune 17888d374b2SPeter Brune #undef __FUNCT__ 1798c40d5fbSBarry Smith #define __FUNCT__ "SNESNCGComputeYtJtF_Private" 18088d374b2SPeter Brune /* 18188d374b2SPeter Brune 18288d374b2SPeter Brune Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian. 18388d374b2SPeter Brune 18488d374b2SPeter Brune */ 18567392de3SBarry Smith PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf) 18667392de3SBarry Smith { 18788d374b2SPeter Brune PetscErrorCode ierr; 18888d374b2SPeter Brune PetscScalar ftf, ftg, fty, h; 18967392de3SBarry Smith 19088d374b2SPeter Brune PetscFunctionBegin; 19188d374b2SPeter Brune ierr = VecDot(F, F, &ftf);CHKERRQ(ierr); 19288d374b2SPeter Brune ierr = VecDot(F, Y, &fty);CHKERRQ(ierr); 19388d374b2SPeter Brune h = 1e-5*fty / fty; 19488d374b2SPeter Brune ierr = VecCopy(X, W);CHKERRQ(ierr); 19588d374b2SPeter Brune ierr = VecAXPY(W, -h, Y);CHKERRQ(ierr); /* this is arbitrary */ 19688d374b2SPeter Brune ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr); 19788d374b2SPeter Brune ierr = VecDot(G, F, &ftg);CHKERRQ(ierr); 19888d374b2SPeter Brune *ytJtf = (ftg - ftf) / h; 19988d374b2SPeter Brune PetscFunctionReturn(0); 20088d374b2SPeter Brune } 20188d374b2SPeter Brune 2020a844d1aSPeter Brune #undef __FUNCT__ 2030a844d1aSPeter Brune #define __FUNCT__ "SNESNCGSetType" 2040a844d1aSPeter Brune /*@ 2050a844d1aSPeter Brune SNESNCGSetType - Sets the conjugate update type for SNESNCG. 2060a844d1aSPeter Brune 2070a844d1aSPeter Brune Logically Collective on SNES 2080a844d1aSPeter Brune 2090a844d1aSPeter Brune Input Parameters: 2100a844d1aSPeter Brune + snes - the iterative context 2110a844d1aSPeter Brune - btype - update type 2120a844d1aSPeter Brune 2130a844d1aSPeter Brune Options Database: 2140a844d1aSPeter Brune . -snes_ncg_type<prp,fr,hs,dy,cd> 2150a844d1aSPeter Brune 2160a844d1aSPeter Brune Level: intermediate 2170a844d1aSPeter Brune 2180a844d1aSPeter Brune SNESNCGSelectTypes: 2190a844d1aSPeter Brune + SNES_NCG_FR - Fletcher-Reeves update 2200a844d1aSPeter Brune . SNES_NCG_PRP - Polak-Ribiere-Polyak update 2210a844d1aSPeter Brune . SNES_NCG_HS - Hestenes-Steifel update 2220a844d1aSPeter Brune . SNES_NCG_DY - Dai-Yuan update 2230a844d1aSPeter Brune - SNES_NCG_CD - Conjugate Descent update 2240a844d1aSPeter Brune 2250a844d1aSPeter Brune Notes: 2260a844d1aSPeter Brune PRP is the default, and the only one that tolerates generalized search directions. 2270a844d1aSPeter Brune 2280a844d1aSPeter Brune .keywords: SNES, SNESNCG, selection, type, set 2290a844d1aSPeter Brune @*/ 2300adebc6cSBarry Smith PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype) 2310adebc6cSBarry Smith { 2320a844d1aSPeter Brune PetscErrorCode ierr; 2336e111a19SKarl Rupp 2340a844d1aSPeter Brune PetscFunctionBegin; 2350a844d1aSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2360a844d1aSPeter Brune ierr = PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));CHKERRQ(ierr); 2370a844d1aSPeter Brune PetscFunctionReturn(0); 2380a844d1aSPeter Brune } 2390a844d1aSPeter Brune 2400a844d1aSPeter Brune #undef __FUNCT__ 2410a844d1aSPeter Brune #define __FUNCT__ "SNESNCGSetType_NCG" 2420adebc6cSBarry Smith PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype) 2430adebc6cSBarry Smith { 2440a844d1aSPeter Brune SNES_NCG *ncg = (SNES_NCG*)snes->data; 2456e111a19SKarl Rupp 2460a844d1aSPeter Brune PetscFunctionBegin; 2470a844d1aSPeter Brune ncg->type = btype; 2480a844d1aSPeter Brune PetscFunctionReturn(0); 2490a844d1aSPeter Brune } 2500a844d1aSPeter Brune 251fef7b6d8SPeter Brune /* 252fef7b6d8SPeter Brune SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method. 253fef7b6d8SPeter Brune 254fef7b6d8SPeter Brune Input Parameters: 255fef7b6d8SPeter Brune . snes - the SNES context 256fef7b6d8SPeter Brune 257fef7b6d8SPeter Brune Output Parameter: 258fef7b6d8SPeter Brune . outits - number of iterations until termination 259fef7b6d8SPeter Brune 260fef7b6d8SPeter Brune Application Interface Routine: SNESSolve() 261fef7b6d8SPeter Brune */ 262fef7b6d8SPeter Brune #undef __FUNCT__ 263fef7b6d8SPeter Brune #define __FUNCT__ "SNESSolve_NCG" 264fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes) 265fef7b6d8SPeter Brune { 266dfb256c7SPeter Brune SNES_NCG *ncg = (SNES_NCG*)snes->data; 26732cc618eSPeter Brune Vec X,dX,lX,F,dXold; 268bbd5d0b3SPeter Brune PetscReal fnorm, ynorm, xnorm, beta = 0.0; 26932cc618eSPeter Brune PetscScalar dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold; 270fef7b6d8SPeter Brune PetscInt maxits, i; 271fef7b6d8SPeter Brune PetscErrorCode ierr; 272fef7b6d8SPeter Brune PetscBool lsSuccess = PETSC_TRUE; 273f1c6b773SPeter Brune SNESLineSearch linesearch; 274b7281c8aSPeter Brune SNESConvergedReason reason; 27588d374b2SPeter Brune 276fef7b6d8SPeter Brune PetscFunctionBegin; 277fef7b6d8SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 278fef7b6d8SPeter Brune 279fef7b6d8SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 280fef7b6d8SPeter Brune X = snes->vec_sol; /* X^n */ 28132cc618eSPeter Brune dXold = snes->work[0]; /* The previous iterate of X */ 28288d374b2SPeter Brune dX = snes->work[1]; /* the preconditioned direction */ 283169e2be8SPeter Brune lX = snes->vec_sol_update; /* the conjugate direction */ 284fef7b6d8SPeter Brune F = snes->vec_func; /* residual vector */ 285fef7b6d8SPeter Brune 2867601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 2879e764e56SPeter Brune 288e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 289fef7b6d8SPeter Brune snes->iter = 0; 290fef7b6d8SPeter Brune snes->norm = 0.; 291e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 292fef7b6d8SPeter Brune 293bbd5d0b3SPeter Brune /* compute the initial function and preconditioned update dX */ 294a71552e2SPeter Brune 295a71552e2SPeter Brune if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) { 296a71552e2SPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,dX);CHKERRQ(ierr); 297b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 298b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 299b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 300b7281c8aSPeter Brune PetscFunctionReturn(0); 301b7281c8aSPeter Brune } 302a71552e2SPeter Brune ierr = VecCopy(dX,F);CHKERRQ(ierr); 303a71552e2SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 304a71552e2SPeter Brune } else { 305e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 306fef7b6d8SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 307fef7b6d8SPeter Brune if (snes->domainerror) { 308fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 309fef7b6d8SPeter Brune PetscFunctionReturn(0); 310fef7b6d8SPeter Brune } 3111aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 312c1c75074SPeter Brune 313fef7b6d8SPeter Brune /* convergence test */ 314a71552e2SPeter Brune ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); 315189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 316189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 317189a9710SBarry Smith PetscFunctionReturn(0); 318189a9710SBarry Smith } 319c1c75074SPeter Brune 320a71552e2SPeter Brune ierr = VecCopy(F,dX);CHKERRQ(ierr); 321a71552e2SPeter Brune } 322a71552e2SPeter Brune if (snes->pc) { 323a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 324a71552e2SPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,dX);CHKERRQ(ierr); 325b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 326b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 327b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 328b7281c8aSPeter Brune PetscFunctionReturn(0); 329b7281c8aSPeter Brune } 330a71552e2SPeter Brune } 331a71552e2SPeter Brune } 332a71552e2SPeter Brune ierr = VecCopy(dX,lX);CHKERRQ(ierr); 33332cc618eSPeter Brune ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 334a71552e2SPeter Brune 335a71552e2SPeter Brune 336e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 337fef7b6d8SPeter Brune snes->norm = fnorm; 338e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 339a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 340fef7b6d8SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 341fef7b6d8SPeter Brune 342fef7b6d8SPeter Brune /* test convergence */ 343fef7b6d8SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 344fef7b6d8SPeter Brune if (snes->reason) PetscFunctionReturn(0); 345fef7b6d8SPeter Brune 346fef7b6d8SPeter Brune /* Call general purpose update function */ 347fef7b6d8SPeter Brune if (snes->ops->update) { 348fef7b6d8SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 349fef7b6d8SPeter Brune } 350fef7b6d8SPeter Brune 351fef7b6d8SPeter Brune /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */ 352bbd5d0b3SPeter Brune 35309c08436SPeter Brune for (i = 1; i < maxits + 1; i++) { 354fef7b6d8SPeter Brune lsSuccess = PETSC_TRUE; 35529c94e34SPeter Brune /* some update types require the old update direction or conjugate direction */ 3560a844d1aSPeter Brune if (ncg->type != SNES_NCG_FR) { 35732cc618eSPeter Brune ierr = VecCopy(dX, dXold);CHKERRQ(ierr); 35888d374b2SPeter Brune } 359f1c6b773SPeter Brune ierr = SNESLineSearchApply(linesearch,X,F,&fnorm,lX);CHKERRQ(ierr); 360f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(linesearch, &lsSuccess);CHKERRQ(ierr); 361fef7b6d8SPeter Brune if (!lsSuccess) { 362fef7b6d8SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 363fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 3645d10c001SPeter Brune PetscFunctionReturn(0); 365fef7b6d8SPeter Brune } 366fef7b6d8SPeter Brune } 367fef7b6d8SPeter Brune if (snes->nfuncs >= snes->max_funcs) { 368fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 3695d10c001SPeter Brune PetscFunctionReturn(0); 370fef7b6d8SPeter Brune } 371fef7b6d8SPeter Brune if (snes->domainerror) { 372fef7b6d8SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 373fef7b6d8SPeter Brune PetscFunctionReturn(0); 374fef7b6d8SPeter Brune } 375f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 376fef7b6d8SPeter Brune /* Monitor convergence */ 377e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 378169e2be8SPeter Brune snes->iter = i; 379fef7b6d8SPeter Brune snes->norm = fnorm; 380e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 381a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 382fef7b6d8SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 383302dbbaeSPeter Brune 384fef7b6d8SPeter Brune /* Test for convergence */ 385d484d688SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3865d10c001SPeter Brune if (snes->reason) PetscFunctionReturn(0); 387302dbbaeSPeter Brune 388302dbbaeSPeter Brune /* Call general purpose update function */ 389302dbbaeSPeter Brune if (snes->ops->update) { 390302dbbaeSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 391302dbbaeSPeter Brune } 392a71552e2SPeter Brune if (snes->pc) { 393a71552e2SPeter Brune if (snes->functype == SNES_FUNCTION_PRECONDITIONED) { 394a71552e2SPeter Brune ierr = SNESApplyPC(snes,X,NULL,NULL,dX);CHKERRQ(ierr); 395b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 396b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 397b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 398b7281c8aSPeter Brune PetscFunctionReturn(0); 399b7281c8aSPeter Brune } 400a71552e2SPeter Brune ierr = VecCopy(dX,F);CHKERRQ(ierr); 401a71552e2SPeter Brune } else { 402a71552e2SPeter Brune ierr = SNESApplyPC(snes,X,F,&fnorm,dX);CHKERRQ(ierr); 403b7281c8aSPeter Brune ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr); 404b7281c8aSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 405b7281c8aSPeter Brune snes->reason = SNES_DIVERGED_INNER; 406b7281c8aSPeter Brune PetscFunctionReturn(0); 407b7281c8aSPeter Brune } 408302dbbaeSPeter Brune } 409a3685310SPeter Brune } else { 410a3685310SPeter Brune ierr = VecCopy(F,dX);CHKERRQ(ierr); 411302dbbaeSPeter Brune } 412302dbbaeSPeter Brune 41329c94e34SPeter Brune /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/ 4140a844d1aSPeter Brune switch (ncg->type) { 4150a844d1aSPeter Brune case SNES_NCG_FR: /* Fletcher-Reeves */ 41632cc618eSPeter Brune dXolddotdXold= dXdotdX; 41732cc618eSPeter Brune ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr); 41832cc618eSPeter Brune beta = PetscRealPart(dXdotdX / dXolddotdXold); 419dfb256c7SPeter Brune break; 4200a844d1aSPeter Brune case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */ 42132cc618eSPeter Brune dXolddotdXold= dXdotdX; 42232cc618eSPeter Brune ierr = VecDotBegin(dX, dX, &dXdotdXold);CHKERRQ(ierr); 42332cc618eSPeter Brune ierr = VecDotBegin(dXold, dX, &dXdotdXold);CHKERRQ(ierr); 42432cc618eSPeter Brune ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 42532cc618eSPeter Brune ierr = VecDotEnd(dXold, dX, &dXdotdXold);CHKERRQ(ierr); 42632cc618eSPeter Brune beta = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold)); 427dfb256c7SPeter Brune if (beta < 0.0) beta = 0.0; /* restart */ 428dfb256c7SPeter Brune break; 4290a844d1aSPeter Brune case SNES_NCG_HS: /* Hestenes-Stiefel */ 43032cc618eSPeter Brune ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 43132cc618eSPeter Brune ierr = VecDotBegin(dX, dXold, &dXdotdXold);CHKERRQ(ierr); 43232cc618eSPeter Brune ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr); 43332cc618eSPeter Brune ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 43432cc618eSPeter Brune ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 43532cc618eSPeter Brune ierr = VecDotEnd(dX, dXold, &dXdotdXold);CHKERRQ(ierr); 43632cc618eSPeter Brune ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr); 43732cc618eSPeter Brune ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 43832cc618eSPeter Brune beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold)); 439dfb256c7SPeter Brune break; 4400a844d1aSPeter Brune case SNES_NCG_DY: /* Dai-Yuan */ 44132cc618eSPeter Brune ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 44232cc618eSPeter Brune ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr); 44332cc618eSPeter Brune ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 44432cc618eSPeter Brune ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 44532cc618eSPeter Brune ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr); 44632cc618eSPeter Brune ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 44732cc618eSPeter Brune beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));CHKERRQ(ierr); 448dfb256c7SPeter Brune break; 4490a844d1aSPeter Brune case SNES_NCG_CD: /* Conjugate Descent */ 45032cc618eSPeter Brune ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr); 45132cc618eSPeter Brune ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 45232cc618eSPeter Brune ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr); 45332cc618eSPeter Brune ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr); 45432cc618eSPeter Brune beta = PetscRealPart(dXdotdX / lXdotdXold);CHKERRQ(ierr); 455dfb256c7SPeter Brune break; 456dfb256c7SPeter Brune } 457dfb256c7SPeter Brune if (ncg->monitor) { 458646217ecSPeter Brune ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);CHKERRQ(ierr); 459dfb256c7SPeter Brune } 460dfb256c7SPeter Brune ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr); 461fef7b6d8SPeter Brune } 462fef7b6d8SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 463fef7b6d8SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 464fef7b6d8SPeter Brune PetscFunctionReturn(0); 465fef7b6d8SPeter Brune } 466fef7b6d8SPeter Brune 4670a844d1aSPeter Brune 4680a844d1aSPeter Brune 469fef7b6d8SPeter Brune /*MC 470fef7b6d8SPeter Brune SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems. 471fef7b6d8SPeter Brune 472fef7b6d8SPeter Brune Level: beginner 473fef7b6d8SPeter Brune 474fef7b6d8SPeter Brune Options Database: 4751867fe5bSPeter Brune + -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter. 47641a8cb50SPeter Brune . -snes_linesearch_type <cp,l2,basic> - Line search type. 4771867fe5bSPeter Brune - -snes_ncg_monitor - Print relevant information about the ncg iteration. 478fef7b6d8SPeter Brune 479fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate 480fef7b6d8SPeter Brune gradient method. This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise 481fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x. 482fef7b6d8SPeter Brune 483fef7b6d8SPeter Brune 48404d7464bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN 485fef7b6d8SPeter Brune M*/ 486fef7b6d8SPeter Brune #undef __FUNCT__ 487fef7b6d8SPeter Brune #define __FUNCT__ "SNESCreate_NCG" 4888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes) 489fef7b6d8SPeter Brune { 490fef7b6d8SPeter Brune PetscErrorCode ierr; 491ea630c6eSPeter Brune SNES_NCG * neP; 492fef7b6d8SPeter Brune 493fef7b6d8SPeter Brune PetscFunctionBegin; 494fef7b6d8SPeter Brune snes->ops->destroy = SNESDestroy_NCG; 495fef7b6d8SPeter Brune snes->ops->setup = SNESSetUp_NCG; 496fef7b6d8SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NCG; 497fef7b6d8SPeter Brune snes->ops->view = SNESView_NCG; 498fef7b6d8SPeter Brune snes->ops->solve = SNESSolve_NCG; 499fef7b6d8SPeter Brune snes->ops->reset = SNESReset_NCG; 500fef7b6d8SPeter Brune 501fef7b6d8SPeter Brune snes->usesksp = PETSC_FALSE; 502fef7b6d8SPeter Brune snes->usespc = PETSC_TRUE; 503a71552e2SPeter Brune snes->pcside = PC_LEFT; 504fef7b6d8SPeter Brune 50588976e71SPeter Brune if (!snes->tolerancesset) { 5060e444f03SPeter Brune snes->max_funcs = 30000; 5070e444f03SPeter Brune snes->max_its = 10000; 508c522fa08SPeter Brune snes->stol = 1e-20; 50988976e71SPeter Brune } 5100e444f03SPeter Brune 511*b00a9115SJed Brown ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 512fef7b6d8SPeter Brune snes->data = (void*) neP; 5130298fd71SBarry Smith neP->monitor = NULL; 5140a844d1aSPeter Brune neP->type = SNES_NCG_PRP; 515bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);CHKERRQ(ierr); 516fef7b6d8SPeter Brune PetscFunctionReturn(0); 517fef7b6d8SPeter Brune } 518