Actual source code: fnutil.c

slepc-3.13.1 2020-04-12
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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:    Utility subroutines common to several impls
 12: */

 14: #include <slepc/private/fnimpl.h>      /*I "slepcfn.h" I*/
 15: #include <slepcblaslapack.h>

 17: /*
 18:    Compute the square root of an upper quasi-triangular matrix T,
 19:    using Higham's algorithm (LAA 88, 1987). T is overwritten with sqrtm(T).
 20:  */
 21: PetscErrorCode SlepcMatDenseSqrt(PetscBLASInt n,PetscScalar *T,PetscBLASInt ld)
 22: {
 23:   PetscScalar  one=1.0,mone=-1.0;
 24:   PetscReal    scal;
 25:   PetscBLASInt i,j,si,sj,r,ione=1,info;
 26: #if !defined(PETSC_USE_COMPLEX)
 27:   PetscReal    alpha,theta,mu,mu2;
 28: #endif

 31:   for (j=0;j<n;j++) {
 32: #if defined(PETSC_USE_COMPLEX)
 33:     sj = 1;
 34:     T[j+j*ld] = PetscSqrtScalar(T[j+j*ld]);
 35: #else
 36:     sj = (j==n-1 || T[j+1+j*ld] == 0.0)? 1: 2;
 37:     if (sj==1) {
 38:       if (T[j+j*ld]<0.0) SETERRQ(PETSC_COMM_SELF,1,"Matrix has a real negative eigenvalue, no real primary square root exists");
 39:       T[j+j*ld] = PetscSqrtReal(T[j+j*ld]);
 40:     } else {
 41:       /* square root of 2x2 block */
 42:       theta = (T[j+j*ld]+T[j+1+(j+1)*ld])/2.0;
 43:       mu    = (T[j+j*ld]-T[j+1+(j+1)*ld])/2.0;
 44:       mu2   = -mu*mu-T[j+1+j*ld]*T[j+(j+1)*ld];
 45:       mu    = PetscSqrtReal(mu2);
 46:       if (theta>0.0) alpha = PetscSqrtReal((theta+PetscSqrtReal(theta*theta+mu2))/2.0);
 47:       else alpha = mu/PetscSqrtReal(2.0*(-theta+PetscSqrtReal(theta*theta+mu2)));
 48:       T[j+j*ld]       /= 2.0*alpha;
 49:       T[j+1+(j+1)*ld] /= 2.0*alpha;
 50:       T[j+(j+1)*ld]   /= 2.0*alpha;
 51:       T[j+1+j*ld]     /= 2.0*alpha;
 52:       T[j+j*ld]       += alpha-theta/(2.0*alpha);
 53:       T[j+1+(j+1)*ld] += alpha-theta/(2.0*alpha);
 54:     }
 55: #endif
 56:     for (i=j-1;i>=0;i--) {
 57: #if defined(PETSC_USE_COMPLEX)
 58:       si = 1;
 59: #else
 60:       si = (i==0 || T[i+(i-1)*ld] == 0.0)? 1: 2;
 61:       if (si==2) i--;
 62: #endif
 63:       /* solve Sylvester equation of order si x sj */
 64:       r = j-i-si;
 65:       if (r) PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&si,&sj,&r,&mone,T+i+(i+si)*ld,&ld,T+i+si+j*ld,&ld,&one,T+i+j*ld,&ld));
 66:       PetscStackCallBLAS("LAPACKtrsyl",LAPACKtrsyl_("N","N",&ione,&si,&sj,T+i+i*ld,&ld,T+j+j*ld,&ld,T+i+j*ld,&ld,&scal,&info));
 67:       SlepcCheckLapackInfo("trsyl",info);
 68:       if (scal!=1.0) SETERRQ1(PETSC_COMM_SELF,1,"Current implementation cannot handle scale factor %g",scal);
 69:     }
 70:     if (sj==2) j++;
 71:   }
 72:   return(0);
 73: }

 75: #define BLOCKSIZE 64

 77: /*
 78:    Schur method for the square root of an upper quasi-triangular matrix T.
 79:    T is overwritten with sqrtm(T).
 80:    If firstonly then only the first column of T will contain relevant values.
 81:  */
 82: PetscErrorCode SlepcSqrtmSchur(PetscBLASInt n,PetscScalar *T,PetscBLASInt ld,PetscBool firstonly)
 83: {
 85:   PetscBLASInt   i,j,k,r,ione=1,sdim,lwork,*s,*p,info,bs=BLOCKSIZE;
 86:   PetscScalar    *wr,*W,*Q,*work,one=1.0,zero=0.0,mone=-1.0;
 87:   PetscInt       m,nblk;
 88:   PetscReal      scal;
 89: #if defined(PETSC_USE_COMPLEX)
 90:   PetscReal      *rwork;
 91: #else
 92:   PetscReal      *wi;
 93: #endif

 96:   m     = n;
 97:   nblk  = (m+bs-1)/bs;
 98:   lwork = 5*n;
 99:   k     = firstonly? 1: n;

101:   /* compute Schur decomposition A*Q = Q*T */
102: #if !defined(PETSC_USE_COMPLEX)
103:   PetscMalloc7(m,&wr,m,&wi,m*k,&W,m*m,&Q,lwork,&work,nblk,&s,nblk,&p);
104:   PetscStackCallBLAS("LAPACKgees",LAPACKgees_("V","N",NULL,&n,T,&ld,&sdim,wr,wi,Q,&ld,work,&lwork,NULL,&info));
105: #else
106:   PetscMalloc7(m,&wr,m,&rwork,m*k,&W,m*m,&Q,lwork,&work,nblk,&s,nblk,&p);
107:   PetscStackCallBLAS("LAPACKgees",LAPACKgees_("V","N",NULL,&n,T,&ld,&sdim,wr,Q,&ld,work,&lwork,rwork,NULL,&info));
108: #endif
109:   SlepcCheckLapackInfo("gees",info);

111:   /* determine block sizes and positions, to avoid cutting 2x2 blocks */
112:   j = 0;
113:   p[j] = 0;
114:   do {
115:     s[j] = PetscMin(bs,n-p[j]);
116: #if !defined(PETSC_USE_COMPLEX)
117:     if (p[j]+s[j]!=n && T[p[j]+s[j]+(p[j]+s[j]-1)*ld]!=0.0) s[j]++;
118: #endif
119:     if (p[j]+s[j]==n) break;
120:     j++;
121:     p[j] = p[j-1]+s[j-1];
122:   } while (1);
123:   nblk = j+1;

125:   for (j=0;j<nblk;j++) {
126:     /* evaluate f(T_jj) */
127:     SlepcMatDenseSqrt(s[j],T+p[j]+p[j]*ld,ld);
128:     for (i=j-1;i>=0;i--) {
129:       /* solve Sylvester equation for block (i,j) */
130:       r = p[j]-p[i]-s[i];
131:       if (r) PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",s+i,s+j,&r,&mone,T+p[i]+(p[i]+s[i])*ld,&ld,T+p[i]+s[i]+p[j]*ld,&ld,&one,T+p[i]+p[j]*ld,&ld));
132:       PetscStackCallBLAS("LAPACKtrsyl",LAPACKtrsyl_("N","N",&ione,s+i,s+j,T+p[i]+p[i]*ld,&ld,T+p[j]+p[j]*ld,&ld,T+p[i]+p[j]*ld,&ld,&scal,&info));
133:       SlepcCheckLapackInfo("trsyl",info);
134:       if (scal!=1.0) SETERRQ1(PETSC_COMM_SELF,1,"Current implementation cannot handle scale factor %g",scal);
135:     }
136:   }

138:   /* backtransform B = Q*T*Q' */
139:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&n,&k,&n,&one,T,&ld,Q,&ld,&zero,W,&ld));
140:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&k,&n,&one,Q,&ld,W,&ld,&zero,T,&ld));

