Actual source code: dspriv.c

slepc-3.10.0 2018-09-18
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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>      /*I "slepcds.h" I*/
 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;

 25:   PetscObjectTypeCompare((PetscObject)ds,DSPEP,&ispep);
 26:   if (ispep) {
 27:     DSPEPGetDegree(ds,&d);
 28:   }
 29:   if (ispep && (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;
 30:   else n = ds->ld;
 31:   switch (m) {
 32:     case DS_MAT_T:
 33:       nelem = 3*ds->ld;
 34:       break;
 35:     case DS_MAT_D:
 36:       nelem = ds->ld;
 37:       break;
 38:     case DS_MAT_X:
 39:       nelem = ds->ld*n;
 40:       break;
 41:     case DS_MAT_Y:
 42:       nelem = ds->ld*n;
 43:       break;
 44:     default:
 45:       nelem = n*n;
 46:   }
 47:   if (isreal) {
 48:     sz = nelem*sizeof(PetscReal);
 49:     if (ds->rmat[m]) {
 50:       PetscFree(ds->rmat[m]);
 51:     } else {
 52:       PetscLogObjectMemory((PetscObject)ds,sz);
 53:     }
 54:     PetscCalloc1(nelem,&ds->rmat[m]);
 55:   } else {
 56:     sz = nelem*sizeof(PetscScalar);
 57:     if (ds->mat[m]) {
 58:       PetscFree(ds->mat[m]);
 59:     } else {
 60:       PetscLogObjectMemory((PetscObject)ds,sz);
 61:     }
 62:     PetscCalloc1(nelem,&ds->mat[m]);
 63:   }
 64:   return(0);
 65: }

 67: PetscErrorCode DSAllocateWork_Private(DS ds,PetscInt s,PetscInt r,PetscInt i)
 68: {

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

 93: /*@C
 94:    DSViewMat - Prints one of the internal DS matrices.

 96:    Collective on DS

 98:    Input Parameters:
 99: +  ds     - the direct solver context
100: .  viewer - visualization context
101: -  m      - matrix to display

103:    Note:
104:    Works only for ascii viewers. Set the viewer in Matlab format if
105:    want to paste into Matlab.

107:    Level: developer

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

124:   DSCheckValidMat(ds,m,3);
125:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ds));
128:   PetscViewerGetFormat(viewer,&format);
129:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
130:   PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
131:   DSMatGetSize(ds,m,&rows,&cols);
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 ",(double)PetscRealPart(*v));
152:       } else {
153:         PetscViewerASCIIPrintf(viewer,"%18.16e%+18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));
154:       }
155: #else
156:       PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*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: }

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

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

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

224: PetscErrorCode DSSortEigenvaluesReal_Private(DS ds,PetscReal *eig,PetscInt *perm)
225: {
227:   PetscScalar    re;
228:   PetscInt       i,j,result,tmp,l,n;

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

248: /*
249:   DSCopyMatrix_Private - Copies the trailing block of a matrix (from
250:   rows/columns l to n).
251: */
252: PetscErrorCode DSCopyMatrix_Private(DS ds,DSMatType dst,DSMatType src)
253: {
255:   PetscInt    j,m,off,ld;
256:   PetscScalar *S,*D;

259:   ld  = ds->ld;
260:   m   = ds->n-ds->l;
261:   off = ds->l+ds->l*ld;
262:   S   = ds->mat[src];
263:   D   = ds->mat[dst];
264:   for (j=0;j<m;j++) {
265:     PetscMemcpy(D+off+j*ld,S+off+j*ld,m*sizeof(PetscScalar));
266:   }
267:   return(0);
268: }

270: PetscErrorCode DSPermuteColumns_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat,PetscInt *perm)
271: {
272:   PetscInt    i,j,k,p,ld;
273:   PetscScalar *Q,rtmp;

276:   ld = ds->ld;
277:   Q  = ds->mat[mat];
278:   for (i=l;i<n;i++) {
279:     p = perm[i];
280:     if (p != i) {
281:       j = i + 1;
282:       while (perm[j] != i) j++;
283:       perm[j] = p; perm[i] = i;
284:       /* swap columns i and j */
285:       for (k=0;k<n;k++) {
286:         rtmp = Q[k+p*ld]; Q[k+p*ld] = Q[k+i*ld]; Q[k+i*ld] = rtmp;
287:       }
288:     }
289:   }
290:   return(0);
291: }

293: PetscErrorCode DSPermuteRows_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat,PetscInt *perm)
294: {
295:   PetscInt    i,j,m=ds->m,k,p,ld;
296:   PetscScalar *Q,rtmp;

299:   if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
300:   ld = ds->ld;
301:   Q  = ds->mat[mat];
302:   for (i=l;i<n;i++) {
303:     p = perm[i];
304:     if (p != i) {
305:       j = i + 1;
306:       while (perm[j] != i) j++;
307:       perm[j] = p; perm[i] = i;
308:       /* swap rows i and j */
309:       for (k=0;k<m;k++) {
310:         rtmp = Q[p+k*ld]; Q[p+k*ld] = Q[i+k*ld]; Q[i+k*ld] = rtmp;
311:       }
312:     }
313:   }
314:   return(0);
315: }

