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