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