Actual source code: pipefcg.c

petsc-3.10.5 2019-03-28
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>

  7: static PetscBool  cited = PETSC_FALSE;
  8: static const char citation[] =
  9:   "@article{SSM2016,\n"
 10:   "  author = {P. Sanan and S.M. Schnepp and D.A. May},\n"
 11:   "  title = {Pipelined, Flexible Krylov Subspace Methods},\n"
 12:   "  journal = {SIAM Journal on Scientific Computing},\n"
 13:   "  volume = {38},\n"
 14:   "  number = {5},\n"
 15:   "  pages = {C441-C470},\n"
 16:   "  year = {2016},\n"
 17:   "  doi = {10.1137/15M1049130},\n"
 18:   "  URL = {http://dx.doi.org/10.1137/15M1049130},\n"
 19:   "  eprint = {http://dx.doi.org/10.1137/15M1049130}\n"
 20:   "}\n";

 22: #define KSPPIPEFCG_DEFAULT_MMAX 15
 23: #define KSPPIPEFCG_DEFAULT_NPREALLOC 5
 24: #define KSPPIPEFCG_DEFAULT_VECB 5
 25: #define KSPPIPEFCG_DEFAULT_TRUNCSTRAT KSP_FCD_TRUNC_TYPE_NOTAY

 27: static PetscErrorCode KSPAllocateVectors_PIPEFCG(KSP ksp, PetscInt nvecsneeded, PetscInt chunksize)
 28: {
 29:   PetscErrorCode  ierr;
 30:   PetscInt        i;
 31:   KSP_PIPEFCG     *pipefcg;
 32:   PetscInt        nnewvecs, nvecsprev;

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

 37:   /* Allocate enough new vectors to add chunksize new vectors, reach nvecsneedtotal, or to reach mmax+1, whichever is smallest */
 38:   if(pipefcg->nvecs < PetscMin(pipefcg->mmax+1,nvecsneeded)){
 39:     nvecsprev = pipefcg->nvecs;
 40:     nnewvecs = PetscMin(PetscMax(nvecsneeded-pipefcg->nvecs,chunksize),pipefcg->mmax+1-pipefcg->nvecs);
 41:     KSPCreateVecs(ksp,nnewvecs,&pipefcg->pQvecs[pipefcg->nchunks],0,NULL);
 42:     PetscLogObjectParents((PetscObject)ksp,nnewvecs,pipefcg->pQvecs[pipefcg->nchunks]);
 43:     KSPCreateVecs(ksp,nnewvecs,&pipefcg->pZETAvecs[pipefcg->nchunks],0,NULL);
 44:     PetscLogObjectParents((PetscObject)ksp,nnewvecs,pipefcg->pZETAvecs[pipefcg->nchunks]);
 45:     KSPCreateVecs(ksp,nnewvecs,&pipefcg->pPvecs[pipefcg->nchunks],0,NULL);
 46:     PetscLogObjectParents((PetscObject)ksp,nnewvecs,pipefcg->pPvecs[pipefcg->nchunks]);
 47:     KSPCreateVecs(ksp,nnewvecs,&pipefcg->pSvecs[pipefcg->nchunks],0,NULL);
 48:     PetscLogObjectParents((PetscObject)ksp,nnewvecs,pipefcg->pSvecs[pipefcg->nchunks]);
 49:     pipefcg->nvecs += nnewvecs;
 50:     for(i=0;i<nnewvecs;++i){
 51:       pipefcg->Qvecs[nvecsprev + i]    = pipefcg->pQvecs[pipefcg->nchunks][i];
 52:       pipefcg->ZETAvecs[nvecsprev + i] = pipefcg->pZETAvecs[pipefcg->nchunks][i];
 53:       pipefcg->Pvecs[nvecsprev + i]    = pipefcg->pPvecs[pipefcg->nchunks][i];
 54:       pipefcg->Svecs[nvecsprev + i]    = pipefcg->pSvecs[pipefcg->nchunks][i];
 55:     }
 56:     pipefcg->chunksizes[pipefcg->nchunks] = nnewvecs;
 57:     ++pipefcg->nchunks;
 58:   }
 59:   return(0);
 60: }

 62: static PetscErrorCode    KSPSetUp_PIPEFCG(KSP ksp)
 63: {
 65:   KSP_PIPEFCG    *pipefcg;
 66:   const PetscInt nworkstd = 5;

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

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

 74:   /* Allocated space for pointers to additional work vectors
 75:    note that mmax is the number of previous directions, so we add 1 for the current direction,
 76:    and an extra 1 for the prealloc (which might be empty) */
 77:   PetscMalloc4(pipefcg->mmax+1,&(pipefcg->Pvecs),pipefcg->mmax+1,&(pipefcg->pPvecs),pipefcg->mmax+1,&(pipefcg->Svecs),pipefcg->mmax+1,&(pipefcg->pSvecs));
 78:   PetscMalloc4(pipefcg->mmax+1,&(pipefcg->Qvecs),pipefcg->mmax+1,&(pipefcg->pQvecs),pipefcg->mmax+1,&(pipefcg->ZETAvecs),pipefcg->mmax+1,&(pipefcg->pZETAvecs));
 79:   PetscMalloc4(pipefcg->mmax+1,&(pipefcg->Pold),pipefcg->mmax+1,&(pipefcg->Sold),pipefcg->mmax+1,&(pipefcg->Qold),pipefcg->mmax+1,&(pipefcg->ZETAold));
 80:   PetscMalloc1(pipefcg->mmax+1,&(pipefcg->chunksizes));
 81:   PetscMalloc3(pipefcg->mmax+2,&(pipefcg->dots),pipefcg->mmax+1,&(pipefcg->etas),pipefcg->mmax+2,&(pipefcg->redux));

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

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

 91:   PetscLogObjectMemory((PetscObject)ksp,(pipefcg->mmax+1)*4*sizeof(Vec*)+(pipefcg->mmax+1)*4*sizeof(Vec**)+(pipefcg->mmax+1)*4*sizeof(Vec*)+
 92:     (pipefcg->mmax+1)*sizeof(PetscInt)+(pipefcg->mmax+2)*sizeof(Vec*)+(pipefcg->mmax+2)*sizeof(PetscScalar)+(pipefcg->mmax+1)*sizeof(PetscReal));
 93:   return(0);
 94: }

 96: static PetscErrorCode KSPSolve_PIPEFCG_cycle(KSP ksp)
 97: {
 99:   PetscInt       i,j,k,idx,kdx,mi;
100:   KSP_PIPEFCG    *pipefcg;
101:   PetscScalar    alpha=0.0,gamma,*betas,*dots;
102:   PetscReal      dp=0.0, delta,*eta,*etas;
103:   Vec            B,R,Z,X,Qcurr,W,ZETAcurr,M,N,Pcurr,Scurr,*redux;
104:   Mat            Amat,Pmat;


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

114: #define VecXDot(x,y,a)         (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecDot       (x,y,a)   : VecTDot       (x,y,a))
115: #define VecXDotBegin(x,y,a)    (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecDotBegin  (x,y,a)   : VecTDotBegin  (x,y,a))
116: #define VecXDotEnd(x,y,a)      (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecDotEnd    (x,y,a)   : VecTDotEnd    (x,y,a))
117: #define VecMXDot(x,n,y,a)      (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecMDot      (x,n,y,a) : VecMTDot      (x,n,y,a))
118: #define VecMXDotBegin(x,n,y,a) (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecMDotBegin (x,n,y,a) : VecMTDotBegin (x,n,y,a))
119: #define VecMXDotEnd(x,n,y,a)   (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecMDotEnd   (x,n,y,a) : VecMTDotEnd   (x,n,y,a))

121:   pipefcg       = (KSP_PIPEFCG*)ksp->data;
122:   X             = ksp->vec_sol;
123:   B             = ksp->vec_rhs;
124:   R             = ksp->work[0];
125:   Z             = ksp->work[1];
126:   W             = ksp->work[2];
127:   M             = ksp->work[3];
128:   N             = ksp->work[4];

130:   redux = pipefcg->redux;
131:   dots  = pipefcg->dots;
132:   etas  = pipefcg->etas;
133:   betas = dots;        /* dots takes the result of all dot products of which the betas are a subset */

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

137:   /* Compute cycle initial residual */
138:   KSP_MatMult(ksp,Amat,X,R);
139:   VecAYPX(R,-1.0,B);                   /* r <- b - Ax */
140:   KSP_PCApply(ksp,R,Z);                /* z <- Br     */

142:   Pcurr = pipefcg->Pvecs[0];
143:   Scurr = pipefcg->Svecs[0];
144:   Qcurr = pipefcg->Qvecs[0];
145:   ZETAcurr = pipefcg->ZETAvecs[0];
146:   VecCopy(Z,Pcurr);
147:   KSP_MatMult(ksp,Amat,Pcurr,Scurr);  /* S = Ap     */
148:   VecCopy(Scurr,W);                   /* w = s = Az */

150:   /* Initial state of pipelining intermediates */
151:   redux[0] = R;
152:   redux[1] = W;
153:   VecMXDotBegin(Z,2,redux,dots);
154:   PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)Z)); /* perform asynchronous reduction */
155:   KSP_PCApply(ksp,W,M);            /* m = B(w) */
156:   KSP_MatMult(ksp,Amat,M,N);       /* n = Am   */
157:   VecCopy(M,Qcurr);                /* q = m    */
158:   VecCopy(N,ZETAcurr);             /* zeta = n */
159:   VecMXDotEnd(Z,2,redux,dots);
160:   gamma    = dots[0];
161:   delta    = PetscRealPart(dots[1]);
162:   etas[0]  = delta;
163:   alpha    = gamma/delta;

