Actual source code: fcg.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  1: /*
  2:     This file implements the FCG (Flexible Conjugate Gradient) method

  4: */

  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*);

 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

 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;

 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:     KSPCreateVecs(ksp,nnewvecs,&fcg->pCvecs[fcg->nchunks],0,NULL);
 30:     PetscLogObjectParents((PetscObject)ksp,nnewvecs,fcg->pCvecs[fcg->nchunks]);
 31:     KSPCreateVecs(ksp,nnewvecs,&fcg->pPvecs[fcg->nchunks],0,NULL);
 32:     PetscLogObjectParents((PetscObject)ksp,nnewvecs,fcg->pPvecs[fcg->nchunks]);
 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:   return(0);
 42: }

 46: PetscErrorCode    KSPSetUp_FCG(KSP ksp)
 47: {
 49:   KSP_FCG        *fcg = (KSP_FCG*)ksp->data;
 50:   PetscInt       maxit = ksp->max_it;
 51:   const PetscInt nworkstd = 2;


 55:   /* Allocate "standard" work vectors (not including the basis and transformed basis vectors) */
 56:   KSPSetWorkVecs(ksp,nworkstd);

 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:   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);
 62:   PetscLogObjectMemory((PetscObject)ksp,2*(fcg->mmax+1)*sizeof(Vec*) + 2*(fcg->mmax + 1)*sizeof(Vec**) + (fcg->mmax + 2)*sizeof(PetscInt));

 64:   /* If the requested number of preallocated vectors is greater than mmax reduce nprealloc */
 65:   if(fcg->nprealloc > fcg->mmax+1){
 66:     PetscInfo2(NULL,"Requested nprealloc=%d is greater than m_max+1=%d. Resetting nprealloc = m_max+1.\n",fcg->nprealloc, fcg->mmax+1);
 67:   }

 69:   /* Preallocate additional work vectors */
 70:   KSPAllocateVectors_FCG(ksp,fcg->nprealloc,fcg->nprealloc);
 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:     PetscMalloc4(maxit,&fcg->e,maxit,&fcg->d,maxit,&fcg->ee,maxit,&fcg->dd);
 78:     PetscLogObjectMemory((PetscObject)ksp,2*(maxit+1)*(sizeof(PetscScalar)+sizeof(PetscReal)));

 80:     ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
 81:     ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_CG;
 82:   }
 83:   return(0);
 84: }

 88: PetscErrorCode KSPSolve_FCG(KSP ksp)
 89: {
 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 */


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))

106:   X             = ksp->vec_sol;
107:   B             = ksp->vec_rhs;
108:   R             = ksp->work[0];
109:   Z             = ksp->work[1];

111:   PCGetOperators(ksp->pc,&Amat,&Pmat);
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:     KSP_MatMult(ksp,Amat,X,R);
117:     VecAYPX(R,-1.0,B);                    /*   r <- b - Ax     */
118:   } else {
119:     VecCopy(B,R);                         /*   r <- b (x is 0) */
120:   }
121:   switch (ksp->normtype) {
122:     case KSP_NORM_PRECONDITIONED:
123:       KSP_PCApply(ksp,R,Z);               /*   z <- Br         */
124:       VecNorm(Z,NORM_2,&dp);              /*   dp <- dqrt(z'*z) = sqrt(e'*A'*B'*B*A*e)     */
125:       break;
126:     case KSP_NORM_UNPRECONDITIONED:
127:       VecNorm(R,NORM_2,&dp);              /*   dp <- sqrt(r'*r) = sqrt(e'*A'*A*e)     */
128:       break;
129:     case KSP_NORM_NATURAL:
130:       KSP_PCApply(ksp,R,Z);               /*   z <- Br         */
131:       VecXDot(R,Z,&s);
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:   }

