Actual source code: pipefcg.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: /*
  2:     Contributed by Patrick Sanan and Sascha M. Schnepp
  3: */

  5: #include <../src/ksp/ksp/impls/fcg/pipefcg/pipefcgimpl.h>       /*I  "petscksp.h"  I*/

  7: #define KSPPIPEFCG_DEFAULT_MMAX 15
  8: #define KSPPIPEFCG_DEFAULT_NPREALLOC 5
  9: #define KSPPIPEFCG_DEFAULT_VECB 5
 10: #define KSPPIPEFCG_DEFAULT_TRUNCSTRAT KSP_FCD_TRUNC_TYPE_NOTAY

 14: static PetscErrorCode KSPAllocateVectors_PIPEFCG(KSP ksp, PetscInt nvecsneeded, PetscInt chunksize)
 15: {
 16:   PetscErrorCode  ierr;
 17:   PetscInt        i;
 18:   KSP_PIPEFCG     *pipefcg;
 19:   PetscInt        nnewvecs, nvecsprev;

 22:   pipefcg = (KSP_PIPEFCG*)ksp->data;

 24:   /* Allocate enough new vectors to add chunksize new vectors, reach nvecsneedtotal, or to reach mmax+1, whichever is smallest */
 25:   if(pipefcg->nvecs < PetscMin(pipefcg->mmax+1,nvecsneeded)){
 26:     nvecsprev = pipefcg->nvecs;
 27:     nnewvecs = PetscMin(PetscMax(nvecsneeded-pipefcg->nvecs,chunksize),pipefcg->mmax+1-pipefcg->nvecs);
 28:     KSPCreateVecs(ksp,nnewvecs,&pipefcg->pQvecs[pipefcg->nchunks],0,NULL);
 29:     PetscLogObjectParents((PetscObject)ksp,nnewvecs,pipefcg->pQvecs[pipefcg->nchunks]);
 30:     KSPCreateVecs(ksp,nnewvecs,&pipefcg->pZETAvecs[pipefcg->nchunks],0,NULL);
 31:     PetscLogObjectParents((PetscObject)ksp,nnewvecs,pipefcg->pZETAvecs[pipefcg->nchunks]);
 32:     KSPCreateVecs(ksp,nnewvecs,&pipefcg->pPvecs[pipefcg->nchunks],0,NULL);
 33:     PetscLogObjectParents((PetscObject)ksp,nnewvecs,pipefcg->pPvecs[pipefcg->nchunks]);
 34:     KSPCreateVecs(ksp,nnewvecs,&pipefcg->pSvecs[pipefcg->nchunks],0,NULL);
 35:     PetscLogObjectParents((PetscObject)ksp,nnewvecs,pipefcg->pSvecs[pipefcg->nchunks]);
 36:     pipefcg->nvecs += nnewvecs;
 37:     for(i=0;i<nnewvecs;++i){
 38:       pipefcg->Qvecs[nvecsprev + i]    = pipefcg->pQvecs[pipefcg->nchunks][i];
 39:       pipefcg->ZETAvecs[nvecsprev + i] = pipefcg->pZETAvecs[pipefcg->nchunks][i];
 40:       pipefcg->Pvecs[nvecsprev + i]    = pipefcg->pPvecs[pipefcg->nchunks][i];
 41:       pipefcg->Svecs[nvecsprev + i]    = pipefcg->pSvecs[pipefcg->nchunks][i];
 42:     }
 43:     pipefcg->chunksizes[pipefcg->nchunks] = nnewvecs;
 44:     ++pipefcg->nchunks;
 45:   }
 46:   return(0);
 47: }

 51: PetscErrorCode    KSPSetUp_PIPEFCG(KSP ksp)
 52: {
 54:   KSP_PIPEFCG    *pipefcg;
 55:   const PetscInt nworkstd = 5;

 58:   pipefcg = (KSP_PIPEFCG*)ksp->data;

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

 63:   /* Allocated space for pointers to additional work vectors
 64:    note that mmax is the number of previous directions, so we add 1 for the current direction,
 65:    and an extra 1 for the prealloc (which might be empty) */
 66:   PetscMalloc4(pipefcg->mmax+1,&(pipefcg->Pvecs),pipefcg->mmax+1,&(pipefcg->pPvecs),pipefcg->mmax+1,&(pipefcg->Svecs),pipefcg->mmax+1,&(pipefcg->pSvecs));
 67:   PetscMalloc4(pipefcg->mmax+1,&(pipefcg->Qvecs),pipefcg->mmax+1,&(pipefcg->pQvecs),pipefcg->mmax+1,&(pipefcg->ZETAvecs),pipefcg->mmax+1,&(pipefcg->pZETAvecs));
 68:   PetscMalloc4(pipefcg->mmax+1,&(pipefcg->Pold),pipefcg->mmax+1,&(pipefcg->Sold),pipefcg->mmax+1,&(pipefcg->Qold),pipefcg->mmax+1,&(pipefcg->ZETAold));
 69:   PetscMalloc1(pipefcg->mmax+1,&(pipefcg->chunksizes));
 70:   PetscMalloc3(pipefcg->mmax+2,&(pipefcg->dots),pipefcg->mmax+1,&(pipefcg->etas),pipefcg->mmax+2,&(pipefcg->redux));

 72:   /* If the requested number of preallocated vectors is greater than mmax reduce nprealloc */
 73:   if(pipefcg->nprealloc > pipefcg->mmax+1){
 74:     PetscInfo2(NULL,"Requested nprealloc=%d is greater than m_max+1=%d. Resetting nprealloc = m_max+1.\n",pipefcg->nprealloc, pipefcg->mmax+1);
 75:   }

 77:   /* Preallocate additional work vectors */
 78:   KSPAllocateVectors_PIPEFCG(ksp,pipefcg->nprealloc,pipefcg->nprealloc);

 80:   PetscLogObjectMemory((PetscObject)ksp,(pipefcg->mmax+1)*4*sizeof(Vec*)+(pipefcg->mmax+1)*4*sizeof(Vec**)+(pipefcg->mmax+1)*4*sizeof(Vec*)+
 81:     (pipefcg->mmax+1)*sizeof(PetscInt)+(pipefcg->mmax+2)*sizeof(Vec*)+(pipefcg->mmax+2)*sizeof(PetscScalar)+(pipefcg->mmax+1)*sizeof(PetscReal));
 82:   return(0);
 83: }

 87: PetscErrorCode KSPSolve_PIPEFCG_cycle(KSP ksp)
 88: {
 90:   PetscInt       i,j,k,idx,kdx,mi;
 91:   KSP_PIPEFCG    *pipefcg;
 92:   PetscScalar    alpha=0.0,gamma,*betas,*dots;
 93:   PetscReal      dp=0.0, delta,*eta,*etas;
 94:   Vec            B,R,Z,X,Qcurr,W,ZETAcurr,M,N,Pcurr,Scurr,*redux;
 95:   Mat            Amat,Pmat;


 99:   /* We have not checked these routines for use with complex numbers. The inner products
100:      are likely not defined correctly for that case */
101: #if (defined(PETSC_USE_COMPLEX) && !defined(PETSC_SKIP_COMPLEX))
102:   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"PIPEFGMRES has not been implemented for use with complex scalars");
103: #endif

