Actual source code: lgmres.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: #include <../src/ksp/ksp/impls/gmres/lgmres/lgmresimpl.h>   /*I petscksp.h I*/

  4: #define LGMRES_DELTA_DIRECTIONS 10
  5: #define LGMRES_DEFAULT_MAXK     30
  6: #define LGMRES_DEFAULT_AUGDIM   2 /*default number of augmentation vectors */
  7: static PetscErrorCode    KSPLGMRESGetNewVectors(KSP,PetscInt);
  8: static PetscErrorCode    KSPLGMRESUpdateHessenberg(KSP,PetscInt,PetscBool,PetscReal*);
  9: static PetscErrorCode    KSPLGMRESBuildSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);

 13: PetscErrorCode  KSPLGMRESSetAugDim(KSP ksp, PetscInt dim)
 14: {

 18:   PetscTryMethod((ksp),"KSPLGMRESSetAugDim_C",(KSP,PetscInt),(ksp,dim));
 19:   return(0);
 20: }

 24: PetscErrorCode  KSPLGMRESSetConstant(KSP ksp)
 25: {

 29:   PetscTryMethod((ksp),"KSPLGMRESSetConstant_C",(KSP),(ksp));
 30:   return(0);
 31: }

 33: /*
 34:     KSPSetUp_LGMRES - Sets up the workspace needed by lgmres.

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

 39: */
 42: PetscErrorCode    KSPSetUp_LGMRES(KSP ksp)
 43: {
 45:   PetscInt       max_k,k, aug_dim;
 46:   KSP_LGMRES     *lgmres = (KSP_LGMRES*)ksp->data;

 49:   max_k   = lgmres->max_k;
 50:   aug_dim = lgmres->aug_dim;
 51:   KSPSetUp_GMRES(ksp);

 53:   /* need array of pointers to augvecs*/
 54:   PetscMalloc1(2*aug_dim + AUG_OFFSET,&lgmres->augvecs);

 56:   lgmres->aug_vecs_allocated = 2 *aug_dim + AUG_OFFSET;

 58:   PetscMalloc1(2*aug_dim + AUG_OFFSET,&lgmres->augvecs_user_work);
 59:   PetscMalloc1(aug_dim,&lgmres->aug_order);
 60:   PetscLogObjectMemory((PetscObject)ksp,(aug_dim)*(4*sizeof(void*) + sizeof(PetscInt)) + AUG_OFFSET*2*sizeof(void*));

 62:   /*  for now we will preallocate the augvecs - because aug_dim << restart
 63:      ... also keep in mind that we need to keep augvecs from cycle to cycle*/
 64:   lgmres->aug_vv_allocated = 2* aug_dim + AUG_OFFSET;
 65:   lgmres->augwork_alloc    =  2* aug_dim + AUG_OFFSET;

 67:   KSPCreateVecs(ksp,lgmres->aug_vv_allocated,&lgmres->augvecs_user_work[0],0,NULL);
 68:   PetscMalloc1(max_k+1,&lgmres->hwork);
 69:   PetscLogObjectParents(ksp,lgmres->aug_vv_allocated,lgmres->augvecs_user_work[0]);
 70:   for (k=0; k<lgmres->aug_vv_allocated; k++) {
 71:     lgmres->augvecs[k] = lgmres->augvecs_user_work[0][k];
 72:   }
 73:   return(0);
 74: }


 77: /*

 79:     KSPLGMRESCycle - Run lgmres, possibly with restart.  Return residual
 80:                   history if requested.

 82:     input parameters:
 83: .        lgmres  - structure containing parameters and work areas

 85:     output parameters:
 86: .        nres    - residuals (from preconditioned system) at each step.
 87:                   If restarting, consider passing nres+it.  If null,
 88:                   ignored
 89: .        itcount - number of iterations used.   nres[0] to nres[itcount]
 90:                   are defined.  If null, ignored.  If null, ignored.
 91: .        converged - 0 if not converged


 94:     Notes:
 95:     On entry, the value in vector VEC_VV(0) should be
 96:     the initial residual.


 99:  */
