Actual source code: lsqr.c
petsc-3.10.5 2019-03-28
2: /* lourens.vanzanen@shell.com contributed the standard error estimates of the solution, Jul 25, 2006 */
3: /* Bas van't Hof contributed the preconditioned aspects Feb 10, 2010 */
5: #define SWAP(a,b,c) { c = a; a = b; b = c; }
7: #include <petsc/private/kspimpl.h>
9: typedef struct {
10: PetscInt nwork_n,nwork_m;
11: Vec *vwork_m; /* work vectors of length m, where the system is size m x n */
12: Vec *vwork_n; /* work vectors of length n */
13: Vec se; /* Optional standard error vector */
14: PetscBool se_flg; /* flag for -ksp_lsqr_set_standard_error */
15: PetscBool exact_norm; /* flag for -ksp_lsqr_exact_mat_norm */
16: PetscReal arnorm; /* Good estimate of norm((A*inv(Pmat))'*r), where r = A*x - b, used in specific stopping criterion */
17: PetscReal anorm; /* Poor estimate of norm(A*inv(Pmat),'fro') used in specific stopping criterion */
18: } KSP_LSQR;
20: static PetscErrorCode VecSquare(Vec v)
21: {
23: PetscScalar *x;
24: PetscInt i, n;
27: VecGetLocalSize(v, &n);
28: VecGetArray(v, &x);
29: for (i = 0; i < n; i++) x[i] *= PetscConj(x[i]);
30: VecRestoreArray(v, &x);
31: return(0);
32: }
34: static PetscErrorCode KSPSetUp_LSQR(KSP ksp)
35: {
37: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
38: PetscBool nopreconditioner;
41: PetscObjectTypeCompare((PetscObject)ksp->pc,PCNONE,&nopreconditioner);
43: if (lsqr->vwork_m) {
44: VecDestroyVecs(lsqr->nwork_m,&lsqr->vwork_m);
45: }
47: if (lsqr->vwork_n) {
48: VecDestroyVecs(lsqr->nwork_n,&lsqr->vwork_n);
49: }
51: lsqr->nwork_m = 2;
52: if (nopreconditioner) lsqr->nwork_n = 4;
53: else lsqr->nwork_n = 5;
54: KSPCreateVecs(ksp,lsqr->nwork_n,&lsqr->vwork_n,lsqr->nwork_m,&lsqr->vwork_m);
56: if (lsqr->se_flg && !lsqr->se) {
57: VecDuplicate(lsqr->vwork_n[0],&lsqr->se);
58: VecSet(lsqr->se,PETSC_INFINITY);
59: } else if (!lsqr->se_flg) {
60: VecDestroy(&lsqr->se);
61: }
62: return(0);
63: }
65: static PetscErrorCode KSPSolve_LSQR(KSP ksp)
66: {
68: PetscInt i,size1,size2;
69: PetscScalar rho,rhobar,phi,phibar,theta,c,s,tmp,tau;
70: PetscReal beta,alpha,rnorm;
71: Vec X,B,V,V1,U,U1,TMP,W,W2,Z = NULL;
72: Mat Amat,Pmat;
73: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
74: PetscBool diagonalscale,nopreconditioner;
77: PCGetDiagonalScale(ksp->pc,&diagonalscale);
78: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
80: PCGetOperators(ksp->pc,&Amat,&Pmat);
81: PetscObjectTypeCompare((PetscObject)ksp->pc,PCNONE,&nopreconditioner);
83: /* vectors of length m, where system size is mxn */
84: B = ksp->vec_rhs;
85: U = lsqr->vwork_m[0];
86: U1 = lsqr->vwork_m[1];
88: /* vectors of length n */
89: X = ksp->vec_sol;
90: W = lsqr->vwork_n[0];
91: V = lsqr->vwork_n[1];
92: V1 = lsqr->vwork_n[2];
93: W2 = lsqr->vwork_n[3];
94: if (!nopreconditioner) Z = lsqr->vwork_n[4];
96: /* standard error vector */
97: if (lsqr->se) {
98: VecSet(lsqr->se,0.0);
99: }
101: /* Compute initial residual, temporarily use work vector u */
102: if (!ksp->guess_zero) {
103: KSP_MatMult(ksp,Amat,X,U); /* u <- b - Ax */
104: VecAYPX(U,-1.0,B);
105: } else {
106: VecCopy(B,U); /* u <- b (x is 0) */
107: }
109: /* Test for nothing to do */
110: VecNorm(U,NORM_2,&rnorm);
111: PetscObjectSAWsTakeAccess((PetscObject)ksp);
112: ksp->its = 0;
113: ksp->rnorm = rnorm;
114: PetscObjectSAWsGrantAccess((PetscObject)ksp);
115: KSPLogResidualHistory(ksp,rnorm);
116: KSPMonitor(ksp,0,rnorm);
117: (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);
118: if (ksp->reason) return(0);
120: beta = rnorm;
121: VecScale(U,1.0/beta);
122: KSP_MatMultTranspose(ksp,Amat,U,V);
123: if (nopreconditioner) {
124: VecNorm(V,NORM_2,&alpha);
125: } else {
126: /* this is an application of the preconditioner for the normal equations; not the operator, see the manual page */
127: PCApply(ksp->pc,V,Z);
128: VecDotRealPart(V,Z,&alpha);
129: if (alpha <= 0.0) {
130: ksp->reason = KSP_DIVERGED_BREAKDOWN;
131: return(0);
132: }
133: alpha = PetscSqrtReal(alpha);
134: VecScale(Z,1.0/alpha);
135: }
136: VecScale(V,1.0/alpha);
138: if (nopreconditioner) {
139: VecCopy(V,W);
140: } else {
141: VecCopy(Z,W);
142: }
144: if (lsqr->exact_norm) {
145: MatNorm(Amat,NORM_FROBENIUS,&lsqr->anorm);
146: } else lsqr->anorm = 0.0;
148: lsqr->arnorm = alpha * beta;
149: phibar = beta;
150: rhobar = alpha;
151: i = 0;
152: do {
153: if (nopreconditioner) {
154: KSP_MatMult(ksp,Amat,V,U1);
155: } else {
156: KSP_MatMult(ksp,Amat,Z,U1);
157: }
158: VecAXPY(U1,-alpha,U);
159: VecNorm(U1,NORM_2,&beta);
160: if (beta > 0.0) {
161: VecScale(U1,1.0/beta); /* beta*U1 = Amat*V - alpha*U */
162: if (!lsqr->exact_norm) {
163: lsqr->anorm = PetscSqrtScalar(PetscSqr(lsqr->anorm) + PetscSqr(alpha) + PetscSqr(beta));
164: }
165: }
167: KSP_MatMultTranspose(ksp,Amat,U1,V1);
168: VecAXPY(V1,-beta,V);
169: if (nopreconditioner) {
170: VecNorm(V1,NORM_2,&alpha);
171: } else {
172: PCApply(ksp->pc,V1,Z);
173: VecDotRealPart(V1,Z,&alpha);
174: if (alpha <= 0.0) {
175: ksp->reason = KSP_DIVERGED_BREAKDOWN;
176: break;
177: }
178: alpha = PetscSqrtReal(alpha);
179: VecScale(Z,1.0/alpha);
180: }
181: VecScale(V1,1.0/alpha); /* alpha*V1 = Amat^T*U1 - beta*V */
182: rho = PetscSqrtScalar(rhobar*rhobar + beta*beta);
183: c = rhobar / rho;
184: s = beta / rho;
185: theta = s * alpha;
186: rhobar = -c * alpha;
187: phi = c * phibar;
188: phibar = s * phibar;
189: tau = s * phi;
191: VecAXPY(X,phi/rho,W); /* x <- x + (phi/rho) w */
193: if (lsqr->se) {
194: VecCopy(W,W2);
195: VecSquare(W2);
196: VecScale(W2,1.0/(rho*rho));
197: VecAXPY(lsqr->se, 1.0, W2); /* lsqr->se <- lsqr->se + (w^2/rho^2) */
198: }
199: if (nopreconditioner) {
200: VecAYPX(W,-theta/rho,V1); /* w <- v - (theta/rho) w */
201: } else {
202: VecAYPX(W,-theta/rho,Z); /* w <- z - (theta/rho) w */
203: }
205: lsqr->arnorm = alpha*PetscAbsScalar(tau);
206: rnorm = PetscRealPart(phibar);
208: PetscObjectSAWsTakeAccess((PetscObject)ksp);
209: ksp->its++;
210: ksp->rnorm = rnorm;
211: PetscObjectSAWsGrantAccess((PetscObject)ksp);
212: KSPLogResidualHistory(ksp,rnorm);
213: KSPMonitor(ksp,i+1,rnorm);
214: (*ksp->converged)(ksp,i+1,rnorm,&ksp->reason,ksp->cnvP);
215: if (ksp->reason) break;
216: SWAP(U1,U,TMP);
217: SWAP(V1,V,TMP);
219: i++;
220: } while (i<ksp->max_it);
221: if (i >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
223: /* Finish off the standard error estimates */
224: if (lsqr->se) {
225: tmp = 1.0;
226: MatGetSize(Amat,&size1,&size2);
227: if (size1 > size2) tmp = size1 - size2;
228: tmp = rnorm / PetscSqrtScalar(tmp);
229: VecSqrtAbs(lsqr->se);
230: VecScale(lsqr->se,tmp);
231: }
232: return(0);
233: }
236: PetscErrorCode KSPDestroy_LSQR(KSP ksp)
237: {
238: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
242: /* Free work vectors */
243: if (lsqr->vwork_n) {
244: VecDestroyVecs(lsqr->nwork_n,&lsqr->vwork_n);
245: }
246: if (lsqr->vwork_m) {
247: VecDestroyVecs(lsqr->nwork_m,&lsqr->vwork_m);
248: }
249: VecDestroy(&lsqr->se);
250: PetscFree(ksp->data);
251: return(0);
252: }
254: /*@
255: KSPLSQRSetComputeStandardErrorVec - Compute vector of standard error estimates during KSPSolve_LSQR().
257: Not Collective
259: Input Parameters:
260: + ksp - iterative context
261: - flg - compute the vector of standard estimates or not
263: Developer notes:
264: Vaclav: I'm not sure whether this vector is useful for anything.
266: Level: intermediate
268: .keywords: KSP, KSPLSQR
270: .seealso: KSPSolve(), KSPLSQR, KSPLSQRGetStandardErrorVec()
271: @*/
272: PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP ksp, PetscBool flg)
273: {
274: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
277: lsqr->se_flg = flg;
278: return(0);
279: }
281: /*@
282: KSPLSQRSetExactMatNorm - Compute exact matrix norm instead of iteratively refined estimate.
284: Not Collective
286: Input Parameters:
287: + ksp - iterative context
288: - flg - compute exact matrix norm or not
290: Notes:
291: By default, flg=PETSC_FALSE. This is usually preferred to avoid possibly expensive computation of the norm.
292: For flg=PETSC_TRUE, we call MatNorm(Amat,NORM_FROBENIUS,&lsqr->anorm) which will work only for some types of explicitly assembled matrices.
293: This can affect convergence rate as KSPLSQRConvergedDefault() assumes different value of ||A|| used in normal equation stopping criterion.
295: Level: intermediate
297: .keywords: KSP, KSPLSQR
299: .seealso: KSPSolve(), KSPLSQR, KSPLSQRGetNorms(), KSPLSQRConvergedDefault()
300: @*/
301: PetscErrorCode KSPLSQRSetExactMatNorm(KSP ksp, PetscBool flg)
302: {
303: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
306: lsqr->exact_norm = flg;
307: return(0);
308: }
310: /*@
311: KSPLSQRGetStandardErrorVec - Get vector of standard error estimates.
312: Only available if -ksp_lsqr_set_standard_error was set to true
313: or KSPLSQRSetComputeStandardErrorVec(ksp, PETSC_TRUE) was called.
314: Otherwise returns NULL.
316: Not Collective
318: Input Parameters:
319: . ksp - iterative context
321: Output Parameters:
322: . se - vector of standard estimates
324: Options Database Keys:
325: . -ksp_lsqr_set_standard_error - set standard error estimates of solution
327: Developer notes:
328: Vaclav: I'm not sure whether this vector is useful for anything.
330: Level: intermediate
332: .keywords: KSP, KSPLSQR
334: .seealso: KSPSolve(), KSPLSQR, KSPLSQRSetComputeStandardErrorVec()
335: @*/
336: PetscErrorCode KSPLSQRGetStandardErrorVec(KSP ksp,Vec *se)
337: {
338: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
341: *se = lsqr->se;
342: return(0);
343: }
345: /*@
346: KSPLSQRGetNorms - Get norm estimates that LSQR computes internally during KSPSolve().
348: Not Collective
350: Input Parameters:
351: . ksp - iterative context
353: Output Parameters:
354: + arnorm - good estimate of norm((A*inv(Pmat))'*r), where r = A*x - b, used in specific stopping criterion
355: - anorm - poor estimate of norm(A*inv(Pmat),'fro') used in specific stopping criterion
357: Notes:
358: Output parameters are meaningful only after KSPSolve().
359: These are the same quantities as normar and norma in MATLAB's lsqr(), whose output lsvec is a vector of normar / norma for all iterations.
360: If -ksp_lsqr_exact_mat_norm is set or KSPLSQRSetExactMatNorm(ksp, PETSC_TRUE) called, then anorm is exact Frobenius norm.
362: Level: intermediate
364: .keywords: KSP, KSPLSQR
366: .seealso: KSPSolve(), KSPLSQR, KSPLSQRSetExactMatNorm()
367: @*/
368: PetscErrorCode KSPLSQRGetNorms(KSP ksp,PetscReal *arnorm, PetscReal *anorm)
369: {
370: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
373: if (arnorm) *arnorm = lsqr->arnorm;
374: if (anorm) *anorm = lsqr->anorm;
375: return(0);
376: }
378: /*@C
379: KSPLSQRMonitorDefault - Print the residual norm at each iteration of the LSQR method,
380: norm of the residual of the normal equations A'*A x = A' b, and estimate of matrix norm ||A||.
382: Collective on KSP
384: Input Parameters:
385: + ksp - iterative context
386: . n - iteration number
387: . rnorm - 2-norm (preconditioned) residual value (may be estimated).
388: - dummy - viewer and format context
390: Level: intermediate
392: .keywords: KSP, KSPLSQR, default, monitor, residual
394: .seealso: KSPLSQR, KSPMonitorSet(), KSPMonitorTrueResidualNorm(), KSPMonitorLGResidualNormCreate(), KSPMonitorDefault()
395: @*/
396: PetscErrorCode KSPLSQRMonitorDefault(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat *dummy)
397: {
399: PetscViewer viewer = dummy->viewer;
400: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
403: PetscViewerPushFormat(viewer,dummy->format);
404: PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
405: if (!n && ((PetscObject)ksp)->prefix) {
406: PetscViewerASCIIPrintf(viewer," Residual norm, norm of normal equations, and matrix norm for %s solve.\n",((PetscObject)ksp)->prefix);
407: }
409: if (!n) {
410: PetscViewerASCIIPrintf(viewer,"%3D KSP resid norm %14.12e\n",n,(double)rnorm);
411: } else {
412: PetscViewerASCIIPrintf(viewer,"%3D KSP resid norm %14.12e normal eq resid norm %14.12e matrix norm %14.12e\n",n,(double)rnorm,(double)lsqr->arnorm,(double)lsqr->anorm);
413: }
415: PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
416: PetscViewerPopFormat(viewer);
417: return(0);
418: }
420: PetscErrorCode KSPSetFromOptions_LSQR(PetscOptionItems *PetscOptionsObject,KSP ksp)
421: {
423: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
426: PetscOptionsHead(PetscOptionsObject,"KSP LSQR Options");
427: PetscOptionsBool("-ksp_lsqr_compute_standard_error","Set Standard Error Estimates of Solution","KSPLSQRSetComputeStandardErrorVec",lsqr->se_flg,&lsqr->se_flg,NULL);
428: PetscOptionsBool("-ksp_lsqr_exact_mat_norm","Compute exact matrix norm instead of iteratively refined estimate","KSPLSQRSetExactMatNorm",lsqr->exact_norm,&lsqr->exact_norm,NULL);
429: KSPMonitorSetFromOptions(ksp,"-ksp_lsqr_monitor","Monitor residual norm and norm of residual of normal equations","KSPMonitorSet",KSPLSQRMonitorDefault);
430: PetscOptionsTail();
431: return(0);
432: }
434: PetscErrorCode KSPView_LSQR(KSP ksp,PetscViewer viewer)
435: {
436: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
438: PetscBool iascii;
441: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
442: if (iascii) {
443: if (lsqr->se) {
444: PetscReal rnorm;
445: VecNorm(lsqr->se,NORM_2,&rnorm);
446: PetscViewerASCIIPrintf(viewer," norm of standard error %g, iterations %d\n",(double)rnorm,ksp->its);
447: } else {
448: PetscViewerASCIIPrintf(viewer," standard error not computed\n");
449: }
450: if (lsqr->exact_norm) {
451: PetscViewerASCIIPrintf(viewer," using exact matrix norm\n");
452: } else {
453: PetscViewerASCIIPrintf(viewer," using inexact matrix norm\n");
454: }
455: }
456: return(0);
457: }
459: /*@C
460: KSPLSQRConvergedDefault - Determines convergence of the LSQR Krylov method.
462: Collective on KSP
464: Input Parameters:
465: + ksp - iterative context
466: . n - iteration number
467: . rnorm - 2-norm residual value (may be estimated)
468: - ctx - convergence context which must be created by KSPConvergedDefaultCreate()
470: reason is set to:
471: + positive - if the iteration has converged;
472: . negative - if residual norm exceeds divergence threshold;
473: - 0 - otherwise.
475: Notes:
476: KSPConvergedDefault() is called first to check for convergence in A*x=b.
477: If that does not determine convergence then checks convergence for the least squares problem, i.e. in min{|b-A*x|}.
478: Possible convergence for the least squares problem (which is based on the residual of the normal equations) are KSP_CONVERGED_RTOL_NORMAL norm and KSP_CONVERGED_ATOL_NORMAL.
479: KSP_CONVERGED_RTOL_NORMAL is returned if ||A'*r|| < rtol * ||A|| * ||r||.
480: Matrix norm ||A|| is iteratively refined estimate, see KSPLSQRGetNorms().
481: This criterion is now largely compatible with that in MATLAB lsqr().
483: Level: intermediate
485: .keywords: KSP, KSPLSQR, default, convergence, residual
487: .seealso: KSPLSQR, KSPSetConvergenceTest(), KSPSetTolerances(), KSPConvergedSkip(), KSPConvergedReason, KSPGetConvergedReason(),
488: KSPConvergedDefaultSetUIRNorm(), KSPConvergedDefaultSetUMIRNorm(), KSPConvergedDefaultCreate(), KSPConvergedDefaultDestroy(), KSPConvergedDefault(), KSPLSQRGetNorms(), KSPLSQRSetExactMatNorm()
489: @*/
490: PetscErrorCode KSPLSQRConvergedDefault(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
491: {
493: KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data;
496: /* check for convergence in A*x=b */
497: KSPConvergedDefault(ksp,n,rnorm,reason,ctx);
498: if (!n || *reason) return(0);
500: /* check for convergence in min{|b-A*x|} */
501: if (lsqr->arnorm < ksp->abstol) {
502: PetscInfo3(ksp,"LSQR solver has converged. Normal equation residual %14.12e is less than absolute tolerance %14.12e at iteration %D\n",(double)lsqr->arnorm,(double)ksp->abstol,n);
503: *reason = KSP_CONVERGED_ATOL_NORMAL;
504: } else if (lsqr->arnorm < ksp->rtol * lsqr->anorm * rnorm) {
505: PetscInfo6(ksp,"LSQR solver has converged. Normal equation residual %14.12e is less than rel. tol. %14.12e times %s Frobenius norm of matrix %14.12e times residual %14.12e at iteration %D\n",(double)lsqr->arnorm,(double)ksp->rtol,lsqr->exact_norm?"exact":"approx.",(double)lsqr->anorm,(double)rnorm,n);
506: *reason = KSP_CONVERGED_RTOL_NORMAL;
507: }
508: return(0);
509: }
511: /*MC
512: KSPLSQR - This implements LSQR
514: Options Database Keys:
515: + -ksp_lsqr_set_standard_error - set standard error estimates of solution, see KSPLSQRSetComputeStandardErrorVec() and KSPLSQRGetStandardErrorVec()
516: . -ksp_lsqr_exact_mat_norm - compute exact matrix norm instead of iteratively refined estimate, see KSPLSQRSetExactMatNorm()
517: - -ksp_lsqr_monitor - monitor residual norm, norm of residual of normal equations A'*A x = A' b, and estimate of matrix norm ||A||
519: Level: beginner
521: Notes:
522: Supports non-square (rectangular) matrices.
524: This varient, when applied with no preconditioning is identical to the original algorithm in exact arithematic; however, in practice, with no preconditioning
525: due to inexact arithematic, it can converge differently. Hence when no preconditioner is used (PCType PCNONE) it automatically reverts to the original algorithm.
527: With the PETSc built-in preconditioners, such as ICC, one should call KSPSetOperators(ksp,A,A'*A)) since the preconditioner needs to work
528: for the normal equations A'*A.
530: Supports only left preconditioning.
532: For least squares problems wit nonzero residual A*x - b, there are additional convergence tests for the residual of the normal equations, A'*(b - Ax), see KSPLSQRConvergedDefault().
534: References:
535: . 1. - The original unpreconditioned algorithm can be found in Paige and Saunders, ACM Transactions on Mathematical Software, Vol 8, 1982.
537: In exact arithmetic the LSQR method (with no preconditioning) is identical to the KSPCG algorithm applied to the normal equations.
538: The preconditioned variant was implemented by Bas van't Hof and is essentially a left preconditioning for the Normal Equations. It appears the implementation with preconditioner
539: track the true norm of the residual and uses that in the convergence test.
541: Developer Notes:
542: How is this related to the KSPCGNE implementation? One difference is that KSPCGNE applies
543: the preconditioner transpose times the preconditioner, so one does not need to pass A'*A as the third argument to KSPSetOperators().
547: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPSolve(), KSPLSQRConvergedDefault(), KSPLSQRSetComputeStandardErrorVec(), KSPLSQRGetStandardErrorVec(), KSPLSQRSetExactMatNorm()
549: M*/
550: PETSC_EXTERN PetscErrorCode KSPCreate_LSQR(KSP ksp)
551: {
552: KSP_LSQR *lsqr;
553: void *ctx;
557: PetscNewLog(ksp,&lsqr);
558: lsqr->se = NULL;
559: lsqr->se_flg = PETSC_FALSE;
560: lsqr->exact_norm = PETSC_FALSE;
561: lsqr->anorm = -1.0;
562: lsqr->arnorm = -1.0;
563: ksp->data = (void*)lsqr;
564: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,3);
566: ksp->ops->setup = KSPSetUp_LSQR;
567: ksp->ops->solve = KSPSolve_LSQR;
568: ksp->ops->destroy = KSPDestroy_LSQR;
569: ksp->ops->setfromoptions = KSPSetFromOptions_LSQR;
570: ksp->ops->view = KSPView_LSQR;
572: KSPConvergedDefaultCreate(&ctx);
573: KSPSetConvergenceTest(ksp,KSPLSQRConvergedDefault,ctx,KSPConvergedDefaultDestroy);
574: return(0);
575: }