105: #define VecXDot(x,y,a)         (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecDot       (x,y,a)   : VecTDot       (x,y,a))
106: #define VecXDotBegin(x,y,a)    (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecDotBegin  (x,y,a)   : VecTDotBegin  (x,y,a))
107: #define VecXDotEnd(x,y,a)      (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecDotEnd    (x,y,a)   : VecTDotEnd    (x,y,a))
108: #define VecMXDot(x,n,y,a)      (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecMDot      (x,n,y,a) : VecMTDot      (x,n,y,a))
109: #define VecMXDotBegin(x,n,y,a) (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecMDotBegin (x,n,y,a) : VecMTDotBegin (x,n,y,a))
110: #define VecMXDotEnd(x,n,y,a)   (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecMDotEnd   (x,n,y,a) : VecMTDotEnd   (x,n,y,a))

112:   pipefcg       = (KSP_PIPEFCG*)ksp->data;
113:   X             = ksp->vec_sol;
114:   B             = ksp->vec_rhs;
115:   R             = ksp->work[0];
116:   Z             = ksp->work[1];
117:   W             = ksp->work[2];
118:   M             = ksp->work[3];
119:   N             = ksp->work[4];

121:   redux = pipefcg->redux;
122:   dots  = pipefcg->dots;
123:   etas  = pipefcg->etas;
124:   betas = dots;        /* dots takes the result of all dot products of which the betas are a subset */

