Actual source code: dspep.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: #include <slepc/private/dsimpl.h>       /*I "slepcds.h" I*/
 23: #include <slepcblaslapack.h>

 25: typedef struct {
 26:   PetscInt d;              /* polynomial degree */
 27: } DS_PEP;

 31: PetscErrorCode DSAllocate_PEP(DS ds,PetscInt ld)
 32: {
 34:   DS_PEP         *ctx = (DS_PEP*)ds->data;
 35:   PetscInt       i;

 38:   if (!ctx->d) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DSPEP requires specifying the polynomial degree via DSPEPSetDegree()");
 39:   DSAllocateMat_Private(ds,DS_MAT_X);
 40:   DSAllocateMat_Private(ds,DS_MAT_Y);
 41:   for (i=0;i<=ctx->d;i++) {
 42:     DSAllocateMat_Private(ds,DSMatExtra[i]);
 43:   }
 44:   PetscFree(ds->perm);
 45:   PetscMalloc1(ld*ctx->d,&ds->perm);
 46:   PetscLogObjectMemory((PetscObject)ds,ld*ctx->d*sizeof(PetscInt));
 47:   return(0);
 48: }

 52: PetscErrorCode DSView_PEP(DS ds,PetscViewer viewer)
 53: {
 54:   PetscErrorCode    ierr;
 55:   DS_PEP            *ctx = (DS_PEP*)ds->data;
 56:   PetscViewerFormat format;
 57:   PetscInt          i;

 60:   PetscViewerGetFormat(viewer,&format);
 61:   PetscViewerASCIIPrintf(viewer,"  polynomial degree: %D\n",ctx->d);
 62:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
 63:   for (i=0;i<=ctx->d;i++) {
 64:     DSViewMat(ds,viewer,DSMatExtra[i]);
 65:   }
 66:   if (ds->state>DS_STATE_INTERMEDIATE) {
 67:     ds->m = ctx->d*ds->n;  /* temporarily set number of columns */
 68:     DSViewMat(ds,viewer,DS_MAT_X);
 69:     ds->m = 0;
 70:   }
 71:   return(0);
 72: }

 76: PetscErrorCode DSVectors_PEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
 77: {
 79:   if (rnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
 80:   switch (mat) {
 81:     case DS_MAT_X:
 82:       break;
 83:     case DS_MAT_Y:
 84:       break;
 85:     default:
 86:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
 87:   }
 88:   return(0);
 89: }

 93: PetscErrorCode DSNormalize_PEP(DS ds,DSMatType mat,PetscInt col)
 94: {
 96:   switch (mat) {
 97:     case DS_MAT_X:
 98:       break;
 99:     case DS_MAT_Y:
100:       break;
101:     default:
102:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
103:   }
104:   return(0);
105: }

109: PetscErrorCode DSSort_PEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *kout)
110: {
112:   DS_PEP         *ctx = (DS_PEP*)ds->data;
113:   PetscInt       n,i,j,k,p,*perm,told,ld;
114:   PetscScalar    *A,*X,*Y,rtmp,rtmp2;

117:   if (!ds->sc) return(0);
118:   n = ds->n*ctx->d;
119:   A  = ds->mat[DS_MAT_A];
120:   perm = ds->perm;
121:   for (i=0;i<n;i++) perm[i] = i;
122:   told = ds->t;
123:   ds->t = n;  /* force the sorting routines to consider d*n eigenvalues */
124:   if (rr) {
125:     DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
126:   } else {
127:     DSSortEigenvalues_Private(ds,wr,wi,perm,PETSC_FALSE);
128:   }
129:   ds->t = told;  /* restore value of t */
130:   for (i=0;i<n;i++) A[i]  = wr[perm[i]];
131:   for (i=0;i<n;i++) wr[i] = A[i];
132:   for (i=0;i<n;i++) A[i]  = wi[perm[i]];
133:   for (i=0;i<n;i++) wi[i] = A[i];
134:   /* cannot use DSPermuteColumns_Private() since matrix is not square */
135:   ld = ds->ld;
136:   X  = ds->mat[DS_MAT_X];
137:   Y  = ds->mat[DS_MAT_Y];
138:   for (i=0;i<n;i++) {
139:     p = perm[i];
140:     if (p != i) {
141:       j = i + 1;
142:       while (perm[j] != i) j++;
143:       perm[j] = p; perm[i] = i;
144:       /* swap columns i and j */
145:       for (k=0;k<ds->n;k++) {
146:         rtmp  = X[k+p*ld]; X[k+p*ld] = X[k+i*ld]; X[k+i*ld] = rtmp;
147:         rtmp2 = Y[k+p*ld]; Y[k+p*ld] = Y[k+i*ld]; Y[k+i*ld] = rtmp2;
148:       }
149:     }
150:   }
151:   return(0);
152: }

