Actual source code: dsnep.c

slepc-3.14.1 2020-12-08
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>
 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:     default:
107:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
108:   }
109:   return(0);
110: }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

322:   if (n<=0) SETERRQ1(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more functions, you have %D",n);
323:   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);
324:   if (ds->ld) { PetscInfo(ds,"DSNEPSetFN() called after DSAllocate()\n"); }
325:   for (i=0;i<n;i++) {
326:     PetscObjectReference((PetscObject)fn[i]);
327:   }
328:   for (i=0;i<ctx->nf;i++) {
329:     FNDestroy(&ctx->f[i]);
330:   }
331:   for (i=0;i<n;i++) ctx->f[i] = fn[i];
332:   ctx->nf = n;
333:   return(0);
334: }

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

340:    Collective on ds

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

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

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

355:    Level: advanced

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

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

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

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

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

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

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

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

398:    Level: advanced

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

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

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

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

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

426:    Not collective

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

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

434:    Level: advanced

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

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

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

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

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

463:    Logically Collective on ds

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

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

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

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

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

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

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

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

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

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

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