2: #include <../src/ksp/ksp/impls/gmres/lgmres/lgmresimpl.h> /*I petscksp.h I*/
4: #define LGMRES_DELTA_DIRECTIONS 10 5: #define LGMRES_DEFAULT_MAXK 30 6: #define LGMRES_DEFAULT_AUGDIM 2 /*default number of augmentation vectors */ 7: static PetscErrorCode KSPLGMRESGetNewVectors(KSP,PetscInt);
8: static PetscErrorCode KSPLGMRESUpdateHessenberg(KSP,PetscInt,PetscBool ,PetscReal *);
9: static PetscErrorCode KSPLGMRESBuildSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);
13: PetscErrorCode KSPLGMRESSetAugDim(KSP ksp, PetscInt dim) 14: {
18: PetscTryMethod((ksp),"KSPLGMRESSetAugDim_C",(KSP,PetscInt),(ksp,dim));
19: return(0);
20: }
24: PetscErrorCode KSPLGMRESSetConstant(KSP ksp) 25: {
29: PetscTryMethod((ksp),"KSPLGMRESSetConstant_C",(KSP),(ksp));
30: return(0);
31: }
33: /*
34: KSPSetUp_LGMRES - Sets up the workspace needed by lgmres.
36: This is called once, usually automatically by KSPSolve() or KSPSetUp(),
37: but can be called directly by KSPSetUp().
39: */
42: PetscErrorCode KSPSetUp_LGMRES(KSP ksp) 43: {
45: PetscInt max_k,k, aug_dim;
46: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
49: max_k = lgmres->max_k;
50: aug_dim = lgmres->aug_dim;
51: KSPSetUp_GMRES(ksp);
53: /* need array of pointers to augvecs*/
54: PetscMalloc((2 * aug_dim + AUG_OFFSET)*sizeof(void*),&lgmres->augvecs);
55: lgmres->aug_vecs_allocated = 2 *aug_dim + AUG_OFFSET;
56: PetscMalloc((2* aug_dim + AUG_OFFSET)*sizeof(void*),&lgmres->augvecs_user_work);
57: PetscMalloc(aug_dim*sizeof(PetscInt),&lgmres->aug_order);
58: PetscLogObjectMemory(ksp,(aug_dim)*(4*sizeof(void*) + sizeof(PetscInt)) + AUG_OFFSET*2*sizeof(void*));
60: /* for now we will preallocate the augvecs - because aug_dim << restart
61: ... also keep in mind that we need to keep augvecs from cycle to cycle*/
62: lgmres->aug_vv_allocated = 2* aug_dim + AUG_OFFSET;
63: lgmres->augwork_alloc = 2* aug_dim + AUG_OFFSET;
64: KSPGetVecs(ksp,lgmres->aug_vv_allocated,&lgmres->augvecs_user_work[0],0,PETSC_NULL);
65: PetscMalloc((max_k+1)*sizeof(PetscScalar),&lgmres->hwork);
66: PetscLogObjectParents(ksp,lgmres->aug_vv_allocated,lgmres->augvecs_user_work[0]);
67: for (k=0; k<lgmres->aug_vv_allocated; k++) {
68: lgmres->augvecs[k] = lgmres->augvecs_user_work[0][k];
69: }
70: return(0);
71: }
74: /*
76: KSPLGMRESCycle - Run lgmres, possibly with restart. Return residual
77: history if requested.
79: input parameters:
80: . lgmres - structure containing parameters and work areas
82: output parameters:
83: . nres - residuals (from preconditioned system) at each step.
84: If restarting, consider passing nres+it. If null,
85: ignored
86: . itcount - number of iterations used. nres[0] to nres[itcount]
87: are defined. If null, ignored. If null, ignored.
88: . converged - 0 if not converged
90: 91: Notes:
92: On entry, the value in vector VEC_VV(0) should be
93: the initial residual.
96: */
99: PetscErrorCode KSPLGMRESCycle(PetscInt *itcount,KSP ksp)100: {
101: KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
102: PetscReal res_norm, res;
103: PetscReal hapbnd, tt;
104: PetscScalar tmp;
105: PetscBool hapend = PETSC_FALSE; /* indicates happy breakdown ending */
107: PetscInt loc_it; /* local count of # of dir. in Krylov space */
108: PetscInt max_k = lgmres->max_k; /* max approx space size */
109: PetscInt max_it = ksp->max_it; /* max # of overall iterations for the method */
110: /* LGMRES_MOD - new variables*/
111: PetscInt aug_dim = lgmres->aug_dim;
112: PetscInt spot = 0;
113: PetscInt order = 0;
114: PetscInt it_arnoldi; /* number of arnoldi steps to take */
115: PetscInt it_total; /* total number of its to take (=approx space size)*/
116: PetscInt ii, jj;
117: PetscReal tmp_norm;
118: PetscScalar inv_tmp_norm;
119: PetscScalar *avec;
122: /* Number of pseudo iterations since last restart is the number
123: of prestart directions */
124: loc_it = 0;
126: /* LGMRES_MOD: determine number of arnoldi steps to take */
127: /* if approx_constant then we keep the space the same size even if
128: we don't have the full number of aug vectors yet*/
129: if (lgmres->approx_constant) {
130: it_arnoldi = max_k - lgmres->aug_ct;
131: } else {
132: it_arnoldi = max_k - aug_dim;
133: }
135: it_total = it_arnoldi + lgmres->aug_ct;
137: /* initial residual is in VEC_VV(0) - compute its norm*/
138: VecNorm(VEC_VV(0),NORM_2,&res_norm);
139: res = res_norm;
140: 141: /* first entry in right-hand-side of hessenberg system is just
142: the initial residual norm */
143: *GRS(0) = res_norm;
145: /* check for the convergence */
146: if (!res) {
147: if (itcount) *itcount = 0;
148: ksp->reason = KSP_CONVERGED_ATOL;
149: PetscInfo(ksp,"Converged due to zero residual norm on entry\n");
150: return(0);
151: }
153: /* scale VEC_VV (the initial residual) */
154: tmp = 1.0/res_norm; VecScale(VEC_VV(0),tmp);
156: ksp->rnorm = res;
159: /* note: (lgmres->it) is always set one less than (loc_it) It is used in
160: KSPBUILDSolution_LGMRES, where it is passed to KSPLGMRESBuildSoln.
161: Note that when KSPLGMRESBuildSoln is called from this function,
162: (loc_it -1) is passed, so the two are equivalent */
163: lgmres->it = (loc_it - 1);
165: 166: /* MAIN ITERATION LOOP BEGINNING*/
169: /* keep iterating until we have converged OR generated the max number
170: of directions OR reached the max number of iterations for the method */
171: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
172: 173: while (!ksp->reason && loc_it < it_total && ksp->its < max_it) { /* LGMRES_MOD: changed to it_total */
174: KSPLogResidualHistory(ksp,res);
175: lgmres->it = (loc_it - 1);
176: KSPMonitor(ksp,ksp->its,res);
178: /* see if more space is needed for work vectors */
179: if (lgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
180: KSPLGMRESGetNewVectors(ksp,loc_it+1);
181: /* (loc_it+1) is passed in as number of the first vector that should
182: be allocated */
183: }
185: /*LGMRES_MOD: decide whether this is an arnoldi step or an aug step */
186: if (loc_it < it_arnoldi) { /* Arnoldi */
187: KSP_PCApplyBAorAB(ksp,VEC_VV(loc_it),VEC_VV(1+loc_it),VEC_TEMP_MATOP);
188: } else { /*aug step */
189: order = loc_it - it_arnoldi + 1; /* which aug step */
190: for (ii=0; ii<aug_dim; ii++) {
191: if (lgmres->aug_order[ii] == order) {
192: spot = ii;
193: break; /* must have this because there will be duplicates before aug_ct = aug_dim */
194: }
195: }
197: VecCopy(A_AUGVEC(spot), VEC_VV(1+loc_it));
198: /*note: an alternate implementation choice would be to only save the AUGVECS and
199: not A_AUGVEC and then apply the PC here to the augvec */
200: }
202: /* update hessenberg matrix and do Gram-Schmidt - new direction is in
203: VEC_VV(1+loc_it)*/
204: (*lgmres->orthog)(ksp,loc_it);
206: /* new entry in hessenburg is the 2-norm of our new direction */
207: VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);
208: *HH(loc_it+1,loc_it) = tt;
209: *HES(loc_it+1,loc_it) = tt;
212: /* check for the happy breakdown */
213: hapbnd = PetscAbsScalar(tt / *GRS(loc_it));/* GRS(loc_it) contains the res_norm from the last iteration */
214: if (hapbnd > lgmres->haptol) hapbnd = lgmres->haptol;
215: if (tt > hapbnd) {
216: tmp = 1.0/tt;
217: VecScale(VEC_VV(loc_it+1),tmp); /* scale new direction by its norm */
218: } else {
219: PetscInfo2(ksp,"Detected happy breakdown, current hapbnd = %G tt = %G\n",hapbnd,tt);
220: hapend = PETSC_TRUE;
221: }
223: /* Now apply rotations to new col of hessenberg (and right side of system),
224: calculate new rotation, and get new residual norm at the same time*/
225: KSPLGMRESUpdateHessenberg(ksp,loc_it,hapend,&res);
226: if (ksp->reason) break;
228: loc_it++;
229: lgmres->it = (loc_it-1); /* Add this here in case it has converged */
230: 231: PetscObjectTakeAccess(ksp);
232: ksp->its++;
233: ksp->rnorm = res;
234: PetscObjectGrantAccess(ksp);
236: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
238: /* Catch error in happy breakdown and signal convergence and break from loop */
239: if (hapend) {
240: if (!ksp->reason) {
241: SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"You reached the happy break down,but convergence was not indicated. Residual norm = %G",res);
242: }
243: break;
244: }
245: }
246: /* END OF ITERATION LOOP */
248: KSPLogResidualHistory(ksp,res);
250: /* Monitor if we know that we will not return for a restart */
251: if (ksp->reason || ksp->its >= max_it) {
252: KSPMonitor(ksp, ksp->its, res);
253: }
255: if (itcount) *itcount = loc_it;
257: /*
258: Down here we have to solve for the "best" coefficients of the Krylov
259: columns, add the solution values together, and possibly unwind the
260: preconditioning from the solution
261: */
262: 263: /* Form the solution (or the solution so far) */
264: /* Note: must pass in (loc_it-1) for iteration count so that KSPLGMRESBuildSoln
265: properly navigates */
267: KSPLGMRESBuildSoln(GRS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);
270: /* LGMRES_MOD collect aug vector and A*augvector for future restarts -
271: only if we will be restarting (i.e. this cycle performed it_total
272: iterations) */
273: if (!ksp->reason && ksp->its < max_it && aug_dim > 0) {
275: /*AUG_TEMP contains the new augmentation vector (assigned in KSPLGMRESBuildSoln) */
276: if (!lgmres->aug_ct) {
277: spot = 0;
278: lgmres->aug_ct++;
279: } else if (lgmres->aug_ct < aug_dim) {
280: spot = lgmres->aug_ct;
281: lgmres->aug_ct++;
282: } else { /* truncate */
283: for (ii=0; ii<aug_dim; ii++) {
284: if (lgmres->aug_order[ii] == aug_dim) {
285: spot = ii;
286: }
287: }
288: }
290: 292: VecCopy(AUG_TEMP, AUGVEC(spot));
293: /*need to normalize */
294: VecNorm(AUGVEC(spot), NORM_2, &tmp_norm);
295: inv_tmp_norm = 1.0/tmp_norm;
296: VecScale(AUGVEC(spot),inv_tmp_norm);
298: /*set new aug vector to order 1 - move all others back one */
299: for (ii=0; ii < aug_dim; ii++) {
300: AUG_ORDER(ii)++;
301: }
302: AUG_ORDER(spot) = 1;
304: /*now add the A*aug vector to A_AUGVEC(spot) - this is independ. of preconditioning type*/
305: /* want V*H*y - y is in GRS, V is in VEC_VV and H is in HES */
307: 308: /* first do H+*y */
309: avec = lgmres->hwork;
310: PetscMemzero(avec,(it_total+1)*sizeof(*avec));
311: for (ii=0; ii < it_total + 1; ii++) {
312: for (jj=0; jj <= ii+1 && jj < it_total+1; jj++) {
313: avec[jj] += *HES(jj ,ii) * *GRS(ii);
314: }
315: }
317: /*now multiply result by V+ */
318: VecSet(VEC_TEMP,0.0);
319: VecMAXPY(VEC_TEMP, it_total+1, avec, &VEC_VV(0)); /*answer is in VEC_TEMP*/
321: /*copy answer to aug location and scale*/
322: VecCopy(VEC_TEMP, A_AUGVEC(spot));
323: VecScale(A_AUGVEC(spot),inv_tmp_norm);
324: }
325: return(0);
326: }
328: /*
329: KSPSolve_LGMRES - This routine applies the LGMRES method.
332: Input Parameter:
333: . ksp - the Krylov space object that was set to use lgmres
335: Output Parameter:
336: . outits - number of iterations used
338: */
342: PetscErrorCode KSPSolve_LGMRES(KSP ksp)343: {
345: PetscInt cycle_its; /* iterations done in a call to KSPLGMRESCycle */
346: PetscInt itcount; /* running total of iterations, incl. those in restarts */
347: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
348: PetscBool guess_zero = ksp->guess_zero;
349: PetscInt ii; /*LGMRES_MOD variable */
352: if (ksp->calc_sings && !lgmres->Rsvd) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ORDER,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
353: 354: PetscObjectTakeAccess(ksp);
355: ksp->its = 0;
356: lgmres->aug_ct = 0;
357: lgmres->matvecs = 0;
358: PetscObjectGrantAccess(ksp);
360: /* initialize */
361: itcount = 0;
362: ksp->reason = KSP_CONVERGED_ITERATING;
363: /*LGMRES_MOD*/
364: for (ii=0; ii<lgmres->aug_dim; ii++) {
365: lgmres->aug_order[ii] = 0;
366: }
368: while (!ksp->reason) {
369: /* calc residual - puts in VEC_VV(0) */
370: KSPInitialResidual(ksp,ksp->vec_sol,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),ksp->vec_rhs);
371: KSPLGMRESCycle(&cycle_its,ksp);
372: itcount += cycle_its;
373: if (itcount >= ksp->max_it) {
374: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
375: break;
376: }
377: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
378: }
379: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
380: return(0);
381: }
383: /*
385: KSPDestroy_LGMRES - Frees all memory space used by the Krylov method.
387: */
390: PetscErrorCode KSPDestroy_LGMRES(KSP ksp)391: {
392: KSP_LGMRES *lgmres = (KSP_LGMRES*)ksp->data;
396: PetscFree(lgmres->augvecs);
397: if (lgmres->augwork_alloc) {
398: VecDestroyVecs(lgmres->augwork_alloc,&lgmres->augvecs_user_work[0]);
399: }
400: PetscFree(lgmres->augvecs_user_work);
401: PetscFree(lgmres->aug_order);
402: PetscFree(lgmres->hwork);
403: KSPDestroy_GMRES(ksp);
404: return(0);
405: }
407: /*
408: KSPLGMRESBuildSoln - create the solution from the starting vector and the
409: current iterates.
411: Input parameters:
412: nrs - work area of size it + 1.
413: vguess - index of initial guess
414: vdest - index of result. Note that vguess may == vdest (replace
415: guess with the solution).
416: it - HH upper triangular part is a block of size (it+1) x (it+1)
418: This is an internal routine that knows about the LGMRES internals.
419: */
422: static PetscErrorCode KSPLGMRESBuildSoln(PetscScalar* nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)423: {
424: PetscScalar tt;
426: PetscInt ii,k,j;
427: KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
428: /*LGMRES_MOD */
429: PetscInt it_arnoldi, it_aug;
430: PetscInt jj, spot = 0;
433: /* Solve for solution vector that minimizes the residual */
435: /* If it is < 0, no lgmres steps have been performed */
436: if (it < 0) {
437: VecCopy(vguess,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
438: return(0);
439: }
441: /* so (it+1) lgmres steps HAVE been performed */
443: /* LGMRES_MOD - determine if we need to use augvecs for the soln - do not assume that
444: this is called after the total its allowed for an approx space */
445: if (lgmres->approx_constant) {
446: it_arnoldi = lgmres->max_k - lgmres->aug_ct;
447: } else {
448: it_arnoldi = lgmres->max_k - lgmres->aug_dim;
449: }
450: if (it_arnoldi >= it +1) {
451: it_aug = 0;
452: it_arnoldi = it+1;
453: } else {
454: it_aug = (it + 1) - it_arnoldi;
455: }
457: /* now it_arnoldi indicates the number of matvecs that took place */
458: lgmres->matvecs += it_arnoldi;
460: 461: /* solve the upper triangular system - GRS is the right side and HH is
462: the upper triangular matrix - put soln in nrs */
463: if (*HH(it,it) == 0.0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"HH(it,it) is identically zero; it = %D GRS(it) = %G",it,PetscAbsScalar(*GRS(it)));
464: if (*HH(it,it) != 0.0) {
465: nrs[it] = *GRS(it) / *HH(it,it);
466: } else {
467: nrs[it] = 0.0;
468: }
470: for (ii=1; ii<=it; ii++) {
471: k = it - ii;
472: tt = *GRS(k);
473: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
474: nrs[k] = tt / *HH(k,k);
475: }
477: /* Accumulate the correction to the soln of the preconditioned prob. in VEC_TEMP */
478: VecSet(VEC_TEMP,0.0); /* set VEC_TEMP components to 0 */
480: /*LGMRES_MOD - if augmenting has happened we need to form the solution
481: using the augvecs */
482: if (!it_aug) { /* all its are from arnoldi */
483: VecMAXPY(VEC_TEMP,it+1,nrs,&VEC_VV(0));
484: } else { /*use aug vecs */
485: /*first do regular krylov directions */
486: VecMAXPY(VEC_TEMP,it_arnoldi,nrs,&VEC_VV(0));
487: /*now add augmented portions - add contribution of aug vectors one at a time*/
490: for (ii=0; ii<it_aug; ii++) {
491: for (jj=0; jj<lgmres->aug_dim; jj++) {
492: if (lgmres->aug_order[jj] == (ii+1)) {
493: spot = jj;
494: break; /* must have this because there will be duplicates before aug_ct = aug_dim */
495: }
496: }
497: VecAXPY(VEC_TEMP,nrs[it_arnoldi+ii],AUGVEC(spot));
498: }
499: }
500: /* now VEC_TEMP is what we want to keep for augmenting purposes - grab before the
501: preconditioner is "unwound" from right-precondtioning*/
502: VecCopy(VEC_TEMP, AUG_TEMP);
504: KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
506: /* add solution to previous solution */
507: /* put updated solution into vdest.*/
508: VecCopy(vguess,vdest);
509: VecAXPY(vdest,1.0,VEC_TEMP);
511: return(0);
512: }
514: /*
516: KSPLGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
517: Return new residual.
519: input parameters:
521: . ksp - Krylov space object
522: . it - plane rotations are applied to the (it+1)th column of the
523: modified hessenberg (i.e. HH(:,it))
524: . hapend - PETSC_FALSE not happy breakdown ending.
526: output parameters:
527: . res - the new residual
528: 529: */
532: static PetscErrorCode KSPLGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool hapend,PetscReal *res)533: {
534: PetscScalar *hh,*cc,*ss,tt;
535: PetscInt j;
536: KSP_LGMRES *lgmres = (KSP_LGMRES *)(ksp->data);
539: hh = HH(0,it); /* pointer to beginning of column to update - so
540: incrementing hh "steps down" the (it+1)th col of HH*/
541: cc = CC(0); /* beginning of cosine rotations */
542: ss = SS(0); /* beginning of sine rotations */
544: /* Apply all the previously computed plane rotations to the new column
545: of the Hessenberg matrix */
546: /* Note: this uses the rotation [conj(c) s ; -s c], c= cos(theta), s= sin(theta) */
548: for (j=1; j<=it; j++) {
549: tt = *hh;
550: #if defined(PETSC_USE_COMPLEX)
551: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
552: #else
553: *hh = *cc * tt + *ss * *(hh+1);
554: #endif
555: hh++;
556: *hh = *cc++ * *hh - (*ss++ * tt);
557: /* hh, cc, and ss have all been incremented one by end of loop */
558: }
560: /*
561: compute the new plane rotation, and apply it to:
562: 1) the right-hand-side of the Hessenberg system (GRS)
563: note: it affects GRS(it) and GRS(it+1)
564: 2) the new column of the Hessenberg matrix
565: note: it affects HH(it,it) which is currently pointed to
566: by hh and HH(it+1, it) (*(hh+1))
567: thus obtaining the updated value of the residual...
568: */
570: /* compute new plane rotation */
572: if (!hapend) {
573: #if defined(PETSC_USE_COMPLEX)
574: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
575: #else
576: tt = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
577: #endif
578: if (tt == 0.0) {
579: ksp->reason = KSP_DIVERGED_NULL;
580: return(0);
581: }
582: *cc = *hh / tt; /* new cosine value */
583: *ss = *(hh+1) / tt; /* new sine value */
585: /* apply to 1) and 2) */
586: *GRS(it+1) = - (*ss * *GRS(it));
587: #if defined(PETSC_USE_COMPLEX)
588: *GRS(it) = PetscConj(*cc) * *GRS(it);
589: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
590: #else
591: *GRS(it) = *cc * *GRS(it);
592: *hh = *cc * *hh + *ss * *(hh+1);
593: #endif
595: /* residual is the last element (it+1) of right-hand side! */
596: *res = PetscAbsScalar(*GRS(it+1));
598: } else { /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
599: another rotation matrix (so RH doesn't change). The new residual is
600: always the new sine term times the residual from last time (GRS(it)),
601: but now the new sine rotation would be zero...so the residual should
602: be zero...so we will multiply "zero" by the last residual. This might
603: not be exactly what we want to do here -could just return "zero". */
604: 605: *res = 0.0;
606: }
607: return(0);
608: }
610: /*
612: KSPLGMRESGetNewVectors - This routine allocates more work vectors, starting from
613: VEC_VV(it)
614: 615: */
618: static PetscErrorCode KSPLGMRESGetNewVectors(KSP ksp,PetscInt it)619: {
620: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
621: PetscInt nwork = lgmres->nwork_alloc; /* number of work vector chunks allocated */
622: PetscInt nalloc; /* number to allocate */
624: PetscInt k;
625: 627: nalloc = lgmres->delta_allocate; /* number of vectors to allocate
628: in a single chunk */
630: /* Adjust the number to allocate to make sure that we don't exceed the
631: number of available slots (lgmres->vecs_allocated)*/
632: if (it + VEC_OFFSET + nalloc >= lgmres->vecs_allocated){
633: nalloc = lgmres->vecs_allocated - it - VEC_OFFSET;
634: }
635: if (!nalloc) return(0);
637: lgmres->vv_allocated += nalloc; /* vv_allocated is the number of vectors allocated */
639: /* work vectors */
640: KSPGetVecs(ksp,nalloc,&lgmres->user_work[nwork],0,PETSC_NULL);
641: PetscLogObjectParents(ksp,nalloc,lgmres->user_work[nwork]);
642: /* specify size of chunk allocated */
643: lgmres->mwork_alloc[nwork] = nalloc;
645: for (k=0; k < nalloc; k++) {
646: lgmres->vecs[it+VEC_OFFSET+k] = lgmres->user_work[nwork][k];
647: }
648: 650: /* LGMRES_MOD - for now we are preallocating the augmentation vectors */
651: 653: /* increment the number of work vector chunks */
654: lgmres->nwork_alloc++;
655: return(0);
656: }
658: /*
660: KSPBuildSolution_LGMRES
662: Input Parameter:
663: . ksp - the Krylov space object
664: . ptr-
666: Output Parameter:
667: . result - the solution
669: Note: this calls KSPLGMRESBuildSoln - the same function that KSPLGMRESCycle
670: calls directly.
672: */
675: PetscErrorCode KSPBuildSolution_LGMRES(KSP ksp,Vec ptr,Vec *result)676: {
677: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
681: if (!ptr) {
682: if (!lgmres->sol_temp) {
683: VecDuplicate(ksp->vec_sol,&lgmres->sol_temp);
684: PetscLogObjectParent(ksp,lgmres->sol_temp);
685: }
686: ptr = lgmres->sol_temp;
687: }
688: if (!lgmres->nrs) {
689: /* allocate the work area */
690: PetscMalloc(lgmres->max_k*sizeof(PetscScalar),&lgmres->nrs);
691: PetscLogObjectMemory(ksp,lgmres->max_k*sizeof(PetscScalar));
692: }
693: 694: KSPLGMRESBuildSoln(lgmres->nrs,ksp->vec_sol,ptr,ksp,lgmres->it);
695: if (result) *result = ptr;
696: 697: return(0);
698: }
702: PetscErrorCode KSPView_LGMRES(KSP ksp,PetscViewer viewer)703: {
704: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
706: PetscBool iascii;
709: KSPView_GMRES(ksp,viewer);
710: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
711: if (iascii) {
712: /*LGMRES_MOD */
713: PetscViewerASCIIPrintf(viewer," LGMRES: aug. dimension=%D\n",lgmres->aug_dim);
714: if (lgmres->approx_constant) {
715: PetscViewerASCIIPrintf(viewer," LGMRES: approx. space size was kept constant.\n");
716: }
717: PetscViewerASCIIPrintf(viewer," LGMRES: number of matvecs=%D\n",lgmres->matvecs);
718: } else {
719: SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for KSP LGMRES",((PetscObject)viewer)->type_name);
720: }
721: return(0);
722: }
726: PetscErrorCode KSPSetFromOptions_LGMRES(KSP ksp)727: {
729: PetscInt aug;
730: KSP_LGMRES *lgmres = (KSP_LGMRES*) ksp->data;
731: PetscBool flg = PETSC_FALSE;
734: KSPSetFromOptions_GMRES(ksp);
735: PetscOptionsHead("KSP LGMRES Options");
736: PetscOptionsBool("-ksp_lgmres_constant","Use constant approx. space size","KSPGMRESSetConstant",flg,&flg,PETSC_NULL);
737: if (flg) { lgmres->approx_constant = 1; }
738: PetscOptionsInt("-ksp_lgmres_augment","Number of error approximations to augment the Krylov space with","KSPLGMRESSetAugDim",lgmres->aug_dim,&aug,&flg);
739: if (flg) { KSPLGMRESSetAugDim(ksp,aug); }
740: PetscOptionsTail();
741: return(0);
742: }
744: /*functions for extra lgmres options here*/
745: EXTERN_C_BEGIN
748: PetscErrorCode KSPLGMRESSetConstant_LGMRES(KSP ksp)749: {
750: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
752: lgmres->approx_constant = 1;
753: return(0);
754: }
755: EXTERN_C_END
757: EXTERN_C_BEGIN
760: PetscErrorCode KSPLGMRESSetAugDim_LGMRES(KSP ksp,PetscInt aug_dim)761: {
762: KSP_LGMRES *lgmres = (KSP_LGMRES *)ksp->data;
765: if (aug_dim < 0) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be positive");
766: if (aug_dim > (lgmres->max_k -1)) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Augmentation dimension must be <= (restart size-1)");
767: lgmres->aug_dim = aug_dim;
768: return(0);
769: }
770: EXTERN_C_END
773: /* end new lgmres functions */
775: /*MC
776: KSPLGMRES - Augments the standard GMRES approximation space with approximations to
777: the error from previous restart cycles.
779: Options Database Keys:
780: + -ksp_gmres_restart <restart> - total approximation space size (Krylov directions + error approximations)
781: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
782: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
783: vectors are allocated as needed)
784: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
785: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
786: . -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
787: stability of the classical Gram-Schmidt orthogonalization.
788: . -ksp_gmres_krylov_monitor - plot the Krylov space generated
789: . -ksp_lgmres_augment <k> - number of error approximations to augment the Krylov space with
790: - -ksp_lgmres_constant - use a constant approx. space size (only affects restart cycles < num. error approx.(k), i.e. the first k restarts)
792: To run LGMRES(m, k) as described in the above paper, use:
793: -ksp_gmres_restart <m+k>
794: -ksp_lgmres_augment <k>
796: Level: beginner
798: Notes: Supports both left and right preconditioning, but not symmetric.
800: References:
801: A. H. Baker, E.R. Jessup, and T.A. Manteuffel. A technique for
802: accelerating the convergence of restarted GMRES. SIAM
803: Journal on Matrix Analysis and Applications, 26 (2005), pp. 962-984.
805: Developer Notes: This object is subclassed off of KSPGMRES807: Contributed by: Allison Baker
809: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPFGMRES, KSPGMRES,
810: KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
811: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
812: KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPLGMRESSetAugDim(),
813: KSPGMRESSetConstant()
815: M*/
817: EXTERN_C_BEGIN
820: PetscErrorCode KSPCreate_LGMRES(KSP ksp)821: {
822: KSP_LGMRES *lgmres;
826: PetscNewLog(ksp,KSP_LGMRES,&lgmres);
827: ksp->data = (void*)lgmres;
828: ksp->ops->buildsolution = KSPBuildSolution_LGMRES;
830: ksp->ops->setup = KSPSetUp_LGMRES;
831: ksp->ops->solve = KSPSolve_LGMRES;
832: ksp->ops->destroy = KSPDestroy_LGMRES;
833: ksp->ops->view = KSPView_LGMRES;
834: ksp->ops->setfromoptions = KSPSetFromOptions_LGMRES;
835: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
836: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
838: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
839: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);
841: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
842: "KSPGMRESSetPreAllocateVectors_GMRES",
843: KSPGMRESSetPreAllocateVectors_GMRES);
844: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
845: "KSPGMRESSetOrthogonalization_GMRES",
846: KSPGMRESSetOrthogonalization_GMRES);
847: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",
848: "KSPGMRESGetOrthogonalization_GMRES",
849: KSPGMRESGetOrthogonalization_GMRES);
850: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
851: "KSPGMRESSetRestart_GMRES",
852: KSPGMRESSetRestart_GMRES);
853: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetRestart_C",
854: "KSPGMRESGetRestart_GMRES",
855: KSPGMRESGetRestart_GMRES);
856: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
857: "KSPGMRESSetHapTol_GMRES",
858: KSPGMRESSetHapTol_GMRES);
859: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
860: "KSPGMRESSetCGSRefinementType_GMRES",
861: KSPGMRESSetCGSRefinementType_GMRES);
862: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",
863: "KSPGMRESGetCGSRefinementType_GMRES",
864: KSPGMRESGetCGSRefinementType_GMRES);
866: /*LGMRES_MOD add extra functions here - like the one to set num of aug vectors */
867: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetConstant_C",
868: "KSPLGMRESSetConstant_LGMRES",
869: KSPLGMRESSetConstant_LGMRES);
871: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPLGMRESSetAugDim_C",
872: "KSPLGMRESSetAugDim_LGMRES",
873: KSPLGMRESSetAugDim_LGMRES);
874: 876: /*defaults */
877: lgmres->haptol = 1.0e-30;
878: lgmres->q_preallocate = 0;
879: lgmres->delta_allocate = LGMRES_DELTA_DIRECTIONS;
880: lgmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
881: lgmres->nrs = 0;
882: lgmres->sol_temp = 0;
883: lgmres->max_k = LGMRES_DEFAULT_MAXK;
884: lgmres->Rsvd = 0;
885: lgmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
886: lgmres->orthogwork = 0;
887: /*LGMRES_MOD - new defaults */
888: lgmres->aug_dim = LGMRES_DEFAULT_AUGDIM;
889: lgmres->aug_ct = 0; /* start with no aug vectors */
890: lgmres->approx_constant = 0;
891: lgmres->matvecs = 0;
893: return(0);
894: }
895: EXTERN_C_END