142:   /* flop count: Schur decomposition, triangular square root, and backtransform */
143:   PetscLogFlops(25.0*n*n*n+n*n*n/3.0+4.0*n*n*k);

145: #if !defined(PETSC_USE_COMPLEX)
146:   PetscFree7(wr,wi,W,Q,work,s,p);
147: #else
148:   PetscFree7(wr,rwork,W,Q,work,s,p);
149: #endif
150:   return(0);
151: }

153: #define DBMAXIT 25

155: /*
156:    Computes the principal square root of the matrix T using the product form
157:    of the Denman-Beavers iteration.
158:    T is overwritten with sqrtm(T) or inv(sqrtm(T)) depending on flag inv.
159:  */
160: PetscErrorCode SlepcSqrtmDenmanBeavers(PetscBLASInt n,PetscScalar *T,PetscBLASInt ld,PetscBool inv)
161: {
162:   PetscScalar        *Told,*M=NULL,*invM,*work,work1,prod,alpha;
163:   PetscScalar        szero=0.0,sone=1.0,smone=-1.0,spfive=0.5,sp25=0.25;
164:   PetscReal          tol,Mres=0.0,detM,g,reldiff,fnormdiff,fnormT,rwork[1];
165:   PetscBLASInt       N,i,it,*piv=NULL,info,query=-1,lwork;
166:   const PetscBLASInt one=1;
167:   PetscBool          converged=PETSC_FALSE,scale=PETSC_FALSE;
168:   PetscErrorCode     ierr;
169:   unsigned int       ftz;

172:   N = n*n;
173:   tol = PetscSqrtReal((PetscReal)n)*PETSC_MACHINE_EPSILON/2;
174:   SlepcSetFlushToZero(&ftz);

176:   /* query work size */
177:   PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,M,&ld,piv,&work1,&query,&info));
178:   PetscBLASIntCast((PetscInt)PetscRealPart(work1),&lwork);
179:   PetscMalloc5(lwork,&work,n,&piv,n*n,&Told,n*n,&M,n*n,&invM);
180:   PetscArraycpy(M,T,n*n);

