Actual source code: dspriv.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    Private DS routines

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <slepc-private/dsimpl.h>      /*I "slepcds.h" I*/
 25: #include <slepcblaslapack.h>

 29: PetscErrorCode DSAllocateMat_Private(DS ds,DSMatType m)
 30: {
 31:   size_t         sz;

 35:   if (m==DS_MAT_T) sz = 3*ds->ld*sizeof(PetscScalar);
 36:   else if (m==DS_MAT_D) sz = ds->ld*sizeof(PetscScalar);
 37:   else sz = ds->ld*ds->ld*sizeof(PetscScalar);
 38:   if (ds->mat[m]) {
 39:     PetscFree(ds->mat[m]);
 40:   } else {
 41:     PetscLogObjectMemory((PetscObject)ds,sz);
 42:   }
 43:   PetscMalloc(sz,&ds->mat[m]);
 44:   PetscMemzero(ds->mat[m],sz);
 45:   return(0);
 46: }

 50: PetscErrorCode DSAllocateMatReal_Private(DS ds,DSMatType m)
 51: {
 52:   size_t         sz;

 56:   if (m==DS_MAT_T) sz = 3*ds->ld*sizeof(PetscReal);
 57:   else if (m==DS_MAT_D) sz = ds->ld*sizeof(PetscReal);
 58:   else sz = ds->ld*ds->ld*sizeof(PetscReal);
 59:   if (!ds->rmat[m]) {
 60:     PetscLogObjectMemory((PetscObject)ds,sz);
 61:     PetscMalloc(sz,&ds->rmat[m]);
 62:   }
 63:   PetscMemzero(ds->rmat[m],sz);
 64:   return(0);
 65: }

 69: PetscErrorCode DSAllocateWork_Private(DS ds,PetscInt s,PetscInt r,PetscInt i)
 70: {

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

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

100:    Collective on DS

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

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

111:    Level: developer

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

126:   PetscViewerGetFormat(viewer,&format);
127:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
128:   PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
129:   if (ds->state==DS_STATE_TRUNCATED && m>=DS_MAT_Q) rows = ds->t;
130:   else rows = (m==DS_MAT_A && ds->extrarow)? ds->n+1: ds->n;
131:   cols = (ds->m!=0)? ds->m: ds->n;
132: #if defined(PETSC_USE_COMPLEX)
133:   /* determine if matrix has all real values */
134:   v = ds->mat[m];
135:   for (i=0;i<rows;i++)
136:     for (j=0;j<cols;j++)
137:       if (PetscImaginaryPart(v[i+j*ds->ld])) { allreal = PETSC_FALSE; break; }
138: #endif
139:   if (format == PETSC_VIEWER_ASCII_MATLAB) {
140:     PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",rows,cols);
141:     PetscViewerASCIIPrintf(viewer,"%s = [\n",DSMatName[m]);
142:   } else {
143:     PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[m]);
144:   }

146:   for (i=0;i<rows;i++) {
147:     v = ds->mat[m]+i;
148:     for (j=0;j<cols;j++) {
149: #if defined(PETSC_USE_COMPLEX)
150:       if (allreal) {
151:         PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));
152:       } else {
153:         PetscViewerASCIIPrintf(viewer,"%18.16e%+18.16ei ",PetscRealPart(*v),PetscImaginaryPart(*v));
154:       }
155: #else
156:       PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);
157: #endif
158:       v += ds->ld;
159:     }
160:     PetscViewerASCIIPrintf(viewer,"\n");
161:   }

163:   if (format == PETSC_VIEWER_ASCII_MATLAB) {
164:     PetscViewerASCIIPrintf(viewer,"];\n");
165:   }
166:   PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
167:   PetscViewerFlush(viewer);
168:   return(0);
169: }

173: PetscErrorCode DSSortEigenvalues_Private(DS ds,PetscScalar *wr,PetscScalar *wi,PetscInt *perm,PetscBool isghiep)
174: {
176:   PetscScalar    re,im,wi0;
177:   PetscInt       n,i,j,result,tmp1,tmp2=0,d=1;

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

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

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

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

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

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

280: PetscErrorCode DSPermuteColumns_Private(DS ds,PetscInt l,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=l;i<n;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: }

305: PetscErrorCode DSPermuteRows_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat,PetscInt *perm)
306: {
307:   PetscInt    i,j,m=ds->m,k,p,ld;
308:   PetscScalar *Q,rtmp;

311:   if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
312:   ld = ds->ld;
313:   Q  = ds->mat[mat];
314:   for (i=l;i<n;i++) {
315:     p = perm[i];
316:     if (p != i) {
317:       j = i + 1;
318:       while (perm[j] != i) j++;
319:       perm[j] = p; perm[i] = i;
320:       /* swap rows i and j */
321:       for (k=0;k<m;k++) {
322:         rtmp = Q[p+k*ld]; Q[p+k*ld] = Q[i+k*ld]; Q[i+k*ld] = rtmp;
323:       }
324:     }
325:   }
326:   return(0);
327: }

