Actual source code: pipefgmres.c

petsc-3.7.2 2016-06-05
Report Typos and Errors
  1: /*
  2:   Contributed by Patrick Sanan and Sascha M. Schnepp
  3: */

  5: #include <../src/ksp/ksp/impls/gmres/pipefgmres/pipefgmresimpl.h>       /*I  "petscksp.h"  I*/

  7: #define PIPEFGMRES_DELTA_DIRECTIONS 10
  8: #define PIPEFGMRES_DEFAULT_MAXK     30

 10: static PetscErrorCode KSPPIPEFGMRESGetNewVectors(KSP,PetscInt);
 11: static PetscErrorCode KSPPIPEFGMRESUpdateHessenberg(KSP,PetscInt,PetscBool*,PetscReal*);
 12: static PetscErrorCode KSPPIPEFGMRESBuildSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);
 13: extern PetscErrorCode KSPReset_PIPEFGMRES(KSP);

 15: /*

 17:     KSPSetUp_PIPEFGMRES - Sets up the workspace needed by pipefgmres.

 19:     This is called once, usually automatically by KSPSolve() or KSPSetUp(),
 20:     but can be called directly by KSPSetUp().

 22: */
 25: static PetscErrorCode KSPSetUp_PIPEFGMRES(KSP ksp)
 26: {
 28:   PetscInt       k;
 29:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data;
 30:   const PetscInt max_k = pipefgmres->max_k;

 33:   KSPSetUp_GMRES(ksp);

 35:   PetscMalloc1((VEC_OFFSET+max_k),&pipefgmres->prevecs);
 36:   PetscMalloc1((VEC_OFFSET+max_k),&pipefgmres->prevecs_user_work);
 37:   PetscLogObjectMemory((PetscObject)ksp,(VEC_OFFSET+max_k)*(2*sizeof(void*)));

 39:   KSPCreateVecs(ksp,pipefgmres->vv_allocated,&pipefgmres->prevecs_user_work[0],0,NULL);
 40:   PetscLogObjectParents(ksp,pipefgmres->vv_allocated,pipefgmres->prevecs_user_work[0]);
 41:   for (k=0; k < pipefgmres->vv_allocated; k++) {
 42:     pipefgmres->prevecs[k] = pipefgmres->prevecs_user_work[0][k];
 43:   }

 45:   PetscMalloc1((VEC_OFFSET+max_k),&pipefgmres->zvecs);
 46:   PetscMalloc1((VEC_OFFSET+max_k),&pipefgmres->zvecs_user_work);
 47:   PetscLogObjectMemory((PetscObject)ksp,(VEC_OFFSET+max_k)*(2*sizeof(void*)));

 49:   PetscMalloc1((VEC_OFFSET+max_k),&pipefgmres->redux);
 50:   PetscLogObjectMemory((PetscObject)ksp,(VEC_OFFSET+max_k)*(sizeof(void*)));

 52:   KSPCreateVecs(ksp,pipefgmres->vv_allocated,&pipefgmres->zvecs_user_work[0],0,NULL);
 53:   PetscLogObjectParents(ksp,pipefgmres->vv_allocated,pipefgmres->zvecs_user_work[0]);
 54:   for (k=0; k < pipefgmres->vv_allocated; k++) {
 55:     pipefgmres->zvecs[k] = pipefgmres->zvecs_user_work[0][k];
 56:   }

 58:   return(0);
 59: }

 61: /*

 63:     KSPPIPEFGMRESCycle - Run pipefgmres, possibly with restart.  Return residual
 64:                   history if requested.

 66:     input parameters:
 67: .        pipefgmres  - structure containing parameters and work areas

 69:     output parameters:
 70: .        itcount - number of iterations used.  If null, ignored.
 71: .        converged - 0 if not converged

 73:     Notes:
 74:     On entry, the value in vector VEC_VV(0) should be
 75:     the initial residual.


 78:  */

 82: static PetscErrorCode KSPPIPEFGMRESCycle(PetscInt *itcount,KSP ksp)
 83: {
 84:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)(ksp->data);
 85:   PetscReal      res_norm;
 86:   PetscReal      hapbnd,tt;
 87:   PetscScalar    *hh,*hes,*lhh,shift = pipefgmres->shift;
 88:   PetscBool      hapend = PETSC_FALSE;  /* indicates happy breakdown ending */
 90:   PetscInt       loc_it;                /* local count of # of dir. in Krylov space */
 91:   PetscInt       max_k = pipefgmres->max_k; /* max # of directions Krylov space */
 92:   PetscInt       i,j,k;
 93:   Mat            Amat,Pmat;
 94:   Vec            Q,W; /* Pipelining vectors */
 95:   Vec            *redux = pipefgmres->redux; /* workspace for single reduction */


 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:   if (itcount) *itcount = 0;

