Actual source code: dspriv.c

slepc-3.16.3 2022-04-11
Report Typos and Errors
  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:    Private DS routines
 12: */

 14: #include <slepc/private/dsimpl.h>
 15: #include <slepcblaslapack.h>

 17: PetscErrorCode DSAllocateMatrix_Private(DS ds,DSMatType m,PetscBool isreal)
 18: {
 19:   size_t         sz;
 20:   PetscInt       n,d,nelem;
 21:   PetscBool      ispep,isnep;

 25:   PetscObjectTypeCompare((PetscObject)ds,DSPEP,&ispep);
 26:   PetscObjectTypeCompare((PetscObject)ds,DSNEP,&isnep);
 27:   if (ispep) {
 28:     DSPEPGetDegree(ds,&d);
 29:   }
 30:   if (isnep) {
 31:     DSNEPGetMinimality(ds,&d);
 32:   }
 33:   if ((ispep || isnep) && (m==DS_MAT_A || m==DS_MAT_B || m==DS_MAT_W || m==DS_MAT_U || m==DS_MAT_X || m==DS_MAT_Y)) n = d*ds->ld;
 34:   else n = ds->ld;

 36:   switch (m) {
 37:     case DS_MAT_T:
 38:       nelem = 3*ds->ld;
 39:       break;
 40:     case DS_MAT_D:
 41:       nelem = ds->ld;
 42:       break;
 43:     case DS_MAT_X:
 44:       nelem = ds->ld*n;
 45:       break;
 46:     case DS_MAT_Y:
 47:       nelem = ds->ld*n;
 48:       break;
 49:     default:
 50:       nelem = n*n;
 51:   }
 52:   if (isreal) {
 53:     sz = nelem*sizeof(PetscReal);
 54:     if (ds->rmat[m]) {
 55:       PetscFree(ds->rmat[m]);
 56:     } else {
 57:       PetscLogObjectMemory((PetscObject)ds,sz);
 58:     }
 59:     PetscCalloc1(nelem,&ds->rmat[m]);
 60:   } else {
 61:     sz = nelem*sizeof(PetscScalar);
 62:     if (ds->mat[m]) {
 63:       PetscFree(ds->mat[m]);
 64:     } else {
 65:       PetscLogObjectMemory((PetscObject)ds,sz);
 66:     }
 67:     PetscCalloc1(nelem,&ds->mat[m]);
 68:   }
 69:   return(0);
 70: }

 72: PetscErrorCode DSAllocateWork_Private(DS ds,PetscInt s,PetscInt r,PetscInt i)
 73: {

 77:   if (s>ds->lwork) {
 78:     PetscFree(ds->work);
 79:     PetscMalloc1(s,&ds->work);
 80:     PetscLogObjectMemory((PetscObject)ds,(s-ds->lwork)*sizeof(PetscScalar));
 81:     ds->lwork = s;
 82:   }
 83:   if (r>ds->lrwork) {
 84:     PetscFree(ds->rwork);
 85:     PetscMalloc1(r,&ds->rwork);
 86:     PetscLogObjectMemory((PetscObject)ds,(r-ds->lrwork)*sizeof(PetscReal));
 87:     ds->lrwork = r;
 88:   }
 89:   if (i>ds->liwork) {
 90:     PetscFree(ds->iwork);
 91:     PetscMalloc1(i,&ds->iwork);
 92:     PetscLogObjectMemory((PetscObject)ds,(i-ds->liwork)*sizeof(PetscBLASInt));
 93:     ds->liwork = i;
 94:   }
 95:   return(0);
 96: }

 98: /*@C
 99:    DSViewMat - Prints one of the internal DS matrices.

101:    Collective on ds

103:    Input Parameters:
104: +  ds     - the direct solver context
105: .  viewer - visualization context
106: -  m      - matrix to display

108:    Note:
109:    Works only for ascii viewers. Set the viewer in Matlab format if
110:    want to paste into Matlab.

112:    Level: developer

114: .seealso: DSView()
115: @*/
116: PetscErrorCode DSViewMat(DS ds,PetscViewer viewer,DSMatType m)
117: {
118:   PetscErrorCode    ierr;
119:   PetscInt          i,j,rows,cols;
120:   PetscScalar       *v;
121:   PetscViewerFormat format;
122: #if defined(PETSC_USE_COMPLEX)
123:   PetscBool         allreal = PETSC_TRUE;
124: #endif

129:   DSCheckValidMat(ds,m,3);
130:   if (!viewer) {
131:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ds),&viewer);
132:   }
135:   PetscViewerGetFormat(viewer,&format);
136:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
137:   PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
138:   DSMatGetSize(ds,m,&rows,&cols);
139: #if defined(PETSC_USE_COMPLEX)
140:   /* determine if matrix has all real values */
141:   v = ds->mat[m];
142:   for (i=0;i<rows;i++)
143:     for (j=0;j<cols;j++)
144:       if (PetscImaginaryPart(v[i+j*ds->ld])) { allreal = PETSC_FALSE; break; }
145: #endif
146:   if (format == PETSC_VIEWER_ASCII_MATLAB) {
147:     PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",rows,cols);
148:     PetscViewerASCIIPrintf(viewer,"%s = [\n",DSMatName[m]);
149:   } else {
150:     PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[m]);
151:   }

