Actual source code: lgmres.c
1: #define PETSCKSP_DLL
3: #include src/ksp/ksp/impls/gmres/lgmres/lgmresp.h
5: #define LGMRES_DELTA_DIRECTIONS 10
6: #define LGMRES_DEFAULT_MAXK 30
7: #define LGMRES_DEFAULT_AUGDIM 2 /*default number of augmentation vectors */
8: static PetscErrorCode LGMRESGetNewVectors(KSP,PetscInt);
9: static PetscErrorCode LGMRESUpdateHessenberg(KSP,PetscInt,PetscTruth,PetscReal *);
10: static PetscErrorCode BuildLgmresSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);
14: PetscErrorCode KSPLGMRESSetAugDim(KSP ksp, PetscInt dim)
15: {
19: PetscTryMethod((ksp),KSPLGMRESSetAugDim_C,(KSP,PetscInt),(ksp,dim));
20: return(0);
21: }
25: PetscErrorCode KSPLGMRESSetConstant(KSP ksp)
26: {
30: PetscTryMethod((ksp),KSPLGMRESSetConstant_C,(KSP),(ksp));
31: return(0);
32: }
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: {
45: PetscInt size,hh,hes,rs,cc;
47: PetscInt max_k,k, aug_dim;
48: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
51: if (ksp->pc_side == PC_SYMMETRIC) {
52: SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPLGMRES");
53: }
54: max_k = lgmres->max_k;
55: aug_dim = lgmres->aug_dim;
56: hh = (max_k + 2) * (max_k + 1);
57: hes = (max_k + 1) * (max_k + 1);
58: rs = (max_k + 2);
59: cc = (max_k + 1); /* SS and CC are the same size */
60: size = (hh + hes + rs + 2*cc) * sizeof(PetscScalar);
62: /* Allocate space and set pointers to beginning */
63: PetscMalloc(size,&lgmres->hh_origin);
64: PetscMemzero(lgmres->hh_origin,size);
65: PetscLogObjectMemory(ksp,size); /* HH - modified (by plane rotations) hessenburg */
66: lgmres->hes_origin = lgmres->hh_origin + hh; /* HES - unmodified hessenburg */
67: lgmres->rs_origin = lgmres->hes_origin + hes; /* RS - the right-hand-side of the
68: Hessenberg system */
69: lgmres->cc_origin = lgmres->rs_origin + rs; /* CC - cosines for rotations */
70: lgmres->ss_origin = lgmres->cc_origin + cc; /* SS - sines for rotations */
72: if (ksp->calc_sings) {
73: /* Allocate workspace to hold Hessenberg matrix needed by Eispack */
74: size = (max_k + 3)*(max_k + 9)*sizeof(PetscScalar);
75: PetscMalloc(size,&lgmres->Rsvd);
76: PetscMalloc(5*(max_k+2)*sizeof(PetscReal),&lgmres->Dsvd);
77: PetscLogObjectMemory(ksp,size+5*(max_k+2)*sizeof(PetscReal));
78: }
80: /* Allocate array to hold pointers to user vectors. Note that we need
81: we need it+1 vectors, and it <= max_k) - vec_offset indicates some initial work vectors*/
82: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&lgmres->vecs);
83: lgmres->vecs_allocated = VEC_OFFSET + 2 + max_k;
84: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&lgmres->user_work);
85: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(PetscInt),&lgmres->mwork_alloc);
86: PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void*)+sizeof(PetscInt)));
88: /* LGMRES_MOD: need array of pointers to augvecs*/
89: PetscMalloc((2 * aug_dim + AUG_OFFSET)*sizeof(void*),&lgmres->augvecs);
90: lgmres->aug_vecs_allocated = 2 *aug_dim + AUG_OFFSET;
91: PetscMalloc((2* aug_dim + AUG_OFFSET)*sizeof(void*),&lgmres->augvecs_user_work);
92: PetscMalloc(aug_dim*sizeof(PetscInt),&lgmres->aug_order);
93: PetscLogObjectMemory(ksp,(aug_dim)*(4*sizeof(void*) + sizeof(PetscInt)) + AUG_OFFSET*2*sizeof(void*));
95:
96: /* if q_preallocate = 0 then only allocate one "chunk" of space (for
97: 5 vectors) - additional will then be allocated from LGMREScycle()
98: as needed. Otherwise, allocate all of the space that could be needed */
99: if (lgmres->q_preallocate) {
100: lgmres->vv_allocated = VEC_OFFSET + 2 + max_k;
101: KSPGetVecs(ksp,lgmres->vv_allocated,&lgmres->user_work[0],0,PETSC_NULL);
102: PetscLogObjectParents(ksp,lgmres->vv_allocated,lgmres->user_work[0]);
103: lgmres->mwork_alloc[0] = lgmres->vv_allocated;
104: lgmres->nwork_alloc = 1;
105: for (k=0; k<lgmres->vv_allocated; k++) {
106: lgmres->vecs[k] = lgmres->user_work[0][k];
107: }
108: } else {
109: lgmres->vv_allocated = 5;
110: KSPGetVecs(ksp,5,&lgmres->user_work[0],0,PETSC_NULL);
111: PetscLogObjectParents(ksp,5,lgmres->user_work[0]);
112: lgmres->mwork_alloc[0] = 5;
113: lgmres->nwork_alloc = 1;
114: for (k=0; k<lgmres->vv_allocated; k++) {
115: lgmres->vecs[k] = lgmres->user_work[0][k];
116: }
117: }
118: /* LGMRES_MOD - for now we will preallocate the augvecs - because aug_dim << restart
119: ... also keep in mind that we need to keep augvecs from cycle to cycle*/
120: lgmres->aug_vv_allocated = 2* aug_dim + AUG_OFFSET;
121: lgmres->augwork_alloc = 2* aug_dim + AUG_OFFSET;
122: KSPGetVecs(ksp,lgmres->aug_vv_allocated,&lgmres->augvecs_user_work[0],0,PETSC_NULL);
123: PetscLogObjectParents(ksp,lgmres->aug_vv_allocated,lgmres->augvecs_user_work[0]);
124: for (k=0; k<lgmres->aug_vv_allocated; k++) {
125: lgmres->augvecs[k] = lgmres->augvecs_user_work[0][k];
126: }
127: return(0);
128: }
131: /*
133: LGMRESCycle - Run lgmres, possibly with restart. Return residual
134: history if requested.
136: input parameters:
137: . lgmres - structure containing parameters and work areas
139: output parameters:
140: . nres - residuals (from preconditioned system) at each step.
141: If restarting, consider passing nres+it. If null,
142: ignored
143: . itcount - number of iterations used. nres[0] to nres[itcount]
144: are defined. If null, ignored. If null, ignored.
145: . converged - 0 if not converged
147:
148: Notes:
149: On entry, the value in vector VEC_VV(0) should be
150: the initial residual.
153: */
156: PetscErrorCode LGMREScycle(PetscInt *itcount,KSP ksp)
157: {
159: KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
160: PetscReal res_norm, res;
161: PetscReal hapbnd, tt;
162: PetscScalar tmp;
163: PetscTruth hapend = PETSC_FALSE; /* indicates happy breakdown ending */
165: PetscInt loc_it; /* local count of # of dir. in Krylov space */
166: PetscInt max_k = lgmres->max_k; /* max approx space size */
167: PetscInt max_it = ksp->max_it; /* max # of overall iterations for the method */
168: /* LGMRES_MOD - new variables*/
169: PetscInt aug_dim = lgmres->aug_dim;
170: PetscInt spot = 0;
171: PetscInt order = 0;
172: PetscInt it_arnoldi; /* number of arnoldi steps to take */
173: PetscInt it_total; /* total number of its to take (=approx space size)*/
174: PetscInt ii, jj;
175: PetscReal tmp_norm;
176: PetscScalar inv_tmp_norm;
177: 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;
195: /* initial residual is in VEC_VV(0) - compute its norm*/
196: VecNorm(VEC_VV(0),NORM_2,&res_norm);
197: res = res_norm;
198:
199: /* first entry in right-hand-side of hessenberg system is just
200: the initial residual norm */
201: *GRS(0) = res_norm;
203: /* check for the convergence */
204: if (!res) {
205: if (itcount) *itcount = 0;
206: ksp->reason = KSP_CONVERGED_ATOL;
207: PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
208: return(0);
209: }
211: /* scale VEC_VV (the initial residual) */
212: tmp = 1.0/res_norm; VecScale(VEC_VV(0),tmp);
214: /* FYI: AMS calls are for memory snooper */
215: PetscObjectTakeAccess(ksp);
216: ksp->rnorm = res;
217: PetscObjectGrantAccess(ksp);
220: /* note: (lgmres->it) is always set one less than (loc_it) It is used in
221: KSPBUILDSolution_LGMRES, where it is passed to BuildLgmresSoln.
222: Note that when BuildLgmresSoln is called from this function,
223: (loc_it -1) is passed, so the two are equivalent */
224: lgmres->it = (loc_it - 1);
226:
227: /* MAIN ITERATION LOOP BEGINNING*/
230: /* keep iterating until we have converged OR generated the max number
231: of directions OR reached the max number of iterations for the method */
232: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
233:
234: while (!ksp->reason && loc_it < it_total && ksp->its < max_it) { /* LGMRES_MOD: changed to it_total */
235: KSPLogResidualHistory(ksp,res);
236: lgmres->it = (loc_it - 1);
237: KSPMonitor(ksp,ksp->its,res);
239: /* see if more space is needed for work vectors */
240: if (lgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
241: LGMRESGetNewVectors(ksp,loc_it+1);
242: /* (loc_it+1) is passed in as number of the first vector that should
243: be allocated */
244: }
246: /*LGMRES_MOD: decide whether this is an arnoldi step or an aug step */
247: if (loc_it < it_arnoldi) { /* Arnoldi */
248: KSP_PCApplyBAorAB(ksp,VEC_VV(loc_it),VEC_VV(1+loc_it),VEC_TEMP_MATOP);
249: } else { /*aug step */
250: order = loc_it - it_arnoldi + 1; /* which aug step */
251: for (ii=0; ii<aug_dim; ii++) {
252: if (lgmres->aug_order[ii] == order) {
253: spot = ii;
254: break; /* must have this because there will be duplicates before aug_ct = aug_dim */
255: }
256: }
258: VecCopy(A_AUGVEC(spot), VEC_VV(1+loc_it));
259: /*note: an alternate implementation choice would be to only save the AUGVECS and
260: not A_AUGVEC and then apply the PC here to the augvec */
261: }
263: /* update hessenberg matrix and do Gram-Schmidt - new direction is in
264: VEC_VV(1+loc_it)*/
265: (*lgmres->orthog)(ksp,loc_it);
267: /* new entry in hessenburg is the 2-norm of our new direction */
268: VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);
269: *HH(loc_it+1,loc_it) = tt;
270: *HES(loc_it+1,loc_it) = tt;
273: /* check for the happy breakdown */
274: hapbnd = PetscAbsScalar(tt / *GRS(loc_it));/* GRS(loc_it) contains the res_norm from the last iteration */
275: if (hapbnd > lgmres->haptol) hapbnd = lgmres->haptol;
276: if (tt > hapbnd) {
277: tmp = 1.0/tt;
278: VecScale(VEC_VV(loc_it+1),tmp); /* scale new direction by its norm */
279: } else {
280: PetscInfo2(ksp,"Detected happy breakdown, current hapbnd = %G tt = %G\n",hapbnd,tt);
281: hapend = PETSC_TRUE;
282: }
284: /* Now apply rotations to new col of hessenberg (and right side of system),
285: calculate new rotation, and get new residual norm at the same time*/
286: LGMRESUpdateHessenberg(ksp,loc_it,hapend,&res);
287: if (ksp->reason) break;
289: loc_it++;
290: lgmres->it = (loc_it-1); /* Add this here in case it has converged */
291:
292: PetscObjectTakeAccess(ksp);
293: ksp->its++;
294: ksp->rnorm = res;
295: PetscObjectGrantAccess(ksp);
297: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
299: /* Catch error in happy breakdown and signal convergence and break from loop */
300: if (hapend) {
301: if (!ksp->reason) {
302: SETERRQ1(0,"You reached the happy break down,but convergence was not indicated. Residual norm = %G",res);
303: }
304: break;
305: }
306: }
307: /* END OF ITERATION LOOP */
309: KSPLogResidualHistory(ksp,res);
311: /* Monitor if we know that we will not return for a restart */
312: if (ksp->reason || ksp->its >= max_it) {
313: KSPMonitor(ksp, ksp->its, res);
314: }
316: if (itcount) *itcount = loc_it;
318: /*
319: Down here we have to solve for the "best" coefficients of the Krylov
320: columns, add the solution values together, and possibly unwind the
321: preconditioning from the solution
322: */
323:
324: /* Form the solution (or the solution so far) */
325: /* Note: must pass in (loc_it-1) for iteration count so that BuildLgmresSoln
326: properly navigates */
328: BuildLgmresSoln(GRS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);
331: /* LGMRES_MOD collect aug vector and A*augvector for future restarts -
332: only if we will be restarting (i.e. this cycle performed it_total
333: iterations) */
334: if (!ksp->reason && ksp->its < max_it && aug_dim > 0) {
336: /*AUG_TEMP contains the new augmentation vector (assigned in BuildLgmresSoln) */
337: if (!lgmres->aug_ct) {
338: spot = 0;
339: lgmres->aug_ct++;
340: } else if (lgmres->aug_ct < aug_dim) {
341: spot = lgmres->aug_ct;
342: lgmres->aug_ct++;
343: } else { /* truncate */
344: for (ii=0; ii<aug_dim; ii++) {
345: if (lgmres->aug_order[ii] == aug_dim) {
346: spot = ii;
347: }
348: }
349: }
351:
353: VecCopy(AUG_TEMP, AUGVEC(spot));
354: /*need to normalize */
355: VecNorm(AUGVEC(spot), NORM_2, &tmp_norm);
356: inv_tmp_norm = 1.0/tmp_norm;
357: VecScale(AUGVEC(spot),inv_tmp_norm);
359: /*set new aug vector to order 1 - move all others back one */
360: for (ii=0; ii < aug_dim; ii++) {
361: AUG_ORDER(ii)++;
362: }
363: AUG_ORDER(spot) = 1;
365: /*now add the A*aug vector to A_AUGVEC(spot) - this is independ. of preconditioning type*/
366: /* want V*H*y - y is in GRS, V is in VEC_VV and H is in HES */
368:
369: /* first do H+*y */
370: VecSet(AUG_TEMP,0.0);
371: VecGetArray(AUG_TEMP, &avec);
372: for (ii=0; ii < it_total + 1; ii++) {
373: for (jj=0; jj <= ii+1; jj++) {
374: avec[jj] += *HES(jj ,ii) * *GRS(ii);
375: }
376: }
378: /*now multiply result by V+ */
379: VecSet(VEC_TEMP,0.0);
380: VecMAXPY(VEC_TEMP, it_total+1, avec, &VEC_VV(0)); /*answer is in VEC_TEMP*/
381: VecRestoreArray(AUG_TEMP, &avec);
382:
383: /*copy answer to aug location and scale*/
384: VecCopy(VEC_TEMP, A_AUGVEC(spot));
385: VecScale(A_AUGVEC(spot),inv_tmp_norm);
388: }
389: return(0);
390: }
392: /*
393: KSPSolve_LGMRES - This routine applies the LGMRES method.
396: Input Parameter:
397: . ksp - the Krylov space object that was set to use lgmres
399: Output Parameter:
400: . outits - number of iterations used
402: */
406: PetscErrorCode KSPSolve_LGMRES(KSP ksp)
407: {
409: PetscInt cycle_its; /* iterations done in a call to LGMREScycle */
410: PetscInt itcount; /* running total of iterations, incl. those in restarts */
411: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
412: PetscTruth guess_zero = ksp->guess_zero;
413: PetscInt ii; /*LGMRES_MOD variable */
416: if (ksp->calc_sings && !lgmres->Rsvd) {
417: SETERRQ(PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
418: }
419: PetscObjectTakeAccess(ksp);
420: ksp->its = 0;
421: lgmres->aug_ct = 0;
422: lgmres->matvecs = 0;
423: PetscObjectGrantAccess(ksp);
425: /* initialize */
426: itcount = 0;
427: ksp->reason = KSP_CONVERGED_ITERATING;
428: /*LGMRES_MOD*/
429: for (ii=0; ii<lgmres->aug_dim; ii++) {
430: lgmres->aug_order[ii] = 0;
431: }
433: while (!ksp->reason) {
434: /* calc residual - puts in VEC_VV(0) */
435: KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
436: LGMREScycle(&cycle_its,ksp);
437: itcount += cycle_its;
438: if (itcount >= ksp->max_it) {
439: ksp->reason = KSP_DIVERGED_ITS;
440: break;
441: }
442: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
443: }
444: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
445: return(0);
446: }
448: /*
450: KSPDestroy_LGMRES - Frees all memory space used by the Krylov method.
452: */
455: PetscErrorCode KSPDestroy_LGMRES(KSP ksp)
456: {
457: KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;
459: PetscInt i;
462: /* Free the Hessenberg matrices */
463: PetscFree(lgmres->hh_origin);
465: /* Free pointers to user variables */
466: PetscFree(lgmres->vecs);
468: /*LGMRES_MOD - free pointers for extra vectors */
469: PetscFree(lgmres->augvecs);
471: /* free work vectors */
472: for (i=0; i < lgmres->nwork_alloc; i++) {
473: VecDestroyVecs(lgmres->user_work[i],lgmres->mwork_alloc[i]);
474: }
475: PetscFree(lgmres->user_work);
477: /*LGMRES_MOD - free aug work vectors also */
478: /*this was all allocated as one "chunk" */
479: if (lgmres->augwork_alloc) {
480: VecDestroyVecs(lgmres->augvecs_user_work[0],lgmres->augwork_alloc);
481: }
482: PetscFree(lgmres->augvecs_user_work);
483: PetscFree(lgmres->aug_order);
484: PetscFree(lgmres->mwork_alloc);
485: PetscFree(lgmres->nrs);
486: if (lgmres->sol_temp) {VecDestroy(lgmres->sol_temp);}
487: PetscFree(lgmres->Rsvd);
488: PetscFree(lgmres->Dsvd);
489: PetscFree(lgmres->orthogwork);
490: PetscFree(ksp->data);
491: /* clear composed functions */
492: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C","",PETSC_NULL);
493: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C","",PETSC_NULL);
494: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C","",PETSC_NULL);
495: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C","",PETSC_NULL);
496: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C","",PETSC_NULL);
497: return(0);
498: }
500: /*
501: BuildLgmresSoln - create the solution from the starting vector and the
502: current iterates.
504: Input parameters:
505: nrs - work area of size it + 1.
506: vguess - index of initial guess
507: vdest - index of result. Note that vguess may == vdest (replace
508: guess with the solution).
509: it - HH upper triangular part is a block of size (it+1) x (it+1)
511: This is an internal routine that knows about the LGMRES internals.
512: */
515: static PetscErrorCode BuildLgmresSoln(PetscScalar* nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
516: {
517: PetscScalar tt;
519: PetscInt ii,k,j;
520: KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
521: /*LGMRES_MOD */
522: PetscInt it_arnoldi, it_aug;
523: PetscInt jj, spot = 0;
526: /* Solve for solution vector that minimizes the residual */
528: /* If it is < 0, no lgmres steps have been performed */
529: if (it < 0) {
530: if (vdest != vguess) {
531: VecCopy(vguess,vdest);
532: }
533: return(0);
534: }
536: /* so (it+1) lgmres steps HAVE been performed */
538: /* LGMRES_MOD - determine if we need to use augvecs for the soln - do not assume that
539: this is called after the total its allowed for an approx space */
540: if (lgmres->approx_constant) {
541: it_arnoldi = lgmres->max_k - lgmres->aug_ct;
542: } else {
543: it_arnoldi = lgmres->max_k - lgmres->aug_dim;
544: }
545: if (it_arnoldi >= it +1) {
546: it_aug = 0;
547: it_arnoldi = it+1;
548: } else {
549: it_aug = (it + 1) - it_arnoldi;
550: }
552: /* now it_arnoldi indicates the number of matvecs that took place */
553: lgmres->matvecs += it_arnoldi;
555:
556: /* solve the upper triangular system - GRS is the right side and HH is
557: the upper triangular matrix - put soln in nrs */
558: if (*HH(it,it) == 0.0) SETERRQ2(PETSC_ERR_CONV_FAILED,"HH(it,it) is identically zero; it = %D GRS(it) = %G",it,PetscAbsScalar(*GRS(it)));
559: if (*HH(it,it) != 0.0) {
560: nrs[it] = *GRS(it) / *HH(it,it);
561: } else {
562: nrs[it] = 0.0;
563: }
565: for (ii=1; ii<=it; ii++) {
566: k = it - ii;
567: tt = *GRS(k);
568: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
569: nrs[k] = tt / *HH(k,k);
570: }
572: /* Accumulate the correction to the soln of the preconditioned prob. in VEC_TEMP */
573: VecSet(VEC_TEMP,0.0); /* set VEC_TEMP components to 0 */
575: /*LGMRES_MOD - if augmenting has happened we need to form the solution
576: using the augvecs */
577: if (!it_aug) { /* all its are from arnoldi */
578: VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
579: } else { /*use aug vecs */
580: /*first do regular krylov directions */
581: VecMAXPY(VEC_TEMP,it_arnoldi,nrs,&VEC_VV(0));
582: /*now add augmented portions - add contribution of aug vectors one at a time*/
585: for (ii=0; ii<it_aug; ii++) {
586: for (jj=0; jj<lgmres->aug_dim; jj++) {
587: if (lgmres->aug_order[jj] == (ii+1)) {
588: spot = jj;
589: break; /* must have this because there will be duplicates before aug_ct = aug_dim */
590: }
591: }
592: VecAXPY(VEC_TEMP,nrs[it_arnoldi+ii],AUGVEC(spot));
593: }
594: }
595: /* now VEC_TEMP is what we want to keep for augmenting purposes - grab before the
596: preconditioner is "unwound" from right-precondtioning*/
597: VecCopy(VEC_TEMP, AUG_TEMP);
599: KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
601: /* add solution to previous solution */
602: /* put updated solution into vdest.*/
603: if (vdest != vguess) {
604: VecCopy(VEC_TEMP,vdest);
605: }
606: VecAXPY(vdest,1.0,VEC_TEMP);
608: return(0);
609: }
611: /*
613: LGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
614: Return new residual.
616: input parameters:
618: . ksp - Krylov space object
619: . it - plane rotations are applied to the (it+1)th column of the
620: modified hessenberg (i.e. HH(:,it))
621: . hapend - PETSC_FALSE not happy breakdown ending.
623: output parameters:
624: . res - the new residual
625:
626: */
629: static PetscErrorCode LGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscTruth hapend,PetscReal *res)
630: {
631: PetscScalar *hh,*cc,*ss,tt;
632: PetscInt j;
633: KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
636: hh = HH(0,it); /* pointer to beginning of column to update - so
637: incrementing hh "steps down" the (it+1)th col of HH*/
638: cc = CC(0); /* beginning of cosine rotations */
639: ss = SS(0); /* beginning of sine rotations */
641: /* Apply all the previously computed plane rotations to the new column
642: of the Hessenberg matrix */
643: /* Note: this uses the rotation [conj(c) s ; -s c], c= cos(theta), s= sin(theta) */
645: for (j=1; j<=it; j++) {
646: tt = *hh;
647: #if defined(PETSC_USE_COMPLEX)
648: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
649: #else
650: *hh = *cc * tt + *ss * *(hh+1);
651: #endif
652: hh++;
653: *hh = *cc++ * *hh - (*ss++ * tt);
654: /* hh, cc, and ss have all been incremented one by end of loop */
655: }
657: /*
658: compute the new plane rotation, and apply it to:
659: 1) the right-hand-side of the Hessenberg system (GRS)
660: note: it affects GRS(it) and GRS(it+1)
661: 2) the new column of the Hessenberg matrix
662: note: it affects HH(it,it) which is currently pointed to
663: by hh and HH(it+1, it) (*(hh+1))
664: thus obtaining the updated value of the residual...
665: */
667: /* compute new plane rotation */
669: if (!hapend) {
670: #if defined(PETSC_USE_COMPLEX)
671: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
672: #else
673: tt = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
674: #endif
675: if (tt == 0.0) {
676: ksp->reason = KSP_DIVERGED_NULL;
677: return(0);
678: }
679: *cc = *hh / tt; /* new cosine value */
680: *ss = *(hh+1) / tt; /* new sine value */
682: /* apply to 1) and 2) */
683: *GRS(it+1) = - (*ss * *GRS(it));
684: #if defined(PETSC_USE_COMPLEX)
685: *GRS(it) = PetscConj(*cc) * *GRS(it);
686: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
687: #else
688: *GRS(it) = *cc * *GRS(it);
689: *hh = *cc * *hh + *ss * *(hh+1);
690: #endif
692: /* residual is the last element (it+1) of right-hand side! */
693: *res = PetscAbsScalar(*GRS(it+1));
695: } else { /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
696: another rotation matrix (so RH doesn't change). The new residual is
697: always the new sine term times the residual from last time (GRS(it)),
698: but now the new sine rotation would be zero...so the residual should
699: be zero...so we will multiply "zero" by the last residual. This might
700: not be exactly what we want to do here -could just return "zero". */
701:
702: *res = 0.0;
703: }
704: return(0);
705: }
707: /*
709: LGMRESGetNewVectors - This routine allocates more work vectors, starting from
710: VEC_VV(it)
711:
712: */
715: static PetscErrorCode LGMRESGetNewVectors(KSP ksp,PetscInt it)
716: {
717: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
718: PetscInt nwork = lgmres->nwork_alloc; /* number of work vector chunks allocated */
719: PetscInt nalloc; /* number to allocate */
721: PetscInt k;
722:
724: nalloc = lgmres->delta_allocate; /* number of vectors to allocate
725: in a single chunk */
727: /* Adjust the number to allocate to make sure that we don't exceed the
728: number of available slots (lgmres->vecs_allocated)*/
729: if (it + VEC_OFFSET + nalloc >= lgmres->vecs_allocated){
730: nalloc = lgmres->vecs_allocated - it - VEC_OFFSET;
731: }
732: if (!nalloc) return(0);
734: lgmres->vv_allocated += nalloc; /* vv_allocated is the number of vectors allocated */
736: /* work vectors */
737: KSPGetVecs(ksp,nalloc,&lgmres->user_work[nwork],0,PETSC_NULL);
738: PetscLogObjectParents(ksp,nalloc,lgmres->user_work[nwork]);
739: /* specify size of chunk allocated */
740: lgmres->mwork_alloc[nwork] = nalloc;
742: for (k=0; k < nalloc; k++) {
743: lgmres->vecs[it+VEC_OFFSET+k] = lgmres->user_work[nwork][k];
744: }
745:
747: /* LGMRES_MOD - for now we are preallocating the augmentation vectors */
748:
750: /* increment the number of work vector chunks */
751: lgmres->nwork_alloc++;
752: return(0);
753: }
755: /*
757: KSPBuildSolution_LGMRES
759: Input Parameter:
760: . ksp - the Krylov space object
761: . ptr-
763: Output Parameter:
764: . result - the solution
766: Note: this calls BuildLgmresSoln - the same function that LGMREScycle
767: calls directly.
769: */
772: PetscErrorCode KSPBuildSolution_LGMRES(KSP ksp,Vec ptr,Vec *result)
773: {
774: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
778: if (!ptr) {
779: if (!lgmres->sol_temp) {
780: VecDuplicate(ksp->vec_sol,&lgmres->sol_temp);
781: PetscLogObjectParent(ksp,lgmres->sol_temp);
782: }
783: ptr = lgmres->sol_temp;
784: }
785: if (!lgmres->nrs) {
786: /* allocate the work area */
787: PetscMalloc(lgmres->max_k*sizeof(PetscScalar),&lgmres->nrs);
788: PetscLogObjectMemory(ksp,lgmres->max_k*sizeof(PetscScalar));
789: }
790:
791: BuildLgmresSoln(lgmres->nrs,ksp->vec_sol,ptr,ksp,lgmres->it);
792: *result = ptr;
793:
794: return(0);
795: }
801: PetscErrorCode KSPView_LGMRES(KSP ksp,PetscViewer viewer)
802: {
803: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
805: PetscTruth iascii;
808: KSPView_GMRES(ksp,viewer);
809: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
810: if (iascii) {
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);
817: } else {
818: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP LGMRES",((PetscObject)viewer)->type_name);
819: }
820: return(0);
821: }
827: PetscErrorCode KSPSetFromOptions_LGMRES(KSP ksp)
828: {
830: PetscInt aug;
831: KSP_LGMRES *lgmres = (KSP_LGMRES*) ksp->data;
832: PetscTruth flg;
835: KSPSetFromOptions_GMRES(ksp);
836: PetscOptionsHead("KSP LGMRES Options");
837: PetscOptionsName("-ksp_lgmres_constant","Use constant approx. space size","KSPGMRESSetConstant",&flg);
838: if (flg) { lgmres->approx_constant = 1; }
839: PetscOptionsInt("-ksp_lgmres_augment","Number of error approximations to augment the Krylov space with","KSPLGMRESSetAugDim",lgmres->aug_dim,&aug,&flg);
840: if (flg) { KSPLGMRESSetAugDim(ksp,aug); }
841: PetscOptionsTail();
842: return(0);
843: }
846: EXTERN PetscErrorCode KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
847: EXTERN PetscErrorCode KSPComputeEigenvalues_GMRES(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);
849: /*functions for extra lgmres options here*/
853: PetscErrorCode KSPLGMRESSetConstant_LGMRES(KSP ksp)
854: {
855: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
857: lgmres->approx_constant = 1;
858: return(0);
859: }
865: PetscErrorCode KSPLGMRESSetAugDim_LGMRES(KSP ksp,PetscInt aug_dim)
866: {
867: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
870: if (aug_dim < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be positive");
871: if (aug_dim > (lgmres->max_k -1)) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be <= (restart size-1)");
872: lgmres->aug_dim = aug_dim;
873: return(0);
874: }
878: /* end new lgmres functions */
881: /* use these options from gmres */
883: EXTERN PetscErrorCode KSPGMRESSetHapTol_GMRES(KSP,double);
884: EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors_GMRES(KSP);
885: EXTERN PetscErrorCode KSPGMRESSetRestart_GMRES(KSP,PetscInt);
886: EXTERN PetscErrorCode KSPGMRESSetOrthogonalization_GMRES(KSP,PetscErrorCode (*)(KSP,PetscInt));
887: EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType_GMRES(KSP,KSPGMRESCGSRefinementType);
890: /*MC
891: KSPLGMRES - Augments the standard GMRES approximation space with approximations to
892: the error from previous restart cycles.
894: Options Database Keys:
895: + -ksp_gmres_restart <restart> - total approximation space size (Krylov directions + error approximations)
896: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
897: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
898: vectors are allocated as needed)
899: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
900: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
901: . -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
902: stability of the classical Gram-Schmidt orthogonalization.
903: . -ksp_gmres_krylov_monitor - plot the Krylov space generated
904: . -ksp_lgmres_augment <k> - number of error approximations to augment the Krylov space with
905: - -ksp_lgmres_constant - use a constant approx. space size (only affects restart cycles < num. error approx.(k), i.e. the first k restarts)
907: Described in:
908: A. H. Baker, E.R. Jessup, and T.A. Manteuffel. A technique for
909: accelerating the convergence of restarted GMRES. SIAM
910: Journal on Matrix Analysis and Applications, 26 (2005), pp. 962-984.
912: To run LGMRES(m, k) as described in the above paper, use:
913: -ksp_gmres_restart <m+k>
914: -ksp_lgmres_augment <k>
916: Level: beginner
918: Notes: This object is subclassed off of KSPGMRES
920: Contributed by: Allison Baker
922: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPGMRES,
923: KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization()
924: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
925: KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPLGMRESSetAugDim(),
926: KSPGMRESSetConstant()
928: M*/
933: PetscErrorCode KSPCreate_LGMRES(KSP ksp)
934: {
935: KSP_LGMRES *lgmres;
939: PetscNew(KSP_LGMRES,&lgmres);
940: PetscLogObjectMemory(ksp,sizeof(KSP_LGMRES));
941: ksp->data = (void*)lgmres;
942: ksp->ops->buildsolution = KSPBuildSolution_LGMRES;
944: ksp->ops->setup = KSPSetUp_LGMRES;
945: ksp->ops->solve = KSPSolve_LGMRES;
946: ksp->ops->destroy = KSPDestroy_LGMRES;
947: ksp->ops->view = KSPView_LGMRES;
948: ksp->ops->setfromoptions = KSPSetFromOptions_LGMRES;
949: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
950: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
952: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
953: "KSPGMRESSetPreAllocateVectors_GMRES",
954: KSPGMRESSetPreAllocateVectors_GMRES);
955: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
956: "KSPGMRESSetOrthogonalization_GMRES",
957: KSPGMRESSetOrthogonalization_GMRES);
958: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
959: "KSPGMRESSetRestart_GMRES",
960: KSPGMRESSetRestart_GMRES);
961: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
962: "KSPGMRESSetHapTol_GMRES",
963: KSPGMRESSetHapTol_GMRES);
964: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
965: "KSPGMRESSetCGSRefinementType_GMRES",
966: KSPGMRESSetCGSRefinementType_GMRES);
968: /*LGMRES_MOD add extra functions here - like the one to set num of aug vectors */
969: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetConstant_C",
970: "KSPLGMRESSetConstant_LGMRES",
971: KSPLGMRESSetConstant_LGMRES);
973: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetAugDim_C",
974: "KSPLGMRESSetAugDim_LGMRES",
975: KSPLGMRESSetAugDim_LGMRES);
976:
978: /*defaults */
979: lgmres->haptol = 1.0e-30;
980: lgmres->q_preallocate = 0;
981: lgmres->delta_allocate = LGMRES_DELTA_DIRECTIONS;
982: lgmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
983: lgmres->nrs = 0;
984: lgmres->sol_temp = 0;
985: lgmres->max_k = LGMRES_DEFAULT_MAXK;
986: lgmres->Rsvd = 0;
987: lgmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
988: lgmres->orthogwork = 0;
989: /*LGMRES_MOD - new defaults */
990: lgmres->aug_dim = LGMRES_DEFAULT_AUGDIM;
991: lgmres->aug_ct = 0; /* start with no aug vectors */
992: lgmres->approx_constant = 0;
993: lgmres->matvecs = 0;
995: return(0);
996: }