140:   /* Initial Convergence Check */
141:   KSPLogResidualHistory(ksp,dp);
142:   KSPMonitor(ksp,0,dp);
143:   ksp->rnorm = dp;
144:   if (ksp->normtype == KSP_NORM_NONE) {
145:     KSPConvergedSkip(ksp,0,dp,&ksp->reason,ksp->cnvP);
146:   } else {
147:     (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
148:   }
149:   if (ksp->reason) return(0);

151:   /* Apply PC if not already done for convergence check */
152:   if(ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->normtype == KSP_NORM_NONE){
153:     KSP_PCApply(ksp,R,Z);               /*   z <- Br         */
154:   }

156:   i = 0;
157:   do {
158:     ksp->its = i+1;

160:     /*  If needbe, allocate a new chunk of vectors in P and C */
161:     KSPAllocateVectors_FCG(ksp,i+1,fcg->vecb);

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];

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:     }

180:     /* Compute a new column of P (Currently does not support modified G-S or iterative refinement)*/
181:     VecCopy(Z,Pcurr);

183:     {
184:       PetscInt l,ndots;

186:       l = PetscMax(0,i-mi);
187:       ndots = i-l;
188:       if (ndots){
189:         PetscInt    j;
190:         Vec         *Pold,  *Cold;
191:         PetscScalar *dots;

193:         PetscMalloc3(ndots,&dots,ndots,&Cold,ndots,&Pold);
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:         VecXMDot(Z,ndots,Cold,dots);
200:         for(k=0;k<ndots;++k){
201:           dots[k] = -dots[k];
202:         }
203:         VecMAXPY(Pcurr,ndots,dots,Pold);
204:         PetscFree3(dots,Cold,Pold);
205:       }
206:     }

208:     /* Update X and R */
209:     betaold = beta;
210:     VecXDot(Pcurr,R,&beta);                 /*  beta <- pi'*r       */
211:     KSP_MatMult(ksp,Amat,Pcurr,Ccurr);      /*  w <- A*pi (stored in ci)   */
212:     VecXDot(Pcurr,Ccurr,&dpi);              /*  dpi <- pi'*w        */
213:     alphaold = alpha;
214:     alpha = beta / dpi;                                          /*  alpha <- beta/dpi    */
215:     VecAXPY(X,alpha,Pcurr);                 /*  x <- x + alpha * pi  */
216:     VecAXPY(R,-alpha,Ccurr);                /*  r <- r - alpha * wi  */

218:     /* Compute norm for convergence check */
219:     switch (ksp->normtype) {
220:       case KSP_NORM_PRECONDITIONED:
221:         KSP_PCApply(ksp,R,Z);               /*   z <- Br             */
222:         VecNorm(Z,NORM_2,&dp);              /*   dp <- sqrt(z'*z) = sqrt(e'*A'*B'*B*A*e)  */
223:         break;
224:       case KSP_NORM_UNPRECONDITIONED:
225:         VecNorm(R,NORM_2,&dp);              /*   dp <- sqrt(r'*r) = sqrt(e'*A'*A*e)   */
226:         break;
227:       case KSP_NORM_NATURAL:
228:         KSP_PCApply(ksp,R,Z);               /*   z <- Br             */
229:         VecXDot(R,Z,&s);
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:     }

238:     /* Check for convergence */
239:     ksp->rnorm = dp;
240:     KSPLogResidualHistory(ksp,dp);
241:     KSPMonitor(ksp,i+1,dp);
242:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
243:     if (ksp->reason) break;

245:     /* Apply PC if not already done for convergence check */
246:     if(ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->normtype == KSP_NORM_NONE){
247:       KSP_PCApply(ksp,R,Z);               /*   z <- Br         */
248:     }

250:     /* Compute current C (which is W/dpi) */
251:     VecScale(Ccurr,1.0/dpi);              /*   w <- ci/dpi   */

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:   return(0);
267: }

271: PetscErrorCode KSPDestroy_FCG(KSP ksp)
272: {
274:   PetscInt       i;
275:   KSP_FCG        *fcg = (KSP_FCG*)ksp->data;


279:   /* Destroy "standard" work vecs */
280:   VecDestroyVecs(ksp->nwork,&ksp->work);

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:       VecDestroyVecs(fcg->chunksizes[i],&fcg->pPvecs[i]);
286:       VecDestroyVecs(fcg->chunksizes[i],&fcg->pCvecs[i]);
287:     }
288:   }
289:   PetscFree5(fcg->Pvecs,fcg->Cvecs,fcg->pPvecs,fcg->pCvecs,fcg->chunksizes);
290:   /* free space used for singular value calculations */
291:   if (ksp->calc_sings) {
292:     PetscFree4(fcg->e,fcg->d,fcg->ee,fcg->dd);
293:   }
294:   KSPDestroyDefault(ksp);
295:   return(0);
296: }