153:   for (i=0;i<rows;i++) {
154:     v = ds->mat[m]+i;
155:     for (j=0;j<cols;j++) {
156: #if defined(PETSC_USE_COMPLEX)
157:       if (allreal) {
158:         PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v));
159:       } else {
160:         PetscViewerASCIIPrintf(viewer,"%18.16e%+18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));
161:       }
162: #else
163:       PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v);
164: #endif
165:       v += ds->ld;
166:     }
167:     PetscViewerASCIIPrintf(viewer,"\n");
168:   }

170:   if (format == PETSC_VIEWER_ASCII_MATLAB) {
171:     PetscViewerASCIIPrintf(viewer,"];\n");
172:   }
173:   PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
174:   PetscViewerFlush(viewer);
175:   return(0);
176: }

178: PetscErrorCode DSSortEigenvalues_Private(DS ds,PetscScalar *wr,PetscScalar *wi,PetscInt *perm,PetscBool isghiep)
179: {
181:   PetscScalar    re,im,wi0;
182:   PetscInt       n,i,j,result,tmp1,tmp2=0,d=1;

185:   n = ds->t;   /* sort only first t pairs if truncated */
186:   /* insertion sort */
187:   i=ds->l+1;
188: #if !defined(PETSC_USE_COMPLEX)
189:   if (wi && wi[perm[i-1]]!=0.0) i++; /* initial value is complex */
190: #else
191:   if (isghiep && PetscImaginaryPart(wr[perm[i-1]])!=0.0) i++;
192: #endif
193:   for (;i<n;i+=d) {
194:     re = wr[perm[i]];
195:     if (wi) im = wi[perm[i]];
196:     else im = 0.0;
197:     tmp1 = perm[i];
198: #if !defined(PETSC_USE_COMPLEX)
199:     if (im!=0.0) { d = 2; tmp2 = perm[i+1]; }
200:     else d = 1;
201: #else
202:     if (isghiep && PetscImaginaryPart(re)!=0.0) { d = 2; tmp2 = perm[i+1]; }
203:     else d = 1;
204: #endif
205:     j = i-1;
206:     if (wi) wi0 = wi[perm[j]];
207:     else wi0 = 0.0;
208:     SlepcSCCompare(ds->sc,re,im,wr[perm[j]],wi0,&result);
209:     while (result<0 && j>=ds->l) {
210:       perm[j+d] = perm[j];
211:       j--;
212: #if !defined(PETSC_USE_COMPLEX)
213:       if (wi && wi[perm[j+1]]!=0)
214: #else
215:       if (isghiep && PetscImaginaryPart(wr[perm[j+1]])!=0)
216: #endif
217:         { perm[j+d] = perm[j]; j--; }

219:       if (j>=ds->l) {
220:         if (wi) wi0 = wi[perm[j]];
221:         else wi0 = 0.0;
222:         SlepcSCCompare(ds->sc,re,im,wr[perm[j]],wi0,&result);
223:       }
224:     }
225:     perm[j+1] = tmp1;
226:     if (d==2) perm[j+2] = tmp2;
227:   }
228:   return(0);
229: }