317: PetscErrorCode DSPermuteBoth_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat1,DSMatType mat2,PetscInt *perm)
318: {
319:   PetscInt    i,j,m=ds->m,k,p,ld;
320:   PetscScalar *U,*VT,rtmp;

323:   if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
324:   ld = ds->ld;
325:   U  = ds->mat[mat1];
326:   VT = ds->mat[mat2];
327:   for (i=l;i<n;i++) {
328:     p = perm[i];
329:     if (p != i) {
330:       j = i + 1;
331:       while (perm[j] != i) j++;
332:       perm[j] = p; perm[i] = i;
333:       /* swap columns i and j of U */
334:       for (k=0;k<n;k++) {
335:         rtmp = U[k+p*ld]; U[k+p*ld] = U[k+i*ld]; U[k+i*ld] = rtmp;
336:       }
337:       /* swap rows i and j of VT */
338:       for (k=0;k<m;k++) {
339:         rtmp = VT[p+k*ld]; VT[p+k*ld] = VT[i+k*ld]; VT[i+k*ld] = rtmp;
340:       }
341:     }
342:   }
343:   return(0);
344: }

346: /*@
347:    DSSetIdentity - Copy the identity (a diagonal matrix with ones) on the
348:    active part of a matrix.

350:    Logically Collective on DS

352:    Input Parameters:
353: +  ds  - the direct solver context
354: -  mat - the matrix to modify

356:    Level: intermediate
357: @*/
358: PetscErrorCode DSSetIdentity(DS ds,DSMatType mat)
359: {
361:   PetscScalar    *x;
362:   PetscInt       i,ld,n,l;

367:   DSCheckValidMat(ds,mat,2);

369:   DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
370:   DSGetLeadingDimension(ds,&ld);
371:   PetscLogEventBegin(DS_Other,ds,0,0,0);
372:   DSGetArray(ds,mat,&x);
373:   PetscMemzero(&x[ld*l],ld*(n-l)*sizeof(PetscScalar));
374:   for (i=l;i<n;i++) x[i+i*ld] = 1.0;
375:   DSRestoreArray(ds,mat,&x);
376:   PetscLogEventEnd(DS_Other,ds,0,0,0);
377:   return(0);
378: }

380: /*@C
381:    DSOrthogonalize - Orthogonalize the columns of a matrix.

383:    Logically Collective on DS

385:    Input Parameters:
386: +  ds   - the direct solver context
387: .  mat  - a matrix
388: -  cols - number of columns to orthogonalize (starting from column zero)

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

393:    Level: developer

395: .seealso: DSPseudoOrthogonalize()
396: @*/
397: PetscErrorCode DSOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscInt *lindcols)
398: {
399: #if defined(PETSC_MISSING_LAPACK_GEQRF) || defined(PETSC_MISSING_LAPACK_ORGQR)
401:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQRF/ORGQR - Lapack routine is unavailable");
402: #else
404:   PetscInt       n,l,ld;
405:   PetscBLASInt   ld_,rA,cA,info,ltau,lw;
406:   PetscScalar    *A,*tau,*w,saux,dummy;

410:   DSCheckAlloc(ds,1);
412:   DSCheckValidMat(ds,mat,2);

415:   DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
416:   DSGetLeadingDimension(ds,&ld);
417:   n = n - l;
418:   if (cols > n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns");
419:   if (n == 0 || cols == 0) return(0);

421:   PetscLogEventBegin(DS_Other,ds,0,0,0);
422:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
423:   DSGetArray(ds,mat,&A);
424:   PetscBLASIntCast(PetscMin(cols,n),&ltau);
425:   PetscBLASIntCast(ld,&ld_);
426:   PetscBLASIntCast(n,&rA);
427:   PetscBLASIntCast(cols,&cA);
428:   lw = -1;
429:   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,A,&ld_,&dummy,&saux,&lw,&info));
430:   SlepcCheckLapackInfo("geqrf",info);
431:   lw = (PetscBLASInt)PetscRealPart(saux);
432:   DSAllocateWork_Private(ds,lw+ltau,0,0);
433:   tau = ds->work;
434:   w = &tau[ltau];
435:   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,&A[ld*l+l],&ld_,tau,w,&lw,&info));
436:   SlepcCheckLapackInfo("geqrf",info);
437:   PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&rA,&ltau,&ltau,&A[ld*l+l],&ld_,tau,w,&lw,&info));
438:   SlepcCheckLapackInfo("orgqr",info);
439:   if (lindcols) *lindcols = ltau;

441:   PetscFPTrapPop();
442:   PetscLogEventEnd(DS_Other,ds,0,0,0);
443:   DSRestoreArray(ds,mat,&A);
444:   PetscObjectStateIncrease((PetscObject)ds);
445:   return(0);
446: #endif
447: }