165:   i = 0;
166:   do {
167:     ksp->its++;

169:     /* Update X, R, Z, W */
170:     VecAXPY(X,+alpha,Pcurr);           /* x <- x + alpha * pi    */
171:     VecAXPY(R,-alpha,Scurr);           /* r <- r - alpha * si    */
172:     VecAXPY(Z,-alpha,Qcurr);           /* z <- z - alpha * qi    */
173:     VecAXPY(W,-alpha,ZETAcurr);        /* w <- w - alpha * zetai */

175:     /* Compute norm for convergence check */
176:     switch (ksp->normtype) {
177:       case KSP_NORM_PRECONDITIONED:
178:         VecNorm(Z,NORM_2,&dp);         /* dp <- sqrt(z'*z) = sqrt(e'*A'*B'*B*A*e) */
179:         break;
180:       case KSP_NORM_UNPRECONDITIONED:
181:         VecNorm(R,NORM_2,&dp);         /* dp <- sqrt(r'*r) = sqrt(e'*A'*A*e)      */
182:         break;
183:       case KSP_NORM_NATURAL:
184:         dp = PetscSqrtReal(PetscAbsScalar(gamma));          /* dp <- sqrt(r'*z) = sqrt(e'*A'*B*A*e)    */
185:         break;
186:       case KSP_NORM_NONE:
187:         dp = 0.0;
188:         break;
189:       default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]);
190:     }

