Actual source code: dspriv.c
slepc-3.6.1 2015-09-03
1: /*
2: Private DS routines
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, 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_X)) 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: default:
54: sz = n*n*sizeof(PetscScalar);
55: }
56: if (ds->mat[m]) {
57: PetscFree(ds->mat[m]);
58: } else {
59: PetscLogObjectMemory((PetscObject)ds,sz);
60: }
61: PetscMalloc(sz,&ds->mat[m]);
62: PetscMemzero(ds->mat[m],sz);
63: return(0);
64: }
68: PetscErrorCode DSAllocateMatReal_Private(DS ds,DSMatType m)
69: {
70: size_t sz;
74: if (m==DS_MAT_T) sz = 3*ds->ld*sizeof(PetscReal);
75: else if (m==DS_MAT_D) sz = ds->ld*sizeof(PetscReal);
76: else sz = ds->ld*ds->ld*sizeof(PetscReal);
77: if (!ds->rmat[m]) {
78: PetscLogObjectMemory((PetscObject)ds,sz);
79: PetscMalloc(sz,&ds->rmat[m]);
80: }
81: PetscMemzero(ds->rmat[m],sz);
82: return(0);
83: }
87: PetscErrorCode DSAllocateWork_Private(DS ds,PetscInt s,PetscInt r,PetscInt i)
88: {
92: if (s>ds->lwork) {
93: PetscFree(ds->work);
94: PetscMalloc1(s,&ds->work);
95: PetscLogObjectMemory((PetscObject)ds,(s-ds->lwork)*sizeof(PetscScalar));
96: ds->lwork = s;
97: }
98: if (r>ds->lrwork) {
99: PetscFree(ds->rwork);
100: PetscMalloc1(r,&ds->rwork);
101: PetscLogObjectMemory((PetscObject)ds,(r-ds->lrwork)*sizeof(PetscReal));
102: ds->lrwork = r;
103: }
104: if (i>ds->liwork) {
105: PetscFree(ds->iwork);
106: PetscMalloc1(i,&ds->iwork);
107: PetscLogObjectMemory((PetscObject)ds,(i-ds->liwork)*sizeof(PetscBLASInt));
108: ds->liwork = i;
109: }
110: return(0);
111: }
115: /*@C
116: DSViewMat - Prints one of the internal DS matrices.
118: Collective on DS
120: Input Parameters:
121: + ds - the direct solver context
122: . viewer - visualization context
123: - m - matrix to display
125: Note:
126: Works only for ascii viewers. Set the viewer in Matlab format if
127: want to paste into Matlab.
129: Level: developer
131: .seealso: DSView()
132: @*/
133: PetscErrorCode DSViewMat(DS ds,PetscViewer viewer,DSMatType m)
134: {
135: PetscErrorCode ierr;
136: PetscInt i,j,rows,cols;
137: PetscScalar *v;
138: PetscViewerFormat format;
139: #if defined(PETSC_USE_COMPLEX)
140: PetscBool allreal = PETSC_TRUE;
141: #endif
144: if (m>=DS_NUM_MAT) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid matrix");
145: if (!ds->mat[m]) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Requested matrix was not created in this DS");
146: PetscViewerGetFormat(viewer,&format);
147: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
148: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
149: if (ds->state==DS_STATE_TRUNCATED && m>=DS_MAT_Q) rows = ds->t;
150: else rows = (m==DS_MAT_A && ds->extrarow)? ds->n+1: ds->n;
151: cols = (ds->m!=0)? ds->m: ds->n;
152: #if defined(PETSC_USE_COMPLEX)
153: /* determine if matrix has all real values */
154: v = ds->mat[m];
155: for (i=0;i<rows;i++)
156: for (j=0;j<cols;j++)
157: if (PetscImaginaryPart(v[i+j*ds->ld])) { allreal = PETSC_FALSE; break; }
158: #endif
159: if (format == PETSC_VIEWER_ASCII_MATLAB) {
160: PetscViewerASCIIPrintf(viewer,"%% Size = %D %D\n",rows,cols);
161: PetscViewerASCIIPrintf(viewer,"%s = [\n",DSMatName[m]);
162: } else {
163: PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[m]);
164: }
166: for (i=0;i<rows;i++) {
167: v = ds->mat[m]+i;
168: for (j=0;j<cols;j++) {
169: #if defined(PETSC_USE_COMPLEX)
170: if (allreal) {
171: PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));
172: } else {
173: PetscViewerASCIIPrintf(viewer,"%18.16e%+18.16ei ",PetscRealPart(*v),PetscImaginaryPart(*v));
174: }
175: #else
176: PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);
177: #endif
178: v += ds->ld;
179: }
180: PetscViewerASCIIPrintf(viewer,"\n");
181: }
183: if (format == PETSC_VIEWER_ASCII_MATLAB) {
184: PetscViewerASCIIPrintf(viewer,"];\n");
185: }
186: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
187: PetscViewerFlush(viewer);
188: return(0);
189: }
193: PetscErrorCode DSSortEigenvalues_Private(DS ds,PetscScalar *wr,PetscScalar *wi,PetscInt *perm,PetscBool isghiep)
194: {
196: PetscScalar re,im,wi0;
197: PetscInt n,i,j,result,tmp1,tmp2=0,d=1;
200: n = ds->t; /* sort only first t pairs if truncated */
201: /* insertion sort */
202: i=ds->l+1;
203: #if !defined(PETSC_USE_COMPLEX)
204: if (wi && wi[perm[i-1]]!=0.0) i++; /* initial value is complex */
205: #else
206: if (isghiep && PetscImaginaryPart(wr[perm[i-1]])!=0.0) i++;
207: #endif
208: for (;i<n;i+=d) {
209: re = wr[perm[i]];
210: if (wi) im = wi[perm[i]];
211: else im = 0.0;
212: tmp1 = perm[i];
213: #if !defined(PETSC_USE_COMPLEX)
214: if (im!=0.0) { d = 2; tmp2 = perm[i+1]; }
215: else d = 1;
216: #else
217: if (isghiep && PetscImaginaryPart(re)!=0.0) { d = 2; tmp2 = perm[i+1]; }
218: else d = 1;
219: #endif
220: j = i-1;
221: if (wi) wi0 = wi[perm[j]];
222: else wi0 = 0.0;
223: SlepcSCCompare(ds->sc,re,im,wr[perm[j]],wi0,&result);
224: while (result<0 && j>=ds->l) {
225: perm[j+d] = perm[j];
226: j--;
227: #if !defined(PETSC_USE_COMPLEX)
228: if (wi && wi[perm[j+1]]!=0)
229: #else
230: if (isghiep && PetscImaginaryPart(wr[perm[j+1]])!=0)
231: #endif
232: { perm[j+d] = perm[j]; j--; }
234: if (j>=ds->l) {
235: if (wi) wi0 = wi[perm[j]];
236: else wi0 = 0.0;
237: SlepcSCCompare(ds->sc,re,im,wr[perm[j]],wi0,&result);
238: }
239: }
240: perm[j+1] = tmp1;
241: if (d==2) perm[j+2] = tmp2;
242: }
243: return(0);
244: }
248: PetscErrorCode DSSortEigenvaluesReal_Private(DS ds,PetscReal *eig,PetscInt *perm)
249: {
251: PetscScalar re;
252: PetscInt i,j,result,tmp,l,n;
255: n = ds->t; /* sort only first t pairs if truncated */
256: l = ds->l;
257: /* insertion sort */
258: for (i=l+1;i<n;i++) {
259: re = eig[perm[i]];
260: j = i-1;
261: SlepcSCCompare(ds->sc,re,0.0,eig[perm[j]],0.0,&result);
262: while (result<0 && j>=l) {
263: tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
264: if (j>=l) {
265: SlepcSCCompare(ds->sc,re,0.0,eig[perm[j]],0.0,&result);
266: }
267: }
268: }
269: return(0);
270: }
274: /*
275: DSCopyMatrix_Private - Copies the trailing block of a matrix (from
276: rows/columns l to n).
277: */
278: PetscErrorCode DSCopyMatrix_Private(DS ds,DSMatType dst,DSMatType src)
279: {
281: PetscInt j,m,off,ld;
282: PetscScalar *S,*D;
285: ld = ds->ld;
286: m = ds->n-ds->l;
287: off = ds->l+ds->l*ld;
288: S = ds->mat[src];
289: D = ds->mat[dst];
290: for (j=0;j<m;j++) {
291: PetscMemcpy(D+off+j*ld,S+off+j*ld,m*sizeof(PetscScalar));
292: }
293: return(0);
294: }
298: PetscErrorCode DSPermuteColumns_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat,PetscInt *perm)
299: {
300: PetscInt i,j,k,p,ld;
301: PetscScalar *Q,rtmp;
304: ld = ds->ld;
305: Q = ds->mat[mat];
306: for (i=l;i<n;i++) {
307: p = perm[i];
308: if (p != i) {
309: j = i + 1;
310: while (perm[j] != i) j++;
311: perm[j] = p; perm[i] = i;
312: /* swap columns i and j */
313: for (k=0;k<n;k++) {
314: rtmp = Q[k+p*ld]; Q[k+p*ld] = Q[k+i*ld]; Q[k+i*ld] = rtmp;
315: }
316: }
317: }
318: return(0);
319: }
323: PetscErrorCode DSPermuteRows_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat,PetscInt *perm)
324: {
325: PetscInt i,j,m=ds->m,k,p,ld;
326: PetscScalar *Q,rtmp;
329: if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
330: ld = ds->ld;
331: Q = ds->mat[mat];
332: for (i=l;i<n;i++) {
333: p = perm[i];
334: if (p != i) {
335: j = i + 1;
336: while (perm[j] != i) j++;
337: perm[j] = p; perm[i] = i;
338: /* swap rows i and j */
339: for (k=0;k<m;k++) {
340: rtmp = Q[p+k*ld]; Q[p+k*ld] = Q[i+k*ld]; Q[i+k*ld] = rtmp;
341: }
342: }
343: }
344: return(0);
345: }
349: PetscErrorCode DSPermuteBoth_Private(DS ds,PetscInt l,PetscInt n,DSMatType mat1,DSMatType mat2,PetscInt *perm)
350: {
351: PetscInt i,j,m=ds->m,k,p,ld;
352: PetscScalar *U,*VT,rtmp;
355: if (!m) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"m was not set");
356: ld = ds->ld;
357: U = ds->mat[mat1];
358: VT = ds->mat[mat2];
359: for (i=l;i<n;i++) {
360: p = perm[i];
361: if (p != i) {
362: j = i + 1;
363: while (perm[j] != i) j++;
364: perm[j] = p; perm[i] = i;
365: /* swap columns i and j of U */
366: for (k=0;k<n;k++) {
367: rtmp = U[k+p*ld]; U[k+p*ld] = U[k+i*ld]; U[k+i*ld] = rtmp;
368: }
369: /* swap rows i and j of VT */
370: for (k=0;k<m;k++) {
371: rtmp = VT[p+k*ld]; VT[p+k*ld] = VT[i+k*ld]; VT[i+k*ld] = rtmp;
372: }
373: }
374: }
375: return(0);
376: }
380: /*@
381: DSSetIdentity - Copy the identity (a diagonal matrix with ones) on the
382: active part of a matrix.
384: Logically Collective on DS
386: Input Parameters:
387: + ds - the direct solver context
388: - mat - the matrix to modify
390: Level: intermediate
391: @*/
392: PetscErrorCode DSSetIdentity(DS ds,DSMatType mat)
393: {
395: PetscScalar *x;
396: PetscInt i,ld,n,l;
399: DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
400: DSGetLeadingDimension(ds,&ld);
401: PetscLogEventBegin(DS_Other,ds,0,0,0);
402: DSGetArray(ds,mat,&x);
403: PetscMemzero(&x[ld*l],ld*(n-l)*sizeof(PetscScalar));
404: for (i=l;i<n;i++) {
405: x[ld*i+i] = 1.0;
406: }
407: DSRestoreArray(ds,mat,&x);
408: PetscLogEventEnd(DS_Other,ds,0,0,0);
409: return(0);
410: }
414: /*
415: DSOrthogonalize - Orthogonalize the columns of a matrix.
417: Input Parameters:
418: + ds - the direct solver context
419: . mat - a matrix
420: - cols - number of columns to orthogonalize (starting from the column zero)
422: Output Parameter:
423: . lindcols - number of linearly independent columns of the matrix (can be NULL)
424: */
425: PetscErrorCode DSOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscInt *lindcols)
426: {
427: #if defined(PETSC_MISSING_LAPACK_GEQRF) || defined(SLEPC_MISSING_LAPACK_ORGQR)
429: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQRF/ORGQR - Lapack routine is unavailable");
430: #else
431: PetscErrorCode ierr;
432: PetscInt n,l,ld;
433: PetscBLASInt ld_,rA,cA,info,ltau,lw;
434: PetscScalar *A,*tau,*w,saux;
440: DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
441: DSGetLeadingDimension(ds,&ld);
442: n = n - l;
443: if (cols > n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns");
444: if (n == 0 || cols == 0) return(0);
445: DSGetArray(ds,mat,&A);
446: PetscBLASIntCast(PetscMin(cols,n),<au);
447: PetscBLASIntCast(ld,&ld_);
448: PetscBLASIntCast(n,&rA);
449: PetscBLASIntCast(cols,&cA);
450: lw = -1;
451: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,A,&ld_,NULL,&saux,&lw,&info));
452: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
453: lw = (PetscBLASInt)PetscRealPart(saux);
454: DSAllocateWork_Private(ds,lw+ltau,0,0);
455: tau = ds->work;
456: w = &tau[ltau];
457: PetscLogEventBegin(DS_Other,ds,0,0,0);
458: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
459: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,&A[ld*l+l],&ld_,tau,w,&lw,&info));
460: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
461: PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&rA,<au,<au,&A[ld*l+l],&ld_,tau,w,&lw,&info));
462: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);
463: PetscFPTrapPop();
464: PetscLogEventEnd(DS_Other,ds,0,0,0);
465: DSRestoreArray(ds,mat,&A);
466: PetscObjectStateIncrease((PetscObject)ds);
467: if (lindcols) *lindcols = ltau;
468: return(0);
469: #endif
470: }
474: /*
475: Compute C <- a*A*B + b*C, where
476: ldC, the leading dimension of C,
477: ldA, the leading dimension of A,
478: rA, cA, rows and columns of A,
479: At, if true use the transpose of A instead,
480: ldB, the leading dimension of B,
481: rB, cB, rows and columns of B,
482: Bt, if true use the transpose of B instead
483: */
484: 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)
485: {
487: PetscInt tmp;
488: PetscBLASInt m, n, k, ldA = _ldA, ldB = _ldB, ldC = _ldC;
489: const char *N = "N", *T = "C", *qA = N, *qB = N;
492: if ((rA == 0) || (cB == 0)) return(0);
497: /* Transpose if needed */
498: if (At) tmp = rA, rA = cA, cA = tmp, qA = T;
499: if (Bt) tmp = rB, rB = cB, cB = tmp, qB = T;
501: /* Check size */
502: if (cA != rB) SETERRQ(PETSC_COMM_SELF,1, "Matrix dimensions do not match");
504: /* Do stub */
505: if ((rA == 1) && (cA == 1) && (cB == 1)) {
506: if (!At && !Bt) *C = *A * *B;
507: else if (At && !Bt) *C = PetscConj(*A) * *B;
508: else if (!At && Bt) *C = *A * PetscConj(*B);
509: else *C = PetscConj(*A) * PetscConj(*B);
510: m = n = k = 1;
511: } else {
512: m = rA; n = cB; k = cA;
513: PetscStackCallBLAS("BLASgemm",BLASgemm_(qA,qB,&m,&n,&k,&a,(PetscScalar*)A,&ldA,(PetscScalar*)B,&ldB,&b,C,&ldC));
514: }
516: PetscLogFlops(m*n*2*k);
517: return(0);
518: }
522: /*
523: DSPseudoOrthogonalize - Orthogonalize the columns of a matrix with Modified
524: Gram-Schmidt in an indefinite inner product space defined by a signature.
526: Input Parameters:
527: + ds - the direct solver context
528: . mat - the matrix
529: . cols - number of columns to orthogonalize (starting from the column zero)
530: - s - the signature that defines the inner product
532: Output Parameter:
533: + lindcols - linear independent columns of the matrix (can be NULL)
534: - ns - the new norm of the vectors (can be NULL)
535: */
536: PetscErrorCode DSPseudoOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscReal *s,PetscInt *lindcols,PetscReal *ns)
537: {
538: PetscErrorCode ierr;
539: PetscInt i,j,k,l,n,ld;
540: PetscBLASInt one=1,rA_;
541: PetscScalar alpha,*A,*A_,*m,*h,nr0;
542: PetscReal nr_o,nr,*ns_;
550: DSGetDimensions(ds,&n,NULL,&l,NULL,NULL);
551: DSGetLeadingDimension(ds,&ld);
552: n = n - l;
553: if (cols > n) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONG,"Invalid number of columns");
554: if (n == 0 || cols == 0) return(0);
555: PetscBLASIntCast(n,&rA_);
556: DSGetArray(ds,mat,&A_);
557: A = &A_[ld*l+l];
558: DSAllocateWork_Private(ds,n+cols,ns?0:cols,0);
559: m = ds->work;
560: h = &m[n];
561: ns_ = ns ? ns : ds->rwork;
562: PetscLogEventBegin(DS_Other,ds,0,0,0);
563: for (i=0; i<cols; i++) {
564: /* m <- diag(s)*A[i] */
565: for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
566: /* nr_o <- mynorm(A[i]'*m), mynorm(x) = sign(x)*sqrt(|x|) */
567: SlepcMatDenseMult(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
568: nr = nr_o = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
569: for (j=0; j<3 && i>0; j++) {
570: /* h <- A[0:i-1]'*m */
571: SlepcMatDenseMult(h,i,0.0,1.0,A,ld,n,i,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
572: /* h <- diag(ns)*h */
573: for (k=0; k<i; k++) h[k] *= ns_[k];
574: /* A[i] <- A[i] - A[0:i-1]*h */
575: SlepcMatDenseMult(&A[ld*i],ld,1.0,-1.0,A,ld,n,i,PETSC_FALSE,h,i,i,1,PETSC_FALSE);
576: /* m <- diag(s)*A[i] */
577: for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
578: /* nr_o <- mynorm(A[i]'*m) */
579: SlepcMatDenseMult(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
580: nr = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
581: if (PetscAbs(nr) < PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_SELF,1,"Linear dependency detected");
582: if (PetscAbs(nr) > 0.7*PetscAbs(nr_o)) break;
583: nr_o = nr;
584: }
585: ns_[i] = PetscSign(nr);
586: /* A[i] <- A[i]/|nr| */
587: alpha = 1.0/PetscAbs(nr);
588: PetscStackCallBLAS("BLASscal",BLASscal_(&rA_,&alpha,&A[i*ld],&one));
589: }
590: PetscLogEventEnd(DS_Other,ds,0,0,0);
591: DSRestoreArray(ds,mat,&A_);
592: PetscObjectStateIncrease((PetscObject)ds);
593: if (lindcols) *lindcols = cols;
594: return(0);
595: }