Actual source code: fgmres.c
1: #define PETSCKSP_DLL
3: /*
4: This file implements FGMRES (a Generalized Minimal Residual) method.
5: Reference: Saad, 1993.
7: Preconditioning: It the preconditioner is constant then this fgmres
8: code is equivalent to RIGHT-PRECONDITIONED GMRES.
10: Restarts: Restarts are basically solves with x0 not equal to zero.
11:
12: Contributed by Allison Baker
14: */
16: #include src/ksp/ksp/impls/gmres/fgmres/fgmresp.h
17: #define FGMRES_DELTA_DIRECTIONS 10
18: #define FGMRES_DEFAULT_MAXK 30
19: static PetscErrorCode FGMRESGetNewVectors(KSP,PetscInt);
20: static PetscErrorCode FGMRESUpdateHessenberg(KSP,PetscInt,PetscTruth,PetscReal *);
21: static PetscErrorCode BuildFgmresSoln(PetscScalar*,Vec,Vec,KSP,PetscInt);
23: EXTERN PetscErrorCode KSPView_GMRES(KSP,PetscViewer);
24: /*
26: KSPSetUp_FGMRES - Sets up the workspace needed by fgmres.
28: This is called once, usually automatically by KSPSolveQ() or KSPSetUp(),
29: but can be called directly by KSPSetUp().
31: */
34: PetscErrorCode KSPSetUp_FGMRES(KSP ksp)
35: {
36: PetscInt size,hh,hes,rs,cc;
38: PetscInt max_k,k;
39: KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;
42: if (ksp->pc_side == PC_SYMMETRIC) {
43: SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPFGMRES");
44: } else if (ksp->pc_side == PC_LEFT) {
45: SETERRQ(PETSC_ERR_SUP,"no left preconditioning for KSPFGMRES");
46: }
47: max_k = fgmres->max_k;
48: hh = (max_k + 2) * (max_k + 1);
49: hes = (max_k + 1) * (max_k + 1);
50: rs = (max_k + 2);
51: cc = (max_k + 1); /* SS and CC are the same size */
52: size = (hh + hes + rs + 2*cc) * sizeof(PetscScalar);
54: /* Allocate space and set pointers to beginning */
55: PetscMalloc(size,&fgmres->hh_origin);
56: PetscMemzero(fgmres->hh_origin,size);
57: PetscLogObjectMemory(ksp,size); /* HH - modified (by plane rotations) hessenburg */
58: fgmres->hes_origin = fgmres->hh_origin + hh; /* HES - unmodified hessenburg */
59: fgmres->rs_origin = fgmres->hes_origin + hes; /* RS - the right-hand-side of the
60: Hessenberg system */
61: fgmres->cc_origin = fgmres->rs_origin + rs; /* CC - cosines for rotations */
62: fgmres->ss_origin = fgmres->cc_origin + cc; /* SS - sines for rotations */
64: if (ksp->calc_sings) {
65: /* Allocate workspace to hold Hessenberg matrix needed by Eispack */
66: size = (max_k + 3)*(max_k + 9)*sizeof(PetscScalar);
67: PetscMalloc(size,&fgmres->Rsvd);
68: PetscMalloc(5*(max_k+2)*sizeof(PetscReal),&fgmres->Dsvd);
69: PetscLogObjectMemory(ksp,size+5*(max_k+2)*sizeof(PetscReal));
70: }
72: /* Allocate array to hold pointers to user vectors. Note that we need
73: 4 + max_k + 1 (since we need it+1 vectors, and it <= max_k) */
74: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&fgmres->vecs);
75: fgmres->vecs_allocated = VEC_OFFSET + 2 + max_k;
76: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&fgmres->user_work);
77: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(PetscInt),&fgmres->mwork_alloc);
78: PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void*)+sizeof(PetscInt)));
80: /* New for FGMRES - Allocate array to hold pointers to preconditioned
81: vectors - same sizes as user vectors above */
82: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&fgmres->prevecs);
83: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void*),&fgmres->prevecs_user_work);
84: PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void*)));
87: /* if q_preallocate = 0 then only allocate one "chunck" of space (for
88: 5 vectors) - additional will then be allocated from FGMREScycle()
89: as needed. Otherwise, allocate all of the space that could be needed */
90: if (fgmres->q_preallocate) {
91: fgmres->vv_allocated = VEC_OFFSET + 2 + max_k;
92: } else {
93: fgmres->vv_allocated = 5;
94: }
96: /* space for work vectors */
97: KSPGetVecs(ksp,fgmres->vv_allocated,&fgmres->user_work[0],0,PETSC_NULL);
98: PetscLogObjectParents(ksp,fgmres->vv_allocated,fgmres->user_work[0]);
99: for (k=0; k < fgmres->vv_allocated; k++) {
100: fgmres->vecs[k] = fgmres->user_work[0][k];
101: }
103: /* space for preconditioned vectors */
104: KSPGetVecs(ksp,fgmres->vv_allocated,&fgmres->prevecs_user_work[0],0,PETSC_NULL);
105: PetscLogObjectParents(ksp,fgmres->vv_allocated,fgmres->prevecs_user_work[0]);
106: for (k=0; k < fgmres->vv_allocated; k++) {
107: fgmres->prevecs[k] = fgmres->prevecs_user_work[0][k];
108: }
110: /* specify how many work vectors have been allocated in this
111: chunck" (the first one) */
112: fgmres->mwork_alloc[0] = fgmres->vv_allocated;
113: fgmres->nwork_alloc = 1;
115: return(0);
116: }
118: /*
119: FGMRESResidual - This routine computes the initial residual (NOT PRECONDITIONED)
120: */
123: static PetscErrorCode FGMRESResidual(KSP ksp)
124: {
125: KSP_FGMRES *fgmres = (KSP_FGMRES *)(ksp->data);
126: Mat Amat,Pmat;
127: MatStructure pflag;
131: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
133: /* put A*x into VEC_TEMP */
134: MatMult(Amat,ksp->vec_sol,VEC_TEMP);
135: /* now put residual (-A*x + f) into vec_vv(0) */
136: VecWAXPY(VEC_VV(0),-1.0,VEC_TEMP,ksp->vec_rhs);
137: return(0);
138: }
140: /*
142: FGMRESCycle - Run fgmres, possibly with restart. Return residual
143: history if requested.
145: input parameters:
146: . fgmres - structure containing parameters and work areas
148: output parameters:
149: . itcount - number of iterations used. If null, ignored.
150: . converged - 0 if not converged
152:
153: Notes:
154: On entry, the value in vector VEC_VV(0) should be
155: the initial residual.
158: */
161: PetscErrorCode FGMREScycle(PetscInt *itcount,KSP ksp)
162: {
164: KSP_FGMRES *fgmres = (KSP_FGMRES *)(ksp->data);
165: PetscReal res_norm;
166: PetscReal hapbnd,tt;
167: PetscTruth hapend = PETSC_FALSE; /* indicates happy breakdown ending */
169: PetscInt loc_it; /* local count of # of dir. in Krylov space */
170: PetscInt max_k = fgmres->max_k; /* max # of directions Krylov space */
171: Mat Amat,Pmat;
172: MatStructure pflag;
176: /* Number of pseudo iterations since last restart is the number
177: of prestart directions */
178: loc_it = 0;
180: /* note: (fgmres->it) is always set one less than (loc_it) It is used in
181: KSPBUILDSolution_FGMRES, where it is passed to BuildFGmresSoln.
182: Note that when BuildFGmresSoln is called from this function,
183: (loc_it -1) is passed, so the two are equivalent */
184: fgmres->it = (loc_it - 1);
186: /* initial residual is in VEC_VV(0) - compute its norm*/
187: VecNorm(VEC_VV(0),NORM_2,&res_norm);
189: /* first entry in right-hand-side of hessenberg system is just
190: the initial residual norm */
191: *RS(0) = res_norm;
193: /* FYI: AMS calls are for memory snooper */
194: PetscObjectTakeAccess(ksp);
195: ksp->rnorm = res_norm;
196: PetscObjectGrantAccess(ksp);
197: KSPLogResidualHistory(ksp,res_norm);
199: /* check for the convergence - maybe the current guess is good enough */
200: (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);
201: if (ksp->reason) {
202: if (itcount) *itcount = 0;
203: return(0);
204: }
206: /* scale VEC_VV (the initial residual) */
207: VecScale(VEC_VV(0),1.0/res_norm);
208:
209: /* MAIN ITERATION LOOP BEGINNING*/
210: /* keep iterating until we have converged OR generated the max number
211: of directions OR reached the max number of iterations for the method */
212: while (!ksp->reason && loc_it < max_k && ksp->its < ksp->max_it) {
213: KSPLogResidualHistory(ksp,res_norm);
214: fgmres->it = (loc_it - 1);
215: KSPMonitor(ksp,ksp->its,res_norm);
217: /* see if more space is needed for work vectors */
218: if (fgmres->vv_allocated <= loc_it + VEC_OFFSET + 1) {
219: FGMRESGetNewVectors(ksp,loc_it+1);
220: /* (loc_it+1) is passed in as number of the first vector that should
221: be allocated */
222: }
224: /* CHANGE THE PRECONDITIONER? */
225: /* ModifyPC is the callback function that can be used to
226: change the PC or its attributes before its applied */
227: (*fgmres->modifypc)(ksp,ksp->its,loc_it,res_norm,fgmres->modifyctx);
228:
229:
230: /* apply PRECONDITIONER to direction vector and store with
231: preconditioned vectors in prevec */
232: KSP_PCApply(ksp,VEC_VV(loc_it),PREVEC(loc_it));
233:
234: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
235: /* Multiply preconditioned vector by operator - put in VEC_VV(loc_it+1) */
236: MatMult(Amat,PREVEC(loc_it),VEC_VV(1+loc_it));
238:
239: /* update hessenberg matrix and do Gram-Schmidt - new direction is in
240: VEC_VV(1+loc_it)*/
241: (*fgmres->orthog)(ksp,loc_it);
243: /* new entry in hessenburg is the 2-norm of our new direction */
244: VecNorm(VEC_VV(loc_it+1),NORM_2,&tt);
245: *HH(loc_it+1,loc_it) = tt;
246: *HES(loc_it+1,loc_it) = tt;
248: /* Happy Breakdown Check */
249: hapbnd = PetscAbsScalar((tt) / *RS(loc_it));
250: /* RS(loc_it) contains the res_norm from the last iteration */
251: hapbnd = PetscMin(fgmres->haptol,hapbnd);
252: if (tt > hapbnd) {
253: /* scale new direction by its norm */
254: VecScale(VEC_VV(loc_it+1),1.0/tt);
255: } else {
256: /* This happens when the solution is exactly reached. */
257: /* So there is no new direction... */
258: VecSet(VEC_TEMP,0.0); /* set VEC_TEMP to 0 */
259: hapend = PETSC_TRUE;
260: }
261: /* note that for FGMRES we could get HES(loc_it+1, loc_it) = 0 and the
262: current solution would not be exact if HES was singular. Note that
263: HH non-singular implies that HES is no singular, and HES is guaranteed
264: to be nonsingular when PREVECS are linearly independent and A is
265: nonsingular (in GMRES, the nonsingularity of A implies the nonsingularity
266: of HES). So we should really add a check to verify that HES is nonsingular.*/
268:
269: /* Now apply rotations to new col of hessenberg (and right side of system),
270: calculate new rotation, and get new residual norm at the same time*/
271: FGMRESUpdateHessenberg(ksp,loc_it,hapend,&res_norm);
272: if (ksp->reason) break;
274: loc_it++;
275: fgmres->it = (loc_it-1); /* Add this here in case it has converged */
276:
277: PetscObjectTakeAccess(ksp);
278: ksp->its++;
279: ksp->rnorm = res_norm;
280: PetscObjectGrantAccess(ksp);
282: (*ksp->converged)(ksp,ksp->its,res_norm,&ksp->reason,ksp->cnvP);
284: /* Catch error in happy breakdown and signal convergence and break from loop */
285: if (hapend) {
286: if (!ksp->reason) {
287: SETERRQ(0,"You reached the happy break down,but convergence was not indicated.");
288: }
289: break;
290: }
291: }
292: /* END OF ITERATION LOOP */
294: KSPLogResidualHistory(ksp,res_norm);
296: /*
297: Monitor if we know that we will not return for a restart */
298: if (ksp->reason || ksp->its >= ksp->max_it) {
299: KSPMonitor(ksp,ksp->its,res_norm);
300: }
302: if (itcount) *itcount = loc_it;
304: /*
305: Down here we have to solve for the "best" coefficients of the Krylov
306: columns, add the solution values together, and possibly unwind the
307: preconditioning from the solution
308: */
309:
310: /* Form the solution (or the solution so far) */
311: /* Note: must pass in (loc_it-1) for iteration count so that BuildFgmresSoln
312: properly navigates */
314: BuildFgmresSoln(RS(0),ksp->vec_sol,ksp->vec_sol,ksp,loc_it-1);
316: return(0);
317: }
319: /*
320: KSPSolve_FGMRES - This routine applies the FGMRES method.
323: Input Parameter:
324: . ksp - the Krylov space object that was set to use fgmres
326: Output Parameter:
327: . outits - number of iterations used
329: */
333: PetscErrorCode KSPSolve_FGMRES(KSP ksp)
334: {
336: PetscInt cycle_its = 0; /* iterations done in a call to FGMREScycle */
337: KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;
338: PetscTruth diagonalscale;
341: PCDiagonalScale(ksp->pc,&diagonalscale);
342: if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",ksp->type_name);
344: PetscObjectTakeAccess(ksp);
345: ksp->its = 0;
346: PetscObjectGrantAccess(ksp);
348: /* Compute the initial (NOT preconditioned) residual */
349: if (!ksp->guess_zero) {
350: FGMRESResidual(ksp);
351: } else { /* guess is 0 so residual is F (which is in ksp->vec_rhs) */
352: VecCopy(ksp->vec_rhs,VEC_VV(0));
353: }
354: /* now the residual is in VEC_VV(0) - which is what
355: FGMREScycle expects... */
356:
357: FGMREScycle(&cycle_its,ksp);
358: while (!ksp->reason) {
359: FGMRESResidual(ksp);
360: if (ksp->its >= ksp->max_it) break;
361: FGMREScycle(&cycle_its,ksp);
362: }
363: /* mark lack of convergence */
364: if (ksp->its >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
366: return(0);
367: }
369: /*
371: KSPDestroy_FGMRES - Frees all memory space used by the Krylov method.
373: */
376: PetscErrorCode KSPDestroy_FGMRES(KSP ksp)
377: {
378: KSP_FGMRES *fgmres = (KSP_FGMRES*)ksp->data;
380: PetscInt i;
383: /* Free the Hessenberg matrices */
384: PetscFree(fgmres->hh_origin);
386: /* Free pointers to user variables */
387: PetscFree(fgmres->vecs);
388: PetscFree (fgmres->prevecs);
390: /* free work vectors */
391: for (i=0; i < fgmres->nwork_alloc; i++) {
392: VecDestroyVecs(fgmres->user_work[i],fgmres->mwork_alloc[i]);
393: }
394: PetscFree(fgmres->user_work);
396: for (i=0; i < fgmres->nwork_alloc; i++) {
397: VecDestroyVecs(fgmres->prevecs_user_work[i],fgmres->mwork_alloc[i]);
398: }
399: PetscFree(fgmres->prevecs_user_work);
401: PetscFree(fgmres->mwork_alloc);
402: PetscFree(fgmres->nrs);
403: if (fgmres->sol_temp) {VecDestroy(fgmres->sol_temp);}
404: PetscFree(fgmres->Rsvd);
405: PetscFree(fgmres->Dsvd);
406: PetscFree(fgmres->orthogwork);
407: if (fgmres->modifydestroy) {
408: (*fgmres->modifydestroy)(fgmres->modifyctx);
409: }
410: PetscFree(fgmres);
411: return(0);
412: }
414: /*
415: BuildFgmresSoln - create the solution from the starting vector and the
416: current iterates.
418: Input parameters:
419: nrs - work area of size it + 1.
420: vguess - index of initial guess
421: vdest - index of result. Note that vguess may == vdest (replace
422: guess with the solution).
423: it - HH upper triangular part is a block of size (it+1) x (it+1)
425: This is an internal routine that knows about the FGMRES internals.
426: */
429: static PetscErrorCode BuildFgmresSoln(PetscScalar* nrs,Vec vguess,Vec vdest,KSP ksp,PetscInt it)
430: {
431: PetscScalar tt;
433: PetscInt ii,k,j;
434: KSP_FGMRES *fgmres = (KSP_FGMRES *)(ksp->data);
437: /* Solve for solution vector that minimizes the residual */
439: /* If it is < 0, no fgmres steps have been performed */
440: if (it < 0) {
441: if (vdest != vguess) {
442: VecCopy(vguess,vdest);
443: }
444: return(0);
445: }
447: /* so fgmres steps HAVE been performed */
449: /* solve the upper triangular system - RS is the right side and HH is
450: the upper triangular matrix - put soln in nrs */
451: nrs[it] = *RS(it) / *HH(it,it);
452: for (ii=1; ii<=it; ii++) {
453: k = it - ii;
454: tt = *RS(k);
455: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
456: nrs[k] = tt / *HH(k,k);
457: }
459: /* Accumulate the correction to the soln of the preconditioned prob. in
460: VEC_TEMP - note that we use the preconditioned vectors */
461: VecSet(VEC_TEMP,0.0); /* set VEC_TEMP components to 0 */
462: VecMAXPY(VEC_TEMP,it+1,nrs,&PREVEC(0));
464: /* put updated solution into vdest.*/
465: if (vdest != vguess) {
466: VecCopy(VEC_TEMP,vdest);
467: VecAXPY(vdest,1.0,vguess);
468: } else {/* replace guess with solution */
469: VecAXPY(vdest,1.0,VEC_TEMP);
470: }
471: return(0);
472: }
474: /*
476: FGMRESUpdateHessenberg - Do the scalar work for the orthogonalization.
477: Return new residual.
479: input parameters:
481: . ksp - Krylov space object
482: . it - plane rotations are applied to the (it+1)th column of the
483: modified hessenberg (i.e. HH(:,it))
484: . hapend - PETSC_FALSE not happy breakdown ending.
486: output parameters:
487: . res - the new residual
488:
489: */
492: static PetscErrorCode FGMRESUpdateHessenberg(KSP ksp,PetscInt it,PetscTruth hapend,PetscReal *res)
493: {
494: PetscScalar *hh,*cc,*ss,tt;
495: PetscInt j;
496: KSP_FGMRES *fgmres = (KSP_FGMRES *)(ksp->data);
499: hh = HH(0,it); /* pointer to beginning of column to update - so
500: incrementing hh "steps down" the (it+1)th col of HH*/
501: cc = CC(0); /* beginning of cosine rotations */
502: ss = SS(0); /* beginning of sine rotations */
504: /* Apply all the previously computed plane rotations to the new column
505: of the Hessenberg matrix */
506: /* Note: this uses the rotation [conj(c) s ; -s c], c= cos(theta), s= sin(theta),
507: and some refs have [c s ; -conj(s) c] (don't be confused!) */
509: for (j=1; j<=it; j++) {
510: tt = *hh;
511: #if defined(PETSC_USE_COMPLEX)
512: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
513: #else
514: *hh = *cc * tt + *ss * *(hh+1);
515: #endif
516: hh++;
517: *hh = *cc++ * *hh - (*ss++ * tt);
518: /* hh, cc, and ss have all been incremented one by end of loop */
519: }
521: /*
522: compute the new plane rotation, and apply it to:
523: 1) the right-hand-side of the Hessenberg system (RS)
524: note: it affects RS(it) and RS(it+1)
525: 2) the new column of the Hessenberg matrix
526: note: it affects HH(it,it) which is currently pointed to
527: by hh and HH(it+1, it) (*(hh+1))
528: thus obtaining the updated value of the residual...
529: */
531: /* compute new plane rotation */
533: if (!hapend) {
534: #if defined(PETSC_USE_COMPLEX)
535: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
536: #else
537: tt = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
538: #endif
539: if (tt == 0.0) {
540: ksp->reason = KSP_DIVERGED_NULL;
541: return(0);
542: }
544: *cc = *hh / tt; /* new cosine value */
545: *ss = *(hh+1) / tt; /* new sine value */
547: /* apply to 1) and 2) */
548: *RS(it+1) = - (*ss * *RS(it));
549: #if defined(PETSC_USE_COMPLEX)
550: *RS(it) = PetscConj(*cc) * *RS(it);
551: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
552: #else
553: *RS(it) = *cc * *RS(it);
554: *hh = *cc * *hh + *ss * *(hh+1);
555: #endif
557: /* residual is the last element (it+1) of right-hand side! */
558: *res = PetscAbsScalar(*RS(it+1));
560: } else { /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
561: another rotation matrix (so RH doesn't change). The new residual is
562: always the new sine term times the residual from last time (RS(it)),
563: but now the new sine rotation would be zero...so the residual should
564: be zero...so we will multiply "zero" by the last residual. This might
565: not be exactly what we want to do here -could just return "zero". */
566:
567: *res = 0.0;
568: }
569: return(0);
570: }
572: /*
574: FGMRESGetNewVectors - This routine allocates more work vectors, starting from
575: VEC_VV(it), and more preconditioned work vectors, starting
576: from PREVEC(i).
578: */
581: static PetscErrorCode FGMRESGetNewVectors(KSP ksp,PetscInt it)
582: {
583: KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;
584: PetscInt nwork = fgmres->nwork_alloc; /* number of work vector chunks allocated */
585: PetscInt nalloc; /* number to allocate */
587: PetscInt k;
588:
590: nalloc = fgmres->delta_allocate; /* number of vectors to allocate
591: in a single chunk */
593: /* Adjust the number to allocate to make sure that we don't exceed the
594: number of available slots (fgmres->vecs_allocated)*/
595: if (it + VEC_OFFSET + nalloc >= fgmres->vecs_allocated){
596: nalloc = fgmres->vecs_allocated - it - VEC_OFFSET;
597: }
598: if (!nalloc) return(0);
600: fgmres->vv_allocated += nalloc; /* vv_allocated is the number of vectors allocated */
602: /* work vectors */
603: KSPGetVecs(ksp,nalloc,&fgmres->user_work[nwork],0,PETSC_NULL);
604: PetscLogObjectParents(ksp,nalloc,fgmres->user_work[nwork]);
605: for (k=0; k < nalloc; k++) {
606: fgmres->vecs[it+VEC_OFFSET+k] = fgmres->user_work[nwork][k];
607: }
608: /* specify size of chunk allocated */
609: fgmres->mwork_alloc[nwork] = nalloc;
611: /* preconditioned vectors */
612: KSPGetVecs(ksp,nalloc,&fgmres->prevecs_user_work[nwork],0,PETSC_NULL);
613: PetscLogObjectParents(ksp,nalloc,fgmres->prevecs_user_work[nwork]);
614: for (k=0; k < nalloc; k++) {
615: fgmres->prevecs[it+VEC_OFFSET+k] = fgmres->prevecs_user_work[nwork][k];
616: }
618: /* increment the number of work vector chunks */
619: fgmres->nwork_alloc++;
620: return(0);
621: }
623: /*
625: KSPBuildSolution_FGMRES
627: Input Parameter:
628: . ksp - the Krylov space object
629: . ptr-
631: Output Parameter:
632: . result - the solution
634: Note: this calls BuildFgmresSoln - the same function that FGMREScycle
635: calls directly.
637: */
640: PetscErrorCode KSPBuildSolution_FGMRES(KSP ksp,Vec ptr,Vec *result)
641: {
642: KSP_FGMRES *fgmres = (KSP_FGMRES *)ksp->data;
646: if (!ptr) {
647: if (!fgmres->sol_temp) {
648: VecDuplicate(ksp->vec_sol,&fgmres->sol_temp);
649: PetscLogObjectParent(ksp,fgmres->sol_temp);
650: }
651: ptr = fgmres->sol_temp;
652: }
653: if (!fgmres->nrs) {
654: /* allocate the work area */
655: PetscMalloc(fgmres->max_k*sizeof(PetscScalar),&fgmres->nrs);
656: PetscLogObjectMemory(ksp,fgmres->max_k*sizeof(PetscScalar));
657: }
658:
659: BuildFgmresSoln(fgmres->nrs,ksp->vec_sol,ptr,ksp,fgmres->it);
660: *result = ptr;
661:
662: return(0);
663: }
669: PetscErrorCode KSPSetFromOptions_FGMRES(KSP ksp)
670: {
672: PetscTruth flg;
675: KSPSetFromOptions_GMRES(ksp);
676: PetscOptionsHead("KSP flexible GMRES Options");
677: PetscOptionsTruthGroupBegin("-ksp_fgmres_modifypcnochange","do not vary the preconditioner","KSPFGMRESSetModifyPC",&flg);
678: if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCNoChange,0,0);}
679: PetscOptionsTruthGroupEnd("-ksp_fgmres_modifypcksp","vary the KSP based preconditioner","KSPFGMRESSetModifyPC",&flg);
680: if (flg) {KSPFGMRESSetModifyPC(ksp,KSPFGMRESModifyPCKSP,0,0);}
681: PetscOptionsTail();
682: return(0);
683: }
685: EXTERN PetscErrorCode KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
686: EXTERN PetscErrorCode KSPComputeEigenvalues_GMRES(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);
689: typedef PetscErrorCode (*FCN2)(void*);
693: PetscErrorCode PETSCKSP_DLLEXPORT KSPFGMRESSetModifyPC_FGMRES(KSP ksp,FCN1 fcn,void *ctx,FCN2 d)
694: {
697: ((KSP_FGMRES *)ksp->data)->modifypc = fcn;
698: ((KSP_FGMRES *)ksp->data)->modifydestroy = d;
699: ((KSP_FGMRES *)ksp->data)->modifyctx = ctx;
700: return(0);
701: }
705: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetPreAllocateVectors_GMRES(KSP);
706: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart_GMRES(KSP,PetscInt);
707: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetOrthogonalization_GMRES(KSP,PetscErrorCode (*)(KSP,PetscInt));
710: EXTERN PetscErrorCode KSPDestroy_GMRES_Internal(KSP);
714: PetscErrorCode KSPDestroy_FGMRES_Internal(KSP ksp)
715: {
716: KSP_FGMRES *gmres = (KSP_FGMRES*)ksp->data;
720: KSPDestroy_GMRES_Internal(ksp);
721: PetscFree (gmres->prevecs);
722: PetscFree(gmres->prevecs_user_work);
723: if (gmres->modifydestroy) {
724: (*gmres->modifydestroy)(gmres->modifyctx);
725: }
726: gmres->prevecs = 0;
727: gmres->prevecs_user_work = 0;
728: return(0);
729: }
734: PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetRestart_FGMRES(KSP ksp,PetscInt max_k)
735: {
736: KSP_FGMRES *gmres = (KSP_FGMRES *)ksp->data;
740: if (max_k < 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Restart must be positive");
741: if (!ksp->setupcalled) {
742: gmres->max_k = max_k;
743: } else if (gmres->max_k != max_k) {
744: gmres->max_k = max_k;
745: ksp->setupcalled = 0;
746: /* free the data structures, then create them again */
747: KSPDestroy_FGMRES_Internal(ksp);
748: }
749: return(0);
750: }
754: EXTERN PetscErrorCode PETSCKSP_DLLEXPORT KSPGMRESSetCGSRefinementType_GMRES(KSP,KSPGMRESCGSRefinementType);
757: /*MC
758: KSPFGMRES - Implements the Flexible Generalized Minimal Residual method.
759: developed by Saad with restart
762: Options Database Keys:
763: + -ksp_gmres_restart <restart> - the number of Krylov directions to orthogonalize against
764: . -ksp_gmres_haptol <tol> - sets the tolerance for "happy ending" (exact convergence)
765: . -ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of
766: vectors are allocated as needed)
767: . -ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
768: . -ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
769: . -ksp_gmres_cgs_refinement_type <never,ifneeded,always> - determine if iterative refinement is used to increase the
770: stability of the classical Gram-Schmidt orthogonalization.
771: . -ksp_gmres_krylov_monitor - plot the Krylov space generated
772: . -ksp_fgmres_modifypcnochange - do not change the preconditioner between iterations
773: - -ksp_fgmres_modifypcksp - modify the preconditioner using KSPFGMRESModifyPCKSP()
775: Level: beginner
777: Notes: See KSPFGMRESSetModifyPC() for how to vary the preconditioner between iterations
778: This object is subclassed off of KSPGMRES
780: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPGMRES, KSPLGMRES,
781: KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization()
782: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(),
783: KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESKrylovMonitor(), KSPFGMRESSetModifyPC(),
784: KSPFGMRESModifyPCKSP()
786: M*/
791: PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate_FGMRES(KSP ksp)
792: {
793: KSP_FGMRES *fgmres;
797: PetscNew(KSP_FGMRES,&fgmres);
798: PetscLogObjectMemory(ksp,sizeof(KSP_FGMRES));
799: ksp->data = (void*)fgmres;
800: ksp->ops->buildsolution = KSPBuildSolution_FGMRES;
802: ksp->ops->setup = KSPSetUp_FGMRES;
803: ksp->ops->solve = KSPSolve_FGMRES;
804: ksp->ops->destroy = KSPDestroy_FGMRES;
805: ksp->ops->view = KSPView_GMRES;
806: ksp->ops->setfromoptions = KSPSetFromOptions_FGMRES;
807: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
808: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
810: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
811: "KSPGMRESSetPreAllocateVectors_GMRES",
812: KSPGMRESSetPreAllocateVectors_GMRES);
813: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
814: "KSPGMRESSetOrthogonalization_GMRES",
815: KSPGMRESSetOrthogonalization_GMRES);
816: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
817: "KSPGMRESSetRestart_FGMRES",
818: KSPGMRESSetRestart_FGMRES);
819: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPFGMRESSetModifyPC_C",
820: "KSPFGMRESSetModifyPC_FGMRES",
821: KSPFGMRESSetModifyPC_FGMRES);
822: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetCGSRefinementType_C",
823: "KSPGMRESSetCGSRefinementType_GMRES",
824: KSPGMRESSetCGSRefinementType_GMRES);
827: fgmres->haptol = 1.0e-30;
828: fgmres->q_preallocate = 0;
829: fgmres->delta_allocate = FGMRES_DELTA_DIRECTIONS;
830: fgmres->orthog = KSPGMRESClassicalGramSchmidtOrthogonalization;
831: fgmres->nrs = 0;
832: fgmres->sol_temp = 0;
833: fgmres->max_k = FGMRES_DEFAULT_MAXK;
834: fgmres->Rsvd = 0;
835: fgmres->orthogwork = 0;
836: fgmres->modifypc = KSPFGMRESModifyPCNoChange;
837: fgmres->modifyctx = PETSC_NULL;
838: fgmres->modifydestroy = PETSC_NULL;
839: fgmres->cgstype = KSP_GMRES_CGS_REFINE_NEVER;
840: /*
841: This is not great since it changes this without explicit request from the user
842: but there is no left preconditioning in the FGMRES
843: */
844: PetscInfo(ksp,"WARNING! Setting PC_SIDE for FGMRES to right!\n");
845: ksp->pc_side = PC_RIGHT;
847: return(0);
848: }