182:   if (inv) {  /* start recurrence with I instead of A */
183:     PetscArrayzero(T,n*n);
184:     for (i=0;i<n;i++) T[i+i*ld] += 1.0;
185:   }

187:   for (it=0;it<DBMAXIT && !converged;it++) {

189:     if (scale) {  /* g = (abs(det(M)))^(-1/(2*n)) */
190:       PetscArraycpy(invM,M,n*n);
191:       PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,invM,&ld,piv,&info));
192:       SlepcCheckLapackInfo("getrf",info);
193:       prod = invM[0];
194:       for (i=1;i<n;i++) prod *= invM[i+i*ld];
195:       detM = PetscAbsScalar(prod);
196:       g = PetscPowReal(detM,-1.0/(2.0*n));
197:       alpha = g;
198:       PetscStackCallBLAS("BLASscal",BLASscal_(&N,&alpha,T,&one));
199:       alpha = g*g;
200:       PetscStackCallBLAS("BLASscal",BLASscal_(&N,&alpha,M,&one));
201:       PetscLogFlops(2.0*n*n*n/3.0+2.0*n*n);
202:     }

204:     PetscArraycpy(Told,T,n*n);
205:     PetscArraycpy(invM,M,n*n);

207:     PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,invM,&ld,piv,&info));
208:     SlepcCheckLapackInfo("getrf",info);
209:     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,invM,&ld,piv,work,&lwork,&info));
210:     SlepcCheckLapackInfo("getri",info);
211:     PetscLogFlops(2.0*n*n*n/3.0+4.0*n*n*n/3.0);

213:     for (i=0;i<n;i++) invM[i+i*ld] += 1.0;
214:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&spfive,Told,&ld,invM,&ld,&szero,T,&ld));
215:     for (i=0;i<n;i++) invM[i+i*ld] -= 1.0;

217:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&N,&sone,invM,&one,M,&one));
218:     PetscStackCallBLAS("BLASscal",BLASscal_(&N,&sp25,M,&one));
219:     for (i=0;i<n;i++) M[i+i*ld] -= 0.5;
220:     PetscLogFlops(2.0*n*n*n+2.0*n*n);