192:     /* Check for convergence */
193:     ksp->rnorm = dp;
194:     KSPLogResidualHistory(ksp,dp);
195:     KSPMonitor(ksp,ksp->its,dp);
196:     (*ksp->converged)(ksp,ksp->its+1,dp,&ksp->reason,ksp->cnvP);
197:     if (ksp->reason) break;

199:     /* Computations of current iteration done */
200:     ++i;

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

205:     /* Note that we wrap around and start clobbering old vectors */
206:     idx = i % (pipefcg->mmax+1);
207:     Pcurr    = pipefcg->Pvecs[idx];
208:     Scurr    = pipefcg->Svecs[idx];
209:     Qcurr    = pipefcg->Qvecs[idx];
210:     ZETAcurr = pipefcg->ZETAvecs[idx];
211:     eta      = pipefcg->etas+idx;

213:     /* number of old directions to orthogonalize against */
214:     switch(pipefcg->truncstrat){
215:       case KSP_FCD_TRUNC_TYPE_STANDARD:
216:         mi = pipefcg->mmax;
217:         break;
218:       case KSP_FCD_TRUNC_TYPE_NOTAY:
219:         mi = ((i-1) % pipefcg->mmax)+1;
220:         break;
221:       default:
222:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unrecognized Truncation Strategy");
223:     }