107:   /* Assign simpler names to these vectors, allocated as pipelining workspace */
108:   Q = VEC_Q;
109:   W = VEC_W;

111:   /* Allocate memory for orthogonalization work (freed in the GMRES Destroy routine)*/
112:   /* Note that we add an extra value here to allow for a single reduction */
113:   if (!pipefgmres->orthogwork) { PetscMalloc1(pipefgmres->max_k + 2 ,&pipefgmres->orthogwork);
114:   }
115:   lhh = pipefgmres->orthogwork;

117:   /* Number of pseudo iterations since last restart is the number
118:      of prestart directions */
119:   loc_it = 0;

121:   /* note: (pipefgmres->it) is always set one less than (loc_it) It is used in
122:      KSPBUILDSolution_PIPEFGMRES, where it is passed to KSPPIPEFGMRESBuildSoln.
123:      Note that when KSPPIPEFGMRESBuildSoln is called from this function,
124:      (loc_it -1) is passed, so the two are equivalent */
125:   pipefgmres->it = (loc_it - 1);

127:   /* initial residual is in VEC_VV(0)  - compute its norm*/
128:   VecNorm(VEC_VV(0),NORM_2,&res_norm);

130:   /* first entry in right-hand-side of hessenberg system is just
131:      the initial residual norm */
132:   *RS(0) = res_norm;

134:   ksp->rnorm = res_norm;
135:   KSPLogResidualHistory(ksp,res_norm);
136:   KSPMonitor(ksp,ksp->its,res_norm);

138:   /* check for the convergence - maybe the current guess is good enough */
139:   (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);
140:   if (ksp->reason) {
141:     if (itcount) *itcount = 0;
142:     return(0);
143:   }

145:   /* scale VEC_VV (the initial residual) */
146:   VecScale(VEC_VV(0),1.0/res_norm);

148:   /* Fill the pipeline */
149:   KSP_PCApply(ksp,VEC_VV(loc_it),PREVEC(loc_it));
150:   PCGetOperators(ksp->pc,&Amat,&Pmat);
151:   KSP_MatMult(ksp,Amat,PREVEC(loc_it),ZVEC(loc_it));
152:   VecAXPY(ZVEC(loc_it),-shift,VEC_VV(loc_it)); /* Note shift */

154:   /* MAIN ITERATION LOOP BEGINNING*/
155:   /* keep iterating until we have converged OR generated the max number
156:      of directions OR reached the max number of iterations for the method */
157:   while (!ksp->reason && loc_it < max_k && ksp->its < ksp->max_it) {
158:     if (loc_it) {
159:       KSPLogResidualHistory(ksp,res_norm);
160:       KSPMonitor(ksp,ksp->its,res_norm);
161:     }
162:     pipefgmres->it = (loc_it - 1);

164:     /* see if more space is needed for work vectors */
165:     if (pipefgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
166:       KSPPIPEFGMRESGetNewVectors(ksp,loc_it+1);
167:       /* (loc_it+1) is passed in as number of the first vector that should
168:          be allocated */
169:     }

171:     /* Note that these inner products are with "Z" now, so
172:        in particular, lhh[loc_it] is the 'barred' or 'shifted' value,
173:        not the value from the equivalent FGMRES run (even in exact arithmetic)
174:        That is, the H we need for the Arnoldi relation is different from the
175:        coefficients we use in the orthogonalization process,because of the shift */

177:     /* Do some local twiddling to allow for a single reduction */
178:     for(i=0;i<loc_it+1;i++){
179:       redux[i] = VEC_VV(i);
180:     }
181:     redux[loc_it+1] = ZVEC(loc_it);

183:     /* note the extra dot product which ends up in lh[loc_it+1], which computes ||z||^2 */
184:     VecMDotBegin(ZVEC(loc_it),loc_it+2,redux,lhh);

186:     /* Start the split reduction (This actually calls the MPI_Iallreduce, otherwise, the reduction is simply delayed until the "end" call)*/
187:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)ZVEC(loc_it)));

