Actual source code: gmres.c
1: /*$Id: gmres.c,v 1.176 2001/08/07 03:03:51 balay 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 eliminate an extra application of B^-1 between
28: restarts as long as we don't require that the solution at the end
29: of an unsuccessful gmres iteration always be the solution x.
30: */
32: #include src/sles/ksp/impls/gmres/gmresp.h
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(PetscScalar*,Vec,Vec,KSP,int);
39: #undef __FUNCT__
41: int KSPSetUp_GMRES(KSP ksp)
42: {
43: unsigned int size,hh,hes,rs,cc;
44: int ierr,max_k,k;
45: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
48: if (ksp->pc_side == PC_SYMMETRIC) {
49: SETERRQ(2,"no symmetric preconditioning for KSPGMRES");
50: } else if (ksp->pc_side == PC_RIGHT) {
51: SETERRQ(2,"no right preconditioning for KSPGMRES");
52: }
54: max_k = gmres->max_k;
55: hh = (max_k + 2) * (max_k + 1);
56: hes = (max_k + 1) * (max_k + 1);
57: rs = (max_k + 2);
58: cc = (max_k + 1);
59: size = (hh + hes + rs + 2*cc) * sizeof(PetscScalar);
61: PetscMalloc(size,&gmres->hh_origin);
62: PetscMemzero(gmres->hh_origin,size);
63: PetscLogObjectMemory(ksp,size);
64: gmres->hes_origin = gmres->hh_origin + hh;
65: gmres->rs_origin = gmres->hes_origin + hes;
66: gmres->cc_origin = gmres->rs_origin + rs;
67: gmres->ss_origin = gmres->cc_origin + cc;
69: if (ksp->calc_sings) {
70: /* Allocate workspace to hold Hessenberg matrix needed by Eispack */
71: size = (max_k + 3)*(max_k + 9)*sizeof(PetscScalar);
72: PetscMalloc(size,&gmres->Rsvd);
73: PetscMalloc(5*(max_k+2)*sizeof(PetscReal),&gmres->Dsvd);
74: PetscLogObjectMemory(ksp,size+5*(max_k+2)*sizeof(PetscReal));
75: }
77: /* Allocate array to hold pointers to user vectors. Note that we need
78: 4 + max_k + 1 (since we need it+1 vectors, and it <= max_k) */
79: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void *),&gmres->vecs);
80: gmres->vecs_allocated = VEC_OFFSET + 2 + max_k;
81: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(void *),&gmres->user_work);
82: PetscMalloc((VEC_OFFSET+2+max_k)*sizeof(int),&gmres->mwork_alloc);
83: PetscLogObjectMemory(ksp,(VEC_OFFSET+2+max_k)*(2*sizeof(void *)+sizeof(int)));
85: if (gmres->q_preallocate) {
86: gmres->vv_allocated = VEC_OFFSET + 2 + max_k;
87: VecDuplicateVecs(VEC_RHS,gmres->vv_allocated,&gmres->user_work[0]);
88: PetscLogObjectParents(ksp,gmres->vv_allocated,gmres->user_work[0]);
89: gmres->mwork_alloc[0] = gmres->vv_allocated;
90: gmres->nwork_alloc = 1;
91: for (k=0; k<gmres->vv_allocated; k++) {
92: gmres->vecs[k] = gmres->user_work[0][k];
93: }
94: } else {
95: gmres->vv_allocated = 5;
96: VecDuplicateVecs(ksp->vec_rhs,5,&gmres->user_work[0]);
97: PetscLogObjectParents(ksp,5,gmres->user_work[0]);
98: gmres->mwork_alloc[0] = 5;
99: gmres->nwork_alloc = 1;
100: for (k=0; k<gmres->vv_allocated; k++) {
101: gmres->vecs[k] = gmres->user_work[0][k];
102: }
103: }
104: return(0);
105: }
107: /*
108: Run gmres, possibly with restart. Return residual history if requested.
109: input parameters:
111: . gmres - structure containing parameters and work areas
113: output parameters:
114: . nres - residuals (from preconditioned system) at each step.
115: If restarting, consider passing nres+it. If null,
116: ignored
117: . itcount - number of iterations used. nres[0] to nres[itcount]
118: are defined. If null, ignored.
119:
120: Notes:
121: On entry, the value in vector VEC_VV(0) should be the initial residual
122: (this allows shortcuts where the initial preconditioned residual is 0).
123: */
124: #undef __FUNCT__
126: int GMREScycle(int *itcount,KSP ksp)
127: {
128: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
129: PetscReal res_norm,res,hapbnd,tt;
130: PetscScalar tmp;
131: int ierr,it = 0, max_k = gmres->max_k,max_it = ksp->max_it;
132: PetscTruth hapend = PETSC_FALSE;
135: ierr = VecNorm(VEC_VV(0),NORM_2,&res_norm);
136: res = res_norm;
137: *GRS(0) = res_norm;
139: /* check for the convergence */
140: if (!res) {
141: if (itcount) *itcount = 0;
142: ksp->reason = KSP_CONVERGED_ATOL;
143: PetscLogInfo(ksp,"GMRESCycle: Converged due to zero residual norm on entryn");
144: return(0);
145: }
147: /* scale VEC_VV (the initial residual) */
148: tmp = 1.0/res_norm; VecScale(&tmp,VEC_VV(0));
150: PetscObjectTakeAccess(ksp);
151: ksp->rnorm = res;
152: PetscObjectGrantAccess(ksp);
153: gmres->it = (it - 1);
154: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
155: while (!ksp->reason && it < max_k && ksp->its < max_it) {
156: KSPLogResidualHistory(ksp,res);
157: gmres->it = (it - 1);
158: KSPMonitor(ksp,ksp->its,res);
159: if (gmres->vv_allocated <= it + VEC_OFFSET + 1) {
160: GMRESGetNewVectors(ksp,it+1);
161: }
162: KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,VEC_VV(it),VEC_VV(1+it),VEC_TEMP_MATOP);
164: /* update hessenberg matrix and do Gram-Schmidt */
165: (*gmres->orthog)(ksp,it);
167: /* vv(i+1) . vv(i+1) */
168: VecNorm(VEC_VV(it+1),NORM_2,&tt);
169: /* save the magnitude */
170: *HH(it+1,it) = tt;
171: *HES(it+1,it) = tt;
173: /* check for the happy breakdown */
174: hapbnd = PetscAbsScalar(tt / *GRS(it));
175: if (hapbnd > gmres->haptol) hapbnd = gmres->haptol;
176: if (tt > hapbnd) {
177: tmp = 1.0/tt; VecScale(&tmp,VEC_VV(it+1));
178: } else {
179: PetscLogInfo(ksp,"Detected happy breakdown, current hapbnd = %g tt = %gn",hapbnd,tt);
180: hapend = PETSC_TRUE;
181: }
182: GMRESUpdateHessenberg(ksp,it,hapend,&res);
183: it++;
184: gmres->it = (it-1); /* For converged */
185: PetscObjectTakeAccess(ksp);
186: ksp->its++;
187: ksp->rnorm = res;
188: PetscObjectGrantAccess(ksp);
190: (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
192: /* Catch error in happy breakdown and signal convergence and break from loop */
193: if (hapend) {
194: if (!ksp->reason) {
195: SETERRQ1(0,"You reached the happy break down, but convergence was not indicated. Residual norm = %g",res);
196: }
197: break;
198: }
199: }
200: KSPLogResidualHistory(ksp,res);
202: /*
203: Monitor if we know that we will not return for a restart */
204: if (ksp->reason || ksp->its >= max_it) {
205: KSPMonitor(ksp,ksp->its,res);
206: }
208: if (itcount) *itcount = it;
211: /*
212: Down here we have to solve for the "best" coefficients of the Krylov
213: columns, add the solution values together, and possibly unwind the
214: preconditioning from the solution
215: */
216: /* Form the solution (or the solution so far) */
217: BuildGmresSoln(GRS(0),VEC_SOLN,VEC_SOLN,ksp,it-1);
219: return(0);
220: }
222: #undef __FUNCT__
224: int KSPSolve_GMRES(KSP ksp,int *outits)
225: {
226: int ierr,its,itcount;
227: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
228: PetscTruth guess_zero = ksp->guess_zero;
231: if (ksp->calc_sings && !gmres->Rsvd) {
232: SETERRQ(1,"Must call KSPSetComputeSingularValues() before KSPSetUp() is called");
233: }
235: ierr = PetscObjectTakeAccess(ksp);
236: ksp->its = 0;
237: ierr = PetscObjectGrantAccess(ksp);
239: itcount = 0;
240: ksp->reason = KSP_CONVERGED_ITERATING;
241: while (!ksp->reason) {
242: ierr = KSPInitialResidual(ksp,VEC_SOLN,VEC_TEMP,VEC_TEMP_MATOP,VEC_VV(0),VEC_RHS);
243: ierr = GMREScycle(&its,ksp);
244: itcount += its;
245: if (itcount >= ksp->max_it) {
246: ksp->reason = KSP_DIVERGED_ITS;
247: break;
248: }
249: ksp->guess_zero = PETSC_FALSE; /* every future call to KSPInitialResidual() will have nonzero guess */
250: }
251: ksp->guess_zero = guess_zero; /* restore if user provided nonzero initial guess */
252: if (outits) *outits = itcount;
253: return(0);
254: }
256: #undef __FUNCT__
258: int KSPDestroy_GMRES_Internal(KSP ksp)
259: {
260: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
261: int i,ierr;
264: /* Free the Hessenberg matrix */
265: if (gmres->hh_origin) {PetscFree(gmres->hh_origin);}
267: /* Free the pointer to user variables */
268: if (gmres->vecs) {PetscFree(gmres->vecs);}
270: /* free work vectors */
271: for (i=0; i<gmres->nwork_alloc; i++) {
272: VecDestroyVecs(gmres->user_work[i],gmres->mwork_alloc[i]);
273: }
274: if (gmres->user_work) {PetscFree(gmres->user_work);}
275: if (gmres->mwork_alloc) {PetscFree(gmres->mwork_alloc);}
276: if (gmres->nrs) {PetscFree(gmres->nrs);}
277: if (gmres->sol_temp) {VecDestroy(gmres->sol_temp);}
278: if (gmres->Rsvd) {PetscFree(gmres->Rsvd);}
279: if (gmres->Dsvd) {PetscFree(gmres->Dsvd);}
281: return(0);
282: }
284: #undef __FUNCT__
286: int KSPDestroy_GMRES(KSP ksp)
287: {
288: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
289: int ierr;
292: KSPDestroy_GMRES_Internal(ksp);
293: PetscFree(gmres);
294: return(0);
295: }
296: /*
297: BuildGmresSoln - create the solution from the starting vector and the
298: current iterates.
300: Input parameters:
301: nrs - work area of size it + 1.
302: vs - index of initial guess
303: vdest - index of result. Note that vs may == vdest (replace
304: guess with the solution).
306: This is an internal routine that knows about the GMRES internals.
307: */
308: #undef __FUNCT__
310: static int BuildGmresSoln(PetscScalar* nrs,Vec vs,Vec vdest,KSP ksp,int it)
311: {
312: PetscScalar tt,zero = 0.0,one = 1.0;
313: int ierr,ii,k,j;
314: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
317: /* Solve for solution vector that minimizes the residual */
319: /* If it is < 0, no gmres steps have been performed */
320: if (it < 0) {
321: if (vdest != vs) {
322: VecCopy(vs,vdest);
323: }
324: return(0);
325: }
326: if (*HH(it,it) == 0.0) SETERRQ2(1,"HH(it,it) is identically zero; it = %d GRS(it) = %g",it,PetscAbsScalar(*GRS(it)));
327: if (*HH(it,it) != 0.0) {
328: nrs[it] = *GRS(it) / *HH(it,it);
329: } else {
330: nrs[it] = 0.0;
331: }
332: for (ii=1; ii<=it; ii++) {
333: k = it - ii;
334: tt = *GRS(k);
335: for (j=k+1; j<=it; j++) tt = tt - *HH(k,j) * nrs[j];
336: nrs[k] = tt / *HH(k,k);
337: }
339: /* Accumulate the correction to the solution of the preconditioned problem in TEMP */
340: VecSet(&zero,VEC_TEMP);
341: VecMAXPY(it+1,nrs,VEC_TEMP,&VEC_VV(0));
343: KSPUnwindPreconditioner(ksp,VEC_TEMP,VEC_TEMP_MATOP);
344: /* add solution to previous solution */
345: if (vdest != vs) {
346: VecCopy(vs,vdest);
347: }
348: VecAXPY(&one,VEC_TEMP,vdest);
349: return(0);
350: }
351: /*
352: Do the scalar work for the orthogonalization. Return new residual.
353: */
354: #undef __FUNCT__
356: static int GMRESUpdateHessenberg(KSP ksp,int it,PetscTruth hapend,PetscReal *res)
357: {
358: PetscScalar *hh,*cc,*ss,tt;
359: int j;
360: KSP_GMRES *gmres = (KSP_GMRES *)(ksp->data);
363: hh = HH(0,it);
364: cc = CC(0);
365: ss = SS(0);
367: /* Apply all the previously computed plane rotations to the new column
368: of the Hessenberg matrix */
369: for (j=1; j<=it; j++) {
370: tt = *hh;
371: #if defined(PETSC_USE_COMPLEX)
372: *hh = PetscConj(*cc) * tt + *ss * *(hh+1);
373: #else
374: *hh = *cc * tt + *ss * *(hh+1);
375: #endif
376: hh++;
377: *hh = *cc++ * *hh - (*ss++ * tt);
378: }
380: /*
381: compute the new plane rotation, and apply it to:
382: 1) the right-hand-side of the Hessenberg system
383: 2) the new column of the Hessenberg matrix
384: thus obtaining the updated value of the residual
385: */
386: if (!hapend) {
387: #if defined(PETSC_USE_COMPLEX)
388: tt = PetscSqrtScalar(PetscConj(*hh) * *hh + PetscConj(*(hh+1)) * *(hh+1));
389: #else
390: tt = PetscSqrtScalar(*hh * *hh + *(hh+1) * *(hh+1));
391: #endif
392: if (tt == 0.0) {SETERRQ(PETSC_ERR_KSP_BRKDWN,"Your matrix or preconditioner is the null operator");}
393: *cc = *hh / tt;
394: *ss = *(hh+1) / tt;
395: *GRS(it+1) = - (*ss * *GRS(it));
396: #if defined(PETSC_USE_COMPLEX)
397: *GRS(it) = PetscConj(*cc) * *GRS(it);
398: *hh = PetscConj(*cc) * *hh + *ss * *(hh+1);
399: #else
400: *GRS(it) = *cc * *GRS(it);
401: *hh = *cc * *hh + *ss * *(hh+1);
402: #endif
403: *res = PetscAbsScalar(*GRS(it+1));
404: } else {
405: /* happy breakdown: HH(it+1, it) = 0, therfore we don't need to apply
406: another rotation matrix (so RH doesn't change). The new residual is
407: always the new sine term times the residual from last time (GRS(it)),
408: but now the new sine rotation would be zero...so the residual should
409: be zero...so we will multiply "zero" by the last residual. This might
410: not be exactly what we want to do here -could just return "zero". */
411:
412: *res = 0.0;
413: }
414: return(0);
415: }
416: /*
417: This routine allocates more work vectors, starting from VEC_VV(it).
418: */
419: #undef __FUNCT__
421: static int GMRESGetNewVectors(KSP ksp,int it)
422: {
423: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
424: int nwork = gmres->nwork_alloc,k,nalloc,ierr;
427: nalloc = gmres->delta_allocate;
428: /* Adjust the number to allocate to make sure that we don't exceed the
429: number of available slots */
430: if (it + VEC_OFFSET + nalloc >= gmres->vecs_allocated){
431: nalloc = gmres->vecs_allocated - it - VEC_OFFSET;
432: }
433: if (!nalloc) return(0);
435: gmres->vv_allocated += nalloc;
436: VecDuplicateVecs(ksp->vec_rhs,nalloc,&gmres->user_work[nwork]);
437: PetscLogObjectParents(ksp,nalloc,gmres->user_work[nwork]);
438: gmres->mwork_alloc[nwork] = nalloc;
439: for (k=0; k<nalloc; k++) {
440: gmres->vecs[it+VEC_OFFSET+k] = gmres->user_work[nwork][k];
441: }
442: gmres->nwork_alloc++;
443: return(0);
444: }
446: #undef __FUNCT__
448: int KSPBuildSolution_GMRES(KSP ksp,Vec ptr,Vec *result)
449: {
450: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
451: int ierr;
454: if (!ptr) {
455: if (!gmres->sol_temp) {
456: VecDuplicate(ksp->vec_sol,&gmres->sol_temp);
457: PetscLogObjectParent(ksp,gmres->sol_temp);
458: }
459: ptr = gmres->sol_temp;
460: }
461: if (!gmres->nrs) {
462: /* allocate the work area */
463: PetscMalloc(gmres->max_k*sizeof(PetscScalar),&gmres->nrs);
464: PetscLogObjectMemory(ksp,gmres->max_k*sizeof(PetscScalar));
465: }
467: BuildGmresSoln(gmres->nrs,VEC_SOLN,ptr,ksp,gmres->it);
468: *result = ptr;
469: return(0);
470: }
472: #undef __FUNCT__
474: int KSPView_GMRES(KSP ksp,PetscViewer viewer)
475: {
476: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
477: char *cstr;
478: int ierr;
479: PetscTruth isascii,isstring;
482: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
483: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
484: if (gmres->orthog == KSPGMRESUnmodifiedGramSchmidtOrthogonalization) {
485: cstr = "Unmodified Gram-Schmidt Orthogonalization";
486: } else if (gmres->orthog == KSPGMRESModifiedGramSchmidtOrthogonalization) {
487: cstr = "Modified Gram-Schmidt Orthogonalization";
488: } else if (gmres->orthog == KSPGMRESIROrthogonalization) {
489: cstr = "Unmodified Gram-Schmidt + 1 step Iterative Refinement Orthogonalization";
490: } else {
491: cstr = "unknown orthogonalization";
492: }
493: if (isascii) {
494: PetscViewerASCIIPrintf(viewer," GMRES: restart=%d, using %sn",gmres->max_k,cstr);
495: PetscViewerASCIIPrintf(viewer," GMRES: happy breakdown tolerance %gn",gmres->haptol);
496: } else if (isstring) {
497: PetscViewerStringSPrintf(viewer,"%s restart %d",cstr,gmres->max_k);
498: } else {
499: SETERRQ1(1,"Viewer type %s not supported for KSP GMRES",((PetscObject)viewer)->type_name);
500: }
501: return(0);
502: }
504: #undef __FUNCT__
506: /*@C
507: KSPGMRESKrylovMonitor - Calls VecView() for each direction in the
508: GMRES accumulated Krylov space.
510: Collective on KSP
512: Input Parameters:
513: + ksp - the KSP context
514: . its - iteration number
515: . fgnorm - 2-norm of residual (or gradient)
516: - a viewers object created with PetscViewersCreate()
518: Level: intermediate
520: .keywords: KSP, nonlinear, vector, monitor, view, Krylov space
522: .seealso: KSPSetMonitor(), KSPDefaultMonitor(), VecView(), PetscViewersCreate(), PetscViewersDestroy()
523: @*/
524: int KSPGMRESKrylovMonitor(KSP ksp,int its,PetscReal fgnorm,void *dummy)
525: {
526: PetscViewers viewers = (PetscViewers)dummy;
527: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
528: int ierr;
529: Vec x;
530: PetscViewer viewer;
533: PetscViewersGetViewer(viewers,gmres->it+1,&viewer);
534: PetscViewerSetType(viewer,PETSC_VIEWER_DRAW);
536: x = VEC_VV(gmres->it+1);
537: ierr = VecView(x,viewer);
539: return(0);
540: }
542: #undef __FUNCT__
544: int KSPSetFromOptions_GMRES(KSP ksp)
545: {
546: int ierr,restart;
547: PetscReal haptol;
548: KSP_GMRES *gmres = (KSP_GMRES*)ksp->data;
549: PetscTruth flg;
552: PetscOptionsHead("KSP GMRES Options");
553: PetscOptionsInt("-ksp_gmres_restart","Number of Krylov search directions","KSPGMRESSetRestart",gmres->max_k,&restart,&flg);
554: if (flg) { KSPGMRESSetRestart(ksp,restart); }
555: PetscOptionsReal("-ksp_gmres_haptol","Tolerance for exact convergence (happy ending)","KSPGMRESSetHapTol",gmres->haptol,&haptol,&flg);
556: if (flg) { KSPGMRESSetHapTol(ksp,haptol); }
557: PetscOptionsName("-ksp_gmres_preallocate","Preallocate Krylov vectors","KSPGMRESSetPreAllocateVectors",&flg);
558: if (flg) {KSPGMRESSetPreAllocateVectors(ksp);}
559: PetscOptionsLogicalGroupBegin("-ksp_gmres_unmodifiedgramschmidt","Classical (unmodified) Gram-Schmidt (fast)","KSPGMRESSetOrthogonalization",&flg);
560: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESUnmodifiedGramSchmidtOrthogonalization);}
561: PetscOptionsLogicalGroup("-ksp_gmres_modifiedgramschmidt","Modified Gram-Schmidt (slow,more stable)","KSPGMRESSetOrthogonalization",&flg);
562: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESModifiedGramSchmidtOrthogonalization);}
563: PetscOptionsLogicalGroupEnd("-ksp_gmres_irorthog","Classical Gram-Schmidt + iterative refinement","KSPGMRESSetOrthogonalization",&flg);
564: if (flg) {KSPGMRESSetOrthogonalization(ksp,KSPGMRESIROrthogonalization);}
565: PetscOptionsName("-ksp_gmres_krylov_monitor","Plot the Krylov directions","KSPSetMonitor",&flg);
566: if (flg) {
567: PetscViewers viewers;
568: PetscViewersCreate(ksp->comm,&viewers);
569: KSPSetMonitor(ksp,KSPGMRESKrylovMonitor,viewers,(int (*)(void*))PetscViewersDestroy);
570: }
571: PetscOptionsTail();
572: return(0);
573: }
575: EXTERN int KSPComputeExtremeSingularValues_GMRES(KSP,PetscReal *,PetscReal *);
576: EXTERN int KSPComputeEigenvalues_GMRES(KSP,int,PetscReal *,PetscReal *,int *);
579: EXTERN_C_BEGIN
580: #undef __FUNCT__
582: int KSPGMRESSetHapTol_GMRES(KSP ksp,double tol)
583: {
584: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
587: if (tol < 0.0) SETERRQ(1,"Tolerance must be non-negative");
588: gmres->haptol = tol;
589: return(0);
590: }
591: EXTERN_C_END
593: EXTERN_C_BEGIN
594: #undef __FUNCT__
596: int KSPGMRESSetRestart_GMRES(KSP ksp,int max_k)
597: {
598: KSP_GMRES *gmres = (KSP_GMRES *)ksp->data;
599: int ierr;
602: if (max_k < 1) SETERRQ(1,"Restart must be positive");
603: if (!ksp->setupcalled) {
604: gmres->max_k = max_k;
605: } else if (gmres->max_k != max_k) {
606: gmres->max_k = max_k;
607: ksp->setupcalled = 0;
608: /* free the data structures, then create them again */
609: KSPDestroy_GMRES_Internal(ksp);
610: }
612: return(0);
613: }
614: EXTERN_C_END
616: EXTERN_C_BEGIN
617: #undef __FUNCT__
619: int KSPGMRESSetOrthogonalization_GMRES(KSP ksp,int (*fcn)(KSP,int))
620: {
623: ((KSP_GMRES *)ksp->data)->orthog = fcn;
624: return(0);
625: }
626: EXTERN_C_END
628: EXTERN_C_BEGIN
629: #undef __FUNCT__
631: int KSPGMRESSetPreAllocateVectors_GMRES(KSP ksp)
632: {
633: KSP_GMRES *gmres;
636: gmres = (KSP_GMRES *)ksp->data;
637: gmres->q_preallocate = 1;
638: return(0);
639: }
640: EXTERN_C_END
642: EXTERN_C_BEGIN
643: #undef __FUNCT__
645: int KSPCreate_GMRES(KSP ksp)
646: {
647: KSP_GMRES *gmres;
648: int ierr;
651: PetscNew(KSP_GMRES,&gmres);
652: ierr = PetscMemzero(gmres,sizeof(KSP_GMRES));
653: PetscLogObjectMemory(ksp,sizeof(KSP_GMRES));
654: ksp->data = (void*)gmres;
655: ksp->ops->buildsolution = KSPBuildSolution_GMRES;
657: ksp->ops->setup = KSPSetUp_GMRES;
658: ksp->ops->solve = KSPSolve_GMRES;
659: ksp->ops->destroy = KSPDestroy_GMRES;
660: ksp->ops->view = KSPView_GMRES;
661: ksp->ops->setfromoptions = KSPSetFromOptions_GMRES;
662: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_GMRES;
663: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_GMRES;
665: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetPreAllocateVectors_C",
666: "KSPGMRESSetPreAllocateVectors_GMRES",
667: KSPGMRESSetPreAllocateVectors_GMRES);
668: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetOrthogonalization_C",
669: "KSPGMRESSetOrthogonalization_GMRES",
670: KSPGMRESSetOrthogonalization_GMRES);
671: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetRestart_C",
672: "KSPGMRESSetRestart_GMRES",
673: KSPGMRESSetRestart_GMRES);
674: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGMRESSetHapTol_C",
675: "KSPGMRESSetHapTol_GMRES",
676: KSPGMRESSetHapTol_GMRES);
678: gmres->haptol = 1.0e-30;
679: gmres->q_preallocate = 0;
680: gmres->delta_allocate = GMRES_DELTA_DIRECTIONS;
681: gmres->orthog = KSPGMRESUnmodifiedGramSchmidtOrthogonalization;
682: gmres->nrs = 0;
683: gmres->sol_temp = 0;
684: gmres->max_k = GMRES_DEFAULT_MAXK;
685: gmres->Rsvd = 0;
686: return(0);
687: }
688: EXTERN_C_END