126:   PCGetOperators(ksp->pc,&Amat,&Pmat);

128:   /* Compute cycle initial residual */
129:   KSP_MatMult(ksp,Amat,X,R);
130:   VecAYPX(R,-1.0,B);                   /* r <- b - Ax */
131:   KSP_PCApply(ksp,R,Z);                /* z <- Br     */

133:   Pcurr = pipefcg->Pvecs[0];
134:   Scurr = pipefcg->Svecs[0];
135:   Qcurr = pipefcg->Qvecs[0];
136:   ZETAcurr = pipefcg->ZETAvecs[0];
137:   VecCopy(Z,Pcurr);
138:   KSP_MatMult(ksp,Amat,Pcurr,Scurr);  /* S = Ap     */
139:   VecCopy(Scurr,W);                   /* w = s = Az */

141:   /* Initial state of pipelining intermediates */
142:   redux[0] = R;
143:   redux[1] = W;
144:   VecMXDotBegin(Z,2,redux,dots);
145:   PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)Z)); /* perform asynchronous reduction */
146:   KSP_PCApply(ksp,W,M);            /* m = B(w) */
147:   KSP_MatMult(ksp,Amat,M,N);       /* n = Am   */
148:   VecCopy(M,Qcurr);                /* q = m    */
149:   VecCopy(N,ZETAcurr);             /* zeta = n */
150:   VecMXDotEnd(Z,2,redux,dots);
151:   gamma    = dots[0];
152:   delta    = PetscRealPart(dots[1]);
153:   etas[0]  = delta;
154:   alpha    = gamma/delta;

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

160:     /* Update X, R, Z, W */
161:     VecAXPY(X,+alpha,Pcurr);           /* x <- x + alpha * pi    */
162:     VecAXPY(R,-alpha,Scurr);           /* r <- r - alpha * si    */
163:     VecAXPY(Z,-alpha,Qcurr);           /* z <- z - alpha * qi    */
164:     VecAXPY(W,-alpha,ZETAcurr);        /* w <- w - alpha * zetai */

166:     /* Compute norm for convergence check */
167:     switch (ksp->normtype) {
168:       case KSP_NORM_PRECONDITIONED:
169:         VecNorm(Z,NORM_2,&dp);         /* dp <- sqrt(z'*z) = sqrt(e'*A'*B'*B*A*e) */
170:         break;
171:       case KSP_NORM_UNPRECONDITIONED:
172:         VecNorm(R,NORM_2,&dp);         /* dp <- sqrt(r'*r) = sqrt(e'*A'*A*e)      */
173:         break;
174:       case KSP_NORM_NATURAL:
175:         dp = PetscSqrtReal(PetscAbsScalar(gamma));          /* dp <- sqrt(r'*z) = sqrt(e'*A'*B*A*e)    */
176:         break;
177:       case KSP_NORM_NONE:
178:         dp = 0.0;
179:         break;
180:       default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]);
181:     }

183:     /* Check for convergence */
184:     ksp->rnorm = dp;
185:     KSPLogResidualHistory(ksp,dp);
186:     KSPMonitor(ksp,ksp->its,dp);
187:     (*ksp->converged)(ksp,ksp->its+1,dp,&ksp->reason,ksp->cnvP);
188:     if (ksp->reason) break;

190:     /* Computations of current iteration done */
191:     ++i;

193:     /* If needbe, allocate a new chunk of vectors in P and C */
194:     KSPAllocateVectors_PIPEFCG(ksp,i+1,pipefcg->vecb);

196:     /* Note that we wrap around and start clobbering old vectors */
197:     idx = i % (pipefcg->mmax+1);
198:     Pcurr    = pipefcg->Pvecs[idx];
199:     Scurr    = pipefcg->Svecs[idx];
200:     Qcurr    = pipefcg->Qvecs[idx];
201:     ZETAcurr = pipefcg->ZETAvecs[idx];
202:     eta      = pipefcg->etas+idx;