189:     /* The work to be overlapped with the inner products follows.
190:        This is application of the preconditioner and the operator
191:        to compute intermediate quantites which will be combined (locally)
192:        with the results of the inner products.
193:        */
194:     KSP_PCApply(ksp,ZVEC(loc_it),Q);
195:     PCGetOperators(ksp->pc,&Amat,&Pmat);
196:     KSP_MatMult(ksp,Amat,Q,W);

198:     /* Compute inner products of the new direction with previous directions,
199:        and the norm of the to-be-orthogonalized direction "Z".
200:        This information is enough to build the required entries
201:        of H. The inner product with VEC_VV(it_loc) is
202:        *different* than in the standard FGMRES and need to be dealt with specially.
203:        That is, for standard FGMRES the orthogonalization coefficients are the same
204:        as the coefficients used in the Arnoldi relation to reconstruct, but here this
205:        is not true (albeit only for the one entry of H which we "unshift" below. */

207:     /* Finish the dot product, retrieving the extra entry */
208:     VecMDotEnd(ZVEC(loc_it),loc_it+2,redux,lhh);
209:     tt = PetscRealPart(lhh[loc_it+1]);

211:     /* Hessenberg entries, and entries for (naive) classical Graham-Schmidt
212:       Note that the Hessenberg entries require a shift, as these are for the
213:       relation AU = VH, which is wrt unshifted basis vectors */
214:     hh = HH(0,loc_it); hes=HES(0,loc_it);
215:     for (j=0; j<loc_it; j++) {
216:       hh[j]  = lhh[j];
217:       hes[j] = lhh[j];
218:     }
219:     hh[loc_it]  = lhh[loc_it] + shift;
220:     hes[loc_it] = lhh[loc_it] + shift;

222:     /* we delay applying the shift here */
223:     for (j=0; j<=loc_it; j++) {
224:       lhh[j]        = -lhh[j]; /* flip sign */
225:     }

227:     /* Compute the norm of the un-normalized new direction using the rearranged formula
228:        Note that these are shifted ("barred") quantities */
229:     for(k=0;k<=loc_it;k++) tt -= ((PetscReal)(PetscAbsScalar(lhh[k]) * PetscAbsScalar(lhh[k])));
230:     if (tt < 0.0) {
231:       /* If we detect square root breakdown in the norm, we must restart the algorithm.
232:          Here this means we simply break the current loop and reconstruct the solution
233:          using the basis we have computed thus far. Note that by breaking immediately,
234:          we do not update the iteration count, so computation done in this iteration
235:          should be disregarded.
236:          */
237:       PetscInfo1(ksp,"Restart due to square root breakdown at it = \n",ksp->its);
238:       break;
239:     } else {
240:       tt = PetscSqrtReal(tt);
241:     }

243:     /* new entry in hessenburg is the 2-norm of our new direction */
244:     hh[loc_it+1]  = tt;
245:     hes[loc_it+1] = tt;

247:     /* The recurred computation for the new direction
248:        The division by tt is delayed to the happy breakdown check later
249:        Note placement BEFORE the unshift
250:        */
251:     VecCopy(ZVEC(loc_it),VEC_VV(loc_it+1));
252:     VecMAXPY(VEC_VV(loc_it+1),loc_it+1,lhh,&VEC_VV(0));
253:     /* (VEC_VV(loc_it+1) is not normalized yet) */

255:     /* The recurred computation for the preconditioned vector (u) */
256:     VecCopy(Q,PREVEC(loc_it+1));
257:     VecMAXPY(PREVEC(loc_it+1),loc_it+1,lhh,&PREVEC(0));
258:     VecScale(PREVEC(loc_it+1),1.0/tt);

260:     /* Unshift an entry in the GS coefficients ("removing the bar") */
261:     lhh[loc_it]         -= shift;

263:     /* The recurred computation for z (Au)
264:        Note placement AFTER the "unshift" */
265:     VecCopy(W,ZVEC(loc_it+1));
266:     VecMAXPY(ZVEC(loc_it+1),loc_it+1,lhh,&ZVEC(0));
267:     VecScale(ZVEC(loc_it+1),1.0/tt);

