Actual source code: dspriv.c
slepc-3.18.0 2022-10-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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 DSAllocateMat_Private(DS ds,DSMatType m)
18: {
19: PetscInt n,d,rows=0,cols;
20: PetscBool ispep,isnep;
22: n = ds->ld;
23: PetscObjectTypeCompare((PetscObject)ds,DSPEP,&ispep);
24: if (ispep) {
25: if (m==DS_MAT_A || m==DS_MAT_B || m==DS_MAT_W || m==DS_MAT_U || m==DS_MAT_X || m==DS_MAT_Y) {
26: DSPEPGetDegree(ds,&d);
27: n = d*ds->ld;
28: }
29: } else {
30: PetscObjectTypeCompare((PetscObject)ds,DSNEP,&isnep);
31: if (isnep) {
32: if (m==DS_MAT_Q || m==DS_MAT_Z || m==DS_MAT_U || m==DS_MAT_V || m==DS_MAT_X || m==DS_MAT_Y) {
33: DSNEPGetMinimality(ds,&d);
34: n = d*ds->ld;
35: }
36: }
37: }
38: cols = n;
40: switch (m) {
41: case DS_MAT_T:
42: cols = PetscDefined(USE_COMPLEX)? 2: 3;
43: rows = n;
44: break;
45: case DS_MAT_D:
46: cols = 1;
47: rows = n;
48: break;
49: case DS_MAT_X:
50: rows = ds->ld;
51: break;
52: case DS_MAT_Y:
53: rows = ds->ld;
54: break;
55: default:
56: rows = n;
57: }
58: if (ds->omat[m]) MatZeroEntries(ds->omat[m]);
59: else {
60: MatCreateSeqDense(PETSC_COMM_SELF,rows,cols,NULL,&ds->omat[m]);
61: }
62: return 0;
63: }
65: PetscErrorCode DSAllocateWork_Private(DS ds,PetscInt s,PetscInt r,PetscInt i)
66: {
67: if (s>ds->lwork) {
68: PetscFree(ds->work);
69: PetscMalloc1(s,&ds->work);
70: ds->lwork = s;
71: }
72: if (r>ds->lrwork) {
73: PetscFree(ds->rwork);
74: PetscMalloc1(r,&ds->rwork);
75: ds->lrwork = r;
76: }
77: if (i>ds->liwork) {
78: PetscFree(ds->iwork);
79: PetscMalloc1(i,&ds->iwork);
80: ds->liwork = i;
81: }
82: return 0;
83: }
85: /*@C
86: DSViewMat - Prints one of the internal DS matrices.
88: Collective on ds
90: Input Parameters:
91: + ds - the direct solver context
92: . viewer - visualization context
93: - m - matrix to display
95: Note:
96: Works only for ascii viewers. Set the viewer in Matlab format if
97: want to paste into Matlab.
99: Level: developer
101: .seealso: DSView()
102: @*/
103: PetscErrorCode DSViewMat(DS ds,PetscViewer viewer,DSMatType m)
104: {
105: PetscInt i,j,rows,cols;
106: const PetscScalar *M=NULL,*v;
107: PetscViewerFormat format;
108: #if defined(PETSC_USE_COMPLEX)
109: PetscBool allreal = PETSC_TRUE;
110: #endif
114: DSCheckValidMat(ds,m,3);
115: if (!viewer) PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ds),&viewer);
118: PetscViewerGetFormat(viewer,&format);
119: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return 0;
120: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
121: DSMatGetSize(ds,m,&rows,&cols);
122: MatDenseGetArrayRead(ds->omat[m],&M);
123: #if defined(PETSC_USE_COMPLEX)
124: /* determine if matrix has all real values */
125: v = M;
126: for (i=0;i<rows;i++)
127: for (j=0;j<cols;j++)
128: if (PetscImaginaryPart(v[i+j*ds->ld])) { allreal = PETSC_FALSE; break; }
129: #endif
130: if (format == PETSC_VIEWER_ASCII_MATLAB) {
131: PetscViewerASCIIPrintf(viewer,"%% Size = %" PetscInt_FMT " %" PetscInt_FMT "\n",rows,cols);
132: PetscViewerASCIIPrintf(viewer,"%s = [\n",DSMatName[m]);
133: } else PetscViewerASCIIPrintf(viewer,"Matrix %s =\n",DSMatName[m]);
135: for (i=0;i<rows;i++) {
136: v = M+i;
137: for (j=0;j<cols;j++) {
138: #if defined(PETSC_USE_COMPLEX)
139: if (allreal) PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v));
140: else PetscViewerASCIIPrintf(viewer,"%18.16e%+18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));
141: #else
142: PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v);
143: #endif
144: v += ds->ld;
145: }
146: PetscViewerASCIIPrintf(viewer,"\n");
147: }
148: MatDenseRestoreArrayRead(ds->omat[m],&M);
150: if (format == PETSC_VIEWER_ASCII_MATLAB) PetscViewerASCIIPrintf(viewer,"];\n");
151: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
152: PetscViewerFlush(viewer);
153: return 0;
154: }
156: PetscErrorCode DSSortEigenvalues_Private(DS ds,PetscScalar *wr,PetscScalar *wi,PetscInt *perm,PetscBool isghiep)
157: {
158: PetscScalar re,im,wi0;
159: PetscInt n,i,j,result,tmp1,tmp2=0,d=1;
161: n = ds->t; /* sort only first t pairs if truncated */
162: /* insertion sort */
163: i=ds->l+1;
164: #if !defined(PETSC_USE_COMPLEX)
165: if (wi && wi[perm[i-1]]!=0.0) i++; /* initial value is complex */
166: #else
167: if (isghiep && PetscImaginaryPart(wr[perm[i-1]])!=0.0) i++;
168: #endif
169: for (;i<n;i+=d) {
170: re = wr[perm[i]];
171: if (wi) im = wi[perm[i]];
172: else im = 0.0;
173: tmp1 = perm[i];
174: #if !defined(PETSC_USE_COMPLEX)
175: if (im!=0.0) { d = 2; tmp2 = perm[i+1]; }
176: else d = 1;
177: #else
178: if (isghiep && PetscImaginaryPart(re)!=0.0) { d = 2; tmp2 = perm[i+1]; }
179: else d = 1;
180: #endif
181: j = i-1;
182: if (wi) wi0 = wi[perm[j]];
183: else wi0 = 0.0;
184: SlepcSCCompare(ds->sc,re,im,wr[perm[j]],wi0,&result);
185: while (result<0 && j>=ds->l) {
186: perm[j+d] = perm[j];
187: j--;
188: #if !defined(PETSC_USE_COMPLEX)
189: if (wi && wi[perm[j+1]]!=0)
190: #else
191: if (isghiep && PetscImaginaryPart(wr[perm[j+1]])!=0)
192: #endif
193: { perm[j+d] = perm[j]; j--; }
195: if (j>=ds->l) {
196: if (wi) wi0 = wi[perm[j]];
197: else wi0 = 0.0;
198: SlepcSCCompare(ds->sc,re,im,wr[perm[j]],wi0,&result);
199: }
200: }
201: perm[j+1] = tmp1;
202: if (d==2) perm[j+2] = tmp2;
203: }
204: return 0;
205: }
207: PetscErrorCode DSSortEigenvaluesReal_Private(DS ds,PetscReal *eig,PetscInt *perm)
208: {
209: PetscScalar re;
210: PetscInt i,j,result,tmp,l,n;
212: n = ds->t; /* sort only first t pairs if truncated */
213: l = ds->l;
214: /* insertion sort */
215: for (i=l+1;i<n;i++) {
216: re = eig[perm[i]];
217: j = i-1;
218: SlepcSCCompare(ds->sc,re,0.0,eig[perm[j]],0.0,&result);
219: while (result<0 && j>=l) {
220: tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
221: if (j>=l) SlepcSCCompare(ds->sc,re,0.0,eig[perm[j]],0.0,&result);
222: }
223: }
224: return 0;
225: }
227: /*
228: Permute comumns [istart..iend-1] of [mat] according to perm. Columns have length n
229: */
230: PetscErrorCode DSPermuteColumns_Private(DS ds,PetscInt istart,PetscInt iend,PetscInt n,DSMatType mat,PetscInt *perm)
231: {
232: PetscInt i,j,k,p,ld;
233: PetscScalar *Q,rtmp;
235: ld = ds->ld;
236: MatDenseGetArray(ds->omat[mat],&Q);
237: for (i=istart;i<iend;i++) {
238: p = perm[i];
239: if (p != i) {
240: j = i + 1;
241: while (perm[j] != i) j++;
242: perm[j] = p; perm[i] = i;
243: /* swap columns i and j */
244: for (k=0;k<n;k++) {
245: rtmp = Q[k+p*ld]; Q[k+p*ld] = Q[k+i*ld]; Q[k+i*ld] = rtmp;
246: }
247: }
248: }
249: MatDenseRestoreArray(ds->omat[mat],&Q);
250: return 0;
251: }
253: /*
254: The same as DSPermuteColumns_Private but for two matrices [mat1] and [mat2]
255: */
256: PetscErrorCode DSPermuteColumnsTwo_Private(DS ds,PetscInt istart,PetscInt iend,PetscInt n,DSMatType mat1,DSMatType mat2,PetscInt *perm)
257: {
258: PetscInt i,j,k,p,ld;
259: PetscScalar *Q,*Z,rtmp,rtmp2;
261: ld = ds->ld;
262: MatDenseGetArray(ds->omat[mat1],&Q);
263: MatDenseGetArray(ds->omat[mat2],&Z);
264: for (i=istart;i<iend;i++) {
265: p = perm[i];
266: if (p != i) {
267: j = i + 1;
268: while (perm[j] != i) j++;
269: perm[j] = p; perm[i] = i;
270: /* swap columns i and j */
271: for (k=0;k<n;k++) {
272: rtmp = Q[k+p*ld]; Q[k+p*ld] = Q[k+i*ld]; Q[k+i*ld] = rtmp;
273: rtmp2 = Z[k+p*ld]; Z[k+p*ld] = Z[k+i*ld]; Z[k+i*ld] = rtmp2;
274: }
275: }
276: }
277: MatDenseRestoreArray(ds->omat[mat1],&Q);
278: MatDenseRestoreArray(ds->omat[mat2],&Z);
279: return 0;
280: }
282: /*
283: Permute rows [istart..iend-1] of [mat] according to perm. Rows have length m
284: */
285: PetscErrorCode DSPermuteRows_Private(DS ds,PetscInt istart,PetscInt iend,PetscInt m,DSMatType mat,PetscInt *perm)
286: {
287: PetscInt i,j,k,p,ld;
288: PetscScalar *Q,rtmp;
290: ld = ds->ld;
291: MatDenseGetArray(ds->omat[mat],&Q);
292: for (i=istart;i<iend;i++) {
293: p = perm[i];
294: if (p != i) {
295: j = i + 1;
296: while (perm[j] != i) j++;
297: perm[j] = p; perm[i] = i;
298: /* swap rows i and j */
299: for (k=0;k<m;k++) {
300: rtmp = Q[p+k*ld]; Q[p+k*ld] = Q[i+k*ld]; Q[i+k*ld] = rtmp;
301: }
302: }
303: }
304: MatDenseRestoreArray(ds->omat[mat],&Q);
305: return 0;
306: }
308: /*
309: Permute columns [istart..iend-1] of [mat1] and [mat2] according to perm.
310: Columns of [mat1] have length n, columns of [mat2] have length m
311: */
312: PetscErrorCode DSPermuteBoth_Private(DS ds,PetscInt istart,PetscInt iend,PetscInt n,PetscInt m,DSMatType mat1,DSMatType mat2,PetscInt *perm)
313: {
314: PetscInt i,j,k,p,ld;
315: PetscScalar *U,*V,rtmp;
317: ld = ds->ld;
318: MatDenseGetArray(ds->omat[mat1],&U);
319: MatDenseGetArray(ds->omat[mat2],&V);
320: for (i=istart;i<iend;i++) {
321: p = perm[i];
322: if (p != i) {
323: j = i + 1;
324: while (perm[j] != i) j++;
325: perm[j] = p; perm[i] = i;
326: /* swap columns i and j of U */
327: for (k=0;k<n;k++) {
328: rtmp = U[k+p*ld]; U[k+p*ld] = U[k+i*ld]; U[k+i*ld] = rtmp;
329: }
330: /* swap columns i and j of V */
331: for (k=0;k<m;k++) {
332: rtmp = V[k+p*ld]; V[k+p*ld] = V[k+i*ld]; V[k+i*ld] = rtmp;
333: }
334: }
335: }
336: MatDenseRestoreArray(ds->omat[mat1],&U);
337: MatDenseRestoreArray(ds->omat[mat2],&V);
338: return 0;
339: }
341: /*@
342: DSSetIdentity - Copy the identity (a diagonal matrix with ones) on the
343: active part of a matrix.
345: Logically Collective on ds
347: Input Parameters:
348: + ds - the direct solver context
349: - mat - the matrix to modify
351: Level: intermediate
353: .seealso: DSGetMat()
354: @*/
355: PetscErrorCode DSSetIdentity(DS ds,DSMatType mat)
356: {
357: PetscScalar *A;
358: PetscInt i,ld,n,l;
362: DSCheckValidMat(ds,mat,2);
364: DSGetDimensions(ds,&n,&l,NULL,NULL);
365: DSGetLeadingDimension(ds,&ld);
366: PetscLogEventBegin(DS_Other,ds,0,0,0);
367: MatDenseGetArray(ds->omat[mat],&A);
368: PetscArrayzero(A+l*ld,ld*(n-l));
369: for (i=l;i<n;i++) A[i+i*ld] = 1.0;
370: MatDenseRestoreArray(ds->omat[mat],&A);
371: PetscLogEventEnd(DS_Other,ds,0,0,0);
372: return 0;
373: }
375: /*@C
376: DSOrthogonalize - Orthogonalize the columns of a matrix.
378: Logically Collective on ds
380: Input Parameters:
381: + ds - the direct solver context
382: . mat - a matrix
383: - cols - number of columns to orthogonalize (starting from column zero)
385: Output Parameter:
386: . lindcols - (optional) number of linearly independent columns of the matrix
388: Level: developer
390: .seealso: DSPseudoOrthogonalize()
391: @*/
392: PetscErrorCode DSOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscInt *lindcols)
393: {
394: PetscInt n,l,ld;
395: PetscBLASInt ld_,rA,cA,info,ltau,lw;
396: PetscScalar *A,*tau,*w,saux,dummy;
399: DSCheckAlloc(ds,1);
401: DSCheckValidMat(ds,mat,2);
404: DSGetDimensions(ds,&n,&l,NULL,NULL);
405: DSGetLeadingDimension(ds,&ld);
406: n = n - l;
408: if (n == 0 || cols == 0) return 0;
410: PetscLogEventBegin(DS_Other,ds,0,0,0);
411: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
412: MatDenseGetArray(ds->omat[mat],&A);
413: PetscBLASIntCast(PetscMin(cols,n),<au);
414: PetscBLASIntCast(ld,&ld_);
415: PetscBLASIntCast(n,&rA);
416: PetscBLASIntCast(cols,&cA);
417: lw = -1;
418: PetscCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,A,&ld_,&dummy,&saux,&lw,&info));
419: SlepcCheckLapackInfo("geqrf",info);
420: lw = (PetscBLASInt)PetscRealPart(saux);
421: DSAllocateWork_Private(ds,lw+ltau,0,0);
422: tau = ds->work;
423: w = &tau[ltau];
424: PetscCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&rA,&cA,&A[ld*l+l],&ld_,tau,w,&lw,&info));
425: SlepcCheckLapackInfo("geqrf",info);
426: PetscCallBLAS("LAPACKorgqr",LAPACKorgqr_(&rA,<au,<au,&A[ld*l+l],&ld_,tau,w,&lw,&info));
427: SlepcCheckLapackInfo("orgqr",info);
428: if (lindcols) *lindcols = ltau;
429: PetscFPTrapPop();
430: MatDenseRestoreArray(ds->omat[mat],&A);
431: PetscLogEventEnd(DS_Other,ds,0,0,0);
432: PetscObjectStateIncrease((PetscObject)ds);
433: return 0;
434: }
436: /*
437: Compute C <- a*A*B + b*C, where
438: ldC, the leading dimension of C,
439: ldA, the leading dimension of A,
440: rA, cA, rows and columns of A,
441: At, if true use the transpose of A instead,
442: ldB, the leading dimension of B,
443: rB, cB, rows and columns of B,
444: Bt, if true use the transpose of B instead
445: */
446: 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)
447: {
448: PetscInt tmp;
449: PetscBLASInt m, n, k, ldA = _ldA, ldB = _ldB, ldC = _ldC;
450: const char *N = "N", *T = "C", *qA = N, *qB = N;
452: if ((rA == 0) || (cB == 0)) return 0;
457: /* Transpose if needed */
458: if (At) tmp = rA, rA = cA, cA = tmp, qA = T;
459: if (Bt) tmp = rB, rB = cB, cB = tmp, qB = T;
461: /* Check size */
464: /* Do stub */
465: if ((rA == 1) && (cA == 1) && (cB == 1)) {
466: if (!At && !Bt) *C = *A * *B;
467: else if (At && !Bt) *C = PetscConj(*A) * *B;
468: else if (!At && Bt) *C = *A * PetscConj(*B);
469: else *C = PetscConj(*A) * PetscConj(*B);
470: m = n = k = 1;
471: } else {
472: m = rA; n = cB; k = cA;
473: PetscCallBLAS("BLASgemm",BLASgemm_(qA,qB,&m,&n,&k,&a,(PetscScalar*)A,&ldA,(PetscScalar*)B,&ldB,&b,C,&ldC));
474: }
476: PetscLogFlops(2.0*m*n*k);
477: return 0;
478: }
480: /*@C
481: DSPseudoOrthogonalize - Orthogonalize the columns of a matrix with Modified
482: Gram-Schmidt in an indefinite inner product space defined by a signature.
484: Logically Collective on ds
486: Input Parameters:
487: + ds - the direct solver context
488: . mat - the matrix
489: . cols - number of columns to orthogonalize (starting from column zero)
490: - s - the signature that defines the inner product
492: Output Parameters:
493: + lindcols - (optional) linearly independent columns of the matrix
494: - ns - (optional) the new signature of the vectors
496: Note:
497: After the call the matrix satisfies A'*s*A = ns.
499: Level: developer
501: .seealso: DSOrthogonalize()
502: @*/
503: PetscErrorCode DSPseudoOrthogonalize(DS ds,DSMatType mat,PetscInt cols,PetscReal *s,PetscInt *lindcols,PetscReal *ns)
504: {
505: PetscInt i,j,k,l,n,ld;
506: PetscBLASInt info,one=1,zero=0,rA_,ld_;
507: PetscScalar *A,*A_,*m,*h,nr0;
508: PetscReal nr_o,nr,nr_abs,*ns_,done=1.0;
511: DSCheckAlloc(ds,1);
513: DSCheckValidMat(ds,mat,2);
516: DSGetDimensions(ds,&n,&l,NULL,NULL);
517: DSGetLeadingDimension(ds,&ld);
518: n = n - l;
520: if (n == 0 || cols == 0) return 0;
521: PetscLogEventBegin(DS_Other,ds,0,0,0);
522: PetscBLASIntCast(n,&rA_);
523: PetscBLASIntCast(ld,&ld_);
524: MatDenseGetArray(ds->omat[mat],&A_);
525: A = &A_[ld*l+l];
526: DSAllocateWork_Private(ds,n+cols,ns?0:cols,0);
527: m = ds->work;
528: h = &m[n];
529: ns_ = ns ? ns : ds->rwork;
530: for (i=0; i<cols; i++) {
531: /* m <- diag(s)*A[i] */
532: for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
533: /* nr_o <- mynorm(A[i]'*m), mynorm(x) = sign(x)*sqrt(|x|) */
534: SlepcMatDenseMult(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
535: nr = nr_o = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
536: for (j=0; j<3 && i>0; j++) {
537: /* h <- A[0:i-1]'*m */
538: SlepcMatDenseMult(h,i,0.0,1.0,A,ld,n,i,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
539: /* h <- diag(ns)*h */
540: for (k=0; k<i; k++) h[k] *= ns_[k];
541: /* A[i] <- A[i] - A[0:i-1]*h */
542: SlepcMatDenseMult(&A[ld*i],ld,1.0,-1.0,A,ld,n,i,PETSC_FALSE,h,i,i,1,PETSC_FALSE);
543: /* m <- diag(s)*A[i] */
544: for (k=0; k<n; k++) m[k] = s[k]*A[k+i*ld];
545: /* nr_o <- mynorm(A[i]'*m) */
546: SlepcMatDenseMult(&nr0,1,0.0,1.0,&A[ld*i],ld,n,1,PETSC_TRUE,m,n,n,1,PETSC_FALSE);
547: nr = PetscSign(PetscRealPart(nr0))*PetscSqrtReal(PetscAbsScalar(nr0));
549: if (PetscAbs(nr) > 0.7*PetscAbs(nr_o)) break;
550: nr_o = nr;
551: }
552: ns_[i] = PetscSign(nr);
553: /* A[i] <- A[i]/|nr| */
554: nr_abs = PetscAbs(nr);
555: PetscCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&nr_abs,&done,&rA_,&one,A+i*ld,&ld_,&info));
556: SlepcCheckLapackInfo("lascl",info);
557: }
558: MatDenseRestoreArray(ds->omat[mat],&A_);
559: PetscLogEventEnd(DS_Other,ds,0,0,0);
560: PetscObjectStateIncrease((PetscObject)ds);
561: if (lindcols) *lindcols = cols;
562: return 0;
563: }