204:     /* number of old directions to orthogonalize against */
205:     switch(pipefcg->truncstrat){
206:       case KSP_FCD_TRUNC_TYPE_STANDARD:
207:         mi = pipefcg->mmax;
208:         break;
209:       case KSP_FCD_TRUNC_TYPE_NOTAY:
210:         mi = ((i-1) % pipefcg->mmax)+1;
211:         break;
212:       default:
213:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unrecognized Truncation Strategy");
214:     }

216:     /* Pick old p,s,q,zeta in a way suitable for VecMDot */
217:     VecCopy(Z,Pcurr);
218:     for(k=PetscMax(0,i-mi),j=0;k<i;++j,++k){
219:       kdx = k % (pipefcg->mmax+1);
220:       pipefcg->Pold[j]    = pipefcg->Pvecs[kdx];
221:       pipefcg->Sold[j]    = pipefcg->Svecs[kdx];
222:       pipefcg->Qold[j]    = pipefcg->Qvecs[kdx];
223:       pipefcg->ZETAold[j] = pipefcg->ZETAvecs[kdx];
224:       redux[j]            = pipefcg->Svecs[kdx];
225:     }
226:     redux[j]   = R;   /* If the above loop is not executed redux contains only R => all beta_k = 0, only gamma, delta != 0 */
227:     redux[j+1] = W;

229:     VecMXDotBegin(Z,j+2,redux,betas);  /* Start split reductions for beta_k = (z,s_k), gamma = (z,r), delta = (z,w) */
230:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)Z)); /* perform asynchronous reduction */
231:     VecWAXPY(N,-1.0,R,W);              /* m = u + B(w-r): (a) ntmp = w-r              */
232:     KSP_PCApply(ksp,N,M);              /* m = u + B(w-r): (b) mtmp = B(ntmp) = B(w-r) */
233:     VecAXPY(M,1.0,Z);                  /* m = u + B(w-r): (c) m = z + mtmp            */
234:     KSP_MatMult(ksp,Amat,M,N);         /* n = Am                                      */
235:     VecMXDotEnd(Z,j+2,redux,betas);    /* Finish split reductions */
236:     gamma = betas[j];
237:     delta = PetscRealPart(betas[j+1]);

239:     *eta = 0.;
240:     for(k=PetscMax(0,i-mi),j=0;k<i;++j,++k){
241:       kdx = k % (pipefcg->mmax+1);
242:       betas[j] /= -etas[kdx];                               /* betak  /= etak */
243:       *eta -= ((PetscReal)(PetscAbsScalar(betas[j])*PetscAbsScalar(betas[j]))) * etas[kdx];
244:                                                             /* etaitmp = -betaik^2 * etak */
245:     }
246:     *eta += delta;                                          /* etai    = delta -betaik^2 * etak */
247:     if(*eta < 0.) {
248:       pipefcg->norm_breakdown = PETSC_TRUE;
249:       PetscInfo1(ksp,"Restart due to square root breakdown at it = \n",ksp->its);
250:       break;
251:     } else {
252:       alpha= gamma/(*eta);                                  /* alpha = gamma/etai */
253:     }

255:     /* project out stored search directions using classical G-S */
256:     VecCopy(Z,Pcurr);
257:     VecCopy(W,Scurr);
258:     VecCopy(M,Qcurr);
259:     VecCopy(N,ZETAcurr);
260:     VecMAXPY(Pcurr   ,j,betas,pipefcg->Pold);    /* pi    <- ui - sum_k beta_k p_k    */
261:     VecMAXPY(Scurr   ,j,betas,pipefcg->Sold);    /* si    <- wi - sum_k beta_k s_k    */
262:     VecMAXPY(Qcurr   ,j,betas,pipefcg->Qold);    /* qi    <- m  - sum_k beta_k q_k    */
263:     VecMAXPY(ZETAcurr,j,betas,pipefcg->ZETAold); /* zetai <- n  - sum_k beta_k zeta_k */

265:   } while (ksp->its < ksp->max_it);
266:   return(0);
267: }

