Actual source code: dspep.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2015, 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:   for (i=0;i<=ctx->d;i++) {
 41:     DSAllocateMat_Private(ds,DSMatExtra[i]);
 42:   }
 43:   PetscFree(ds->perm);
 44:   PetscMalloc1(ld*ctx->d,&ds->perm);
 45:   PetscLogObjectMemory((PetscObject)ds,ld*ctx->d*sizeof(PetscInt));
 46:   return(0);
 47: }

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

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

 75: PetscErrorCode DSVectors_PEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
 76: {
 78:   if (rnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
 79:   switch (mat) {
 80:     case DS_MAT_X:
 81:       break;
 82:     case DS_MAT_Y:
 83:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
 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:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented yet");
101:       break;
102:     default:
103:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
104:   }
105:   return(0);
106: }

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

118:   if (!ds->sc) return(0);
119:   n = ds->n*ctx->d;
120:   A  = ds->mat[DS_MAT_A];
121:   perm = ds->perm;
122:   for (i=0;i<n;i++) perm[i] = i;
123:   told = ds->t;
124:   ds->t = n;  /* force the sorting routines to consider d*n eigenvalues */
125:   if (rr) {
126:     DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
127:   } else {
128:     DSSortEigenvalues_Private(ds,wr,wi,perm,PETSC_FALSE);
129:   }
130:   ds->t = told;  /* restore value of t */
131:   for (i=0;i<n;i++) A[i]  = wr[perm[i]];
132:   for (i=0;i<n;i++) wr[i] = A[i];
133:   for (i=0;i<n;i++) A[i]  = wi[perm[i]];
134:   for (i=0;i<n;i++) wi[i] = A[i];
135:   /* cannot use DSPermuteColumns_Private() since matrix is not square */
136:   ld = ds->ld;
137:   X  = ds->mat[DS_MAT_X];
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:       }
148:     }
149:   }
150:   return(0);
151: }

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

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

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

217:   /* solve generalized eigenproblem */
218: #if defined(PETSC_USE_COMPLEX)
219:   rwork = ds->rwork;
220:   PetscStackCallBLAS("LAPACKggev",LAPACKggev_("N","V",&nd,A,&ldd,B,&ldd,wr,beta,NULL,&ldd,W,&ldd,work,&lwork,rwork,&info));
221:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack ZGGEV %d",info);
222: #else
223:   PetscStackCallBLAS("LAPACKggev",LAPACKggev_("N","V",&nd,A,&ldd,B,&ldd,wr,wi,beta,NULL,&ldd,W,&ldd,work,&lwork,&info));
224:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DGGEV %d",info);
225: #endif

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

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

265: static PetscErrorCode DSPEPSetDegree_PEP(DS ds,PetscInt d)
266: {
267:   DS_PEP *ctx = (DS_PEP*)ds->data;

270:   if (d<0) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"The degree must be a non-negative integer");
271:   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);
272:   ctx->d = d;
273:   return(0);
274: }

278: /*@
279:    DSPEPSetDegree - Sets the polynomial degree for a DSPEP.

281:    Logically Collective on DS

283:    Input Parameters:
284: +  ds - the direct solver context
285: -  d  - the degree

287:    Level: intermediate

289: .seealso: DSPEPGetDegree()
290: @*/
291: PetscErrorCode DSPEPSetDegree(DS ds,PetscInt d)
292: {

298:   PetscTryMethod(ds,"DSPEPSetDegree_C",(DS,PetscInt),(ds,d));
299:   return(0);
300: }

304: static PetscErrorCode DSPEPGetDegree_PEP(DS ds,PetscInt *d)
305: {
306:   DS_PEP *ctx = (DS_PEP*)ds->data;

309:   *d = ctx->d;
310:   return(0);
311: }

315: /*@
316:    DSPEPGetDegree - Returns the polynomial degree for a DSPEP.

318:    Not collective

320:    Input Parameter:
321: .  ds - the direct solver context

323:    Output Parameters:
324: .  d - the degree

326:    Level: intermediate

328: .seealso: DSPEPSetDegree()
329: @*/
330: PetscErrorCode DSPEPGetDegree(DS ds,PetscInt *d)
331: {

337:   PetscTryMethod(ds,"DSPEPGetDegree_C",(DS,PetscInt*),(ds,d));
338:   return(0);
339: }

343: PetscErrorCode DSDestroy_PEP(DS ds)
344: {

348:   PetscFree(ds->data);
349:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetDegree_C",NULL);
350:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetDegree_C",NULL);
351:   return(0);
352: }

356: PETSC_EXTERN PetscErrorCode DSCreate_PEP(DS ds)
357: {
358:   DS_PEP         *ctx;

362:   PetscNewLog(ds,&ctx);
363:   ds->data = (void*)ctx;

365:   ds->ops->allocate      = DSAllocate_PEP;
366:   ds->ops->view          = DSView_PEP;
367:   ds->ops->vectors       = DSVectors_PEP;
368:   ds->ops->solve[0]      = DSSolve_PEP_QZ;
369:   ds->ops->sort          = DSSort_PEP;
370:   ds->ops->normalize     = DSNormalize_PEP;
371:   ds->ops->destroy       = DSDestroy_PEP;
372:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPSetDegree_C",DSPEPSetDegree_PEP);
373:   PetscObjectComposeFunction((PetscObject)ds,"DSPEPGetDegree_C",DSPEPGetDegree_PEP);
374:   return(0);
375: }