156: PetscErrorCode DSSolve_PEP_QZ(DS ds,PetscScalar *wr,PetscScalar *wi)
157: {
158: #if defined(SLEPC_MISSING_LAPACK_GGEV)
160:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GGEV - Lapack routine is unavailable");
161: #else
163:   DS_PEP         *ctx = (DS_PEP*)ds->data;
164:   PetscInt       i,j,off;
165:   PetscScalar    *A,*B,*W,*X,*U,*Y,*E,*work,*beta,norm;
166:   PetscBLASInt   info,n,ldd,nd,lrwork=0,lwork,one=1;
167: #if defined(PETSC_USE_COMPLEX)
168:   PetscReal      *rwork;
169: #else
170:   PetscScalar    norm0;
171: #endif

174:   if (!ds->mat[DS_MAT_A]) {
175:     DSAllocateMat_Private(ds,DS_MAT_A);
176:   }
177:   if (!ds->mat[DS_MAT_B]) {
178:     DSAllocateMat_Private(ds,DS_MAT_B);
179:   }
180:   if (!ds->mat[DS_MAT_W]) {
181:     DSAllocateMat_Private(ds,DS_MAT_W);
182:   }
183:   if (!ds->mat[DS_MAT_U]) {
184:     DSAllocateMat_Private(ds,DS_MAT_U);
185:   }
186:   PetscBLASIntCast(ds->n*ctx->d,&nd);
187:   PetscBLASIntCast(ds->n,&n);
188:   PetscBLASIntCast(ds->ld*ctx->d,&ldd);
189: #if defined(PETSC_USE_COMPLEX)
190:   PetscBLASIntCast(nd+2*nd,&lwork);
191:   PetscBLASIntCast(8*nd,&lrwork);
192: #else
193:   PetscBLASIntCast(nd+8*nd,&lwork);
194: #endif
195:   DSAllocateWork_Private(ds,lwork,lrwork,0);
196:   beta = ds->work;
197:   work = ds->work + nd;
198:   lwork -= nd;
199:   A = ds->mat[DS_MAT_A];
200:   B = ds->mat[DS_MAT_B];
201:   W = ds->mat[DS_MAT_W];
202:   U = ds->mat[DS_MAT_U];
203:   X = ds->mat[DS_MAT_X];
204:   Y = ds->mat[DS_MAT_Y];
205:   E = ds->mat[DSMatExtra[ctx->d]];

207:   /* build matrices A and B of the linearization */
208:   PetscMemzero(A,ldd*ldd*sizeof(PetscScalar));
209:   for (i=0;i<nd-ds->n;i++) A[i+(i+ds->n)*ldd] = -1.0;
210:   for (i=0;i<ctx->d;i++) {
211:     off = i*ds->n*ldd+(ctx->d-1)*ds->n;
212:     for (j=0;j<ds->n;j++) {
213:       PetscMemcpy(A+off+j*ldd,ds->mat[DSMatExtra[i]]+j*ds->ld,ds->n*sizeof(PetscScalar));
214:     }
215:   }
216:   PetscMemzero(B,ldd*ldd*sizeof(PetscScalar));
217:   for (i=0;i<nd-ds->n;i++) B[i+i*ldd] = -1.0;
218:   off = (ctx->d-1)*ds->n*(ldd+1);
219:   for (j=0;j<ds->n;j++) {
220:     for (i=0;i<ds->n;i++) B[off+i+j*ldd] = -E[i+j*ds->ld];
221:   }

223:   /* solve generalized eigenproblem */
224: #if defined(PETSC_USE_COMPLEX)
225:   rwork = ds->rwork;
226:   PetscStackCallBLAS("LAPACKggev",LAPACKggev_("V","V",&nd,A,&ldd,B,&ldd,wr,beta,U,&ldd,W,&ldd,work,&lwork,rwork,&info));
227:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack ZGGEV %d",info);
228: #else
229:   PetscStackCallBLAS("LAPACKggev",LAPACKggev_("V","V",&nd,A,&ldd,B,&ldd,wr,wi,beta,U,&ldd,W,&ldd,work,&lwork,&info));
230:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DGGEV %d",info);
231: #endif

233:   /* copy eigenvalues */
234:   for (i=0;i<nd;i++) {
235:     if (beta[i]==0.0) wr[i] = (PetscRealPart(wr[i])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
236:     else wr[i] /= beta[i];
237: #if !defined(PETSC_USE_COMPLEX)
238:     if (beta[i]==0.0) wi[i] = (wi[i]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
239:     else wi[i] /= beta[i];
240: #else
241:     if (wi) wi[i] = 0.0;
242: #endif
243:   }

245:   /* copy and normalize eigenvectors */
246:   for (j=0;j<nd;j++) {
247:     PetscMemcpy(X+j*ds->ld,W+j*ldd,ds->n*sizeof(PetscScalar));
248:     PetscMemcpy(Y+j*ds->ld,U+ds->n*(ctx->d-1)+j*ldd,ds->n*sizeof(PetscScalar));
249:   }
250:   for (j=0;j<nd;j++) {
251: #if !defined(PETSC_USE_COMPLEX)
252:     if (wi[j] != 0.0) {
253:       norm = BLASnrm2_(&n,X+j*ds->ld,&one);
254:       norm0 = BLASnrm2_(&n,X+(j+1)*ds->ld,&one);
255:       norm = 1.0/SlepcAbsEigenvalue(norm,norm0);
256:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,X+j*ds->ld,&one));
257:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,X+(j+1)*ds->ld,&one));
258:       norm = BLASnrm2_(&n,Y+j*ds->ld,&one);
259:       norm0 = BLASnrm2_(&n,Y+(j+1)*ds->ld,&one);
260:       norm = 1.0/SlepcAbsEigenvalue(norm,norm0);
261:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,Y+j*ds->ld,&one));
262:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,Y+(j+1)*ds->ld,&one));
263:       j++;
264:     } else
265: #endif
266:     {
267:       norm = 1.0/BLASnrm2_(&n,X+j*ds->ld,&one);
268:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,X+j*ds->ld,&one));
269:       norm = 1.0/BLASnrm2_(&n,Y+j*ds->ld,&one);
270:       PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,Y+j*ds->ld,&one));
271:     }
272:   }
273:   return(0);
274: #endif
275: }

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

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

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: }