102: PetscErrorCode KSPLGMRESCycle(PetscInt *itcount,KSP ksp)
103: {
104:   KSP_LGMRES     *lgmres = (KSP_LGMRES*)(ksp->data);
105:   PetscReal      res_norm, res;
106:   PetscReal      hapbnd, tt;
107:   PetscScalar    tmp;
108:   PetscBool      hapend = PETSC_FALSE;  /* indicates happy breakdown ending */
110:   PetscInt       loc_it;                /* local count of # of dir. in Krylov space */
111:   PetscInt       max_k  = lgmres->max_k; /* max approx space size */
112:   PetscInt       max_it = ksp->max_it;  /* max # of overall iterations for the method */

114:   /* LGMRES_MOD - new variables*/
115:   PetscInt    aug_dim = lgmres->aug_dim;
116:   PetscInt    spot    = 0;
117:   PetscInt    order   = 0;
118:   PetscInt    it_arnoldi;                /* number of arnoldi steps to take */
119:   PetscInt    it_total;                  /* total number of its to take (=approx space size)*/
120:   PetscInt    ii, jj;
121:   PetscReal   tmp_norm;
122:   PetscScalar inv_tmp_norm;
123:   PetscScalar *avec;

126:   /* Number of pseudo iterations since last restart is the number
127:      of prestart directions */
128:   loc_it = 0;

130:   /* LGMRES_MOD: determine number of arnoldi steps to take */
131:   /* if approx_constant then we keep the space the same size even if
132:      we don't have the full number of aug vectors yet*/
133:   if (lgmres->approx_constant) it_arnoldi = max_k - lgmres->aug_ct;
134:   else it_arnoldi = max_k - aug_dim;

136:   it_total =  it_arnoldi + lgmres->aug_ct;

138:   /* initial residual is in VEC_VV(0)  - compute its norm*/
139:   VecNorm(VEC_VV(0),NORM_2,&res_norm);
140:   KSPCheckNorm(ksp,res_norm);
141:   res  = res_norm;

143:   /* first entry in right-hand-side of hessenberg system is just
144:      the initial residual norm */
145:   *GRS(0) = res_norm;

147:   /* check for the convergence */
148:   if (!res) {
149:     if (itcount) *itcount = 0;
150:     ksp->reason = KSP_CONVERGED_ATOL;
151:     PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
152:     return(0);
153:   }

155:   /* scale VEC_VV (the initial residual) */
156:   tmp = 1.0/res_norm; VecScale(VEC_VV(0),tmp);

158:   ksp->rnorm = res;


161:   /* note: (lgmres->it) is always set one less than (loc_it) It is used in
162:      KSPBUILDSolution_LGMRES, where it is passed to KSPLGMRESBuildSoln.
163:      Note that when KSPLGMRESBuildSoln is called from this function,
164:      (loc_it -1) is passed, so the two are equivalent */
165:   lgmres->it = (loc_it - 1);


168:   /* MAIN ITERATION LOOP BEGINNING*/


171:   /* keep iterating until we have converged OR generated the max number
172:      of directions OR reached the max number of iterations for the method */
173:   (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);

175:   while (!ksp->reason && loc_it < it_total && ksp->its < max_it) { /* LGMRES_MOD: changed to it_total */
176:     KSPLogResidualHistory(ksp,res);
177:     lgmres->it = (loc_it - 1);
178:     KSPMonitor(ksp,ksp->its,res);

180:     /* see if more space is needed for work vectors */
181:     if (lgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
182:       KSPLGMRESGetNewVectors(ksp,loc_it+1);
183:       /* (loc_it+1) is passed in as number of the first vector that should
184:           be allocated */
185:     }

187:     /*LGMRES_MOD: decide whether this is an arnoldi step or an aug step */
188:     if (loc_it < it_arnoldi) { /* Arnoldi */
189:       KSP_PCApplyBAorAB(ksp,VEC_VV(loc_it),VEC_VV(1+loc_it),VEC_TEMP_MATOP);
190:     } else { /*aug step */
191:       order = loc_it - it_arnoldi + 1; /* which aug step */
192:       for (ii=0; ii<aug_dim; ii++) {
193:         if (lgmres->aug_order[ii] == order) {
194:           spot = ii;
195:           break; /* must have this because there will be duplicates before aug_ct = aug_dim */
196:         }
197:       }

199:       VecCopy(A_AUGVEC(spot), VEC_VV(1+loc_it));
200:       /*note: an alternate implementation choice would be to only save the AUGVECS and
201:         not A_AUGVEC and then apply the PC here to the augvec */
202:     }

204:     /* update hessenberg matrix and do Gram-Schmidt - new direction is in
205:        VEC_VV(1+loc_it)*/
206:     (*lgmres->orthog)(ksp,loc_it);

208:     /* new entry in hessenburg is the 2-norm of our new direction */
209:     VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);

211:     *HH(loc_it+1,loc_it)  = tt;
212:     *HES(loc_it+1,loc_it) = tt;


215:     /* check for the happy breakdown */
216:     hapbnd = PetscAbsScalar(tt / *GRS(loc_it)); /* GRS(loc_it) contains the res_norm from the last iteration  */
217:     if (hapbnd > lgmres->haptol) hapbnd = lgmres->haptol;
218:     if (tt > hapbnd) {
219:       tmp  = 1.0/tt;
220:       VecScale(VEC_VV(loc_it+1),tmp); /* scale new direction by its norm */
221:     } else {
222:       PetscInfo2(ksp,"Detected happy breakdown, current hapbnd = %g tt = %g\n",(double)hapbnd,(double)tt);
223:       hapend = PETSC_TRUE;
224:     }

226:     /* Now apply rotations to new col of hessenberg (and right side of system),
227:        calculate new rotation, and get new residual norm at the same time*/
228:     KSPLGMRESUpdateHessenberg(ksp,loc_it,hapend,&res);
229:     if (ksp->reason) break;

231:     loc_it++;
232:     lgmres->it = (loc_it-1);   /* Add this here in case it has converged */

234:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
235:     ksp->its++;
236:     ksp->rnorm = res;
237:     PetscObjectSAWsGrantAccess((PetscObject)ksp);

239:     (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);

241:     /* Catch error in happy breakdown and signal convergence and break from loop */
242:     if (hapend) {
243:       if (!ksp->reason) {
244:         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);
245:         else {
246:           ksp->reason = KSP_DIVERGED_BREAKDOWN;
247:           break;
248:         }
249:       }
250:     }
251:   }
252:   /* END OF ITERATION LOOP */
253:   KSPLogResidualHistory(ksp,res);

