Actual source code: fgmres.c
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.
13: Contributed by Allison Baker
15: */
17: #include <../src/ksp/ksp/impls/gmres/fgmres/fgmresimpl.h>
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: */
32: PetscErrorCode KSPSetUp_FGMRES(KSP ksp)
33: {
35: PetscInt max_k,k;
36: KSP_FGMRES *fgmres = (KSP_FGMRES*)ksp->data;
39: max_k = fgmres->max_k;
41: KSPSetUp_GMRES(ksp);
43: PetscMalloc1(max_k+2,&fgmres->prevecs);
44: PetscMalloc1(max_k+2,&fgmres->prevecs_user_work);
45: PetscLogObjectMemory((PetscObject)ksp,(max_k+2)*(2*sizeof(void*)));
47: /* fgmres->vv_allocated includes extra work vectors, which are not used in the additional
48: block of vectors used to store the preconditioned directions, hence the -VEC_OFFSET
49: term for this first allocation of vectors holding preconditioned directions */
50: KSPCreateVecs(ksp,fgmres->vv_allocated-VEC_OFFSET,&fgmres->prevecs_user_work[0],0,NULL);
51: PetscLogObjectParents(ksp,fgmres->vv_allocated-VEC_OFFSET,fgmres->prevecs_user_work[0]);
52: for (k=0; k < fgmres->vv_allocated - VEC_OFFSET ; k++) {
53: fgmres->prevecs[k] = fgmres->prevecs_user_work[0][k];
54: }
55: return(0);
56: }
58: /*
59: KSPFGMRESResidual - This routine computes the initial residual (NOT PRECONDITIONED)
60: */
61: static PetscErrorCode KSPFGMRESResidual(KSP ksp)
62: {
63: KSP_FGMRES *fgmres = (KSP_FGMRES*)(ksp->data);
64: Mat Amat,Pmat;
68: PCGetOperators(ksp->pc,&Amat,&Pmat);
70: /* put A*x into VEC_TEMP */
71: KSP_MatMult(ksp,Amat,ksp->vec_sol,VEC_TEMP);
72: /* now put residual (-A*x + f) into vec_vv(0) */
73: VecWAXPY(VEC_VV(0),-1.0,VEC_TEMP,ksp->vec_rhs);
74: return(0);
75: }
77: /*
79: KSPFGMRESCycle - Run fgmres, possibly with restart. Return residual
80: history if requested.
82: input parameters:
83: . fgmres - structure containing parameters and work areas
85: output parameters:
86: . itcount - number of iterations used. If null, ignored.
87: . converged - 0 if not converged
89: Notes:
90: On entry, the value in vector VEC_VV(0) should be
91: the initial residual.
93: */
94: PetscErrorCode KSPFGMRESCycle(PetscInt *itcount,KSP ksp)
95: {
97: KSP_FGMRES *fgmres = (KSP_FGMRES*)(ksp->data);
98: PetscReal res_norm;
99: PetscReal hapbnd,tt;
100: PetscBool hapend = PETSC_FALSE; /* indicates happy breakdown ending */
102: PetscInt loc_it; /* local count of # of dir. in Krylov space */
103: PetscInt max_k = fgmres->max_k; /* max # of directions Krylov space */
104: Mat Amat,Pmat;
107: /* Number of pseudo iterations since last restart is the number
108: of prestart directions */
109: loc_it = 0;
111: /* note: (fgmres->it) is always set one less than (loc_it) It is used in
112: KSPBUILDSolution_FGMRES, where it is passed to KSPFGMRESBuildSoln.
113: Note that when KSPFGMRESBuildSoln is called from this function,
114: (loc_it -1) is passed, so the two are equivalent */
115: fgmres->it = (loc_it - 1);
117: /* initial residual is in VEC_VV(0) - compute its norm*/
118: VecNorm(VEC_VV(0),NORM_2,&res_norm);
119: KSPCheckNorm(ksp,res_norm);
121: /* first entry in right-hand-side of hessenberg system is just
122: the initial residual norm */
123: *RS(0) = res_norm;
125: ksp->rnorm = res_norm;
126: KSPLogResidualHistory(ksp,res_norm);
127: KSPMonitor(ksp,ksp->its,res_norm);
129: /* check for the convergence - maybe the current guess is good enough */
130: (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);
131: if (ksp->reason) {
132: if (itcount) *itcount = 0;
133: return(0);
134: }
136: /* scale VEC_VV (the initial residual) */
137: VecScale(VEC_VV(0),1.0/res_norm);
139: /* MAIN ITERATION LOOP BEGINNING*/
140: /* keep iterating until we have converged OR generated the max number
141: of directions OR reached the max number of iterations for the method */
142: while (!ksp->reason && loc_it < max_k && ksp->its < ksp->max_it) {
143: if (loc_it) {
144: KSPLogResidualHistory(ksp,res_norm);
145: KSPMonitor(ksp,ksp->its,res_norm);
146: }
147: fgmres->it = (loc_it - 1);
149: /* see if more space is needed for work vectors */
150: if (fgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
151: KSPFGMRESGetNewVectors(ksp,loc_it+1);
152: /* (loc_it+1) is passed in as number of the first vector that should
153: be allocated */
154: }
156: /* CHANGE THE PRECONDITIONER? */
157: /* ModifyPC is the callback function that can be used to
158: change the PC or its attributes before its applied */
159: (*fgmres->modifypc)(ksp,ksp->its,loc_it,res_norm,fgmres->modifyctx);
161: /* apply PRECONDITIONER to direction vector and store with
162: preconditioned vectors in prevec */
163: KSP_PCApply(ksp,VEC_VV(loc_it),PREVEC(loc_it));
165: PCGetOperators(ksp->pc,&Amat,&Pmat);
166: /* Multiply preconditioned vector by operator - put in VEC_VV(loc_it+1) */
167: KSP_MatMult(ksp,Amat,PREVEC(loc_it),VEC_VV(1+loc_it));
169: /* update hessenberg matrix and do Gram-Schmidt - new direction is in
170: VEC_VV(1+loc_it)*/
171: (*fgmres->orthog)(ksp,loc_it);
173: /* new entry in hessenburg is the 2-norm of our new direction */
174: VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);
176: *HH(loc_it+1,loc_it) = tt;
177: *HES(loc_it+1,loc_it) = tt;
179: /* Happy Breakdown Check */
180: hapbnd = PetscAbsScalar((tt) / *RS(loc_it));
181: /* RS(loc_it) contains the res_norm from the last iteration */
182: hapbnd = PetscMin(fgmres->haptol,hapbnd);
183: if (tt > hapbnd) {
184: /* scale new direction by its norm */
185: VecScale(VEC_VV(loc_it+1),1.0/tt);
186: } else {
187: /* This happens when the solution is exactly reached. */
188: /* So there is no new direction... */
189: VecSet(VEC_TEMP,0.0); /* set VEC_TEMP to 0 */
190: hapend = PETSC_TRUE;
191: }
192: /* note that for FGMRES we could get HES(loc_it+1, loc_it) = 0 and the
193: current solution would not be exact if HES was singular. Note that
194: HH non-singular implies that HES is no singular, and HES is guaranteed
195: to be nonsingular when PREVECS are linearly independent and A is
196: nonsingular (in GMRES, the nonsingularity of A implies the nonsingularity
197: of HES). So we should really add a check to verify that HES is nonsingular.*/
199: /* Now apply rotations to new col of hessenberg (and right side of system),
200: calculate new rotation, and get new residual norm at the same time*/
201: KSPFGMRESUpdateHessenberg(ksp,loc_it,hapend,&res_norm);
202: if (ksp->reason) break;
204: loc_it++;
205: fgmres->it = (loc_it-1); /* Add this here in case it has converged */
207: PetscObjectSAWsTakeAccess((PetscObject)ksp);
208: ksp->its++;
209: ksp->rnorm = res_norm;
210: PetscObjectSAWsGrantAccess((PetscObject)ksp);
212: (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);
214: /* Catch error in happy breakdown and signal convergence and break from loop */
215: if (hapend) {
216: if (!ksp->reason) {
217: if (ksp->errorifnotconverged) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"You reached the happy break down, but convergence was not indicated. Residual norm = %g",(double)res_norm);
218: else {
219: ksp->reason = KSP_DIVERGED_BREAKDOWN;
220: break;
221: }
222: }
223: }
224: }
225: /* END OF ITERATION LOOP */
226: KSPLogResidualHistory(ksp,res_norm);
228: /*
229: Monitor if we know that we will not return for a restart */
230: if (loc_it && (ksp->reason || ksp->its >= ksp->max_it)) {
231: KSPMonitor(ksp,ksp->its,res_norm);
232: }
234: if (itcount) *itcount = loc_it;
236: /*
237: Down here we have to solve for the "best" coefficients of the Krylov
238: columns, add the solution values together, and possibly unwind the
239: preconditioning from the solution
240: */
242: /* Form the solution (or the solution so far) */
243: /* Note: must pass in (loc_it-1) for iteration count so that KSPFGMRESBuildSoln
244: properly navigates */
246: KSPFGMRESBuildSoln(RS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);
247: return(0);
248: }
250: /*
251: KSPSolve_FGMRES - This routine applies the FGMRES method.
253: Input Parameter:
254: . ksp - the Krylov space object that was set to use fgmres
256: Output Parameter:
257: . outits - number of iterations used
259: */
261: PetscErrorCode KSPSolve_FGMRES(KSP ksp)
262: {
264: PetscInt cycle_its = 0; /* iterations done in a call to KSPFGMRESCycle */
265: KSP_FGMRES *fgmres = (KSP_FGMRES*)ksp->data;
266: PetscBool diagonalscale;
269: PCGetDiagonalScale(ksp->pc,&diagonalscale);
270: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
272: PetscObjectSAWsTakeAccess((PetscObject)ksp);
273: ksp->its = 0;
274: PetscObjectSAWsGrantAccess((PetscObject)ksp);
276: /* Compute the initial (NOT preconditioned) residual */
277: if (!ksp->guess_zero) {
278: KSPFGMRESResidual(ksp);
279: } else { /* guess is 0 so residual is F (which is in ksp->vec_rhs) */
280: VecCopy(ksp->vec_rhs,VEC_VV(0));
281: }
282: /* now the residual is in VEC_VV(0) - which is what
283: KSPFGMRESCycle expects... */
285: KSPFGMRESCycle(&cycle_its,ksp);
286: while (!ksp->reason) {
287: KSPFGMRESResidual(ksp);
288: if (ksp->its >= ksp->max_it) break;
289: KSPFGMRESCycle(&cycle_its,ksp);
290: }
291: /* mark lack of convergence */
292: if (ksp->its >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
293: return(0);
294: }
296: extern PetscErrorCode KSPReset_FGMRES(KSP);
297: /*
299: KSPDestroy_FGMRES - Frees all memory space used by the Krylov method.
301: */
302: PetscErrorCode KSPDestroy_FGMRES(KSP ksp)
303: {
307: KSPReset_FGMRES(ksp);
308: PetscObjectComposeFunction((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",NULL);
309: KSPDestroy_GMRES(ksp);
310: return(0);
311: }
313: /*
314: KSPFGMRESBuildSoln - create the solution from the starting vector and the
315: current iterates.
317: Input parameters:
318: nrs - work area of size it + 1.
319: vguess - index of initial guess
320: vdest - index of result. Note that vguess may == vdest (replace
321: guess with the solution).
322: it - HH upper triangular part is a block of size (it+1) x (it+1)
324: This is an internal routine that knows about the FGMRES internals.
325: */
326: static PetscErrorCode KSPFGMRESBuildSoln(PetscScalar *nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
327: {
328: PetscScalar tt;
330: PetscInt ii,k,j;
331: KSP_FGMRES *fgmres = (KSP_FGMRES*)(ksp->data);
334: /* Solve for solution vector that minimizes the residual */
336: /* If it is < 0, no fgmres steps have been performed */
337: if (it < 0) {
338: VecCopy(vguess,vdest); /* VecCopy() is smart, exists immediately if vguess == vdest */
339: return(0);
340: }
342: /* so fgmres steps HAVE been performed */
344: /* solve the upper triangular system - RS is the right side and HH is
345: the upper triangular matrix - put soln in nrs */
346: if (*HH(it,it) != 0.0) {
347: nrs[it] = *RS(it) / *HH(it,it);
348: } else {
349: nrs[it] = 0.0;
350: }
351: for (ii=1; ii<=it; ii++) {
352: k = it - ii;
353: tt = *RS(k);
354: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
355: nrs[k] = tt / *HH(k,k);
356: }
358: /* Accumulate the correction to the soln of the preconditioned prob. in
359: VEC_TEMP - note that we use the preconditioned vectors */
360: VecSet(VEC_TEMP,0.0); /* set VEC_TEMP components to 0 */
361: VecMAXPY(VEC_TEMP,it+1,nrs,&PREVEC(0));
363: /* put updated solution into vdest.*/
364: if (vdest != vguess) {
365: VecCopy(VEC_TEMP,vdest);
366: VecAXPY(vdest,1.0,vguess);
367: } else { /* replace guess with solution */
368: VecAXPY(vdest,1.0,VEC_TEMP);
369: }
370: return(0);
371: }
373: /*
375: KSPFGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
376: Return new residual.
378: input parameters:
380: . ksp - Krylov space object
381: . it - plane rotations are applied to the (it+1)th column of the
382: modified hessenberg (i.e. HH(:,it))
383: . hapend - PETSC_FALSE not happy breakdown ending.
385: output parameters:
386: . res - the new residual
388: */
389: static PetscErrorCode KSPFGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscBool hapend,PetscReal *res)
390: {
391: PetscScalar *hh,*cc,*ss,tt;
392: PetscInt j;
393: KSP_FGMRES *fgmres = (KSP_FGMRES*)(ksp->data);
396: hh = HH(0,it); /* pointer to beginning of column to update - so
397: incrementing hh "steps down" the (it+1)th col of HH*/
398: cc = CC(0); /* beginning of cosine rotations */
399: ss = SS(0); /* beginning of sine rotations */
401: /* Apply all the previously computed plane rotations to the new column
402: of the Hessenberg matrix */
403: /* Note: this uses the rotation [conj(c) s ; -s c], c= cos(theta), s= sin(theta),
404: and some refs have [c s ; -conj(s) c] (don't be confused!) */
406: for (j=1; j<=it; j++) {
407: tt = *hh;
408: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
409: hh++;
410: *hh = *cc++ * *hh - (*ss++ * tt);
411: /* hh, cc, and ss have all been incremented one by end of loop */
412: }
414: /*
415: compute the new plane rotation, and apply it to:
416: 1) the right-hand-side of the Hessenberg system (RS)
417: note: it affects RS(it) and RS(it+1)
418: 2) the new column of the Hessenberg matrix
419: note: it affects HH(it,it) which is currently pointed to
420: by hh and HH(it+1, it) (*(hh+1))
421: thus obtaining the updated value of the residual...
422: */
424: /* compute new plane rotation */
426: if (!hapend) {
427: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
428: if (tt == 0.0) {
429: ksp->reason = KSP_DIVERGED_NULL;
430: return(0);
431: }
433: *cc = *hh / tt; /* new cosine value */
434: *ss = *(hh+1) / tt; /* new sine value */
436: /* apply to 1) and 2) */
437: *RS(it+1) = -(*ss * *RS(it));
438: *RS(it) = PetscConj(*cc) * *RS(it);
439: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
441: /* residual is the last element (it+1) of right-hand side! */
442: *res = PetscAbsScalar(*RS(it+1));
444: } else { /* happy breakdown: HH(it+1, it) = 0, therefore we don't need to apply
445: another rotation matrix (so RH doesn't change). The new residual is
446: always the new sine term times the residual from last time (RS(it)),
447: but now the new sine rotation would be zero...so the residual should
448: be zero...so we will multiply "zero" by the last residual. This might
449: not be exactly what we want to do here -could just return "zero". */
451: *res = 0.0;
452: }
453: return(0);
454: }
456: /*
458: KSPFGMRESGetNewVectors - This routine allocates more work vectors, starting from
459: VEC_VV(it), and more preconditioned work vectors, starting
460: from PREVEC(i).
462: */
463: static PetscErrorCode KSPFGMRESGetNewVectors(KSP ksp,PetscInt it)
464: {
465: KSP_FGMRES *fgmres = (KSP_FGMRES*)ksp->data;
466: PetscInt nwork = fgmres->nwork_alloc; /* number of work vector chunks allocated */
467: PetscInt nalloc; /* number to allocate */
469: PetscInt k;
472: nalloc = fgmres->delta_allocate; /* number of vectors to allocate
473: in a single chunk */
475: /* Adjust the number to allocate to make sure that we don't exceed the
476: number of available slots (fgmres->vecs_allocated)*/
477: if (it + VEC_OFFSET + nalloc >= fgmres->vecs_allocated) {
478: nalloc = fgmres->vecs_allocated - it - VEC_OFFSET;
479: }
480: if (!nalloc) return(0);
482: fgmres->vv_allocated += nalloc; /* vv_allocated is the number of vectors allocated */
484: /* work vectors */
485: KSPCreateVecs(ksp,nalloc,&fgmres->user_work[nwork],0,NULL);
486: PetscLogObjectParents(ksp,nalloc,fgmres->user_work[nwork]);
487: for (k=0; k < nalloc; k++) {
488: fgmres->vecs[it+VEC_OFFSET+k] = fgmres->user_work[nwork][k];
489: }
490: /* specify size of chunk allocated */
491: fgmres->mwork_alloc[nwork] = nalloc;
493: /* preconditioned vectors */
494: KSPCreateVecs(ksp,nalloc,&fgmres->prevecs_user_work[nwork],0,NULL);
495: PetscLogObjectParents(ksp,nalloc,fgmres->prevecs_user_work[nwork]);
496: for (k=0; k < nalloc; k++) {
497: fgmres->prevecs[it+k] = fgmres->prevecs_user_work[nwork][k];
498: }
500: /* increment the number of work vector chunks */
501: fgmres->nwork_alloc++;
502: return(0);
503: }
505: /*
507: KSPBuildSolution_FGMRES
509: Input Parameter:
510: . ksp - the Krylov space object
511: . ptr-
513: Output Parameter:
514: . result - the solution
516: Note: this calls KSPFGMRESBuildSoln - the same function that KSPFGMRESCycle
517: calls directly.
519: */
520: PetscErrorCode KSPBuildSolution_FGMRES(KSP ksp,Vec ptr,Vec *result)
521: {
522: KSP_FGMRES *fgmres = (KSP_FGMRES*)ksp->data;
526: if (!ptr) {
527: if (!fgmres->sol_temp) {
528: VecDuplicate(ksp->vec_sol,&fgmres->sol_temp);
529: PetscLogObjectParent((PetscObject)ksp,(PetscObject)fgmres->sol_temp);
530: }
531: ptr = fgmres->sol_temp;
532: }
533: if (!fgmres->nrs) {
534: /* allocate the work area */
535: PetscMalloc1(fgmres->max_k,&fgmres->nrs);
536: PetscLogObjectMemory((PetscObject)ksp,fgmres->max_k*sizeof(PetscScalar));
537: }
539: KSPFGMRESBuildSoln(fgmres->nrs,ksp->vec_sol,ptr,ksp,fgmres->it);
540: if (result) *result = ptr;
541: return(0);
542: }
544: PetscErrorCode KSPSetFromOptions_FGMRES(PetscOptionItems *PetscOptionsObject,KSP ksp)
545: {
547: PetscBool flg;
550: KSPSetFromOptions_GMRES(PetscOptionsObject,ksp);
551: PetscOptionsHead(PetscOptionsObject,"KSP flexible GMRES Options");
552: PetscOptionsBoolGroupBegin("-ksp_fgmres_modifypcnochange","do not vary the preconditioner","KSPFGMRESSetModifyPC",&flg);
553: if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCNoChange,NULL,NULL);}
554: PetscOptionsBoolGroupEnd("-ksp_fgmres_modifypcksp","vary the KSP based preconditioner","KSPFGMRESSetModifyPC",&flg);
555: if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCKSP,NULL,NULL);}
556: PetscOptionsTail();
557: return(0);
558: }
560: typedef PetscErrorCode (*FCN1)(KSP,PetscInt,PetscInt,PetscReal,void*); /* force argument to next function to not be extern C*/
561: typedef PetscErrorCode (*FCN2)(void*);
563: static PetscErrorCode KSPFGMRESSetModifyPC_FGMRES(KSP ksp,FCN1 fcn,void *ctx,FCN2 d)
564: {
567: ((KSP_FGMRES*)ksp->data)->modifypc = fcn;
568: ((KSP_FGMRES*)ksp->data)->modifydestroy = d;
569: ((KSP_FGMRES*)ksp->data)->modifyctx = ctx;
570: return(0);
571: }
573: PetscErrorCode KSPReset_FGMRES(KSP ksp)
574: {
575: KSP_FGMRES *fgmres = (KSP_FGMRES*)ksp->data;
577: PetscInt i;
580: PetscFree (fgmres->prevecs);
581: if (fgmres->nwork_alloc>0) {
582: i=0;
583: /* In the first allocation we allocated VEC_OFFSET fewer vectors in prevecs */
584: VecDestroyVecs(fgmres->mwork_alloc[i]-VEC_OFFSET,&fgmres->prevecs_user_work[i]);
585: for (i=1; i<fgmres->nwork_alloc; i++) {
586: VecDestroyVecs(fgmres->mwork_alloc[i],&fgmres->prevecs_user_work[i]);
587: }
588: }
589: PetscFree(fgmres->prevecs_user_work);
590: if (fgmres->modifydestroy) {
591: (*fgmres->modifydestroy)(fgmres->modifyctx);
592: }
593: KSPReset_GMRES(ksp);
594: return(0);
595: }
597: PetscErrorCode KSPGMRESSetRestart_FGMRES(KSP ksp,PetscInt max_k)
598: {
599: KSP_FGMRES *gmres = (KSP_FGMRES*)ksp->data;
603: if (max_k < 1) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
604: if (!ksp->setupstage) {
605: gmres->max_k = max_k;
606: } else if (gmres->max_k != max_k) {
607: gmres->max_k = max_k;
608: ksp->setupstage = KSP_SETUP_NEW;
609: /* free the data structures, then create them again */
610: KSPReset_FGMRES(ksp);
611: }
612: return(0);
613: }
615: PetscErrorCode KSPGMRESGetRestart_FGMRES(KSP ksp,PetscInt *max_k)
616: {
617: KSP_FGMRES *gmres = (KSP_FGMRES*)ksp->data;
620: *max_k = gmres->max_k;
621: return(0);
622: }
624: /*MC
625: KSPFGMRES - Implements the Flexible Generalized Minimal Residual method.
626: developed by Saad with restart
628: Options Database Keys:
629: + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
630: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
631: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
632: vectors are allocated as needed)
633: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
634: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
635: . -ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the
636: stability of the classical Gram-Schmidt orthogonalization.
637: . -ksp_gmres_krylov_monitor - plot the Krylov space generated
638: . -ksp_fgmres_modifypcnochange - do not change the preconditioner between iterations
639: - -ksp_fgmres_modifypcksp - modify the preconditioner using KSPFGMRESModifyPCKSP()
641: Level: beginner
643: Notes:
644: See KSPFGMRESSetModifyPC() for how to vary the preconditioner between iterations
645: Only right preconditioning is supported.
647: Notes:
648: 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)
649: be bi-CG-stab with a preconditioner of Jacobi.
651: Developer Notes:
652: This object is subclassed off of KSPGMRES
654: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES,
655: KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(),
656: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
657: KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPFGMRESSetModifyPC(),
658: KSPFGMRESModifyPCKSP()
660: M*/
662: PETSC_EXTERN PetscErrorCode KSPCreate_FGMRES(KSP ksp)
663: {
664: KSP_FGMRES *fgmres;
668: PetscNewLog(ksp,&fgmres);
670: ksp->data = (void*)fgmres;
671: ksp->ops->buildsolution = KSPBuildSolution_FGMRES;
672: ksp->ops->setup = KSPSetUp_FGMRES;
673: ksp->ops->solve = KSPSolve_FGMRES;
674: ksp->ops->reset = KSPReset_FGMRES;
675: ksp->ops->destroy = KSPDestroy_FGMRES;
676: ksp->ops->view = KSPView_GMRES;
677: ksp->ops->setfromoptions = KSPSetFromOptions_FGMRES;
678: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
679: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
681: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,3);
682: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
684: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",KSPGMRESSetPreAllocateVectors_GMRES);
685: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",KSPGMRESSetOrthogonalization_GMRES);
686: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetOrthogonalization_C",KSPGMRESGetOrthogonalization_GMRES);
687: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetRestart_C",KSPGMRESSetRestart_FGMRES);
688: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetRestart_C",KSPGMRESGetRestart_FGMRES);
689: PetscObjectComposeFunction((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",KSPFGMRESSetModifyPC_FGMRES);
690: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",KSPGMRESSetCGSRefinementType_GMRES);
691: PetscObjectComposeFunction((PetscObject)ksp,"KSPGMRESGetCGSRefinementType_C",KSPGMRESGetCGSRefinementType_GMRES);
693: fgmres->haptol = 1.0e-30;
694: fgmres->q_preallocate = 0;
695: fgmres->delta_allocate = FGMRES_DELTA_DIRECTIONS;
696: fgmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
697: fgmres->nrs = NULL;
698: fgmres->sol_temp = NULL;
699: fgmres->max_k = FGMRES_DEFAULT_MAXK;
700: fgmres->Rsvd = NULL;
701: fgmres->orthogwork = NULL;
702: fgmres->modifypc = KSPFGMRESModifyPCNoChange;
703: fgmres->modifyctx = NULL;
704: fgmres->modifydestroy = NULL;
705: fgmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
706: return(0);
707: }