269:     /* Happy Breakdown Check */
270:     hapbnd = PetscAbsScalar((tt) / *RS(loc_it));
271:     /* RS(loc_it) contains the res_norm from the last iteration  */
272:     hapbnd = PetscMin(pipefgmres->haptol,hapbnd);
273:     if (tt > hapbnd) {
274:       /* scale new direction by its norm  */
275:       VecScale(VEC_VV(loc_it+1),1.0/tt);
276:     } else {
277:       /* This happens when the solution is exactly reached. */
278:       /* So there is no new direction... */
279:       VecSet(VEC_TEMP,0.0);     /* set VEC_TEMP to 0 */
280:       hapend = PETSC_TRUE;
281:     }
282:     /* note that for pipefgmres we could get HES(loc_it+1, loc_it)  = 0 and the
283:        current solution would not be exact if HES was singular.  Note that
284:        HH non-singular implies that HES is not singular, and HES is guaranteed
285:        to be nonsingular when PREVECS are linearly independent and A is
286:        nonsingular (in GMRES, the nonsingularity of A implies the nonsingularity
287:        of HES). So we should really add a check to verify that HES is nonsingular.*/

289:     /* Note that to be thorough, in debug mode, one could call a LAPACK routine
290:        here to check that the hessenberg matrix is indeed non-singular (since
291:        FGMRES does not guarantee this) */

293:     /* Now apply rotations to new col of hessenberg (and right side of system),
294:        calculate new rotation, and get new residual norm at the same time*/
295:     KSPPIPEFGMRESUpdateHessenberg(ksp,loc_it,&hapend,&res_norm);
296:     if (ksp->reason) break;

298:     loc_it++;
299:     pipefgmres->it = (loc_it-1);   /* Add this here in case it has converged */

301:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
302:     ksp->its++;
303:     ksp->rnorm = res_norm;
304:     PetscObjectSAWsGrantAccess((PetscObject)ksp);

306:     (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);

308:     /* Catch error in happy breakdown and signal convergence and break from loop */
309:     if (hapend) {
310:       if (!ksp->reason) {
311:         if (ksp->errorifnotconverged) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"You reached the happy break down, but convergence was not indicated. Residual norm = %g",(double)res_norm);
312:         else {
313:           ksp->reason = KSP_DIVERGED_BREAKDOWN;
314:           break;
315:         }
316:       }
317:     }
318:   }
319:   /* END OF ITERATION LOOP */
320:   KSPLogResidualHistory(ksp,res_norm);

322:   /*
323:      Monitor if we know that we will not return for a restart */
324:   if (loc_it && (ksp->reason || ksp->its >= ksp->max_it)) {
325:     KSPMonitor(ksp,ksp->its,res_norm);
326:   }

328:   if (itcount) *itcount = loc_it;

330:   /*
331:     Down here we have to solve for the "best" coefficients of the Krylov
332:     columns, add the solution values together, and possibly unwind the
333:     preconditioning from the solution
334:    */

336:   /* Form the solution (or the solution so far) */
337:   /* Note: must pass in (loc_it-1) for iteration count so that KSPPIPEGMRESIIBuildSoln
338:      properly navigates */

340:   KSPPIPEFGMRESBuildSoln(RS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);

342:   return(0);
343: }

345: /*
346:     KSPSolve_PIPEFGMRES - This routine applies the PIPEFGMRES method.


349:    Input Parameter:
350: .     ksp - the Krylov space object that was set to use pipefgmres

352:    Output Parameter:
353: .     outits - number of iterations used

355: */
358: static PetscErrorCode KSPSolve_PIPEFGMRES(KSP ksp)
359: {
361:   PetscInt       its,itcount;
362:   KSP_PIPEFGMRES *pipefgmres    = (KSP_PIPEFGMRES*)ksp->data;
363:   PetscBool      guess_zero = ksp->guess_zero;

366:   if (ksp->calc_sings && !pipefgmres->Rsvd) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
367:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
368:   ksp->its = 0;
369:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

371:   itcount     = 0;
372:   ksp->reason = KSP_CONVERGED_ITERATING;
373:   while (!ksp->reason) {
374:     KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
375:     KSPPIPEFGMRESCycle(&its,ksp);
376:     itcount += its;
377:     if (itcount >= ksp->max_it) {
378:       if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
379:       break;
380:     }
381:     ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
382:   }
383:   ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
384:   return(0);
385: }

