Actual source code: dspep.c
slepc-3.18.2 2023-01-26
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: */
11: #include <slepc/private/dsimpl.h>
12: #include <slepcblaslapack.h>
14: typedef struct {
15: PetscInt d; /* polynomial degree */
16: PetscReal *pbc; /* polynomial basis coefficients */
17: } DS_PEP;
19: PetscErrorCode DSAllocate_PEP(DS ds,PetscInt ld)
20: {
21: DS_PEP *ctx = (DS_PEP*)ds->data;
22: PetscInt i;
25: DSAllocateMat_Private(ds,DS_MAT_X);
26: DSAllocateMat_Private(ds,DS_MAT_Y);
27: for (i=0;i<=ctx->d;i++) DSAllocateMat_Private(ds,DSMatExtra[i]);
28: PetscFree(ds->perm);
29: PetscMalloc1(ld*ctx->d,&ds->perm);
30: return 0;
31: }
33: PetscErrorCode DSView_PEP(DS ds,PetscViewer viewer)
34: {
35: DS_PEP *ctx = (DS_PEP*)ds->data;
36: PetscViewerFormat format;
37: PetscInt i;
39: PetscViewerGetFormat(viewer,&format);
40: if (format == PETSC_VIEWER_ASCII_INFO) return 0;
41: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
42: PetscViewerASCIIPrintf(viewer,"polynomial degree: %" PetscInt_FMT "\n",ctx->d);
43: return 0;
44: }
45: for (i=0;i<=ctx->d;i++) DSViewMat(ds,viewer,DSMatExtra[i]);
46: if (ds->state>DS_STATE_INTERMEDIATE) DSViewMat(ds,viewer,DS_MAT_X);
47: return 0;
48: }
50: PetscErrorCode DSVectors_PEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
51: {
53: switch (mat) {
54: case DS_MAT_X:
55: break;
56: case DS_MAT_Y:
57: break;
58: default:
59: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
60: }
61: return 0;
62: }
64: PetscErrorCode DSSort_PEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *kout)
65: {
66: DS_PEP *ctx = (DS_PEP*)ds->data;
67: PetscInt n,i,*perm,told;
68: PetscScalar *A;
70: if (!ds->sc) return 0;
71: n = ds->n*ctx->d;
72: perm = ds->perm;
73: for (i=0;i<n;i++) perm[i] = i;
74: told = ds->t;
75: ds->t = n; /* force the sorting routines to consider d*n eigenvalues */
76: if (rr) DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
77: else DSSortEigenvalues_Private(ds,wr,wi,perm,PETSC_FALSE);
78: ds->t = told; /* restore value of t */
79: MatDenseGetArray(ds->omat[DS_MAT_A],&A);
80: for (i=0;i<n;i++) A[i] = wr[perm[i]];
81: for (i=0;i<n;i++) wr[i] = A[i];
82: for (i=0;i<n;i++) A[i] = wi[perm[i]];
83: for (i=0;i<n;i++) wi[i] = A[i];
84: MatDenseRestoreArray(ds->omat[DS_MAT_A],&A);
85: DSPermuteColumnsTwo_Private(ds,0,n,ds->n,DS_MAT_X,DS_MAT_Y,perm);
86: return 0;
87: }
89: PetscErrorCode DSSolve_PEP_QZ(DS ds,PetscScalar *wr,PetscScalar *wi)
90: {
91: DS_PEP *ctx = (DS_PEP*)ds->data;
92: PetscInt i,j,k,off;
93: PetscScalar *A,*B,*W,*X,*U,*Y,*work,*beta;
94: const PetscScalar *Ed,*Ei;
95: PetscReal *ca,*cb,*cg,norm,done=1.0;
96: PetscBLASInt info,n,ld,ldd,nd,lrwork=0,lwork,one=1,zero=0,cols;
97: #if defined(PETSC_USE_COMPLEX)
98: PetscReal *rwork;
99: #endif
101: PetscBLASIntCast(ds->n*ctx->d,&nd);
102: PetscBLASIntCast(ds->n,&n);
103: PetscBLASIntCast(ds->ld,&ld);
104: PetscBLASIntCast(ds->ld*ctx->d,&ldd);
105: #if defined(PETSC_USE_COMPLEX)
106: PetscBLASIntCast(nd+2*nd,&lwork);
107: PetscBLASIntCast(8*nd,&lrwork);
108: #else
109: PetscBLASIntCast(nd+8*nd,&lwork);
110: #endif
111: DSAllocateWork_Private(ds,lwork,lrwork,0);
112: beta = ds->work;
113: work = ds->work + nd;
114: lwork -= nd;
115: DSAllocateMat_Private(ds,DS_MAT_A);
116: DSAllocateMat_Private(ds,DS_MAT_B);
117: DSAllocateMat_Private(ds,DS_MAT_W);
118: DSAllocateMat_Private(ds,DS_MAT_U);
119: MatDenseGetArray(ds->omat[DS_MAT_A],&A);
120: MatDenseGetArray(ds->omat[DS_MAT_B],&B);
122: /* build matrices A and B of the linearization */
123: MatDenseGetArrayRead(ds->omat[DSMatExtra[ctx->d]],&Ed);
124: PetscArrayzero(A,ldd*ldd);
125: if (!ctx->pbc) { /* monomial basis */
126: for (i=0;i<nd-ds->n;i++) A[i+(i+ds->n)*ldd] = 1.0;
127: for (i=0;i<ctx->d;i++) {
128: MatDenseGetArrayRead(ds->omat[DSMatExtra[i]],&Ei);
129: off = i*ds->n*ldd+(ctx->d-1)*ds->n;
130: for (j=0;j<ds->n;j++) PetscArraycpy(A+off+j*ldd,Ei+j*ds->ld,ds->n);
131: MatDenseRestoreArrayRead(ds->omat[DSMatExtra[i]],&Ei);
132: }
133: } else {
134: ca = ctx->pbc;
135: cb = ca+ctx->d+1;
136: cg = cb+ctx->d+1;
137: for (i=0;i<ds->n;i++) {
138: A[i+(i+ds->n)*ldd] = ca[0];
139: A[i+i*ldd] = cb[0];
140: }
141: for (;i<nd-ds->n;i++) {
142: j = i/ds->n;
143: A[i+(i+ds->n)*ldd] = ca[j];
144: A[i+i*ldd] = cb[j];
145: A[i+(i-ds->n)*ldd] = cg[j];
146: }
147: for (i=0;i<ctx->d-2;i++) {
148: MatDenseGetArrayRead(ds->omat[DSMatExtra[i]],&Ei);
149: off = i*ds->n*ldd+(ctx->d-1)*ds->n;
150: for (j=0;j<ds->n;j++)
151: for (k=0;k<ds->n;k++)
152: A[off+j*ldd+k] = Ei[j*ds->ld+k]*ca[ctx->d-1];
153: MatDenseRestoreArrayRead(ds->omat[DSMatExtra[i]],&Ei);
154: }
155: MatDenseGetArrayRead(ds->omat[DSMatExtra[i]],&Ei);
156: off = i*ds->n*ldd+(ctx->d-1)*ds->n;
157: for (j=0;j<ds->n;j++)
158: for (k=0;k<ds->n;k++)
159: A[off+j*ldd+k] = Ei[j*ds->ld+k]*ca[ctx->d-1]-Ed[j*ds->ld+k]*cg[ctx->d-1];
160: MatDenseRestoreArrayRead(ds->omat[DSMatExtra[i]],&Ei);
161: i++;
162: MatDenseGetArrayRead(ds->omat[DSMatExtra[i]],&Ei);
163: off = i*ds->n*ldd+(ctx->d-1)*ds->n;
164: for (j=0;j<ds->n;j++)
165: for (k=0;k<ds->n;k++)
166: A[off+j*ldd+k] = Ei[j*ds->ld+k]*ca[ctx->d-1]-Ed[j*ds->ld+k]*cb[ctx->d-1];
167: MatDenseRestoreArrayRead(ds->omat[DSMatExtra[i]],&Ei);
168: }
169: PetscArrayzero(B,ldd*ldd);
170: for (i=0;i<nd-ds->n;i++) B[i+i*ldd] = 1.0;
171: off = (ctx->d-1)*ds->n*(ldd+1);
172: for (j=0;j<ds->n;j++) {
173: for (i=0;i<ds->n;i++) B[off+i+j*ldd] = -Ed[i+j*ds->ld];
174: }
175: MatDenseRestoreArrayRead(ds->omat[DSMatExtra[ctx->d]],&Ed);
177: /* solve generalized eigenproblem */
178: MatDenseGetArray(ds->omat[DS_MAT_W],&W);
179: MatDenseGetArray(ds->omat[DS_MAT_U],&U);
180: #if defined(PETSC_USE_COMPLEX)
181: rwork = ds->rwork;
182: PetscCallBLAS("LAPACKggev",LAPACKggev_("V","V",&nd,A,&ldd,B,&ldd,wr,beta,U,&ldd,W,&ldd,work,&lwork,rwork,&info));
183: #else
184: PetscCallBLAS("LAPACKggev",LAPACKggev_("V","V",&nd,A,&ldd,B,&ldd,wr,wi,beta,U,&ldd,W,&ldd,work,&lwork,&info));
185: #endif
186: SlepcCheckLapackInfo("ggev",info);
187: MatDenseRestoreArray(ds->omat[DS_MAT_A],&A);
188: MatDenseRestoreArray(ds->omat[DS_MAT_B],&B);
190: /* copy eigenvalues */
191: for (i=0;i<nd;i++) {
192: if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
193: else wr[i] /= beta[i];
194: #if !defined(PETSC_USE_COMPLEX)
195: if (beta[i]==0.0) wi[i] = 0.0;
196: else wi[i] /= beta[i];
197: #else
198: if (wi) wi[i] = 0.0;
199: #endif
200: }
202: /* copy and normalize eigenvectors */
203: MatDenseGetArray(ds->omat[DS_MAT_X],&X);
204: MatDenseGetArray(ds->omat[DS_MAT_Y],&Y);
205: for (j=0;j<nd;j++) {
206: PetscArraycpy(X+j*ds->ld,W+j*ldd,ds->n);
207: PetscArraycpy(Y+j*ds->ld,U+ds->n*(ctx->d-1)+j*ldd,ds->n);
208: }
209: MatDenseRestoreArray(ds->omat[DS_MAT_W],&W);
210: MatDenseRestoreArray(ds->omat[DS_MAT_U],&U);
211: for (j=0;j<nd;j++) {
212: cols = 1;
213: norm = BLASnrm2_(&n,X+j*ds->ld,&one);
214: #if !defined(PETSC_USE_COMPLEX)
215: if (wi[j] != 0.0) {
216: norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,X+(j+1)*ds->ld,&one));
217: cols = 2;
218: }
219: #endif
220: PetscCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,X+j*ds->ld,&ld,&info));
221: SlepcCheckLapackInfo("lascl",info);
222: norm = BLASnrm2_(&n,Y+j*ds->ld,&one);
223: #if !defined(PETSC_USE_COMPLEX)
224: if (wi[j] != 0.0) norm = SlepcAbsEigenvalue(norm,BLASnrm2_(&n,Y+(j+1)*ds->ld,&one));
225: #endif
226: PetscCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&cols,Y+j*ds->ld,&ld,&info));
227: SlepcCheckLapackInfo("lascl",info);
228: #if !defined(PETSC_USE_COMPLEX)
229: if (wi[j] != 0.0) j++;
230: #endif
231: }
232: MatDenseRestoreArray(ds->omat[DS_MAT_X],&X);
233: MatDenseRestoreArray(ds->omat[DS_MAT_Y],&Y);
234: return 0;
235: }
237: #if !defined(PETSC_HAVE_MPIUNI)
238: PetscErrorCode DSSynchronize_PEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
239: {
240: DS_PEP *ctx = (DS_PEP*)ds->data;
241: PetscInt ld=ds->ld,k=0;
242: PetscMPIInt ldnd,rank,off=0,size,dn;
243: PetscScalar *X,*Y;
245: if (ds->state>=DS_STATE_CONDENSED) k += 2*ctx->d*ds->n*ld;
246: if (eigr) k += ctx->d*ds->n;
247: if (eigi) k += ctx->d*ds->n;
248: DSAllocateWork_Private(ds,k,0,0);
249: PetscMPIIntCast(k*sizeof(PetscScalar),&size);
250: PetscMPIIntCast(ds->n*ctx->d*ld,&ldnd);
251: PetscMPIIntCast(ctx->d*ds->n,&dn);
252: if (ds->state>=DS_STATE_CONDENSED) {
253: MatDenseGetArray(ds->omat[DS_MAT_X],&X);
254: MatDenseGetArray(ds->omat[DS_MAT_Y],&Y);
255: }
256: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
257: if (!rank) {
258: if (ds->state>=DS_STATE_CONDENSED) {
259: MPI_Pack(X,ldnd,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
260: MPI_Pack(Y,ldnd,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
261: }
262: if (eigr) MPI_Pack(eigr,dn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
263: #if !defined(PETSC_USE_COMPLEX)
264: if (eigi) MPI_Pack(eigi,dn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
265: #endif
266: }
267: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
268: if (rank) {
269: if (ds->state>=DS_STATE_CONDENSED) {
270: MPI_Unpack(ds->work,size,&off,X,ldnd,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
271: MPI_Unpack(ds->work,size,&off,Y,ldnd,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
272: }
273: if (eigr) MPI_Unpack(ds->work,size,&off,eigr,dn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
274: #if !defined(PETSC_USE_COMPLEX)
275: if (eigi) MPI_Unpack(ds->work,size,&off,eigi,dn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
276: #endif
277: }
278: if (ds->state>=DS_STATE_CONDENSED) {
279: MatDenseRestoreArray(ds->omat[DS_MAT_X],&X);
280: MatDenseRestoreArray(ds->omat[DS_MAT_Y],&Y);
281: }
282: return 0;
283: }
284: #endif
286: static PetscErrorCode DSPEPSetDegree_PEP(DS ds,PetscInt d)
287: {
288: DS_PEP *ctx = (DS_PEP*)ds->data;
292: ctx->d = d;
293: return 0;
294: }
296: /*@
297: DSPEPSetDegree - Sets the polynomial degree for a DSPEP.
299: Logically Collective on ds
301: Input Parameters:
302: + ds - the direct solver context
303: - d - the degree
305: Level: intermediate
307: .seealso: DSPEPGetDegree()
308: @*/
309: PetscErrorCode DSPEPSetDegree(DS ds,PetscInt d)
310: {
313: PetscTryMethod(ds,"DSPEPSetDegree_C",(DS,PetscInt),(ds,d));
314: return 0;
315: }
317: static PetscErrorCode DSPEPGetDegree_PEP(DS ds,PetscInt *d)
318: {
319: DS_PEP *ctx = (DS_PEP*)ds->data;
321: *d = ctx->d;
322: return 0;
323: }
325: /*@
326: DSPEPGetDegree - Returns the polynomial degree for a DSPEP.
328: Not collective
330: Input Parameter:
331: . ds - the direct solver context
333: Output Parameters:
334: . d - the degree
336: Level: intermediate
338: .seealso: DSPEPSetDegree()
339: @*/
340: PetscErrorCode DSPEPGetDegree(DS ds,PetscInt *d)
341: {
344: PetscUseMethod(ds,"DSPEPGetDegree_C",(DS,PetscInt*),(ds,d));
345: return 0;
346: }
348: static PetscErrorCode DSPEPSetCoefficients_PEP(DS ds,PetscReal *pbc)
349: {
350: DS_PEP *ctx = (DS_PEP*)ds->data;
351: PetscInt i;
354: if (ctx->pbc) PetscFree(ctx->pbc);
355: PetscMalloc1(3*(ctx->d+1),&ctx->pbc);
356: for (i=0;i<3*(ctx->d+1);i++) ctx->pbc[i] = pbc[i];
357: ds->state = DS_STATE_RAW;
358: return 0;
359: }
361: /*@C
362: DSPEPSetCoefficients - Sets the polynomial basis coefficients for a DSPEP.
364: Logically Collective on ds
366: Input Parameters:
367: + ds - the direct solver context
368: - pbc - the polynomial basis coefficients
370: Notes:
371: This function is required only in the case of a polynomial specified in a
372: non-monomial basis, to provide the coefficients that will be used
373: during the linearization, multiplying the identity blocks on the three main
374: diagonal blocks. Depending on the polynomial basis (Chebyshev, Legendre, ...)
375: the coefficients must be different.
377: There must be a total of 3*(d+1) coefficients, where d is the degree of the
378: polynomial. The coefficients are arranged in three groups, alpha, beta, and
379: gamma, according to the definition of the three-term recurrence. In the case
380: of the monomial basis, alpha=1 and beta=gamma=0, in which case it is not
381: necessary to invoke this function.
383: Level: advanced
385: .seealso: DSPEPGetCoefficients(), DSPEPSetDegree()
386: @*/
387: PetscErrorCode DSPEPSetCoefficients(DS ds,PetscReal *pbc)
388: {
390: PetscTryMethod(ds,"DSPEPSetCoefficients_C",(DS,PetscReal*),(ds,pbc));
391: return 0;
392: }
394: static PetscErrorCode DSPEPGetCoefficients_PEP(DS ds,PetscReal **pbc)
395: {
396: DS_PEP *ctx = (DS_PEP*)ds->data;
397: PetscInt i;
400: PetscCalloc1(3*(ctx->d+1),pbc);
401: if (ctx->pbc) for (i=0;i<3*(ctx->d+1);i++) (*pbc)[i] = ctx->pbc[i];
402: else for (i=0;i<ctx->d+1;i++) (*pbc)[i] = 1.0;
403: return 0;
404: }
406: /*@C
407: DSPEPGetCoefficients - Returns the polynomial basis coefficients for a DSPEP.
409: Not collective
411: Input Parameter:
412: . ds - the direct solver context
414: Output Parameters:
415: . pbc - the polynomial basis coefficients
417: Note:
418: The returned array has length 3*(d+1) and should be freed by the user.
420: Fortran Notes:
421: The calling sequence from Fortran is
422: .vb
423: DSPEPGetCoefficients(eps,pbc,ierr)
424: double precision pbc(d+1) output
425: .ve
427: Level: advanced
429: .seealso: DSPEPSetCoefficients()
430: @*/
431: PetscErrorCode DSPEPGetCoefficients(DS ds,PetscReal **pbc)
432: {
435: PetscUseMethod(ds,"DSPEPGetCoefficients_C",(DS,PetscReal**),(ds,pbc));
436: return 0;
437: }
439: PetscErrorCode DSDestroy_PEP(DS ds)
440: {
441: DS_PEP *ctx = (DS_PEP*)ds->data;
443: if (ctx->pbc) PetscFree(ctx->pbc);
444: PetscFree(ds->data);
445: PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetDegree_C",NULL);
446: PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetDegree_C",NULL);
447: PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetCoefficients_C",NULL);
448: PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetCoefficients_C",NULL);
449: return 0;
450: }
452: PetscErrorCode DSMatGetSize_PEP(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
453: {
454: DS_PEP *ctx = (DS_PEP*)ds->data;
457: *rows = ds->n;
458: if (t==DS_MAT_A || t==DS_MAT_B || t==DS_MAT_W || t==DS_MAT_U) *rows *= ctx->d;
459: *cols = ds->n;
460: if (t==DS_MAT_A || t==DS_MAT_B || t==DS_MAT_W || t==DS_MAT_U || t==DS_MAT_X || t==DS_MAT_Y) *cols *= ctx->d;
461: return 0;
462: }
464: /*MC
465: DSPEP - Dense Polynomial Eigenvalue Problem.
467: Level: beginner
469: Notes:
470: The problem is expressed as P(lambda)*x = 0, where P(.) is a matrix
471: polynomial of degree d. The eigenvalues lambda are the arguments
472: returned by DSSolve().
474: The degree of the polynomial, d, can be set with DSPEPSetDegree(), with
475: the first d+1 extra matrices of the DS defining the matrix polynomial. By
476: default, the polynomial is expressed in the monomial basis, but a
477: different basis can be used by setting the corresponding coefficients
478: via DSPEPSetCoefficients().
480: The problem is solved via linearization, by building a pencil (A,B) of
481: size p*n and solving the corresponding GNHEP.
483: Used DS matrices:
484: + DS_MAT_Ex - coefficients of the matrix polynomial
485: . DS_MAT_A - (workspace) first matrix of the linearization
486: . DS_MAT_B - (workspace) second matrix of the linearization
487: . DS_MAT_W - (workspace) right eigenvectors of the linearization
488: - DS_MAT_U - (workspace) left eigenvectors of the linearization
490: Implemented methods:
491: . 0 - QZ iteration on the linearization (_ggev)
493: .seealso: DSCreate(), DSSetType(), DSType, DSPEPSetDegree(), DSPEPSetCoefficients()
494: M*/
495: SLEPC_EXTERN PetscErrorCode DSCreate_PEP(DS ds)
496: {
497: DS_PEP *ctx;
499: PetscNew(&ctx);
500: ds->data = (void*)ctx;
502: ds->ops->allocate = DSAllocate_PEP;
503: ds->ops->view = DSView_PEP;
504: ds->ops->vectors = DSVectors_PEP;
505: ds->ops->solve[0] = DSSolve_PEP_QZ;
506: ds->ops->sort = DSSort_PEP;
507: #if !defined(PETSC_HAVE_MPIUNI)
508: ds->ops->synchronize = DSSynchronize_PEP;
509: #endif
510: ds->ops->destroy = DSDestroy_PEP;
511: ds->ops->matgetsize = DSMatGetSize_PEP;
512: PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetDegree_C",DSPEPSetDegree_PEP);
513: PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetDegree_C",DSPEPGetDegree_PEP);
514: PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetCoefficients_C",DSPEPSetCoefficients_PEP);
515: PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetCoefficients_C",DSPEPGetCoefficients_PEP);
516: return 0;
517: }