231: PetscErrorCode DSSortEigenvaluesReal_Private(DS ds,PetscReal *eig,PetscInt *perm)
232: {
234:   PetscScalar    re;
235:   PetscInt       i,j,result,tmp,l,n;

238:   n = ds->t;   /* sort only first t pairs if truncated */
239:   l = ds->l;
240:   /* insertion sort */
241:   for (i=l+1;i<n;i++) {
242:     re = eig[perm[i]];
243:     j = i-1;
244:     SlepcSCCompare(ds->sc,re,0.0,eig[perm[j]],0.0,&result);
245:     while (result<0 && j>=l) {
246:       tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
247:       if (j>=l) {
248:         SlepcSCCompare(ds->sc,re,0.0,eig[perm[j]],0.0,&result);
249:       }
250:     }
251:   }
252:   return(0);
253: }

255: /*
256:   DSCopyMatrix_Private - Copies the trailing block of a matrix (from
257:   rows/columns l to n).
258: */
259: PetscErrorCode DSCopyMatrix_Private(DS ds,DSMatType dst,DSMatType src)
260: {
262:   PetscInt    j,m,off,ld;
263:   PetscScalar *S,*D;

266:   ld  = ds->ld;
267:   m   = ds->n-ds->l;
268:   off = ds->l+ds->l*ld;
269:   S   = ds->mat[src];
270:   D   = ds->mat[dst];
271:   for (j=0;j<m;j++) {
272:     PetscArraycpy(D+off+j*ld,S+off+j*ld,m);
273:   }
274:   return(0);
275: }

277: /*
278:   Permute comumns [istart..iend-1] of [mat] according to perm. Columns have length n
279:  */
280: PetscErrorCode DSPermuteColumns_Private(DS ds,PetscInt istart,PetscInt iend,PetscInt n,DSMatType mat,PetscInt *perm)
281: {
282:   PetscInt    i,j,k,p,ld;
283:   PetscScalar *Q,rtmp;

286:   ld = ds->ld;
287:   Q  = ds->mat[mat];
288:   for (i=istart;i<iend;i++) {
289:     p = perm[i];
290:     if (p != i) {
291:       j = i + 1;
292:       while (perm[j] != i) j++;
293:       perm[j] = p; perm[i] = i;
294:       /* swap columns i and j */
295:       for (k=0;k<n;k++) {
296:         rtmp = Q[k+p*ld]; Q[k+p*ld] = Q[k+i*ld]; Q[k+i*ld] = rtmp;
297:       }
298:     }
299:   }
300:   return(0);
301: }

303: /*
304:   The same as DSPermuteColumns_Private but for two matrices [mat1] and [mat2]
305:  */
306: PetscErrorCode DSPermuteColumnsTwo_Private(DS ds,PetscInt istart,PetscInt iend,PetscInt n,DSMatType mat1,DSMatType mat2,PetscInt *perm)
307: {
308:   PetscInt    i,j,k,p,ld;
309:   PetscScalar *Q,*Z,rtmp,rtmp2;

312:   ld = ds->ld;
313:   Q  = ds->mat[mat1];
314:   Z  = ds->mat[mat2];
315:   for (i=istart;i<iend;i++) {
316:     p = perm[i];
317:     if (p != i) {
318:       j = i + 1;
319:       while (perm[j] != i) j++;
320:       perm[j] = p; perm[i] = i;
321:       /* swap columns i and j */
322:       for (k=0;k<n;k++) {
323:         rtmp  = Q[k+p*ld]; Q[k+p*ld] = Q[k+i*ld]; Q[k+i*ld] = rtmp;
324:         rtmp2 = Z[k+p*ld]; Z[k+p*ld] = Z[k+i*ld]; Z[k+i*ld] = rtmp2;
325:       }
326:     }
327:   }
328:   return(0);
329: }