389: static PetscErrorCode KSPDestroy_PIPEFGMRES(KSP ksp)
390: {

394:   KSPReset_PIPEFGMRES(ksp);
395:   KSPDestroy_GMRES(ksp);
396:   return(0);
397: }

399: /*
400:     KSPPIPEFGMRESBuildSoln - create the solution from the starting vector and the
401:                       current iterates.

403:     Input parameters:
404:         nrs - work area of size it + 1.
405:         vguess  - index of initial guess
406:         vdest - index of result.  Note that vguess may == vdest (replace
407:                 guess with the solution).
408:         it - HH upper triangular part is a block of size (it+1) x (it+1)

410:      This is an internal routine that knows about the PIPEFGMRES internals.
411:  */
414: static PetscErrorCode KSPPIPEFGMRESBuildSoln(PetscScalar *nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
415: {
416:   PetscScalar    tt;
418:   PetscInt       k,j;
419:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)(ksp->data);

422:   /* Solve for solution vector that minimizes the residual */

424:   if (it < 0) {                                 /* no pipefgmres steps have been performed */
425:     VecCopy(vguess,vdest); /* VecCopy() is smart, exits immediately if vguess == vdest */
426:     return(0);
427:   }

429:   /* solve the upper triangular system - RS is the right side and HH is
430:      the upper triangular matrix  - put soln in nrs */
431:   if (*HH(it,it) != 0.0) nrs[it] = *RS(it) / *HH(it,it);
432:   else nrs[it] = 0.0;

434:   for (k=it-1; k>=0; k--) {
435:     tt = *RS(k);
436:     for (j=k+1; j<=it; j++) tt -= *HH(k,j) * nrs[j];
437:     nrs[k] = tt / *HH(k,k);
438:   }

440:   /* Accumulate the correction to the solution of the preconditioned problem in VEC_TEMP */
441:   VecZeroEntries(VEC_TEMP);
442:   VecMAXPY(VEC_TEMP,it+1,nrs,&PREVEC(0));

444:   /* add solution to previous solution */
445:   if (vdest == vguess) {
446:     VecAXPY(vdest,1.0,VEC_TEMP);
447:   } else {
448:     VecWAXPY(vdest,1.0,VEC_TEMP,vguess);
449:   }
450:   return(0);
451: }

453: /*

455:     KSPPIPEFGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
456:                             Return new residual.

458:     input parameters:

460: .        ksp -    Krylov space object
461: .        it  -    plane rotations are applied to the (it+1)th column of the
462:                   modified hessenberg (i.e. HH(:,it))
463: .        hapend - PETSC_FALSE not happy breakdown ending.

465:     output parameters:
466: .        res - the new residual

468:  */
471: /*
472: .  it - column of the Hessenberg that is complete, PIPEFGMRES is actually computing two columns ahead of this
473:  */
474: static PetscErrorCode KSPPIPEFGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool *hapend,PetscReal *res)
475: {
476:   PetscScalar    *hh,*cc,*ss,*rs;
477:   PetscInt       j;
478:   PetscReal      hapbnd;
479:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)(ksp->data);

483:   hh = HH(0,it);   /* pointer to beginning of column to update */
484:   cc = CC(0);      /* beginning of cosine rotations */
485:   ss = SS(0);      /* beginning of sine rotations */
486:   rs = RS(0);      /* right hand side of least squares system */

488:   /* The Hessenberg matrix is now correct through column it, save that form for possible spectral analysis */
489:   for (j=0; j<=it+1; j++) *HES(j,it) = hh[j];

491:   /* check for the happy breakdown */
492:   hapbnd = PetscMin(PetscAbsScalar(hh[it+1] / rs[it]),pipefgmres->haptol);
493:   if (PetscAbsScalar(hh[it+1]) < hapbnd) {
494:     PetscInfo4(ksp,"Detected happy breakdown, current hapbnd = %14.12e H(%D,%D) = %14.12e\n",(double)hapbnd,it+1,it,(double)PetscAbsScalar(*HH(it+1,it)));
495:     *hapend = PETSC_TRUE;
496:   }

