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