222:     Mres = LAPACKlange_("F",&n,&n,M,&n,rwork);
223:     for (i=0;i<n;i++) M[i+i*ld] += 1.0;

225:     if (scale) {
226:       /* reldiff = norm(T - Told,'fro')/norm(T,'fro') */
227:       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&N,&smone,T,&one,Told,&one));
228:       fnormdiff = LAPACKlange_("F",&n,&n,Told,&n,rwork);
229:       fnormT = LAPACKlange_("F",&n,&n,T,&n,rwork);
230:       PetscLogFlops(7.0*n*n);
231:       reldiff = fnormdiff/fnormT;
232:       PetscInfo4(NULL,"it: %D reldiff: %g scale: %g tol*scale: %g\n",it,(double)reldiff,(double)g,(double)tol*g);
233:       if (reldiff<1e-2) scale = PETSC_FALSE;  /* Switch off scaling */
234:     }

236:     if (Mres<=tol) converged = PETSC_TRUE;
237:   }

239:   if (Mres>tol) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"SQRTM not converged after %d iterations",DBMAXIT);
240:   PetscFree5(work,piv,Told,M,invM);
241:   SlepcResetFlushToZero(&ftz);
242:   return(0);
243: }

245: #define NSMAXIT 50

247: /*
248:    Computes the principal square root of the matrix A using the Newton-Schulz iteration.
249:    T is overwritten with sqrtm(T) or inv(sqrtm(T)) depending on flag inv.
250:  */
251: PetscErrorCode SlepcSqrtmNewtonSchulz(PetscBLASInt n,PetscScalar *A,PetscBLASInt ld,PetscBool inv)
252: {
253:   PetscScalar        *Y=A,*Yold,*Z,*Zold,*M,alpha,sqrtnrm;
254:   PetscScalar        szero=0.0,sone=1.0,smone=-1.0,spfive=0.5,sthree=3.0;
255:   PetscReal          tol,Yres=0.0,nrm,rwork[1];
256:   PetscBLASInt       i,it,N;
257:   const PetscBLASInt one=1;
258:   PetscBool          converged=PETSC_FALSE;
259:   PetscErrorCode     ierr;
260:   unsigned int       ftz;

263:   N = n*n;
264:   tol = PetscSqrtReal((PetscReal)n)*PETSC_MACHINE_EPSILON/2;
265:   SlepcSetFlushToZero(&ftz);

267:   PetscMalloc4(N,&Yold,N,&Z,N,&Zold,N,&M);

269:   /* scale A so that ||I-A|| < 1 */
270:   PetscArraycpy(Z,A,N);
271:   for (i=0;i<n;i++) Z[i+i*ld] -= 1.0;
272:   nrm = LAPACKlange_("fro",&n,&n,Z,&n,rwork);
273:   sqrtnrm = PetscSqrtReal(nrm);
274:   alpha = 1.0/nrm;
275:   PetscStackCallBLAS("BLASscal",BLASscal_(&N,&alpha,A,&one));
276:   tol *= nrm;
277:   PetscInfo2(NULL,"||I-A||_F = %g, new tol: %g\n",(double)nrm,(double)tol);
278:   PetscLogFlops(2.0*n*n);

280:   /* Z = I */
281:   PetscArrayzero(Z,N);
282:   for (i=0;i<n;i++) Z[i+i*ld] = 1.0;

284:   for (it=0;it<NSMAXIT && !converged;it++) {
285:     /* Yold = Y, Zold = Z */
286:     PetscArraycpy(Yold,Y,N);
287:     PetscArraycpy(Zold,Z,N);

289:     /* M = (3*I-Zold*Yold) */
290:     PetscArrayzero(M,N);
291:     for (i=0;i<n;i++) M[i+i*ld] = sthree;
292:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&smone,Zold,&ld,Yold,&ld,&sone,M,&ld));

