2: /*
3: This file implements FGMRES (a Generalized Minimal Residual) method.
4: Reference: Saad, 1993.
6: Preconditioning: If the preconditioner is constant then this fgmres
7: code is equivalent to RIGHT-PRECONDITIONED GMRES.
8: FGMRES is a modification of gmres that allows the preconditioner to change
9: at each iteration.
11: Restarts: Restarts are basically solves with x0 not equal to zero.
12: 13: Contributed by Allison Baker
15: */
17: #include <../src/ksp/ksp/impls/gmres/fgmres/fgmresimpl.h> /*I "petscksp.h" I*/
18: #define FGMRES_DELTA_DIRECTIONS 10 19: #define FGMRES_DEFAULT_MAXK 30 20: static PetscErrorCode KSPFGMRESGetNewVectors(KSP,PetscInt);
21: static PetscErrorCode KSPFGMRESUpdateHessenberg(KSP,PetscInt,PetscBool ,PetscReal *);
22: static PetscErrorCode KSPFGMRESBuildSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);
24: /*
26: KSPSetUp_FGMRES - Sets up the workspace needed by fgmres.
28: This is called once, usually automatically by KSPSolve() or KSPSetUp(),
29: but can be called directly by KSPSetUp().
31: */
34: PetscErrorCode KSPSetUp_FGMRES(KSP ksp) 35: {
37: PetscInt max_k,k;
38: KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;
41: max_k = fgmres->max_k;
43: KSPSetUp_GMRES(ksp);
45: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&fgmres->prevecs);
46: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&fgmres->prevecs_user_work);
47: PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void*)));
49: KSPGetVecs(ksp,fgmres->vv_allocated,&fgmres->prevecs_user_work[0],0,PETSC_NULL);
50: PetscLogObjectParents(ksp,fgmres->vv_allocated,fgmres->prevecs_user_work[0]);
51: for (k=0; k < fgmres->vv_allocated; k++) {
52: fgmres->prevecs[k] = fgmres->prevecs_user_work[0][k];
53: }
55: return(0);
56: }
58: /*
59: KSPFGMRESResidual - This routine computes the initial residual (NOT PRECONDITIONED)
60: */
63: static PetscErrorCode KSPFGMRESResidual(KSP ksp) 64: {
65: KSP_FGMRES *fgmres = (KSP_FGMRES *)(ksp->data);
66: Mat Amat,Pmat;
67: MatStructure pflag;
71: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
73: /* put A*x into VEC_TEMP */
74: MatMult(Amat,ksp->vec_sol,VEC_TEMP);
75: /* now put residual (-A*x + f) into vec_vv(0) */
76: VecWAXPY(VEC_VV(0),-1.0,VEC_TEMP,ksp->vec_rhs);
77: return(0);
78: }
80: /*
82: KSPFGMRESCycle - Run fgmres, possibly with restart. Return residual
83: history if requested.
85: input parameters:
86: . fgmres - structure containing parameters and work areas
88: output parameters:
89: . itcount - number of iterations used. If null, ignored.
90: . converged - 0 if not converged
92: 93: Notes:
94: On entry, the value in vector VEC_VV(0) should be
95: the initial residual.
98: */
101: PetscErrorCode KSPFGMRESCycle(PetscInt *itcount,KSP ksp)102: {
104: KSP_FGMRES *fgmres = (KSP_FGMRES *)(ksp->data);
105: PetscReal res_norm;
106: PetscReal hapbnd,tt;
107: PetscBool hapend = PETSC_FALSE; /* indicates happy breakdown ending */
109: PetscInt loc_it; /* local count of # of dir. in Krylov space */
110: PetscInt max_k = fgmres->max_k; /* max # of directions Krylov space */
111: Mat Amat,Pmat;
112: MatStructure pflag;
116: /* Number of pseudo iterations since last restart is the number
117: of prestart directions */
118: loc_it = 0;
120: /* note: (fgmres->it) is always set one less than (loc_it) It is used in
121: KSPBUILDSolution_FGMRES, where it is passed to KSPFGMRESBuildSoln.
122: Note that when KSPFGMRESBuildSoln is called from this function,
123: (loc_it -1) is passed, so the two are equivalent */
124: fgmres->it = (loc_it - 1);
126: /* initial residual is in VEC_VV(0) - compute its norm*/
127: VecNorm(VEC_VV(0),NORM_2,&res_norm);
129: /* first entry in right-hand-side of hessenberg system is just
130: the initial residual norm */
131: *RS(0) = res_norm;
133: ksp->rnorm = res_norm;
134: KSPLogResidualHistory(ksp,res_norm);
136: /* check for the convergence - maybe the current guess is good enough */
137: (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);
138: if (ksp->reason) {
139: if (itcount) *itcount = 0;
140: return(0);
141: }
143: /* scale VEC_VV (the initial residual) */
144: VecScale(VEC_VV(0),1.0/res_norm);
145: 146: /* MAIN ITERATION LOOP BEGINNING*/
147: /* keep iterating until we have converged OR generated the max number
148: of directions OR reached the max number of iterations for the method */
149: while (!ksp->reason && loc_it < max_k && ksp->its < ksp->max_it) {
150: if (loc_it) KSPLogResidualHistory(ksp,res_norm);
151: fgmres->it = (loc_it - 1);
152: KSPMonitor(ksp,ksp->its,res_norm);
154: /* see if more space is needed for work vectors */
155: if (fgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
156: KSPFGMRESGetNewVectors(ksp,loc_it+1);
157: /* (loc_it+1) is passed in as number of the first vector that should
158: be allocated */
159: }
161: /* CHANGE THE PRECONDITIONER? */
162: /* ModifyPC is the callback function that can be used to
163: change the PC or its attributes before its applied */
164: (*fgmres->modifypc)(ksp,ksp->its,loc_it,res_norm,fgmres->modifyctx);
165: 166: 167: /* apply PRECONDITIONER to direction vector and store with
168: preconditioned vectors in prevec */
169: KSP_PCApply(ksp,VEC_VV(loc_it),PREVEC(loc_it));
170: 171: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
172: /* Multiply preconditioned vector by operator - put in VEC_VV(loc_it+1) */
173: MatMult(Amat,PREVEC(loc_it),VEC_VV(1+loc_it));
175: 176: /* update hessenberg matrix and do Gram-Schmidt - new direction is in
177: VEC_VV(1+loc_it)*/
178: (*fgmres->orthog)(ksp,loc_it);
180: /* new entry in hessenburg is the 2-norm of our new direction */
181: VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);
182: *HH(loc_it+1,loc_it) = tt;
183: *HES(loc_it+1,loc_it) = tt;
185: /* Happy Breakdown Check */
186: hapbnd = PetscAbsScalar((tt) / *RS(loc_it));
187: /* RS(loc_it) contains the res_norm from the last iteration */
188: hapbnd = PetscMin(fgmres->haptol,hapbnd);
189: if (tt > hapbnd) {
190: /* scale new direction by its norm */
191: VecScale(VEC_VV(loc_it+1),1.0/tt);
192: } else {
193: /* This happens when the solution is exactly reached. */
194: /* So there is no new direction... */
195: VecSet(VEC_TEMP,0.0); /* set VEC_TEMP to 0 */
196: hapend = PETSC_TRUE;
197: }
198: /* note that for FGMRES we could get HES(loc_it+1, loc_it) = 0 and the
199: current solution would not be exact if HES was singular. Note that
200: HH non-singular implies that HES is no singular, and HES is guaranteed
201: to be nonsingular when PREVECS are linearly independent and A is
202: nonsingular (in GMRES, the nonsingularity of A implies the nonsingularity
203: of HES). So we should really add a check to verify that HES is nonsingular.*/
205: 206: /* Now apply rotations to new col of hessenberg (and right side of system),
207: calculate new rotation, and get new residual norm at the same time*/
208: KSPFGMRESUpdateHessenberg(ksp,loc_it,hapend,&res_norm);
209: if (ksp->reason) break;
211: loc_it++;
212: fgmres->it = (loc_it-1); /* Add this here in case it has converged */
213: 214: PetscObjectTakeAccess(ksp);
215: ksp->its++;
216: ksp->rnorm = res_norm;
217: PetscObjectGrantAccess(ksp);
219: (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);
221: /* Catch error in happy breakdown and signal convergence and break from loop */
222: if (hapend) {
223: if (!ksp->reason) {
224: SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"You reached the happy break down,but convergence was not indicated.");
225: }
226: break;
227: }
228: }
229: /* END OF ITERATION LOOP */
231: KSPLogResidualHistory(ksp,res_norm);
233: /*
234: Monitor if we know that we will not return for a restart */
235: if (ksp->reason || ksp->its >= ksp->max_it) {
236: KSPMonitor(ksp,ksp->its,res_norm);
237: }
239: if (itcount) *itcount = loc_it;
241: /*
242: Down here we have to solve for the "best" coefficients of the Krylov
243: columns, add the solution values together, and possibly unwind the
244: preconditioning from the solution
245: */
246: 247: /* Form the solution (or the solution so far) */
248: /* Note: must pass in (loc_it-1) for iteration count so that KSPFGMRESBuildSoln
249: properly navigates */
251: KSPFGMRESBuildSoln(RS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);
253: return(0);
254: }
256: /*
257: KSPSolve_FGMRES - This routine applies the FGMRES method.
260: Input Parameter:
261: . ksp - the Krylov space object that was set to use fgmres
263: Output Parameter:
264: . outits - number of iterations used
266: */
270: PetscErrorCode KSPSolve_FGMRES(KSP ksp)271: {
273: PetscInt cycle_its = 0; /* iterations done in a call to KSPFGMRESCycle */
274: KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;
275: PetscBool diagonalscale;
278: PCGetDiagonalScale(ksp->pc,&diagonalscale);
279: if (diagonalscale) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
281: PetscObjectTakeAccess(ksp);
282: ksp->its = 0;
283: PetscObjectGrantAccess(ksp);
285: /* Compute the initial (NOT preconditioned) residual */
286: if (!ksp->guess_zero) {
287: KSPFGMRESResidual(ksp);
288: } else { /* guess is 0 so residual is F (which is in ksp->vec_rhs) */
289: VecCopy(ksp->vec_rhs,VEC_VV(0));
290: }
291: /* now the residual is in VEC_VV(0) - which is what
292: KSPFGMRESCycle expects... */
293: 294: KSPFGMRESCycle(&cycle_its,ksp);
295: while (!ksp->reason) {
296: KSPFGMRESResidual(ksp);
297: if (ksp->its >= ksp->max_it) break;
298: KSPFGMRESCycle(&cycle_its,ksp);
299: }
300: /* mark lack of convergence */
301: if (ksp->its >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
303: return(0);
304: }
306: extern PetscErrorCode KSPReset_FGMRES(KSP);
307: /*
309: KSPDestroy_FGMRES - Frees all memory space used by the Krylov method.
311: */
314: PetscErrorCode KSPDestroy_FGMRES(KSP ksp)315: {
319: KSPReset_FGMRES(ksp);
320: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPFGMRESSetModifyPC_C","",PETSC_NULL);
321: KSPDestroy_GMRES(ksp);
322: return(0);
323: }
325: /*
326: KSPFGMRESBuildSoln - create the solution from the starting vector and the
327: current iterates.
329: Input parameters:
330: nrs - work area of size it + 1.
331: vguess - index of initial guess
332: vdest - index of result. Note that vguess may == vdest (replace
333: guess with the solution).
334: it - HH upper triangular part is a block of size (it+1) x (it+1)
336: This is an internal routine that knows about the FGMRES internals.
337: */
340: static PetscErrorCode KSPFGMRESBuildSoln(PetscScalar* nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)341: {
342: PetscScalar tt;
344: PetscInt ii,k,j;
345: KSP_FGMRES *fgmres = (KSP_FGMRES *)(ksp->data);
348: /* Solve for solution vector that minimizes the residual */
350: /* If it is < 0, no fgmres steps have been performed */
351: if (it < 0) {
352: VecCopy(vguess,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
353: return(0);
354: }
356: /* so fgmres steps HAVE been performed */
358: /* solve the upper triangular system - RS is the right side and HH is
359: the upper triangular matrix - put soln in nrs */
360: if (*HH(it,it) != 0.0) {
361: nrs[it] = *RS(it) / *HH(it,it);
362: } else {
363: nrs[it] = 0.0;
364: }
365: for (ii=1; ii<=it; ii++) {
366: k = it - ii;
367: tt = *RS(k);
368: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
369: nrs[k] = tt / *HH(k,k);
370: }
372: /* Accumulate the correction to the soln of the preconditioned prob. in
373: VEC_TEMP - note that we use the preconditioned vectors */
374: VecSet(VEC_TEMP,0.0); /* set VEC_TEMP components to 0 */
375: VecMAXPY(VEC_TEMP,it+1,nrs,&PREVEC(0));
377: /* put updated solution into vdest.*/
378: if (vdest != vguess) {
379: VecCopy(VEC_TEMP,vdest);
380: VecAXPY(vdest,1.0,vguess);
381: } else {/* replace guess with solution */
382: VecAXPY(vdest,1.0,VEC_TEMP);
383: }
384: return(0);
385: }
387: /*
389: KSPFGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
390: Return new residual.
392: input parameters:
394: . ksp - Krylov space object
395: . it - plane rotations are applied to the (it+1)th column of the
396: modified hessenberg (i.e. HH(:,it))
397: . hapend - PETSC_FALSE not happy breakdown ending.
399: output parameters:
400: . res - the new residual
401: 402: */
405: static PetscErrorCode KSPFGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool hapend,PetscReal *res)406: {
407: PetscScalar *hh,*cc,*ss,tt;
408: PetscInt j;
409: KSP_FGMRES *fgmres = (KSP_FGMRES *)(ksp->data);
412: hh = HH(0,it); /* pointer to beginning of column to update - so
413: incrementing hh "steps down" the (it+1)th col of HH*/
414: cc = CC(0); /* beginning of cosine rotations */
415: ss = SS(0); /* beginning of sine rotations */
417: /* Apply all the previously computed plane rotations to the new column
418: of the Hessenberg matrix */
419: /* Note: this uses the rotation [conj(c) s ; -s c], c= cos(theta), s= sin(theta),
420: and some refs have [c s ; -conj(s) c] (don't be confused!) */
422: for (j=1; j<=it; j++) {
423: tt = *hh;
424: #if defined(PETSC_USE_COMPLEX)
425: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
426: #else
427: *hh = *cc * tt + *ss * *(hh+1);
428: #endif
429: hh++;
430: *hh = *cc++ * *hh - (*ss++ * tt);
431: /* hh, cc, and ss have all been incremented one by end of loop */
432: }
434: /*
435: compute the new plane rotation, and apply it to:
436: 1) the right-hand-side of the Hessenberg system (RS)
437: note: it affects RS(it) and RS(it+1)
438: 2) the new column of the Hessenberg matrix
439: note: it affects HH(it,it) which is currently pointed to
440: by hh and HH(it+1, it) (*(hh+1))
441: thus obtaining the updated value of the residual...
442: */
444: /* compute new plane rotation */
446: if (!hapend) {
447: #if defined(PETSC_USE_COMPLEX)
448: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
449: #else
450: tt = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
451: #endif
452: if (tt == 0.0) {
453: ksp->reason = KSP_DIVERGED_NULL;
454: return(0);
455: }
457: *cc = *hh / tt; /* new cosine value */
458: *ss = *(hh+1) / tt; /* new sine value */
460: /* apply to 1) and 2) */
461: *RS(it+1) = - (*ss * *RS(it));
462: #if defined(PETSC_USE_COMPLEX)
463: *RS(it) = PetscConj(*cc) * *RS(it);
464: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
465: #else
466: *RS(it) = *cc * *RS(it);
467: *hh = *cc * *hh + *ss * *(hh+1);
468: #endif
470: /* residual is the last element (it+1) of right-hand side! */
471: *res = PetscAbsScalar(*RS(it+1));
473: } else { /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
474: another rotation matrix (so RH doesn't change). The new residual is
475: always the new sine term times the residual from last time (RS(it)),
476: but now the new sine rotation would be zero...so the residual should
477: be zero...so we will multiply "zero" by the last residual. This might
478: not be exactly what we want to do here -could just return "zero". */
479: 480: *res = 0.0;
481: }
482: return(0);
483: }
485: /*
487: KSPFGMRESGetNewVectors - This routine allocates more work vectors, starting from
488: VEC_VV(it), and more preconditioned work vectors, starting
489: from PREVEC(i).
491: */
494: static PetscErrorCode KSPFGMRESGetNewVectors(KSP ksp,PetscInt it)495: {
496: KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;
497: PetscInt nwork = fgmres->nwork_alloc; /* number of work vector chunks allocated */
498: PetscInt nalloc; /* number to allocate */
500: PetscInt k;
501: 503: nalloc = fgmres->delta_allocate; /* number of vectors to allocate
504: in a single chunk */
506: /* Adjust the number to allocate to make sure that we don't exceed the
507: number of available slots (fgmres->vecs_allocated)*/
508: if (it + VEC_OFFSET + nalloc >= fgmres->vecs_allocated){
509: nalloc = fgmres->vecs_allocated - it - VEC_OFFSET;
510: }
511: if (!nalloc) return(0);
513: fgmres->vv_allocated += nalloc; /* vv_allocated is the number of vectors allocated */
515: /* work vectors */
516: KSPGetVecs(ksp,nalloc,&fgmres->user_work[nwork],0,PETSC_NULL);
517: PetscLogObjectParents(ksp,nalloc,fgmres->user_work[nwork]);
518: for (k=0; k < nalloc; k++) {
519: fgmres->vecs[it+VEC_OFFSET+k] = fgmres->user_work[nwork][k];
520: }
521: /* specify size of chunk allocated */
522: fgmres->mwork_alloc[nwork] = nalloc;
524: /* preconditioned vectors */
525: KSPGetVecs(ksp,nalloc,&fgmres->prevecs_user_work[nwork],0,PETSC_NULL);
526: PetscLogObjectParents(ksp,nalloc,fgmres->prevecs_user_work[nwork]);
527: for (k=0; k < nalloc; k++) {
528: fgmres->prevecs[it+VEC_OFFSET+k] = fgmres->prevecs_user_work[nwork][k];
529: }
531: /* increment the number of work vector chunks */
532: fgmres->nwork_alloc++;
533: return(0);
534: }
536: /*
538: KSPBuildSolution_FGMRES
540: Input Parameter:
541: . ksp - the Krylov space object
542: . ptr-
544: Output Parameter:
545: . result - the solution
547: Note: this calls KSPFGMRESBuildSoln - the same function that KSPFGMRESCycle
548: calls directly.
550: */
553: PetscErrorCode KSPBuildSolution_FGMRES(KSP ksp,Vec ptr,Vec *result)554: {
555: KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;
559: if (!ptr) {
560: if (!fgmres->sol_temp) {
561: VecDuplicate(ksp->vec_sol,&fgmres->sol_temp);
562: PetscLogObjectParent(ksp,fgmres->sol_temp);
563: }
564: ptr = fgmres->sol_temp;
565: }
566: if (!fgmres->nrs) {
567: /* allocate the work area */
568: PetscMalloc(fgmres->max_k*sizeof(PetscScalar),&fgmres->nrs);
569: PetscLogObjectMemory(ksp,fgmres->max_k*sizeof(PetscScalar));
570: }
571: 572: KSPFGMRESBuildSoln(fgmres->nrs,ksp->vec_sol,ptr,ksp,fgmres->it);
573: if (result) *result = ptr;
574: 575: return(0);
576: }
580: PetscErrorCode KSPSetFromOptions_FGMRES(KSP ksp)581: {
583: PetscBool flg;
586: KSPSetFromOptions_GMRES(ksp);
587: PetscOptionsHead("KSP flexible GMRES Options");
588: PetscOptionsBoolGroupBegin("-ksp_fgmres_modifypcnochange","do not vary the preconditioner","KSPFGMRESSetModifyPC",&flg);
589: if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCNoChange,0,0);}
590: PetscOptionsBoolGroupEnd("-ksp_fgmres_modifypcksp","vary the KSP based preconditioner","KSPFGMRESSetModifyPC",&flg);
591: if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCKSP,0,0);}
592: PetscOptionsTail();
593: return(0);
594: }
596: typedef PetscErrorCode (*FCN1)(KSP,PetscInt,PetscInt,PetscReal,void*); /* force argument to next function to not be extern C*/
597: typedef PetscErrorCode (*FCN2)(void*);
598: EXTERN_C_BEGIN
601: PetscErrorCode KSPFGMRESSetModifyPC_FGMRES(KSP ksp,FCN1 fcn,void *ctx,FCN2 d)602: {
605: ((KSP_FGMRES *)ksp->data)->modifypc = fcn;
606: ((KSP_FGMRES *)ksp->data)->modifydestroy = d;
607: ((KSP_FGMRES *)ksp->data)->modifyctx = ctx;
608: return(0);
609: }
610: EXTERN_C_END
614: PetscErrorCode KSPReset_FGMRES(KSP ksp)615: {
616: KSP_FGMRES *fgmres = (KSP_FGMRES*)ksp->data;
618: PetscInt i;
621: PetscFree (fgmres->prevecs);
622: for (i=0; i<fgmres->nwork_alloc; i++) {
623: VecDestroyVecs(fgmres->mwork_alloc[i],&fgmres->prevecs_user_work[i]);
624: }
625: PetscFree(fgmres->prevecs_user_work);
626: if (fgmres->modifydestroy) {
627: (*fgmres->modifydestroy)(fgmres->modifyctx);
628: }
629: KSPReset_GMRES(ksp);
630: return(0);
631: }
633: EXTERN_C_BEGIN
636: PetscErrorCode KSPGMRESSetRestart_FGMRES(KSP ksp,PetscInt max_k)637: {
638: KSP_FGMRES *gmres = (KSP_FGMRES *)ksp->data;
642: if (max_k < 1) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
643: if (!ksp->setupstage) {
644: gmres->max_k = max_k;
645: } else if (gmres->max_k != max_k) {
646: gmres->max_k = max_k;
647: ksp->setupstage = KSP_SETUP_NEW;
648: /* free the data structures, then create them again */
649: KSPReset_FGMRES(ksp);
650: }
651: return(0);
652: }
653: EXTERN_C_END
655: EXTERN_C_BEGIN
658: PetscErrorCode KSPGMRESGetRestart_FGMRES(KSP ksp,PetscInt *max_k)659: {
660: KSP_FGMRES *gmres = (KSP_FGMRES *)ksp->data;
663: *max_k = gmres->max_k;
664: return(0);
665: }
666: EXTERN_C_END
668: /*MC
669: KSPFGMRES - Implements the Flexible Generalized Minimal Residual method.
670: developed by Saad with restart
673: Options Database Keys:
674: + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
675: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
676: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
677: vectors are allocated as needed)
678: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
679: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
680: . -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
681: stability of the classical Gram-Schmidt orthogonalization.
682: . -ksp_gmres_krylov_monitor - plot the Krylov space generated
683: . -ksp_fgmres_modifypcnochange - do not change the preconditioner between iterations
684: - -ksp_fgmres_modifypcksp - modify the preconditioner using KSPFGMRESModifyPCKSP()
686: Level: beginner
688: Notes: See KSPFGMRESSetModifyPC() for how to vary the preconditioner between iterations
689: Only right preconditioning is supported.
691: Notes: The following options -ksp_type fgmres -pc_type ksp -ksp_ksp_type bcgs -ksp_view -ksp_pc_type jacobi make the preconditioner (or inner solver)
692: be bi-CG-stab with a preconditioner of Jacobi.
694: Developer Notes: This object is subclassed off of KSPGMRES696: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES,
697: KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
698: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
699: KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPFGMRESSetModifyPC(),
700: KSPFGMRESModifyPCKSP()
702: M*/
704: EXTERN_C_BEGIN
707: PetscErrorCode KSPCreate_FGMRES(KSP ksp)708: {
709: KSP_FGMRES *fgmres;
713: PetscNewLog(ksp,KSP_FGMRES,&fgmres);
714: ksp->data = (void*)fgmres;
715: ksp->ops->buildsolution = KSPBuildSolution_FGMRES;
716: ksp->ops->setup = KSPSetUp_FGMRES;
717: ksp->ops->solve = KSPSolve_FGMRES;
718: ksp->ops->reset = KSPReset_FGMRES;
719: ksp->ops->destroy = KSPDestroy_FGMRES;
720: ksp->ops->view = KSPView_GMRES;
721: ksp->ops->setfromoptions = KSPSetFromOptions_FGMRES;
722: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
723: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
725: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
726: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,0);
728: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
729: "KSPGMRESSetPreAllocateVectors_GMRES",
730: KSPGMRESSetPreAllocateVectors_GMRES);
731: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
732: "KSPGMRESSetOrthogonalization_GMRES",
733: KSPGMRESSetOrthogonalization_GMRES);
734: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",
735: "KSPGMRESGetOrthogonalization_GMRES",
736: KSPGMRESGetOrthogonalization_GMRES);
737: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
738: "KSPGMRESSetRestart_FGMRES",
739: KSPGMRESSetRestart_FGMRES);
740: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetRestart_C",
741: "KSPGMRESGetRestart_FGMRES",
742: KSPGMRESGetRestart_FGMRES);
743: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",
744: "KSPFGMRESSetModifyPC_FGMRES",
745: KSPFGMRESSetModifyPC_FGMRES);
746: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
747: "KSPGMRESSetCGSRefinementType_GMRES",
748: KSPGMRESSetCGSRefinementType_GMRES);
749: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",
750: "KSPGMRESGetCGSRefinementType_GMRES",
751: KSPGMRESGetCGSRefinementType_GMRES);
754: fgmres->haptol = 1.0e-30;
755: fgmres->q_preallocate = 0;
756: fgmres->delta_allocate = FGMRES_DELTA_DIRECTIONS;
757: fgmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
758: fgmres->nrs = 0;
759: fgmres->sol_temp = 0;
760: fgmres->max_k = FGMRES_DEFAULT_MAXK;
761: fgmres->Rsvd = 0;
762: fgmres->orthogwork = 0;
763: fgmres->modifypc = KSPFGMRESModifyPCNoChange;
764: fgmres->modifyctx = PETSC_NULL;
765: fgmres->modifydestroy = PETSC_NULL;
766: fgmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
767: return(0);
768: }
769: EXTERN_C_END