300: PetscErrorCode KSPView_FCG(KSP ksp,PetscViewer viewer)
301: {
302:   KSP_FCG        *fcg = (KSP_FCG*)ksp->data;
304:   PetscBool      iascii,isstring;
305:   const char     *truncstr;

308:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
309:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);

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");

315:   if (iascii) {
316:     PetscViewerASCIIPrintf(viewer,"  FCG: m_max=%D\n",fcg->mmax);
317:     PetscViewerASCIIPrintf(viewer,"  FCG: preallocated %D directions\n",PetscMin(fcg->nprealloc,fcg->mmax+1));
318:     PetscViewerASCIIPrintf(viewer,"  FCG: %s\n",truncstr);
319:   } else if (isstring) {
320:     PetscViewerStringSPrintf(viewer,"m_max %D nprealloc %D %s",fcg->mmax,fcg->nprealloc,truncstr);
321:   }
322:   return(0);
323: }

327: /*@
328:   KSPFCGSetMmax - set the maximum number of previous directions FCG will store for orthogonalization

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())

334:   Logically Collective on KSP

336:   Input Parameters:
337: +  ksp - the Krylov space context
338: -  mmax - the maximum number of previous directions to orthogonalize againt

340:   Level: intermediate

342:   Options Database:
343: . -ksp_fcg_mmax <N>

345: .seealso: KSPFCG, KSPFCGGetTruncationType(), KSPFCGGetNprealloc()
346: @*/
347: PetscErrorCode KSPFCGSetMmax(KSP ksp,PetscInt mmax)
348: {
349:   KSP_FCG *fcg = (KSP_FCG*)ksp->data;

354:   fcg->mmax = mmax;
355:   return(0);
356: }

360: /*@
361:   KSPFCGGetMmax - get the maximum number of previous directions FCG will store

363:   Note: FCG stores mmax+1 directions at most (mmax previous ones, and one current one)

365:    Not Collective

367:    Input Parameter:
368: .  ksp - the Krylov space context

370:    Output Parameter:
371: .  mmax - the maximum number of previous directons allowed for orthogonalization

373:   Options Database:
374: . -ksp_fcg_mmax <N>

376:    Level: intermediate

378: .keywords: KSP, FCG, truncation

380: .seealso: KSPFCG, KSPFCGGetTruncationType(), KSPFCGGetNprealloc(), KSPFCGSetMmax()
381: @*/

383: PetscErrorCode KSPFCGGetMmax(KSP ksp,PetscInt *mmax)
384: {
385:   KSP_FCG *fcg=(KSP_FCG*)ksp->data;

389:   *mmax = fcg->mmax;
390:   return(0);
391: }

395: /*@
396:   KSPFCGSetNprealloc - set the number of directions to preallocate with FCG

398:   Logically Collective on KSP

400:   Input Parameters:
401: +  ksp - the Krylov space context
402: -  nprealloc - the number of vectors to preallocate

404:   Level: advanced

406:   Options Database:
407: . -ksp_fcg_nprealloc <N> - number of directions to preallocate

409: .seealso: KSPFCG, KSPFCGGetTruncationType(), KSPFCGGetNprealloc()
410: @*/
411: PetscErrorCode KSPFCGSetNprealloc(KSP ksp,PetscInt nprealloc)
412: {
413:   KSP_FCG *fcg=(KSP_FCG*)ksp->data;

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:   return(0);
421: }

425: /*@
426:   KSPFCGGetNprealloc - get the number of directions preallocate by FCG

428:    Not Collective

430:    Input Parameter:
431: .  ksp - the Krylov space context

433:    Output Parameter:
434: .  nprealloc - the number of directions preallocated

436:    Level: advanced

438: .keywords: KSP, FCG, truncation

440: .seealso: KSPFCG, KSPFCGGetTruncationType(), KSPFCGSetNprealloc()
441: @*/
442: PetscErrorCode KSPFCGGetNprealloc(KSP ksp,PetscInt *nprealloc)
443: {
444:   KSP_FCG *fcg=(KSP_FCG*)ksp->data;

448:   *nprealloc = fcg->nprealloc;
449:   return(0);
450: }

