Actual source code: dspep.c

slepc-3.12.0 2019-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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>       /*I "slepcds.h" I*/
 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: {
 22:   DS_PEP         *ctx = (DS_PEP*)ds->data;
 23:   PetscInt       i;

 26:   if (!ctx->d) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"DSPEP requires specifying the polynomial degree via DSPEPSetDegree()");
 27:   DSAllocateMat_Private(ds,DS_MAT_X);
 28:   DSAllocateMat_Private(ds,DS_MAT_Y);
 29:   for (i=0;i<=ctx->d;i++) {
 30:     DSAllocateMat_Private(ds,DSMatExtra[i]);
 31:   }
 32:   PetscFree(ds->perm);
 33:   PetscMalloc1(ld*ctx->d,&ds->perm);
 34:   PetscLogObjectMemory((PetscObject)ds,ld*ctx->d*sizeof(PetscInt));
 35:   return(0);
 36: }

 38: PetscErrorCode DSView_PEP(DS ds,PetscViewer viewer)
 39: {
 40:   PetscErrorCode    ierr;
 41:   DS_PEP            *ctx = (DS_PEP*)ds->data;
 42:   PetscViewerFormat format;
 43:   PetscInt          i;

 46:   PetscViewerGetFormat(viewer,&format);
 47:   PetscViewerASCIIPrintf(viewer,"polynomial degree: %D\n",ctx->d);
 48:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
 49:   for (i=0;i<=ctx->d;i++) {
 50:     DSViewMat(ds,viewer,DSMatExtra[i]);
 51:   }
 52:   if (ds->state>DS_STATE_INTERMEDIATE) { DSViewMat(ds,viewer,DS_MAT_X); }
 53:   return(0);
 54: }

 56: PetscErrorCode DSVectors_PEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
 57: {
 59:   if (rnorm) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
 60:   switch (mat) {
 61:     case DS_MAT_X:
 62:       break;
 63:     case DS_MAT_Y:
 64:       break;
 65:     default:
 66:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
 67:   }
 68:   return(0);
 69: }

 71: PetscErrorCode DSSort_PEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *kout)
 72: {
 74:   DS_PEP         *ctx = (DS_PEP*)ds->data;
 75:   PetscInt       n,i,j,k,p,*perm,told,ld;
 76:   PetscScalar    *A,*X,*Y,rtmp,rtmp2;

 79:   if (!ds->sc) return(0);
 80:   n = ds->n*ctx->d;
 81:   A  = ds->mat[DS_MAT_A];
 82:   perm = ds->perm;
 83:   for (i=0;i<n;i++) perm[i] = i;
 84:   told = ds->t;
 85:   ds->t = n;  /* force the sorting routines to consider d*n eigenvalues */
 86:   if (rr) {
 87:     DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
 88:   } else {
 89:     DSSortEigenvalues_Private(ds,wr,wi,perm,PETSC_FALSE);
 90:   }
 91:   ds->t = told;  /* restore value of t */
 92:   for (i=0;i<n;i++) A[i]  = wr[perm[i]];
 93:   for (i=0;i<n;i++) wr[i] = A[i];
 94:   for (i=0;i<n;i++) A[i]  = wi[perm[i]];
 95:   for (i=0;i<n;i++) wi[i] = A[i];
 96:   /* cannot use DSPermuteColumns_Private() since matrix is not square */
 97:   ld = ds->ld;
 98:   X  = ds->mat[DS_MAT_X];
 99:   Y  = ds->mat[DS_MAT_Y];
100:   for (i=0;i<n;i++) {
101:     p = perm[i];
102:     if (p != i) {
103:       j = i + 1;
104:       while (perm[j] != i) j++;
105:       perm[j] = p; perm[i] = i;
106:       /* swap columns i and j */
107:       for (k=0;k<ds->n;k++) {
108:         rtmp  = X[k+p*ld]; X[k+p*ld] = X[k+i*ld]; X[k+i*ld] = rtmp;
109:         rtmp2 = Y[k+p*ld]; Y[k+p*ld] = Y[k+i*ld]; Y[k+i*ld] = rtmp2;
110:       }
111:     }
112:   }
113:   return(0);
114: }