255:   /* Monitor if we know that we will not return for a restart */
256:   if (ksp->reason || ksp->its >= max_it) {
257:     KSPMonitor(ksp, ksp->its, res);
258:   }

260:   if (itcount) *itcount = loc_it;

262:   /*
263:     Down here we have to solve for the "best" coefficients of the Krylov
264:     columns, add the solution values together, and possibly unwind the
265:     preconditioning from the solution
266:    */

268:   /* Form the solution (or the solution so far) */
269:   /* Note: must pass in (loc_it-1) for iteration count so that KSPLGMRESBuildSoln
270:      properly navigates */

272:   KSPLGMRESBuildSoln(GRS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);


275:   /* LGMRES_MOD collect aug vector and A*augvector for future restarts -
276:      only if we will be restarting (i.e. this cycle performed it_total
277:      iterations)  */
278:   if (!ksp->reason && ksp->its < max_it && aug_dim > 0) {

280:     /*AUG_TEMP contains the new augmentation vector (assigned in  KSPLGMRESBuildSoln) */
281:     if (!lgmres->aug_ct) {
282:       spot = 0;
283:       lgmres->aug_ct++;
284:     } else if (lgmres->aug_ct < aug_dim) {
285:       spot = lgmres->aug_ct;
286:       lgmres->aug_ct++;
287:     } else { /* truncate */
288:       for (ii=0; ii<aug_dim; ii++) {
289:         if (lgmres->aug_order[ii] == aug_dim) spot = ii;
290:       }
291:     }



295:     VecCopy(AUG_TEMP, AUGVEC(spot));
296:     /*need to normalize */
297:     VecNorm(AUGVEC(spot), NORM_2, &tmp_norm);

299:     inv_tmp_norm = 1.0/tmp_norm;

301:     VecScale(AUGVEC(spot),inv_tmp_norm);

303:     /*set new aug vector to order 1  - move all others back one */
304:     for (ii=0; ii < aug_dim; ii++) AUG_ORDER(ii)++;
305:     AUG_ORDER(spot) = 1;

307:     /*now add the A*aug vector to A_AUGVEC(spot)  - this is independ. of preconditioning type*/
308:     /* want V*H*y - y is in GRS, V is in VEC_VV and H is in HES */


311:     /* first do H+*y */
312:     avec = lgmres->hwork;
313:     PetscMemzero(avec,(it_total+1)*sizeof(*avec));
314:     for (ii=0; ii < it_total + 1; ii++) {
315:       for (jj=0; jj <= ii+1 && jj < it_total+1; jj++) {
316:         avec[jj] += *HES(jj ,ii) * *GRS(ii);
317:       }
318:     }

320:     /*now multiply result by V+ */
321:     VecSet(VEC_TEMP,0.0);
322:     VecMAXPY(VEC_TEMP, it_total+1, avec, &VEC_VV(0)); /*answer is in VEC_TEMP*/

324:     /*copy answer to aug location  and scale*/
325:     VecCopy(VEC_TEMP,  A_AUGVEC(spot));
326:     VecScale(A_AUGVEC(spot),inv_tmp_norm);
327:   }
328:   return(0);
329: }