225:     /* Pick old p,s,q,zeta in a way suitable for VecMDot */
226:     VecCopy(Z,Pcurr);
227:     for(k=PetscMax(0,i-mi),j=0;k<i;++j,++k){
228:       kdx = k % (pipefcg->mmax+1);
229:       pipefcg->Pold[j]    = pipefcg->Pvecs[kdx];
230:       pipefcg->Sold[j]    = pipefcg->Svecs[kdx];
231:       pipefcg->Qold[j]    = pipefcg->Qvecs[kdx];
232:       pipefcg->ZETAold[j] = pipefcg->ZETAvecs[kdx];
233:       redux[j]            = pipefcg->Svecs[kdx];
234:     }
235:     redux[j]   = R;   /* If the above loop is not executed redux contains only R => all beta_k = 0, only gamma, delta != 0 */
236:     redux[j+1] = W;

238:     VecMXDotBegin(Z,j+2,redux,betas);  /* Start split reductions for beta_k = (z,s_k), gamma = (z,r), delta = (z,w) */
239:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)Z)); /* perform asynchronous reduction */
240:     VecWAXPY(N,-1.0,R,W);              /* m = u + B(w-r): (a) ntmp = w-r              */
241:     KSP_PCApply(ksp,N,M);              /* m = u + B(w-r): (b) mtmp = B(ntmp) = B(w-r) */
242:     VecAXPY(M,1.0,Z);                  /* m = u + B(w-r): (c) m = z + mtmp            */
243:     KSP_MatMult(ksp,Amat,M,N);         /* n = Am                                      */
244:     VecMXDotEnd(Z,j+2,redux,betas);    /* Finish split reductions */
245:     gamma = betas[j];
246:     delta = PetscRealPart(betas[j+1]);

248:     *eta = 0.;
249:     for(k=PetscMax(0,i-mi),j=0;k<i;++j,++k){
250:       kdx = k % (pipefcg->mmax+1);
251:       betas[j] /= -etas[kdx];                               /* betak  /= etak */
252:       *eta -= ((PetscReal)(PetscAbsScalar(betas[j])*PetscAbsScalar(betas[j]))) * etas[kdx];
253:                                                             /* etaitmp = -betaik^2 * etak */
254:     }
255:     *eta += delta;                                          /* etai    = delta -betaik^2 * etak */
256:     if(*eta < 0.) {
257:       pipefcg->norm_breakdown = PETSC_TRUE;
258:       PetscInfo1(ksp,"Restart due to square root breakdown at it = \n",ksp->its);
259:       break;
260:     } else {
261:       alpha= gamma/(*eta);                                  /* alpha = gamma/etai */
262:     }

264:     /* project out stored search directions using classical G-S */
265:     VecCopy(Z,Pcurr);
266:     VecCopy(W,Scurr);
267:     VecCopy(M,Qcurr);
268:     VecCopy(N,ZETAcurr);
269:     VecMAXPY(Pcurr   ,j,betas,pipefcg->Pold);    /* pi    <- ui - sum_k beta_k p_k    */
270:     VecMAXPY(Scurr   ,j,betas,pipefcg->Sold);    /* si    <- wi - sum_k beta_k s_k    */
271:     VecMAXPY(Qcurr   ,j,betas,pipefcg->Qold);    /* qi    <- m  - sum_k beta_k q_k    */
272:     VecMAXPY(ZETAcurr,j,betas,pipefcg->ZETAold); /* zetai <- n  - sum_k beta_k zeta_k */

274:   } while (ksp->its < ksp->max_it);
275:   return(0);
276: }