271: PetscErrorCode KSPSolve_PIPEFCG(KSP ksp)
272: {
274:   KSP_PIPEFCG    *pipefcg;
275:   PetscScalar    gamma;
276:   PetscReal      dp=0.0;
277:   Vec            B,R,Z,X;
278:   Mat            Amat,Pmat;

280: #define VecXDot(x,y,a)         (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecDot       (x,y,a)   : VecTDot       (x,y,a))

283:   pipefcg       = (KSP_PIPEFCG*)ksp->data;
284:   X             = ksp->vec_sol;
285:   B             = ksp->vec_rhs;
286:   R             = ksp->work[0];
287:   Z             = ksp->work[1];

289:   PCGetOperators(ksp->pc,&Amat,&Pmat);

291:   /* Compute initial residual needed for convergence check*/
292:   ksp->its = 0;
293:   if (!ksp->guess_zero) {
294:     KSP_MatMult(ksp,Amat,X,R);
295:     VecAYPX(R,-1.0,B);                 /* r <- b - Ax                             */
296:   } else {
297:     VecCopy(B,R);                      /* r <- b (x is 0)                         */
298:   }
299:   switch (ksp->normtype) {
300:     case KSP_NORM_PRECONDITIONED:
301:       KSP_PCApply(ksp,R,Z);            /* z <- Br                                 */
302:       VecNorm(Z,NORM_2,&dp);           /* dp <- dqrt(z'*z) = sqrt(e'*A'*B'*B*A*e) */
303:       break;
304:     case KSP_NORM_UNPRECONDITIONED:
305:       VecNorm(R,NORM_2,&dp);           /* dp <- sqrt(r'*r) = sqrt(e'*A'*A*e)      */
306:       break;
307:     case KSP_NORM_NATURAL:
308:       KSP_PCApply(ksp,R,Z);            /* z <- Br                                 */
309:       VecXDot(Z,R,&gamma);
310:       dp = PetscSqrtReal(PetscAbsScalar(gamma));            /* dp <- sqrt(r'*z) = sqrt(e'*A'*B*A*e)    */
311:       break;
312:     case KSP_NORM_NONE:
313:       dp = 0.0;
314:       break;
315:     default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]);
316:   }