454: /*@
455:   KSPFCGSetTruncationType - specify how many of its stored previous directions FCG uses during orthoganalization

457:   Logically Collective on KSP

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,..

462:   Input Parameters:
463: +  ksp - the Krylov space context
464: -  truncstrat - the choice of strategy

466:   Level: intermediate

468:   Options Database:
469: . -ksp_fcg_truncation_type <standard, notay> - specify how many of its stored previous directions FCG uses during orthoganalization

471:   .seealso: KSPFCDTruncationType, KSPFCGGetTruncationType
472: @*/
473: PetscErrorCode KSPFCGSetTruncationType(KSP ksp,KSPFCDTruncationType truncstrat)
474: {
475:   KSP_FCG *fcg=(KSP_FCG*)ksp->data;

480:   fcg->truncstrat=truncstrat;
481:   return(0);
482: }

486: /*@
487:   KSPFCGGetTruncationType - get the truncation strategy employed by FCG

489:    Not Collective

491:    Input Parameter:
492: .  ksp - the Krylov space context

494:    Output Parameter:
495: .  truncstrat - the strategy type

497:    Level: intermediate

499: .keywords: KSP, FCG, truncation

501: .seealso: KSPFCG, KSPFCGSetTruncationType, KSPFCDTruncationType
502: @*/
503: PetscErrorCode KSPFCGGetTruncationType(KSP ksp,KSPFCDTruncationType *truncstrat)
504: {
505:   KSP_FCG *fcg=(KSP_FCG*)ksp->data;

509:   *truncstrat=fcg->truncstrat;
510:   return(0);
511: }

515: PetscErrorCode KSPSetFromOptions_FCG(PetscOptionItems *PetscOptionsObject,KSP ksp)
516: {
518:   KSP_FCG        *fcg=(KSP_FCG*)ksp->data;
519:   PetscInt       mmax,nprealloc;
520:   PetscBool      flg;

523:   PetscOptionsHead(PetscOptionsObject,"KSP FCG Options");
524:   PetscOptionsInt("-ksp_fcg_mmax","Maximum number of search directions to store","KSPFCGSetMmax",fcg->mmax,&mmax,&flg);
525:   if (flg) {
526:     KSPFCGSetMmax(ksp,mmax);
527:   }
528:   PetscOptionsInt("-ksp_fcg_nprealloc","Number of directions to preallocate","KSPFCGSetNprealloc",fcg->nprealloc,&nprealloc,&flg);
529:   if (flg) {
530:     KSPFCGSetNprealloc(ksp,nprealloc);
531:   }
532:   PetscOptionsEnum("-ksp_fcg_truncation_type","Truncation approach for directions","KSPFCGSetTruncationType",KSPFCDTruncationTypes,(PetscEnum)fcg->truncstrat,(PetscEnum*)&fcg->truncstrat,NULL);
533:   PetscOptionsTail();
534:   return(0);
535: }

537: /*MC
538:       KSPFCG - Implements the Flexible Conjugate Gradient method (FCG)

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

545:     Contributed by Patrick Sanan

547:    Notes:
548:    Supports left preconditioning only.

550:    Level: beginner

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

557:  .seealso : KSPGCR, KSPFGMRES, KSPCG, KSPFCGSetMmax(), KSPFCGGetMmax(), KSPFCGSetNprealloc(), KSPFCGGetNprealloc(), KSPFCGSetTruncationType(), KSPFCGGetTruncationType()

559: M*/
562: PETSC_EXTERN PetscErrorCode KSPCreate_FCG(KSP ksp)
563: {
565:   KSP_FCG        *fcg;

568:   PetscNewLog(ksp,&fcg);
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;

581:   ksp->data = (void*)fcg;

583:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
584:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,1);
585:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,1);

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:   return(0);
595: }