331: /*
332:   Permute rows [istart..iend-1] of [mat] according to perm. Rows have length m
333:  */
334: PetscErrorCode DSPermuteRows_Private(DS ds,PetscInt istart,PetscInt iend,PetscInt m,DSMatType mat,PetscInt *perm)
335: {
336:   PetscInt    i,j,k,p,ld;
337:   PetscScalar *Q,rtmp;

340:   ld = ds->ld;
341:   Q  = ds->mat[mat];
342:   for (i=istart;i<iend;i++) {
343:     p = perm[i];
344:     if (p != i) {
345:       j = i + 1;
346:       while (perm[j] != i) j++;
347:       perm[j] = p; perm[i] = i;
348:       /* swap rows i and j */
349:       for (k=0;k<m;k++) {
350:         rtmp = Q[p+k*ld]; Q[p+k*ld] = Q[i+k*ld]; Q[i+k*ld] = rtmp;
351:       }
352:     }
353:   }
354:   return(0);
355: }

357: /*
358:   Permute columns [istart..iend-1] of [mat1] and [mat2] according to perm.
359:   Columns of [mat1] have length n, columns of [mat2] have length m
360:  */
361: PetscErrorCode DSPermuteBoth_Private(DS ds,PetscInt istart,PetscInt iend,PetscInt n,PetscInt m,DSMatType mat1,DSMatType mat2,PetscInt *perm)
362: {
363:   PetscInt    i,j,k,p,ld;
364:   PetscScalar *U,*V,rtmp;

367:   ld = ds->ld;
368:   U  = ds->mat[mat1];
369:   V  = ds->mat[mat2];
370:   for (i=istart;i<iend;i++) {
371:     p = perm[i];
372:     if (p != i) {
373:       j = i + 1;
374:       while (perm[j] != i) j++;
375:       perm[j] = p; perm[i] = i;
376:       /* swap columns i and j of U */
377:       for (k=0;k<n;k++) {
378:         rtmp = U[k+p*ld]; U[k+p*ld] = U[k+i*ld]; U[k+i*ld] = rtmp;
379:       }
380:       /* swap columns i and j of V */
381:       for (k=0;k<m;k++) {
382:         rtmp = V[k+p*ld]; V[k+p*ld] = V[k+i*ld]; V[k+i*ld] = rtmp;
383:       }
384:     }
385:   }
386:   return(0);
387: }

389: /*@
390:    DSSetIdentity - Copy the identity (a diagonal matrix with ones) on the
391:    active part of a matrix.

393:    Logically Collective on ds

395:    Input Parameters:
396: +  ds  - the direct solver context
397: -  mat - the matrix to modify

399:    Level: intermediate
400: @*/
401: PetscErrorCode DSSetIdentity(DS ds,DSMatType mat)
402: {
404:   PetscScalar    *x;
405:   PetscInt       i,ld,n,l;

410:   DSCheckValidMat(ds,mat,2);

412:   DSGetDimensions(ds,&n,&l,NULL,NULL);
413:   DSGetLeadingDimension(ds,&ld);
414:   PetscLogEventBegin(DS_Other,ds,0,0,0);
415:   DSGetArray(ds,mat,&x);
416:   PetscArrayzero(&x[ld*l],ld*(n-l));
417:   for (i=l;i<n;i++) x[i+i*ld] = 1.0;
418:   DSRestoreArray(ds,mat,&x);
419:   PetscLogEventEnd(DS_Other,ds,0,0,0);
420:   return(0);
421: }