116: PetscErrorCode DSSolve_PEP_QZ(DS ds,PetscScalar *wr,PetscScalar *wi)
117: {
118: #if defined(SLEPC_MISSING_LAPACK_GGEV)
120:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GGEV - Lapack routine is unavailable");
121: #else
123:   DS_PEP         *ctx = (DS_PEP*)ds->data;
124:   PetscInt       i,j,k,off;
125:   PetscScalar    *A,*B,*W,*X,*U,*Y,*E,*work,*beta,norm;
126:   PetscReal      *ca,*cb,*cg;
127:   PetscBLASInt   info,n,ldd,nd,lrwork=0,lwork,one=1;
128: #if defined(PETSC_USE_COMPLEX)
129:   PetscReal      *rwork;
130: #else
131:   PetscScalar    norm0;
132: #endif

135:   if (!ds->mat[DS_MAT_A]) {
136:     DSAllocateMat_Private(ds,DS_MAT_A);
137:   }
138:   if (!ds->mat[DS_MAT_B]) {
139:     DSAllocateMat_Private(ds,DS_MAT_B);
140:   }
141:   if (!ds->mat[DS_MAT_W]) {
142:     DSAllocateMat_Private(ds,DS_MAT_W);
143:   }
144:   if (!ds->mat[DS_MAT_U]) {
145:     DSAllocateMat_Private(ds,DS_MAT_U);
146:   }
147:   PetscBLASIntCast(ds->n*ctx->d,&nd);
148:   PetscBLASIntCast(ds->n,&n);
149:   PetscBLASIntCast(ds->ld*ctx->d,&ldd);
150: #if defined(PETSC_USE_COMPLEX)
151:   PetscBLASIntCast(nd+2*nd,&lwork);
152:   PetscBLASIntCast(8*nd,&lrwork);
153: #else
154:   PetscBLASIntCast(nd+8*nd,&lwork);
155: #endif
156:   DSAllocateWork_Private(ds,lwork,lrwork,0);
157:   beta = ds->work;
158:   work = ds->work + nd;
159:   lwork -= nd;
160:   A = ds->mat[DS_MAT_A];
161:   B = ds->mat[DS_MAT_B];
162:   W = ds->mat[DS_MAT_W];
163:   U = ds->mat[DS_MAT_U];
164:   X = ds->mat[DS_MAT_X];
165:   Y = ds->mat[DS_MAT_Y];
166:   E = ds->mat[DSMatExtra[ctx->d]];

168:   /* build matrices A and B of the linearization */
169:   PetscArrayzero(A,ldd*ldd);
170:   if (!ctx->pbc) { /* monomial basis */
171:     for (i=0;i<nd-ds->n;i++) A[i+(i+ds->n)*ldd] = 1.0;
172:     for (i=0;i<ctx->d;i++) {
173:       off = i*ds->n*ldd+(ctx->d-1)*ds->n;
174:       for (j=0;j<ds->n;j++) {
175:         PetscArraycpy(A+off+j*ldd,ds->mat[DSMatExtra[i]]+j*ds->ld,ds->n);
176:       }
177:     }
178:   } else {
179:     ca = ctx->pbc;
180:     cb = ca+ctx->d+1;
181:     cg = cb+ctx->d+1;
182:     for (i=0;i<ds->n;i++) {
183:       A[i+(i+ds->n)*ldd] = ca[0];
184:       A[i+i*ldd] = cb[0];
185:     }
186:     for (;i<nd-ds->n;i++) {
187:       j = i/ds->n;
188:       A[i+(i+ds->n)*ldd] = ca[j];
189:       A[i+i*ldd] = cb[j];
190:       A[i+(i-ds->n)*ldd] = cg[j];
191:     }
192:     for (i=0;i<ctx->d-2;i++) {
193:       off = i*ds->n*ldd+(ctx->d-1)*ds->n;
194:       for (j=0;j<ds->n;j++)
195:         for (k=0;k<ds->n;k++)
196:           *(A+off+j*ldd+k) = *(ds->mat[DSMatExtra[i]]+j*ds->ld+k)*ca[ctx->d-1];
197:     }
198:     off = i*ds->n*ldd+(ctx->d-1)*ds->n;
199:     for (j=0;j<ds->n;j++)
200:       for (k=0;k<ds->n;k++)
201:         *(A+off+j*ldd+k) = *(ds->mat[DSMatExtra[i]]+j*ds->ld+k)*ca[ctx->d-1]-E[j*ds->ld+k]*cg[ctx->d-1];
202:     off = (++i)*ds->n*ldd+(ctx->d-1)*ds->n;
203:     for (j=0;j<ds->n;j++)
204:       for (k=0;k<ds->n;k++)
205:         *(A+off+j*ldd+k) = *(ds->mat[DSMatExtra[i]]+j*ds->ld+k)*ca[ctx->d-1]-E[j*ds->ld+k]*cb[ctx->d-1];
206:   }
207:   PetscArrayzero(B,ldd*ldd);
208:   for (i=0;i<nd-ds->n;i++) B[i+i*ldd] = 1.0;
209:   off = (ctx->d-1)*ds->n*(ldd+1);
210:   for (j=0;j<ds->n;j++) {
211:     for (i=0;i<ds->n;i++) B[off+i+j*ldd] = -E[i+j*ds->ld];
212:   }

214:   /* solve generalized eigenproblem */
215: #if defined(PETSC_USE_COMPLEX)
216:   rwork = ds->rwork;
217:   PetscStackCallBLAS("LAPACKggev",LAPACKggev_("V","V",&nd,A,&ldd,B,&ldd,wr,beta,U,&ldd,W,&ldd,work,&lwork,rwork,&info));
218: #else
219:   PetscStackCallBLAS("LAPACKggev",LAPACKggev_("V","V",&nd,A,&ldd,B,&ldd,wr,wi,beta,U,&ldd,W,&ldd,work,&lwork,&info));
220: #endif
221:   SlepcCheckLapackInfo("ggev",info);

223:   /* copy eigenvalues */
224:   for (i=0;i<nd;i++) {
225:     if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
226:     else wr[i] /= beta[i];
227: #if !defined(PETSC_USE_COMPLEX)
228:     if (beta[i]==0.0) wi[i] = 0.0;
229:     else wi[i] /= beta[i];
230: #else
231:     if (wi) wi[i] = 0.0;
232: #endif
233:   }

235:   /* copy and normalize eigenvectors */
236:   for (j=0;j<nd;j++) {
237:     PetscArraycpy(X+j*ds->ld,W+j*ldd,ds->n);
238:     PetscArraycpy(Y+j*ds->ld,U+ds->n*(ctx->d-1)+j*ldd,ds->n);
239:   }
240:   for (j=0;j<nd;j++) {
241: #if !defined(PETSC_USE_COMPLEX)
242:     if (wi[j] != 0.0) {
243:       norm = BLASnrm2_(&n,X+j*ds->ld,&one);
244:       norm0 = BLASnrm2_(&n,X+(j+1)*ds->ld,&one);
245:       norm = 1.0/SlepcAbsEigenvalue(norm,norm0);
246:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,X+j*ds->ld,&one));
247:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,X+(j+1)*ds->ld,&one));
248:       norm = BLASnrm2_(&n,Y+j*ds->ld,&one);
249:       norm0 = BLASnrm2_(&n,Y+(j+1)*ds->ld,&one);
250:       norm = 1.0/SlepcAbsEigenvalue(norm,norm0);
251:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,Y+j*ds->ld,&one));
252:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,Y+(j+1)*ds->ld,&one));
253:       j++;
254:     } else
255: #endif
256:     {
257:       norm = 1.0/BLASnrm2_(&n,X+j*ds->ld,&one);
258:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,X+j*ds->ld,&one));
259:       norm = 1.0/BLASnrm2_(&n,Y+j*ds->ld,&one);
260:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,Y+j*ds->ld,&one));
261:     }
262:   }
263:   return(0);
264: #endif
265: }

