Actual source code: lgmres.c
1: /* The LGMRES method
3: Contributed by: Allison Baker
5: Augments the standard GMRES approximation space with approximation to
6: the error from previous restart cycles.
8: Can be combined with left or right preconditioning.
10: Described in:
11: A. H. Baker, E.R. Jessup, and T.A. Manteuffel. A technique for
12: accelerating the convergence of restarted GMRES. Submitted to SIAM
13: Journal on Matrix Analysis and Applications. Also available as
14: Technical Report #CU-CS-945-03, University of Colorado, Department of
15: Computer Science, January, 2003.
17: */
19: #include lgmresp.h
21: #define LGMRES_DELTA_DIRECTIONS 10
22: #define LGMRES_DEFAULT_MAXK 30
23: #define LGMRES_DEFAULT_AUGDIM 2 /*default number of augmentation vectors */
24: static int LGMRESGetNewVectors(KSP,int);
25: static int LGMRESUpdateHessenberg(KSP,int,PetscTruth,PetscReal *);
26: static int BuildLgmresSoln(PetscScalar*,Vec,Vec,KSP,int);
28: /*
29: KSPSetUp_LGMRES - Sets up the workspace needed by lgmres.
31: This is called once, usually automatically by SLESSolve() or SLESSetUp(),
32: but can be called directly by KSPSetUp().
34: */
37: int KSPSetUp_LGMRES(KSP ksp)
38: {
39: unsigned int size,hh,hes,rs,cc;
40: int ierr,max_k,k, aug_dim;
41: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
44: if (ksp->pc_side == PC_SYMMETRIC) {
45: SETERRQ(2,"no symmetric preconditioning for KSPLGMRES");
46: }
50: max_k = lgmres->max_k;
51: aug_dim = lgmres->aug_dim;
52: hh = (max_k + 2) * (max_k + 1);
53: hes = (max_k + 1) * (max_k + 1);
54: rs = (max_k + 2);
55: cc = (max_k + 1); /* SS and CC are the same size */
56: size = (hh + hes + rs + 2*cc) * sizeof(PetscScalar);
58: /* Allocate space and set pointers to beginning */
59: PetscMalloc(size,&lgmres->hh_origin);
60: PetscMemzero(lgmres->hh_origin,size);
61: PetscLogObjectMemory(ksp,size); /* HH - modified (by plane
62: rotations) hessenburg */
63: lgmres->hes_origin = lgmres->hh_origin + hh; /* HES - unmodified hessenburg */
64: lgmres->rs_origin = lgmres->hes_origin + hes; /* RS - the right-hand-side of the
65: Hessenberg system */
66: lgmres->cc_origin = lgmres->rs_origin + rs; /* CC - cosines for rotations */
67: lgmres->ss_origin = lgmres->cc_origin + cc; /* SS - sines for rotations */
69: if (ksp->calc_sings) {
70: /* Allocate workspace to hold Hessenberg matrix needed by Eispack */
71: size = (max_k + 3)*(max_k + 9)*sizeof(PetscScalar);
72: PetscMalloc(size,&lgmres->Rsvd);
73: PetscMalloc(5*(max_k+2)*sizeof(PetscReal),&lgmres->Dsvd);
74: PetscLogObjectMemory(ksp,size+5*(max_k+2)*sizeof(PetscReal));
75: }
77: /* Allocate array to hold pointers to user vectors. Note that we need
78: we need it+1 vectors, and it <= max_k) - vec_offset indicates some initial work vectors*/
79: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void *),&lgmres->vecs);
80: lgmres->vecs_allocated = VEC_OFFSET + 2 + max_k;
81: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void *),&lgmres->user_work);
82: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(int),&lgmres->mwork_alloc);
83: PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void *)+sizeof(int)));
85: /* LGMRES_MOD: need array of pointers to augvecs*/
86: PetscMalloc((2 * aug_dim + AUG_OFFSET)*sizeof(void *),&lgmres->augvecs);
87: lgmres->aug_vecs_allocated = 2 *aug_dim + AUG_OFFSET;
88: PetscMalloc((2* aug_dim + AUG_OFFSET)*sizeof(void *),&lgmres->augvecs_user_work);
89: PetscMalloc(aug_dim*sizeof(int),&lgmres->aug_order);
90: PetscLogObjectMemory(ksp,(aug_dim)*(4*sizeof(void *) + sizeof(int)) + AUG_OFFSET*2*sizeof(void *) );
92:
93: /* if q_preallocate = 0 then only allocate one "chunk" of space (for
94: 5 vectors) - additional will then be allocated from LGMREScycle()
95: as needed. Otherwise, allocate all of the space that could be needed */
96: if (lgmres->q_preallocate) {
97: lgmres->vv_allocated = VEC_OFFSET + 2 + max_k;
98: VecDuplicateVecs(VEC_RHS,lgmres->vv_allocated,&lgmres->user_work[0]);
99: PetscLogObjectParents(ksp,lgmres->vv_allocated,lgmres->user_work[0]);
100: lgmres->mwork_alloc[0] = lgmres->vv_allocated;
101: lgmres->nwork_alloc = 1;
102: for (k=0; k<lgmres->vv_allocated; k++) {
103: lgmres->vecs[k] = lgmres->user_work[0][k];
104: }
105: } else {
106: lgmres->vv_allocated = 5;
107: VecDuplicateVecs(ksp->vec_rhs,5,&lgmres->user_work[0]);
108: PetscLogObjectParents(ksp,5,lgmres->user_work[0]);
109: lgmres->mwork_alloc[0] = 5;
110: lgmres->nwork_alloc = 1;
111: for (k=0; k<lgmres->vv_allocated; k++) {
112: lgmres->vecs[k] = lgmres->user_work[0][k];
113: }
114: }
115: /* LGMRES_MOD - for now we will preallocate the augvecs - because aug_dim << restart
116: ... also keep in mind that we need to keep augvecs from cycle to cycle*/
117: lgmres->aug_vv_allocated = 2* aug_dim + AUG_OFFSET;
118: lgmres->augwork_alloc = 2* aug_dim + AUG_OFFSET;
119: VecDuplicateVecs(VEC_RHS,lgmres->aug_vv_allocated,&lgmres->augvecs_user_work[0]);
120: PetscLogObjectParents(ksp,lgmres->aug_vv_allocated,lgmres->augvecs_user_work[0]);
121: for (k=0; k<lgmres->aug_vv_allocated; k++) {
122: lgmres->augvecs[k] = lgmres->augvecs_user_work[0][k];
123: }
125: return(0);
126: }
129: /*
131: LGMRESCycle - Run lgmres, possibly with restart. Return residual
132: history if requested.
134: input parameters:
135: . lgmres - structure containing parameters and work areas
137: output parameters:
138: . nres - residuals (from preconditioned system) at each step.
139: If restarting, consider passing nres+it. If null,
140: ignored
141: . itcount - number of iterations used. nres[0] to nres[itcount]
142: are defined. If null, ignored. If null, ignored.
143: . converged - 0 if not converged
145:
146: Notes:
147: On entry, the value in vector VEC_VV(0) should be
148: the initial residual.
151: */
154: int LGMREScycle(int *itcount,KSP ksp)
155: {
157: KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
158: PetscReal res_norm, res;
159: PetscReal hapbnd, tt;
160: PetscScalar zero = 0.0;
161: PetscScalar tmp;
162: PetscTruth hapend = PETSC_FALSE; /* indicates happy breakdown ending */
163: int ierr;
164: int loc_it; /* local count of # of dir. in Krylov space */
165: int max_k = lgmres->max_k; /* max approx space size */
166: int max_it = ksp->max_it; /* max # of overall iterations for the method */
167: /* LGMRES_MOD - new variables*/
168: int aug_dim = lgmres->aug_dim;
169: int spot = 0;
170: int order = 0;
171: int it_arnoldi; /* number of arnoldi steps to take */
172: int it_total; /* total number of its to take (=approx space size)*/
173: int ii, jj;
174: PetscReal tmp_norm;
175: PetscScalar inv_tmp_norm;
176: PetscScalar *avec;
180: /* Number of pseudo iterations since last restart is the number
181: of prestart directions */
182: loc_it = 0;
184: /* LGMRES_MOD: determine number of arnoldi steps to take */
185: /* if approx_constant then we keep the space the same size even if
186: we don't have the full number of aug vectors yet*/
187: if (lgmres->approx_constant) {
188: it_arnoldi = max_k - lgmres->aug_ct;
189: } else {
190: it_arnoldi = max_k - aug_dim;
191: }
193: it_total = it_arnoldi + lgmres->aug_ct;
196: /* initial residual is in VEC_VV(0) - compute its norm*/
197: VecNorm(VEC_VV(0),NORM_2,&res_norm);
198: res = res_norm;
199:
200: /* first entry in right-hand-side of hessenberg system is just
201: the initial residual norm */
202: *GRS(0) = res_norm;
204: /* check for the convergence */
205: if (!res) {
206: if (itcount) *itcount = 0;
207: ksp->reason = KSP_CONVERGED_ATOL;
208: PetscLogInfo(ksp,"GMRESCycle: Converged due to zero residual norm on entry\n");
209: return(0);
210: }
212: /* scale VEC_VV (the initial residual) */
213: tmp = 1.0/res_norm; VecScale(&tmp,VEC_VV(0));
215: /* FYI: AMS calls are for memory snooper */
216: PetscObjectTakeAccess(ksp);
217: ksp->rnorm = res;
218: PetscObjectGrantAccess(ksp);
221: /* note: (lgmres->it) is always set one less than (loc_it) It is used in
222: KSPBUILDSolution_LGMRES, where it is passed to BuildLgmresSoln.
223: Note that when BuildLgmresSoln is called from this function,
224: (loc_it -1) is passed, so the two are equivalent */
225: lgmres->it = (loc_it - 1);
227:
228: /* MAIN ITERATION LOOP BEGINNING*/
231: /* keep iterating until we have converged OR generated the max number
232: of directions OR reached the max number of iterations for the method */
233: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
234:
235: while (!ksp->reason && loc_it < it_total && ksp->its < max_it) { /* LGMRES_MOD: changed to it_total */
236: KSPLogResidualHistory(ksp,res);
237: lgmres->it = (loc_it - 1);
238: KSPMonitor(ksp,ksp->its,res);
240: /* see if more space is needed for work vectors */
241: if (lgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
242: LGMRESGetNewVectors(ksp,loc_it+1);
243: /* (loc_it+1) is passed in as number of the first vector that should
244: be allocated */
245: }
247: /*LGMRES_MOD: decide whether this is an arnoldi step or an aug step */
248: if (loc_it < it_arnoldi) { /* arnoldi */
249: KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,VEC_VV(loc_it),VEC_VV(1+loc_it),VEC_TEMP_MATOP);
250: } else { /*aug step */
251: order = loc_it - it_arnoldi + 1; /* which aug step */
252: for (ii=0; ii<aug_dim; ii++) {
253: if (lgmres->aug_order[ii] == order) {
254: spot = ii;
255: break; /* must have this because there will be duplicates before aug_ct = aug_dim */
256: }
257: }
259: VecCopy(A_AUGVEC(spot), VEC_VV(1+loc_it));
260: /*note: an alternate implementation choice would be to only save the AUGVECS and
261: not A_AUGVEC and then apply the PC here to the augvec */
262: }
264: /* update hessenberg matrix and do Gram-Schmidt - new direction is in
265: VEC_VV(1+loc_it)*/
266: (*lgmres->orthog)(ksp,loc_it);
268: /* new entry in hessenburg is the 2-norm of our new direction */
269: VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);
270: *HH(loc_it+1,loc_it) = tt;
271: *HES(loc_it+1,loc_it) = tt;
274: /* check for the happy breakdown */
275: hapbnd = PetscAbsScalar(tt / *GRS(loc_it));/* GRS(loc_it) contains the res_norm from the last iteration */
276: if (hapbnd > lgmres->haptol) hapbnd = lgmres->haptol;
277: if (tt > hapbnd) {
278: tmp = 1.0/tt;
279: VecScale(&tmp,VEC_VV(loc_it+1)); /* scale new direction by its norm */
280: } else {
281: PetscLogInfo(ksp,"Detected happy breakdown, current hapbnd = %g tt = %g\n",hapbnd,tt);
282: hapend = PETSC_TRUE;
283: }
285: /* Now apply rotations to new col of hessenberg (and right side of system),
286: calculate new rotation, and get new residual norm at the same time*/
287: LGMRESUpdateHessenberg(ksp,loc_it,hapend,&res);
288: loc_it++;
289: lgmres->it = (loc_it-1); /* Add this here in case it has converged */
290:
291: PetscObjectTakeAccess(ksp);
292: ksp->its++;
293: ksp->rnorm = res;
294: PetscObjectGrantAccess(ksp);
296: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
298: /* Catch error in happy breakdown and signal convergence and break from loop */
299: if (hapend) {
300: if (!ksp->reason) {
301: SETERRQ1(0,"You reached the happy break down,but convergence was not indicated. Residual norm = %g",res);
302: }
303: break;
304: }
305: }
306: /* END OF ITERATION LOOP */
308: KSPLogResidualHistory(ksp,res);
310: /* Monitor if we know that we will not return for a restart */
311: if (ksp->reason || ksp->its >= max_it) {
312: KSPMonitor(ksp, ksp->its, res);
313: }
315: if (itcount) *itcount = loc_it;
317: /*
318: Down here we have to solve for the "best" coefficients of the Krylov
319: columns, add the solution values together, and possibly unwind the
320: preconditioning from the solution
321: */
322:
323: /* Form the solution (or the solution so far) */
324: /* Note: must pass in (loc_it-1) for iteration count so that BuildLgmresSoln
325: properly navigates */
327: BuildLgmresSoln(GRS(0),VEC_SOLN,VEC_SOLN,ksp,loc_it-1);
330: /* LGMRES_MOD collect aug vector and A*augvector for future restarts -
331: only if we will be restarting (i.e. this cycle performed it_total
332: iterations) */
333: if (!ksp->reason && ksp->its < max_it && aug_dim > 0) {
335: /*AUG_TEMP contains the new augmentation vector (assigned in BuildLgmresSoln) */
336: if (lgmres->aug_ct == 0) {
337: spot = 0;
338: lgmres->aug_ct++;
339: } else if (lgmres->aug_ct < aug_dim) {
340: spot = lgmres->aug_ct;
341: lgmres->aug_ct++;
342: } else { /* truncate */
343: for (ii=0; ii<aug_dim; ii++) {
344: if (lgmres->aug_order[ii] == aug_dim) {
345: spot = ii;
346: }
347: }
348: }
350:
352: VecCopy(AUG_TEMP, AUGVEC(spot));
353: /*need to normalize */
354: VecNorm(AUGVEC(spot), NORM_2, &tmp_norm);
355: inv_tmp_norm = 1.0/tmp_norm;
356: VecScale(&inv_tmp_norm, AUGVEC(spot));
358: /*set new aug vector to order 1 - move all others back one */
359: for (ii=0; ii < aug_dim; ii++) {
360: AUG_ORDER(ii)++;
361: }
362: AUG_ORDER(spot) = 1;
364: /*now add the A*aug vector to A_AUGVEC(spot) - this is independ. of preconditioning type*/
365: /* want V*H*y - y is in GRS, V is in VEC_VV and H is in HES */
367:
368: /* first do H+*y */
369: VecSet(&zero, AUG_TEMP);
370: VecGetArray(AUG_TEMP, &avec);
371: for (ii=0; ii < it_total + 1; ii++) {
372: for (jj=0; jj <= ii+1; jj++) {
373: avec[jj] += *HES(jj ,ii) * *GRS(ii);
374: }
375: }
377: /*now multiply result by V+ */
378: VecSet(&zero, VEC_TEMP);
379: VecMAXPY(it_total+1, avec, VEC_TEMP, &VEC_VV(0)); /*answer is in VEC_TEMP*/
380: VecRestoreArray(AUG_TEMP, &avec);
381:
382: /*copy answer to aug location and scale*/
383: VecCopy(VEC_TEMP, A_AUGVEC(spot));
384: VecScale(&inv_tmp_norm, A_AUGVEC(spot));
387: }
388: return(0);
389: }
391: /*
392: KSPSolve_LGMRES - This routine applies the LGMRES method.
395: Input Parameter:
396: . ksp - the Krylov space object that was set to use lgmres
398: Output Parameter:
399: . outits - number of iterations used
401: */
405: int KSPSolve_LGMRES(KSP ksp,int *outits)
406: {
407: int ierr;
408: int cycle_its; /* iterations done in a call to LGMREScycle */
409: int itcount; /* running total of iterations, incl. those in restarts */
410: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
411: PetscTruth guess_zero = ksp->guess_zero;
412: /*LGMRES_MOD variable */
413: int ii;
416: if (ksp->calc_sings && !lgmres->Rsvd) {
417: SETERRQ(1,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
418: }
419: PetscObjectTakeAccess(ksp);
420: ksp->its = 0;
421: PetscObjectGrantAccess(ksp);
423: /* initialize */
424: itcount = 0;
425: ksp->reason = KSP_CONVERGED_ITERATING;
426: /*LGMRES_MOD*/
427: for (ii=0; ii<lgmres->aug_dim; ii++) {
428: lgmres->aug_order[ii] = 0;
429: }
431: while (!ksp->reason) {
432: /* calc residual - puts in VEC_VV(0) */
433: KSPInitialResidual(ksp,VEC_SOLN,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),VEC_RHS);
434: LGMREScycle(&cycle_its,ksp);
435: itcount += cycle_its;
436: if (itcount >= ksp->max_it) {
437: ksp->reason = KSP_DIVERGED_ITS;
438: break;
439: }
440: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
441: }
442: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
443: if (outits) *outits = itcount;
444: return(0);
445: }
447: /*
449: KSPDestroy_LGMRES - Frees all memory space used by the Krylov method.
451: */
454: int KSPDestroy_LGMRES(KSP ksp)
455: {
456: KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;
457: int i,ierr;
460: /* Free the Hessenberg matrices */
461: if (lgmres->hh_origin) {PetscFree(lgmres->hh_origin);}
463: /* Free pointers to user variables */
464: if (lgmres->vecs) {PetscFree(lgmres->vecs);}
466: /*LGMRES_MOD - free pointers for extra vectors */
467: if (lgmres->augvecs) {PetscFree(lgmres->augvecs);}
469: /* free work vectors */
470: for (i=0; i < lgmres->nwork_alloc; i++) {
471: VecDestroyVecs(lgmres->user_work[i],lgmres->mwork_alloc[i]);
472: }
473: if (lgmres->user_work) {PetscFree(lgmres->user_work);}
475: /*LGMRES_MOD - free aug work vectors also */
476: /*this was all allocated as one "chunk" */
477: VecDestroyVecs(lgmres->augvecs_user_work[0],lgmres->augwork_alloc);
478: if (lgmres->augvecs_user_work) {PetscFree(lgmres->augvecs_user_work);}
479: if (lgmres->aug_order) {PetscFree(lgmres->aug_order);}
481: if (lgmres->mwork_alloc) {PetscFree(lgmres->mwork_alloc);}
482: if (lgmres->nrs) {PetscFree(lgmres->nrs);}
483: if (lgmres->sol_temp) {VecDestroy(lgmres->sol_temp);}
484: if (lgmres->Rsvd) {PetscFree(lgmres->Rsvd);}
485: if (lgmres->Dsvd) {PetscFree(lgmres->Dsvd);}
486: PetscFree(lgmres);
487: return(0);
488: }
490: /*
491: BuildLgmresSoln - create the solution from the starting vector and the
492: current iterates.
494: Input parameters:
495: nrs - work area of size it + 1.
496: vguess - index of initial guess
497: vdest - index of result. Note that vguess may == vdest (replace
498: guess with the solution).
499: it - HH upper triangular part is a block of size (it+1) x (it+1)
501: This is an internal routine that knows about the LGMRES internals.
502: */
505: static int BuildLgmresSoln(PetscScalar* nrs,Vec vguess,Vec vdest,KSP ksp,int it)
506: {
507: PetscScalar tt,zero = 0.0,one = 1.0;
508: int ierr,ii,k,j;
509: KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
510: /*LGMRES_MOD */
511: int it_arnoldi, it_aug;
512: int jj, spot = 0;
515: /* Solve for solution vector that minimizes the residual */
517: /* If it is < 0, no lgmres steps have been performed */
518: if (it < 0) {
519: if (vdest != vguess) {
520: VecCopy(vguess,vdest);
521: }
522: return(0);
523: }
525: /* so (it+1) lgmres steps HAVE been performed */
527: /* LGMRES_MOD - determine if we need to use augvecs for the soln - do not assume that
528: this is called after the total its allowed for an approx space */
529: if (lgmres->approx_constant) {
530: it_arnoldi = lgmres->max_k - lgmres->aug_ct;
531: } else {
532: it_arnoldi = lgmres->max_k - lgmres->aug_dim;
533: }
534: if (it_arnoldi >= it +1) {
535: it_aug = 0;
536: it_arnoldi = it+1;
537: } else {
538: it_aug = (it + 1) - it_arnoldi;
539: }
541: /* now it_arnoldi indicates the number of matvecs that took place */
542: lgmres->matvecs += it_arnoldi;
544:
545: /* solve the upper triangular system - GRS is the right side and HH is
546: the upper triangular matrix - put soln in nrs */
547: if (*HH(it,it) == 0.0) SETERRQ2(1,"HH(it,it) is identically zero; it = %d GRS(it) = %g",it,PetscAbsScalar(*GRS(it)));
548: if (*HH(it,it) != 0.0) {
549: nrs[it] = *GRS(it) / *HH(it,it);
550: } else {
551: nrs[it] = 0.0;
552: }
554: for (ii=1; ii<=it; ii++) {
555: k = it - ii;
556: tt = *GRS(k);
557: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
558: nrs[k] = tt / *HH(k,k);
559: }
561: /* Accumulate the correction to the soln of the preconditioned prob. in VEC_TEMP */
562: VecSet(&zero,VEC_TEMP); /* set VEC_TEMP components to 0 */
564: /*LGMRES_MOD - if augmenting has happened we need to form the solution
565: using the augvecs */
566: if (it_aug == 0) { /* all its are from arnoldi */
567: VecMAXPY(it+1,nrs,VEC_TEMP,&VEC_VV(0));
568: } else { /*use aug vecs */
569: /*first do regular krylov directions */
570: VecMAXPY(it_arnoldi,nrs,VEC_TEMP,&VEC_VV(0));
571: /*now add augmented portions - add contribution of aug vectors one at a time*/
574: for (ii=0; ii<it_aug; ii++) {
575: for (jj=0; jj<lgmres->aug_dim; jj++) {
576: if (lgmres->aug_order[jj] == (ii+1)) {
577: spot = jj;
578: break; /* must have this because there will be duplicates before aug_ct = aug_dim */
579: }
580: }
581: VecAXPY(&nrs[it_arnoldi+ii],AUGVEC(spot),VEC_TEMP);
582: }
583: }
584: /* now VEC_TEMP is what we want to keep for augmenting purposes - grab before the
585: preconditioner is "unwound" from right-precondtioning*/
586: VecCopy(VEC_TEMP, AUG_TEMP);
588: KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
590: /* add solution to previous solution */
591: /* put updated solution into vdest.*/
592: if (vdest != vguess) {
593: VecCopy(VEC_TEMP,vdest);
594: }
595: VecAXPY(&one,VEC_TEMP,vdest);
597: return(0);
598: }
600: /*
602: LGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
603: Return new residual.
605: input parameters:
607: . ksp - Krylov space object
608: . it - plane rotations are applied to the (it+1)th column of the
609: modified hessenberg (i.e. HH(:,it))
610: . hapend - PETSC_FALSE not happy breakdown ending.
612: output parameters:
613: . res - the new residual
614:
615: */
618: static int LGMRESUpdateHessenberg(KSP ksp,int it,PetscTruth hapend,PetscReal *res)
619: {
620: PetscScalar *hh,*cc,*ss,tt;
621: int j;
622: KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
625: hh = HH(0,it); /* pointer to beginning of column to update - so
626: incrementing hh "steps down" the (it+1)th col of HH*/
627: cc = CC(0); /* beginning of cosine rotations */
628: ss = SS(0); /* beginning of sine rotations */
630: /* Apply all the previously computed plane rotations to the new column
631: of the Hessenberg matrix */
632: /* Note: this uses the rotation [conj(c) s ; -s c], c= cos(theta), s= sin(theta) */
634: for (j=1; j<=it; j++) {
635: tt = *hh;
636: #if defined(PETSC_USE_COMPLEX)
637: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
638: #else
639: *hh = *cc * tt + *ss * *(hh+1);
640: #endif
641: hh++;
642: *hh = *cc++ * *hh - (*ss++ * tt);
643: /* hh, cc, and ss have all been incremented one by end of loop */
644: }
646: /*
647: compute the new plane rotation, and apply it to:
648: 1) the right-hand-side of the Hessenberg system (GRS)
649: note: it affects GRS(it) and GRS(it+1)
650: 2) the new column of the Hessenberg matrix
651: note: it affects HH(it,it) which is currently pointed to
652: by hh and HH(it+1, it) (*(hh+1))
653: thus obtaining the updated value of the residual...
654: */
656: /* compute new plane rotation */
658: if (!hapend) {
659: #if defined(PETSC_USE_COMPLEX)
660: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
661: #else
662: tt = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
663: #endif
664: if (tt == 0.0) {SETERRQ(PETSC_ERR_KSP_BRKDWN,"Your matrix or preconditioner is the null operator");}
665: *cc = *hh / tt; /* new cosine value */
666: *ss = *(hh+1) / tt; /* new sine value */
668: /* apply to 1) and 2) */
669: *GRS(it+1) = - (*ss * *GRS(it));
670: #if defined(PETSC_USE_COMPLEX)
671: *GRS(it) = PetscConj(*cc) * *GRS(it);
672: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
673: #else
674: *GRS(it) = *cc * *GRS(it);
675: *hh = *cc * *hh + *ss * *(hh+1);
676: #endif
678: /* residual is the last element (it+1) of right-hand side! */
679: *res = PetscAbsScalar(*GRS(it+1));
681: } else { /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
682: another rotation matrix (so RH doesn't change). The new residual is
683: always the new sine term times the residual from last time (GRS(it)),
684: but now the new sine rotation would be zero...so the residual should
685: be zero...so we will multiply "zero" by the last residual. This might
686: not be exactly what we want to do here -could just return "zero". */
687:
688: *res = 0.0;
689: }
690: return(0);
691: }
693: /*
695: LGMRESGetNewVectors - This routine allocates more work vectors, starting from
696: VEC_VV(it)
697:
698: */
701: static int LGMRESGetNewVectors(KSP ksp,int it)
702: {
703: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
704: int nwork = lgmres->nwork_alloc; /* number of work vector chunks allocated */
705: int nalloc; /* number to allocate */
706: int k,ierr;
707:
709: nalloc = lgmres->delta_allocate; /* number of vectors to allocate
710: in a single chunk */
712: /* Adjust the number to allocate to make sure that we don't exceed the
713: number of available slots (lgmres->vecs_allocated)*/
714: if (it + VEC_OFFSET + nalloc >= lgmres->vecs_allocated){
715: nalloc = lgmres->vecs_allocated - it - VEC_OFFSET;
716: }
717: if (!nalloc) return(0);
719: lgmres->vv_allocated += nalloc; /* vv_allocated is the number of vectors allocated */
721: /* work vectors */
722: VecDuplicateVecs(ksp->vec_rhs,nalloc,&lgmres->user_work[nwork]);
723: PetscLogObjectParents(ksp,nalloc,lgmres->user_work[nwork]);
724: /* specify size of chunk allocated */
725: lgmres->mwork_alloc[nwork] = nalloc;
727: for (k=0; k < nalloc; k++) {
728: lgmres->vecs[it+VEC_OFFSET+k] = lgmres->user_work[nwork][k];
729: }
730:
732: /* LGMRES_MOD - for now we are preallocating the augmentation vectors */
733:
735: /* increment the number of work vector chunks */
736: lgmres->nwork_alloc++;
737: return(0);
738: }
740: /*
742: KSPBuildSolution_LGMRES
744: Input Parameter:
745: . ksp - the Krylov space object
746: . ptr-
748: Output Parameter:
749: . result - the solution
751: Note: this calls BuildLgmresSoln - the same function that LGMREScycle
752: calls directly.
754: */
757: int KSPBuildSolution_LGMRES(KSP ksp,Vec ptr,Vec *result)
758: {
759: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
760: int ierr;
763: if (!ptr) {
764: if (!lgmres->sol_temp) {
765: VecDuplicate(ksp->vec_sol,&lgmres->sol_temp);
766: PetscLogObjectParent(ksp,lgmres->sol_temp);
767: }
768: ptr = lgmres->sol_temp;
769: }
770: if (!lgmres->nrs) {
771: /* allocate the work area */
772: PetscMalloc(lgmres->max_k*sizeof(PetscScalar),&lgmres->nrs);
773: PetscLogObjectMemory(ksp,lgmres->max_k*sizeof(PetscScalar));
774: }
775:
776: BuildLgmresSoln(lgmres->nrs,VEC_SOLN,ptr,ksp,lgmres->it);
777: *result = ptr;
778:
779: return(0);
780: }
782: /*
784: KSPView_LGMRES -Prints information about the current Krylov method
785: being used.
787: */
790: int KSPView_LGMRES(KSP ksp,PetscViewer viewer)
791: {
792: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
793: char *cstr;
794: int ierr;
795: PetscTruth isascii,isstring;
798: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
799: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
800: if (lgmres->orthog == KSPGMRESUnmodifiedGramSchmidtOrthogonalization) {
801: cstr = "Unmodified Gram-Schmidt Orthogonalization";
802: } else if (lgmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
803: cstr = "Modified Gram-Schmidt Orthogonalization";
804: } else if (lgmres->orthog == KSPGMRESIROrthogonalization) {
805: cstr = "Unmodified Gram-Schmidt + 1 step Iterative Refinement Orthogonalization";
806: } else {
807: cstr = "unknown orthogonalization";
808: }
809: if (isascii) {
810: PetscViewerASCIIPrintf(viewer," LGMRES: restart=%d, using %s\n",lgmres->max_k,cstr);
811: /*LGMRES_MOD */
812: PetscViewerASCIIPrintf(viewer," LGMRES: aug. dimension=%d\n",lgmres->aug_dim);
813: if (lgmres->approx_constant) {
814: PetscViewerASCIIPrintf(viewer," LGMRES: approx. space size was kept constant.\n");
815: }
816: PetscViewerASCIIPrintf(viewer," LGMRES: number of matvecs=%d\n",lgmres->matvecs);
818: PetscViewerASCIIPrintf(viewer," LGMRES: happy breakdown tolerance %g\n",lgmres->haptol);
819: } else if (isstring) {
820: PetscViewerStringSPrintf(viewer,"%s restart %d",cstr,lgmres->max_k);
821: } else {
822: SETERRQ1(1,"Viewer type %s not supported for KSP LGMRES",((PetscObject)viewer)->type_name);
823: }
824: return(0);
827: }
831: int KSPSetFromOptions_LGMRES(KSP ksp)
832: {
833: int ierr, restart, aug;
834: PetscReal haptol;
835: KSP_LGMRES *lgmres = (KSP_LGMRES*) ksp->data;
836: PetscTruth flg;
839: PetscOptionsHead("KSP LGMRES Options");
841:
843: PetscOptionsInt("-ksp_gmres_restart","For LGMRES, this is the maximum size of the approximation space","KSPGMRESSetRestart",lgmres->max_k,&restart,&flg);
844: if (flg) { KSPGMRESSetRestart(ksp,restart); }
845: PetscOptionsReal("-ksp_gmres_haptol","Tolerance for declaring exact convergence (happy ending)","KSPGMRESSetHapTol",lgmres->haptol,&haptol,&flg);
846: if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
847: PetscOptionsName("-ksp_gmres_preallocate","Preallocate all Krylov vectors","KSPGMRESSetPreAllocateVectors",&flg);
848: if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
849: PetscOptionsLogicalGroupBegin("-ksp_gmres_unmodifiedgramschmidt","Use classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
850: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESUnmodifiedGramSchmidtOrthogonalization);}
851: PetscOptionsLogicalGroup("-ksp_gmres_modifiedgramschmidt","Use modified Gram-Schmidt (slow but more stable)","KSPGMRESSetOrthogonalization",&flg);
852: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
853: PetscOptionsLogicalGroupEnd("-ksp_gmres_irorthog","Use classical Gram-Schmidt with iterative refinement","KSPGMRESSetOrthogonalization",&flg);
854: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESIROrthogonalization);}
855: PetscOptionsName("-ksp_gmres_krylov_monitor","Graphically plot the Krylov directions","KSPSetMonitor",&flg);
856: if (flg) {
857: PetscViewers viewers;
858: PetscViewersCreate(ksp->comm,&viewers);
859: KSPSetMonitor(ksp,KSPGMRESKrylovMonitor,viewers,(int (*)(void*))PetscViewersDestroy);
860: }
862: /* LGMRES_MOD - specify number of augmented vectors and whether the space should be a constant size*/
863: PetscOptionsName("-ksp_lgmres_constant","Use constant approx. space size","KSPGMRESSetConstant",&flg);
864: /*if (flg) {KSPGMRESSetConstant(ksp);}*/ /*<--doesn't like this */
865: if (flg) { lgmres->approx_constant = 1; } /* in favor of this line....*/
867: PetscOptionsInt("-ksp_lgmres_augment","Number of error approximations to augment the Krylov space with","KSPLGMRESSetAugDim",lgmres->aug_dim,&aug,&flg);
868: if (flg) { KSPLGMRESSetAugDim(ksp,aug); }
872: PetscOptionsTail();
873: return(0);
874: }
877: EXTERN int KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
878: EXTERN int KSPComputeEigenvalues_GMRES(KSP,int,PetscReal *,PetscReal *,int *);
880: /*functions for extra lgmres options here*/
881: EXTERN_C_BEGIN
884: int KSPLGMRESSetConstant_LGMRES(KSP ksp)
885: {
886: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
888: lgmres->approx_constant = 1;
889:
890: return(0);
891: }
892: EXTERN_C_END
894: EXTERN_C_BEGIN
897: int KSPLGMRESSetAugDim_LGMRES(KSP ksp,int aug_dim)
898: {
899: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
903: if (aug_dim < 0) SETERRQ(1,"Augmentation dimension must be positive");
904: if (aug_dim > (lgmres->max_k -1)) SETERRQ(1,"Augmentation dimension must be <= (restart size-1)");
906: lgmres->aug_dim = aug_dim;
908: return(0);
909: }
910: EXTERN_C_END
913: /* end new lgmres functions */
916: /* use these options from gmres */
917: EXTERN_C_BEGIN
918: EXTERN int KSPGMRESSetHapTol_GMRES(KSP,double);
919: EXTERN int KSPGMRESSetPreAllocateVectors_GMRES(KSP);
920: EXTERN int KSPGMRESSetRestart_GMRES(KSP,int);
921: EXTERN int KSPGMRESSetOrthogonalization_GMRES(KSP,int (*)(KSP,int));
922: EXTERN_C_END
924: EXTERN_C_BEGIN
927: int KSPCreate_LGMRES(KSP ksp)
928: {
929: KSP_LGMRES *lgmres;
930: int ierr;
933: PetscNew(KSP_LGMRES,&lgmres);
934: PetscMemzero(lgmres,sizeof(KSP_LGMRES));
935: PetscLogObjectMemory(ksp,sizeof(KSP_LGMRES));
936: ksp->data = (void*)lgmres;
937: ksp->ops->buildsolution = KSPBuildSolution_LGMRES;
939: ksp->ops->setup = KSPSetUp_LGMRES;
940: ksp->ops->solve = KSPSolve_LGMRES;
941: ksp->ops->destroy = KSPDestroy_LGMRES;
942: ksp->ops->view = KSPView_LGMRES;
943: ksp->ops->setfromoptions = KSPSetFromOptions_LGMRES;
944: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
945: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
947: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
948: "KSPGMRESSetPreAllocateVectors_GMRES",
949: KSPGMRESSetPreAllocateVectors_GMRES);
950: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
951: "KSPGMRESSetOrthogonalization_GMRES",
952: KSPGMRESSetOrthogonalization_GMRES);
953: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
954: "KSPGMRESSetRestart_GMRES",
955: KSPGMRESSetRestart_GMRES);
956: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
957: "KSPGMRESSetHapTol_GMRES",
958: KSPGMRESSetHapTol_GMRES);
960: /*LGMRES_MOD add extra functions here - like the one to set num of aug vectors */
961: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetConstant_C",
962: "KSPLGMRESSetConstant_LGMRES",
963: KSPLGMRESSetConstant_LGMRES);
965: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetAugDim_C",
966: "KSPLGMRESSetAugDim_LGMRES",
967: KSPLGMRESSetAugDim_LGMRES);
968:
970: /*defaults */
971: lgmres->haptol = 1.0e-30;
972: lgmres->q_preallocate = 0;
973: lgmres->delta_allocate = LGMRES_DELTA_DIRECTIONS;
974: lgmres->orthog = KSPGMRESUnmodifiedGramSchmidtOrthogonalization;
975: lgmres->nrs = 0;
976: lgmres->sol_temp = 0;
977: lgmres->max_k = LGMRES_DEFAULT_MAXK;
978: lgmres->Rsvd = 0;
980: /*LGMRES_MOD - new defaults */
981: lgmres->aug_dim = LGMRES_DEFAULT_AUGDIM;
982: lgmres->aug_ct = 0; /* start with no aug vectors */
983: lgmres->approx_constant = 0;
984: lgmres->matvecs = 0;
986: return(0);
987: }
988: EXTERN_C_END