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: BV developer functions needed in contour integral methods
12: */
14: #include <slepc/private/bvimpl.h> 15: #include <slepcblaslapack.h> 17: #define p_id(i) (i*subcomm->n + subcomm->color) 19: /*@
20: BVScatter - Scatters the columns of a BV to another BV created in a
21: subcommunicator.
23: Collective on Vin
25: Input Parameters:
26: + Vin - input basis vectors (defined on the whole communicator)
27: . scat - VecScatter object that contains the info for the communication
28: - xdup - an auxiliary vector
30: Output Parameter:
31: . Vout - output basis vectors (defined on the subcommunicator)
33: Notes:
34: Currently implemented as a loop for each the active column, where each
35: column is scattered independently. The vector xdup is defined on the
36: contiguous parent communicator and have enough space to store one
37: duplicate of the original vector per each subcommunicator.
39: Level: developer
41: .seealso: BVGetColumn()
42: @*/
43: PetscErrorCode BVScatter(BV Vin,BV Vout,VecScatter scat,Vec xdup) 44: {
45: PetscErrorCode ierr;
46: PetscInt i;
47: Vec v;
48: const PetscScalar *array;
55: for (i=Vin->l;i<Vin->k;i++) {
56: BVGetColumn(Vin,i,&v);
57: VecScatterBegin(scat,v,xdup,INSERT_VALUES,SCATTER_FORWARD);
58: VecScatterEnd(scat,v,xdup,INSERT_VALUES,SCATTER_FORWARD);
59: BVRestoreColumn(Vin,i,&v);
60: VecGetArrayRead(xdup,&array);
61: VecPlaceArray(Vout->t,array);
62: BVInsertVec(Vout,i,Vout->t);
63: VecResetArray(Vout->t);
64: VecRestoreArrayRead(xdup,&array);
65: }
66: return(0);
67: }
69: /*@
70: BVSumQuadrature - Computes the sum of terms required in the quadrature
71: rule to approximate the contour integral.
73: Collective on S
75: Input Parameters:
76: + Y - input basis vectors
77: . M - number of moments
78: . L - block size
79: . L_max - maximum block size
80: . w - quadrature weights
81: . zn - normalized quadrature points
82: . scat - (optional) VecScatter object to communicate between subcommunicators
83: . subcomm - subcommunicator layout
84: . npoints - number of points to process by the subcommunicator
85: - useconj - whether conjugate points can be used or not
87: Output Parameter:
88: . S - output basis vectors
90: Notes:
91: This is a generalization of BVMult(). The resulting matrix S consists of M
92: panels of L columns, and the following formula is computed for each panel:
93: S_k = sum_j w_j*zn_j^k*Y_j, where Y_j is the j-th panel of Y containing
94: the result of solving T(z_j)^{-1}*X for each integration point j. L_max is
95: the width of the panels in Y.
97: When using subcommunicators, Y is stored in the subcommunicators for a subset
98: of integration points. In that case, the computation is done in the subcomm
99: and then scattered to the whole communicator in S using the VecScatter scat.
100: The value npoints is the number of points to be processed in this subcomm
101: and the flag useconj indicates whether symmetric points can be reused.
103: Level: developer
105: .seealso: BVMult(), BVScatter(), BVDotQuadrature(), RGComputeQuadrature(), RGCanUseConjugates()
106: @*/
107: PetscErrorCode BVSumQuadrature(BV S,BV Y,PetscInt M,PetscInt L,PetscInt L_max,PetscScalar *w,PetscScalar *zn,VecScatter scat,PetscSubcomm subcomm,PetscInt npoints,PetscBool useconj)108: {
110: PetscInt i,j,k,nloc;
111: Vec v,sj;
112: PetscScalar *ppk,*pv,one=1.0;
119: BVGetSizes(Y,&nloc,NULL,NULL);
120: PetscMalloc1(npoints,&ppk);
121: for (i=0;i<npoints;i++) ppk[i] = 1.0;
122: BVCreateVec(Y,&v);
123: for (k=0;k<M;k++) {
124: for (j=0;j<L;j++) {
125: VecSet(v,0.0);
126: for (i=0;i<npoints;i++) {
127: BVSetActiveColumns(Y,i*L_max+j,i*L_max+j+1);
128: BVMultVec(Y,ppk[i]*w[p_id(i)],1.0,v,&one);
129: }
130: if (useconj) {
131: VecGetArray(v,&pv);
132: for (i=0;i<nloc;i++) pv[i] = 2.0*PetscRealPart(pv[i]);
133: VecRestoreArray(v,&pv);
134: }
135: BVGetColumn(S,k*L+j,&sj);
136: if (scat) {
137: VecScatterBegin(scat,v,sj,ADD_VALUES,SCATTER_REVERSE);
138: VecScatterEnd(scat,v,sj,ADD_VALUES,SCATTER_REVERSE);
139: } else {
140: VecCopy(v,sj);
141: }
142: BVRestoreColumn(S,k*L+j,&sj);
143: }
144: for (i=0;i<npoints;i++) ppk[i] *= zn[p_id(i)];
145: }
146: PetscFree(ppk);
147: VecDestroy(&v);
148: return(0);
149: }
151: /*@
152: BVDotQuadrature - Computes the projection terms required in the quadrature
153: rule to approximate the contour integral.
155: Collective on Y
157: Input Parameters:
158: + Y - first basis vectors
159: . V - second basis vectors
160: . M - number of moments
161: . L - block size
162: . L_max - maximum block size
163: . w - quadrature weights
164: . zn - normalized quadrature points
165: . subcomm - subcommunicator layout
166: . npoints - number of points to process by the subcommunicator
167: - useconj - whether conjugate points can be used or not
169: Output Parameter:
170: . Mu - computed result
172: Notes:
173: This is a generalization of BVDot(). The resulting matrix Mu consists of M
174: blocks of size LxL (placed horizontally), each of them computed as:
175: Mu_k = sum_j w_j*zn_j^k*V'*Y_j, where Y_j is the j-th panel of Y containing
176: the result of solving T(z_j)^{-1}*X for each integration point j. L_max is
177: the width of the panels in Y.
179: When using subcommunicators, Y is stored in the subcommunicators for a subset
180: of integration points. In that case, the computation is done in the subcomm
181: and then the final result is combined via reduction.
182: The value npoints is the number of points to be processed in this subcomm
183: and the flag useconj indicates whether symmetric points can be reused.
185: Level: developer
187: .seealso: BVDot(), BVScatter(), BVSumQuadrature(), RGComputeQuadrature(), RGCanUseConjugates()
188: @*/
189: PetscErrorCode BVDotQuadrature(BV Y,BV V,PetscScalar *Mu,PetscInt M,PetscInt L,PetscInt L_max,PetscScalar *w,PetscScalar *zn,PetscSubcomm subcomm,PetscInt npoints,PetscBool useconj)190: {
192: PetscMPIInt sub_size,count;
193: PetscInt i,j,k,s;
194: PetscScalar *temp,*temp2,*ppk,alp;
195: Mat H;
196: const PetscScalar *pH;
202: MPI_Comm_size(PetscSubcommChild(subcomm),&sub_size);
203: PetscMalloc3(npoints*L*(L+1),&temp,2*M*L*L,&temp2,npoints,&ppk);
204: MatCreateSeqDense(PETSC_COMM_SELF,L,L_max*npoints,NULL,&H);
205: PetscArrayzero(temp2,2*M*L*L);
206: BVSetActiveColumns(Y,0,L_max*npoints);
207: BVSetActiveColumns(V,0,L);
208: BVDot(Y,V,H);
209: MatDenseGetArrayRead(H,&pH);
210: for (i=0;i<npoints;i++) {
211: for (j=0;j<L;j++) {
212: for (k=0;k<L;k++) {
213: temp[k+j*L+i*L*L] = pH[k+j*L+i*L*L_max];
214: }
215: }
216: }
217: MatDenseRestoreArrayRead(H,&pH);
218: for (i=0;i<npoints;i++) ppk[i] = 1;
219: for (k=0;k<2*M;k++) {
220: for (j=0;j<L;j++) {
221: for (i=0;i<npoints;i++) {
222: alp = ppk[i]*w[p_id(i)];
223: for (s=0;s<L;s++) {
224: if (useconj) temp2[s+(j+k*L)*L] += 2.0*PetscRealPart(alp*temp[s+(j+i*L)*L]);
225: else temp2[s+(j+k*L)*L] += alp*temp[s+(j+i*L)*L];
226: }
227: }
228: }
229: for (i=0;i<npoints;i++) ppk[i] *= zn[p_id(i)];
230: }
231: for (i=0;i<2*M*L*L;i++) temp2[i] /= sub_size;
232: PetscMPIIntCast(2*M*L*L,&count);
233: MPIU_Allreduce(temp2,Mu,count,MPIU_SCALAR,MPIU_SUM,PetscSubcommParent(subcomm));
234: PetscFree3(temp,temp2,ppk);
235: MatDestroy(&H);
236: return(0);
237: }
239: /*@
240: BVTraceQuadrature - Computes an estimate of the number of eigenvalues
241: inside a region via quantities computed in the quadrature rule of
242: contour integral methods.
244: Collective on Y
246: Input Parameters:
247: + Y - first basis vectors
248: . V - second basis vectors
249: . L - block size
250: . L_max - maximum block size
251: . w - quadrature weights
252: . scat - (optional) VecScatter object to communicate between subcommunicators
253: . subcomm - subcommunicator layout
254: . npoints - number of points to process by the subcommunicator
255: - useconj - whether conjugate points can be used or not
257: Output Parameter:
258: . est_eig - estimated eigenvalue count
260: Notes:
261: This function returns an estimation of the number of eigenvalues in the
262: region, computed as trace(V'*S_0), where S_0 is the first panel of S
263: computed by BVSumQuadrature().
265: When using subcommunicators, Y is stored in the subcommunicators for a subset
266: of integration points. In that case, the computation is done in the subcomm
267: and then scattered to the whole communicator in S using the VecScatter scat.
268: The value npoints is the number of points to be processed in this subcomm
269: and the flag useconj indicates whether symmetric points can be reused.
271: Level: developer
273: .seealso: BVScatter(), BVDotQuadrature(), BVSumQuadrature(), RGComputeQuadrature(), RGCanUseConjugates()
274: @*/
275: PetscErrorCode BVTraceQuadrature(BV Y,BV V,PetscInt L,PetscInt L_max,PetscScalar *w,VecScatter scat,PetscSubcomm subcomm,PetscInt npoints,PetscBool useconj,PetscReal *est_eig)276: {
278: PetscInt i,j;
279: Vec y,yall,vj;
280: PetscScalar dot,sum=0.0,one=1.0;
287: BVCreateVec(Y,&y);
288: BVCreateVec(V,&yall);
289: for (j=0;j<L;j++) {
290: VecSet(y,0.0);
291: for (i=0;i<npoints;i++) {
292: BVSetActiveColumns(Y,i*L_max+j,i*L_max+j+1);
293: BVMultVec(Y,w[p_id(i)],1.0,y,&one);
294: }
295: BVGetColumn(V,j,&vj);
296: if (scat) {
297: VecScatterBegin(scat,y,yall,ADD_VALUES,SCATTER_REVERSE);
298: VecScatterEnd(scat,y,yall,ADD_VALUES,SCATTER_REVERSE);
299: VecDot(vj,yall,&dot);
300: } else {
301: VecDot(vj,y,&dot);
302: }
303: BVRestoreColumn(V,j,&vj);
304: if (useconj) sum += 2.0*PetscRealPart(dot);
305: else sum += dot;
306: }
307: *est_eig = PetscAbsScalar(sum)/(PetscReal)L;
308: VecDestroy(&y);
309: VecDestroy(&yall);
310: return(0);
311: }
313: PetscErrorCode BVSVDAndRank_Refine(BV S,PetscReal delta,PetscScalar *pA,PetscReal *sigma,PetscInt *rank)314: {
316: PetscInt i,j,k,ml=S->k;
317: PetscMPIInt len;
318: PetscScalar *work,*B,*tempB,*sarray,*Q1,*Q2,*temp2,alpha=1.0,beta=0.0;
319: PetscBLASInt l,m,n,lda,ldu,ldvt,lwork,info,ldb,ldc;
320: #if defined(PETSC_USE_COMPLEX)
321: PetscReal *rwork;
322: #endif
325: BVGetArray(S,&sarray);
326: PetscMalloc6(ml*ml,&temp2,S->n*ml,&Q1,S->n*ml,&Q2,ml*ml,&B,ml*ml,&tempB,5*ml,&work);
327: #if defined(PETSC_USE_COMPLEX)
328: PetscMalloc1(5*ml,&rwork);
329: #endif
330: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
332: PetscArrayzero(B,ml*ml);
333: for (i=0;i<ml;i++) B[i*ml+i]=1;
335: for (k=0;k<2;k++) {
336: PetscBLASIntCast(S->n,&m);
337: PetscBLASIntCast(ml,&l);
338: n = l; lda = m; ldb = m; ldc = l;
339: if (!k) {
340: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&l,&n,&m,&alpha,sarray,&lda,sarray,&ldb,&beta,pA,&ldc));
341: } else {
342: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&l,&n,&m,&alpha,Q1,&lda,Q1,&ldb,&beta,pA,&ldc));
343: }
344: PetscArrayzero(temp2,ml*ml);
345: PetscMPIIntCast(ml*ml,&len);
346: MPIU_Allreduce(pA,temp2,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)S));
348: PetscBLASIntCast(ml,&m);
349: n = m; lda = m; lwork = 5*m, ldu = 1; ldvt = 1;
350: #if defined(PETSC_USE_COMPLEX)
351: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&m,&n,temp2,&lda,sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,rwork,&info));
352: #else
353: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&m,&n,temp2,&lda,sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,&info));
354: #endif
355: SlepcCheckLapackInfo("gesvd",info);
357: PetscBLASIntCast(S->n,&l);
358: PetscBLASIntCast(ml,&n);
359: m = n; lda = l; ldb = m; ldc = l;
360: if (!k) {
361: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&l,&n,&m,&alpha,sarray,&lda,temp2,&ldb,&beta,Q1,&ldc));
362: } else {
363: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&l,&n,&m,&alpha,Q1,&lda,temp2,&ldb,&beta,Q2,&ldc));
364: }
366: PetscBLASIntCast(ml,&l);
367: m = l; n = l; lda = l; ldb = m; ldc = l;
368: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&l,&n,&m,&alpha,B,&lda,temp2,&ldb,&beta,tempB,&ldc));
369: for (i=0;i<ml;i++) {
370: sigma[i] = PetscSqrtReal(sigma[i]);
371: for (j=0;j<S->n;j++) {
372: if (k%2) Q2[j+i*S->n] /= sigma[i];
373: else Q1[j+i*S->n] /= sigma[i];
374: }
375: for (j=0;j<ml;j++) B[j+i*ml] = tempB[j+i*ml]*sigma[i];
376: }
377: }
379: PetscBLASIntCast(ml,&m);
380: n = m; lda = m; ldu=1; ldvt=1;
381: #if defined(PETSC_USE_COMPLEX)
382: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","O",&m,&n,B,&lda,sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,rwork,&info));
383: #else
384: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","O",&m,&n,B,&lda,sigma,NULL,&ldu,NULL,&ldvt,work,&lwork,&info));
385: #endif
386: SlepcCheckLapackInfo("gesvd",info);
388: PetscBLASIntCast(S->n,&l);
389: PetscBLASIntCast(ml,&n);
390: m = n; lda = l; ldb = m; ldc = l;
391: if (k%2) {
392: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&l,&n,&m,&alpha,Q1,&lda,B,&ldb,&beta,sarray,&ldc));
393: } else {
394: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&l,&n,&m,&alpha,Q2,&lda,B,&ldb,&beta,sarray,&ldc));
395: }
397: PetscFPTrapPop();
398: BVRestoreArray(S,&sarray);
400: if (rank) {
401: (*rank) = 0;
402: for (i=0;i<ml;i++) {
403: if (sigma[i]/PetscMax(sigma[0],1.0)>delta) (*rank)++;
404: }
405: }
406: PetscFree6(temp2,Q1,Q2,B,tempB,work);
407: #if defined(PETSC_USE_COMPLEX)
408: PetscFree(rwork);
409: #endif
410: return(0);
411: }
413: PetscErrorCode BVSVDAndRank_QR(BV S,PetscReal delta,PetscScalar *pA,PetscReal *sigma,PetscInt *rank)414: {
416: PetscInt i,n,ml=S->k;
417: PetscBLASInt m,lda,lwork,info;
418: PetscScalar *work;
419: PetscReal *rwork;
420: Mat A;
421: Vec v;
424: /* Compute QR factorizaton of S */
425: BVGetSizes(S,NULL,&n,NULL);
426: n = PetscMin(n,ml);
427: BVSetActiveColumns(S,0,n);
428: PetscArrayzero(pA,ml*n);
429: MatCreateDense(PETSC_COMM_SELF,n,n,PETSC_DECIDE,PETSC_DECIDE,pA,&A);
430: BVOrthogonalize(S,A);
431: if (n<ml) {
432: /* the rest of the factorization */
433: for (i=n;i<ml;i++) {
434: BVGetColumn(S,i,&v);
435: BVOrthogonalizeVec(S,v,pA+i*n,NULL,NULL);
436: BVRestoreColumn(S,i,&v);
437: }
438: }
439: PetscBLASIntCast(n,&lda);
440: PetscBLASIntCast(ml,&m);
441: PetscMalloc2(5*ml,&work,5*ml,&rwork);
442: lwork = 5*m;
443: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
444: #if !defined (PETSC_USE_COMPLEX)
445: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&lda,&m,pA,&lda,sigma,NULL,&lda,NULL,&lda,work,&lwork,&info));
446: #else
447: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("O","N",&lda,&m,pA,&lda,sigma,NULL,&lda,NULL,&lda,work,&lwork,rwork,&info));
448: #endif
449: SlepcCheckLapackInfo("gesvd",info);
450: PetscFPTrapPop();
451: *rank = 0;
452: for (i=0;i<n;i++) {
453: if (sigma[i]/PetscMax(sigma[0],1)>delta) (*rank)++;
454: }
455: /* n first columns of A have the left singular vectors */
456: BVMultInPlace(S,A,0,*rank);
457: PetscFree2(work,rwork);
458: MatDestroy(&A);
459: return(0);
460: }
462: PetscErrorCode BVSVDAndRank_QR_CAA(BV S,PetscInt M,PetscInt L,PetscReal delta,PetscScalar *pA,PetscReal *sigma,PetscInt *rank)463: {
465: PetscInt i,j,n,ml=S->k;
466: PetscBLASInt m,k_,lda,lwork,info;
467: PetscScalar *work,*T,*U,*R,sone=1.0,zero=0.0;
468: PetscReal *rwork;
469: Mat A;
472: /* Compute QR factorizaton of S */
473: BVGetSizes(S,NULL,&n,NULL);
474: if (n<ml) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_SUP,"The QR_CAA method does not support problem size n < m*L");
475: BVSetActiveColumns(S,0,ml);
476: PetscArrayzero(pA,ml*ml);
477: MatCreateDense(PETSC_COMM_SELF,ml,ml,PETSC_DECIDE,PETSC_DECIDE,pA,&A);
478: BVOrthogonalize(S,A);
479: MatDestroy(&A);
481: /* SVD of first (M-1)*L diagonal block */
482: PetscBLASIntCast((M-1)*L,&m);
483: PetscMalloc5(m*m,&T,m*m,&R,m*m,&U,5*ml,&work,5*ml,&rwork);
484: for (j=0;j<m;j++) {
485: PetscArraycpy(R+j*m,pA+j*ml,m);
486: }
487: lwork = 5*m;
488: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
489: #if !defined (PETSC_USE_COMPLEX)
490: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","O",&m,&m,R,&m,sigma,U,&m,NULL,&m,work,&lwork,&info));
491: #else
492: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","O",&m,&m,R,&m,sigma,U,&m,NULL,&m,work,&lwork,rwork,&info));
493: #endif
494: SlepcCheckLapackInfo("gesvd",info);
495: PetscFPTrapPop();
496: *rank = 0;
497: for (i=0;i<m;i++) {
498: if (sigma[i]/PetscMax(sigma[0],1)>delta) (*rank)++;
499: }
500: MatCreateDense(PETSC_COMM_SELF,m,m,PETSC_DECIDE,PETSC_DECIDE,U,&A);
501: BVSetActiveColumns(S,0,m);
502: BVMultInPlace(S,A,0,*rank);
503: MatDestroy(&A);
504: /* Projected linear system */
505: /* m first columns of A have the right singular vectors */
506: PetscBLASIntCast(*rank,&k_);
507: PetscBLASIntCast(ml,&lda);
508: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&m,&k_,&m,&sone,pA+L*lda,&lda,R,&m,&zero,T,&m));
509: PetscArrayzero(pA,ml*ml);
510: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&k_,&k_,&m,&sone,U,&m,T,&m,&zero,pA,&k_));
511: for (j=0;j<k_;j++) for (i=0;i<k_;i++) pA[j*k_+i] /= sigma[j];
512: PetscFree5(T,R,U,work,rwork);
513: return(0);
514: }
516: /*@
517: BVSVDAndRank - Compute the SVD (left singular vectors only, and singular
518: values) and determine the numerical rank according to a tolerance.
520: Collective on S
522: Input Parameters:
523: + S - the basis vectors
524: . m - the moment degree
525: . l - the block size
526: . delta - the tolerance used to determine the rank
527: - meth - the method to be used
529: Output Parameters:
530: + A - workspace, on output contains relevant values in the CAA method
531: . sigma - computed singular values
532: - rank - estimated rank (optional)
534: Notes:
535: This function computes [U,Sigma,V] = svd(S) and replaces S with U.
536: The current implementation computes this via S'*S, and it may include
537: some kind of iterative refinement to improve accuracy in some cases.
539: The parameters m and l refer to the moment and block size of contour
540: integral methods. All columns up to m*l are modified, and the active
541: columns are set to 0..m*l.
543: The method is one of BV_SVD_METHOD_REFINE, BV_SVD_METHOD_QR, BV_SVD_METHOD_QR_CAA.
545: The A workspace should be m*l*m*l in size.
547: Once the decomposition is computed, the numerical rank is estimated
548: by counting the number of singular values that are larger than the
549: tolerance delta, relative to the first singular value.
551: Level: developer
553: .seealso: BVSetActiveColumns()
554: @*/
555: PetscErrorCode BVSVDAndRank(BV S,PetscInt m,PetscInt l,PetscReal delta,BVSVDMethod meth,PetscScalar *A,PetscReal *sigma,PetscInt *rank)556: {
569: PetscLogEventBegin(BV_SVDAndRank,S,0,0,0);
570: BVSetActiveColumns(S,0,m*l);
571: switch (meth) {
572: case BV_SVD_METHOD_REFINE:
573: BVSVDAndRank_Refine(S,delta,A,sigma,rank);
574: break;
575: case BV_SVD_METHOD_QR:
576: BVSVDAndRank_QR(S,delta,A,sigma,rank);
577: break;
578: case BV_SVD_METHOD_QR_CAA:
579: BVSVDAndRank_QR_CAA(S,m,l,delta,A,sigma,rank);
580: break;
581: }
582: PetscLogEventEnd(BV_SVDAndRank,S,0,0,0);
583: return(0);
584: }