xref: /petsc/src/snes/impls/ncg/snesncg.c (revision 8cc058d9cd56c1ccb3be12a47760ddfc446aaffc)
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 
44*8cc058d9SJed 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;
53bbd5d0b3SPeter Brune   ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr);
54bbd5d0b3SPeter Brune   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
550298fd71SBarry Smith   ierr = SNESLineSearchRegisterDynamic(SNESLINESEARCHNCGLINEAR, NULL,"SNESLineSearchCreate_NCGLinear", SNESLineSearchCreate_NCGLinear);CHKERRQ(ierr);
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) {
86f1c6b773SPeter Brune     ierr = SNESGetSNESLineSearch(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   }
93fef7b6d8SPeter Brune   PetscFunctionReturn(0);
94fef7b6d8SPeter Brune }
95fef7b6d8SPeter Brune 
96fef7b6d8SPeter Brune /*
97fef7b6d8SPeter Brune   SNESView_NCG - Prints info from the SNESNCG data structure.
98fef7b6d8SPeter Brune 
99fef7b6d8SPeter Brune   Input Parameters:
100fef7b6d8SPeter Brune + SNES - the SNES context
101fef7b6d8SPeter Brune - viewer - visualization context
102fef7b6d8SPeter Brune 
103fef7b6d8SPeter Brune   Application Interface Routine: SNESView()
104fef7b6d8SPeter Brune */
105fef7b6d8SPeter Brune #undef __FUNCT__
106fef7b6d8SPeter Brune #define __FUNCT__ "SNESView_NCG"
107fef7b6d8SPeter Brune static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
108fef7b6d8SPeter Brune {
109fef7b6d8SPeter Brune   PetscBool      iascii;
110fef7b6d8SPeter Brune   PetscErrorCode ierr;
111fef7b6d8SPeter Brune 
112fef7b6d8SPeter Brune   PetscFunctionBegin;
113251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
114fef7b6d8SPeter Brune   if (iascii) {
115fef7b6d8SPeter Brune   }
116fef7b6d8SPeter Brune   PetscFunctionReturn(0);
117fef7b6d8SPeter Brune }
118fef7b6d8SPeter Brune 
119fef7b6d8SPeter Brune #undef __FUNCT__
120f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchApply_NCGLinear"
121f1c6b773SPeter Brune PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
122ea630c6eSPeter Brune {
123ea630c6eSPeter Brune   PetscScalar    alpha, ptAp;
124bbd5d0b3SPeter Brune   Vec            X, Y, F, W;
125bbd5d0b3SPeter Brune   SNES           snes;
126ea630c6eSPeter Brune   PetscErrorCode ierr;
127bbd5d0b3SPeter Brune   PetscReal      *fnorm, *xnorm, *ynorm;
128ea630c6eSPeter Brune   MatStructure   flg = DIFFERENT_NONZERO_PATTERN;
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 
140ea630c6eSPeter Brune   /*
141ea630c6eSPeter Brune 
142d2140ca3SPeter Brune    The exact step size for unpreconditioned linear CG is just:
143ea630c6eSPeter Brune    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
144ea630c6eSPeter Brune    */
145a5892487SPeter Brune   ierr  = SNESComputeJacobian(snes, X, &snes->jacobian, &snes->jacobian_pre, &flg);CHKERRQ(ierr);
146ea630c6eSPeter Brune   ierr  = VecDot(F, F, &alpha);CHKERRQ(ierr);
147ea630c6eSPeter Brune   ierr  = MatMult(snes->jacobian, Y, W);CHKERRQ(ierr);
148ea630c6eSPeter Brune   ierr  = VecDot(Y, W, &ptAp);CHKERRQ(ierr);
149ea630c6eSPeter Brune   alpha = alpha / ptAp;
150ce94432eSBarry Smith   ierr  = PetscPrintf(PetscObjectComm((PetscObject)snes), "alpha: %G\n", PetscRealPart(alpha));CHKERRQ(ierr);
151bbd5d0b3SPeter Brune   ierr  = VecAXPY(X, alpha, Y);CHKERRQ(ierr);
152bbd5d0b3SPeter Brune   ierr  = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
153bbd5d0b3SPeter Brune 
154bbd5d0b3SPeter Brune   ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr);
155bbd5d0b3SPeter Brune   ierr = VecNorm(X, NORM_2, xnorm);CHKERRQ(ierr);
156bbd5d0b3SPeter Brune   ierr = VecNorm(Y, NORM_2, ynorm);CHKERRQ(ierr);
157ea630c6eSPeter Brune   PetscFunctionReturn(0);
158ea630c6eSPeter Brune }
159ea630c6eSPeter Brune 
160bbd5d0b3SPeter Brune #undef __FUNCT__
161f1c6b773SPeter Brune #define __FUNCT__ "SNESLineSearchCreate_NCGLinear"
162*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
163bbd5d0b3SPeter Brune {
164bbd5d0b3SPeter Brune   PetscFunctionBegin;
165f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
1660298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
1670298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
1680298fd71SBarry Smith   linesearch->ops->reset          = NULL;
1690298fd71SBarry Smith   linesearch->ops->view           = NULL;
1700298fd71SBarry Smith   linesearch->ops->setup          = NULL;
171bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
172bbd5d0b3SPeter Brune }
173bbd5d0b3SPeter Brune 
17488d374b2SPeter Brune #undef __FUNCT__
1758c40d5fbSBarry Smith #define __FUNCT__ "SNESNCGComputeYtJtF_Private"
17688d374b2SPeter Brune /*
17788d374b2SPeter Brune 
17888d374b2SPeter Brune  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.
17988d374b2SPeter Brune 
18088d374b2SPeter Brune  */
18167392de3SBarry Smith PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
18267392de3SBarry Smith {
18388d374b2SPeter Brune   PetscErrorCode ierr;
18488d374b2SPeter Brune   PetscScalar    ftf, ftg, fty, h;
18567392de3SBarry Smith 
18688d374b2SPeter Brune   PetscFunctionBegin;
18788d374b2SPeter Brune   ierr   = VecDot(F, F, &ftf);CHKERRQ(ierr);
18888d374b2SPeter Brune   ierr   = VecDot(F, Y, &fty);CHKERRQ(ierr);
18988d374b2SPeter Brune   h      = 1e-5*fty / fty;
19088d374b2SPeter Brune   ierr   = VecCopy(X, W);CHKERRQ(ierr);
19188d374b2SPeter Brune   ierr   = VecAXPY(W, -h, Y);CHKERRQ(ierr);          /* this is arbitrary */
19288d374b2SPeter Brune   ierr   = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
19388d374b2SPeter Brune   ierr   = VecDot(G, F, &ftg);CHKERRQ(ierr);
19488d374b2SPeter Brune   *ytJtf = (ftg - ftf) / h;
19588d374b2SPeter Brune   PetscFunctionReturn(0);
19688d374b2SPeter Brune }
19788d374b2SPeter Brune 
1980a844d1aSPeter Brune #undef __FUNCT__
1990a844d1aSPeter Brune #define __FUNCT__ "SNESNCGSetType"
2000a844d1aSPeter Brune /*@
2010a844d1aSPeter Brune     SNESNCGSetType - Sets the conjugate update type for SNESNCG.
2020a844d1aSPeter Brune 
2030a844d1aSPeter Brune     Logically Collective on SNES
2040a844d1aSPeter Brune 
2050a844d1aSPeter Brune     Input Parameters:
2060a844d1aSPeter Brune +   snes - the iterative context
2070a844d1aSPeter Brune -   btype - update type
2080a844d1aSPeter Brune 
2090a844d1aSPeter Brune     Options Database:
2100a844d1aSPeter Brune .   -snes_ncg_type<prp,fr,hs,dy,cd>
2110a844d1aSPeter Brune 
2120a844d1aSPeter Brune     Level: intermediate
2130a844d1aSPeter Brune 
2140a844d1aSPeter Brune     SNESNCGSelectTypes:
2150a844d1aSPeter Brune +   SNES_NCG_FR - Fletcher-Reeves update
2160a844d1aSPeter Brune .   SNES_NCG_PRP - Polak-Ribiere-Polyak update
2170a844d1aSPeter Brune .   SNES_NCG_HS - Hestenes-Steifel update
2180a844d1aSPeter Brune .   SNES_NCG_DY - Dai-Yuan update
2190a844d1aSPeter Brune -   SNES_NCG_CD - Conjugate Descent update
2200a844d1aSPeter Brune 
2210a844d1aSPeter Brune    Notes:
2220a844d1aSPeter Brune    PRP is the default, and the only one that tolerates generalized search directions.
2230a844d1aSPeter Brune 
2240a844d1aSPeter Brune .keywords: SNES, SNESNCG, selection, type, set
2250a844d1aSPeter Brune @*/
2260adebc6cSBarry Smith PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
2270adebc6cSBarry Smith {
2280a844d1aSPeter Brune   PetscErrorCode ierr;
2296e111a19SKarl Rupp 
2300a844d1aSPeter Brune   PetscFunctionBegin;
2310a844d1aSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2320a844d1aSPeter Brune   ierr = PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));CHKERRQ(ierr);
2330a844d1aSPeter Brune   PetscFunctionReturn(0);
2340a844d1aSPeter Brune }
2350a844d1aSPeter Brune 
2360a844d1aSPeter Brune #undef __FUNCT__
2370a844d1aSPeter Brune #define __FUNCT__ "SNESNCGSetType_NCG"
2380adebc6cSBarry Smith PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
2390adebc6cSBarry Smith {
2400a844d1aSPeter Brune   SNES_NCG *ncg = (SNES_NCG*)snes->data;
2416e111a19SKarl Rupp 
2420a844d1aSPeter Brune   PetscFunctionBegin;
2430a844d1aSPeter Brune   ncg->type = btype;
2440a844d1aSPeter Brune   PetscFunctionReturn(0);
2450a844d1aSPeter Brune }
2460a844d1aSPeter Brune 
247fef7b6d8SPeter Brune /*
248fef7b6d8SPeter Brune   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.
249fef7b6d8SPeter Brune 
250fef7b6d8SPeter Brune   Input Parameters:
251fef7b6d8SPeter Brune . snes - the SNES context
252fef7b6d8SPeter Brune 
253fef7b6d8SPeter Brune   Output Parameter:
254fef7b6d8SPeter Brune . outits - number of iterations until termination
255fef7b6d8SPeter Brune 
256fef7b6d8SPeter Brune   Application Interface Routine: SNESSolve()
257fef7b6d8SPeter Brune */
258fef7b6d8SPeter Brune #undef __FUNCT__
259fef7b6d8SPeter Brune #define __FUNCT__ "SNESSolve_NCG"
260fef7b6d8SPeter Brune PetscErrorCode SNESSolve_NCG(SNES snes)
261fef7b6d8SPeter Brune {
262dfb256c7SPeter Brune   SNES_NCG            *ncg = (SNES_NCG*)snes->data;
263bbd5d0b3SPeter Brune   Vec                 X, dX, lX, F, B, Fold;
264bbd5d0b3SPeter Brune   PetscReal           fnorm, ynorm, xnorm, beta = 0.0;
2655d115551SPeter Brune   PetscScalar         dXdotF, dXolddotFold, dXdotFold, lXdotF, lXdotFold;
266fef7b6d8SPeter Brune   PetscInt            maxits, i;
267fef7b6d8SPeter Brune   PetscErrorCode      ierr;
268fef7b6d8SPeter Brune   SNESConvergedReason reason;
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 */
281302dbbaeSPeter Brune   B      = snes->vec_rhs;            /* the right hand side */
282fef7b6d8SPeter Brune 
283f1c6b773SPeter Brune   ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
2849e764e56SPeter Brune 
285fef7b6d8SPeter Brune   ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
286fef7b6d8SPeter Brune   snes->iter = 0;
287fef7b6d8SPeter Brune   snes->norm = 0.;
288fef7b6d8SPeter Brune   ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
289fef7b6d8SPeter Brune 
290bbd5d0b3SPeter Brune   /* compute the initial function and preconditioned update dX */
291e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
292fef7b6d8SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
293fef7b6d8SPeter Brune     if (snes->domainerror) {
294fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
295fef7b6d8SPeter Brune       PetscFunctionReturn(0);
296fef7b6d8SPeter Brune     }
2971aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
2981aa26658SKarl Rupp 
299e4ed7901SPeter Brune   if (!snes->norm_init_set) {
300fef7b6d8SPeter Brune     /* convergence test */
301fef7b6d8SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
302189a9710SBarry Smith     if (PetscIsInfOrNanReal(fnorm)) {
303189a9710SBarry Smith       snes->reason = SNES_DIVERGED_FNORM_NAN;
304189a9710SBarry Smith       PetscFunctionReturn(0);
305189a9710SBarry Smith     }
306e4ed7901SPeter Brune   } else {
307e4ed7901SPeter Brune     fnorm               = snes->norm_init;
308e4ed7901SPeter Brune     snes->norm_init_set = PETSC_FALSE;
309e4ed7901SPeter Brune   }
310fef7b6d8SPeter Brune   ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
311fef7b6d8SPeter Brune   snes->norm = fnorm;
312fef7b6d8SPeter Brune   ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
313a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
314fef7b6d8SPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
315fef7b6d8SPeter Brune 
316fef7b6d8SPeter Brune   /* set parameter for default relative tolerance convergence test */
317fef7b6d8SPeter Brune   snes->ttol = fnorm*snes->rtol;
318fef7b6d8SPeter Brune   /* test convergence */
319fef7b6d8SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
320fef7b6d8SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
321fef7b6d8SPeter Brune 
322fef7b6d8SPeter Brune   /* Call general purpose update function */
323fef7b6d8SPeter Brune   if (snes->ops->update) {
324fef7b6d8SPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
325fef7b6d8SPeter Brune   }
326fef7b6d8SPeter Brune 
327fef7b6d8SPeter Brune   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */
328bbd5d0b3SPeter Brune 
329c40d0f55SPeter Brune   if (snes->pc && snes->pcside == PC_RIGHT) {
33029c94e34SPeter Brune     ierr = VecCopy(X, dX);CHKERRQ(ierr);
331217b9c2eSPeter Brune     ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
332217b9c2eSPeter Brune     ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
33363e7833aSPeter Brune 
33463e7833aSPeter Brune     ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,dX,B,0);CHKERRQ(ierr);
33529c94e34SPeter Brune     ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
33663e7833aSPeter Brune     ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,dX,B,0);CHKERRQ(ierr);
33763e7833aSPeter Brune 
338fef7b6d8SPeter Brune     ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
33951e86f29SPeter Brune     if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
340fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
341fef7b6d8SPeter Brune       PetscFunctionReturn(0);
342fef7b6d8SPeter Brune     }
34329c94e34SPeter Brune     ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
3445d10c001SPeter Brune   } else {
34529c94e34SPeter Brune     ierr = VecCopy(F, dX);CHKERRQ(ierr);
346fef7b6d8SPeter Brune   }
34729c94e34SPeter Brune   ierr = VecCopy(dX, lX);CHKERRQ(ierr);
348bbd5d0b3SPeter Brune   ierr = VecDot(F, dX, &dXdotF);CHKERRQ(ierr);
349bbd5d0b3SPeter Brune   /*
35088d374b2SPeter Brune   } else {
3518c40d5fbSBarry Smith     ierr = SNESNCGComputeYtJtF_Private(snes, X, F, dX, W, G, &dXdotF);CHKERRQ(ierr);
35288d374b2SPeter Brune   }
353bbd5d0b3SPeter Brune    */
35409c08436SPeter Brune   for (i = 1; i < maxits + 1; i++) {
355fef7b6d8SPeter Brune     lsSuccess = PETSC_TRUE;
35629c94e34SPeter Brune     /* some update types require the old update direction or conjugate direction */
3570a844d1aSPeter Brune     if (ncg->type != SNES_NCG_FR) {
3585d115551SPeter Brune       ierr = VecCopy(F, Fold);CHKERRQ(ierr);
35988d374b2SPeter Brune     }
360f1c6b773SPeter Brune     ierr = SNESLineSearchApply(linesearch, X, F, &fnorm, lX);CHKERRQ(ierr);
361f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(linesearch, &lsSuccess);CHKERRQ(ierr);
362fef7b6d8SPeter Brune     if (!lsSuccess) {
363fef7b6d8SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
364fef7b6d8SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
3655d10c001SPeter Brune         PetscFunctionReturn(0);
366fef7b6d8SPeter Brune       }
367fef7b6d8SPeter Brune     }
368fef7b6d8SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
369fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
3705d10c001SPeter Brune       PetscFunctionReturn(0);
371fef7b6d8SPeter Brune     }
372fef7b6d8SPeter Brune     if (snes->domainerror) {
373fef7b6d8SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
374fef7b6d8SPeter Brune       PetscFunctionReturn(0);
375fef7b6d8SPeter Brune     }
376f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
377fef7b6d8SPeter Brune     /* Monitor convergence */
378fef7b6d8SPeter Brune     ierr       = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
379169e2be8SPeter Brune     snes->iter = i;
380fef7b6d8SPeter Brune     snes->norm = fnorm;
381fef7b6d8SPeter Brune     ierr       = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
382a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
383fef7b6d8SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
384302dbbaeSPeter Brune 
385fef7b6d8SPeter Brune     /* Test for convergence */
386d484d688SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3875d10c001SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
388302dbbaeSPeter Brune 
389302dbbaeSPeter Brune     /* Call general purpose update function */
390302dbbaeSPeter Brune     if (snes->ops->update) {
391302dbbaeSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
392302dbbaeSPeter Brune     }
393c40d0f55SPeter Brune     if (snes->pc && snes->pcside == PC_RIGHT) {
394302dbbaeSPeter Brune       ierr = VecCopy(X,dX);CHKERRQ(ierr);
395217b9c2eSPeter Brune       ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
396217b9c2eSPeter Brune       ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
397cf949ffaSPeter Brune       ierr = SNESSolve(snes->pc, B, dX);CHKERRQ(ierr);
398302dbbaeSPeter Brune       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
39951e86f29SPeter Brune       if (reason < 0 && (reason != SNES_DIVERGED_MAX_IT)) {
400302dbbaeSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
401302dbbaeSPeter Brune         PetscFunctionReturn(0);
402302dbbaeSPeter Brune       }
403eb1825c3SPeter Brune       ierr = VecAYPX(dX,-1.0,X);CHKERRQ(ierr);
404a3685310SPeter Brune     } else {
405a3685310SPeter Brune       ierr = VecCopy(F, dX);CHKERRQ(ierr);
406302dbbaeSPeter Brune     }
407302dbbaeSPeter Brune 
40829c94e34SPeter Brune     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
4090a844d1aSPeter Brune     switch (ncg->type) {
4100a844d1aSPeter Brune     case SNES_NCG_FR: /* Fletcher-Reeves */
4115d115551SPeter Brune       dXolddotFold = dXdotF;
412bbd5d0b3SPeter Brune       ierr         = VecDot(dX, dX, &dXdotF);CHKERRQ(ierr);
4135d115551SPeter Brune       beta         = PetscRealPart(dXdotF / dXolddotFold);
414dfb256c7SPeter Brune       break;
4150a844d1aSPeter Brune     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
4165d115551SPeter Brune       dXolddotFold = dXdotF;
41715f0db4aSPeter Brune       ierr         = VecDotBegin(F, dX, &dXdotF);CHKERRQ(ierr);
41815f0db4aSPeter Brune       ierr         = VecDotBegin(Fold, dX, &dXdotFold);CHKERRQ(ierr);
41915f0db4aSPeter Brune       ierr         = VecDotEnd(F, dX, &dXdotF);CHKERRQ(ierr);
42015f0db4aSPeter Brune       ierr         = VecDotEnd(Fold, dX, &dXdotFold);CHKERRQ(ierr);
4215d115551SPeter Brune       beta         = PetscRealPart(((dXdotF - dXdotFold) / dXolddotFold));
422dfb256c7SPeter Brune       if (beta < 0.0) beta = 0.0; /* restart */
423dfb256c7SPeter Brune       break;
4240a844d1aSPeter Brune     case SNES_NCG_HS: /* Hestenes-Stiefel */
42515f0db4aSPeter Brune       ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr);
42615f0db4aSPeter Brune       ierr = VecDotBegin(dX, Fold, &dXdotFold);CHKERRQ(ierr);
42715f0db4aSPeter Brune       ierr = VecDotBegin(lX, F, &lXdotF);CHKERRQ(ierr);
42815f0db4aSPeter Brune       ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr);
42915f0db4aSPeter Brune       ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr);
43015f0db4aSPeter Brune       ierr = VecDotEnd(dX, Fold, &dXdotFold);CHKERRQ(ierr);
43115f0db4aSPeter Brune       ierr = VecDotEnd(lX, F, &lXdotF);CHKERRQ(ierr);
43215f0db4aSPeter Brune       ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr);
4335d115551SPeter Brune       beta = PetscRealPart((dXdotF - dXdotFold) / (lXdotF - lXdotFold));
434dfb256c7SPeter Brune       break;
4350a844d1aSPeter Brune     case SNES_NCG_DY: /* Dai-Yuan */
43615f0db4aSPeter Brune       ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr);
43715f0db4aSPeter Brune       ierr = VecDotBegin(lX, F, &lXdotF);CHKERRQ(ierr);
43815f0db4aSPeter Brune       ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr);
43915f0db4aSPeter Brune       ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr);
44015f0db4aSPeter Brune       ierr = VecDotEnd(lX, F, &lXdotF);CHKERRQ(ierr);
44115f0db4aSPeter Brune       ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr);
4425d115551SPeter Brune       beta = PetscRealPart(dXdotF / (lXdotFold - lXdotF));CHKERRQ(ierr);
443dfb256c7SPeter Brune       break;
4440a844d1aSPeter Brune     case SNES_NCG_CD: /* Conjugate Descent */
44515f0db4aSPeter Brune       ierr = VecDotBegin(dX, F, &dXdotF);CHKERRQ(ierr);
44615f0db4aSPeter Brune       ierr = VecDotBegin(lX, Fold, &lXdotFold);CHKERRQ(ierr);
44715f0db4aSPeter Brune       ierr = VecDotEnd(dX, F, &dXdotF);CHKERRQ(ierr);
44815f0db4aSPeter Brune       ierr = VecDotEnd(lX, Fold, &lXdotFold);CHKERRQ(ierr);
4495d115551SPeter Brune       beta = PetscRealPart(dXdotF / lXdotFold);CHKERRQ(ierr);
450dfb256c7SPeter Brune       break;
451dfb256c7SPeter Brune     }
452dfb256c7SPeter Brune     if (ncg->monitor) {
453646217ecSPeter Brune       ierr = PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", beta);CHKERRQ(ierr);
454dfb256c7SPeter Brune     }
455dfb256c7SPeter Brune     ierr = VecAYPX(lX, beta, dX);CHKERRQ(ierr);
456fef7b6d8SPeter Brune   }
457fef7b6d8SPeter Brune   ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
458fef7b6d8SPeter Brune   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
459fef7b6d8SPeter Brune   PetscFunctionReturn(0);
460fef7b6d8SPeter Brune }
461fef7b6d8SPeter Brune 
4620a844d1aSPeter Brune 
4630a844d1aSPeter Brune 
464fef7b6d8SPeter Brune /*MC
465fef7b6d8SPeter Brune   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.
466fef7b6d8SPeter Brune 
467fef7b6d8SPeter Brune   Level: beginner
468fef7b6d8SPeter Brune 
469fef7b6d8SPeter Brune   Options Database:
4701867fe5bSPeter Brune +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter.
47141a8cb50SPeter Brune .   -snes_linesearch_type <cp,l2,basic> - Line search type.
4721867fe5bSPeter Brune -   -snes_ncg_monitor - Print relevant information about the ncg iteration.
473fef7b6d8SPeter Brune 
474fef7b6d8SPeter Brune Notes: This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
475fef7b6d8SPeter Brune gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
476fef7b6d8SPeter Brune chooses the initial search direction as F(x) for the initial guess x.
477fef7b6d8SPeter Brune 
478fef7b6d8SPeter Brune 
47904d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN
480fef7b6d8SPeter Brune M*/
481fef7b6d8SPeter Brune #undef __FUNCT__
482fef7b6d8SPeter Brune #define __FUNCT__ "SNESCreate_NCG"
483*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
484fef7b6d8SPeter Brune {
485fef7b6d8SPeter Brune   PetscErrorCode ierr;
486ea630c6eSPeter Brune   SNES_NCG       * neP;
487fef7b6d8SPeter Brune 
488fef7b6d8SPeter Brune   PetscFunctionBegin;
489fef7b6d8SPeter Brune   snes->ops->destroy        = SNESDestroy_NCG;
490fef7b6d8SPeter Brune   snes->ops->setup          = SNESSetUp_NCG;
491fef7b6d8SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
492fef7b6d8SPeter Brune   snes->ops->view           = SNESView_NCG;
493fef7b6d8SPeter Brune   snes->ops->solve          = SNESSolve_NCG;
494fef7b6d8SPeter Brune   snes->ops->reset          = SNESReset_NCG;
495fef7b6d8SPeter Brune 
496fef7b6d8SPeter Brune   snes->usesksp = PETSC_FALSE;
497fef7b6d8SPeter Brune   snes->usespc  = PETSC_TRUE;
498fef7b6d8SPeter Brune 
49988976e71SPeter Brune   if (!snes->tolerancesset) {
5000e444f03SPeter Brune     snes->max_funcs = 30000;
5010e444f03SPeter Brune     snes->max_its   = 10000;
502c522fa08SPeter Brune     snes->stol      = 1e-20;
50388976e71SPeter Brune   }
5040e444f03SPeter Brune 
505fef7b6d8SPeter Brune   ierr         = PetscNewLog(snes, SNES_NCG, &neP);CHKERRQ(ierr);
506fef7b6d8SPeter Brune   snes->data   = (void*) neP;
5070298fd71SBarry Smith   neP->monitor = NULL;
5080a844d1aSPeter Brune   neP->type    = SNES_NCG_PRP;
50900de8ff0SBarry Smith   ierr         = PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C","SNESNCGSetType_NCG", SNESNCGSetType_NCG);CHKERRQ(ierr);
510fef7b6d8SPeter Brune   PetscFunctionReturn(0);
511fef7b6d8SPeter Brune }
512