Actual source code: lmedense.c
slepc-3.16.1 2021-11-17
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: Routines for solving dense matrix equations, in some cases calling SLICOT
12: */
14: #include <slepc/private/lmeimpl.h>
15: #include <slepcblaslapack.h>
17: /*
18: LMEDenseRankSVD - given a square matrix A, compute its SVD U*S*V', and determine the
19: numerical rank. On exit, U contains U*S and A is overwritten with V'
20: */
21: PetscErrorCode LMEDenseRankSVD(LME lme,PetscInt n,PetscScalar *A,PetscInt lda,PetscScalar *U,PetscInt ldu,PetscInt *rank)
22: {
24: PetscInt i,j,rk=0;
25: PetscScalar *work;
26: PetscReal tol,*sg,*rwork;
27: PetscBLASInt n_,lda_,ldu_,info,lw_;
30: PetscCalloc3(n,&sg,10*n,&work,5*n,&rwork);
31: PetscBLASIntCast(n,&n_);
32: PetscBLASIntCast(lda,&lda_);
33: PetscBLASIntCast(ldu,&ldu_);
34: lw_ = 10*n_;
35: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
36: #if !defined (PETSC_USE_COMPLEX)
37: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","O",&n_,&n_,A,&lda_,sg,U,&ldu_,NULL,&n_,work,&lw_,&info));
38: #else
39: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","O",&n_,&n_,A,&lda_,sg,U,&ldu_,NULL,&n_,work,&lw_,rwork,&info));
40: #endif
41: PetscFPTrapPop();
42: SlepcCheckLapackInfo("gesvd",info);
43: tol = 10*PETSC_MACHINE_EPSILON*n*sg[0];
44: for (j=0;j<n;j++) {
45: if (sg[j]>tol) {
46: for (i=0;i<n;i++) U[i+j*n] *= sg[j];
47: rk++;
48: } else break;
49: }
50: *rank = rk;
51: PetscFree3(sg,work,rwork);
52: return(0);
53: }
55: #if defined(PETSC_USE_INFO)
56: /*
57: LyapunovCholResidual - compute the residual norm ||A*U'*U+U'*U*A'+B*B'||
58: */
59: static PetscErrorCode LyapunovCholResidual(PetscInt m,PetscScalar *A,PetscInt lda,PetscInt k,PetscScalar *B,PetscInt ldb,PetscScalar *U,PetscInt ldu,PetscReal *res)
60: {
62: PetscBLASInt n,kk,la,lb,lu;
63: PetscScalar *M,*R,zero=0.0,done=1.0;
66: *res = 0;
67: PetscBLASIntCast(lda,&la);
68: PetscBLASIntCast(ldb,&lb);
69: PetscBLASIntCast(ldu,&lu);
70: PetscBLASIntCast(m,&n);
71: PetscBLASIntCast(k,&kk);
72: PetscMalloc2(m*m,&M,m*m,&R);
74: /* R = B*B' */
75: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&n,&n,&kk,&done,B,&lb,B,&lb,&zero,R,&n));
76: /* M = A*U' */
77: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&n,&n,&n,&done,A,&la,U,&lu,&zero,M,&n));
78: /* R = R+M*U */
79: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&done,M,&n,U,&lu,&done,R,&n));
80: /* R = R+U'*M' */
81: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","C",&n,&n,&n,&done,U,&lu,M,&n,&done,R,&n));
83: *res = LAPACKlange_("F",&n,&n,R,&n,NULL);
84: PetscFree2(M,R);
85: return(0);
86: }
88: /*
89: LyapunovResidual - compute the residual norm ||A*X+X*A'+B||
90: */
91: static PetscErrorCode LyapunovResidual(PetscInt m,PetscScalar *A,PetscInt lda,PetscScalar *B,PetscInt ldb,PetscScalar *X,PetscInt ldx,PetscReal *res)
92: {
94: PetscInt i;
95: PetscBLASInt n,la,lb,lx;
96: PetscScalar *R,done=1.0;
99: *res = 0;
100: PetscBLASIntCast(lda,&la);
101: PetscBLASIntCast(ldb,&lb);
102: PetscBLASIntCast(ldx,&lx);
103: PetscBLASIntCast(m,&n);
104: PetscMalloc1(m*m,&R);
106: /* R = B+A*X */
107: for (i=0;i<m;i++) {
108: PetscArraycpy(R+i*m,B+i*ldb,m);
109: }
110: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&done,A,&la,X,&lx,&done,R,&n));
111: /* R = R+X*A' */
112: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&n,&n,&n,&done,X,&lx,A,&la,&done,R,&n));
114: *res = LAPACKlange_("F",&n,&n,R,&n,NULL);
115: PetscFree(R);
116: return(0);
117: }
118: #endif
120: #if defined(SLEPC_HAVE_SLICOT)
121: /*
122: HessLyapunovChol_SLICOT - implementation used when SLICOT is available
123: */
124: static PetscErrorCode HessLyapunovChol_SLICOT(PetscInt m,PetscScalar *H,PetscInt ldh,PetscInt k,PetscScalar *B,PetscInt ldb,PetscScalar *U,PetscInt ldu,PetscReal *res)
125: {
127: PetscBLASInt lwork,info,n,kk,lu,ione=1,sdim;
128: PetscInt i,j;
129: PetscReal scal;
130: PetscScalar *Q,*W,*wr,*wi,*work;
133: PetscBLASIntCast(ldu,&lu);
134: PetscBLASIntCast(m,&n);
135: PetscBLASIntCast(k,&kk);
136: PetscBLASIntCast(6*m,&lwork);
137: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
139: /* transpose W = H' */
140: PetscMalloc5(m*m,&W,m*m,&Q,m,&wr,m,&wi,lwork,&work);
141: for (j=0;j<m;j++) {
142: for (i=0;i<m;i++) W[i+j*m] = H[j+i*ldh];
143: }
145: /* compute the real Schur form of W */
146: PetscStackCallBLAS("LAPACKgees",LAPACKgees_("V","N",NULL,&n,W,&n,&sdim,wr,wi,Q,&n,work,&lwork,NULL,&info));
147: SlepcCheckLapackInfo("gees",info);
148: #if defined(PETSC_USE_DEBUG)
149: for (i=0;i<m;i++) if (PetscRealPart(wr[i])>=0.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER_INPUT,"Eigenvalue with non-negative real part, the coefficient matrix is not stable");
150: #endif
152: /* copy B' into first rows of U */
153: for (i=0;i<k;i++) {
154: for (j=0;j<m;j++) U[i+j*ldu] = B[j+i*ldb];
155: }
157: /* solve Lyapunov equation (Hammarling) */
158: PetscStackCallBLAS("SLICOTsb03od",SLICOTsb03od_("C","F","N",&n,&kk,W,&n,Q,&n,U,&lu,&scal,wr,wi,work,&lwork,&info));
159: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SLICOT subroutine SB03OD: info=%d",(int)info);
160: if (scal!=1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Current implementation cannot handle scale factor %g",scal);
162: /* resnorm = norm(H(m+1,:)*U'*U), use Q(:,1) = U'*U(:,m) */
163: if (res) {
164: for (j=0;j<m;j++) Q[j] = U[j+(m-1)*ldu];
165: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,U,&lu,Q,&ione));
166: if (k!=1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Residual error is intended for k=1 only, but you set k=%d",(int)k);
167: *res *= BLASnrm2_(&n,Q,&ione);
168: }
170: PetscFPTrapPop();
171: PetscFree5(W,Q,wr,wi,work);
172: return(0);
173: }
175: #else
177: /*
178: Compute the upper Cholesky factor of A
179: */
180: static PetscErrorCode CholeskyFactor(PetscInt m,PetscScalar *A,PetscInt lda)
181: {
183: PetscInt i;
184: PetscScalar *S;
185: PetscBLASInt info,n,ld;
188: PetscBLASIntCast(m,&n);
189: PetscBLASIntCast(lda,&ld);
190: PetscMalloc1(m*m,&S);
191: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
193: /* save a copy of matrix in S */
194: for (i=0;i<m;i++) {
195: PetscArraycpy(S+i*m,A+i*lda,m);
196: }
198: /* compute upper Cholesky factor in R */
199: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("U",&n,A,&ld,&info));
200: PetscLogFlops((1.0*n*n*n)/3.0);
202: if (info) {
203: for (i=0;i<m;i++) {
204: PetscArraycpy(A+i*lda,S+i*m,m);
205: A[i+i*lda] += 50.0*PETSC_MACHINE_EPSILON;
206: }
207: PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("U",&n,A,&ld,&info));
208: SlepcCheckLapackInfo("potrf",info);
209: PetscLogFlops((1.0*n*n*n)/3.0);
210: }
212: /* Zero out entries below the diagonal */
213: for (i=0;i<m-1;i++) {
214: PetscArrayzero(A+i*lda+i+1,m-i-1);
215: }
216: PetscFPTrapPop();
217: PetscFree(S);
218: return(0);
219: }
221: /*
222: HessLyapunovChol_LAPACK - alternative implementation when SLICOT is not available
223: */
224: static PetscErrorCode HessLyapunovChol_LAPACK(PetscInt m,PetscScalar *H,PetscInt ldh,PetscInt k,PetscScalar *B,PetscInt ldb,PetscScalar *U,PetscInt ldu,PetscReal *res)
225: {
227: PetscBLASInt ilo=1,lwork,info,n,kk,lu,lb,ione=1;
228: PetscInt i,j;
229: PetscReal scal;
230: PetscScalar *Q,*C,*W,*Z,*wr,*work,zero=0.0,done=1.0,dmone=-1.0;
231: #if !defined(PETSC_USE_COMPLEX)
232: PetscScalar *wi;
233: #endif
236: PetscBLASIntCast(ldb,&lb);
237: PetscBLASIntCast(ldu,&lu);
238: PetscBLASIntCast(m,&n);
239: PetscBLASIntCast(k,&kk);
240: PetscBLASIntCast(6*m,&lwork);
241: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
242: C = U;
244: #if !defined(PETSC_USE_COMPLEX)
245: PetscMalloc6(m*m,&Q,m*m,&W,m*k,&Z,m,&wr,m,&wi,lwork,&work);
246: #else
247: PetscMalloc5(m*m,&Q,m*m,&W,m*k,&Z,m,&wr,lwork,&work);
248: #endif
250: /* save a copy W = H */
251: for (j=0;j<m;j++) {
252: for (i=0;i<m;i++) W[i+j*m] = H[i+j*ldh];
253: }
255: /* compute the (real) Schur form of W */
256: #if !defined(PETSC_USE_COMPLEX)
257: PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("S","I",&n,&ilo,&n,W,&n,wr,wi,Q,&n,work,&lwork,&info));
258: #else
259: PetscStackCallBLAS("LAPACKhseqr",LAPACKhseqr_("S","I",&n,&ilo,&n,W,&n,wr,Q,&n,work,&lwork,&info));
260: #endif
261: SlepcCheckLapackInfo("hseqr",info);
262: #if defined(PETSC_USE_DEBUG)
263: for (i=0;i<m;i++) if (PetscRealPart(wr[i])>=0.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER_INPUT,"Eigenvalue with non-negative real part %g, the coefficient matrix is not stable",PetscRealPart(wr[i]));
264: #endif
266: /* C = -Z*Z', Z = Q'*B */
267: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n,&kk,&n,&done,Q,&n,B,&lb,&zero,Z,&n));
268: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&n,&n,&kk,&dmone,Z,&n,Z,&n,&zero,C,&lu));
270: /* solve triangular Sylvester equation */
271: PetscStackCallBLAS("LAPACKtrsyl",LAPACKtrsyl_("N","C",&ione,&n,&n,W,&n,W,&n,C,&lu,&scal,&info));
272: SlepcCheckLapackInfo("trsyl",info);
273: if (scal!=1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Current implementation cannot handle scale factor %g",scal);
275: /* back-transform C = Q*C*Q' */
276: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&done,Q,&n,C,&n,&zero,W,&n));
277: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&n,&n,&n,&done,W,&n,Q,&n,&zero,C,&lu));
279: /* resnorm = norm(H(m+1,:)*Y) */
280: if (res) {
281: if (k!=1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Residual error is intended for k=1 only, but you set k=%d",(int)k);
282: *res *= BLASnrm2_(&n,C+m-1,&n);
283: }
285: /* U = chol(C) */
286: CholeskyFactor(m,C,ldu);
288: PetscFPTrapPop();
289: #if !defined(PETSC_USE_COMPLEX)
290: PetscFree6(Q,W,Z,wr,wi,work);
291: #else
292: PetscFree5(Q,W,Z,wr,work);
293: #endif
294: return(0);
295: }
297: #endif /* SLEPC_HAVE_SLICOT */
299: /*@C
300: LMEDenseHessLyapunovChol - Computes the Cholesky factor of the solution of a
301: dense Lyapunov equation with an upper Hessenberg coefficient matrix.
303: Logically Collective on lme
305: Input Parameters:
306: + lme - linear matrix equation solver context
307: . m - number of rows and columns of H
308: . H - coefficient matrix
309: . ldh - leading dimension of H
310: . k - number of columns of B
311: . B - right-hand side matrix
312: . ldb - leading dimension of B
313: - ldu - leading dimension of U
315: Output Parameter:
316: . U - Cholesky factor of the solution
318: Input/Output Parameter:
319: . res - (optional) residual norm, on input it should contain H(m+1,m)
321: Note:
322: The Lyapunov equation has the form H*X + X*H' = -B*B', where H is an mxm
323: upper Hessenberg matrix, B is an mxk matrix and the solution is expressed
324: as X = U'*U, where U is upper triangular. H is supposed to be stable.
326: When k=1 and the res argument is provided, the last row of X is used to
327: compute the residual norm of a Lyapunov equation projected via Arnoldi.
329: Level: developer
331: .seealso: LMEDenseLyapunov(), LMESolve()
332: @*/
333: PetscErrorCode LMEDenseHessLyapunovChol(LME lme,PetscInt m,PetscScalar *H,PetscInt ldh,PetscInt k,PetscScalar *B,PetscInt ldb,PetscScalar *U,PetscInt ldu,PetscReal *res)
334: {
336: #if defined(PETSC_USE_INFO)
337: PetscReal error;
338: #endif
352: #if defined(SLEPC_HAVE_SLICOT)
353: HessLyapunovChol_SLICOT(m,H,ldh,k,B,ldb,U,ldu,res);
354: #else
355: HessLyapunovChol_LAPACK(m,H,ldh,k,B,ldb,U,ldu,res);
356: #endif
358: #if defined(PETSC_USE_INFO)
359: if (PetscLogPrintInfo) {
360: LyapunovCholResidual(m,H,ldh,k,B,ldb,U,ldu,&error);
361: PetscInfo1(lme,"Residual norm of dense Lyapunov equation = %g\n",error);
362: }
363: #endif
364: return(0);
365: }
367: #if defined(SLEPC_HAVE_SLICOT)
368: /*
369: Lyapunov_SLICOT - implementation used when SLICOT is available
370: */
371: static PetscErrorCode Lyapunov_SLICOT(PetscInt m,PetscScalar *H,PetscInt ldh,PetscScalar *B,PetscInt ldb,PetscScalar *X,PetscInt ldx)
372: {
374: PetscBLASInt sdim,lwork,info,n,lx,*iwork;
375: PetscInt i,j;
376: PetscReal scal,sep,ferr,*work;
377: PetscScalar *Q,*W,*wr,*wi;
380: PetscBLASIntCast(ldx,&lx);
381: PetscBLASIntCast(m,&n);
382: PetscBLASIntCast(PetscMax(20,m*m),&lwork);
383: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
385: /* transpose W = H' */
386: PetscMalloc6(m*m,&W,m*m,&Q,m,&wr,m,&wi,m*m,&iwork,lwork,&work);
387: for (j=0;j<m;j++) {
388: for (i=0;i<m;i++) W[i+j*m] = H[j+i*ldh];
389: }
391: /* compute the real Schur form of W */
392: PetscStackCallBLAS("LAPACKgees",LAPACKgees_("V","N",NULL,&n,W,&n,&sdim,wr,wi,Q,&n,work,&lwork,NULL,&info));
393: SlepcCheckLapackInfo("gees",info);
395: /* copy -B into X */
396: for (i=0;i<m;i++) {
397: for (j=0;j<m;j++) X[i+j*ldx] = -B[i+j*ldb];
398: }
400: /* solve Lyapunov equation (Hammarling) */
401: PetscStackCallBLAS("SLICOTsb03md",SLICOTsb03md_("C","X","F","N",&n,W,&n,Q,&n,X,&lx,&scal,&sep,&ferr,wr,wi,iwork,work,&lwork,&info));
402: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SLICOT subroutine SB03OD: info=%d",(int)info);
403: if (scal!=1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Current implementation cannot handle scale factor %g",scal);
405: PetscFPTrapPop();
406: PetscFree6(W,Q,wr,wi,iwork,work);
407: return(0);
408: }
410: #else
412: /*
413: Lyapunov_LAPACK - alternative implementation when SLICOT is not available
414: */
415: static PetscErrorCode Lyapunov_LAPACK(PetscInt m,PetscScalar *A,PetscInt lda,PetscScalar *B,PetscInt ldb,PetscScalar *X,PetscInt ldx)
416: {
418: PetscBLASInt sdim,lwork,info,n,lx,lb,ione=1;
419: PetscInt i,j;
420: PetscReal scal;
421: PetscScalar *Q,*W,*Z,*wr,*work,zero=0.0,done=1.0,dmone=-1.0;
422: #if defined(PETSC_USE_COMPLEX)
423: PetscReal *rwork;
424: #else
425: PetscScalar *wi;
426: #endif
429: PetscBLASIntCast(ldb,&lb);
430: PetscBLASIntCast(ldx,&lx);
431: PetscBLASIntCast(m,&n);
432: PetscBLASIntCast(6*m,&lwork);
433: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
435: #if !defined(PETSC_USE_COMPLEX)
436: PetscMalloc6(m*m,&Q,m*m,&W,m*m,&Z,m,&wr,m,&wi,lwork,&work);
437: #else
438: PetscMalloc6(m*m,&Q,m*m,&W,m*m,&Z,m,&wr,lwork,&work,m,&rwork);
439: #endif
441: /* save a copy W = A */
442: for (j=0;j<m;j++) {
443: for (i=0;i<m;i++) W[i+j*m] = A[i+j*lda];
444: }
446: /* compute the (real) Schur form of W */
447: #if !defined(PETSC_USE_COMPLEX)
448: PetscStackCallBLAS("LAPACKgees",LAPACKgees_("V","N",NULL,&n,W,&n,&sdim,wr,wi,Q,&n,work,&lwork,NULL,&info));
449: #else
450: PetscStackCallBLAS("LAPACKgees",LAPACKgees_("V","N",NULL,&n,W,&n,&sdim,wr,Q,&n,work,&lwork,rwork,NULL,&info));
451: #endif
452: SlepcCheckLapackInfo("gees",info);
454: /* X = -Q'*B*Q */
455: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n,&n,&n,&done,Q,&n,B,&lb,&zero,Z,&n));
456: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&dmone,Z,&n,Q,&n,&zero,X,&lx));
458: /* solve triangular Sylvester equation */
459: PetscStackCallBLAS("LAPACKtrsyl",LAPACKtrsyl_("N","C",&ione,&n,&n,W,&n,W,&n,X,&lx,&scal,&info));
460: SlepcCheckLapackInfo("trsyl",info);
461: if (scal!=1.0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Current implementation cannot handle scale factor %g",scal);
463: /* back-transform X = Q*X*Q' */
464: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&done,Q,&n,X,&n,&zero,W,&n));
465: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&n,&n,&n,&done,W,&n,Q,&n,&zero,X,&lx));
467: PetscFPTrapPop();
468: #if !defined(PETSC_USE_COMPLEX)
469: PetscFree6(Q,W,Z,wr,wi,work);
470: #else
471: PetscFree6(Q,W,Z,wr,work,rwork);
472: #endif
473: return(0);
474: }
476: #endif /* SLEPC_HAVE_SLICOT */
478: /*@C
479: LMEDenseLyapunov - Computes the solution of a dense continuous-time Lyapunov
480: equation.
482: Logically Collective on lme
484: Input Parameters:
485: + lme - linear matrix equation solver context
486: . m - number of rows and columns of A
487: . A - coefficient matrix
488: . lda - leading dimension of A
489: . B - right-hand side matrix
490: . ldb - leading dimension of B
491: - ldx - leading dimension of X
493: Output Parameter:
494: . X - the solution
496: Note:
497: The Lyapunov equation has the form A*X + X*A' = -B, where all are mxm
498: matrices, and B is symmetric.
500: Level: developer
502: .seealso: LMEDenseHessLyapunovChol(), LMESolve()
503: @*/
504: PetscErrorCode LMEDenseLyapunov(LME lme,PetscInt m,PetscScalar *A,PetscInt lda,PetscScalar *B,PetscInt ldb,PetscScalar *X,PetscInt ldx)
505: {
507: #if defined(PETSC_USE_INFO)
508: PetscReal error;
509: #endif
521: #if defined(SLEPC_HAVE_SLICOT)
522: Lyapunov_SLICOT(m,A,lda,B,ldb,X,ldx);
523: #else
524: Lyapunov_LAPACK(m,A,lda,B,ldb,X,ldx);
525: #endif
527: #if defined(PETSC_USE_INFO)
528: if (PetscLogPrintInfo) {
529: LyapunovResidual(m,A,lda,B,ldb,X,ldx,&error);
530: PetscInfo1(lme,"Residual norm of dense Lyapunov equation = %g\n",error);
531: }
532: #endif
533: return(0);
534: }