Actual source code: lgmres.c
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 LGMRESGetNewVectors(KSP,PetscInt);
8: static PetscErrorCode LGMRESUpdateHessenberg(KSP,PetscInt,PetscBool ,PetscReal *);
9: static PetscErrorCode BuildLgmresSoln(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: }
34: /*
35: KSPSetUp_LGMRES - Sets up the workspace needed by lgmres.
37: This is called once, usually automatically by KSPSolve() or KSPSetUp(),
38: but can be called directly by KSPSetUp().
40: */
43: PetscErrorCode KSPSetUp_LGMRES(KSP ksp)
44: {
46: PetscInt max_k,k, aug_dim;
47: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
50: max_k = lgmres->max_k;
51: aug_dim = lgmres->aug_dim;
52: KSPSetUp_GMRES(ksp);
54: /* need array of pointers to augvecs*/
55: PetscMalloc((2 * aug_dim + AUG_OFFSET)*sizeof(void*),&lgmres->augvecs);
56: lgmres->aug_vecs_allocated = 2 *aug_dim + AUG_OFFSET;
57: PetscMalloc((2* aug_dim + AUG_OFFSET)*sizeof(void*),&lgmres->augvecs_user_work);
58: PetscMalloc(aug_dim*sizeof(PetscInt),&lgmres->aug_order);
59: PetscLogObjectMemory(ksp,(aug_dim)*(4*sizeof(void*) + sizeof(PetscInt)) + AUG_OFFSET*2*sizeof(void*));
61: /* for now we will preallocate the augvecs - because aug_dim << restart
62: ... also keep in mind that we need to keep augvecs from cycle to cycle*/
63: lgmres->aug_vv_allocated = 2* aug_dim + AUG_OFFSET;
64: lgmres->augwork_alloc = 2* aug_dim + AUG_OFFSET;
65: KSPGetVecs(ksp,lgmres->aug_vv_allocated,&lgmres->augvecs_user_work[0],0,PETSC_NULL);
66: PetscMalloc((max_k+1)*sizeof(PetscScalar),&lgmres->hwork);
67: PetscLogObjectParents(ksp,lgmres->aug_vv_allocated,lgmres->augvecs_user_work[0]);
68: for (k=0; k<lgmres->aug_vv_allocated; k++) {
69: lgmres->augvecs[k] = lgmres->augvecs_user_work[0][k];
70: }
71: return(0);
72: }
75: /*
77: LGMRESCycle - Run lgmres, possibly with restart. Return residual
78: history if requested.
80: input parameters:
81: . lgmres - structure containing parameters and work areas
83: output parameters:
84: . nres - residuals (from preconditioned system) at each step.
85: If restarting, consider passing nres+it. If null,
86: ignored
87: . itcount - number of iterations used. nres[0] to nres[itcount]
88: are defined. If null, ignored. If null, ignored.
89: . converged - 0 if not converged
91:
92: Notes:
93: On entry, the value in vector VEC_VV(0) should be
94: the initial residual.
97: */
100: PetscErrorCode LGMREScycle(PetscInt *itcount,KSP ksp)
101: {
102: KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
103: PetscReal res_norm, res;
104: PetscReal hapbnd, tt;
105: PetscScalar tmp;
106: PetscBool hapend = PETSC_FALSE; /* indicates happy breakdown ending */
108: PetscInt loc_it; /* local count of # of dir. in Krylov space */
109: PetscInt max_k = lgmres->max_k; /* max approx space size */
110: PetscInt max_it = ksp->max_it; /* max # of overall iterations for the method */
111: /* LGMRES_MOD - new variables*/
112: PetscInt aug_dim = lgmres->aug_dim;
113: PetscInt spot = 0;
114: PetscInt order = 0;
115: PetscInt it_arnoldi; /* number of arnoldi steps to take */
116: PetscInt it_total; /* total number of its to take (=approx space size)*/
117: PetscInt ii, jj;
118: PetscReal tmp_norm;
119: PetscScalar inv_tmp_norm;
120: PetscScalar *avec;
123: /* Number of pseudo iterations since last restart is the number
124: of prestart directions */
125: loc_it = 0;
127: /* LGMRES_MOD: determine number of arnoldi steps to take */
128: /* if approx_constant then we keep the space the same size even if
129: we don't have the full number of aug vectors yet*/
130: if (lgmres->approx_constant) {
131: it_arnoldi = max_k - lgmres->aug_ct;
132: } else {
133: it_arnoldi = max_k - aug_dim;
134: }
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: res = res_norm;
141:
142: /* first entry in right-hand-side of hessenberg system is just
143: the initial residual norm */
144: *GRS(0) = res_norm;
146: /* check for the convergence */
147: if (!res) {
148: if (itcount) *itcount = 0;
149: ksp->reason = KSP_CONVERGED_ATOL;
150: PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
151: return(0);
152: }
154: /* scale VEC_VV (the initial residual) */
155: tmp = 1.0/res_norm; VecScale(VEC_VV(0),tmp);
157: ksp->rnorm = res;
160: /* note: (lgmres->it) is always set one less than (loc_it) It is used in
161: KSPBUILDSolution_LGMRES, where it is passed to BuildLgmresSoln.
162: Note that when BuildLgmresSoln is called from this function,
163: (loc_it -1) is passed, so the two are equivalent */
164: lgmres->it = (loc_it - 1);
166:
167: /* MAIN ITERATION LOOP BEGINNING*/
170: /* keep iterating until we have converged OR generated the max number
171: of directions OR reached the max number of iterations for the method */
172: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
173:
174: while (!ksp->reason && loc_it < it_total && ksp->its < max_it) { /* LGMRES_MOD: changed to it_total */
175: KSPLogResidualHistory(ksp,res);
176: lgmres->it = (loc_it - 1);
177: KSPMonitor(ksp,ksp->its,res);
179: /* see if more space is needed for work vectors */
180: if (lgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
181: LGMRESGetNewVectors(ksp,loc_it+1);
182: /* (loc_it+1) is passed in as number of the first vector that should
183: be allocated */
184: }
186: /*LGMRES_MOD: decide whether this is an arnoldi step or an aug step */
187: if (loc_it < it_arnoldi) { /* Arnoldi */
188: KSP_PCApplyBAorAB(ksp,VEC_VV(loc_it),VEC_VV(1+loc_it),VEC_TEMP_MATOP);
189: } else { /*aug step */
190: order = loc_it - it_arnoldi + 1; /* which aug step */
191: for (ii=0; ii<aug_dim; ii++) {
192: if (lgmres->aug_order[ii] == order) {
193: spot = ii;
194: break; /* must have this because there will be duplicates before aug_ct = aug_dim */
195: }
196: }
198: VecCopy(A_AUGVEC(spot), VEC_VV(1+loc_it));
199: /*note: an alternate implementation choice would be to only save the AUGVECS and
200: not A_AUGVEC and then apply the PC here to the augvec */
201: }
203: /* update hessenberg matrix and do Gram-Schmidt - new direction is in
204: VEC_VV(1+loc_it)*/
205: (*lgmres->orthog)(ksp,loc_it);
207: /* new entry in hessenburg is the 2-norm of our new direction */
208: VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);
209: *HH(loc_it+1,loc_it) = tt;
210: *HES(loc_it+1,loc_it) = tt;
213: /* check for the happy breakdown */
214: hapbnd = PetscAbsScalar(tt / *GRS(loc_it));/* GRS(loc_it) contains the res_norm from the last iteration */
215: if (hapbnd > lgmres->haptol) hapbnd = lgmres->haptol;
216: if (tt > hapbnd) {
217: tmp = 1.0/tt;
218: VecScale(VEC_VV(loc_it+1),tmp); /* scale new direction by its norm */
219: } else {
220: PetscInfo2(ksp,"Detected happy breakdown, current hapbnd = %G tt = %G\n",hapbnd,tt);
221: hapend = PETSC_TRUE;
222: }
224: /* Now apply rotations to new col of hessenberg (and right side of system),
225: calculate new rotation, and get new residual norm at the same time*/
226: LGMRESUpdateHessenberg(ksp,loc_it,hapend,&res);
227: if (ksp->reason) break;
229: loc_it++;
230: lgmres->it = (loc_it-1); /* Add this here in case it has converged */
231:
232: PetscObjectTakeAccess(ksp);
233: ksp->its++;
234: ksp->rnorm = res;
235: PetscObjectGrantAccess(ksp);
237: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
239: /* Catch error in happy breakdown and signal convergence and break from loop */
240: if (hapend) {
241: if (!ksp->reason) {
242: SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"You reached the happy break down,but convergence was not indicated. Residual norm = %G",res);
243: }
244: break;
245: }
246: }
247: /* END OF ITERATION LOOP */
249: KSPLogResidualHistory(ksp,res);
251: /* Monitor if we know that we will not return for a restart */
252: if (ksp->reason || ksp->its >= max_it) {
253: KSPMonitor(ksp, ksp->its, res);
254: }
256: if (itcount) *itcount = loc_it;
258: /*
259: Down here we have to solve for the "best" coefficients of the Krylov
260: columns, add the solution values together, and possibly unwind the
261: preconditioning from the solution
262: */
263:
264: /* Form the solution (or the solution so far) */
265: /* Note: must pass in (loc_it-1) for iteration count so that BuildLgmresSoln
266: properly navigates */
268: BuildLgmresSoln(GRS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);
271: /* LGMRES_MOD collect aug vector and A*augvector for future restarts -
272: only if we will be restarting (i.e. this cycle performed it_total
273: iterations) */
274: if (!ksp->reason && ksp->its < max_it && aug_dim > 0) {
276: /*AUG_TEMP contains the new augmentation vector (assigned in BuildLgmresSoln) */
277: if (!lgmres->aug_ct) {
278: spot = 0;
279: lgmres->aug_ct++;
280: } else if (lgmres->aug_ct < aug_dim) {
281: spot = lgmres->aug_ct;
282: lgmres->aug_ct++;
283: } else { /* truncate */
284: for (ii=0; ii<aug_dim; ii++) {
285: if (lgmres->aug_order[ii] == aug_dim) {
286: spot = ii;
287: }
288: }
289: }
291:
293: VecCopy(AUG_TEMP, AUGVEC(spot));
294: /*need to normalize */
295: VecNorm(AUGVEC(spot), NORM_2, &tmp_norm);
296: inv_tmp_norm = 1.0/tmp_norm;
297: VecScale(AUGVEC(spot),inv_tmp_norm);
299: /*set new aug vector to order 1 - move all others back one */
300: for (ii=0; ii < aug_dim; ii++) {
301: AUG_ORDER(ii)++;
302: }
303: AUG_ORDER(spot) = 1;
305: /*now add the A*aug vector to A_AUGVEC(spot) - this is independ. of preconditioning type*/
306: /* want V*H*y - y is in GRS, V is in VEC_VV and H is in HES */
308:
309: /* first do H+*y */
310: avec = lgmres->hwork;
311: PetscMemzero(avec,(it_total+1)*sizeof(*avec));
312: for (ii=0; ii < it_total + 1; ii++) {
313: for (jj=0; jj <= ii+1 && jj < it_total+1; jj++) {
314: avec[jj] += *HES(jj ,ii) * *GRS(ii);
315: }
316: }
318: /*now multiply result by V+ */
319: VecSet(VEC_TEMP,0.0);
320: VecMAXPY(VEC_TEMP, it_total+1, avec, &VEC_VV(0)); /*answer is in VEC_TEMP*/
322: /*copy answer to aug location and scale*/
323: VecCopy(VEC_TEMP, A_AUGVEC(spot));
324: VecScale(A_AUGVEC(spot),inv_tmp_norm);
325: }
326: return(0);
327: }
329: /*
330: KSPSolve_LGMRES - This routine applies the LGMRES method.
333: Input Parameter:
334: . ksp - the Krylov space object that was set to use lgmres
336: Output Parameter:
337: . outits - number of iterations used
339: */
343: PetscErrorCode KSPSolve_LGMRES(KSP ksp)
344: {
346: PetscInt cycle_its; /* iterations done in a call to LGMREScycle */
347: PetscInt itcount; /* running total of iterations, incl. those in restarts */
348: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
349: PetscBool guess_zero = ksp->guess_zero;
350: PetscInt ii; /*LGMRES_MOD variable */
353: if (ksp->calc_sings && !lgmres->Rsvd) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
354:
355: PetscObjectTakeAccess(ksp);
356: ksp->its = 0;
357: lgmres->aug_ct = 0;
358: lgmres->matvecs = 0;
359: PetscObjectGrantAccess(ksp);
361: /* initialize */
362: itcount = 0;
363: ksp->reason = KSP_CONVERGED_ITERATING;
364: /*LGMRES_MOD*/
365: for (ii=0; ii<lgmres->aug_dim; ii++) {
366: lgmres->aug_order[ii] = 0;
367: }
369: while (!ksp->reason) {
370: /* calc residual - puts in VEC_VV(0) */
371: KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
372: LGMREScycle(&cycle_its,ksp);
373: itcount += cycle_its;
374: if (itcount >= ksp->max_it) {
375: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
376: break;
377: }
378: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
379: }
380: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
381: return(0);
382: }
385: /*
387: KSPDestroy_LGMRES - Frees all memory space used by the Krylov method.
389: */
392: PetscErrorCode KSPDestroy_LGMRES(KSP ksp)
393: {
394: KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;
398: PetscFree(lgmres->augvecs);
399: if (lgmres->augwork_alloc) {
400: VecDestroyVecs(lgmres->augwork_alloc,&lgmres->augvecs_user_work[0]);
401: }
402: PetscFree(lgmres->augvecs_user_work);
403: PetscFree(lgmres->aug_order);
404: PetscFree(lgmres->hwork);
405: KSPDestroy_GMRES(ksp);
406: return(0);
407: }
409: /*
410: BuildLgmresSoln - create the solution from the starting vector and the
411: current iterates.
413: Input parameters:
414: nrs - work area of size it + 1.
415: vguess - index of initial guess
416: vdest - index of result. Note that vguess may == vdest (replace
417: guess with the solution).
418: it - HH upper triangular part is a block of size (it+1) x (it+1)
420: This is an internal routine that knows about the LGMRES internals.
421: */
424: static PetscErrorCode BuildLgmresSoln(PetscScalar* nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
425: {
426: PetscScalar tt;
428: PetscInt ii,k,j;
429: KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
430: /*LGMRES_MOD */
431: PetscInt it_arnoldi, it_aug;
432: PetscInt jj, spot = 0;
435: /* Solve for solution vector that minimizes the residual */
437: /* If it is < 0, no lgmres steps have been performed */
438: if (it < 0) {
439: VecCopy(vguess,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
440: return(0);
441: }
443: /* so (it+1) lgmres steps HAVE been performed */
445: /* LGMRES_MOD - determine if we need to use augvecs for the soln - do not assume that
446: this is called after the total its allowed for an approx space */
447: if (lgmres->approx_constant) {
448: it_arnoldi = lgmres->max_k - lgmres->aug_ct;
449: } else {
450: it_arnoldi = lgmres->max_k - lgmres->aug_dim;
451: }
452: if (it_arnoldi >= it +1) {
453: it_aug = 0;
454: it_arnoldi = it+1;
455: } else {
456: it_aug = (it + 1) - it_arnoldi;
457: }
459: /* now it_arnoldi indicates the number of matvecs that took place */
460: lgmres->matvecs += it_arnoldi;
462:
463: /* solve the upper triangular system - GRS is the right side and HH is
464: the upper triangular matrix - put soln in nrs */
465: 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,PetscAbsScalar(*GRS(it)));
466: if (*HH(it,it) != 0.0) {
467: nrs[it] = *GRS(it) / *HH(it,it);
468: } else {
469: nrs[it] = 0.0;
470: }
472: for (ii=1; ii<=it; ii++) {
473: k = it - ii;
474: tt = *GRS(k);
475: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
476: nrs[k] = tt / *HH(k,k);
477: }
479: /* Accumulate the correction to the soln of the preconditioned prob. in VEC_TEMP */
480: VecSet(VEC_TEMP,0.0); /* set VEC_TEMP components to 0 */
482: /*LGMRES_MOD - if augmenting has happened we need to form the solution
483: using the augvecs */
484: if (!it_aug) { /* all its are from arnoldi */
485: VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
486: } else { /*use aug vecs */
487: /*first do regular krylov directions */
488: VecMAXPY(VEC_TEMP,it_arnoldi,nrs,&VEC_VV(0));
489: /*now add augmented portions - add contribution of aug vectors one at a time*/
492: for (ii=0; ii<it_aug; ii++) {
493: for (jj=0; jj<lgmres->aug_dim; jj++) {
494: if (lgmres->aug_order[jj] == (ii+1)) {
495: spot = jj;
496: break; /* must have this because there will be duplicates before aug_ct = aug_dim */
497: }
498: }
499: VecAXPY(VEC_TEMP,nrs[it_arnoldi+ii],AUGVEC(spot));
500: }
501: }
502: /* now VEC_TEMP is what we want to keep for augmenting purposes - grab before the
503: preconditioner is "unwound" from right-precondtioning*/
504: VecCopy(VEC_TEMP, AUG_TEMP);
506: KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
508: /* add solution to previous solution */
509: /* put updated solution into vdest.*/
510: VecCopy(vguess,vdest);
511: VecAXPY(vdest,1.0,VEC_TEMP);
513: return(0);
514: }
516: /*
518: LGMRESUpdateHessenberg - 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
530:
531: */
534: static PetscErrorCode LGMRESUpdateHessenberg(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: #if defined(PETSC_USE_COMPLEX)
553: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
554: #else
555: *hh = *cc * tt + *ss * *(hh+1);
556: #endif
557: hh++;
558: *hh = *cc++ * *hh - (*ss++ * tt);
559: /* hh, cc, and ss have all been incremented one by end of loop */
560: }
562: /*
563: compute the new plane rotation, and apply it to:
564: 1) the right-hand-side of the Hessenberg system (GRS)
565: note: it affects GRS(it) and GRS(it+1)
566: 2) the new column of the Hessenberg matrix
567: note: it affects HH(it,it) which is currently pointed to
568: by hh and HH(it+1, it) (*(hh+1))
569: thus obtaining the updated value of the residual...
570: */
572: /* compute new plane rotation */
574: if (!hapend) {
575: #if defined(PETSC_USE_COMPLEX)
576: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
577: #else
578: tt = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
579: #endif
580: if (tt == 0.0) {
581: ksp->reason = KSP_DIVERGED_NULL;
582: return(0);
583: }
584: *cc = *hh / tt; /* new cosine value */
585: *ss = *(hh+1) / tt; /* new sine value */
587: /* apply to 1) and 2) */
588: *GRS(it+1) = - (*ss * *GRS(it));
589: #if defined(PETSC_USE_COMPLEX)
590: *GRS(it) = PetscConj(*cc) * *GRS(it);
591: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
592: #else
593: *GRS(it) = *cc * *GRS(it);
594: *hh = *cc * *hh + *ss * *(hh+1);
595: #endif
597: /* residual is the last element (it+1) of right-hand side! */
598: *res = PetscAbsScalar(*GRS(it+1));
600: } else { /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
601: another rotation matrix (so RH doesn't change). The new residual is
602: always the new sine term times the residual from last time (GRS(it)),
603: but now the new sine rotation would be zero...so the residual should
604: be zero...so we will multiply "zero" by the last residual. This might
605: not be exactly what we want to do here -could just return "zero". */
606:
607: *res = 0.0;
608: }
609: return(0);
610: }
612: /*
614: LGMRESGetNewVectors - This routine allocates more work vectors, starting from
615: VEC_VV(it)
616:
617: */
620: static PetscErrorCode LGMRESGetNewVectors(KSP ksp,PetscInt it)
621: {
622: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
623: PetscInt nwork = lgmres->nwork_alloc; /* number of work vector chunks allocated */
624: PetscInt nalloc; /* number to allocate */
626: PetscInt k;
627:
629: nalloc = lgmres->delta_allocate; /* number of vectors to allocate
630: in a single chunk */
632: /* Adjust the number to allocate to make sure that we don't exceed the
633: number of available slots (lgmres->vecs_allocated)*/
634: if (it + VEC_OFFSET + nalloc >= lgmres->vecs_allocated){
635: nalloc = lgmres->vecs_allocated - it - VEC_OFFSET;
636: }
637: if (!nalloc) return(0);
639: lgmres->vv_allocated += nalloc; /* vv_allocated is the number of vectors allocated */
641: /* work vectors */
642: KSPGetVecs(ksp,nalloc,&lgmres->user_work[nwork],0,PETSC_NULL);
643: PetscLogObjectParents(ksp,nalloc,lgmres->user_work[nwork]);
644: /* specify size of chunk allocated */
645: lgmres->mwork_alloc[nwork] = nalloc;
647: for (k=0; k < nalloc; k++) {
648: lgmres->vecs[it+VEC_OFFSET+k] = lgmres->user_work[nwork][k];
649: }
650:
652: /* LGMRES_MOD - for now we are preallocating the augmentation vectors */
653:
655: /* increment the number of work vector chunks */
656: lgmres->nwork_alloc++;
657: return(0);
658: }
660: /*
662: KSPBuildSolution_LGMRES
664: Input Parameter:
665: . ksp - the Krylov space object
666: . ptr-
668: Output Parameter:
669: . result - the solution
671: Note: this calls BuildLgmresSoln - the same function that LGMREScycle
672: calls directly.
674: */
677: PetscErrorCode KSPBuildSolution_LGMRES(KSP ksp,Vec ptr,Vec *result)
678: {
679: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
683: if (!ptr) {
684: if (!lgmres->sol_temp) {
685: VecDuplicate(ksp->vec_sol,&lgmres->sol_temp);
686: PetscLogObjectParent(ksp,lgmres->sol_temp);
687: }
688: ptr = lgmres->sol_temp;
689: }
690: if (!lgmres->nrs) {
691: /* allocate the work area */
692: PetscMalloc(lgmres->max_k*sizeof(PetscScalar),&lgmres->nrs);
693: PetscLogObjectMemory(ksp,lgmres->max_k*sizeof(PetscScalar));
694: }
695:
696: BuildLgmresSoln(lgmres->nrs,ksp->vec_sol,ptr,ksp,lgmres->it);
697: if (result) *result = ptr;
698:
699: return(0);
700: }
706: PetscErrorCode KSPView_LGMRES(KSP ksp,PetscViewer viewer)
707: {
708: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
710: PetscBool iascii;
713: KSPView_GMRES(ksp,viewer);
714: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
715: if (iascii) {
716: /*LGMRES_MOD */
717: PetscViewerASCIIPrintf(viewer," LGMRES: aug. dimension=%D\n",lgmres->aug_dim);
718: if (lgmres->approx_constant) {
719: PetscViewerASCIIPrintf(viewer," LGMRES: approx. space size was kept constant.\n");
720: }
721: PetscViewerASCIIPrintf(viewer," LGMRES: number of matvecs=%D\n",lgmres->matvecs);
722: } else {
723: SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for KSP LGMRES",((PetscObject)viewer)->type_name);
724: }
725: return(0);
726: }
732: PetscErrorCode KSPSetFromOptions_LGMRES(KSP ksp)
733: {
735: PetscInt aug;
736: KSP_LGMRES *lgmres = (KSP_LGMRES*) ksp->data;
737: PetscBool flg = PETSC_FALSE;
740: KSPSetFromOptions_GMRES(ksp);
741: PetscOptionsHead("KSP LGMRES Options");
742: PetscOptionsBool("-ksp_lgmres_constant","Use constant approx. space size","KSPGMRESSetConstant",flg,&flg,PETSC_NULL);
743: if (flg) { lgmres->approx_constant = 1; }
744: PetscOptionsInt("-ksp_lgmres_augment","Number of error approximations to augment the Krylov space with","KSPLGMRESSetAugDim",lgmres->aug_dim,&aug,&flg);
745: if (flg) { KSPLGMRESSetAugDim(ksp,aug); }
746: PetscOptionsTail();
747: return(0);
748: }
754: /*functions for extra lgmres options here*/
758: PetscErrorCode KSPLGMRESSetConstant_LGMRES(KSP ksp)
759: {
760: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
762: lgmres->approx_constant = 1;
763: return(0);
764: }
770: PetscErrorCode KSPLGMRESSetAugDim_LGMRES(KSP ksp,PetscInt aug_dim)
771: {
772: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
775: if (aug_dim < 0) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be positive");
776: if (aug_dim > (lgmres->max_k -1)) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be <= (restart size-1)");
777: lgmres->aug_dim = aug_dim;
778: return(0);
779: }
783: /* end new lgmres functions */
786: /* use these options from gmres */
798: /*MC
799: KSPLGMRES - Augments the standard GMRES approximation space with approximations to
800: the error from previous restart cycles.
802: Options Database Keys:
803: + -ksp_gmres_restart <restart> - total approximation space size (Krylov directions + error approximations)
804: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
805: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
806: vectors are allocated as needed)
807: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
808: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
809: . -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
810: stability of the classical Gram-Schmidt orthogonalization.
811: . -ksp_gmres_krylov_monitor - plot the Krylov space generated
812: . -ksp_lgmres_augment <k> - number of error approximations to augment the Krylov space with
813: - -ksp_lgmres_constant - use a constant approx. space size (only affects restart cycles < num. error approx.(k), i.e. the first k restarts)
815: To run LGMRES(m, k) as described in the above paper, use:
816: -ksp_gmres_restart <m+k>
817: -ksp_lgmres_augment <k>
819: Level: beginner
821: Notes: Supports both left and right preconditioning, but not symmetric.
823: References:
824: A. H. Baker, E.R. Jessup, and T.A. Manteuffel. A technique for
825: accelerating the convergence of restarted GMRES. SIAM
826: Journal on Matrix Analysis and Applications, 26 (2005), pp. 962-984.
828: Developer Notes: This object is subclassed off of KSPGMRES
830: Contributed by: Allison Baker
832: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPGMRES,
833: KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
834: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
835: KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPLGMRESSetAugDim(),
836: KSPGMRESSetConstant()
838: M*/
843: PetscErrorCode KSPCreate_LGMRES(KSP ksp)
844: {
845: KSP_LGMRES *lgmres;
849: PetscNewLog(ksp,KSP_LGMRES,&lgmres);
850: ksp->data = (void*)lgmres;
851: ksp->ops->buildsolution = KSPBuildSolution_LGMRES;
853: ksp->ops->setup = KSPSetUp_LGMRES;
854: ksp->ops->solve = KSPSolve_LGMRES;
855: ksp->ops->destroy = KSPDestroy_LGMRES;
856: ksp->ops->view = KSPView_LGMRES;
857: ksp->ops->setfromoptions = KSPSetFromOptions_LGMRES;
858: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
859: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
861: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
862: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);
864: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
865: "KSPGMRESSetPreAllocateVectors_GMRES",
866: KSPGMRESSetPreAllocateVectors_GMRES);
867: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
868: "KSPGMRESSetOrthogonalization_GMRES",
869: KSPGMRESSetOrthogonalization_GMRES);
870: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",
871: "KSPGMRESGetOrthogonalization_GMRES",
872: KSPGMRESGetOrthogonalization_GMRES);
873: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
874: "KSPGMRESSetRestart_GMRES",
875: KSPGMRESSetRestart_GMRES);
876: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetRestart_C",
877: "KSPGMRESGetRestart_GMRES",
878: KSPGMRESGetRestart_GMRES);
879: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
880: "KSPGMRESSetHapTol_GMRES",
881: KSPGMRESSetHapTol_GMRES);
882: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
883: "KSPGMRESSetCGSRefinementType_GMRES",
884: KSPGMRESSetCGSRefinementType_GMRES);
885: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",
886: "KSPGMRESGetCGSRefinementType_GMRES",
887: KSPGMRESGetCGSRefinementType_GMRES);
889: /*LGMRES_MOD add extra functions here - like the one to set num of aug vectors */
890: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetConstant_C",
891: "KSPLGMRESSetConstant_LGMRES",
892: KSPLGMRESSetConstant_LGMRES);
894: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetAugDim_C",
895: "KSPLGMRESSetAugDim_LGMRES",
896: KSPLGMRESSetAugDim_LGMRES);
897:
899: /*defaults */
900: lgmres->haptol = 1.0e-30;
901: lgmres->q_preallocate = 0;
902: lgmres->delta_allocate = LGMRES_DELTA_DIRECTIONS;
903: lgmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
904: lgmres->nrs = 0;
905: lgmres->sol_temp = 0;
906: lgmres->max_k = LGMRES_DEFAULT_MAXK;
907: lgmres->Rsvd = 0;
908: lgmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
909: lgmres->orthogwork = 0;
910: /*LGMRES_MOD - new defaults */
911: lgmres->aug_dim = LGMRES_DEFAULT_AUGDIM;
912: lgmres->aug_ct = 0; /* start with no aug vectors */
913: lgmres->approx_constant = 0;
914: lgmres->matvecs = 0;
916: return(0);
917: }