331: PetscErrorCode DSPermuteBoth_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat1,DSMatType mat2,PetscInt *perm)
332: {
333:   PetscInt    i,j,m=ds->m,k,p,ld;
334:   PetscScalar *U,*VT,rtmp;

337:   if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
338:   ld = ds->ld;
339:   U  = ds->mat[mat1];
340:   VT = ds->mat[mat2];
341:   for (i=l;i<n;i++) {
342:     p = perm[i];
343:     if (p != i) {
344:       j = i + 1;
345:       while (perm[j] != i) j++;
346:       perm[j] = p; perm[i] = i;
347:       /* swap columns i and j of U */
348:       for (k=0;k<n;k++) {
349:         rtmp = U[k+p*ld]; U[k+p*ld] = U[k+i*ld]; U[k+i*ld] = rtmp;
350:       }
351:       /* swap rows i and j of VT */
352:       for (k=0;k<m;k++) {
353:         rtmp = VT[p+k*ld]; VT[p+k*ld] = VT[i+k*ld]; VT[i+k*ld] = rtmp;
354:       }
355:     }
356:   }
357:   return(0);
358: }

362: /*
363:    DSSetIdentity - Copy the identity (a diagonal matrix with ones) on the
364:    active part of a matrix.
365: */
366: PetscErrorCode DSSetIdentity(DS ds,DSMatType mat)
367: {
369:   PetscScalar    *x;
370:   PetscInt       i,ld,n,l;

373:   DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
374:   DSGetLeadingDimension(ds,&ld);
375:   PetscLogEventBegin(DS_Other,ds,0,0,0);
376:   DSGetArray(ds,mat,&x);
377:   PetscMemzero(&x[ld*l],ld*(n-l)*sizeof(PetscScalar));
378:   for (i=l;i<n;i++) {
379:     x[ld*i+i] = 1.0;
380:   }
381:   DSRestoreArray(ds,mat,&x);
382:   PetscLogEventEnd(DS_Other,ds,0,0,0);
383:   return(0);
384: }

388: /*
389:    DSComputeMatrix - Build the matrix associated to a nonlinear operator
390:    T(lambda) or its derivative T'(lambda), given the parameter lambda, where
391:    T(lambda) = sum_i E_i*f_i(lambda). The result is written in mat.
392: */
393: PetscErrorCode DSComputeMatrix(DS ds,PetscScalar lambda,PetscBool deriv,DSMatType mat)
394: {
396:   PetscScalar    *T,*E,alpha;
397:   PetscInt       i,ld,n;
398:   PetscBLASInt   k,inc=1;

401:   DSGetDimensions(ds,&n,NULL,NULL,NULL,NULL);
402:   DSGetLeadingDimension(ds,&ld);
403:   PetscBLASIntCast(ld*n,&k);
404:   PetscLogEventBegin(DS_Other,ds,0,0,0);
405:   DSGetArray(ds,mat,&T);
406:   PetscMemzero(T,k*sizeof(PetscScalar));
407:   for (i=0;i<ds->nf;i++) {
408:     if (deriv) {
409:       FNEvaluateDerivative(ds->f[i],lambda,&alpha);
410:     } else {
411:       FNEvaluateFunction(ds->f[i],lambda,&alpha);
412:     }
413:     E = ds->mat[DSMatExtra[i]];
414:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k,&alpha,E,&inc,T,&inc));
415:   }
416:   DSRestoreArray(ds,mat,&T);
417:   PetscLogEventEnd(DS_Other,ds,0,0,0);
418:   return(0);
419: }

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

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

431:    Output Parameter:
432: .  lindcols - number of linearly independent columns of the matrix (can be NULL)
433: */
434: PetscErrorCode DSOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscInt *lindcols)
435: {
436: #if defined(PETSC_MISSING_LAPACK_GEQRF) || defined(SLEPC_MISSING_LAPACK_ORGQR)
438:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQRF/ORGQR - Lapack routine is unavailable");
439: #else
440:   PetscErrorCode  ierr;
441:   PetscInt        n,l,ld;
442:   PetscBLASInt    ld_,rA,cA,info,ltau,lw;
443:   PetscScalar     *A,*tau,*w,saux;

449:   DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
450:   DSGetLeadingDimension(ds,&ld);
451:   n = n - l;
452:   if (cols > n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns");
453:   if (n == 0 || cols == 0) return(0);
454:   DSGetArray(ds,mat,&A);
455:   PetscBLASIntCast(PetscMin(cols,n),&ltau);
456:   PetscBLASIntCast(ld,&ld_);
457:   PetscBLASIntCast(n,&rA);
458:   PetscBLASIntCast(cols,&cA);
459:   lw = -1;
460:   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,A,&ld_,NULL,&saux,&lw,&info));
461:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
462:   lw = (PetscBLASInt)PetscRealPart(saux);
463:   DSAllocateWork_Private(ds,lw+ltau,0,0);
464:   tau = ds->work;
465:   w = &tau[ltau];
466:   PetscLogEventBegin(DS_Other,ds,0,0,0);
467:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
468:   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,&A[ld*l+l],&ld_,tau,w,&lw,&info));
469:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
470:   PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&rA,&ltau,&ltau,&A[ld*l+l],&ld_,tau,w,&lw,&info));
471:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);
472:   PetscFPTrapPop();
473:   PetscLogEventEnd(DS_Other,ds,0,0,0);
474:   DSRestoreArray(ds,mat,&A);
475:   PetscObjectStateIncrease((PetscObject)ds);
476:   if (lindcols) *lindcols = ltau;
477:   return(0);
478: #endif
479: }