294:     /* Y = (1/2)*Yold*M, Z = (1/2)*M*Zold */
295:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&spfive,Yold,&ld,M,&ld,&szero,Y,&ld));
296:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&spfive,M,&ld,Zold,&ld,&szero,Z,&ld));

298:     /* reldiff = norm(Y-Yold,'fro')/norm(Y,'fro') */
299:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&N,&smone,Y,&one,Yold,&one));
300:     Yres = LAPACKlange_("fro",&n,&n,Yold,&n,rwork);
301:     PetscIsNanReal(Yres);
302:     if (Yres<=tol) converged = PETSC_TRUE;
303:     PetscInfo2(NULL,"it: %D res: %g\n",it,(double)Yres);

305:     PetscLogFlops(6.0*n*n*n+2.0*n*n);
306:   }

308:   if (Yres>tol) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"SQRTM not converged after %d iterations",NSMAXIT);

310:   /* undo scaling */
311:   if (inv) {
312:     PetscArraycpy(A,Z,N);
313:     sqrtnrm = 1.0/sqrtnrm;
314:     PetscStackCallBLAS("BLASscal",BLASscal_(&N,&sqrtnrm,A,&one));
315:   } else PetscStackCallBLAS("BLASscal",BLASscal_(&N,&sqrtnrm,A,&one));

317:   PetscFree4(Yold,Z,Zold,M);
318:   SlepcResetFlushToZero(&ftz);
319:   return(0);
320: }

322: #define ITMAX 5
323: #define SWAP(a,b,t) {t=a;a=b;b=t;}

325: /*
326:    Estimate norm(A^m,1) by block 1-norm power method (required workspace is 11*n)
327: */
328: static PetscErrorCode SlepcNormEst1(PetscBLASInt n,PetscScalar *A,PetscInt m,PetscScalar *work,PetscRandom rand,PetscReal *nrm)
329: {
330:   PetscScalar    *X,*Y,*Z,*S,*S_old,*aux,val,sone=1.0,szero=0.0;
331:   PetscReal      est=0.0,est_old,vals[2]={0.0,0.0},*zvals,maxzval[2],raux;
332:   PetscBLASInt   i,j,t=2,it=0,ind[2],est_j=0,m1;

336:   X = work;
337:   Y = work + 2*n;
338:   Z = work + 4*n;
339:   S = work + 6*n;
340:   S_old = work + 8*n;
341:   zvals = (PetscReal*)(work + 10*n);

343:   for (i=0;i<n;i++) {  /* X has columns of unit 1-norm */
344:     X[i] = 1.0/n;
345:     PetscRandomGetValue(rand,&val);
346:     if (PetscRealPart(val) < 0.5) X[i+n] = -1.0/n;
347:     else X[i+n] = 1.0/n;
348:   }
349:   for (i=0;i<t*n;i++) S[i] = 0.0;
350:   ind[0] = 0; ind[1] = 0;
351:   est_old = 0;
352:   while (1) {
353:     it++;
354:     for (j=0;j<m;j++) {  /* Y = A^m*X */
355:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&t,&n,&sone,A,&n,X,&n,&szero,Y,&n));
356:       if (j<m-1) SWAP(X,Y,aux);
357:     }
358:     for (j=0;j<t;j++) {  /* vals[j] = norm(Y(:,j),1) */
359:       vals[j] = 0.0;
360:       for (i=0;i<n;i++) vals[j] += PetscAbsScalar(Y[i+j*n]);
361:     }
362:     if (vals[0]<vals[1]) {
363:       SWAP(vals[0],vals[1],raux);
364:       m1 = 1;
365:     } else m1 = 0;
366:     est = vals[0];
367:     if (est>est_old || it==2) est_j = ind[m1];
368:     if (it>=2 && est<=est_old) {
369:       est = est_old;
370:       break;
371:     }
372:     est_old = est;
373:     if (it>ITMAX) break;
374:     SWAP(S,S_old,aux);
375:     for (i=0;i<t*n;i++) {  /* S = sign(Y) */
376:       S[i] = (PetscRealPart(Y[i]) < 0.0)? -1.0: 1.0;
377:     }
378:     for (j=0;j<m;j++) {  /* Z = (A^T)^m*S */
379:       PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n,&t,&n,&sone,A,&n,S,&n,&szero,Z,&n));
380:       if (j<m-1) SWAP(S,Z,aux);
381:     }
382:     maxzval[0] = -1; maxzval[1] = -1;
383:     ind[0] = 0; ind[1] = 0;
384:     for (i=0;i<n;i++) {  /* zvals[i] = norm(Z(i,:),inf) */
385:       zvals[i] = PetscMax(PetscAbsScalar(Z[i+0*n]),PetscAbsScalar(Z[i+1*n]));
386:       if (zvals[i]>maxzval[0]) {
387:         maxzval[0] = zvals[i];
388:         ind[0] = i;
389:       } else if (zvals[i]>maxzval[1]) {
390:         maxzval[1] = zvals[i];
391:         ind[1] = i;
392:       }
393:     }
394:     if (it>=2 && maxzval[0]==zvals[est_j]) break;
395:     for (i=0;i<t*n;i++) X[i] = 0.0;
396:     for (j=0;j<t;j++) X[ind[j]+j*n] = 1.0;
397:   }
398:   *nrm = est;
399:   /* Flop count is roughly (it * 2*m * t*gemv) = 4*its*m*t*n*n */
400:   PetscLogFlops(4.0*it*m*t*n*n);
401:   return(0);
402: }