423: /*@C
424:    DSOrthogonalize - Orthogonalize the columns of a matrix.

426:    Logically Collective on ds

428:    Input Parameters:
429: +  ds   - the direct solver context
430: .  mat  - a matrix
431: -  cols - number of columns to orthogonalize (starting from column zero)

433:    Output Parameter:
434: .  lindcols - (optional) number of linearly independent columns of the matrix

436:    Level: developer

438: .seealso: DSPseudoOrthogonalize()
439: @*/
440: PetscErrorCode DSOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscInt *lindcols)
441: {
443:   PetscInt       n,l,ld;
444:   PetscBLASInt   ld_,rA,cA,info,ltau,lw;
445:   PetscScalar    *A,*tau,*w,saux,dummy;

449:   DSCheckAlloc(ds,1);
451:   DSCheckValidMat(ds,mat,2);

454:   DSGetDimensions(ds,&n,&l,NULL,NULL);
455:   DSGetLeadingDimension(ds,&ld);
456:   n = n - l;
457:   if (cols > n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns");
458:   if (n == 0 || cols == 0) return(0);

460:   PetscLogEventBegin(DS_Other,ds,0,0,0);
461:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
462:   DSGetArray(ds,mat,&A);
463:   PetscBLASIntCast(PetscMin(cols,n),&ltau);
464:   PetscBLASIntCast(ld,&ld_);
465:   PetscBLASIntCast(n,&rA);
466:   PetscBLASIntCast(cols,&cA);
467:   lw = -1;
468:   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,A,&ld_,&dummy,&saux,&lw,&info));
469:   SlepcCheckLapackInfo("geqrf",info);
470:   lw = (PetscBLASInt)PetscRealPart(saux);
471:   DSAllocateWork_Private(ds,lw+ltau,0,0);
472:   tau = ds->work;
473:   w = &tau[ltau];
474:   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,&A[ld*l+l],&ld_,tau,w,&lw,&info));
475:   SlepcCheckLapackInfo("geqrf",info);
476:   PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&rA,&ltau,&ltau,&A[ld*l+l],&ld_,tau,w,&lw,&info));
477:   SlepcCheckLapackInfo("orgqr",info);
478:   if (lindcols) *lindcols = ltau;

480:   PetscFPTrapPop();
481:   PetscLogEventEnd(DS_Other,ds,0,0,0);
482:   DSRestoreArray(ds,mat,&A);
483:   PetscObjectStateIncrease((PetscObject)ds);
484:   return(0);
485: }

487: /*
488:   Compute C <- a*A*B + b*C, where
489:     ldC, the leading dimension of C,
490:     ldA, the leading dimension of A,
491:     rA, cA, rows and columns of A,
492:     At, if true use the transpose of A instead,
493:     ldB, the leading dimension of B,
494:     rB, cB, rows and columns of B,
495:     Bt, if true use the transpose of B instead
496: */
497: static PetscErrorCode SlepcMatDenseMult(PetscScalar *C,PetscInt _ldC,PetscScalar b,PetscScalar a,const PetscScalar *A,PetscInt _ldA,PetscInt rA,PetscInt cA,PetscBool At,const PetscScalar *B,PetscInt _ldB,PetscInt rB,PetscInt cB,PetscBool Bt)
498: {
500:   PetscInt       tmp;
501:   PetscBLASInt   m, n, k, ldA = _ldA, ldB = _ldB, ldC = _ldC;
502:   const char     *N = "N", *T = "C", *qA = N, *qB = N;

505:   if ((rA == 0) || (cB == 0)) return(0);

510:   /* Transpose if needed */
511:   if (At) tmp = rA, rA = cA, cA = tmp, qA = T;
512:   if (Bt) tmp = rB, rB = cB, cB = tmp, qB = T;

514:   /* Check size */
515:   if (cA != rB) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Matrix dimensions do not match");

517:   /* Do stub */
518:   if ((rA == 1) && (cA == 1) && (cB == 1)) {
519:     if (!At && !Bt) *C = *A * *B;
520:     else if (At && !Bt) *C = PetscConj(*A) * *B;
521:     else if (!At && Bt) *C = *A * PetscConj(*B);
522:     else *C = PetscConj(*A) * PetscConj(*B);
523:     m = n = k = 1;
524:   } else {
525:     m = rA; n = cB; k = cA;
526:     PetscStackCallBLAS("BLASgemm",BLASgemm_(qA,qB,&m,&n,&k,&a,(PetscScalar*)A,&ldA,(PetscScalar*)B,&ldB,&b,C,&ldC));
527:   }

529:   PetscLogFlops(2.0*m*n*k);
530:   return(0);
531: }