483: /*
484:   Compute C <- a*A*B + b*C, where
485:     ldC, the leading dimension of C,
486:     ldA, the leading dimension of A,
487:     rA, cA, rows and columns of A,
488:     At, if true use the transpose of A instead,
489:     ldB, the leading dimension of B,
490:     rB, cB, rows and columns of B,
491:     Bt, if true use the transpose of B instead
492: */
493: static PetscErrorCode SlepcDenseMatProd(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)
494: {
496:   PetscInt       tmp;
497:   PetscBLASInt   m, n, k, ldA = _ldA, ldB = _ldB, ldC = _ldC;
498:   const char     *N = "N", *T = "C", *qA = N, *qB = N;

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

506:   /* Transpose if needed */
507:   if (At) tmp = rA, rA = cA, cA = tmp, qA = T;
508:   if (Bt) tmp = rB, rB = cB, cB = tmp, qB = T;

510:   /* Check size */
511:   if (cA != rB) SETERRQ(PETSC_COMM_SELF,1, "Matrix dimensions do not match");

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

525:   PetscLogFlops(m*n*2*k);
526:   return(0);
527: }

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

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

541:    Output Parameter:
542: +  lindcols - linear independent columns of the matrix (can be NULL)
543: -  ns - the new norm of the vectors (can be NULL)
544: */
545: PetscErrorCode DSPseudoOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscReal *s,PetscInt *lindcols,PetscReal *ns)
546: {
547:   PetscErrorCode  ierr;
548:   PetscInt        i,j,k,l,n,ld;
549:   PetscBLASInt    one=1,rA_;
550:   PetscScalar     alpha,*A,*A_,*m,*h,nr0;
551:   PetscReal       nr_o,nr,*ns_;

559:   DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
560:   DSGetLeadingDimension(ds,&ld);
561:   n = n - l;
562:   if (cols > n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns");
563:   if (n == 0 || cols == 0) return(0);
564:   PetscBLASIntCast(n,&rA_);
565:   DSGetArray(ds,mat,&A_);
566:   A = &A_[ld*l+l];
567:   DSAllocateWork_Private(ds,n+cols,ns?0:cols,0);
568:   m = ds->work;
569:   h = &m[n];
570:   ns_ = ns ? ns : ds->rwork;
571:   PetscLogEventBegin(DS_Other,ds,0,0,0);
572:   for (i=0; i<cols; i++) {
573:     /* m <- diag(s)*A[i] */
574:     for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
575:     /* nr_o <- mynorm(A[i]'*m), mynorm(x) = sign(x)*sqrt(|x|) */
576:     SlepcDenseMatProd(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
577:     nr = nr_o = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
578:     for (j=0; j<3 && i>0; j++) {
579:       /* h <- A[0:i-1]'*m */
580:       SlepcDenseMatProd(h,i,0.0,1.0,A,ld,n,i,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
581:       /* h <- diag(ns)*h */
582:       for (k=0; k<i; k++) h[k] *= ns_[k];
583:       /* A[i] <- A[i] - A[0:i-1]*h */
584:       SlepcDenseMatProd(&A[ld*i],ld,1.0,-1.0,A,ld,n,i,PETSC_FALSE,h,i,i,1,PETSC_FALSE);
585:       /* m <- diag(s)*A[i] */
586:       for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
587:       /* nr_o <- mynorm(A[i]'*m) */
588:       SlepcDenseMatProd(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
589:       nr = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
590:       if (PetscAbs(nr) < PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_SELF,1,"Linear dependency detected");
591:       if (PetscAbs(nr) > 0.7*PetscAbs(nr_o)) break;
592:       nr_o = nr;
593:     }
594:     ns_[i] = PetscSign(nr);
595:     /* A[i] <- A[i]/|nr| */
596:     alpha = 1.0/PetscAbs(nr);
597:     PetscStackCallBLAS("BLASscal",BLASscal_(&rA_,&alpha,&A[i*ld],&one));
598:   }
599:   PetscLogEventEnd(DS_Other,ds,0,0,0);
600:   DSRestoreArray(ds,mat,&A_);
601:   PetscObjectStateIncrease((PetscObject)ds);
602:   if (lindcols) *lindcols = cols;
603:   return(0);
604: }