404: #define SMALLN 100

406: /*
407:    Estimate norm(A^m,1) (required workspace is 2*n*n)
408: */
409: PetscErrorCode SlepcNormAm(PetscBLASInt n,PetscScalar *A,PetscInt m,PetscScalar *work,PetscRandom rand,PetscReal *nrm)
410: {
411:   PetscScalar    *v=work,*w=work+n*n,*aux,sone=1.0,szero=0.0;
412:   PetscReal      rwork[1],tmp;
413:   PetscBLASInt   i,j,one=1;
414:   PetscBool      isrealpos=PETSC_TRUE;

418:   if (n<SMALLN) {   /* compute matrix power explicitly */
419:     if (m==1) {
420:       *nrm = LAPACKlange_("O",&n,&n,A,&n,rwork);
421:       PetscLogFlops(1.0*n*n);
422:     } else {  /* m>=2 */
423:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&sone,A,&n,A,&n,&szero,v,&n));
424:       for (j=0;j<m-2;j++) {
425:         PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&n,&n,&sone,A,&n,v,&n,&szero,w,&n));
426:         SWAP(v,w,aux);
427:       }
428:       *nrm = LAPACKlange_("O",&n,&n,v,&n,rwork);
429:       PetscLogFlops(2.0*n*n*n*(m-1)+1.0*n*n);
430:     }
431:   } else {
432:     for (i=0;i<n;i++)
433:       for (j=0;j<n;j++)
434: #if defined(PETSC_USE_COMPLEX)
435:         if (PetscRealPart(A[i+j*n])<0.0 || PetscImaginaryPart(A[i+j*n])!=0.0) { isrealpos = PETSC_FALSE; break; }
436: #else
437:         if (A[i+j*n]<0.0) { isrealpos = PETSC_FALSE; break; }
438: #endif
439:     if (isrealpos) {   /* for positive matrices only */
440:       for (i=0;i<n;i++) v[i] = 1.0;
441:       for (j=0;j<m;j++) {  /* w = A'*v */
442:         PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&sone,A,&n,v,&one,&szero,w,&one));
443:         SWAP(v,w,aux);
444:       }
445:       PetscLogFlops(2.0*n*n*m);
446:       *nrm = 0.0;
447:       for (i=0;i<n;i++) if ((tmp = PetscAbsScalar(v[i])) > *nrm) *nrm = tmp;   /* norm(v,inf) */
448:     } else {
449:       SlepcNormEst1(n,A,m,work,rand,nrm);
450:     }
451:   }
452:   return(0);
453: }