498:   /* Apply all the previously computed plane rotations to the new column
499:      of the Hessenberg matrix */
500:   /* Note: this uses the rotation [conj(c)  s ; -s   c], c= cos(theta), s= sin(theta),
501:      and some refs have [c   s ; -conj(s)  c] (don't be confused!) */

503:   for (j=0; j<it; j++) {
504:     PetscScalar hhj = hh[j];
505:     hh[j]   = PetscConj(cc[j])*hhj + ss[j]*hh[j+1];
506:     hh[j+1] =          -ss[j] *hhj + cc[j]*hh[j+1];
507:   }

509:   /*
510:     compute the new plane rotation, and apply it to:
511:      1) the right-hand-side of the Hessenberg system (RS)
512:         note: it affects RS(it) and RS(it+1)
513:      2) the new column of the Hessenberg matrix
514:         note: it affects HH(it,it) which is currently pointed to
515:         by hh and HH(it+1, it) (*(hh+1))
516:     thus obtaining the updated value of the residual...
517:   */

519:   /* compute new plane rotation */

521:   if (!*hapend) {
522:     PetscReal delta = PetscSqrtReal(PetscSqr(PetscAbsScalar(hh[it])) + PetscSqr(PetscAbsScalar(hh[it+1])));
523:     if (delta == 0.0) {
524:       ksp->reason = KSP_DIVERGED_NULL;
525:       return(0);
526:     }

528:     cc[it] = hh[it] / delta;    /* new cosine value */
529:     ss[it] = hh[it+1] / delta;  /* new sine value */

531:     hh[it]   = PetscConj(cc[it])*hh[it] + ss[it]*hh[it+1];
532:     rs[it+1] = -ss[it]*rs[it];
533:     rs[it]   = PetscConj(cc[it])*rs[it];
534:     *res     = PetscAbsScalar(rs[it+1]);
535:   } else { /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply
536:             another rotation matrix (so RH doesn't change).  The new residual is
537:             always the new sine term times the residual from last time (RS(it)),
538:             but now the new sine rotation would be zero...so the residual should
539:             be zero...so we will multiply "zero" by the last residual.  This might
540:             not be exactly what we want to do here -could just return "zero". */

542:     *res = 0.0;
543:   }
544:   return(0);
545: }

547: /*
548:    KSPBuildSolution_PIPEFGMRES

550:      Input Parameter:
551: .     ksp - the Krylov space object
552: .     ptr-

554:    Output Parameter:
555: .     result - the solution

557:    Note: this calls KSPPIPEFGMRESBuildSoln - the same function that KSPPIPEFGMRESCycle
558:    calls directly.

560: */
563: PetscErrorCode KSPBuildSolution_PIPEFGMRES(KSP ksp,Vec ptr,Vec *result)
564: {
565:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data;

569:   if (!ptr) {
570:     if (!pipefgmres->sol_temp) {
571:       VecDuplicate(ksp->vec_sol,&pipefgmres->sol_temp);
572:       PetscLogObjectParent((PetscObject)ksp,(PetscObject)pipefgmres->sol_temp);
573:     }
574:     ptr = pipefgmres->sol_temp;
575:   }
576:   if (!pipefgmres->nrs) {
577:     /* allocate the work area */
578:     PetscMalloc1(pipefgmres->max_k,&pipefgmres->nrs);
579:     PetscLogObjectMemory((PetscObject)ksp,pipefgmres->max_k*sizeof(PetscScalar));
580:   }

582:   KSPPIPEFGMRESBuildSoln(pipefgmres->nrs,ksp->vec_sol,ptr,ksp,pipefgmres->it);
583:   if (result) *result = ptr;
584:   return(0);
585: }

589: PetscErrorCode KSPSetFromOptions_PIPEFGMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
590: {
592:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data;
593:   PetscBool      flg;
594:   PetscScalar    shift;

597:   KSPSetFromOptions_GMRES(PetscOptionsObject,ksp);
598:   PetscOptionsHead(PetscOptionsObject,"KSP pipelined FGMRES Options");
599:   PetscOptionsScalar("-ksp_pipefgmres_shift","shift parameter","KSPPIPEFGMRESSetShift",pipefgmres->shift,&shift,&flg);
600:   if (flg) { KSPPIPEFGMRESSetShift(ksp,shift); }
601:   PetscOptionsTail();
602:   return(0);
603: }