533: /*@C
534:    DSPseudoOrthogonalize - Orthogonalize the columns of a matrix with Modified
535:    Gram-Schmidt in an indefinite inner product space defined by a signature.

537:    Logically Collective on ds

539:    Input Parameters:
540: +  ds   - the direct solver context
541: .  mat  - the matrix
542: .  cols - number of columns to orthogonalize (starting from column zero)
543: -  s    - the signature that defines the inner product

545:    Output Parameters:
546: +  lindcols - (optional) linearly independent columns of the matrix
547: -  ns   - (optional) the new signature of the vectors

549:    Note:
550:    After the call the matrix satisfies A'*s*A = ns.

552:    Level: developer

554: .seealso: DSOrthogonalize()
555: @*/
556: PetscErrorCode DSPseudoOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscReal *s,PetscInt *lindcols,PetscReal *ns)
557: {
559:   PetscInt       i,j,k,l,n,ld;
560:   PetscBLASInt   info,one=1,zero=0,rA_,ld_;
561:   PetscScalar    *A,*A_,*m,*h,nr0;
562:   PetscReal      nr_o,nr,nr_abs,*ns_,done=1.0;

566:   DSCheckAlloc(ds,1);
568:   DSCheckValidMat(ds,mat,2);
571:   DSGetDimensions(ds,&n,&l,NULL,NULL);
572:   DSGetLeadingDimension(ds,&ld);
573:   n = n - l;
574:   if (cols > n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns");
575:   if (n == 0 || cols == 0) return(0);
576:   PetscBLASIntCast(n,&rA_);
577:   PetscBLASIntCast(ld,&ld_);
578:   DSGetArray(ds,mat,&A_);
579:   A = &A_[ld*l+l];
580:   DSAllocateWork_Private(ds,n+cols,ns?0:cols,0);
581:   m = ds->work;
582:   h = &m[n];
583:   ns_ = ns ? ns : ds->rwork;
584:   PetscLogEventBegin(DS_Other,ds,0,0,0);
585:   for (i=0; i<cols; i++) {
586:     /* m <- diag(s)*A[i] */
587:     for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
588:     /* nr_o <- mynorm(A[i]'*m), mynorm(x) = sign(x)*sqrt(|x|) */
589:     SlepcMatDenseMult(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
590:     nr = nr_o = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
591:     for (j=0; j<3 && i>0; j++) {
592:       /* h <- A[0:i-1]'*m */
593:       SlepcMatDenseMult(h,i,0.0,1.0,A,ld,n,i,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
594:       /* h <- diag(ns)*h */
595:       for (k=0; k<i; k++) h[k] *= ns_[k];
596:       /* A[i] <- A[i] - A[0:i-1]*h */
597:       SlepcMatDenseMult(&A[ld*i],ld,1.0,-1.0,A,ld,n,i,PETSC_FALSE,h,i,i,1,PETSC_FALSE);
598:       /* m <- diag(s)*A[i] */
599:       for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
600:       /* nr_o <- mynorm(A[i]'*m) */
601:       SlepcMatDenseMult(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
602:       nr = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
603:       if (PetscAbs(nr) < PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Linear dependency detected");
604:       if (PetscAbs(nr) > 0.7*PetscAbs(nr_o)) break;
605:       nr_o = nr;
606:     }
607:     ns_[i] = PetscSign(nr);
608:     /* A[i] <- A[i]/|nr| */
609:     nr_abs = PetscAbs(nr);
610:     PetscStackCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&nr_abs,&done,&rA_,&one,A+i*ld,&ld_,&info));
611:     SlepcCheckLapackInfo("lascl",info);
612:   }
613:   PetscLogEventEnd(DS_Other,ds,0,0,0);
614:   DSRestoreArray(ds,mat,&A_);
615:   PetscObjectStateIncrease((PetscObject)ds);
616:   if (lindcols) *lindcols = cols;
617:   return(0);
618: }