318:   /* Initial Convergence Check */
319:   KSPLogResidualHistory(ksp,dp);
320:   KSPMonitor(ksp,0,dp);
321:   ksp->rnorm = dp;
322:   if (ksp->normtype == KSP_NORM_NONE) {
323:     KSPConvergedSkip (ksp,0,dp,&ksp->reason,ksp->cnvP);
324:   } else {
325:     (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
326:   }
327:   if (ksp->reason) return(0);

329:   do {
330:     /* A cycle is broken only if a norm breakdown occurs. If not the entire solve happens in a single cycle.
331:        This is coded this way to allow both truncation and truncation-restart strategy
332:        (see KSPFCDGetNumOldDirections()) */
333:     KSPSolve_PIPEFCG_cycle(ksp);
334:     if (ksp->reason) break;
335:     if (pipefcg->norm_breakdown) {
336:       pipefcg->n_restarts++;
337:       pipefcg->norm_breakdown = PETSC_FALSE;
338:     }
339:   } while (ksp->its < ksp->max_it);

341:   if (ksp->its >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
342:   return(0);
343: }

347: PetscErrorCode KSPDestroy_PIPEFCG(KSP ksp)
348: {
350:   PetscInt       i;
351:   KSP_PIPEFCG    *pipefcg;

354:   pipefcg = (KSP_PIPEFCG*)ksp->data;

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

359:   /* Destroy vectors of old directions and the arrays that manage pointers to them */
360:   if(pipefcg->nvecs){
361:     for(i=0;i<pipefcg->nchunks;++i){
362:       VecDestroyVecs(pipefcg->chunksizes[i],&pipefcg->pPvecs[i]);
363:       VecDestroyVecs(pipefcg->chunksizes[i],&pipefcg->pSvecs[i]);
364:       VecDestroyVecs(pipefcg->chunksizes[i],&pipefcg->pQvecs[i]);
365:       VecDestroyVecs(pipefcg->chunksizes[i],&pipefcg->pZETAvecs[i]);
366:     }
367:   }
368:   PetscFree4(pipefcg->Pvecs,pipefcg->Svecs,pipefcg->pPvecs,pipefcg->pSvecs);
369:   PetscFree4(pipefcg->Qvecs,pipefcg->ZETAvecs,pipefcg->pQvecs,pipefcg->pZETAvecs);
370:   PetscFree4(pipefcg->Pold,pipefcg->Sold,pipefcg->Qold,pipefcg->ZETAold);
371:   PetscFree(pipefcg->chunksizes);
372:   PetscFree3(pipefcg->dots,pipefcg->etas,pipefcg->redux);
373:   KSPDestroyDefault(ksp);
374:   return(0);
375: }

379: PetscErrorCode KSPView_PIPEFCG(KSP ksp,PetscViewer viewer)
380: {
381:   KSP_PIPEFCG    *pipefcg = (KSP_PIPEFCG*)ksp->data;
383:   PetscBool      iascii,isstring;
384:   const char     *truncstr;

387:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
388:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);

390:   if(pipefcg->truncstrat == KSP_FCD_TRUNC_TYPE_STANDARD){
391:     truncstr = "Using standard truncation strategy";
392:   } else if(pipefcg->truncstrat == KSP_FCD_TRUNC_TYPE_NOTAY){
393:     truncstr = "Using Notay's truncation strategy";
394:   } else {
395:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Undefined FCD truncation strategy");
396:   }

398:   if (iascii) {
399:     PetscViewerASCIIPrintf(viewer,"  PIPEFCG: max previous directions = %D\n",pipefcg->mmax);
400:     PetscViewerASCIIPrintf(viewer,"  PIPEFCG: preallocated %D directions\n",PetscMin(pipefcg->nprealloc,pipefcg->mmax+1));
401:     PetscViewerASCIIPrintf(viewer,"  PIPEFCG: %s\n",truncstr);
402:     PetscViewerASCIIPrintf(viewer,"  PIPEFCG: restarts performed = %D \n", pipefcg->n_restarts);
403:   } else if (isstring) {
404:     PetscViewerStringSPrintf(viewer,
405:       "max previous directions = %D, preallocated %D directions, %s truncation strategy",
406:       pipefcg->mmax,pipefcg->nprealloc,truncstr);
407:   }
408:   return(0);
409: }

413: /*@
414:   KSPPIPEFCGSetMmax - set the maximum number of previous directions PIPEFCG will store for orthogonalization

416:   Note: mmax + 1 directions are stored (mmax previous ones along with the current one)
417:   and whether all are used in each iteration also depends on the truncation strategy
418:   (see KSPPIPEFCGSetTruncationType)

420:   Logically Collective on KSP

422:   Input Parameters:
423: +  ksp - the Krylov space context
424: -  mmax - the maximum number of previous directions to orthogonalize against

426:   Level: intermediate

428:   Options Database:
429: . -ksp_pipefcg_mmax <N>

431: .seealso: KSPPIPEFCG, KSPPIPEFCGSetTruncationType(), KSPPIPEFCGSetNprealloc()
432: @*/
433: PetscErrorCode KSPPIPEFCGSetMmax(KSP ksp,PetscInt mmax)
434: {
435:   KSP_PIPEFCG *pipefcg=(KSP_PIPEFCG*)ksp->data;;

440:   pipefcg->mmax=mmax;
441:   return(0);
442: }

446: /*@
447:   KSPPIPEFCGGetMmax - get the maximum number of previous directions PIPEFCG will store

449:   Note: PIPEFCG stores mmax+1 directions at most (mmax previous ones, and the current one)

451:    Not Collective

453:    Input Parameter:
454: .  ksp - the Krylov space context

456:    Output Parameter:
457: .  mmax - the maximum number of previous directons allowed for orthogonalization

459:   Options Database:
460: . -ksp_pipefcg_mmax <N>

462:    Level: intermediate

464: .keywords: KSP, PIPEFCG, truncation

466: .seealso: KSPPIPEFCG, KSPPIPEFCGGetTruncationType(), KSPPIPEFCGGetNprealloc(), KSPPIPEFCGSetMmax()
467: @*/
468: PetscErrorCode KSPPIPEFCGGetMmax(KSP ksp,PetscInt *mmax)
469: {
470:   KSP_PIPEFCG *pipefcg=(KSP_PIPEFCG*)ksp->data;;

474:   *mmax=pipefcg->mmax;
475:   return(0);
476: }

