Actual source code: dspep.c

slepc-3.8.0 2017-10-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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: } DS_PEP;

 18: PetscErrorCode DSAllocate_PEP(DS ds,PetscInt ld)
 19: {
 21:   DS_PEP         *ctx = (DS_PEP*)ds->data;
 22:   PetscInt       i;

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

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

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

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

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

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

117: PetscErrorCode DSSolve_PEP_QZ(DS ds,PetscScalar *wr,PetscScalar *wi)
118: {
119: #if defined(SLEPC_MISSING_LAPACK_GGEV)
121:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GGEV - Lapack routine is unavailable");
122: #else
124:   DS_PEP         *ctx = (DS_PEP*)ds->data;
125:   PetscInt       i,j,off;
126:   PetscScalar    *A,*B,*W,*X,*U,*Y,*E,*work,*beta,norm;
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:   PetscMemzero(A,ldd*ldd*sizeof(PetscScalar));
170:   for (i=0;i<nd-ds->n;i++) A[i+(i+ds->n)*ldd] = -1.0;
171:   for (i=0;i<ctx->d;i++) {
172:     off = i*ds->n*ldd+(ctx->d-1)*ds->n;
173:     for (j=0;j<ds->n;j++) {
174:       PetscMemcpy(A+off+j*ldd,ds->mat[DSMatExtra[i]]+j*ds->ld,ds->n*sizeof(PetscScalar));
175:     }
176:   }
177:   PetscMemzero(B,ldd*ldd*sizeof(PetscScalar));
178:   for (i=0;i<nd-ds->n;i++) B[i+i*ldd] = -1.0;
179:   off = (ctx->d-1)*ds->n*(ldd+1);
180:   for (j=0;j<ds->n;j++) {
181:     for (i=0;i<ds->n;i++) B[off+i+j*ldd] = -E[i+j*ds->ld];
182:   }

184:   /* solve generalized eigenproblem */
185: #if defined(PETSC_USE_COMPLEX)
186:   rwork = ds->rwork;
187:   PetscStackCallBLAS("LAPACKggev",LAPACKggev_("V","V",&nd,A,&ldd,B,&ldd,wr,beta,U,&ldd,W,&ldd,work,&lwork,rwork,&info));
188: #else
189:   PetscStackCallBLAS("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);

193:   /* copy eigenvalues */
194:   for (i=0;i<nd;i++) {
195:     if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
196:     else wr[i] /= beta[i];
197: #if !defined(PETSC_USE_COMPLEX)
198:     if (beta[i]==0.0) wi[i] = 0.0;
199:     else wi[i] /= beta[i];
200: #else
201:     if (wi) wi[i] = 0.0;
202: #endif
203:   }

205:   /* copy and normalize eigenvectors */
206:   for (j=0;j<nd;j++) {
207:     PetscMemcpy(X+j*ds->ld,W+j*ldd,ds->n*sizeof(PetscScalar));
208:     PetscMemcpy(Y+j*ds->ld,U+ds->n*(ctx->d-1)+j*ldd,ds->n*sizeof(PetscScalar));
209:   }
210:   for (j=0;j<nd;j++) {
211: #if !defined(PETSC_USE_COMPLEX)
212:     if (wi[j] != 0.0) {
213:       norm = BLASnrm2_(&n,X+j*ds->ld,&one);
214:       norm0 = BLASnrm2_(&n,X+(j+1)*ds->ld,&one);
215:       norm = 1.0/SlepcAbsEigenvalue(norm,norm0);
216:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,X+j*ds->ld,&one));
217:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,X+(j+1)*ds->ld,&one));
218:       norm = BLASnrm2_(&n,Y+j*ds->ld,&one);
219:       norm0 = BLASnrm2_(&n,Y+(j+1)*ds->ld,&one);
220:       norm = 1.0/SlepcAbsEigenvalue(norm,norm0);
221:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,Y+j*ds->ld,&one));
222:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,Y+(j+1)*ds->ld,&one));
223:       j++;
224:     } else
225: #endif
226:     {
227:       norm = 1.0/BLASnrm2_(&n,X+j*ds->ld,&one);
228:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,X+j*ds->ld,&one));
229:       norm = 1.0/BLASnrm2_(&n,Y+j*ds->ld,&one);
230:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,Y+j*ds->ld,&one));
231:     }
232:   }
233:   return(0);
234: #endif
235: }