331: /*
332:     KSPSolve_LGMRES - This routine applies the LGMRES method.


335:    Input Parameter:
336: .     ksp - the Krylov space object that was set to use lgmres

338:    Output Parameter:
339: .     outits - number of iterations used

341: */

345: PetscErrorCode KSPSolve_LGMRES(KSP ksp)
346: {
348:   PetscInt       cycle_its; /* iterations done in a call to KSPLGMRESCycle */
349:   PetscInt       itcount;   /* running total of iterations, incl. those in restarts */
350:   KSP_LGMRES     *lgmres    = (KSP_LGMRES*)ksp->data;
351:   PetscBool      guess_zero = ksp->guess_zero;
352:   PetscInt       ii;        /*LGMRES_MOD variable */

355:   if (ksp->calc_sings && !lgmres->Rsvd) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");

357:   PetscObjectSAWsTakeAccess((PetscObject)ksp);

359:   ksp->its        = 0;
360:   lgmres->aug_ct  = 0;
361:   lgmres->matvecs = 0;

363:   PetscObjectSAWsGrantAccess((PetscObject)ksp);

365:   /* initialize */
366:   itcount     = 0;
367:   ksp->reason = KSP_CONVERGED_ITERATING;
368:   /*LGMRES_MOD*/
369:   for (ii=0; ii<lgmres->aug_dim; ii++) lgmres->aug_order[ii] = 0;

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

386: /*

388:    KSPDestroy_LGMRES - Frees all memory space used by the Krylov method.

390: */
393: PetscErrorCode KSPDestroy_LGMRES(KSP ksp)
394: {
395:   KSP_LGMRES     *lgmres = (KSP_LGMRES*)ksp->data;

399:   PetscFree(lgmres->augvecs);
400:   if (lgmres->augwork_alloc) {
401:     VecDestroyVecs(lgmres->augwork_alloc,&lgmres->augvecs_user_work[0]);
402:   }
403:   PetscFree(lgmres->augvecs_user_work);
404:   PetscFree(lgmres->aug_order);
405:   PetscFree(lgmres->hwork);
406:   KSPDestroy_GMRES(ksp);
407:   return(0);
408: }

410: /*
411:     KSPLGMRESBuildSoln - create the solution from the starting vector and the
412:                       current iterates.

414:     Input parameters:
415:         nrs - work area of size it + 1.
416:         vguess  - index of initial guess
417:         vdest - index of result.  Note that vguess may == vdest (replace
418:                 guess with the solution).
419:         it - HH upper triangular part is a block of size (it+1) x (it+1)

421:      This is an internal routine that knows about the LGMRES internals.
422:  */
425: static PetscErrorCode KSPLGMRESBuildSoln(PetscScalar *nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
426: {
427:   PetscScalar    tt;
429:   PetscInt       ii,k,j;
430:   KSP_LGMRES     *lgmres = (KSP_LGMRES*)(ksp->data);
431:   /*LGMRES_MOD */
432:   PetscInt it_arnoldi, it_aug;
433:   PetscInt jj, spot = 0;

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

438:   /* If it is < 0, no lgmres steps have been performed */
439:   if (it < 0) {
440:     VecCopy(vguess,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
441:     return(0);
442:   }

444:   /* so (it+1) lgmres steps HAVE been performed */

446:   /* LGMRES_MOD - determine if we need to use augvecs for the soln  - do not assume that
447:      this is called after the total its allowed for an approx space */
448:   if (lgmres->approx_constant) {
449:     it_arnoldi = lgmres->max_k - lgmres->aug_ct;
450:   } else {
451:     it_arnoldi = lgmres->max_k - lgmres->aug_dim;
452:   }
453:   if (it_arnoldi >= it +1) {
454:     it_aug     = 0;
455:     it_arnoldi = it+1;
456:   } else {
457:     it_aug = (it + 1) - it_arnoldi;
458:   }

460:   /* now it_arnoldi indicates the number of matvecs that took place */
461:   lgmres->matvecs += it_arnoldi;


464:   /* solve the upper triangular system - GRS is the right side and HH is
465:      the upper triangular matrix  - put soln in nrs */
466:   if (*HH(it,it) == 0.0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"HH(it,it) is identically zero; it = %D GRS(it) = %g",it,(double)PetscAbsScalar(*GRS(it)));
467:   if (*HH(it,it) != 0.0) {
468:     nrs[it] = *GRS(it) / *HH(it,it);
469:   } else {
470:     nrs[it] = 0.0;
471:   }

473:   for (ii=1; ii<=it; ii++) {
474:     k  = it - ii;
475:     tt = *GRS(k);
476:     for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
477:     nrs[k] = tt / *HH(k,k);
478:   }

480:   /* Accumulate the correction to the soln of the preconditioned prob. in VEC_TEMP */
481:   VecSet(VEC_TEMP,0.0); /* set VEC_TEMP components to 0 */

483:   /*LGMRES_MOD - if augmenting has happened we need to form the solution
484:     using the augvecs */
485:   if (!it_aug) { /* all its are from arnoldi */
486:     VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
487:   } else { /*use aug vecs */
488:     /*first do regular krylov directions */
489:     VecMAXPY(VEC_TEMP,it_arnoldi,nrs,&VEC_VV(0));
490:     /*now add augmented portions - add contribution of aug vectors one at a time*/


493:     for (ii=0; ii<it_aug; ii++) {
494:       for (jj=0; jj<lgmres->aug_dim; jj++) {
495:         if (lgmres->aug_order[jj] == (ii+1)) {
496:           spot = jj;
497:           break; /* must have this because there will be duplicates before aug_ct = aug_dim */
498:         }
499:       }
500:       VecAXPY(VEC_TEMP,nrs[it_arnoldi+ii],AUGVEC(spot));
501:     }
502:   }
503:   /* now VEC_TEMP is what we want to keep for augmenting purposes - grab before the
504:      preconditioner is "unwound" from right-precondtioning*/
505:   VecCopy(VEC_TEMP, AUG_TEMP);

507:   KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);

509:   /* add solution to previous solution */
510:   /* put updated solution into vdest.*/
511:   VecCopy(vguess,vdest);
512:   VecAXPY(vdest,1.0,VEC_TEMP);
513:   return(0);
514: }

516: /*

518:     KSPLGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
519:                             Return new residual.

521:     input parameters:

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

528:     output parameters:
529: .        res - the new residual

531:  */
534: static PetscErrorCode KSPLGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool hapend,PetscReal *res)
535: {
536:   PetscScalar *hh,*cc,*ss,tt;
537:   PetscInt    j;
538:   KSP_LGMRES  *lgmres = (KSP_LGMRES*)(ksp->data);

541:   hh = HH(0,it);   /* pointer to beginning of column to update - so
542:                       incrementing hh "steps down" the (it+1)th col of HH*/
543:   cc = CC(0);      /* beginning of cosine rotations */
544:   ss = SS(0);      /* beginning of sine rotations */

546:   /* Apply all the previously computed plane rotations to the new column
547:      of the Hessenberg matrix */
548:   /* Note: this uses the rotation [conj(c)  s ; -s   c], c= cos(theta), s= sin(theta) */

550:   for (j=1; j<=it; j++) {
551:     tt  = *hh;
552:     *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
553:     hh++;
554:     *hh = *cc++ * *hh - (*ss++ * tt);
555:     /* hh, cc, and ss have all been incremented one by end of loop */
556:   }

558:   /*
559:     compute the new plane rotation, and apply it to:
560:      1) the right-hand-side of the Hessenberg system (GRS)
561:         note: it affects GRS(it) and GRS(it+1)
562:      2) the new column of the Hessenberg matrix
563:         note: it affects HH(it,it) which is currently pointed to
564:         by hh and HH(it+1, it) (*(hh+1))
565:     thus obtaining the updated value of the residual...
566:   */

568:   /* compute new plane rotation */

570:   if (!hapend) {
571:     tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
572:     if (tt == 0.0) {
573:       ksp->reason = KSP_DIVERGED_NULL;
574:       return(0);
575:     }
576:     *cc = *hh / tt;         /* new cosine value */
577:     *ss = *(hh+1) / tt;        /* new sine value */

579:     /* apply to 1) and 2) */
580:     *GRS(it+1) = -(*ss * *GRS(it));
581:     *GRS(it)   = PetscConj(*cc) * *GRS(it);
582:     *hh        = PetscConj(*cc) * *hh + *ss * *(hh+1);

584:     /* residual is the last element (it+1) of right-hand side! */
585:     *res = PetscAbsScalar(*GRS(it+1));

587:   } else { /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
588:             another rotation matrix (so RH doesn't change).  The new residual is
589:             always the new sine term times the residual from last time (GRS(it)),
590:             but now the new sine rotation would be zero...so the residual should
591:             be zero...so we will multiply "zero" by the last residual.  This might
592:             not be exactly what we want to do here -could just return "zero". */

594:     *res = 0.0;
595:   }
596:   return(0);
597: }

599: /*

601:    KSPLGMRESGetNewVectors - This routine allocates more work vectors, starting from
602:                          VEC_VV(it)

604: */
607: static PetscErrorCode KSPLGMRESGetNewVectors(KSP ksp,PetscInt it)
608: {
609:   KSP_LGMRES     *lgmres = (KSP_LGMRES*)ksp->data;
610:   PetscInt       nwork   = lgmres->nwork_alloc; /* number of work vector chunks allocated */
611:   PetscInt       nalloc;                      /* number to allocate */
613:   PetscInt       k;

616:   nalloc = lgmres->delta_allocate; /* number of vectors to allocate
617:                                       in a single chunk */

619:   /* Adjust the number to allocate to make sure that we don't exceed the
620:      number of available slots (lgmres->vecs_allocated)*/
621:   if (it + VEC_OFFSET + nalloc >= lgmres->vecs_allocated) {
622:     nalloc = lgmres->vecs_allocated - it - VEC_OFFSET;
623:   }
624:   if (!nalloc) return(0);

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

628:   /* work vectors */
629:   KSPCreateVecs(ksp,nalloc,&lgmres->user_work[nwork],0,NULL);
630:   PetscLogObjectParents(ksp,nalloc,lgmres->user_work[nwork]);
631:   /* specify size of chunk allocated */
632:   lgmres->mwork_alloc[nwork] = nalloc;

634:   for (k=0; k < nalloc; k++) {
635:     lgmres->vecs[it+VEC_OFFSET+k] = lgmres->user_work[nwork][k];
636:   }


639:   /* LGMRES_MOD - for now we are preallocating the augmentation vectors */


642:   /* increment the number of work vector chunks */
643:   lgmres->nwork_alloc++;
644:   return(0);
645: }

647: /*

649:    KSPBuildSolution_LGMRES

651:      Input Parameter:
652: .     ksp - the Krylov space object
653: .     ptr-

655:    Output Parameter:
656: .     result - the solution

658:    Note: this calls KSPLGMRESBuildSoln - the same function that KSPLGMRESCycle
659:    calls directly.

661: */
664: PetscErrorCode KSPBuildSolution_LGMRES(KSP ksp,Vec ptr,Vec *result)
665: {
666:   KSP_LGMRES     *lgmres = (KSP_LGMRES*)ksp->data;

670:   if (!ptr) {
671:     if (!lgmres->sol_temp) {
672:       VecDuplicate(ksp->vec_sol,&lgmres->sol_temp);
673:       PetscLogObjectParent((PetscObject)ksp,(PetscObject)lgmres->sol_temp);
674:     }
675:     ptr = lgmres->sol_temp;
676:   }
677:   if (!lgmres->nrs) {
678:     /* allocate the work area */
679:     PetscMalloc1(lgmres->max_k,&lgmres->nrs);
680:     PetscLogObjectMemory((PetscObject)ksp,lgmres->max_k*sizeof(PetscScalar));
681:   }

683:   KSPLGMRESBuildSoln(lgmres->nrs,ksp->vec_sol,ptr,ksp,lgmres->it);
684:   if (result) *result = ptr;
685:   return(0);
686: }

690: PetscErrorCode KSPView_LGMRES(KSP ksp,PetscViewer viewer)
691: {
692:   KSP_LGMRES     *lgmres = (KSP_LGMRES*)ksp->data;
694:   PetscBool      iascii;

697:   KSPView_GMRES(ksp,viewer);
698:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
699:   if (iascii) {
700:     /*LGMRES_MOD */
701:     PetscViewerASCIIPrintf(viewer,"  LGMRES: aug. dimension=%D\n",lgmres->aug_dim);
702:     if (lgmres->approx_constant) {
703:       PetscViewerASCIIPrintf(viewer,"  LGMRES: approx. space size was kept constant.\n");
704:     }
705:     PetscViewerASCIIPrintf(viewer,"  LGMRES: number of matvecs=%D\n",lgmres->matvecs);
706:   }
707:   return(0);
708: }

712: PetscErrorCode KSPSetFromOptions_LGMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
713: {
715:   PetscInt       aug;
716:   KSP_LGMRES     *lgmres = (KSP_LGMRES*) ksp->data;
717:   PetscBool      flg     = PETSC_FALSE;

720:   KSPSetFromOptions_GMRES(PetscOptionsObject,ksp);
721:   PetscOptionsHead(PetscOptionsObject,"KSP LGMRES Options");
722:   PetscOptionsBool("-ksp_lgmres_constant","Use constant approx. space size","KSPGMRESSetConstant",lgmres->approx_constant,&lgmres->approx_constant,NULL);
723:   PetscOptionsInt("-ksp_lgmres_augment","Number of error approximations to augment the Krylov space with","KSPLGMRESSetAugDim",lgmres->aug_dim,&aug,&flg);
724:   if (flg) { KSPLGMRESSetAugDim(ksp,aug); }
725:   PetscOptionsTail();
726:   return(0);
727: }

729: /*functions for extra lgmres options here*/
732: static PetscErrorCode  KSPLGMRESSetConstant_LGMRES(KSP ksp)
733: {
734:   KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;

737:   lgmres->approx_constant = PETSC_TRUE;
738:   return(0);
739: }

743: static PetscErrorCode  KSPLGMRESSetAugDim_LGMRES(KSP ksp,PetscInt aug_dim)
744: {
745:   KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;

748:   if (aug_dim < 0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be positive");
749:   if (aug_dim > (lgmres->max_k -1)) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be <= (restart size-1)");
750:   lgmres->aug_dim = aug_dim;
751:   return(0);
752: }

754: /* end new lgmres functions */

756: /*MC
757:     KSPLGMRES - Augments the standard GMRES approximation space with approximations to
758:                 the error from previous restart cycles.

760:   Options Database Keys:
761: +   -ksp_gmres_restart <restart> - total approximation space size (Krylov directions + error approximations)
762: .   -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
763: .   -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
764:                             vectors are allocated as needed)
765: .   -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
766: .   -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
767: .   -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
768:                                   stability of the classical Gram-Schmidt  orthogonalization.
769: .   -ksp_gmres_krylov_monitor - plot the Krylov space generated
770: .   -ksp_lgmres_augment <k> - number of error approximations to augment the Krylov space with
771: -   -ksp_lgmres_constant - use a constant approx. space size (only affects restart cycles < num. error approx.(k), i.e. the first k restarts)

773:     To run LGMRES(m, k) as described in the above paper, use:
774:        -ksp_gmres_restart <m+k>
775:        -ksp_lgmres_augment <k>

777:   Level: beginner

779:    Notes: Supports both left and right preconditioning, but not symmetric.

781:    References:
782: .    1. - A. H. Baker, E.R. Jessup, and T.A. Manteuffel. A technique for accelerating the convergence of restarted GMRES. SIAM Journal on Matrix Analysis and Applications, 26 (2005).

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

786:   Contributed by: Allison Baker

788: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPGMRES,
789:           KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
790:           KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
791:           KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPLGMRESSetAugDim(),
792:           KSPGMRESSetConstant()

794: M*/

798: PETSC_EXTERN PetscErrorCode KSPCreate_LGMRES(KSP ksp)
799: {
800:   KSP_LGMRES     *lgmres;

804:   PetscNewLog(ksp,&lgmres);

806:   ksp->data               = (void*)lgmres;
807:   ksp->ops->buildsolution = KSPBuildSolution_LGMRES;

809:   ksp->ops->setup                        = KSPSetUp_LGMRES;
810:   ksp->ops->solve                        = KSPSolve_LGMRES;
811:   ksp->ops->destroy                      = KSPDestroy_LGMRES;
812:   ksp->ops->view                         = KSPView_LGMRES;
813:   ksp->ops->setfromoptions               = KSPSetFromOptions_LGMRES;
814:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
815:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_GMRES;

817:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
818:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);

820:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
821:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
822:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
823:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);
824:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);
825:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetHapTol_C",KSPGMRESSetHapTol_GMRES);
826:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
827:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);

