Actual source code: dsnep.c

slepc-3.13.1 2020-04-12
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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       nf;                 /* number of functions in f[] */
 16:   FN             f[DS_NUM_EXTRA];    /* functions defining the nonlinear operator */
 17:   PetscInt       neig;               /* number of available eigenpairs */
 18:   void           *computematrixctx;
 19:   PetscErrorCode (*computematrix)(DS,PetscScalar,PetscBool,DSMatType,void*);
 20: } DS_NEP;

 22: /*
 23:    DSNEPComputeMatrix - Build the matrix associated with a nonlinear operator
 24:    T(lambda) or its derivative T'(lambda), given the parameter lambda, where
 25:    T(lambda) = sum_i E_i*f_i(lambda). The result is written in mat.
 26: */
 27: static PetscErrorCode DSNEPComputeMatrix(DS ds,PetscScalar lambda,PetscBool deriv,DSMatType mat)
 28: {
 30:   DS_NEP         *ctx = (DS_NEP*)ds->data;
 31:   PetscScalar    *T,*E,alpha;
 32:   PetscInt       i,ld,n;
 33:   PetscBLASInt   k,inc=1;

 36:   PetscLogEventBegin(DS_Other,ds,0,0,0);
 37:   if (ctx->computematrix) {
 38:     (*ctx->computematrix)(ds,lambda,deriv,mat,ctx->computematrixctx);
 39:   } else {
 40:     DSGetDimensions(ds,&n,NULL,NULL,NULL,NULL);
 41:     DSGetLeadingDimension(ds,&ld);
 42:     PetscBLASIntCast(ld*n,&k);
 43:     DSGetArray(ds,mat,&T);
 44:     PetscArrayzero(T,k);
 45:     for (i=0;i<ctx->nf;i++) {
 46:       if (deriv) {
 47:         FNEvaluateDerivative(ctx->f[i],lambda,&alpha);
 48:       } else {
 49:         FNEvaluateFunction(ctx->f[i],lambda,&alpha);
 50:       }
 51:       E = ds->mat[DSMatExtra[i]];
 52:       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k,&alpha,E,&inc,T,&inc));
 53:     }
 54:     DSRestoreArray(ds,mat,&T);
 55:   }
 56:   PetscLogEventEnd(DS_Other,ds,0,0,0);
 57:   return(0);
 58: }

 60: PetscErrorCode DSAllocate_NEP(DS ds,PetscInt ld)
 61: {
 63:   DS_NEP         *ctx = (DS_NEP*)ds->data;
 64:   PetscInt       i;

 67:   if (!ctx->nf) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_WRONGSTATE,"DSNEP requires passing some functions via DSSetFN()");
 68:   DSAllocateMat_Private(ds,DS_MAT_X);
 69:   for (i=0;i<ctx->nf;i++) {
 70:     DSAllocateMat_Private(ds,DSMatExtra[i]);
 71:   }
 72:   PetscFree(ds->perm);
 73:   PetscMalloc1(ld,&ds->perm);
 74:   PetscLogObjectMemory((PetscObject)ds,ld*sizeof(PetscInt));
 75:   return(0);
 76: }

 78: PetscErrorCode DSView_NEP(DS ds,PetscViewer viewer)
 79: {
 80:   PetscErrorCode    ierr;
 81:   DS_NEP            *ctx = (DS_NEP*)ds->data;
 82:   PetscViewerFormat format;
 83:   PetscInt          i;

 86:   PetscViewerGetFormat(viewer,&format);
 87:   PetscViewerASCIIPrintf(viewer,"number of functions: %D\n",ctx->nf);
 88:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);
 89:   for (i=0;i<ctx->nf;i++) {
 90:     FNView(ctx->f[i],viewer);
 91:     DSViewMat(ds,viewer,DSMatExtra[i]);
 92:   }
 93:   if (ds->state>DS_STATE_INTERMEDIATE) { DSViewMat(ds,viewer,DS_MAT_X); }
 94:   return(0);
 95: }

 97: PetscErrorCode DSVectors_NEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
 98: {
100:   if (rnorm) SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
101:   switch (mat) {
102:     case DS_MAT_X:
103:       break;
104:     case DS_MAT_Y:
105:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
106:       break;
107:     default:
108:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
109:   }
110:   return(0);
111: }

