xref: /petsc/src/snes/impls/ncg/snesncg.c (revision 422a814ef4a731b8529cf3a6428db526d183e312)
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"
688c34d3f5SBarry Smith static PetscErrorCode SNESSetFromOptions_NCG(PetscOptions *PetscOptionsObject,SNES snes)
69fef7b6d8SPeter Brune {
70dfb256c7SPeter Brune   SNES_NCG       *ncg = (SNES_NCG*)snes->data;
71fef7b6d8SPeter Brune   PetscErrorCode ierr;
7294ae4db5SBarry Smith   PetscBool      debug = PETSC_FALSE;
73f1c6b773SPeter Brune   SNESLineSearch linesearch;
74a11d90f7SPeter Brune   SNESNCGType    ncgtype=ncg->type;
75fef7b6d8SPeter Brune 
76fef7b6d8SPeter Brune   PetscFunctionBegin;
77e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"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);
79dfb256c7SPeter Brune   if (debug) {
80ce94432eSBarry Smith     ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));CHKERRQ(ierr);
81dfb256c7SPeter Brune   }
8294ae4db5SBarry Smith   ierr = PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);CHKERRQ(ierr);
8394ae4db5SBarry Smith   ierr = SNESNCGSetType(snes, ncgtype);CHKERRQ(ierr);
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 
130ea630c6eSPeter Brune   PetscFunctionBegin;
131f1c6b773SPeter Brune   ierr  = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
132bbd5d0b3SPeter Brune   X     = linesearch->vec_sol;
133bbd5d0b3SPeter Brune   W     = linesearch->vec_sol_new;
134bbd5d0b3SPeter Brune   F     = linesearch->vec_func;
135bbd5d0b3SPeter Brune   Y     = linesearch->vec_update;
136bbd5d0b3SPeter Brune   fnorm = &linesearch->fnorm;
137bbd5d0b3SPeter Brune   xnorm = &linesearch->xnorm;
138bbd5d0b3SPeter Brune   ynorm = &linesearch->ynorm;
139bbd5d0b3SPeter Brune 
1409105210eSPeter Brune   if (!snes->jacobian) {
1419105210eSPeter Brune     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
1429105210eSPeter Brune   }
1439105210eSPeter Brune 
144ea630c6eSPeter Brune   /*
145ea630c6eSPeter Brune 
146d2140ca3SPeter Brune    The exact step size for unpreconditioned linear CG is just:
147ea630c6eSPeter Brune    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
148ea630c6eSPeter Brune    */
149d1e9a80fSBarry Smith   ierr  = SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre);CHKERRQ(ierr);
150ea630c6eSPeter Brune   ierr  = VecDot(F, F, &alpha);CHKERRQ(ierr);
151ea630c6eSPeter Brune   ierr  = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr);
152ea630c6eSPeter Brune   ierr  = VecDot(Y, W, &ptAp);CHKERRQ(ierr);
153ea630c6eSPeter Brune   alpha = alpha / ptAp;
1549105210eSPeter Brune   ierr  = VecAXPY(X,-alpha,Y);CHKERRQ(ierr);
155bbd5d0b3SPeter Brune   ierr  = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
156bbd5d0b3SPeter Brune 
157bbd5d0b3SPeter Brune   ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);
158bbd5d0b3SPeter Brune   ierr = VecNorm(X,NORM_2,xnorm);CHKERRQ(ierr);
159bbd5d0b3SPeter Brune   ierr = VecNorm(Y,NORM_2,ynorm);CHKERRQ(ierr);
160ea630c6eSPeter Brune   PetscFunctionReturn(0);
161ea630c6eSPeter Brune }
162ea630c6eSPeter Brune 
163bbd5d0b3SPeter Brune #undef __FUNCT__
164f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_NCGLinear"
1658cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
166bbd5d0b3SPeter Brune {
167bbd5d0b3SPeter Brune   PetscFunctionBegin;
168f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
1690298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
1700298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
1710298fd71SBarry Smith   linesearch->ops->reset          = NULL;
1720298fd71SBarry Smith   linesearch->ops->view           = NULL;
1730298fd71SBarry Smith   linesearch->ops->setup          = NULL;
174bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
175bbd5d0b3SPeter Brune }
176bbd5d0b3SPeter Brune 
17788d374b2SPeter Brune #undef __FUNCT__
1788c40d5fbSBarry Smith #define __FUNCT__ "SNESNCGComputeYtJtF_Private"
17988d374b2SPeter Brune /*
18088d374b2SPeter Brune 
18188d374b2SPeter Brune  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
18288d374b2SPeter Brune 
18388d374b2SPeter Brune  */
18467392de3SBarry Smith PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
18567392de3SBarry Smith {
18688d374b2SPeter Brune   PetscErrorCode ierr;
18788d374b2SPeter Brune   PetscScalar    ftf, ftg, fty, h;
18867392de3SBarry Smith 
18988d374b2SPeter Brune   PetscFunctionBegin;
19088d374b2SPeter Brune   ierr   = VecDot(F, F, &ftf);CHKERRQ(ierr);
19188d374b2SPeter Brune   ierr   = VecDot(F, Y, &fty);CHKERRQ(ierr);
19288d374b2SPeter Brune   h      = 1e-5*fty / fty;
19388d374b2SPeter Brune   ierr   = VecCopy(X, W);CHKERRQ(ierr);
19488d374b2SPeter Brune   ierr   = VecAXPY(W, -h, Y);CHKERRQ(ierr);          /* this is arbitrary */
19588d374b2SPeter Brune   ierr   = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
19688d374b2SPeter Brune   ierr   = VecDot(G, F, &ftg);CHKERRQ(ierr);
19788d374b2SPeter Brune   *ytJtf = (ftg - ftf) / h;
19888d374b2SPeter Brune   PetscFunctionReturn(0);
19988d374b2SPeter Brune }
20088d374b2SPeter Brune 
2010a844d1aSPeter Brune #undef __FUNCT__
2020a844d1aSPeter Brune #define __FUNCT__ "SNESNCGSetType"
2030a844d1aSPeter Brune /*@
2040a844d1aSPeter Brune     SNESNCGSetType - Sets the conjugate update type for SNESNCG.
2050a844d1aSPeter Brune 
2060a844d1aSPeter Brune     Logically Collective on SNES
2070a844d1aSPeter Brune 
2080a844d1aSPeter Brune     Input Parameters:
2090a844d1aSPeter Brune +   snes - the iterative context
2100a844d1aSPeter Brune -   btype - update type
2110a844d1aSPeter Brune 
2120a844d1aSPeter Brune     Options Database:
2130a844d1aSPeter Brune .   -snes_ncg_type<prp,fr,hs,dy,cd>
2140a844d1aSPeter Brune 
2150a844d1aSPeter Brune     Level: intermediate
2160a844d1aSPeter Brune 
2170a844d1aSPeter Brune     SNESNCGSelectTypes:
2180a844d1aSPeter Brune +   SNES_NCG_FR - Fletcher-Reeves update
2190a844d1aSPeter Brune .   SNES_NCG_PRP - Polak-Ribiere-Polyak update
2200a844d1aSPeter Brune .   SNES_NCG_HS - Hestenes-Steifel update
2210a844d1aSPeter Brune .   SNES_NCG_DY - Dai-Yuan update
2220a844d1aSPeter Brune -   SNES_NCG_CD - Conjugate Descent update
2230a844d1aSPeter Brune 
2240a844d1aSPeter Brune    Notes:
2250a844d1aSPeter Brune    PRP is the default, and the only one that tolerates generalized search directions.
2260a844d1aSPeter Brune 
2270a844d1aSPeter Brune .keywords: SNES, SNESNCG, selection, type, set
2280a844d1aSPeter Brune @*/
2290adebc6cSBarry Smith PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
2300adebc6cSBarry Smith {
2310a844d1aSPeter Brune   PetscErrorCode ierr;
2326e111a19SKarl Rupp 
2330a844d1aSPeter Brune   PetscFunctionBegin;
2340a844d1aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2350a844d1aSPeter Brune   ierr = PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));CHKERRQ(ierr);
2360a844d1aSPeter Brune   PetscFunctionReturn(0);
2370a844d1aSPeter Brune }
2380a844d1aSPeter Brune 
2390a844d1aSPeter Brune #undef __FUNCT__
2400a844d1aSPeter Brune #define __FUNCT__ "SNESNCGSetType_NCG"
2410adebc6cSBarry Smith PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
2420adebc6cSBarry Smith {
2430a844d1aSPeter Brune   SNES_NCG *ncg = (SNES_NCG*)snes->data;
2446e111a19SKarl Rupp 
2450a844d1aSPeter Brune   PetscFunctionBegin;
2460a844d1aSPeter Brune   ncg->type = btype;
2470a844d1aSPeter Brune   PetscFunctionReturn(0);
2480a844d1aSPeter Brune }
2490a844d1aSPeter Brune 
250fef7b6d8SPeter Brune /*
251fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
252fef7b6d8SPeter Brune 
253fef7b6d8SPeter Brune   Input Parameters:
254fef7b6d8SPeter Brune . snes - the SNES context
255fef7b6d8SPeter Brune 
256fef7b6d8SPeter Brune   Output Parameter:
257fef7b6d8SPeter Brune . outits - number of iterations until termination
258fef7b6d8SPeter Brune 
259fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
260fef7b6d8SPeter Brune */
261fef7b6d8SPeter Brune #undef __FUNCT__
262fef7b6d8SPeter Brune #define __FUNCT__ "SNESSolve_NCG"
263fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes)
264fef7b6d8SPeter Brune {
265dfb256c7SPeter Brune   SNES_NCG             *ncg = (SNES_NCG*)snes->data;
26632cc618eSPeter Brune   Vec                  X,dX,lX,F,dXold;
267bbd5d0b3SPeter Brune   PetscReal            fnorm, ynorm, xnorm, beta = 0.0;
26832cc618eSPeter Brune   PetscScalar          dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
269fef7b6d8SPeter Brune   PetscInt             maxits, i;
270fef7b6d8SPeter Brune   PetscErrorCode       ierr;
271*422a814eSBarry Smith   SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
272f1c6b773SPeter Brune   SNESLineSearch       linesearch;
273b7281c8aSPeter Brune   SNESConvergedReason  reason;
27488d374b2SPeter Brune 
275fef7b6d8SPeter Brune   PetscFunctionBegin;
276c579b300SPatrick Farrell 
277c579b300SPatrick Farrell   if (snes->xl || snes->xu || snes->ops->computevariablebounds) {
278c579b300SPatrick Farrell     SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
279c579b300SPatrick Farrell   }
280c579b300SPatrick Farrell 
281fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
282fef7b6d8SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
283fef7b6d8SPeter Brune 
284fef7b6d8SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
285fef7b6d8SPeter Brune   X      = snes->vec_sol;            /* X^n */
28632cc618eSPeter Brune   dXold  = snes->work[0];            /* The previous iterate of X */
28788d374b2SPeter Brune   dX     = snes->work[1];            /* the preconditioned direction */
288169e2be8SPeter Brune   lX     = snes->vec_sol_update;     /* the conjugate direction */
289fef7b6d8SPeter Brune   F      = snes->vec_func;           /* residual vector */
290fef7b6d8SPeter Brune 
2917601faf0SJed Brown   ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
2929e764e56SPeter Brune 
293e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
294fef7b6d8SPeter Brune   snes->iter = 0;
295fef7b6d8SPeter Brune   snes->norm = 0.;
296e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
297fef7b6d8SPeter Brune 
298bbd5d0b3SPeter Brune   /* compute the initial function and preconditioned update dX */
299a71552e2SPeter Brune 
300a71552e2SPeter Brune   if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
301be95d8f1SBarry Smith     ierr = SNESApplyNPC(snes,X,NULL,dX);CHKERRQ(ierr);
302b7281c8aSPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
303b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
304b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
305b7281c8aSPeter Brune       PetscFunctionReturn(0);
306b7281c8aSPeter Brune     }
307a71552e2SPeter Brune     ierr = VecCopy(dX,F);CHKERRQ(ierr);
308a71552e2SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
309a71552e2SPeter Brune   } else {
310e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
311fef7b6d8SPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
3121aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
313c1c75074SPeter Brune 
314fef7b6d8SPeter Brune     /* convergence test */
315a71552e2SPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
316*422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
317a71552e2SPeter Brune     ierr = VecCopy(F,dX);CHKERRQ(ierr);
318a71552e2SPeter Brune   }
319a71552e2SPeter Brune   if (snes->pc) {
320a71552e2SPeter Brune     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
321be95d8f1SBarry Smith       ierr = SNESApplyNPC(snes,X,F,dX);CHKERRQ(ierr);
322b7281c8aSPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
323b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
324b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
325b7281c8aSPeter Brune         PetscFunctionReturn(0);
326b7281c8aSPeter Brune       }
327a71552e2SPeter Brune     }
328a71552e2SPeter Brune   }
329a71552e2SPeter Brune   ierr = VecCopy(dX,lX);CHKERRQ(ierr);
33032cc618eSPeter Brune   ierr = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr);
331a71552e2SPeter Brune 
332a71552e2SPeter Brune 
333e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
334fef7b6d8SPeter Brune   snes->norm = fnorm;
335e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
336a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
337fef7b6d8SPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
338fef7b6d8SPeter Brune 
339fef7b6d8SPeter Brune   /* test convergence */
340fef7b6d8SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
341fef7b6d8SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
342fef7b6d8SPeter Brune 
343fef7b6d8SPeter Brune   /* Call general purpose update function */
344fef7b6d8SPeter Brune   if (snes->ops->update) {
345fef7b6d8SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
346fef7b6d8SPeter Brune   }
347fef7b6d8SPeter Brune 
348fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
349bbd5d0b3SPeter Brune 
35009c08436SPeter Brune   for (i = 1; i < maxits + 1; i++) {
35129c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
3520a844d1aSPeter Brune     if (ncg->type != SNES_NCG_FR) {
35332cc618eSPeter Brune       ierr = VecCopy(dX, dXold);CHKERRQ(ierr);
35488d374b2SPeter Brune     }
355f1c6b773SPeter Brune     ierr = SNESLineSearchApply(linesearch,X,F,&fnorm,lX);CHKERRQ(ierr);
356*422a814eSBarry Smith     ierr = SNESLineSearchGetReason(linesearch, &lsresult);CHKERRQ(ierr);
357*422a814eSBarry Smith     ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
358*422a814eSBarry Smith     if (lsresult) {
359fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
360fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3615d10c001SPeter Brune         PetscFunctionReturn(0);
362fef7b6d8SPeter Brune       }
363fef7b6d8SPeter Brune     }
364fef7b6d8SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
365fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
3665d10c001SPeter Brune       PetscFunctionReturn(0);
367fef7b6d8SPeter Brune     }
368fef7b6d8SPeter Brune     /* Monitor convergence */
369e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
370169e2be8SPeter Brune     snes->iter = i;
371fef7b6d8SPeter Brune     snes->norm = fnorm;
372e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
373a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
374fef7b6d8SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
375302dbbaeSPeter Brune 
376fef7b6d8SPeter Brune     /* Test for convergence */
377d484d688SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3785d10c001SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
379302dbbaeSPeter Brune 
380302dbbaeSPeter Brune     /* Call general purpose update function */
381302dbbaeSPeter Brune     if (snes->ops->update) {
382302dbbaeSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
383302dbbaeSPeter Brune     }
384a71552e2SPeter Brune     if (snes->pc) {
385a71552e2SPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
386be95d8f1SBarry Smith         ierr = SNESApplyNPC(snes,X,NULL,dX);CHKERRQ(ierr);
387b7281c8aSPeter Brune         ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
388b7281c8aSPeter Brune         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
389b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
390b7281c8aSPeter Brune           PetscFunctionReturn(0);
391b7281c8aSPeter Brune         }
392a71552e2SPeter Brune         ierr = VecCopy(dX,F);CHKERRQ(ierr);
393a71552e2SPeter Brune       } else {
394be95d8f1SBarry Smith         ierr = SNESApplyNPC(snes,X,F,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         }
400302dbbaeSPeter Brune       }
401a3685310SPeter Brune     } else {
402a3685310SPeter Brune       ierr = VecCopy(F,dX);CHKERRQ(ierr);
403302dbbaeSPeter Brune     }
404302dbbaeSPeter Brune 
40529c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
4060a844d1aSPeter Brune     switch (ncg->type) {
4070a844d1aSPeter Brune     case SNES_NCG_FR: /* Fletcher-Reeves */
40832cc618eSPeter Brune       dXolddotdXold= dXdotdX;
40932cc618eSPeter Brune       ierr         = VecDot(dX, dX, &dXdotdX);CHKERRQ(ierr);
41032cc618eSPeter Brune       beta         = PetscRealPart(dXdotdX / dXolddotdXold);
411dfb256c7SPeter Brune       break;
4120a844d1aSPeter Brune     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
41332cc618eSPeter Brune       dXolddotdXold= dXdotdX;
41432cc618eSPeter Brune       ierr         = VecDotBegin(dX, dX, &dXdotdXold);CHKERRQ(ierr);
41532cc618eSPeter Brune       ierr         = VecDotBegin(dXold, dX, &dXdotdXold);CHKERRQ(ierr);
41632cc618eSPeter Brune       ierr         = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr);
41732cc618eSPeter Brune       ierr         = VecDotEnd(dXold, dX, &dXdotdXold);CHKERRQ(ierr);
41832cc618eSPeter Brune       beta         = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
419dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
420dfb256c7SPeter Brune       break;
4210a844d1aSPeter Brune     case SNES_NCG_HS: /* Hestenes-Stiefel */
42232cc618eSPeter Brune       ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr);
42332cc618eSPeter Brune       ierr = VecDotBegin(dX, dXold, &dXdotdXold);CHKERRQ(ierr);
42432cc618eSPeter Brune       ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr);
42532cc618eSPeter Brune       ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr);
42632cc618eSPeter Brune       ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr);
42732cc618eSPeter Brune       ierr = VecDotEnd(dX, dXold, &dXdotdXold);CHKERRQ(ierr);
42832cc618eSPeter Brune       ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr);
42932cc618eSPeter Brune       ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr);
43032cc618eSPeter Brune       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
431dfb256c7SPeter Brune       break;
4320a844d1aSPeter Brune     case SNES_NCG_DY: /* Dai-Yuan */
43332cc618eSPeter Brune       ierr = VecDotBegin(dX, dX, &dXdotdX);CHKERRQ(ierr);
43432cc618eSPeter Brune       ierr = VecDotBegin(lX, dX, &lXdotdX);CHKERRQ(ierr);
43532cc618eSPeter Brune       ierr = VecDotBegin(lX, dXold, &lXdotdXold);CHKERRQ(ierr);
43632cc618eSPeter Brune       ierr = VecDotEnd(dX, dX, &dXdotdX);CHKERRQ(ierr);
43732cc618eSPeter Brune       ierr = VecDotEnd(lX, dX, &lXdotdX);CHKERRQ(ierr);
43832cc618eSPeter Brune       ierr = VecDotEnd(lX, dXold, &lXdotdXold);CHKERRQ(ierr);
43932cc618eSPeter Brune       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));CHKERRQ(ierr);
440dfb256c7SPeter Brune       break;
4410a844d1aSPeter Brune     case SNES_NCG_CD: /* Conjugate Descent */
44232cc618eSPeter Brune       ierr = VecDotBegin(dX, dX, &dXdotdX);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, dXold, &lXdotdXold);CHKERRQ(ierr);
44632cc618eSPeter Brune       beta = PetscRealPart(dXdotdX / lXdotdXold);CHKERRQ(ierr);
447dfb256c7SPeter Brune       break;
448dfb256c7SPeter Brune     }
449dfb256c7SPeter Brune     if (ncg->monitor) {
4507904a332SBarry Smith       ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta);CHKERRQ(ierr);
451dfb256c7SPeter Brune     }
452dfb256c7SPeter Brune     ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr);
453fef7b6d8SPeter Brune   }
454fef7b6d8SPeter Brune   ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
455fef7b6d8SPeter Brune   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
456fef7b6d8SPeter Brune   PetscFunctionReturn(0);
457fef7b6d8SPeter Brune }
458fef7b6d8SPeter Brune 
4590a844d1aSPeter Brune 
4600a844d1aSPeter Brune 
461fef7b6d8SPeter Brune /*MC
462fef7b6d8SPeter Brune   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
463fef7b6d8SPeter Brune 
464fef7b6d8SPeter Brune   Level: beginner
465fef7b6d8SPeter Brune 
466fef7b6d8SPeter Brune   Options Database:
4671867fe5bSPeter Brune +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter.
46841a8cb50SPeter Brune .   -snes_linesearch_type <cp,l2,basic> - Line search type.
4691867fe5bSPeter Brune -   -snes_ncg_monitor - Print relevant information about the ncg iteration.
470fef7b6d8SPeter Brune 
471fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
472fef7b6d8SPeter Brune gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
473fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x.
474fef7b6d8SPeter Brune 
475fef7b6d8SPeter Brune 
47604d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN
477fef7b6d8SPeter Brune M*/
478fef7b6d8SPeter Brune #undef __FUNCT__
479fef7b6d8SPeter Brune #define __FUNCT__ "SNESCreate_NCG"
4808cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
481fef7b6d8SPeter Brune {
482fef7b6d8SPeter Brune   PetscErrorCode ierr;
483ea630c6eSPeter Brune   SNES_NCG       * neP;
484fef7b6d8SPeter Brune 
485fef7b6d8SPeter Brune   PetscFunctionBegin;
486fef7b6d8SPeter Brune   snes->ops->destroy        = SNESDestroy_NCG;
487fef7b6d8SPeter Brune   snes->ops->setup          = SNESSetUp_NCG;
488fef7b6d8SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
489fef7b6d8SPeter Brune   snes->ops->view           = SNESView_NCG;
490fef7b6d8SPeter Brune   snes->ops->solve          = SNESSolve_NCG;
491fef7b6d8SPeter Brune   snes->ops->reset          = SNESReset_NCG;
492fef7b6d8SPeter Brune 
493fef7b6d8SPeter Brune   snes->usesksp = PETSC_FALSE;
494fef7b6d8SPeter Brune   snes->usespc  = PETSC_TRUE;
495a71552e2SPeter Brune   snes->pcside  = PC_LEFT;
496fef7b6d8SPeter Brune 
49788976e71SPeter Brune   if (!snes->tolerancesset) {
4980e444f03SPeter Brune     snes->max_funcs = 30000;
4990e444f03SPeter Brune     snes->max_its   = 10000;
500c522fa08SPeter Brune     snes->stol      = 1e-20;
50188976e71SPeter Brune   }
5020e444f03SPeter Brune 
503b00a9115SJed Brown   ierr         = PetscNewLog(snes,&neP);CHKERRQ(ierr);
504fef7b6d8SPeter Brune   snes->data   = (void*) neP;
5050298fd71SBarry Smith   neP->monitor = NULL;
5060a844d1aSPeter Brune   neP->type    = SNES_NCG_PRP;
507bdf89e91SBarry Smith   ierr         = PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);CHKERRQ(ierr);
508fef7b6d8SPeter Brune   PetscFunctionReturn(0);
509fef7b6d8SPeter Brune }
510