829:   /*LGMRES_MOD add extra functions here - like the one to set num of aug vectors */
830:   PetscObjectComposeFunction((PetscObject)ksp,"KSPLGMRESSetConstant_C",KSPLGMRESSetConstant_LGMRES);
831:   PetscObjectComposeFunction((PetscObject)ksp,"KSPLGMRESSetAugDim_C",KSPLGMRESSetAugDim_LGMRES);


834:   /*defaults */
835:   lgmres->haptol         = 1.0e-30;
836:   lgmres->q_preallocate  = 0;
837:   lgmres->delta_allocate = LGMRES_DELTA_DIRECTIONS;
838:   lgmres->orthog         = KSPGMRESClassicalGramSchmidtOrthogonalization;
839:   lgmres->nrs            = 0;
840:   lgmres->sol_temp       = 0;
841:   lgmres->max_k          = LGMRES_DEFAULT_MAXK;
842:   lgmres->Rsvd           = 0;
843:   lgmres->cgstype        = KSP_GMRES_CGS_REFINE_NEVER;
844:   lgmres->orthogwork     = 0;

846:   /*LGMRES_MOD - new defaults */
847:   lgmres->aug_dim         = LGMRES_DEFAULT_AUGDIM;
848:   lgmres->aug_ct          = 0;     /* start with no aug vectors */
849:   lgmres->approx_constant = PETSC_FALSE;
850:   lgmres->matvecs         = 0;
851:   return(0);
852: }