607: PetscErrorCode KSPView_PIPEFGMRES(KSP ksp,PetscViewer viewer)
608: {
609:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data;
611:   PetscBool      iascii,isstring;

614:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
615:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);

617:   if (iascii) {
618:     PetscViewerASCIIPrintf(viewer,"  GMRES: restart=%D\n",pipefgmres->max_k);
619:     PetscViewerASCIIPrintf(viewer,"  GMRES: happy breakdown tolerance %g\n",(double)pipefgmres->haptol);
620: #if defined(PETSC_USE_COMPLEX)
621:     PetscViewerASCIIPrintf(viewer," PIPEFGMRES: shift=%g+%gi\n",PetscRealPart(pipefgmres->shift),PetscImaginaryPart(pipefgmres->shift));
622: #else
623:     PetscViewerASCIIPrintf(viewer," PIPEFGMRES: shift=%g\n",pipefgmres->shift);
624: #endif
625:   } else if (isstring) {
626:     PetscViewerStringSPrintf(viewer,"restart %D",pipefgmres->max_k);
627: #if defined(PETSC_USE_COMPLEX)
628:     PetscViewerStringSPrintf(viewer,"   PIPEFGMRES: shift=%g+%gi\n",PetscRealPart(pipefgmres->shift),PetscImaginaryPart(pipefgmres->shift));
629: #else
630:     PetscViewerStringSPrintf(viewer,"   PIPEFGMRES: shift=%g\n",pipefgmres->shift);
631: #endif
632:   }
633:   return(0);
634: }

638: PetscErrorCode KSPReset_PIPEFGMRES(KSP ksp)
639: {
640:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data;
641:   PetscErrorCode   ierr;
642:   PetscInt         i;

645:   PetscFree(pipefgmres->prevecs);
646:   PetscFree(pipefgmres->zvecs);
647:   for (i=0; i<pipefgmres->nwork_alloc; i++) {
648:     VecDestroyVecs(pipefgmres->mwork_alloc[i],&pipefgmres->prevecs_user_work[i]);
649:     VecDestroyVecs(pipefgmres->mwork_alloc[i],&pipefgmres->zvecs_user_work[i]);
650:   }
651:   PetscFree(pipefgmres->prevecs_user_work);
652:   PetscFree(pipefgmres->zvecs_user_work);
653:   PetscFree(pipefgmres->redux);
654:   KSPReset_GMRES(ksp);
655:   return(0);
656: }

658: /*MC
659:    KSPPIPEFGMRES - Implements the Pipelined Generalized Minimal Residual method.

661:    A Flexible, 1-stage pipelined variant of GMRES

663:    This variant is not "explicitly normalized" like PGMRES, and requires a shift parameter.

665:    A heuristic for choosing the shift parameter is the largest eigenvalue of the preconditioned operator.

667:    Only right preconditioning is supported (but this preconditioner may be nonlinear, as with FGMRES)

669:    Options Database Keys:
670: +   -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
671: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
672: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
673: .   -ksp_pipefgmres_shift - the shift to use (defaults to 1. See KSPPIPEFGMRESSetShift()
674:                              vectors are allocated as needed)
675: -   -ksp_gmres_krylov_monitor - plot the Krylov space generated


678:    Level: beginner

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

684:    Developer Notes: This object is subclassed off of KSPGMRES

686: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPLGMRES, KSPPIPECG, KSPPIPECR, KSPPGMRES, KSPFGMRES
687:            KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESMonitorKrylov(), KSPPIPEGMRESSetShift()
688: M*/

692: PETSC_EXTERN PetscErrorCode KSPCreate_PIPEFGMRES(KSP ksp)
693: {
694:   KSP_PIPEFGMRES *pipefgmres;

698:   PetscNewLog(ksp,&pipefgmres);

700:   ksp->data                              = (void*)pipefgmres;
701:   ksp->ops->buildsolution                = KSPBuildSolution_PIPEFGMRES;
702:   ksp->ops->setup                        = KSPSetUp_PIPEFGMRES;
703:   ksp->ops->solve                        = KSPSolve_PIPEFGMRES;
704:   ksp->ops->reset                        = KSPReset_PIPEFGMRES;
705:   ksp->ops->destroy                      = KSPDestroy_PIPEFGMRES;
706:   ksp->ops->view                         = KSPView_PIPEFGMRES;
707:   ksp->ops->setfromoptions               = KSPSetFromOptions_PIPEFGMRES;
708:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
709:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

711:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,3);

713:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
714:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);
715:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);

717:   pipefgmres->nextra_vecs    = 1;
718:   pipefgmres->haptol         = 1.0e-30;
719:   pipefgmres->q_preallocate  = 0;
720:   pipefgmres->delta_allocate = PIPEFGMRES_DELTA_DIRECTIONS;
721:   pipefgmres->orthog         = 0;
722:   pipefgmres->nrs            = 0;
723:   pipefgmres->sol_temp       = 0;
724:   pipefgmres->max_k          = PIPEFGMRES_DEFAULT_MAXK;
725:   pipefgmres->Rsvd           = 0;
726:   pipefgmres->orthogwork     = 0;
727:   pipefgmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
728:   pipefgmres->shift          = 1.0;
729:   return(0);
730: }

734: static PetscErrorCode KSPPIPEFGMRESGetNewVectors(KSP ksp,PetscInt it)
735: {
736:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data;
737:   PetscInt       nwork   = pipefgmres->nwork_alloc; /* number of work vector chunks allocated */
738:   PetscInt       nalloc;                            /* number to allocate */
740:   PetscInt       k;

743:   nalloc = pipefgmres->delta_allocate; /* number of vectors to allocate
744:                                       in a single chunk */

746:   /* Adjust the number to allocate to make sure that we don't exceed the
747:      number of available slots (pipefgmres->vecs_allocated)*/
748:   if (it + VEC_OFFSET + nalloc >= pipefgmres->vecs_allocated) {
749:     nalloc = pipefgmres->vecs_allocated - it - VEC_OFFSET;
750:   }
751:   if (!nalloc) return(0);

753:   pipefgmres->vv_allocated += nalloc; /* vv_allocated is the number of vectors allocated */

755:   /* work vectors */
756:   KSPCreateVecs(ksp,nalloc,&pipefgmres->user_work[nwork],0,NULL);
757:   PetscLogObjectParents(ksp,nalloc,pipefgmres->user_work[nwork]);
758:   for (k=0; k < nalloc; k++) {
759:     pipefgmres->vecs[it+VEC_OFFSET+k] = pipefgmres->user_work[nwork][k];
760:   }
761:   /* specify size of chunk allocated */
762:   pipefgmres->mwork_alloc[nwork] = nalloc;

764:   /* preconditioned vectors (note we don't use VEC_OFFSET) */
765:   KSPCreateVecs(ksp,nalloc,&pipefgmres->prevecs_user_work[nwork],0,NULL);
766:   PetscLogObjectParents(ksp,nalloc,pipefgmres->prevecs_user_work[nwork]);
767:   for (k=0; k < nalloc; k++) {
768:     pipefgmres->prevecs[it+k] = pipefgmres->prevecs_user_work[nwork][k];
769:   }

771:   KSPCreateVecs(ksp,nalloc,&pipefgmres->zvecs_user_work[nwork],0,NULL);
772:   PetscLogObjectParents(ksp,nalloc,pipefgmres->zvecs_user_work[nwork]);
773:   for (k=0; k < nalloc; k++) {
774:     pipefgmres->zvecs[it+k] = pipefgmres->zvecs_user_work[nwork][k];
775:   }

777:   /* increment the number of work vector chunks */
778:   pipefgmres->nwork_alloc++;
779:   return(0);
780: }
783: /*@
784:   KSPPIPEFGMRESSetShift - Set the shift parameter for the flexible, pipelined GMRES solver.

786:   A heuristic is to set this to be comparable to the largest eigenvalue of the preconditioned operator. This can be acheived with PETSc itself by using a few iterations of a Krylov method. See KSPComputeEigenvalues (and note the caveats there).

788: Logically Collective on KSP

790: Input Parameters:
791: +  ksp - the Krylov space context
792: -  shift - the shift

794: Level: intermediate

796: Options Database:
797: . -ksp_pipefgmres_shift <shift>

799: .seealso: KSPComputeEigenvalues()
800: @*/
801: PetscErrorCode KSPPIPEFGMRESSetShift(KSP ksp,PetscScalar shift)
802: {
803:   KSP_PIPEFGMRES *pipefgmres = (KSP_PIPEFGMRES*)ksp->data;

808:   pipefgmres->shift = shift;
809:   return(0);
810: }