113: PetscErrorCode DSSort_NEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *dummy)
114: {
116:   DS_NEP         *ctx = (DS_NEP*)ds->data;
117:   PetscInt       n,l,i,j,k,p,*perm,told,ld=ds->ld;
118:   PetscScalar    *A,*X,rtmp;

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

154: PetscErrorCode DSSolve_NEP_SLP(DS ds,PetscScalar *wr,PetscScalar *wi)
155: {
157:   DS_NEP         *ctx = (DS_NEP*)ds->data;
158:   PetscScalar    *A,*B,*W,*X,*work,*alpha,*beta;
159:   PetscScalar    norm,sigma,lambda,mu,re,re2,sone=1.0,zero=0.0;
160:   PetscBLASInt   info,n,ld,lrwork=0,lwork,one=1;
161:   PetscInt       it,pos,j,maxit=100,result;
162:   PetscReal      tol;
163: #if defined(PETSC_USE_COMPLEX)
164:   PetscReal      *rwork;
165: #else
166:   PetscReal      *alphai,im,im2;
167: #endif

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

203:   sigma = 0.0;
204:   if (ds->sc->comparison==SlepcCompareTargetMagnitude || ds->sc->comparison==SlepcCompareTargetReal) sigma = *(PetscScalar*)ds->sc->comparisonctx;
205:   lambda = sigma;
206:   tol = 1000*n*PETSC_MACHINE_EPSILON;

208:   for (it=0;it<maxit;it++) {

210:     /* evaluate T and T' */
211:     DSNEPComputeMatrix(ds,lambda,PETSC_FALSE,DS_MAT_A);
212:     if (it) {
213:       PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,A,&ld,X,&one,&zero,X+ld,&one));
214:       norm = BLASnrm2_(&n,X+ld,&one);
215:       if (PetscRealPart(norm)/PetscAbsScalar(lambda)<=tol) break;
216:     }
217:     DSNEPComputeMatrix(ds,lambda,PETSC_TRUE,DS_MAT_B);

219:     /* compute eigenvalue correction mu and eigenvector u */
220: #if defined(PETSC_USE_COMPLEX)
221:     rwork = ds->rwork;
222:     PetscStackCallBLAS("LAPACKggev",LAPACKggev_("N","V",&n,A,&ld,B,&ld,alpha,beta,NULL,&ld,W,&ld,work,&lwork,rwork,&info));
223: #else
224:     PetscStackCallBLAS("LAPACKggev",LAPACKggev_("N","V",&n,A,&ld,B,&ld,alpha,alphai,beta,NULL,&ld,W,&ld,work,&lwork,&info));
225: #endif
226:     SlepcCheckLapackInfo("ggev",info);

228:     /* find smallest eigenvalue */
229:     j = 0;
230:     if (beta[j]==0.0) re = (PetscRealPart(alpha[j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
231:     else re = alpha[j]/beta[j];
232: #if !defined(PETSC_USE_COMPLEX)
233:     if (beta[j]==0.0) im = (alphai[j]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
234:     else im = alphai[j]/beta[j];
235: #endif
236:     pos = 0;
237:     for (j=1;j<n;j++) {
238:       if (beta[j]==0.0) re2 = (PetscRealPart(alpha[j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
239:       else re2 = alpha[j]/beta[j];
240: #if !defined(PETSC_USE_COMPLEX)
241:       if (beta[j]==0.0) im2 = (alphai[j]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
242:       else im2 = alphai[j]/beta[j];
243:       SlepcCompareSmallestMagnitude(re,im,re2,im2,&result,NULL);
244: #else
245:       SlepcCompareSmallestMagnitude(re,0.0,re2,0.0,&result,NULL);
246: #endif
247:       if (result > 0) {
248:         re = re2;
249: #if !defined(PETSC_USE_COMPLEX)
250:         im = im2;
251: #endif
252:         pos = j;
253:       }
254:     }

256: #if !defined(PETSC_USE_COMPLEX)
257:     if (im!=0.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DSNEP found a complex eigenvalue; try rerunning with complex scalars");
258: #endif
259:     mu = alpha[pos]/beta[pos];
260:     PetscArraycpy(X,W+pos*ld,n);
261:     norm = BLASnrm2_(&n,X,&one);
262:     norm = 1.0/norm;
263:     PetscStackCallBLAS("BLASscal",BLASscal_(&n,&norm,X,&one));

265:     /* correct eigenvalue approximation */
266:     lambda = lambda - mu;
267:   }

269:   if (it==maxit) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"DSNEP did not converge");
270:   ctx->neig = 1;
271:   wr[0] = lambda;
272:   if (wi) wi[0] = 0.0;
273:   return(0);
274: }

276: PetscErrorCode DSSynchronize_NEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
277: {
279:   PetscInt       k=0;
280:   PetscMPIInt    n,rank,size,off=0;

283:   if (ds->state>=DS_STATE_CONDENSED) k += ds->n;
284:   if (eigr) k += 1;
285:   if (eigi) k += 1;
286:   DSAllocateWork_Private(ds,k,0,0);
287:   PetscMPIIntCast(k*sizeof(PetscScalar),&size);
288:   PetscMPIIntCast(ds->n,&n);
289:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
290:   if (!rank) {
291:     if (ds->state>=DS_STATE_CONDENSED) {
292:       MPI_Pack(ds->mat[DS_MAT_X],n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
293:     }
294:     if (eigr) {
295:       MPI_Pack(eigr,1,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
296:     }
297:     if (eigi) {
298:       MPI_Pack(eigi,1,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
299:     }
300:   }
301:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
302:   if (rank) {
303:     if (ds->state>=DS_STATE_CONDENSED) {
304:       MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_X],n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
305:     }
306:     if (eigr) {
307:       MPI_Unpack(ds->work,size,&off,eigr,1,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
308:     }
309:     if (eigi) {
310:       MPI_Unpack(ds->work,size,&off,eigi,1,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
311:     }
312:   }
313:   return(0);
314: }

316: static PetscErrorCode DSNEPSetFN_NEP(DS ds,PetscInt n,FN fn[])
317: {
319:   DS_NEP         *ctx = (DS_NEP*)ds->data;
320:   PetscInt       i;

323:   if (n<=0) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more functions, you have %D",n);
324:   if (n>DS_NUM_EXTRA) SETERRQ2(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Too many functions, you specified %D but the limit is %D",n,DS_NUM_EXTRA);
325:   if (ds->ld) { PetscInfo(ds,"DSNEPSetFN() called after DSAllocate()\n"); }
326:   for (i=0;i<n;i++) {
327:     PetscObjectReference((PetscObject)fn[i]);
328:   }
329:   for (i=0;i<ctx->nf;i++) {
330:     FNDestroy(&ctx->f[i]);
331:   }
332:   for (i=0;i<n;i++) ctx->f[i] = fn[i];
333:   ctx->nf = n;
334:   return(0);
335: }

337: /*@
338:    DSNEPSetFN - Sets a number of functions that define the nonlinear
339:    eigenproblem.

341:    Collective on ds

343:    Input Parameters:
344: +  ds - the direct solver context
345: .  n  - number of functions
346: -  fn - array of functions

348:    Notes:
349:    The nonlinear eigenproblem is defined in terms of the split nonlinear
350:    operator T(lambda) = sum_i A_i*f_i(lambda).

352:    This function must be called before DSAllocate(). Then DSAllocate()
353:    will allocate an extra matrix A_i per each function, that can be
354:    filled in the usual way.

356:    Level: advanced

358: .seealso: DSNEPGetFN(), DSAllocate()
359:  @*/
360: PetscErrorCode DSNEPSetFN(DS ds,PetscInt n,FN fn[])
361: {
362:   PetscInt       i;

369:   for (i=0;i<n;i++) {
372:   }
373:   PetscTryMethod(ds,"DSNEPSetFN_C",(DS,PetscInt,FN[]),(ds,n,fn));
374:   return(0);
375: }

377: static PetscErrorCode DSNEPGetFN_NEP(DS ds,PetscInt k,FN *fn)
378: {
379:   DS_NEP *ctx = (DS_NEP*)ds->data;

382:   if (k<0 || k>=ctx->nf) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %D",ctx->nf-1);
383:   *fn = ctx->f[k];
384:   return(0);
385: }

387: /*@
388:    DSNEPGetFN - Gets the functions associated with the nonlinear DS.

390:    Not collective, though parallel FNs are returned if the DS is parallel

392:    Input Parameter:
393: +  ds - the direct solver context
394: -  k  - the index of the requested function (starting in 0)

396:    Output Parameter:
397: .  fn - the function

399:    Level: advanced

401: .seealso: DSNEPSetFN()
402: @*/
403: PetscErrorCode DSNEPGetFN(DS ds,PetscInt k,FN *fn)
404: {

410:   PetscUseMethod(ds,"DSNEPGetFN_C",(DS,PetscInt,FN*),(ds,k,fn));
411:   return(0);
412: }

414: static PetscErrorCode DSNEPGetNumFN_NEP(DS ds,PetscInt *n)
415: {
416:   DS_NEP *ctx = (DS_NEP*)ds->data;

419:   *n = ctx->nf;
420:   return(0);
421: }

423: /*@
424:    DSNEPGetNumFN - Returns the number of functions stored internally by
425:    the DS.

427:    Not collective

429:    Input Parameter:
430: .  ds - the direct solver context

432:    Output Parameters:
433: .  n - the number of functions passed in DSNEPSetFN()

435:    Level: advanced

437: .seealso: DSNEPSetFN()
438: @*/
439: PetscErrorCode DSNEPGetNumFN(DS ds,PetscInt *n)
440: {

446:   PetscUseMethod(ds,"DSNEPGetNumFN_C",(DS,PetscInt*),(ds,n));
447:   return(0);
448: }

450: static PetscErrorCode DSNEPSetComputeMatrixFunction_NEP(DS ds,PetscErrorCode (*fun)(DS,PetscScalar,PetscBool,DSMatType,void*),void *ctx)
451: {
452:   DS_NEP *dsctx = (DS_NEP*)ds->data;

455:   dsctx->computematrix    = fun;
456:   dsctx->computematrixctx = ctx;
457:   return(0);
458: }

460: /*@
461:    DSNEPSetComputeMatrixFunction - Sets a user-provided subroutine to compute
462:    the matrices T(lambda) or T'(lambda).

464:    Logically Collective on ds

466:    Input Parameters:
467: +  ds  - the direct solver context
468: .  fun - a pointer to the user function
469: -  ctx - a context pointer (the last parameter to the user function)

471:    Calling Sequence of fun:
472: $   fun(DS ds,PetscScalar lambda,PetscBool deriv,DSMatType mat,void *ctx)

474: +   ds     - the direct solver object
475: .   lambda - point where T(lambda) or T'(lambda) must be evaluated
476: .   deriv  - if true compute T'(lambda), otherwise compute T(lambda)
477: .   mat    - the DS matrix where the result must be stored
478: -   ctx    - optional context, as set by DSNEPSetComputeMatrixFunction()

480:    Note:
481:    The result is computed as T(lambda) = sum_i E_i*f_i(lambda), and similarly
482:    for the derivative.

484:    Level: developer
485: @*/
486: PetscErrorCode DSNEPSetComputeMatrixFunction(DS ds,PetscErrorCode (*fun)(DS,PetscScalar,PetscBool,DSMatType,void*),void *ctx)
487: {

492:   PetscTryMethod(ds,"DSNEPSetComputeMatrixFunction_C",(DS,PetscErrorCode (*)(DS,PetscScalar,PetscBool,DSMatType,void*),void*),(ds,fun,ctx));
493:   return(0);
494: }

496: PetscErrorCode DSDestroy_NEP(DS ds)
497: {
499:   DS_NEP         *ctx = (DS_NEP*)ds->data;
500:   PetscInt       i;

503:   for (i=0;i<ctx->nf;i++) {
504:     FNDestroy(&ctx->f[i]);
505:   }
506:   PetscFree(ds->data);
507:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetFN_C",NULL);
508:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetFN_C",NULL);
509:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetNumFN_C",NULL);
510:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetComputeMatrixFunction_C",NULL);
511:   return(0);
512: }

514: SLEPC_EXTERN PetscErrorCode DSCreate_NEP(DS ds)
515: {
516:   DS_NEP         *ctx;

520:   PetscNewLog(ds,&ctx);
521:   ds->data = (void*)ctx;

523:   ds->ops->allocate      = DSAllocate_NEP;
524:   ds->ops->view          = DSView_NEP;
525:   ds->ops->vectors       = DSVectors_NEP;
526:   ds->ops->solve[0]      = DSSolve_NEP_SLP;
527:   ds->ops->sort          = DSSort_NEP;
528:   ds->ops->synchronize   = DSSynchronize_NEP;
529:   ds->ops->destroy       = DSDestroy_NEP;
530:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetFN_C",DSNEPSetFN_NEP);
531:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetFN_C",DSNEPGetFN_NEP);
532:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetNumFN_C",DSNEPGetNumFN_NEP);
533:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetComputeMatrixFunction_C",DSNEPSetComputeMatrixFunction_NEP);
534:   return(0);
535: }