449: /*
450:   Compute C <- a*A*B + b*C, where
451:     ldC, the leading dimension of C,
452:     ldA, the leading dimension of A,
453:     rA, cA, rows and columns of A,
454:     At, if true use the transpose of A instead,
455:     ldB, the leading dimension of B,
456:     rB, cB, rows and columns of B,
457:     Bt, if true use the transpose of B instead
458: */
459: 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)
460: {
462:   PetscInt       tmp;
463:   PetscBLASInt   m, n, k, ldA = _ldA, ldB = _ldB, ldC = _ldC;
464:   const char     *N = "N", *T = "C", *qA = N, *qB = N;

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

472:   /* Transpose if needed */
473:   if (At) tmp = rA, rA = cA, cA = tmp, qA = T;
474:   if (Bt) tmp = rB, rB = cB, cB = tmp, qB = T;

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

479:   /* Do stub */
480:   if ((rA == 1) && (cA == 1) && (cB == 1)) {
481:     if (!At && !Bt) *C = *A * *B;
482:     else if (At && !Bt) *C = PetscConj(*A) * *B;
483:     else if (!At && Bt) *C = *A * PetscConj(*B);
484:     else *C = PetscConj(*A) * PetscConj(*B);
485:     m = n = k = 1;
486:   } else {
487:     m = rA; n = cB; k = cA;
488:     PetscStackCallBLAS("BLASgemm",BLASgemm_(qA,qB,&m,&n,&k,&a,(PetscScalar*)A,&ldA,(PetscScalar*)B,&ldB,&b,C,&ldC));
489:   }

491:   PetscLogFlops(2.0*m*n*k);
492:   return(0);
493: }

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

499:    Logically Collective on DS

501:    Input Parameters:
502: +  ds   - the direct solver context
503: .  mat  - the matrix
504: .  cols - number of columns to orthogonalize (starting from column zero)
505: -  s    - the signature that defines the inner product

507:    Output Parameters:
508: +  lindcols - (optional) linearly independent columns of the matrix
509: -  ns   - (optional) the new signature of the vectors

511:    Note:
512:    After the call the matrix satisfies A'*s*A = ns.

514:    Level: developer

516: .seealso: DSOrthogonalize()
517: @*/
518: PetscErrorCode DSPseudoOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscReal *s,PetscInt *lindcols,PetscReal *ns)
519: {
521:   PetscInt       i,j,k,l,n,ld;
522:   PetscBLASInt   one=1,rA_;
523:   PetscScalar    alpha,*A,*A_,*m,*h,nr0;
524:   PetscReal      nr_o,nr,*ns_;

528:   DSCheckAlloc(ds,1);
530:   DSCheckValidMat(ds,mat,2);
533:   DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
534:   DSGetLeadingDimension(ds,&ld);
535:   n = n - l;
536:   if (cols > n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns");
537:   if (n == 0 || cols == 0) return(0);
538:   PetscBLASIntCast(n,&rA_);
539:   DSGetArray(ds,mat,&A_);
540:   A = &A_[ld*l+l];
541:   DSAllocateWork_Private(ds,n+cols,ns?0:cols,0);
542:   m = ds->work;
543:   h = &m[n];
544:   ns_ = ns ? ns : ds->rwork;
545:   PetscLogEventBegin(DS_Other,ds,0,0,0);
546:   for (i=0; i<cols; i++) {
547:     /* m <- diag(s)*A[i] */
548:     for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
549:     /* nr_o <- mynorm(A[i]'*m), mynorm(x) = sign(x)*sqrt(|x|) */
550:     SlepcMatDenseMult(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
551:     nr = nr_o = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
552:     for (j=0; j<3 && i>0; j++) {
553:       /* h <- A[0:i-1]'*m */
554:       SlepcMatDenseMult(h,i,0.0,1.0,A,ld,n,i,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
555:       /* h <- diag(ns)*h */
556:       for (k=0; k<i; k++) h[k] *= ns_[k];
557:       /* A[i] <- A[i] - A[0:i-1]*h */
558:       SlepcMatDenseMult(&A[ld*i],ld,1.0,-1.0,A,ld,n,i,PETSC_FALSE,h,i,i,1,PETSC_FALSE);
559:       /* m <- diag(s)*A[i] */
560:       for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
561:       /* nr_o <- mynorm(A[i]'*m) */
562:       SlepcMatDenseMult(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
563:       nr = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
564:       if (PetscAbs(nr) < PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_SELF,1,"Linear dependency detected");
565:       if (PetscAbs(nr) > 0.7*PetscAbs(nr_o)) break;
566:       nr_o = nr;
567:     }
568:     ns_[i] = PetscSign(nr);
569:     /* A[i] <- A[i]/|nr| */
570:     alpha = 1.0/PetscAbs(nr);
571:     PetscStackCallBLAS("BLASscal",BLASscal_(&rA_,&alpha,&A[i*ld],&one));
572:   }
573:   PetscLogEventEnd(DS_Other,ds,0,0,0);
574:   DSRestoreArray(ds,mat,&A_);
575:   PetscObjectStateIncrease((PetscObject)ds);
576:   if (lindcols) *lindcols = cols;
577:   return(0);
578: }