237: PetscErrorCode DSSynchronize_PEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
238: {
240:   DS_PEP         *ctx = (DS_PEP*)ds->data;  
241:   PetscInt       ld=ds->ld,k=0;
242:   PetscMPIInt    ldnd,rank,off=0,size,dn;

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:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
253:   if (!rank) {
254:     if (ds->state>=DS_STATE_CONDENSED) {
255:       MPI_Pack(ds->mat[DS_MAT_X],ldnd,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
256:       MPI_Pack(ds->mat[DS_MAT_Y],ldnd,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
257:     }
258:     if (eigr) {
259:       MPI_Pack(eigr,dn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
260:     }
261:     if (eigi) {
262:       MPI_Pack(eigi,dn,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
263:     }
264:   }
265:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
266:   if (rank) {
267:     if (ds->state>=DS_STATE_CONDENSED) {
268:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_X],ldnd,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
269:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_Y],ldnd,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
270:     }
271:     if (eigr) {
272:       MPI_Unpack(ds->work,size,&off,eigr,dn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
273:     }
274:     if (eigi) {
275:       MPI_Unpack(ds->work,size,&off,eigi,dn,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
276:     }
277:   }
278:   return(0);
279: }

281: static PetscErrorCode DSPEPSetDegree_PEP(DS ds,PetscInt d)
282: {
283:   DS_PEP *ctx = (DS_PEP*)ds->data;

286:   if (d<0) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"The degree must be a non-negative integer");
287:   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);
288:   ctx->d = d;
289:   return(0);
290: }

292: /*@
293:    DSPEPSetDegree - Sets the polynomial degree for a DSPEP.

295:    Logically Collective on DS

297:    Input Parameters:
298: +  ds - the direct solver context
299: -  d  - the degree

301:    Level: intermediate

303: .seealso: DSPEPGetDegree()
304: @*/
305: PetscErrorCode DSPEPSetDegree(DS ds,PetscInt d)
306: {

312:   PetscTryMethod(ds,"DSPEPSetDegree_C",(DS,PetscInt),(ds,d));
313:   return(0);
314: }

316: static PetscErrorCode DSPEPGetDegree_PEP(DS ds,PetscInt *d)
317: {
318:   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: {

347:   PetscUseMethod(ds,"DSPEPGetDegree_C",(DS,PetscInt*),(ds,d));
348:   return(0);
349: }

351: PetscErrorCode DSDestroy_PEP(DS ds)
352: {

356:   PetscFree(ds->data);
357:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetDegree_C",NULL);
358:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetDegree_C",NULL);
359:   return(0);
360: }

362: PetscErrorCode DSMatGetSize_PEP(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
363: {
364:   DS_PEP *ctx = (DS_PEP*)ds->data;

367:   if (!ctx->d) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"DSPEP requires specifying the polynomial degree via DSPEPSetDegree()");
368:   *rows = ds->n;
369:   if (t==DS_MAT_A || t==DS_MAT_B || t==DS_MAT_W || t==DS_MAT_U) *rows *= ctx->d;
370:   *cols = ds->n;
371:   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;
372:   return(0);
373: }

375: PETSC_EXTERN PetscErrorCode DSCreate_PEP(DS ds)
376: {
377:   DS_PEP         *ctx;

381:   PetscNewLog(ds,&ctx);
382:   ds->data = (void*)ctx;

384:   ds->ops->allocate      = DSAllocate_PEP;
385:   ds->ops->view          = DSView_PEP;
386:   ds->ops->vectors       = DSVectors_PEP;
387:   ds->ops->solve[0]      = DSSolve_PEP_QZ;
388:   ds->ops->sort          = DSSort_PEP;
389:   ds->ops->synchronize   = DSSynchronize_PEP;
390:   ds->ops->destroy       = DSDestroy_PEP;
391:   ds->ops->matgetsize    = DSMatGetSize_PEP;
392:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetDegree_C",DSPEPSetDegree_PEP);
393:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetDegree_C",DSPEPGetDegree_PEP);
394:   return(0);
395: }