267: PetscErrorCode DSSynchronize_PEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
268: {
270:   DS_PEP         *ctx = (DS_PEP*)ds->data;
271:   PetscInt       ld=ds->ld,k=0;
272:   PetscMPIInt    ldnd,rank,off=0,size,dn;

275:   if (ds->state>=DS_STATE_CONDENSED) k += 2*ctx->d*ds->n*ld;
276:   if (eigr) k += ctx->d*ds->n;
277:   if (eigi) k += ctx->d*ds->n;
278:   DSAllocateWork_Private(ds,k,0,0);
279:   PetscMPIIntCast(k*sizeof(PetscScalar),&size);
280:   PetscMPIIntCast(ds->n*ctx->d*ld,&ldnd);
281:   PetscMPIIntCast(ctx->d*ds->n,&dn);
282:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
283:   if (!rank) {
284:     if (ds->state>=DS_STATE_CONDENSED) {
285:       MPI_Pack(ds->mat[DS_MAT_X],ldnd,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
286:       MPI_Pack(ds->mat[DS_MAT_Y],ldnd,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
287:     }
288:     if (eigr) {
289:       MPI_Pack(eigr,dn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
290:     }
291:     if (eigi) {
292:       MPI_Pack(eigi,dn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
293:     }
294:   }
295:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
296:   if (rank) {
297:     if (ds->state>=DS_STATE_CONDENSED) {
298:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_X],ldnd,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
299:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Y],ldnd,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
300:     }
301:     if (eigr) {
302:       MPI_Unpack(ds->work,size,&off,eigr,dn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
303:     }
304:     if (eigi) {
305:       MPI_Unpack(ds->work,size,&off,eigi,dn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
306:     }
307:   }
308:   return(0);
309: }

311: static PetscErrorCode DSPEPSetDegree_PEP(DS ds,PetscInt d)
312: {
313:   DS_PEP *ctx = (DS_PEP*)ds->data;

316:   if (d<0) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"The degree must be a non-negative integer");
317:   if (d>=DS_NUM_EXTRA) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Only implemented for polynomials of degree at most %D",DS_NUM_EXTRA-1);
318:   ctx->d = d;
319:   return(0);
320: }

322: /*@
323:    DSPEPSetDegree - Sets the polynomial degree for a DSPEP.

325:    Logically Collective on ds

327:    Input Parameters:
328: +  ds - the direct solver context
329: -  d  - the degree

331:    Level: intermediate

333: .seealso: DSPEPGetDegree()
334: @*/
335: PetscErrorCode DSPEPSetDegree(DS ds,PetscInt d)
336: {

342:   PetscTryMethod(ds,"DSPEPSetDegree_C",(DS,PetscInt),(ds,d));
343:   return(0);
344: }

346: static PetscErrorCode DSPEPGetDegree_PEP(DS ds,PetscInt *d)
347: {
348:   DS_PEP *ctx = (DS_PEP*)ds->data;

351:   *d = ctx->d;
352:   return(0);
353: }

355: /*@
356:    DSPEPGetDegree - Returns the polynomial degree for a DSPEP.

358:    Not collective

360:    Input Parameter:
361: .  ds - the direct solver context

363:    Output Parameters:
364: .  d - the degree

366:    Level: intermediate

368: .seealso: DSPEPSetDegree()
369: @*/
370: PetscErrorCode DSPEPGetDegree(DS ds,PetscInt *d)
371: {

377:   PetscUseMethod(ds,"DSPEPGetDegree_C",(DS,PetscInt*),(ds,d));
378:   return(0);
379: }

381: static PetscErrorCode DSPEPSetCoefficients_PEP(DS ds,PetscReal *pbc)
382: {
384:   DS_PEP         *ctx = (DS_PEP*)ds->data;
385:   PetscInt       i;

388:   if (!ctx->d) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Must first specify the polynomial degree via DSPEPSetDegree()");
389:   if (ctx->pbc) { PetscFree(ctx->pbc); }
390:   PetscMalloc1(3*(ctx->d+1),&ctx->pbc);
391:   for (i=0;i<3*(ctx->d+1);i++) ctx->pbc[i] = pbc[i];
392:   ds->state = DS_STATE_RAW;
393:   return(0);
394: }

396: /*@C
397:    DSPEPSetCoefficients - Sets the polynomial basis coefficients for a DSPEP.

399:    Logically Collective on ds

401:    Input Parameters:
402: +  ds  - the direct solver context
403: -  pbc - the polynomial basis coefficients

405:    Notes:
406:    This function is required only in the case of a polynomial specified in a
407:    non-monomial basis, to provide the coefficients that will be used
408:    during the linearization, multiplying the identity blocks on the three main
409:    diagonal blocks. Depending on the polynomial basis (Chebyshev, Legendre, ...)
410:    the coefficients must be different.

412:    There must be a total of 3*(d+1) coefficients, where d is the degree of the
413:    polynomial. The coefficients are arranged in three groups: alpha, beta, and
414:    gamma, according to the definition of the three-term recurrence. In the case
415:    of the monomial basis, alpha=1 and beta=gamma=0, in which case it is not
416:    necessary to invoke this function.

418:    Level: advanced

420: .seealso: DSPEPGetCoefficients(), DSPEPSetDegree()
421: @*/
422: PetscErrorCode DSPEPSetCoefficients(DS ds,PetscReal *pbc)
423: {

428:   PetscTryMethod(ds,"DSPEPSetCoefficients_C",(DS,PetscReal*),(ds,pbc));
429:   return(0);
430: }

432: static PetscErrorCode DSPEPGetCoefficients_PEP(DS ds,PetscReal **pbc)
433: {
435:   DS_PEP         *ctx = (DS_PEP*)ds->data;
436:   PetscInt       i;

439:   if (!ctx->d) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"Must first specify the polynomial degree via DSPEPSetDegree()");
440:   PetscCalloc1(3*(ctx->d+1),pbc);
441:   if (ctx->pbc) for (i=0;i<3*(ctx->d+1);i++) (*pbc)[i] = ctx->pbc[i];
442:   else for (i=0;i<ctx->d+1;i++) (*pbc)[i] = 1.0;
443:   return(0);
444: }

446: /*@C
447:    DSPEPGetCoefficients - Returns the polynomial basis coefficients for a DSPEP.

449:    Not collective

451:    Input Parameter:
452: .  ds - the direct solver context

454:    Output Parameters:
455: .  pbc - the polynomial basis coefficients

457:    Note:
458:    The returned array has length 3*(d+1) and should be freed by the user.

460:    Fortran Note:
461:    The calling sequence from Fortran is
462: .vb
463:    DSPEPGetCoefficients(eps,pbc,ierr)
464:    double precision pbc(d+1) output
465: .ve

467:    Level: advanced

469: .seealso: DSPEPSetCoefficients()
470: @*/
471: PetscErrorCode DSPEPGetCoefficients(DS ds,PetscReal **pbc)
472: {

478:   PetscUseMethod(ds,"DSPEPGetCoefficients_C",(DS,PetscReal**),(ds,pbc));
479:   return(0);
480: }

482: PetscErrorCode DSDestroy_PEP(DS ds)
483: {
485:   DS_PEP         *ctx = (DS_PEP*)ds->data;

488:   if (ctx->pbc) { PetscFree(ctx->pbc); }
489:   PetscFree(ds->data);
490:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetDegree_C",NULL);
491:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetDegree_C",NULL);
492:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetCoefficients_C",NULL);
493:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetCoefficients_C",NULL);
494:   return(0);
495: }

497: PetscErrorCode DSMatGetSize_PEP(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
498: {
499:   DS_PEP *ctx = (DS_PEP*)ds->data;

502:   if (!ctx->d) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"DSPEP requires specifying the polynomial degree via DSPEPSetDegree()");
503:   *rows = ds->n;
504:   if (t==DS_MAT_A || t==DS_MAT_B || t==DS_MAT_W || t==DS_MAT_U) *rows *= ctx->d;
505:   *cols = ds->n;
506:   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;
507:   return(0);
508: }

510: SLEPC_EXTERN PetscErrorCode DSCreate_PEP(DS ds)
511: {
512:   DS_PEP         *ctx;

516:   PetscNewLog(ds,&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:   ds->ops->synchronize   = DSSynchronize_PEP;
525:   ds->ops->destroy       = DSDestroy_PEP;
526:   ds->ops->matgetsize    = DSMatGetSize_PEP;
527:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetDegree_C",DSPEPSetDegree_PEP);
528:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetDegree_C",DSPEPGetDegree_PEP);
529:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetCoefficients_C",DSPEPSetCoefficients_PEP);
530:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetCoefficients_C",DSPEPGetCoefficients_PEP);
531:   return(0);
532: }