xref: /petsc/src/ksp/ksp/impls/fcg/fcg.c (revision e0f5bfbec699682fa3e8b8532b1176849ea4e12a)
1 /*
2     This file implements the FCG (Flexible Conjugate Gradient) method
3 
4 */
5 
6 #include <../src/ksp/ksp/impls/fcg/fcgimpl.h> /*I  "petscksp.h"  I*/
7 extern PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP, PetscReal *, PetscReal *);
8 extern PetscErrorCode KSPComputeEigenvalues_CG(KSP, PetscInt, PetscReal *, PetscReal *, PetscInt *);
9 
10 #define KSPFCG_DEFAULT_MMAX       30 /* maximum number of search directions to keep */
11 #define KSPFCG_DEFAULT_NPREALLOC  10 /* number of search directions to preallocate */
12 #define KSPFCG_DEFAULT_VECB       5  /* number of search directions to allocate each time new direction vectors are needed */
13 #define KSPFCG_DEFAULT_TRUNCSTRAT KSP_FCD_TRUNC_TYPE_NOTAY
14 
15 static PetscErrorCode KSPAllocateVectors_FCG(KSP ksp, PetscInt nvecsneeded, PetscInt chunksize) {
16   PetscInt i;
17   KSP_FCG *fcg = (KSP_FCG *)ksp->data;
18   PetscInt nnewvecs, nvecsprev;
19 
20   PetscFunctionBegin;
21   /* Allocate enough new vectors to add chunksize new vectors, reach nvecsneedtotal, or to reach mmax+1, whichever is smallest */
22   if (fcg->nvecs < PetscMin(fcg->mmax + 1, nvecsneeded)) {
23     nvecsprev = fcg->nvecs;
24     nnewvecs  = PetscMin(PetscMax(nvecsneeded - fcg->nvecs, chunksize), fcg->mmax + 1 - fcg->nvecs);
25     PetscCall(KSPCreateVecs(ksp, nnewvecs, &fcg->pCvecs[fcg->nchunks], 0, NULL));
26     PetscCall(KSPCreateVecs(ksp, nnewvecs, &fcg->pPvecs[fcg->nchunks], 0, NULL));
27     fcg->nvecs += nnewvecs;
28     for (i = 0; i < nnewvecs; ++i) {
29       fcg->Cvecs[nvecsprev + i] = fcg->pCvecs[fcg->nchunks][i];
30       fcg->Pvecs[nvecsprev + i] = fcg->pPvecs[fcg->nchunks][i];
31     }
32     fcg->chunksizes[fcg->nchunks] = nnewvecs;
33     ++fcg->nchunks;
34   }
35   PetscFunctionReturn(0);
36 }
37 
38 static PetscErrorCode KSPSetUp_FCG(KSP ksp) {
39   KSP_FCG       *fcg      = (KSP_FCG *)ksp->data;
40   PetscInt       maxit    = ksp->max_it;
41   const PetscInt nworkstd = 2;
42 
43   PetscFunctionBegin;
44 
45   /* Allocate "standard" work vectors (not including the basis and transformed basis vectors) */
46   PetscCall(KSPSetWorkVecs(ksp, nworkstd));
47 
48   /* Allocated space for pointers to additional work vectors
49    note that mmax is the number of previous directions, so we add 1 for the current direction,
50    and an extra 1 for the prealloc (which might be empty) */
51   PetscCall(PetscMalloc5(fcg->mmax + 1, &fcg->Pvecs, fcg->mmax + 1, &fcg->Cvecs, fcg->mmax + 1, &fcg->pPvecs, fcg->mmax + 1, &fcg->pCvecs, fcg->mmax + 2, &fcg->chunksizes));
52 
53   /* If the requested number of preallocated vectors is greater than mmax reduce nprealloc */
54   if (fcg->nprealloc > fcg->mmax + 1) PetscCall(PetscInfo(NULL, "Requested nprealloc=%" PetscInt_FMT " is greater than m_max+1=%" PetscInt_FMT ". Resetting nprealloc = m_max+1.\n", fcg->nprealloc, fcg->mmax + 1));
55 
56   /* Preallocate additional work vectors */
57   PetscCall(KSPAllocateVectors_FCG(ksp, fcg->nprealloc, fcg->nprealloc));
58   /*
59   If user requested computations of eigenvalues then allocate work
60   work space needed
61   */
62   if (ksp->calc_sings) {
63     /* get space to store tridiagonal matrix for Lanczos */
64     PetscCall(PetscMalloc4(maxit, &fcg->e, maxit, &fcg->d, maxit, &fcg->ee, maxit, &fcg->dd));
65 
66     ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
67     ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_CG;
68   }
69   PetscFunctionReturn(0);
70 }
71 
72 static PetscErrorCode KSPSolve_FCG(KSP ksp) {
73   PetscInt    i, k, idx, mi;
74   KSP_FCG    *fcg   = (KSP_FCG *)ksp->data;
75   PetscScalar alpha = 0.0, beta = 0.0, dpi, s;
76   PetscReal   dp = 0.0;
77   Vec         B, R, Z, X, Pcurr, Ccurr;
78   Mat         Amat, Pmat;
79   PetscInt    eigs          = ksp->calc_sings; /* Variables for eigen estimation - START*/
80   PetscInt    stored_max_it = ksp->max_it;
81   PetscScalar alphaold = 0, betaold = 1.0, *e = NULL, *d = NULL; /* Variables for eigen estimation  - FINISH */
82 
83   PetscFunctionBegin;
84 
85 #define VecXDot(x, y, a)     (((fcg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x, y, a) : VecTDot(x, y, a))
86 #define VecXMDot(a, b, c, d) (((fcg->type) == (KSP_CG_HERMITIAN)) ? VecMDot(a, b, c, d) : VecMTDot(a, b, c, d))
87 
88   X = ksp->vec_sol;
89   B = ksp->vec_rhs;
90   R = ksp->work[0];
91   Z = ksp->work[1];
92 
93   PetscCall(PCGetOperators(ksp->pc, &Amat, &Pmat));
94   if (eigs) {
95     e    = fcg->e;
96     d    = fcg->d;
97     e[0] = 0.0;
98   }
99   /* Compute initial residual needed for convergence check*/
100   ksp->its = 0;
101   if (!ksp->guess_zero) {
102     PetscCall(KSP_MatMult(ksp, Amat, X, R));
103     PetscCall(VecAYPX(R, -1.0, B)); /*   r <- b - Ax     */
104   } else {
105     PetscCall(VecCopy(B, R)); /*   r <- b (x is 0) */
106   }
107   switch (ksp->normtype) {
108   case KSP_NORM_PRECONDITIONED:
109     PetscCall(KSP_PCApply(ksp, R, Z));  /*   z <- Br         */
110     PetscCall(VecNorm(Z, NORM_2, &dp)); /*   dp <- dqrt(z'*z) = sqrt(e'*A'*B'*B*A*e)     */
111     KSPCheckNorm(ksp, dp);
112     break;
113   case KSP_NORM_UNPRECONDITIONED:
114     PetscCall(VecNorm(R, NORM_2, &dp)); /*   dp <- sqrt(r'*r) = sqrt(e'*A'*A*e)     */
115     KSPCheckNorm(ksp, dp);
116     break;
117   case KSP_NORM_NATURAL:
118     PetscCall(KSP_PCApply(ksp, R, Z)); /*   z <- Br         */
119     PetscCall(VecXDot(R, Z, &s));
120     KSPCheckDot(ksp, s);
121     dp = PetscSqrtReal(PetscAbsScalar(s)); /*   dp <- sqrt(r'*z) = sqrt(e'*A'*B*A*e)  */
122     break;
123   case KSP_NORM_NONE: dp = 0.0; break;
124   default: SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s", KSPNormTypes[ksp->normtype]);
125   }
126 
127   /* Initial Convergence Check */
128   PetscCall(KSPLogResidualHistory(ksp, dp));
129   PetscCall(KSPMonitor(ksp, 0, dp));
130   ksp->rnorm = dp;
131   if (ksp->normtype == KSP_NORM_NONE) {
132     PetscCall(KSPConvergedSkip(ksp, 0, dp, &ksp->reason, ksp->cnvP));
133   } else {
134     PetscCall((*ksp->converged)(ksp, 0, dp, &ksp->reason, ksp->cnvP));
135   }
136   if (ksp->reason) PetscFunctionReturn(0);
137 
138   /* Apply PC if not already done for convergence check */
139   if (ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->normtype == KSP_NORM_NONE) { PetscCall(KSP_PCApply(ksp, R, Z)); /*   z <- Br         */ }
140 
141   i = 0;
142   do {
143     ksp->its = i + 1;
144 
145     /*  If needbe, allocate a new chunk of vectors in P and C */
146     PetscCall(KSPAllocateVectors_FCG(ksp, i + 1, fcg->vecb));
147 
148     /* Note that we wrap around and start clobbering old vectors */
149     idx   = i % (fcg->mmax + 1);
150     Pcurr = fcg->Pvecs[idx];
151     Ccurr = fcg->Cvecs[idx];
152 
153     /* number of old directions to orthogonalize against */
154     switch (fcg->truncstrat) {
155     case KSP_FCD_TRUNC_TYPE_STANDARD: mi = fcg->mmax; break;
156     case KSP_FCD_TRUNC_TYPE_NOTAY: mi = ((i - 1) % fcg->mmax) + 1; break;
157     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unrecognized Truncation Strategy");
158     }
159 
160     /* Compute a new column of P (Currently does not support modified G-S or iterative refinement)*/
161     PetscCall(VecCopy(Z, Pcurr));
162 
163     {
164       PetscInt l, ndots;
165 
166       l     = PetscMax(0, i - mi);
167       ndots = i - l;
168       if (ndots) {
169         PetscInt     j;
170         Vec         *Pold, *Cold;
171         PetscScalar *dots;
172 
173         PetscCall(PetscMalloc3(ndots, &dots, ndots, &Cold, ndots, &Pold));
174         for (k = l, j = 0; j < ndots; ++k, ++j) {
175           idx     = k % (fcg->mmax + 1);
176           Cold[j] = fcg->Cvecs[idx];
177           Pold[j] = fcg->Pvecs[idx];
178         }
179         PetscCall(VecXMDot(Z, ndots, Cold, dots));
180         for (k = 0; k < ndots; ++k) dots[k] = -dots[k];
181         PetscCall(VecMAXPY(Pcurr, ndots, dots, Pold));
182         PetscCall(PetscFree3(dots, Cold, Pold));
183       }
184     }
185 
186     /* Update X and R */
187     betaold = beta;
188     PetscCall(VecXDot(Pcurr, R, &beta)); /*  beta <- pi'*r       */
189     KSPCheckDot(ksp, beta);
190     PetscCall(KSP_MatMult(ksp, Amat, Pcurr, Ccurr)); /*  w <- A*pi (stored in ci)   */
191     PetscCall(VecXDot(Pcurr, Ccurr, &dpi));          /*  dpi <- pi'*w        */
192     alphaold = alpha;
193     alpha    = beta / dpi;                /*  alpha <- beta/dpi    */
194     PetscCall(VecAXPY(X, alpha, Pcurr));  /*  x <- x + alpha * pi  */
195     PetscCall(VecAXPY(R, -alpha, Ccurr)); /*  r <- r - alpha * wi  */
196 
197     /* Compute norm for convergence check */
198     switch (ksp->normtype) {
199     case KSP_NORM_PRECONDITIONED:
200       PetscCall(KSP_PCApply(ksp, R, Z));  /*   z <- Br             */
201       PetscCall(VecNorm(Z, NORM_2, &dp)); /*   dp <- sqrt(z'*z) = sqrt(e'*A'*B'*B*A*e)  */
202       KSPCheckNorm(ksp, dp);
203       break;
204     case KSP_NORM_UNPRECONDITIONED:
205       PetscCall(VecNorm(R, NORM_2, &dp)); /*   dp <- sqrt(r'*r) = sqrt(e'*A'*A*e)   */
206       KSPCheckNorm(ksp, dp);
207       break;
208     case KSP_NORM_NATURAL:
209       PetscCall(KSP_PCApply(ksp, R, Z)); /*   z <- Br             */
210       PetscCall(VecXDot(R, Z, &s));
211       KSPCheckDot(ksp, s);
212       dp = PetscSqrtReal(PetscAbsScalar(s)); /*   dp <- sqrt(r'*z) = sqrt(e'*A'*B*A*e)  */
213       break;
214     case KSP_NORM_NONE: dp = 0.0; break;
215     default: SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s", KSPNormTypes[ksp->normtype]);
216     }
217 
218     /* Check for convergence */
219     ksp->rnorm = dp;
220     PetscCall(KSPLogResidualHistory(ksp, dp));
221     PetscCall(KSPMonitor(ksp, i + 1, dp));
222     PetscCall((*ksp->converged)(ksp, i + 1, dp, &ksp->reason, ksp->cnvP));
223     if (ksp->reason) break;
224 
225     /* Apply PC if not already done for convergence check */
226     if (ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->normtype == KSP_NORM_NONE) { PetscCall(KSP_PCApply(ksp, R, Z)); /*   z <- Br         */ }
227 
228     /* Compute current C (which is W/dpi) */
229     PetscCall(VecScale(Ccurr, 1.0 / dpi)); /*   w <- ci/dpi   */
230 
231     if (eigs) {
232       if (i > 0) {
233         PetscCheck(ksp->max_it == stored_max_it, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Can not change maxit AND calculate eigenvalues");
234         e[i] = PetscSqrtReal(PetscAbsScalar(beta / betaold)) / alphaold;
235         d[i] = PetscSqrtReal(PetscAbsScalar(beta / betaold)) * e[i] + 1.0 / alpha;
236       } else {
237         d[i] = PetscSqrtReal(PetscAbsScalar(beta)) * e[i] + 1.0 / alpha;
238       }
239       fcg->ned = ksp->its - 1;
240     }
241     ++i;
242   } while (i < ksp->max_it);
243   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
244   PetscFunctionReturn(0);
245 }
246 
247 static PetscErrorCode KSPDestroy_FCG(KSP ksp) {
248   PetscInt i;
249   KSP_FCG *fcg = (KSP_FCG *)ksp->data;
250 
251   PetscFunctionBegin;
252 
253   /* Destroy "standard" work vecs */
254   VecDestroyVecs(ksp->nwork, &ksp->work);
255 
256   /* Destroy P and C vectors and the arrays that manage pointers to them */
257   if (fcg->nvecs) {
258     for (i = 0; i < fcg->nchunks; ++i) {
259       PetscCall(VecDestroyVecs(fcg->chunksizes[i], &fcg->pPvecs[i]));
260       PetscCall(VecDestroyVecs(fcg->chunksizes[i], &fcg->pCvecs[i]));
261     }
262   }
263   PetscCall(PetscFree5(fcg->Pvecs, fcg->Cvecs, fcg->pPvecs, fcg->pCvecs, fcg->chunksizes));
264   /* free space used for singular value calculations */
265   if (ksp->calc_sings) PetscCall(PetscFree4(fcg->e, fcg->d, fcg->ee, fcg->dd));
266   PetscCall(KSPDestroyDefault(ksp));
267   PetscFunctionReturn(0);
268 }
269 
270 static PetscErrorCode KSPView_FCG(KSP ksp, PetscViewer viewer) {
271   KSP_FCG    *fcg = (KSP_FCG *)ksp->data;
272   PetscBool   iascii, isstring;
273   const char *truncstr;
274 
275   PetscFunctionBegin;
276   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
277   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
278 
279   if (fcg->truncstrat == KSP_FCD_TRUNC_TYPE_STANDARD) truncstr = "Using standard truncation strategy";
280   else if (fcg->truncstrat == KSP_FCD_TRUNC_TYPE_NOTAY) truncstr = "Using Notay's truncation strategy";
281   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Undefined FCG truncation strategy");
282 
283   if (iascii) {
284     PetscCall(PetscViewerASCIIPrintf(viewer, "  m_max=%" PetscInt_FMT "\n", fcg->mmax));
285     PetscCall(PetscViewerASCIIPrintf(viewer, "  preallocated %" PetscInt_FMT " directions\n", PetscMin(fcg->nprealloc, fcg->mmax + 1)));
286     PetscCall(PetscViewerASCIIPrintf(viewer, "  %s\n", truncstr));
287   } else if (isstring) {
288     PetscCall(PetscViewerStringSPrintf(viewer, "m_max %" PetscInt_FMT " nprealloc %" PetscInt_FMT " %s", fcg->mmax, fcg->nprealloc, truncstr));
289   }
290   PetscFunctionReturn(0);
291 }
292 
293 /*@
294   KSPFCGSetMmax - set the maximum number of previous directions FCG will store for orthogonalization
295 
296   Note: mmax + 1 directions are stored (mmax previous ones along with a current one)
297   and whether all are used in each iteration also depends on the truncation strategy
298   (see KSPFCGSetTruncationType())
299 
300   Logically Collective on ksp
301 
302   Input Parameters:
303 +  ksp - the Krylov space context
304 -  mmax - the maximum number of previous directions to orthogonalize againt
305 
306   Level: intermediate
307 
308 .seealso: `KSPFCG`, `KSPFCGGetTruncationType()`, `KSPFCGGetNprealloc()`
309 @*/
310 PetscErrorCode KSPFCGSetMmax(KSP ksp, PetscInt mmax) {
311   KSP_FCG *fcg = (KSP_FCG *)ksp->data;
312 
313   PetscFunctionBegin;
314   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
315   PetscValidLogicalCollectiveInt(ksp, mmax, 2);
316   fcg->mmax = mmax;
317   PetscFunctionReturn(0);
318 }
319 
320 /*@
321   KSPFCGGetMmax - get the maximum number of previous directions FCG will store
322 
323   Note: FCG stores mmax+1 directions at most (mmax previous ones, and one current one)
324 
325    Not Collective
326 
327    Input Parameter:
328 .  ksp - the Krylov space context
329 
330    Output Parameter:
331 .  mmax - the maximum number of previous directions allowed for orthogonalization
332 
333    Level: intermediate
334 
335 .seealso: `KSPFCG`, `KSPFCGGetTruncationType()`, `KSPFCGGetNprealloc()`, `KSPFCGSetMmax()`
336 @*/
337 
338 PetscErrorCode KSPFCGGetMmax(KSP ksp, PetscInt *mmax) {
339   KSP_FCG *fcg = (KSP_FCG *)ksp->data;
340 
341   PetscFunctionBegin;
342   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
343   *mmax = fcg->mmax;
344   PetscFunctionReturn(0);
345 }
346 
347 /*@
348   KSPFCGSetNprealloc - set the number of directions to preallocate with FCG
349 
350   Logically Collective on ksp
351 
352   Input Parameters:
353 +  ksp - the Krylov space context
354 -  nprealloc - the number of vectors to preallocate
355 
356   Level: advanced
357 
358   Options Database:
359 . -ksp_fcg_nprealloc <N> - number of directions to preallocate
360 
361 .seealso: `KSPFCG`, `KSPFCGGetTruncationType()`, `KSPFCGGetNprealloc()`
362 @*/
363 PetscErrorCode KSPFCGSetNprealloc(KSP ksp, PetscInt nprealloc) {
364   KSP_FCG *fcg = (KSP_FCG *)ksp->data;
365 
366   PetscFunctionBegin;
367   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
368   PetscValidLogicalCollectiveInt(ksp, nprealloc, 2);
369   PetscCheck(nprealloc <= fcg->mmax + 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot preallocate more than m_max+1 vectors");
370   fcg->nprealloc = nprealloc;
371   PetscFunctionReturn(0);
372 }
373 
374 /*@
375   KSPFCGGetNprealloc - get the number of directions preallocate by FCG
376 
377    Not Collective
378 
379    Input Parameter:
380 .  ksp - the Krylov space context
381 
382    Output Parameter:
383 .  nprealloc - the number of directions preallocated
384 
385    Level: advanced
386 
387 .seealso: `KSPFCG`, `KSPFCGGetTruncationType()`, `KSPFCGSetNprealloc()`
388 @*/
389 PetscErrorCode KSPFCGGetNprealloc(KSP ksp, PetscInt *nprealloc) {
390   KSP_FCG *fcg = (KSP_FCG *)ksp->data;
391 
392   PetscFunctionBegin;
393   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
394   *nprealloc = fcg->nprealloc;
395   PetscFunctionReturn(0);
396 }
397 
398 /*@
399   KSPFCGSetTruncationType - specify how many of its stored previous directions FCG uses during orthoganalization
400 
401   Logically Collective on ksp
402 
403   KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
404   KSP_FCD_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1,..
405 
406   Input Parameters:
407 +  ksp - the Krylov space context
408 -  truncstrat - the choice of strategy
409 
410   Level: intermediate
411 
412   Options Database:
413 . -ksp_fcg_truncation_type <standard, notay> - specify how many of its stored previous directions FCG uses during orthoganalization
414 
415   .seealso: `KSPFCDTruncationType`, `KSPFCGGetTruncationType`
416 @*/
417 PetscErrorCode KSPFCGSetTruncationType(KSP ksp, KSPFCDTruncationType truncstrat) {
418   KSP_FCG *fcg = (KSP_FCG *)ksp->data;
419 
420   PetscFunctionBegin;
421   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
422   PetscValidLogicalCollectiveEnum(ksp, truncstrat, 2);
423   fcg->truncstrat = truncstrat;
424   PetscFunctionReturn(0);
425 }
426 
427 /*@
428   KSPFCGGetTruncationType - get the truncation strategy employed by FCG
429 
430    Not Collective
431 
432    Input Parameter:
433 .  ksp - the Krylov space context
434 
435    Output Parameter:
436 .  truncstrat - the strategy type
437 
438    Level: intermediate
439 
440 .seealso: `KSPFCG`, `KSPFCGSetTruncationType`, `KSPFCDTruncationType`
441 @*/
442 PetscErrorCode KSPFCGGetTruncationType(KSP ksp, KSPFCDTruncationType *truncstrat) {
443   KSP_FCG *fcg = (KSP_FCG *)ksp->data;
444 
445   PetscFunctionBegin;
446   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
447   *truncstrat = fcg->truncstrat;
448   PetscFunctionReturn(0);
449 }
450 
451 static PetscErrorCode KSPSetFromOptions_FCG(KSP ksp, PetscOptionItems *PetscOptionsObject) {
452   KSP_FCG  *fcg = (KSP_FCG *)ksp->data;
453   PetscInt  mmax, nprealloc;
454   PetscBool flg;
455 
456   PetscFunctionBegin;
457   PetscOptionsHeadBegin(PetscOptionsObject, "KSP FCG Options");
458   PetscCall(PetscOptionsInt("-ksp_fcg_mmax", "Maximum number of search directions to store", "KSPFCGSetMmax", fcg->mmax, &mmax, &flg));
459   if (flg) PetscCall(KSPFCGSetMmax(ksp, mmax));
460   PetscCall(PetscOptionsInt("-ksp_fcg_nprealloc", "Number of directions to preallocate", "KSPFCGSetNprealloc", fcg->nprealloc, &nprealloc, &flg));
461   if (flg) PetscCall(KSPFCGSetNprealloc(ksp, nprealloc));
462   PetscCall(PetscOptionsEnum("-ksp_fcg_truncation_type", "Truncation approach for directions", "KSPFCGSetTruncationType", KSPFCDTruncationTypes, (PetscEnum)fcg->truncstrat, (PetscEnum *)&fcg->truncstrat, NULL));
463   PetscOptionsHeadEnd();
464   PetscFunctionReturn(0);
465 }
466 
467 /*MC
468       KSPFCG - Implements the Flexible Conjugate Gradient method (FCG)
469 
470   Options Database Keys:
471 +   -ksp_fcg_mmax <N>  - maximum number of search directions
472 .   -ksp_fcg_nprealloc <N> - number of directions to preallocate
473 -   -ksp_fcg_truncation_type <standard,notay> - truncation approach for directions
474 
475     Contributed by Patrick Sanan
476 
477    Notes:
478    Supports left preconditioning only.
479 
480    Level: beginner
481 
482   References:
483 + * - Notay, Y."Flexible Conjugate Gradients", SIAM J. Sci. Comput. 22:4, 2000
484 - * - Axelsson, O. and Vassilevski, P. S. "A Black Box Generalized Conjugate Gradient Solver with Inner Iterations and Variable step Preconditioning",
485     SIAM J. Matrix Anal. Appl. 12:4, 1991
486 
487  .seealso `:` `KSPGCR`, `KSPFGMRES`, `KSPCG`, `KSPFCGSetMmax()`, `KSPFCGGetMmax()`, `KSPFCGSetNprealloc()`, `KSPFCGGetNprealloc()`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`
488 
489 M*/
490 PETSC_EXTERN PetscErrorCode KSPCreate_FCG(KSP ksp) {
491   KSP_FCG *fcg;
492 
493   PetscFunctionBegin;
494   PetscCall(PetscNew(&fcg));
495 #if !defined(PETSC_USE_COMPLEX)
496   fcg->type = KSP_CG_SYMMETRIC;
497 #else
498   fcg->type = KSP_CG_HERMITIAN;
499 #endif
500   fcg->mmax       = KSPFCG_DEFAULT_MMAX;
501   fcg->nprealloc  = KSPFCG_DEFAULT_NPREALLOC;
502   fcg->nvecs      = 0;
503   fcg->vecb       = KSPFCG_DEFAULT_VECB;
504   fcg->nchunks    = 0;
505   fcg->truncstrat = KSPFCG_DEFAULT_TRUNCSTRAT;
506 
507   ksp->data = (void *)fcg;
508 
509   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 2));
510   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_LEFT, 1));
511   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NATURAL, PC_LEFT, 1));
512   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1));
513 
514   ksp->ops->setup          = KSPSetUp_FCG;
515   ksp->ops->solve          = KSPSolve_FCG;
516   ksp->ops->destroy        = KSPDestroy_FCG;
517   ksp->ops->view           = KSPView_FCG;
518   ksp->ops->setfromoptions = KSPSetFromOptions_FCG;
519   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
520   ksp->ops->buildresidual  = KSPBuildResidualDefault;
521   PetscFunctionReturn(0);
522 }
523