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