Actual source code: dspriv.c
slepc-3.7.3 2016-09-29
1: /*
2: Private DS routines
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, 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;
32: PetscInt n,d;
33: PetscBool ispep;
37: PetscObjectTypeCompare((PetscObject)ds,DSPEP,&ispep);
38: if (ispep) {
39: DSPEPGetDegree(ds,&d);
40: }
41: 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;
42: else n = ds->ld;
43: switch (m) {
44: case DS_MAT_T:
45: sz = 3*ds->ld*sizeof(PetscScalar);
46: break;
47: case DS_MAT_D:
48: sz = ds->ld*sizeof(PetscScalar);
49: break;
50: case DS_MAT_X:
51: sz = ds->ld*n*sizeof(PetscScalar);
52: break;
53: case DS_MAT_Y:
54: sz = ds->ld*n*sizeof(PetscScalar);
55: break;
56: default:
57: sz = n*n*sizeof(PetscScalar);
58: }
59: if (ds->mat[m]) {
60: PetscFree(ds->mat[m]);
61: } else {
62: PetscLogObjectMemory((PetscObject)ds,sz);
63: }
64: PetscMalloc(sz,&ds->mat[m]);
65: PetscMemzero(ds->mat[m],sz);
66: return(0);
67: }
71: PetscErrorCode DSAllocateMatReal_Private(DS ds,DSMatType m)
72: {
73: size_t sz;
77: if (m==DS_MAT_T) sz = 3*ds->ld*sizeof(PetscReal);
78: else if (m==DS_MAT_D) sz = ds->ld*sizeof(PetscReal);
79: else sz = ds->ld*ds->ld*sizeof(PetscReal);
80: if (!ds->rmat[m]) {
81: PetscLogObjectMemory((PetscObject)ds,sz);
82: PetscMalloc(sz,&ds->rmat[m]);
83: }
84: PetscMemzero(ds->rmat[m],sz);
85: return(0);
86: }
90: PetscErrorCode DSAllocateWork_Private(DS ds,PetscInt s,PetscInt r,PetscInt i)
91: {
95: if (s>ds->lwork) {
96: PetscFree(ds->work);
97: PetscMalloc1(s,&ds->work);
98: PetscLogObjectMemory((PetscObject)ds,(s-ds->lwork)*sizeof(PetscScalar));
99: ds->lwork = s;
100: }
101: if (r>ds->lrwork) {
102: PetscFree(ds->rwork);
103: PetscMalloc1(r,&ds->rwork);
104: PetscLogObjectMemory((PetscObject)ds,(r-ds->lrwork)*sizeof(PetscReal));
105: ds->lrwork = r;
106: }
107: if (i>ds->liwork) {
108: PetscFree(ds->iwork);
109: PetscMalloc1(i,&ds->iwork);
110: PetscLogObjectMemory((PetscObject)ds,(i-ds->liwork)*sizeof(PetscBLASInt));
111: ds->liwork = i;
112: }
113: return(0);
114: }
118: /*@C
119: DSViewMat - Prints one of the internal DS matrices.
121: Collective on DS
123: Input Parameters:
124: + ds - the direct solver context
125: . viewer - visualization context
126: - m - matrix to display
128: Note:
129: Works only for ascii viewers. Set the viewer in Matlab format if
130: want to paste into Matlab.
132: Level: developer
134: .seealso: DSView()
135: @*/
136: PetscErrorCode DSViewMat(DS ds,PetscViewer viewer,DSMatType m)
137: {
138: PetscErrorCode ierr;
139: PetscInt i,j,rows,cols,d;
140: PetscScalar *v;
141: PetscViewerFormat format;
142: PetscBool ispep;
143: #if defined(PETSC_USE_COMPLEX)
144: PetscBool allreal = PETSC_TRUE;
145: #endif
148: if (m>=DS_NUM_MAT) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid matrix");
149: if (!ds->mat[m]) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Requested matrix was not created in this DS");
150: PetscViewerGetFormat(viewer,&format);
151: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
152: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
153: if (ds->state==DS_STATE_TRUNCATED && m>=DS_MAT_Q) rows = ds->t;
154: else rows = (m==DS_MAT_A && ds->extrarow)? ds->n+1: ds->n;
155: cols = (ds->m!=0)? ds->m: ds->n;
156: PetscObjectTypeCompare((PetscObject)ds,DSPEP,&ispep);
157: if (ispep) {
158: DSPEPGetDegree(ds,&d);
159: }
160: if (ispep && (m==DS_MAT_X || m==DS_MAT_Y)) cols = d*ds->n;
161: #if defined(PETSC_USE_COMPLEX)
162: /* determine if matrix has all real values */
163: v = ds->mat[m];
164: for (i=0;i<rows;i++)
165: for (j=0;j<cols;j++)
166: if (PetscImaginaryPart(v[i+j*ds->ld])) { allreal = PETSC_FALSE; break; }
167: #endif
168: if (format == PETSC_VIEWER_ASCII_MATLAB) {
169: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",rows,cols);
170: PetscViewerASCIIPrintf(viewer,"%s = [\n",DSMatName[m]);
171: } else {
172: PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[m]);
173: }
175: for (i=0;i<rows;i++) {
176: v = ds->mat[m]+i;
177: for (j=0;j<cols;j++) {
178: #if defined(PETSC_USE_COMPLEX)
179: if (allreal) {
180: PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));
181: } else {
182: PetscViewerASCIIPrintf(viewer,"%18.16e%+18.16ei ",PetscRealPart(*v),PetscImaginaryPart(*v));
183: }
184: #else
185: PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);
186: #endif
187: v += ds->ld;
188: }
189: PetscViewerASCIIPrintf(viewer,"\n");
190: }
192: if (format == PETSC_VIEWER_ASCII_MATLAB) {
193: PetscViewerASCIIPrintf(viewer,"];\n");
194: }
195: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
196: PetscViewerFlush(viewer);
197: return(0);
198: }
202: PetscErrorCode DSSortEigenvalues_Private(DS ds,PetscScalar *wr,PetscScalar *wi,PetscInt *perm,PetscBool isghiep)
203: {
205: PetscScalar re,im,wi0;
206: PetscInt n,i,j,result,tmp1,tmp2=0,d=1;
209: n = ds->t; /* sort only first t pairs if truncated */
210: /* insertion sort */
211: i=ds->l+1;
212: #if !defined(PETSC_USE_COMPLEX)
213: if (wi && wi[perm[i-1]]!=0.0) i++; /* initial value is complex */
214: #else
215: if (isghiep && PetscImaginaryPart(wr[perm[i-1]])!=0.0) i++;
216: #endif
217: for (;i<n;i+=d) {
218: re = wr[perm[i]];
219: if (wi) im = wi[perm[i]];
220: else im = 0.0;
221: tmp1 = perm[i];
222: #if !defined(PETSC_USE_COMPLEX)
223: if (im!=0.0) { d = 2; tmp2 = perm[i+1]; }
224: else d = 1;
225: #else
226: if (isghiep && PetscImaginaryPart(re)!=0.0) { d = 2; tmp2 = perm[i+1]; }
227: else d = 1;
228: #endif
229: j = i-1;
230: if (wi) wi0 = wi[perm[j]];
231: else wi0 = 0.0;
232: SlepcSCCompare(ds->sc,re,im,wr[perm[j]],wi0,&result);
233: while (result<0 && j>=ds->l) {
234: perm[j+d] = perm[j];
235: j--;
236: #if !defined(PETSC_USE_COMPLEX)
237: if (wi && wi[perm[j+1]]!=0)
238: #else
239: if (isghiep && PetscImaginaryPart(wr[perm[j+1]])!=0)
240: #endif
241: { perm[j+d] = perm[j]; j--; }
243: if (j>=ds->l) {
244: if (wi) wi0 = wi[perm[j]];
245: else wi0 = 0.0;
246: SlepcSCCompare(ds->sc,re,im,wr[perm[j]],wi0,&result);
247: }
248: }
249: perm[j+1] = tmp1;
250: if (d==2) perm[j+2] = tmp2;
251: }
252: return(0);
253: }
257: PetscErrorCode DSSortEigenvaluesReal_Private(DS ds,PetscReal *eig,PetscInt *perm)
258: {
260: PetscScalar re;
261: PetscInt i,j,result,tmp,l,n;
264: n = ds->t; /* sort only first t pairs if truncated */
265: l = ds->l;
266: /* insertion sort */
267: for (i=l+1;i<n;i++) {
268: re = eig[perm[i]];
269: j = i-1;
270: SlepcSCCompare(ds->sc,re,0.0,eig[perm[j]],0.0,&result);
271: while (result<0 && j>=l) {
272: tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
273: if (j>=l) {
274: SlepcSCCompare(ds->sc,re,0.0,eig[perm[j]],0.0,&result);
275: }
276: }
277: }
278: return(0);
279: }
283: /*
284: DSCopyMatrix_Private - Copies the trailing block of a matrix (from
285: rows/columns l to n).
286: */
287: PetscErrorCode DSCopyMatrix_Private(DS ds,DSMatType dst,DSMatType src)
288: {
290: PetscInt j,m,off,ld;
291: PetscScalar *S,*D;
294: ld = ds->ld;
295: m = ds->n-ds->l;
296: off = ds->l+ds->l*ld;
297: S = ds->mat[src];
298: D = ds->mat[dst];
299: for (j=0;j<m;j++) {
300: PetscMemcpy(D+off+j*ld,S+off+j*ld,m*sizeof(PetscScalar));
301: }
302: return(0);
303: }
307: PetscErrorCode DSPermuteColumns_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat,PetscInt *perm)
308: {
309: PetscInt i,j,k,p,ld;
310: PetscScalar *Q,rtmp;
313: ld = ds->ld;
314: Q = ds->mat[mat];
315: for (i=l;i<n;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: }
325: }
326: }
327: return(0);
328: }
332: PetscErrorCode DSPermuteRows_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat,PetscInt *perm)
333: {
334: PetscInt i,j,m=ds->m,k,p,ld;
335: PetscScalar *Q,rtmp;
338: if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
339: ld = ds->ld;
340: Q = ds->mat[mat];
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 rows i and j */
348: for (k=0;k<m;k++) {
349: rtmp = Q[p+k*ld]; Q[p+k*ld] = Q[i+k*ld]; Q[i+k*ld] = rtmp;
350: }
351: }
352: }
353: return(0);
354: }
358: PetscErrorCode DSPermuteBoth_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat1,DSMatType mat2,PetscInt *perm)
359: {
360: PetscInt i,j,m=ds->m,k,p,ld;
361: PetscScalar *U,*VT,rtmp;
364: if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
365: ld = ds->ld;
366: U = ds->mat[mat1];
367: VT = ds->mat[mat2];
368: for (i=l;i<n;i++) {
369: p = perm[i];
370: if (p != i) {
371: j = i + 1;
372: while (perm[j] != i) j++;
373: perm[j] = p; perm[i] = i;
374: /* swap columns i and j of U */
375: for (k=0;k<n;k++) {
376: rtmp = U[k+p*ld]; U[k+p*ld] = U[k+i*ld]; U[k+i*ld] = rtmp;
377: }
378: /* swap rows i and j of VT */
379: for (k=0;k<m;k++) {
380: rtmp = VT[p+k*ld]; VT[p+k*ld] = VT[i+k*ld]; VT[i+k*ld] = rtmp;
381: }
382: }
383: }
384: return(0);
385: }
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;
408: DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
409: DSGetLeadingDimension(ds,&ld);
410: PetscLogEventBegin(DS_Other,ds,0,0,0);
411: DSGetArray(ds,mat,&x);
412: PetscMemzero(&x[ld*l],ld*(n-l)*sizeof(PetscScalar));
413: for (i=l;i<n;i++) {
414: x[ld*i+i] = 1.0;
415: }
416: DSRestoreArray(ds,mat,&x);
417: PetscLogEventEnd(DS_Other,ds,0,0,0);
418: return(0);
419: }
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: {
442: #if defined(PETSC_MISSING_LAPACK_GEQRF) || defined(PETSC_MISSING_LAPACK_ORGQR)
444: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQRF/ORGQR - Lapack routine is unavailable");
445: #else
447: PetscInt n,l,ld;
448: PetscBLASInt ld_,rA,cA,info,ltau,lw;
449: PetscScalar *A,*tau,*w,saux,dummy;
453: DSCheckAlloc(ds,1);
457: DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
458: DSGetLeadingDimension(ds,&ld);
459: n = n - l;
460: if (cols > n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns");
461: if (n == 0 || cols == 0) return(0);
463: PetscLogEventBegin(DS_Other,ds,0,0,0);
464: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
465: DSGetArray(ds,mat,&A);
466: PetscBLASIntCast(PetscMin(cols,n),<au);
467: PetscBLASIntCast(ld,&ld_);
468: PetscBLASIntCast(n,&rA);
469: PetscBLASIntCast(cols,&cA);
470: lw = -1;
471: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,A,&ld_,&dummy,&saux,&lw,&info));
472: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
473: lw = (PetscBLASInt)PetscRealPart(saux);
474: DSAllocateWork_Private(ds,lw+ltau,0,0);
475: tau = ds->work;
476: w = &tau[ltau];
477: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,&A[ld*l+l],&ld_,tau,w,&lw,&info));
478: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
479: PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&rA,<au,<au,&A[ld*l+l],&ld_,tau,w,&lw,&info));
480: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);
481: if (lindcols) *lindcols = ltau;
483: PetscFPTrapPop();
484: PetscLogEventEnd(DS_Other,ds,0,0,0);
485: DSRestoreArray(ds,mat,&A);
486: PetscObjectStateIncrease((PetscObject)ds);
487: return(0);
488: #endif
489: }
493: /*
494: Compute C <- a*A*B + b*C, where
495: ldC, the leading dimension of C,
496: ldA, the leading dimension of A,
497: rA, cA, rows and columns of A,
498: At, if true use the transpose of A instead,
499: ldB, the leading dimension of B,
500: rB, cB, rows and columns of B,
501: Bt, if true use the transpose of B instead
502: */
503: 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)
504: {
506: PetscInt tmp;
507: PetscBLASInt m, n, k, ldA = _ldA, ldB = _ldB, ldC = _ldC;
508: const char *N = "N", *T = "C", *qA = N, *qB = N;
511: if ((rA == 0) || (cB == 0)) return(0);
516: /* Transpose if needed */
517: if (At) tmp = rA, rA = cA, cA = tmp, qA = T;
518: if (Bt) tmp = rB, rB = cB, cB = tmp, qB = T;
520: /* Check size */
521: if (cA != rB) SETERRQ(PETSC_COMM_SELF,1, "Matrix dimensions do not match");
523: /* Do stub */
524: if ((rA == 1) && (cA == 1) && (cB == 1)) {
525: if (!At && !Bt) *C = *A * *B;
526: else if (At && !Bt) *C = PetscConj(*A) * *B;
527: else if (!At && Bt) *C = *A * PetscConj(*B);
528: else *C = PetscConj(*A) * PetscConj(*B);
529: m = n = k = 1;
530: } else {
531: m = rA; n = cB; k = cA;
532: PetscStackCallBLAS("BLASgemm",BLASgemm_(qA,qB,&m,&n,&k,&a,(PetscScalar*)A,&ldA,(PetscScalar*)B,&ldB,&b,C,&ldC));
533: }
535: PetscLogFlops(m*n*2*k);
536: return(0);
537: }
541: /*@C
542: DSPseudoOrthogonalize - Orthogonalize the columns of a matrix with Modified
543: Gram-Schmidt in an indefinite inner product space defined by a signature.
545: Logically Collective on DS
547: Input Parameters:
548: + ds - the direct solver context
549: . mat - the matrix
550: . cols - number of columns to orthogonalize (starting from column zero)
551: - s - the signature that defines the inner product
553: Output Parameters:
554: + lindcols - (optional) linearly independent columns of the matrix
555: - ns - (optional) the new norm of the vectors
557: Level: developer
559: .seealso: DSOrthogonalize()
560: @*/
561: PetscErrorCode DSPseudoOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscReal *s,PetscInt *lindcols,PetscReal *ns)
562: {
564: PetscInt i,j,k,l,n,ld;
565: PetscBLASInt one=1,rA_;
566: PetscScalar alpha,*A,*A_,*m,*h,nr0;
567: PetscReal nr_o,nr,*ns_;
571: DSCheckAlloc(ds,1);
575: DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
576: DSGetLeadingDimension(ds,&ld);
577: n = n - l;
578: if (cols > n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns");
579: if (n == 0 || cols == 0) return(0);
580: PetscBLASIntCast(n,&rA_);
581: DSGetArray(ds,mat,&A_);
582: A = &A_[ld*l+l];
583: DSAllocateWork_Private(ds,n+cols,ns?0:cols,0);
584: m = ds->work;
585: h = &m[n];
586: ns_ = ns ? ns : ds->rwork;
587: PetscLogEventBegin(DS_Other,ds,0,0,0);
588: for (i=0; i<cols; i++) {
589: /* m <- diag(s)*A[i] */
590: for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
591: /* nr_o <- mynorm(A[i]'*m), mynorm(x) = sign(x)*sqrt(|x|) */
592: SlepcMatDenseMult(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
593: nr = nr_o = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
594: for (j=0; j<3 && i>0; j++) {
595: /* h <- A[0:i-1]'*m */
596: SlepcMatDenseMult(h,i,0.0,1.0,A,ld,n,i,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
597: /* h <- diag(ns)*h */
598: for (k=0; k<i; k++) h[k] *= ns_[k];
599: /* A[i] <- A[i] - A[0:i-1]*h */
600: SlepcMatDenseMult(&A[ld*i],ld,1.0,-1.0,A,ld,n,i,PETSC_FALSE,h,i,i,1,PETSC_FALSE);
601: /* m <- diag(s)*A[i] */
602: for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
603: /* nr_o <- mynorm(A[i]'*m) */
604: SlepcMatDenseMult(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
605: nr = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
606: if (PetscAbs(nr) < PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_SELF,1,"Linear dependency detected");
607: if (PetscAbs(nr) > 0.7*PetscAbs(nr_o)) break;
608: nr_o = nr;
609: }
610: ns_[i] = PetscSign(nr);
611: /* A[i] <- A[i]/|nr| */
612: alpha = 1.0/PetscAbs(nr);
613: PetscStackCallBLAS("BLASscal",BLASscal_(&rA_,&alpha,&A[i*ld],&one));
614: }
615: PetscLogEventEnd(DS_Other,ds,0,0,0);
616: DSRestoreArray(ds,mat,&A_);
617: PetscObjectStateIncrease((PetscObject)ds);
618: if (lindcols) *lindcols = cols;
619: return(0);
620: }