480: /*@
481:   KSPPIPEFCGSetNprealloc - set the number of directions to preallocate with PIPEFCG

483:   Logically Collective on KSP

485:   Input Parameters:
486: +  ksp - the Krylov space context
487: -  nprealloc - the number of vectors to preallocate

489:   Level: advanced

491:   Options Database:
492: . -ksp_pipefcg_nprealloc <N>

494: .seealso: KSPPIPEFCG, KSPPIPEFCGSetTruncationType(), KSPPIPEFCGGetNprealloc()
495: @*/
496: PetscErrorCode KSPPIPEFCGSetNprealloc(KSP ksp,PetscInt nprealloc)
497: {
498:   KSP_PIPEFCG *pipefcg=(KSP_PIPEFCG*)ksp->data;;

503:   pipefcg->nprealloc = nprealloc;
504:   return(0);
505: }

509: /*@
510:   KSPPIPEFCGGetNprealloc - get the number of directions to preallocate by PIPEFCG

512:    Not Collective

514:    Input Parameter:
515: .  ksp - the Krylov space context

517:    Output Parameter:
518: .  nprealloc - the number of directions preallocated

520:   Options Database:
521: . -ksp_pipefcg_nprealloc <N>

523:    Level: advanced

525: .keywords: KSP, PIPEFCG, truncation

527: .seealso: KSPPIPEFCG, KSPPIPEFCGGetTruncationType(), KSPPIPEFCGSetNprealloc()
528: @*/
529: PetscErrorCode KSPPIPEFCGGetNprealloc(KSP ksp,PetscInt *nprealloc)
530: {
531:   KSP_PIPEFCG *pipefcg=(KSP_PIPEFCG*)ksp->data;;

535:   *nprealloc = pipefcg->nprealloc;
536:   return(0);
537: }

541: /*@
542:   KSPPIPEFCGSetTruncationType - specify how many of its stored previous directions PIPEFCG uses during orthoganalization

544:   Logically Collective on KSP

546:   KSP_FCD_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions
547:   KSP_FCD_TRUNC_TYPE_NOTAY uses max(1,mod(i,mmax)) stored directions at iteration i=0,1,..

549:   Input Parameters:
550: +  ksp - the Krylov space context
551: -  truncstrat - the choice of strategy

553:   Level: intermediate

555:   Options Database:
556: . -ksp_pipefcg_truncation, -ksp_pipefcg_truncation_restart

558: .seealso: KSPPIPEFCG, KSPPIPEFCGGetTruncationType, KSPFCDTruncationType
559: @*/
560: PetscErrorCode KSPPIPEFCGSetTruncationType(KSP ksp,KSPFCDTruncationType truncstrat)
561: {
562:   KSP_PIPEFCG *pipefcg=(KSP_PIPEFCG*)ksp->data;;

567:   pipefcg->truncstrat=truncstrat;
568:   return(0);
569: }

573: /*@
574:   KSPPIPEFCGGetTruncationType - get the truncation strategy employed by PIPEFCG

576:    Not Collective

578:    Input Parameter:
579: .  ksp - the Krylov space context

581:    Output Parameter:
582: .  truncstrat - the strategy type

584:   Options Database:
585: . -ksp_pipefcg_truncation, -ksp_pipefcg_truncation_restart

587:    Level: intermediate

589: .keywords: KSP, PIPEFCG, truncation

591: .seealso: KSPPIPEFCG, KSPPIPEFCGSetTruncationType, KSPFCDTruncationType
592: @*/
593: PetscErrorCode KSPPIPEFCGGetTruncationType(KSP ksp,KSPFCDTruncationType *truncstrat)
594: {
595:   KSP_PIPEFCG *pipefcg=(KSP_PIPEFCG*)ksp->data;;

599:   *truncstrat=pipefcg->truncstrat;
600:   return(0);
601: }