318: static PetscErrorCode DSPEPGetDegree_PEP(DS ds,PetscInt *d)
319: {
320:   DS_PEP *ctx = (DS_PEP*)ds->data;

323:   *d = ctx->d;
324:   return(0);
325: }

329: /*@
330:    DSPEPGetDegree - Returns the polynomial degree for a DSPEP.

332:    Not collective

334:    Input Parameter:
335: .  ds - the direct solver context

337:    Output Parameters:
338: .  d - the degree

340:    Level: intermediate

342: .seealso: DSPEPSetDegree()
343: @*/
344: PetscErrorCode DSPEPGetDegree(DS ds,PetscInt *d)
345: {

351:   PetscUseMethod(ds,"DSPEPGetDegree_C",(DS,PetscInt*),(ds,d));
352:   return(0);
353: }

357: PetscErrorCode DSDestroy_PEP(DS ds)
358: {

362:   PetscFree(ds->data);
363:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetDegree_C",NULL);
364:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetDegree_C",NULL);
365:   return(0);
366: }

370: PETSC_EXTERN PetscErrorCode DSCreate_PEP(DS ds)
371: {
372:   DS_PEP         *ctx;

376:   PetscNewLog(ds,&ctx);
377:   ds->data = (void*)ctx;

379:   ds->ops->allocate      = DSAllocate_PEP;
380:   ds->ops->view          = DSView_PEP;
381:   ds->ops->vectors       = DSVectors_PEP;
382:   ds->ops->solve[0]      = DSSolve_PEP_QZ;
383:   ds->ops->sort          = DSSort_PEP;
384:   ds->ops->normalize     = DSNormalize_PEP;
385:   ds->ops->destroy       = DSDestroy_PEP;
386:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetDegree_C",DSPEPSetDegree_PEP);
387:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetDegree_C",DSPEPGetDegree_PEP);
388:   return(0);
389: }