Actual source code: gmres.c
1: /*$Id: gmres.c,v 1.172 2001/04/10 19:36:32 bsmith Exp $*/
3: /*
4: This file implements GMRES (a Generalized Minimal Residual) method.
5: Reference: Saad and Schultz, 1986.
8: Some comments on left vs. right preconditioning, and restarts.
9: Left and right preconditioning.
10: If right preconditioning is chosen, then the problem being solved
11: by gmres is actually
12: My = AB^-1 y = f
13: so the initial residual is
14: r = f - Mx
15: Note that B^-1 y = x or y = B x, and if x is non-zero, the initial
16: residual is
17: r = f - A x
18: The final solution is then
19: x = B^-1 y
21: If left preconditioning is chosen, then the problem being solved is
22: My = B^-1 A x = B^-1 f,
23: and the initial residual is
24: r = B^-1(f - Ax)
26: Restarts: Restarts are basically solves with x0 not equal to zero.
27: Note that we can elliminate an extra application of B^-1 between
28: restarts as long as we don't require that the solution at the end
29: of a unsuccessful gmres iteration always be the solution x.
30: */
32: #include "src/sles/ksp/impls/gmres/gmresp.h" /*I "petscksp.h" I*/
33: #define GMRES_DELTA_DIRECTIONS 10
34: #define GMRES_DEFAULT_MAXK 30
35: static int GMRESGetNewVectors(KSP,int);
36: static int GMRESUpdateHessenberg(KSP,int,PetscTruth,PetscReal*);
37: static int BuildGmresSoln(Scalar*,Vec,Vec,KSP,int);
39: int KSPSetUp_GMRES(KSP ksp)
40: {
41: unsigned int size,hh,hes,rs,cc;
42: int ierr,max_k,k;
43: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
46: if (ksp->pc_side == PC_SYMMETRIC) {
47: SETERRQ(2,"no symmetric preconditioning for KSPGMRES");
48: }
49: max_k = gmres->max_k;
50: hh = (max_k + 2) * (max_k + 1);
51: hes = (max_k + 1) * (max_k + 1);
52: rs = (max_k + 2);
53: cc = (max_k + 1);
54: size = (hh + hes + rs + 2*cc) * sizeof(Scalar);
56: PetscMalloc(size,&gmres->hh_origin);
57: PetscMemzero(gmres->hh_origin,size);
58: PetscLogObjectMemory(ksp,size);
59: gmres->hes_origin = gmres->hh_origin + hh;
60: gmres->rs_origin = gmres->hes_origin + hes;
61: gmres->cc_origin = gmres->rs_origin + rs;
62: gmres->ss_origin = gmres->cc_origin + cc;
64: if (ksp->calc_sings) {
65: /* Allocate workspace to hold Hessenberg matrix needed by Eispack */
66: size = (max_k + 3)*(max_k + 9)*sizeof(Scalar);
67: PetscMalloc(size,&gmres->Rsvd);
68: PetscMalloc(5*(max_k+2)*sizeof(PetscReal),&gmres->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 *),&gmres->vecs);
75: gmres->vecs_allocated = VEC_OFFSET + 2 + max_k;
76: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void *),&gmres->user_work);
77: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(int),&gmres->mwork_alloc);
78: PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void *)+sizeof(int)));
80: if (gmres->q_preallocate) {
81: gmres->vv_allocated = VEC_OFFSET + 2 + max_k;
82: VecDuplicateVecs(VEC_RHS,gmres->vv_allocated,&gmres->user_work[0]);
83: PetscLogObjectParents(ksp,gmres->vv_allocated,gmres->user_work[0]);
84: gmres->mwork_alloc[0] = gmres->vv_allocated;
85: gmres->nwork_alloc = 1;
86: for (k=0; k<gmres->vv_allocated; k++) {
87: gmres->vecs[k] = gmres->user_work[0][k];
88: }
89: } else {
90: gmres->vv_allocated = 5;
91: VecDuplicateVecs(ksp->vec_rhs,5,&gmres->user_work[0]);
92: PetscLogObjectParents(ksp,5,gmres->user_work[0]);
93: gmres->mwork_alloc[0] = 5;
94: gmres->nwork_alloc = 1;
95: for (k=0; k<gmres->vv_allocated; k++) {
96: gmres->vecs[k] = gmres->user_work[0][k];
97: }
98: }
99: return(0);
100: }
102: /*
103: Run gmres, possibly with restart. Return residual history if requested.
104: input parameters:
106: . gmres - structure containing parameters and work areas
108: output parameters:
109: . nres - residuals (from preconditioned system) at each step.
110: If restarting, consider passing nres+it. If null,
111: ignored
112: . itcount - number of iterations used. nres[0] to nres[itcount]
113: are defined. If null, ignored.
114:
115: Notes:
116: On entry, the value in vector VEC_VV(0) should be the initial residual
117: (this allows shortcuts where the initial preconditioned residual is 0).
118: */
119: int GMREScycle(int *itcount,KSP ksp)
120: {
121: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
122: PetscReal res_norm,res,hapbnd,tt;
123: Scalar tmp;
124: int ierr,it = 0, max_k = gmres->max_k,max_it = ksp->max_it;
125: PetscTruth hapend = PETSC_FALSE;
128: ierr = VecNorm(VEC_VV(0),NORM_2,&res_norm);
129: res = res_norm;
130: *GRS(0) = res_norm;
132: /* check for the convergence */
133: if (!res) {
134: if (itcount) *itcount = 0;
135: ksp->reason = KSP_CONVERGED_ATOL;
136: PetscLogInfo(ksp,"GMRESCycle: Converged due to zero residual norm on entryn");
137: return(0);
138: }
140: /* scale VEC_VV (the initial residual) */
141: tmp = 1.0/res_norm; VecScale(&tmp,VEC_VV(0));
143: PetscObjectTakeAccess(ksp);
144: ksp->rnorm = res;
145: PetscObjectGrantAccess(ksp);
146: gmres->it = (it - 1);
147: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
148: while (!ksp->reason && it < max_k && ksp->its < max_it) {
149: KSPLogResidualHistory(ksp,res);
150: gmres->it = (it - 1);
151: KSPMonitor(ksp,ksp->its,res);
152: if (gmres->vv_allocated <= it + VEC_OFFSET + 1) {
153: GMRESGetNewVectors(ksp,it+1);
154: }
155: KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,VEC_VV(it),VEC_VV(1+it),VEC_TEMP_MATOP);
157: /* update hessenberg matrix and do Gram-Schmidt */
158: (*gmres->orthog)(ksp,it);
160: /* vv(i+1) . vv(i+1) */
161: VecNorm(VEC_VV(it+1),NORM_2,&tt);
162: /* save the magnitude */
163: *HH(it+1,it) = tt;
164: *HES(it+1,it) = tt;
166: /* check for the happy breakdown */
167: hapbnd = PetscAbsScalar(tt / *GRS(it));
168: if (hapbnd > gmres->haptol) hapbnd = gmres->haptol;
169: if (tt > hapbnd) {
170: tmp = 1.0/tt; VecScale(&tmp,VEC_VV(it+1));
171: } else {
172: PetscLogInfo(ksp,"Detected happy breakdown, current hapbnd = %g tt = %gn",hapbnd,tt);
173: hapend = PETSC_TRUE;
174: }
175: GMRESUpdateHessenberg(ksp,it,hapend,&res);
176: it++;
177: gmres->it = (it-1); /* For converged */
178: PetscObjectTakeAccess(ksp);
179: ksp->its++;
180: ksp->rnorm = res;
181: PetscObjectGrantAccess(ksp);
183: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
185: /* Catch error in happy breakdown and signal convergence and break from loop */
186: if (hapend) {
187: if (!ksp->reason) {
188: SETERRQ1(0,"You reached the happy break down, but convergence was not indicated. Residual norm = %",res);
189: }
190: break;
191: }
192: }
193: KSPLogResidualHistory(ksp,res);
195: /*
196: Monitor if we know that we will not return for a restart */
197: if (ksp->reason || ksp->its >= max_it) {
198: KSPMonitor(ksp,ksp->its,res);
199: }
201: if (itcount) *itcount = it;
204: /*
205: Down here we have to solve for the "best" coefficients of the Krylov
206: columns, add the solution values together, and possibly unwind the
207: preconditioning from the solution
208: */
209: /* Form the solution (or the solution so far) */
210: BuildGmresSoln(GRS(0),VEC_SOLN,VEC_SOLN,ksp,it-1);
212: return(0);
213: }
215: int KSPSolve_GMRES(KSP ksp,int *outits)
216: {
217: int ierr,its,itcount;
218: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
219: PetscTruth guess_zero = ksp->guess_zero;
222: if (ksp->calc_sings && !gmres->Rsvd) {
223: SETERRQ(1,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
224: }
226: ierr = PetscObjectTakeAccess(ksp);
227: ksp->its = 0;
228: ierr = PetscObjectGrantAccess(ksp);
230: itcount = 0;
231: ksp->reason = KSP_CONVERGED_ITERATING;
232: while (!ksp->reason) {
233: ierr = KSPInitialResidual(ksp,VEC_SOLN,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),VEC_BINVF,VEC_RHS);
234: ierr = GMREScycle(&its,ksp);
235: itcount += its;
236: if (itcount >= ksp->max_it) {
237: ksp->reason = KSP_DIVERGED_ITS;
238: break;
239: }
240: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
241: }
242: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
243: *outits = itcount;
244: return(0);
245: }
247: int KSPDestroy_GMRES(KSP ksp)
248: {
249: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
250: int i,ierr;
253: /* Free the Hessenberg matrix */
254: if (gmres->hh_origin) {PetscFree(gmres->hh_origin);}
256: /* Free the pointer to user variables */
257: if (gmres->vecs) {PetscFree(gmres->vecs);}
259: /* free work vectors */
260: for (i=0; i<gmres->nwork_alloc; i++) {
261: VecDestroyVecs(gmres->user_work[i],gmres->mwork_alloc[i]);
262: }
263: if (gmres->user_work) {PetscFree(gmres->user_work);}
264: if (gmres->mwork_alloc) {PetscFree(gmres->mwork_alloc);}
265: if (gmres->nrs) {PetscFree(gmres->nrs);}
266: if (gmres->sol_temp) {VecDestroy(gmres->sol_temp);}
267: if (gmres->Rsvd) {PetscFree(gmres->Rsvd);}
268: if (gmres->Dsvd) {PetscFree(gmres->Dsvd);}
269: PetscFree(gmres);
271: return(0);
272: }
273: /*
274: BuildGmresSoln - create the solution from the starting vector and the
275: current iterates.
277: Input parameters:
278: nrs - work area of size it + 1.
279: vs - index of initial guess
280: vdest - index of result. Note that vs may == vdest (replace
281: guess with the solution).
283: This is an internal routine that knows about the GMRES internals.
284: */
285: static int BuildGmresSoln(Scalar* nrs,Vec vs,Vec vdest,KSP ksp,int it)
286: {
287: Scalar tt,zero = 0.0,one = 1.0;
288: int ierr,ii,k,j;
289: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
292: /* Solve for solution vector that minimizes the residual */
294: /* If it is < 0, no gmres steps have been performed */
295: if (it < 0) {
296: if (vdest != vs) {
297: VecCopy(vs,vdest);
298: }
299: return(0);
300: }
301: if (*HH(it,it) == 0.0) SETERRQ2(1,"HH(it,it) is identically zero; it = %d GRS(it) = %g",it,PetscAbsScalar(*GRS(it)));
302: if (*HH(it,it) != 0.0) {
303: nrs[it] = *GRS(it) / *HH(it,it);
304: } else {
305: nrs[it] = 0.0;
306: }
307: for (ii=1; ii<=it; ii++) {
308: k = it - ii;
309: tt = *GRS(k);
310: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
311: nrs[k] = tt / *HH(k,k);
312: }
314: /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
315: VecSet(&zero,VEC_TEMP);
316: VecMAXPY(it+1,nrs,VEC_TEMP,&VEC_VV(0));
318: KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
319: /* add solution to previous solution */
320: if (vdest != vs) {
321: VecCopy(vs,vdest);
322: }
323: VecAXPY(&one,VEC_TEMP,vdest);
324: return(0);
325: }
326: /*
327: Do the scalar work for the orthogonalization. Return new residual.
328: */
329: static int GMRESUpdateHessenberg(KSP ksp,int it,PetscTruth hapend,PetscReal *res)
330: {
331: Scalar *hh,*cc,*ss,tt;
332: int j;
333: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
336: hh = HH(0,it);
337: cc = CC(0);
338: ss = SS(0);
340: /* Apply all the previously computed plane rotations to the new column
341: of the Hessenberg matrix */
342: for (j=1; j<=it; j++) {
343: tt = *hh;
344: #if defined(PETSC_USE_COMPLEX)
345: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
346: #else
347: *hh = *cc * tt + *ss * *(hh+1);
348: #endif
349: hh++;
350: *hh = *cc++ * *hh - (*ss++ * tt);
351: }
353: /*
354: compute the new plane rotation, and apply it to:
355: 1) the right-hand-side of the Hessenberg system
356: 2) the new column of the Hessenberg matrix
357: thus obtaining the updated value of the residual
358: */
359: if (!hapend) {
360: #if defined(PETSC_USE_COMPLEX)
361: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
362: #else
363: tt = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
364: #endif
365: if (tt == 0.0) {SETERRQ(PETSC_ERR_KSP_BRKDWN,"Your matrix or preconditioner is the null operator");}
366: *cc = *hh / tt;
367: *ss = *(hh+1) / tt;
368: *GRS(it+1) = - (*ss * *GRS(it));
369: #if defined(PETSC_USE_COMPLEX)
370: *GRS(it) = PetscConj(*cc) * *GRS(it);
371: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
372: #else
373: *GRS(it) = *cc * *GRS(it);
374: *hh = *cc * *hh + *ss * *(hh+1);
375: #endif
376: *res = PetscAbsScalar(*GRS(it+1));
377: } else {
378: /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
379: another rotation matrix (so RH doesn't change). The new residual is
380: always the new sine term times the residual from last time (GRS(it)),
381: but now the new sine rotation would be zero...so the residual should
382: be zero...so we will multiply "zero" by the last residual. This might
383: not be exactly what we want to do here -could just return "zero". */
384:
385: *res = 0.0;
386: }
387: return(0);
388: }
389: /*
390: This routine allocates more work vectors, starting from VEC_VV(it).
391: */
392: static int GMRESGetNewVectors(KSP ksp,int it)
393: {
394: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
395: int nwork = gmres->nwork_alloc,k,nalloc,ierr;
398: nalloc = gmres->delta_allocate;
399: /* Adjust the number to allocate to make sure that we don't exceed the
400: number of available slots */
401: if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated){
402: nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
403: }
404: if (!nalloc) return(0);
406: gmres->vv_allocated += nalloc;
407: VecDuplicateVecs(ksp->vec_rhs,nalloc,&gmres->user_work[nwork]);
408: PetscLogObjectParents(ksp,nalloc,gmres->user_work[nwork]);
409: gmres->mwork_alloc[nwork] = nalloc;
410: for (k=0; k<nalloc; k++) {
411: gmres->vecs[it+VEC_OFFSET+k] = gmres->user_work[nwork][k];
412: }
413: gmres->nwork_alloc++;
414: return(0);
415: }
417: int KSPBuildSolution_GMRES(KSP ksp,Vec ptr,Vec *result)
418: {
419: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
420: int ierr;
423: if (!ptr) {
424: if (!gmres->sol_temp) {
425: VecDuplicate(ksp->vec_sol,&gmres->sol_temp);
426: PetscLogObjectParent(ksp,gmres->sol_temp);
427: }
428: ptr = gmres->sol_temp;
429: }
430: if (!gmres->nrs) {
431: /* allocate the work area */
432: PetscMalloc(gmres->max_k*sizeof(Scalar),&gmres->nrs);
433: PetscLogObjectMemory(ksp,gmres->max_k*sizeof(Scalar));
434: }
436: BuildGmresSoln(gmres->nrs,VEC_SOLN,ptr,ksp,gmres->it);
437: *result = ptr;
438: return(0);
439: }
441: int KSPView_GMRES(KSP ksp,PetscViewer viewer)
442: {
443: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
444: char *cstr;
445: int ierr;
446: PetscTruth isascii,isstring;
449: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
450: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
451: if (gmres->orthog == KSPGMRESUnmodifiedGramSchmidtOrthogonalization) {
452: cstr = "Unmodified Gram-Schmidt Orthogonalization";
453: } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
454: cstr = "Modified Gram-Schmidt Orthogonalization";
455: } else if (gmres->orthog == KSPGMRESIROrthogonalization) {
456: cstr = "Unmodified Gram-Schmidt + 1 step Iterative Refinement Orthogonalization";
457: } else {
458: cstr = "unknown orthogonalization";
459: }
460: if (isascii) {
461: PetscViewerASCIIPrintf(viewer," GMRES: restart=%d, using %sn",gmres->max_k,cstr);
462: PetscViewerASCIIPrintf(viewer," GMRES: happy breakdown tolerance %gn",gmres->haptol);
463: } else if (isstring) {
464: PetscViewerStringSPrintf(viewer,"%s restart %d",cstr,gmres->max_k);
465: } else {
466: SETERRQ1(1,"Viewer type %s not supported for KSP GMRES",((PetscObject)viewer)->type_name);
467: }
468: return(0);
469: }
471: /*@C
472: KSPGMRESKrylovMonitor - Calls VecView() for each direction in the
473: GMRES accumulated Krylov space.
475: Collective on KSP
477: Input Parameters:
478: + ksp - the KSP context
479: . its - iteration number
480: . fgnorm - 2-norm of residual (or gradient)
481: - a viewers object created with PetscViewersCreate()
483: Level: intermediate
485: .keywords: KSP, nonlinear, vector, monitor, view, Krylov space
487: .seealso: KSPSetMonitor(), KSPDefaultMonitor(), VecView(), PetscViewersCreate(), PetscViewersDestroy()
488: @*/
489: int KSPGMRESKrylovMonitor(KSP ksp,int its,PetscReal fgnorm,void *dummy)
490: {
491: PetscViewers viewers = (PetscViewers)dummy;
492: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
493: int ierr;
494: Vec x;
495: PetscViewer viewer;
498: PetscViewersGetViewer(viewers,gmres->it+1,&viewer);
499: PetscViewerSetType(viewer,PETSC_VIEWER_DRAW);
501: x = VEC_VV(gmres->it+1);
502: ierr = VecView(x,viewer);
504: return(0);
505: }
507: int KSPSetFromOptions_GMRES(KSP ksp)
508: {
509: int ierr,restart;
510: double haptol;
511: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
512: PetscTruth flg;
515: PetscOptionsHead("KSP GMRES Options");
516: PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
517: if (flg) { KSPGMRESSetRestart(ksp,restart); }
518: PetscOptionsDouble("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
519: if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
520: PetscOptionsName("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",&flg);
521: if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
522: PetscOptionsLogicalGroupBegin("-ksp_gmres_unmodifiedgramschmidt","Classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
523: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESUnmodifiedGramSchmidtOrthogonalization);}
524: PetscOptionsLogicalGroup("-ksp_gmres_modifiedgramschmidt","Modified Gram-Schmidt (slow,more stable)","KSPGMRESSetOrthogonalization",&flg);
525: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
526: PetscOptionsLogicalGroupEnd("-ksp_gmres_irorthog","Classical Gram-Schmidt + iterative refinement","KSPGMRESSetOrthogonalization",&flg);
527: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESIROrthogonalization);}
528: PetscOptionsName("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPSetMonitor",&flg);
529: if (flg) {
530: PetscViewers viewers;
531: PetscViewersCreate(ksp->comm,&viewers);
532: KSPSetMonitor(ksp,KSPGMRESKrylovMonitor,viewers,(int (*)(void*))PetscViewersDestroy);
533: }
534: PetscOptionsTail();
535: return(0);
536: }
538: EXTERN int KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
539: EXTERN int KSPComputeEigenvalues_GMRES(KSP,int,PetscReal *,PetscReal *,int *);
542: EXTERN_C_BEGIN
543: int KSPGMRESSetHapTol_GMRES(KSP ksp,double tol)
544: {
545: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
548: gmres->haptol = tol;
549: return(0);
550: }
551: EXTERN_C_END
553: EXTERN_C_BEGIN
554: int KSPGMRESSetRestart_GMRES(KSP ksp,int max_k)
555: {
556: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
559: gmres->max_k = max_k;
560: return(0);
561: }
562: EXTERN_C_END
564: EXTERN_C_BEGIN
565: int KSPGMRESSetOrthogonalization_GMRES(KSP ksp,int (*fcn)(KSP,int))
566: {
569: ((KSP_GMRES *)ksp->data)->orthog = fcn;
570: return(0);
571: }
572: EXTERN_C_END
574: EXTERN_C_BEGIN
575: int KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
576: {
577: KSP_GMRES *gmres;
580: gmres = (KSP_GMRES *)ksp->data;
581: gmres->q_preallocate = 1;
582: return(0);
583: }
584: EXTERN_C_END
586: EXTERN_C_BEGIN
587: int KSPCreate_GMRES(KSP ksp)
588: {
589: KSP_GMRES *gmres;
590: int ierr;
593: PetscNew(KSP_GMRES,&gmres);
594: ierr = PetscMemzero(gmres,sizeof(KSP_GMRES));
595: PetscLogObjectMemory(ksp,sizeof(KSP_GMRES));
596: ksp->data = (void*)gmres;
597: ksp->ops->buildsolution = KSPBuildSolution_GMRES;
599: ksp->ops->setup = KSPSetUp_GMRES;
600: ksp->ops->solve = KSPSolve_GMRES;
601: ksp->ops->destroy = KSPDestroy_GMRES;
602: ksp->ops->view = KSPView_GMRES;
603: ksp->ops->setfromoptions = KSPSetFromOptions_GMRES;
604: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
605: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
607: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
608: "KSPGMRESSetPreAllocateVectors_GMRES",
609: KSPGMRESSetPreAllocateVectors_GMRES);
610: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
611: "KSPGMRESSetOrthogonalization_GMRES",
612: KSPGMRESSetOrthogonalization_GMRES);
613: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
614: "KSPGMRESSetRestart_GMRES",
615: KSPGMRESSetRestart_GMRES);
616: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
617: "KSPGMRESSetHapTol_GMRES",
618: KSPGMRESSetHapTol_GMRES);
620: gmres->haptol = 1.0e-30;
621: gmres->q_preallocate = 0;
622: gmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
623: gmres->orthog = KSPGMRESUnmodifiedGramSchmidtOrthogonalization;
624: gmres->nrs = 0;
625: gmres->sol_temp = 0;
626: gmres->max_k = GMRES_DEFAULT_MAXK;
627: gmres->Rsvd = 0;
628: return(0);
629: }
630: EXTERN_C_END