605: PetscErrorCode KSPSetFromOptions_PIPEFCG(PetscOptionItems *PetscOptionsObject,KSP ksp)
606: {
608:   KSP_PIPEFCG    *pipefcg=(KSP_PIPEFCG*)ksp->data;
609:   PetscInt       mmax,nprealloc;
610:   PetscBool      flg;

613:   PetscOptionsHead(PetscOptionsObject,"KSP PIPEFCG options");
614:   PetscOptionsInt("-ksp_pipefcg_mmax","Number of search directions to storue","KSPPIPEFCGSetMmax",pipefcg->mmax,&mmax,&flg);
615:   if (flg) KSPPIPEFCGSetMmax(ksp,mmax);
616:   PetscOptionsInt("-ksp_pipefcg_nprealloc","Number of directions to preallocate","KSPPIPEFCGSetNprealloc",pipefcg->nprealloc,&nprealloc,&flg);
617:   if (flg) { KSPPIPEFCGSetNprealloc(ksp,nprealloc); }
618:   PetscOptionsEnum("-ksp_pipefcg_truncation_type","Truncation approach for directions","KSPFCGSetTruncationType",KSPFCDTruncationTypes,(PetscEnum)pipefcg->truncstrat,(PetscEnum*)&pipefcg->truncstrat,NULL);
619:   PetscOptionsTail();
620:   return(0);
621: }

623: /*MC

625:  KSPPIPEFCG - A Pipelined, Flexible Conjugate Gradient method

627:     The natural norm for this method is (u,Au). This norm is available at no computational costs. Choosing norm types preconditioned or unpreconditioned involves an extra blocking global reduction, thus removing any benefit from pipelining.

629:  Supports left preconditioning only.
630:   Reference:
631:     Pipelined, Flexible Krylov Subspace Methods
632:     Patrick Sanan, Sascha M. Schnepp, Dave A. May

634:  Options Database Keys:
635: + -ksp_pipefcg_mmax <N>
636: . -ksp_pipefcg_nprealloc <N>
637: . -ksp_pipefcg_truncation
638: - -ksp_pipefcg_trancation_restart

640:   Level: intermediate

642: .seealso : KSPFCG, KSPPIPECG, KSPPIPECR, KSPGCR, KSPPIPEGCR, KSPFGMRES, KSPCG, KSPPIPEFCGSetMmax(), KSPPIPEFCGGetMmax(), KSPPIPEFCGSetNprealloc(), KSPPIPEFCGGetNprealloc(), KSPPIPEFCGSetTruncationType(), KSPPIPEFCGGetTruncationType()

644: M*/
647: PETSC_EXTERN PetscErrorCode KSPCreate_PIPEFCG(KSP ksp)
648: {
650:   KSP_PIPEFCG    *pipefcg;

653:   PetscNewLog(ksp,&pipefcg);
654: #if !defined(PETSC_USE_COMPLEX)
655:   pipefcg->type       = KSP_CG_SYMMETRIC;
656: #else
657:   pipefcg->type       = KSP_CG_HERMITIAN;
658: #endif
659:   pipefcg->mmax       = KSPPIPEFCG_DEFAULT_MMAX;
660:   pipefcg->nprealloc  = KSPPIPEFCG_DEFAULT_NPREALLOC;
661:   pipefcg->nvecs      = 0;
662:   pipefcg->vecb       = KSPPIPEFCG_DEFAULT_VECB;
663:   pipefcg->nchunks    = 0;
664:   pipefcg->truncstrat = KSPPIPEFCG_DEFAULT_TRUNCSTRAT;
665:   pipefcg->n_restarts = 0;

667:   ksp->data = (void*)pipefcg;

669:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
670:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,1);
671:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,1);

673:   ksp->ops->setup          = KSPSetUp_PIPEFCG;
674:   ksp->ops->solve          = KSPSolve_PIPEFCG;
675:   ksp->ops->destroy        = KSPDestroy_PIPEFCG;
676:   ksp->ops->view           = KSPView_PIPEFCG;
677:   ksp->ops->setfromoptions = KSPSetFromOptions_PIPEFCG;
678:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
679:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
680:   return(0);
681: }