Actual source code: lgmres.c
2: #include <../src/ksp/ksp/impls/gmres/lgmres/lgmresimpl.h>
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);
11: PetscErrorCode KSPLGMRESSetAugDim(KSP ksp, PetscInt dim)
12: {
16: PetscTryMethod((ksp),"KSPLGMRESSetAugDim_C",(KSP,PetscInt),(ksp,dim));
17: return(0);
18: }
20: PetscErrorCode KSPLGMRESSetConstant(KSP ksp)
21: {
25: PetscTryMethod((ksp),"KSPLGMRESSetConstant_C",(KSP),(ksp));
26: return(0);
27: }
29: /*
30: KSPSetUp_LGMRES - Sets up the workspace needed by lgmres.
32: This is called once, usually automatically by KSPSolve() or KSPSetUp(),
33: but can be called directly by KSPSetUp().
35: */
36: PetscErrorCode KSPSetUp_LGMRES(KSP ksp)
37: {
39: PetscInt max_k,k, aug_dim;
40: KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;
43: max_k = lgmres->max_k;
44: aug_dim = lgmres->aug_dim;
45: KSPSetUp_GMRES(ksp);
47: /* need array of pointers to augvecs*/
48: PetscMalloc1(2*aug_dim + AUG_OFFSET,&lgmres->augvecs);
50: lgmres->aug_vecs_allocated = 2 *aug_dim + AUG_OFFSET;
52: PetscMalloc1(2*aug_dim + AUG_OFFSET,&lgmres->augvecs_user_work);
53: PetscMalloc1(aug_dim,&lgmres->aug_order);
54: PetscLogObjectMemory((PetscObject)ksp,(aug_dim)*(4*sizeof(void*) + sizeof(PetscInt)) + AUG_OFFSET*2*sizeof(void*));
56: /* for now we will preallocate the augvecs - because aug_dim << restart
57: ... also keep in mind that we need to keep augvecs from cycle to cycle*/
58: lgmres->aug_vv_allocated = 2* aug_dim + AUG_OFFSET;
59: lgmres->augwork_alloc = 2* aug_dim + AUG_OFFSET;
61: KSPCreateVecs(ksp,lgmres->aug_vv_allocated,&lgmres->augvecs_user_work[0],0,NULL);
62: PetscMalloc1(max_k+1,&lgmres->hwork);
63: PetscLogObjectParents(ksp,lgmres->aug_vv_allocated,lgmres->augvecs_user_work[0]);
64: for (k=0; k<lgmres->aug_vv_allocated; k++) {
65: lgmres->augvecs[k] = lgmres->augvecs_user_work[0][k];
66: }
67: return(0);
68: }
70: /*
72: KSPLGMRESCycle - Run lgmres, possibly with restart. Return residual
73: history if requested.
75: input parameters:
76: . lgmres - structure containing parameters and work areas
78: output parameters:
79: . nres - residuals (from preconditioned system) at each step.
80: If restarting, consider passing nres+it. If null,
81: ignored
82: . itcount - number of iterations used. nres[0] to nres[itcount]
83: are defined. If null, ignored. If null, ignored.
84: . converged - 0 if not converged
86: Notes:
87: On entry, the value in vector VEC_VV(0) should be
88: the initial residual.
90: */
91: PetscErrorCode KSPLGMRESCycle(PetscInt *itcount,KSP ksp)
92: {
93: KSP_LGMRES *lgmres = (KSP_LGMRES*)(ksp->data);
94: PetscReal res_norm, res;
95: PetscReal hapbnd, tt;
96: PetscScalar tmp;
97: PetscBool hapend = PETSC_FALSE; /* indicates happy breakdown ending */
99: PetscInt loc_it; /* local count of # of dir. in Krylov space */
100: PetscInt max_k = lgmres->max_k; /* max approx space size */
101: PetscInt max_it = ksp->max_it; /* max # of overall iterations for the method */
103: /* LGMRES_MOD - new variables*/
104: PetscInt aug_dim = lgmres->aug_dim;
105: PetscInt spot = 0;
106: PetscInt order = 0;
107: PetscInt it_arnoldi; /* number of arnoldi steps to take */
108: PetscInt it_total; /* total number of its to take (=approx space size)*/
109: PetscInt ii, jj;
110: PetscReal tmp_norm;
111: PetscScalar inv_tmp_norm;
112: PetscScalar *avec;
115: /* Number of pseudo iterations since last restart is the number
116: of prestart directions */
117: loc_it = 0;
119: /* LGMRES_MOD: determine number of arnoldi steps to take */
120: /* if approx_constant then we keep the space the same size even if
121: we don't have the full number of aug vectors yet*/
122: if (lgmres->approx_constant) it_arnoldi = max_k - lgmres->aug_ct;
123: else it_arnoldi = max_k - aug_dim;
125: it_total = it_arnoldi + lgmres->aug_ct;
127: /* initial residual is in VEC_VV(0) - compute its norm*/
128: VecNorm(VEC_VV(0),NORM_2,&res_norm);
129: KSPCheckNorm(ksp,res_norm);
130: res = res_norm;
132: /* first entry in right-hand-side of hessenberg system is just
133: the initial residual norm */
134: *GRS(0) = res_norm;
136: /* check for the convergence */
137: if (!res) {
138: if (itcount) *itcount = 0;
139: ksp->reason = KSP_CONVERGED_ATOL;
140: PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
141: return(0);
142: }
144: /* scale VEC_VV (the initial residual) */
145: tmp = 1.0/res_norm; VecScale(VEC_VV(0),tmp);
147: if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = res;
148: else ksp->rnorm = 0.0;
150: /* note: (lgmres->it) is always set one less than (loc_it) It is used in
151: KSPBUILDSolution_LGMRES, where it is passed to KSPLGMRESBuildSoln.
152: Note that when KSPLGMRESBuildSoln is called from this function,
153: (loc_it -1) is passed, so the two are equivalent */
154: lgmres->it = (loc_it - 1);
156: /* MAIN ITERATION LOOP BEGINNING*/
158: /* keep iterating until we have converged OR generated the max number
159: of directions OR reached the max number of iterations for the method */
160: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
162: while (!ksp->reason && loc_it < it_total && ksp->its < max_it) { /* LGMRES_MOD: changed to it_total */
163: KSPLogResidualHistory(ksp,res);
164: lgmres->it = (loc_it - 1);
165: KSPMonitor(ksp,ksp->its,res);
167: /* see if more space is needed for work vectors */
168: if (lgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
169: KSPLGMRESGetNewVectors(ksp,loc_it+1);
170: /* (loc_it+1) is passed in as number of the first vector that should
171: be allocated */
172: }
174: /*LGMRES_MOD: decide whether this is an arnoldi step or an aug step */
175: if (loc_it < it_arnoldi) { /* Arnoldi */
176: KSP_PCApplyBAorAB(ksp,VEC_VV(loc_it),VEC_VV(1+loc_it),VEC_TEMP_MATOP);
177: } else { /*aug step */
178: order = loc_it - it_arnoldi + 1; /* which aug step */
179: for (ii=0; ii<aug_dim; ii++) {
180: if (lgmres->aug_order[ii] == order) {
181: spot = ii;
182: break; /* must have this because there will be duplicates before aug_ct = aug_dim */
183: }
184: }
186: VecCopy(A_AUGVEC(spot), VEC_VV(1+loc_it));
187: /*note: an alternate implementation choice would be to only save the AUGVECS and
188: not A_AUGVEC and then apply the PC here to the augvec */
189: }
191: /* update hessenberg matrix and do Gram-Schmidt - new direction is in
192: VEC_VV(1+loc_it)*/
193: (*lgmres->orthog)(ksp,loc_it);
195: /* new entry in hessenburg is the 2-norm of our new direction */
196: VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);
198: *HH(loc_it+1,loc_it) = tt;
199: *HES(loc_it+1,loc_it) = tt;
201: /* check for the happy breakdown */
202: hapbnd = PetscAbsScalar(tt / *GRS(loc_it)); /* GRS(loc_it) contains the res_norm from the last iteration */
203: if (hapbnd > lgmres->haptol) hapbnd = lgmres->haptol;
204: if (tt > hapbnd) {
205: tmp = 1.0/tt;
206: VecScale(VEC_VV(loc_it+1),tmp); /* scale new direction by its norm */
207: } else {
208: PetscInfo2(ksp,"Detected happy breakdown, current hapbnd = %g tt = %g\n",(double)hapbnd,(double)tt);
209: hapend = PETSC_TRUE;
210: }
212: /* Now apply rotations to new col of hessenberg (and right side of system),
213: calculate new rotation, and get new residual norm at the same time*/
214: KSPLGMRESUpdateHessenberg(ksp,loc_it,hapend,&res);
215: if (ksp->reason) break;
217: loc_it++;
218: lgmres->it = (loc_it-1); /* Add this here in case it has converged */
220: PetscObjectSAWsTakeAccess((PetscObject)ksp);
221: ksp->its++;
222: if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = res;
223: else ksp->rnorm = 0.0;
224: PetscObjectSAWsGrantAccess((PetscObject)ksp);
226: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
228: /* Catch error in happy breakdown and signal convergence and break from loop */
229: if (hapend) {
230: if (!ksp->reason) {
231: 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);
232: else {
233: ksp->reason = KSP_DIVERGED_BREAKDOWN;
234: break;
235: }
236: }
237: }
238: }
239: /* END OF ITERATION LOOP */
240: KSPLogResidualHistory(ksp,res);
242: /* Monitor if we know that we will not return for a restart */
243: if (ksp->reason || ksp->its >= max_it) {
244: KSPMonitor(ksp, ksp->its, res);
245: }
247: if (itcount) *itcount = loc_it;
249: /*
250: Down here we have to solve for the "best" coefficients of the Krylov
251: columns, add the solution values together, and possibly unwind the
252: preconditioning from the solution
253: */
255: /* Form the solution (or the solution so far) */
256: /* Note: must pass in (loc_it-1) for iteration count so that KSPLGMRESBuildSoln
257: properly navigates */
259: KSPLGMRESBuildSoln(GRS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);
261: /* LGMRES_MOD collect aug vector and A*augvector for future restarts -
262: only if we will be restarting (i.e. this cycle performed it_total
263: iterations) */
264: if (!ksp->reason && ksp->its < max_it && aug_dim > 0) {
266: /*AUG_TEMP contains the new augmentation vector (assigned in KSPLGMRESBuildSoln) */
267: if (!lgmres->aug_ct) {
268: spot = 0;
269: lgmres->aug_ct++;
270: } else if (lgmres->aug_ct < aug_dim) {
271: spot = lgmres->aug_ct;
272: lgmres->aug_ct++;
273: } else { /* truncate */
274: for (ii=0; ii<aug_dim; ii++) {
275: if (lgmres->aug_order[ii] == aug_dim) spot = ii;
276: }
277: }
279: VecCopy(AUG_TEMP, AUGVEC(spot));
280: /*need to normalize */
281: VecNorm(AUGVEC(spot), NORM_2, &tmp_norm);
283: inv_tmp_norm = 1.0/tmp_norm;
285: VecScale(AUGVEC(spot),inv_tmp_norm);
287: /*set new aug vector to order 1 - move all others back one */
288: for (ii=0; ii < aug_dim; ii++) AUG_ORDER(ii)++;
289: AUG_ORDER(spot) = 1;
291: /*now add the A*aug vector to A_AUGVEC(spot) - this is independ. of preconditioning type*/
292: /* want V*H*y - y is in GRS, V is in VEC_VV and H is in HES */
294: /* first do H+*y */
295: avec = lgmres->hwork;
296: PetscArrayzero(avec,it_total+1);
297: for (ii=0; ii < it_total + 1; ii++) {
298: for (jj=0; jj <= ii+1 && jj < it_total+1; jj++) {
299: avec[jj] += *HES(jj ,ii) * *GRS(ii);
300: }
301: }
303: /*now multiply result by V+ */
304: VecSet(VEC_TEMP,0.0);
305: VecMAXPY(VEC_TEMP, it_total+1, avec, &VEC_VV(0)); /*answer is in VEC_TEMP*/
307: /*copy answer to aug location and scale*/
308: VecCopy(VEC_TEMP, A_AUGVEC(spot));
309: VecScale(A_AUGVEC(spot),inv_tmp_norm);
310: }
311: return(0);
312: }
314: /*
315: KSPSolve_LGMRES - This routine applies the LGMRES method.
317: Input Parameter:
318: . ksp - the Krylov space object that was set to use lgmres
320: Output Parameter:
321: . outits - number of iterations used
323: */
325: PetscErrorCode KSPSolve_LGMRES(KSP ksp)
326: {
328: PetscInt cycle_its; /* iterations done in a call to KSPLGMRESCycle */
329: PetscInt itcount; /* running total of iterations, incl. those in restarts */
330: KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;
331: PetscBool guess_zero = ksp->guess_zero;
332: PetscInt ii; /*LGMRES_MOD variable */
335: if (ksp->calc_sings && !lgmres->Rsvd) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
337: PetscObjectSAWsTakeAccess((PetscObject)ksp);
339: ksp->its = 0;
340: lgmres->aug_ct = 0;
341: lgmres->matvecs = 0;
343: PetscObjectSAWsGrantAccess((PetscObject)ksp);
345: /* initialize */
346: itcount = 0;
347: ksp->reason = KSP_CONVERGED_ITERATING;
348: /*LGMRES_MOD*/
349: for (ii=0; ii<lgmres->aug_dim; ii++) lgmres->aug_order[ii] = 0;
351: while (!ksp->reason) {
352: /* calc residual - puts in VEC_VV(0) */
353: KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
354: KSPLGMRESCycle(&cycle_its,ksp);
355: itcount += cycle_its;
356: if (itcount >= ksp->max_it) {
357: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
358: break;
359: }
360: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
361: }
362: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
363: return(0);
364: }
366: /*
368: KSPDestroy_LGMRES - Frees all memory space used by the Krylov method.
370: */
371: PetscErrorCode KSPDestroy_LGMRES(KSP ksp)
372: {
373: KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;
377: PetscFree(lgmres->augvecs);
378: if (lgmres->augwork_alloc) {
379: VecDestroyVecs(lgmres->augwork_alloc,&lgmres->augvecs_user_work[0]);
380: }
381: PetscFree(lgmres->augvecs_user_work);
382: PetscFree(lgmres->aug_order);
383: PetscFree(lgmres->hwork);
384: KSPDestroy_GMRES(ksp);
385: return(0);
386: }
388: /*
389: KSPLGMRESBuildSoln - create the solution from the starting vector and the
390: current iterates.
392: Input parameters:
393: nrs - work area of size it + 1.
394: vguess - index of initial guess
395: vdest - index of result. Note that vguess may == vdest (replace
396: guess with the solution).
397: it - HH upper triangular part is a block of size (it+1) x (it+1)
399: This is an internal routine that knows about the LGMRES internals.
400: */
401: static PetscErrorCode KSPLGMRESBuildSoln(PetscScalar *nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
402: {
403: PetscScalar tt;
405: PetscInt ii,k,j;
406: KSP_LGMRES *lgmres = (KSP_LGMRES*)(ksp->data);
407: /*LGMRES_MOD */
408: PetscInt it_arnoldi, it_aug;
409: PetscInt jj, spot = 0;
412: /* Solve for solution vector that minimizes the residual */
414: /* If it is < 0, no lgmres steps have been performed */
415: if (it < 0) {
416: VecCopy(vguess,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
417: return(0);
418: }
420: /* so (it+1) lgmres steps HAVE been performed */
422: /* LGMRES_MOD - determine if we need to use augvecs for the soln - do not assume that
423: this is called after the total its allowed for an approx space */
424: if (lgmres->approx_constant) {
425: it_arnoldi = lgmres->max_k - lgmres->aug_ct;
426: } else {
427: it_arnoldi = lgmres->max_k - lgmres->aug_dim;
428: }
429: if (it_arnoldi >= it +1) {
430: it_aug = 0;
431: it_arnoldi = it+1;
432: } else {
433: it_aug = (it + 1) - it_arnoldi;
434: }
436: /* now it_arnoldi indicates the number of matvecs that took place */
437: lgmres->matvecs += it_arnoldi;
439: /* solve the upper triangular system - GRS is the right side and HH is
440: the upper triangular matrix - put soln in nrs */
441: 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)));
442: if (*HH(it,it) != 0.0) {
443: nrs[it] = *GRS(it) / *HH(it,it);
444: } else {
445: nrs[it] = 0.0;
446: }
448: for (ii=1; ii<=it; ii++) {
449: k = it - ii;
450: tt = *GRS(k);
451: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
452: nrs[k] = tt / *HH(k,k);
453: }
455: /* Accumulate the correction to the soln of the preconditioned prob. in VEC_TEMP */
456: VecSet(VEC_TEMP,0.0); /* set VEC_TEMP components to 0 */
458: /*LGMRES_MOD - if augmenting has happened we need to form the solution
459: using the augvecs */
460: if (!it_aug) { /* all its are from arnoldi */
461: VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
462: } else { /*use aug vecs */
463: /*first do regular krylov directions */
464: VecMAXPY(VEC_TEMP,it_arnoldi,nrs,&VEC_VV(0));
465: /*now add augmented portions - add contribution of aug vectors one at a time*/
467: for (ii=0; ii<it_aug; ii++) {
468: for (jj=0; jj<lgmres->aug_dim; jj++) {
469: if (lgmres->aug_order[jj] == (ii+1)) {
470: spot = jj;
471: break; /* must have this because there will be duplicates before aug_ct = aug_dim */
472: }
473: }
474: VecAXPY(VEC_TEMP,nrs[it_arnoldi+ii],AUGVEC(spot));
475: }
476: }
477: /* now VEC_TEMP is what we want to keep for augmenting purposes - grab before the
478: preconditioner is "unwound" from right-precondtioning*/
479: VecCopy(VEC_TEMP, AUG_TEMP);
481: KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
483: /* add solution to previous solution */
484: /* put updated solution into vdest.*/
485: VecCopy(vguess,vdest);
486: VecAXPY(vdest,1.0,VEC_TEMP);
487: return(0);
488: }
490: /*
492: KSPLGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
493: Return new residual.
495: input parameters:
497: . ksp - Krylov space object
498: . it - plane rotations are applied to the (it+1)th column of the
499: modified hessenberg (i.e. HH(:,it))
500: . hapend - PETSC_FALSE not happy breakdown ending.
502: output parameters:
503: . res - the new residual
505: */
506: static PetscErrorCode KSPLGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool hapend,PetscReal *res)
507: {
508: PetscScalar *hh,*cc,*ss,tt;
509: PetscInt j;
510: KSP_LGMRES *lgmres = (KSP_LGMRES*)(ksp->data);
513: hh = HH(0,it); /* pointer to beginning of column to update - so
514: incrementing hh "steps down" the (it+1)th col of HH*/
515: cc = CC(0); /* beginning of cosine rotations */
516: ss = SS(0); /* beginning of sine rotations */
518: /* Apply all the previously computed plane rotations to the new column
519: of the Hessenberg matrix */
520: /* Note: this uses the rotation [conj(c) s ; -s c], c= cos(theta), s= sin(theta) */
522: for (j=1; j<=it; j++) {
523: tt = *hh;
524: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
525: hh++;
526: *hh = *cc++ * *hh - (*ss++ * tt);
527: /* hh, cc, and ss have all been incremented one by end of loop */
528: }
530: /*
531: compute the new plane rotation, and apply it to:
532: 1) the right-hand-side of the Hessenberg system (GRS)
533: note: it affects GRS(it) and GRS(it+1)
534: 2) the new column of the Hessenberg matrix
535: note: it affects HH(it,it) which is currently pointed to
536: by hh and HH(it+1, it) (*(hh+1))
537: thus obtaining the updated value of the residual...
538: */
540: /* compute new plane rotation */
542: if (!hapend) {
543: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
544: if (tt == 0.0) {
545: ksp->reason = KSP_DIVERGED_NULL;
546: return(0);
547: }
548: *cc = *hh / tt; /* new cosine value */
549: *ss = *(hh+1) / tt; /* new sine value */
551: /* apply to 1) and 2) */
552: *GRS(it+1) = -(*ss * *GRS(it));
553: *GRS(it) = PetscConj(*cc) * *GRS(it);
554: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
556: /* residual is the last element (it+1) of right-hand side! */
557: *res = PetscAbsScalar(*GRS(it+1));
559: } else { /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply
560: another rotation matrix (so RH doesn't change). The new residual is
561: always the new sine term times the residual from last time (GRS(it)),
562: but now the new sine rotation would be zero...so the residual should
563: be zero...so we will multiply "zero" by the last residual. This might
564: not be exactly what we want to do here -could just return "zero". */
566: *res = 0.0;
567: }
568: return(0);
569: }
571: /*
573: KSPLGMRESGetNewVectors - This routine allocates more work vectors, starting from
574: VEC_VV(it)
576: */
577: static PetscErrorCode KSPLGMRESGetNewVectors(KSP ksp,PetscInt it)
578: {
579: KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;
580: PetscInt nwork = lgmres->nwork_alloc; /* number of work vector chunks allocated */
581: PetscInt nalloc; /* number to allocate */
583: PetscInt k;
586: nalloc = lgmres->delta_allocate; /* number of vectors to allocate
587: in a single chunk */
589: /* Adjust the number to allocate to make sure that we don't exceed the
590: number of available slots (lgmres->vecs_allocated)*/
591: if (it + VEC_OFFSET + nalloc >= lgmres->vecs_allocated) {
592: nalloc = lgmres->vecs_allocated - it - VEC_OFFSET;
593: }
594: if (!nalloc) return(0);
596: lgmres->vv_allocated += nalloc; /* vv_allocated is the number of vectors allocated */
598: /* work vectors */
599: KSPCreateVecs(ksp,nalloc,&lgmres->user_work[nwork],0,NULL);
600: PetscLogObjectParents(ksp,nalloc,lgmres->user_work[nwork]);
601: /* specify size of chunk allocated */
602: lgmres->mwork_alloc[nwork] = nalloc;
604: for (k=0; k < nalloc; k++) {
605: lgmres->vecs[it+VEC_OFFSET+k] = lgmres->user_work[nwork][k];
606: }
608: /* LGMRES_MOD - for now we are preallocating the augmentation vectors */
610: /* increment the number of work vector chunks */
611: lgmres->nwork_alloc++;
612: return(0);
613: }
615: /*
617: KSPBuildSolution_LGMRES
619: Input Parameter:
620: . ksp - the Krylov space object
621: . ptr-
623: Output Parameter:
624: . result - the solution
626: Note: this calls KSPLGMRESBuildSoln - the same function that KSPLGMRESCycle
627: calls directly.
629: */
630: PetscErrorCode KSPBuildSolution_LGMRES(KSP ksp,Vec ptr,Vec *result)
631: {
632: KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;
636: if (!ptr) {
637: if (!lgmres->sol_temp) {
638: VecDuplicate(ksp->vec_sol,&lgmres->sol_temp);
639: PetscLogObjectParent((PetscObject)ksp,(PetscObject)lgmres->sol_temp);
640: }
641: ptr = lgmres->sol_temp;
642: }
643: if (!lgmres->nrs) {
644: /* allocate the work area */
645: PetscMalloc1(lgmres->max_k,&lgmres->nrs);
646: PetscLogObjectMemory((PetscObject)ksp,lgmres->max_k*sizeof(PetscScalar));
647: }
649: KSPLGMRESBuildSoln(lgmres->nrs,ksp->vec_sol,ptr,ksp,lgmres->it);
650: if (result) *result = ptr;
651: return(0);
652: }
654: PetscErrorCode KSPView_LGMRES(KSP ksp,PetscViewer viewer)
655: {
656: KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;
658: PetscBool iascii;
661: KSPView_GMRES(ksp,viewer);
662: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
663: if (iascii) {
664: /*LGMRES_MOD */
665: PetscViewerASCIIPrintf(viewer," aug. dimension=%D\n",lgmres->aug_dim);
666: if (lgmres->approx_constant) {
667: PetscViewerASCIIPrintf(viewer," approx. space size was kept constant.\n");
668: }
669: PetscViewerASCIIPrintf(viewer," number of matvecs=%D\n",lgmres->matvecs);
670: }
671: return(0);
672: }
674: PetscErrorCode KSPSetFromOptions_LGMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
675: {
677: PetscInt aug;
678: KSP_LGMRES *lgmres = (KSP_LGMRES*) ksp->data;
679: PetscBool flg = PETSC_FALSE;
682: KSPSetFromOptions_GMRES(PetscOptionsObject,ksp);
683: PetscOptionsHead(PetscOptionsObject,"KSP LGMRES Options");
684: PetscOptionsBool("-ksp_lgmres_constant","Use constant approx. space size","KSPGMRESSetConstant",lgmres->approx_constant,&lgmres->approx_constant,NULL);
685: PetscOptionsInt("-ksp_lgmres_augment","Number of error approximations to augment the Krylov space with","KSPLGMRESSetAugDim",lgmres->aug_dim,&aug,&flg);
686: if (flg) { KSPLGMRESSetAugDim(ksp,aug); }
687: PetscOptionsTail();
688: return(0);
689: }
691: /*functions for extra lgmres options here*/
692: static PetscErrorCode KSPLGMRESSetConstant_LGMRES(KSP ksp)
693: {
694: KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;
697: lgmres->approx_constant = PETSC_TRUE;
698: return(0);
699: }
701: static PetscErrorCode KSPLGMRESSetAugDim_LGMRES(KSP ksp,PetscInt aug_dim)
702: {
703: KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;
706: if (aug_dim < 0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be positive");
707: if (aug_dim > (lgmres->max_k -1)) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be <= (restart size-1)");
708: lgmres->aug_dim = aug_dim;
709: return(0);
710: }
712: /* end new lgmres functions */
714: /*MC
715: KSPLGMRES - Augments the standard GMRES approximation space with approximations to
716: the error from previous restart cycles.
718: Options Database Keys:
719: + -ksp_gmres_restart <restart> - total approximation space size (Krylov directions + error approximations)
720: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
721: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
722: vectors are allocated as needed)
723: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
724: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
725: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the
726: stability of the classical Gram-Schmidt orthogonalization.
727: . -ksp_gmres_krylov_monitor - plot the Krylov space generated
728: . -ksp_lgmres_augment <k> - number of error approximations to augment the Krylov space with
729: - -ksp_lgmres_constant - use a constant approx. space size (only affects restart cycles < num. error approx.(k), i.e. the first k restarts)
731: To run LGMRES(m, k) as described in the above paper, use:
732: -ksp_gmres_restart <m+k>
733: -ksp_lgmres_augment <k>
735: Level: beginner
737: Notes:
738: Supports both left and right preconditioning, but not symmetric.
740: References:
741: . 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).
743: Developer Notes:
744: This object is subclassed off of KSPGMRES
746: Contributed by: Allison Baker
748: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPGMRES,
749: KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
750: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
751: KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPLGMRESSetAugDim(),
752: KSPGMRESSetConstant()
754: M*/
756: PETSC_EXTERN PetscErrorCode KSPCreate_LGMRES(KSP ksp)
757: {
758: KSP_LGMRES *lgmres;
762: PetscNewLog(ksp,&lgmres);
764: ksp->data = (void*)lgmres;
765: ksp->ops->buildsolution = KSPBuildSolution_LGMRES;
767: ksp->ops->setup = KSPSetUp_LGMRES;
768: ksp->ops->solve = KSPSolve_LGMRES;
769: ksp->ops->destroy = KSPDestroy_LGMRES;
770: ksp->ops->view = KSPView_LGMRES;
771: ksp->ops->setfromoptions = KSPSetFromOptions_LGMRES;
772: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
773: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
775: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
776: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
777: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
779: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
780: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
781: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
782: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_GMRES);
783: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_GMRES);
784: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetHapTol_C",KSPGMRESSetHapTol_GMRES);
785: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
786: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);
788: /*LGMRES_MOD add extra functions here - like the one to set num of aug vectors */
789: PetscObjectComposeFunction((PetscObject)ksp,"KSPLGMRESSetConstant_C",KSPLGMRESSetConstant_LGMRES);
790: PetscObjectComposeFunction((PetscObject)ksp,"KSPLGMRESSetAugDim_C",KSPLGMRESSetAugDim_LGMRES);
792: /*defaults */
793: lgmres->haptol = 1.0e-30;
794: lgmres->q_preallocate = 0;
795: lgmres->delta_allocate = LGMRES_DELTA_DIRECTIONS;
796: lgmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
797: lgmres->nrs = NULL;
798: lgmres->sol_temp = NULL;
799: lgmres->max_k = LGMRES_DEFAULT_MAXK;
800: lgmres->Rsvd = NULL;
801: lgmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
802: lgmres->orthogwork = NULL;
804: /*LGMRES_MOD - new defaults */
805: lgmres->aug_dim = LGMRES_DEFAULT_AUGDIM;
806: lgmres->aug_ct = 0; /* start with no aug vectors */
807: lgmres->approx_constant = PETSC_FALSE;
808: lgmres->matvecs = 0;
809: return(0);
810: }