278: static PetscErrorCode KSPSolve_PIPEFCG(KSP ksp)
279: {
281:   KSP_PIPEFCG    *pipefcg;
282:   PetscScalar    gamma;
283:   PetscReal      dp=0.0;
284:   Vec            B,R,Z,X;
285:   Mat            Amat,Pmat;

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

290:   PetscCitationsRegister(citation,&cited);

292:   pipefcg       = (KSP_PIPEFCG*)ksp->data;
293:   X             = ksp->vec_sol;
294:   B             = ksp->vec_rhs;
295:   R             = ksp->work[0];
296:   Z             = ksp->work[1];

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

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

327:   /* Initial Convergence Check */
328:   KSPLogResidualHistory(ksp,dp);
329:   KSPMonitor(ksp,0,dp);
330:   ksp->rnorm = dp;
331:   if (ksp->normtype == KSP_NORM_NONE) {
332:     KSPConvergedSkip (ksp,0,dp,&ksp->reason,ksp->cnvP);
333:   } else {
334:     (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
335:   }
336:   if (ksp->reason) return(0);

338:   do {
339:     /* A cycle is broken only if a norm breakdown occurs. If not the entire solve happens in a single cycle.
340:        This is coded this way to allow both truncation and truncation-restart strategy
341:        (see KSPFCDGetNumOldDirections()) */
342:     KSPSolve_PIPEFCG_cycle(ksp);
343:     if (ksp->reason) break;
344:     if (pipefcg->norm_breakdown) {
345:       pipefcg->n_restarts++;
346:       pipefcg->norm_breakdown = PETSC_FALSE;
347:     }
348:   } while (ksp->its < ksp->max_it);

350:   if (ksp->its >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
351:   return(0);
352: }

354: static PetscErrorCode KSPDestroy_PIPEFCG(KSP ksp)
355: {
357:   PetscInt       i;
358:   KSP_PIPEFCG    *pipefcg;

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

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

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

384: static PetscErrorCode KSPView_PIPEFCG(KSP ksp,PetscViewer viewer)
385: {
386:   KSP_PIPEFCG    *pipefcg = (KSP_PIPEFCG*)ksp->data;
388:   PetscBool      iascii,isstring;
389:   const char     *truncstr;

392:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
393:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);

395:   if(pipefcg->truncstrat == KSP_FCD_TRUNC_TYPE_STANDARD){
396:     truncstr = "Using standard truncation strategy";
397:   } else if(pipefcg->truncstrat == KSP_FCD_TRUNC_TYPE_NOTAY){
398:     truncstr = "Using Notay's truncation strategy";
399:   } else {
400:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Undefined FCD truncation strategy");
401:   }

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

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

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

423:   Logically Collective on KSP

425:   Input Parameters:
426: +  ksp - the Krylov space context
427: -  mmax - the maximum number of previous directions to orthogonalize against

429:   Level: intermediate

431:   Options Database:
432: . -ksp_pipefcg_mmax <N>

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

443:   pipefcg->mmax=mmax;
444:   return(0);
445: }

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

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

452:    Not Collective

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

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

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

463:    Level: intermediate

465: .keywords: KSP, PIPEFCG, truncation

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

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

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

482:   Logically Collective on KSP

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

488:   Level: advanced

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

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

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

506: /*@
507:   KSPPIPEFCGGetNprealloc - get the number of directions to preallocate by PIPEFCG

509:    Not Collective

511:    Input Parameter:
512: .  ksp - the Krylov space context

514:    Output Parameter:
515: .  nprealloc - the number of directions preallocated

517:   Options Database:
518: . -ksp_pipefcg_nprealloc <N>

520:    Level: advanced

522: .keywords: KSP, PIPEFCG, truncation

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

532:   *nprealloc = pipefcg->nprealloc;
533:   return(0);
534: }

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

539:   Logically Collective on KSP

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

544:   Input Parameters:
545: +  ksp - the Krylov space context
546: -  truncstrat - the choice of strategy

548:   Level: intermediate

550:   Options Database:
551: .  -ksp_pipefcg_truncation_type <standard,notay> - which stored search directions to orthogonalize against

553: .seealso: KSPPIPEFCG, KSPPIPEFCGGetTruncationType, KSPFCDTruncationType
554: @*/
555: PetscErrorCode KSPPIPEFCGSetTruncationType(KSP ksp,KSPFCDTruncationType truncstrat)
556: {
557:   KSP_PIPEFCG *pipefcg=(KSP_PIPEFCG*)ksp->data;;

562:   pipefcg->truncstrat=truncstrat;
563:   return(0);
564: }

566: /*@
567:   KSPPIPEFCGGetTruncationType - get the truncation strategy employed by PIPEFCG

569:    Not Collective

571:    Input Parameter:
572: .  ksp - the Krylov space context

574:    Output Parameter:
575: .  truncstrat - the strategy type

577:   Options Database:
578: . -ksp_pipefcg_truncation_type <standard,notay> - which stored basis vectors to orthogonalize against

580:    Level: intermediate

582: .keywords: KSP, PIPEFCG, truncation

584: .seealso: KSPPIPEFCG, KSPPIPEFCGSetTruncationType, KSPFCDTruncationType
585: @*/
586: PetscErrorCode KSPPIPEFCGGetTruncationType(KSP ksp,KSPFCDTruncationType *truncstrat)
587: {
588:   KSP_PIPEFCG *pipefcg=(KSP_PIPEFCG*)ksp->data;;

592:   *truncstrat=pipefcg->truncstrat;
593:   return(0);
594: }

596: static PetscErrorCode KSPSetFromOptions_PIPEFCG(PetscOptionItems *PetscOptionsObject,KSP ksp)
597: {
599:   KSP_PIPEFCG    *pipefcg=(KSP_PIPEFCG*)ksp->data;
600:   PetscInt       mmax,nprealloc;
601:   PetscBool      flg;

604:   PetscOptionsHead(PetscOptionsObject,"KSP PIPEFCG options");
605:   PetscOptionsInt("-ksp_pipefcg_mmax","Number of search directions to storue","KSPPIPEFCGSetMmax",pipefcg->mmax,&mmax,&flg);
606:   if (flg) KSPPIPEFCGSetMmax(ksp,mmax);
607:   PetscOptionsInt("-ksp_pipefcg_nprealloc","Number of directions to preallocate","KSPPIPEFCGSetNprealloc",pipefcg->nprealloc,&nprealloc,&flg);
608:   if (flg) { KSPPIPEFCGSetNprealloc(ksp,nprealloc); }
609:   PetscOptionsEnum("-ksp_pipefcg_truncation_type","Truncation approach for directions","KSPFCGSetTruncationType",KSPFCDTruncationTypes,(PetscEnum)pipefcg->truncstrat,(PetscEnum*)&pipefcg->truncstrat,NULL);
610:   PetscOptionsTail();
611:   return(0);
612: }

614: /*MC

616:   KSPPIPEFCG - Implements a Pipelined, Flexible Conjugate Gradient method.

618:   Options Database Keys:
619: .   -ksp_pipefcg_mmax <N> - The number of previous search directions to store
620: .   -ksp_pipefcg_nprealloc <N> - The number of previous search directions to preallocate
621: .   -ksp_pipefcg_truncation_type <standard,notay> - which stored search directions to orthogonalize against

623:   Notes:
624:    Supports left preconditioning only.

626:    The natural "norm" for this method is (u,Au), where u is the preconditioned residual. As with standard CG, this norm is available at no additional computational cost. Choosing preconditioned or unpreconditioned norms involve an extra blocking global reduction, thus removing any benefit from pipelining.

628:    MPI configuration may be necessary for reductions to make asynchronous progress, which is important for performance of pipelined methods.
629:    See the FAQ on the PETSc website for details.

631:   Reference:
632:     P. Sanan, S.M. Schnepp, and D.A. May,
633:     "Pipelined, Flexible Krylov Subspace Methods,"
634:     SIAM Journal on Scientific Computing 2016 38:5, C441-C470,
635:     DOI: 10.1137/15M1049130

637:   Level: intermediate

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

641: M*/
642: PETSC_EXTERN PetscErrorCode KSPCreate_PIPEFCG(KSP ksp)
643: {
645:   KSP_PIPEFCG    *pipefcg;

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

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

664:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
665:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,1);
666:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,1);
667:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

669:   ksp->ops->setup          = KSPSetUp_PIPEFCG;
670:   ksp->ops->solve          = KSPSolve_PIPEFCG;
671:   ksp->ops->destroy        = KSPDestroy_PIPEFCG;
672:   ksp->ops->view           = KSPView_PIPEFCG;
673:   ksp->ops->setfromoptions = KSPSetFromOptions_PIPEFCG;
674:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
675:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
676:   return(0);
677: }