Actual source code: dsnep.c

slepc-3.18.0 2022-10-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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 <slepc/private/rgimpl.h>
 13: #include <slepcblaslapack.h>

 15: typedef struct {
 16:   PetscInt       nf;                 /* number of functions in f[] */
 17:   FN             f[DS_NUM_EXTRA];    /* functions defining the nonlinear operator */
 18:   PetscInt       max_mid;            /* maximum minimality index */
 19:   PetscInt       nnod;               /* number of nodes for quadrature rules */
 20:   PetscInt       spls;               /* number of sampling columns for quadrature rules */
 21:   PetscInt       Nit;                /* number of refinement iterations */
 22:   PetscReal      rtol;               /* tolerance of Newton refinement */
 23:   RG             rg;                 /* region for contour integral */
 24:   PetscLayout    map;                /* used to distribute work among MPI processes */
 25:   void           *computematrixctx;
 26:   PetscErrorCode (*computematrix)(DS,PetscScalar,PetscBool,DSMatType,void*);
 27: } DS_NEP;

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

 42:   PetscLogEventBegin(DS_Other,ds,0,0,0);
 43:   if (ctx->computematrix) (*ctx->computematrix)(ds,lambda,deriv,mat,ctx->computematrixctx);
 44:   else {
 45:     DSGetDimensions(ds,&n,NULL,NULL,NULL);
 46:     DSGetLeadingDimension(ds,&ld);
 47:     PetscBLASIntCast(ld*n,&k);
 48:     MatDenseGetArray(ds->omat[mat],&T);
 49:     PetscArrayzero(T,k);
 50:     for (i=0;i<ctx->nf;i++) {
 51:       if (deriv) FNEvaluateDerivative(ctx->f[i],lambda,&alpha);
 52:       else FNEvaluateFunction(ctx->f[i],lambda,&alpha);
 53:       MatDenseGetArrayRead(ds->omat[DSMatExtra[i]],&E);
 54:       PetscCallBLAS("BLASaxpy",BLASaxpy_(&k,&alpha,E,&inc,T,&inc));
 55:       MatDenseRestoreArrayRead(ds->omat[DSMatExtra[i]],&E);
 56:     }
 57:     MatDenseRestoreArray(ds->omat[mat],&T);
 58:   }
 59:   PetscLogEventEnd(DS_Other,ds,0,0,0);
 60:   return 0;
 61: }

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

 68:   DSAllocateMat_Private(ds,DS_MAT_X);
 69:   for (i=0;i<ctx->nf;i++) DSAllocateMat_Private(ds,DSMatExtra[i]);
 70:   PetscFree(ds->perm);
 71:   PetscMalloc1(ld*ctx->max_mid,&ds->perm);
 72:   return 0;
 73: }

 75: PetscErrorCode DSView_NEP(DS ds,PetscViewer viewer)
 76: {
 77:   DS_NEP            *ctx = (DS_NEP*)ds->data;
 78:   PetscViewerFormat format;
 79:   PetscInt          i;
 80:   const char        *methodname[] = {
 81:                      "Successive Linear Problems",
 82:                      "Contour Integral"
 83:   };
 84:   const int         nmeth=PETSC_STATIC_ARRAY_LENGTH(methodname);

 86:   PetscViewerGetFormat(viewer,&format);
 87:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
 88:     if (ds->method<nmeth) PetscViewerASCIIPrintf(viewer,"solving the problem with: %s\n",methodname[ds->method]);
 89: #if defined(PETSC_USE_COMPLEX)
 90:     if (ds->method==1) {  /* contour integral method */
 91:       PetscViewerASCIIPrintf(viewer,"number of integration points: %" PetscInt_FMT "\n",ctx->nnod);
 92:       PetscViewerASCIIPrintf(viewer,"maximum minimality index: %" PetscInt_FMT "\n",ctx->max_mid);
 93:       if (ctx->spls) PetscViewerASCIIPrintf(viewer,"number of sampling columns for quadrature: %" PetscInt_FMT "\n",ctx->spls);
 94:       if (ctx->Nit) PetscViewerASCIIPrintf(viewer,"doing iterative refinement (%" PetscInt_FMT " its, tolerance %g)\n",ctx->Nit,(double)ctx->rtol);
 95:       RGView(ctx->rg,viewer);
 96:     }
 97: #endif
 98:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscViewerASCIIPrintf(viewer,"number of functions: %" PetscInt_FMT "\n",ctx->nf);
 99:     return 0;
100:   }
101:   for (i=0;i<ctx->nf;i++) {
102:     FNView(ctx->f[i],viewer);
103:     DSViewMat(ds,viewer,DSMatExtra[i]);
104:   }
105:   if (ds->state>DS_STATE_INTERMEDIATE) DSViewMat(ds,viewer,DS_MAT_X);
106:   return 0;
107: }

109: PetscErrorCode DSVectors_NEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
110: {
112:   switch (mat) {
113:     case DS_MAT_X:
114:       break;
115:     case DS_MAT_Y:
116:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
117:     default:
118:       SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
119:   }
120:   return 0;
121: }

123: PetscErrorCode DSSort_NEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *dummy)
124: {
125:   DS_NEP         *ctx = (DS_NEP*)ds->data;
126:   PetscInt       n,l,i,*perm,lds;
127:   PetscScalar    *Q;

129:   if (!ds->sc) return 0;
130:   if (!ds->method) return 0;  /* SLP computes just one eigenvalue */
131:   n = ds->n*ctx->max_mid;
132:   lds = ds->ld*ctx->max_mid;
133:   l = ds->l;
134:   perm = ds->perm;
135:   for (i=0;i<n;i++) perm[i] = i;
136:   if (rr) DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
137:   else DSSortEigenvalues_Private(ds,wr,NULL,perm,PETSC_FALSE);
138:   MatDenseGetArray(ds->omat[DS_MAT_Q],&Q);
139:   for (i=l;i<ds->t;i++) Q[i+i*lds] = wr[perm[i]];
140:   for (i=l;i<ds->t;i++) wr[i] = Q[i+i*lds];
141:   MatDenseRestoreArray(ds->omat[DS_MAT_Q],&Q);
142:   /* n != ds->n */
143:   DSPermuteColumns_Private(ds,0,ds->t,ds->n,DS_MAT_X,perm);
144:   return 0;
145: }

147: PetscErrorCode DSSolve_NEP_SLP(DS ds,PetscScalar *wr,PetscScalar *wi)
148: {
149:   PetscScalar    *A,*B,*W,*X,*work,*alpha,*beta;
150:   PetscScalar    sigma,lambda,mu,re,re2,sone=1.0,szero=0.0;
151:   PetscBLASInt   info,n,ld,lrwork=0,lwork,one=1,zero=0;
152:   PetscInt       it,pos,j,maxit=100,result;
153:   PetscReal      norm,tol,done=1.0;
154: #if defined(PETSC_USE_COMPLEX)
155:   PetscReal      *rwork;
156: #else
157:   PetscReal      *alphai,im,im2;
158: #endif

160:   PetscBLASIntCast(ds->n,&n);
161:   PetscBLASIntCast(ds->ld,&ld);
162: #if defined(PETSC_USE_COMPLEX)
163:   PetscBLASIntCast(2*ds->n+2*ds->n,&lwork);
164:   PetscBLASIntCast(8*ds->n,&lrwork);
165: #else
166:   PetscBLASIntCast(3*ds->n+8*ds->n,&lwork);
167: #endif
168:   DSAllocateWork_Private(ds,lwork,lrwork,0);
169:   alpha = ds->work;
170:   beta = ds->work + ds->n;
171: #if defined(PETSC_USE_COMPLEX)
172:   work = ds->work + 2*ds->n;
173:   lwork -= 2*ds->n;
174: #else
175:   alphai = ds->work + 2*ds->n;
176:   work = ds->work + 3*ds->n;
177:   lwork -= 3*ds->n;
178: #endif
179:   DSAllocateMat_Private(ds,DS_MAT_A);
180:   DSAllocateMat_Private(ds,DS_MAT_B);
181:   DSAllocateMat_Private(ds,DS_MAT_W);
182:   MatDenseGetArray(ds->omat[DS_MAT_A],&A);
183:   MatDenseGetArray(ds->omat[DS_MAT_B],&B);
184:   MatDenseGetArray(ds->omat[DS_MAT_W],&W);
185:   MatDenseGetArray(ds->omat[DS_MAT_X],&X);

187:   sigma = 0.0;
188:   if (ds->sc->comparison==SlepcCompareTargetMagnitude || ds->sc->comparison==SlepcCompareTargetReal) sigma = *(PetscScalar*)ds->sc->comparisonctx;
189:   lambda = sigma;
190:   tol = n*PETSC_MACHINE_EPSILON/PetscSqrtReal(PETSC_SQRT_MACHINE_EPSILON);

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

194:     /* evaluate T and T' */
195:     DSNEPComputeMatrix(ds,lambda,PETSC_FALSE,DS_MAT_A);
196:     if (it) {
197:       PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,A,&ld,X,&one,&szero,X+ld,&one));
198:       norm = BLASnrm2_(&n,X+ld,&one);
199:       if (norm/PetscAbsScalar(lambda)<=tol) break;
200:     }
201:     DSNEPComputeMatrix(ds,lambda,PETSC_TRUE,DS_MAT_B);

203:     /* compute eigenvalue correction mu and eigenvector u */
204: #if defined(PETSC_USE_COMPLEX)
205:     rwork = ds->rwork;
206:     PetscCallBLAS("LAPACKggev",LAPACKggev_("N","V",&n,A,&ld,B,&ld,alpha,beta,NULL,&ld,W,&ld,work,&lwork,rwork,&info));
207: #else
208:     PetscCallBLAS("LAPACKggev",LAPACKggev_("N","V",&n,A,&ld,B,&ld,alpha,alphai,beta,NULL,&ld,W,&ld,work,&lwork,&info));
209: #endif
210:     SlepcCheckLapackInfo("ggev",info);

212:     /* find smallest eigenvalue */
213:     j = 0;
214:     if (beta[j]==0.0) re = (PetscRealPart(alpha[j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
215:     else re = alpha[j]/beta[j];
216: #if !defined(PETSC_USE_COMPLEX)
217:     if (beta[j]==0.0) im = (alphai[j]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
218:     else im = alphai[j]/beta[j];
219: #endif
220:     pos = 0;
221:     for (j=1;j<n;j++) {
222:       if (beta[j]==0.0) re2 = (PetscRealPart(alpha[j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
223:       else re2 = alpha[j]/beta[j];
224: #if !defined(PETSC_USE_COMPLEX)
225:       if (beta[j]==0.0) im2 = (alphai[j]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
226:       else im2 = alphai[j]/beta[j];
227:       SlepcCompareSmallestMagnitude(re,im,re2,im2,&result,NULL);
228: #else
229:       SlepcCompareSmallestMagnitude(re,0.0,re2,0.0,&result,NULL);
230: #endif
231:       if (result > 0) {
232:         re = re2;
233: #if !defined(PETSC_USE_COMPLEX)
234:         im = im2;
235: #endif
236:         pos = j;
237:       }
238:     }

240: #if !defined(PETSC_USE_COMPLEX)
242: #endif
243:     mu = alpha[pos]/beta[pos];
244:     PetscArraycpy(X,W+pos*ld,n);
245:     norm = BLASnrm2_(&n,X,&one);
246:     PetscCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&one,X,&n,&info));
247:     SlepcCheckLapackInfo("lascl",info);

249:     /* correct eigenvalue approximation */
250:     lambda = lambda - mu;
251:   }
252:   MatDenseRestoreArray(ds->omat[DS_MAT_A],&A);
253:   MatDenseRestoreArray(ds->omat[DS_MAT_B],&B);
254:   MatDenseRestoreArray(ds->omat[DS_MAT_W],&W);
255:   MatDenseRestoreArray(ds->omat[DS_MAT_X],&X);

258:   ds->t = 1;
259:   wr[0] = lambda;
260:   if (wi) wi[0] = 0.0;
261:   return 0;
262: }

264: #if defined(PETSC_USE_COMPLEX)
265: /*
266:   Newton refinement for eigenpairs computed with contour integral.
267:   k  - number of eigenpairs to refine
268:   wr - eigenvalues (eigenvectors are stored in DS_MAT_X)
269: */
270: static PetscErrorCode DSNEPNewtonRefine(DS ds,PetscInt k,PetscScalar *wr)
271: {
272:   DS_NEP         *ctx = (DS_NEP*)ds->data;
273:   PetscScalar    *X,*W,*U,*R,sone=1.0,szero=0.0;
274:   PetscReal      norm;
275:   PetscInt       i,j,ii,nwu=0,*p,jstart=0,jend=k;
276:   const PetscInt *range;
277:   PetscBLASInt   n,*perm,info,ld,one=1,n1;
278:   PetscMPIInt    len,size,root;
279:   PetscLayout    map;

281:   MatDenseGetArray(ds->omat[DS_MAT_X],&X);
282:   PetscBLASIntCast(ds->n,&n);
283:   PetscBLASIntCast(ds->ld,&ld);
284:   n1 = n+1;
285:   p  = ds->perm;
286:   PetscArrayzero(p,k);
287:   DSAllocateWork_Private(ds,(n+2)*(n+1),0,n+1);
288:   U    = ds->work+nwu;    nwu += (n+1)*(n+1);
289:   R    = ds->work+nwu;    /*nwu += n+1;*/
290:   perm = ds->iwork;
291:   if (ds->pmode==DS_PARALLEL_DISTRIBUTED) {
292:     PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)ds),PETSC_DECIDE,k,1,&map);
293:     PetscLayoutGetRange(map,&jstart,&jend);
294:   }
295:   for (ii=0;ii<ctx->Nit;ii++) {
296:     for (j=jstart;j<jend;j++) {
297:       if (p[j]<2) {
298:         DSNEPComputeMatrix(ds,wr[j],PETSC_FALSE,DS_MAT_W);
299:         MatDenseGetArray(ds->omat[DS_MAT_W],&W);
300:         PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,W,&ld,X+ld*j,&one,&szero,R,&one));
301:         MatDenseRestoreArray(ds->omat[DS_MAT_W],&W);
302:         norm = BLASnrm2_(&n,R,&one);
303:         if (norm/PetscAbsScalar(wr[j]) > ctx->rtol) {
304:           PetscInfo(NULL,"Refining eigenpair %" PetscInt_FMT ", residual=%g\n",j,(double)(norm/PetscAbsScalar(wr[j])));
305:           p[j] = 1;
306:           R[n] = 0.0;
307:           MatDenseGetArray(ds->omat[DS_MAT_W],&W);
308:           for (i=0;i<n;i++) {
309:             PetscArraycpy(U+i*n1,W+i*ld,n);
310:             U[n+i*n1] = PetscConj(X[j*ld+i]);
311:           }
312:           MatDenseRestoreArray(ds->omat[DS_MAT_W],&W);
313:           U[n+n*n1] = 0.0;
314:           DSNEPComputeMatrix(ds,wr[j],PETSC_TRUE,DS_MAT_W);
315:           MatDenseGetArray(ds->omat[DS_MAT_W],&W);
316:           PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,W,&ld,X+ld*j,&one,&szero,U+n*(n+1),&one));
317:           MatDenseRestoreArray(ds->omat[DS_MAT_W],&W);
318:           /* solve system  */
319:           PetscFPTrapPush(PETSC_FP_TRAP_OFF);
320:           PetscCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n1,&n1,U,&n1,perm,&info));
321:           SlepcCheckLapackInfo("getrf",info);
322:           PetscCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&n1,&one,U,&n1,perm,R,&n1,&info));
323:           SlepcCheckLapackInfo("getrs",info);
324:           PetscFPTrapPop();
325:           wr[j] -= R[n];
326:           for (i=0;i<n;i++) X[j*ld+i] -= R[i];
327:           /* normalization */
328:           norm = BLASnrm2_(&n,X+ld*j,&one);
329:           for (i=0;i<n;i++) X[ld*j+i] /= norm;
330:         } else p[j] = 2;
331:       }
332:     }
333:   }
334:   if (ds->pmode==DS_PARALLEL_DISTRIBUTED) {  /* communicate results */
335:     PetscMPIIntCast(k,&len);
336:     MPIU_Allreduce(MPI_IN_PLACE,p,len,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)ds));
337:     MPI_Comm_size(PetscObjectComm((PetscObject)ds),&size);
338:     PetscLayoutGetRanges(map,&range);
339:     for (j=0;j<k;j++) {
340:       if (p[j]) {  /* j-th eigenpair has been refined */
341:         for (root=0;root<size;root++) if (range[root+1]>j) break;
342:         PetscMPIIntCast(1,&len);
343:         MPI_Bcast(wr+j,len,MPIU_SCALAR,root,PetscObjectComm((PetscObject)ds));
344:         PetscMPIIntCast(n,&len);
345:         MPI_Bcast(X+ld*j,len,MPIU_SCALAR,root,PetscObjectComm((PetscObject)ds));
346:       }
347:     }
348:     PetscLayoutDestroy(&map);
349:   }
350:   MatDenseRestoreArray(ds->omat[DS_MAT_X],&X);
351:   return 0;
352: }

354: PetscErrorCode DSSolve_NEP_Contour(DS ds,PetscScalar *wr,PetscScalar *wi)
355: {
356:   DS_NEP         *ctx = (DS_NEP*)ds->data;
357:   PetscScalar    *alpha,*beta,*Q,*Z,*X,*U,*V,*W,*work,*Rc,*R,*w,*z,*zn,*S;
358:   PetscScalar    sone=1.0,szero=0.0,center;
359:   PetscReal      *rwork,norm,radius,vscale,rgscale,*sigma;
360:   PetscBLASInt   info,n,*perm,p,pp,ld,lwork,k_,rk_,colA,rowA,one=1;
361:   PetscInt       mid,lds,nnod=ctx->nnod,k,i,ii,jj,j,s,off,rk,nwu=0,nw,lrwork,*inside,kstart=0,kend=nnod;
362:   PetscMPIInt    len;
363:   PetscBool      isellipse;
364:   PetscRandom    rand;

367:   /* Contour parameters */
368:   PetscObjectTypeCompare((PetscObject)ctx->rg,RGELLIPSE,&isellipse);
370:   RGEllipseGetParameters(ctx->rg,&center,&radius,&vscale);
371:   RGGetScale(ctx->rg,&rgscale);
372:   if (ds->pmode==DS_PARALLEL_DISTRIBUTED) {
373:     if (!ctx->map) PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)ds),PETSC_DECIDE,ctx->nnod,1,&ctx->map);
374:     PetscLayoutGetRange(ctx->map,&kstart,&kend);
375:   }

377:   DSAllocateMat_Private(ds,DS_MAT_W); /* size n */
378:   DSAllocateMat_Private(ds,DS_MAT_Q); /* size mid*n */
379:   DSAllocateMat_Private(ds,DS_MAT_Z); /* size mid*n */
380:   DSAllocateMat_Private(ds,DS_MAT_U); /* size mid*n */
381:   DSAllocateMat_Private(ds,DS_MAT_V); /* size mid*n */
382:   MatDenseGetArray(ds->omat[DS_MAT_Q],&Q);
383:   MatDenseGetArray(ds->omat[DS_MAT_Z],&Z);
384:   MatDenseGetArray(ds->omat[DS_MAT_U],&U);
385:   MatDenseGetArray(ds->omat[DS_MAT_V],&V);
386:   mid  = ctx->max_mid;
387:   PetscBLASIntCast(ds->n,&n);
388:   p    = n;   /* maximum number of columns for the probing matrix */
389:   PetscBLASIntCast(ds->ld,&ld);
390:   PetscBLASIntCast(mid*n,&rowA);
391:   PetscBLASIntCast(5*rowA,&lwork);
392:   nw   = n*(2*p+7*mid)+3*nnod+2*mid*n*p;
393:   lrwork = mid*n*6+8*n;
394:   DSAllocateWork_Private(ds,nw,lrwork,n+1);

396:   sigma = ds->rwork;
397:   rwork = ds->rwork+mid*n;
398:   perm  = ds->iwork;
399:   z     = ds->work+nwu;    nwu += nnod;         /* quadrature points */
400:   zn    = ds->work+nwu;    nwu += nnod;         /* normalized quadrature points */
401:   w     = ds->work+nwu;    nwu += nnod;         /* quadrature weights */
402:   Rc    = ds->work+nwu;    nwu += n*p;
403:   R     = ds->work+nwu;    nwu += n*p;
404:   alpha = ds->work+nwu;    nwu += mid*n;
405:   beta  = ds->work+nwu;    nwu += mid*n;
406:   S     = ds->work+nwu;    nwu += 2*mid*n*p;
407:   work  = ds->work+nwu;    /*nwu += mid*n*5;*/

409:   /* Compute quadrature parameters */
410:   RGComputeQuadrature(ctx->rg,RG_QUADRULE_TRAPEZOIDAL,nnod,z,zn,w);

412:   /* Set random matrix */
413:   PetscRandomCreate(PetscObjectComm((PetscObject)ds),&rand);
414:   PetscRandomSetSeed(rand,0x12345678);
415:   PetscRandomSeed(rand);
416:   for (j=0;j<p;j++)
417:     for (i=0;i<n;i++) PetscRandomGetValue(rand,Rc+i+j*n);
418:   PetscArrayzero(S,2*mid*n*p);
419:   /* Loop of integration points */
420:   for (k=kstart;k<kend;k++) {
421:     PetscInfo(NULL,"Solving integration point %" PetscInt_FMT "\n",k);
422:     PetscArraycpy(R,Rc,p*n);
423:     DSNEPComputeMatrix(ds,z[k],PETSC_FALSE,DS_MAT_W);

425:     /* LU factorization */
426:     MatDenseGetArray(ds->omat[DS_MAT_W],&W);
427:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
428:     PetscCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,W,&ld,perm,&info));
429:     SlepcCheckLapackInfo("getrf",info);
430:     PetscCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&n,&p,W,&ld,perm,R,&n,&info));
431:     SlepcCheckLapackInfo("getrs",info);
432:     PetscFPTrapPop();
433:     MatDenseRestoreArray(ds->omat[DS_MAT_W],&W);

435:     /* Moments computation */
436:     for (s=0;s<2*ctx->max_mid;s++) {
437:       off = s*n*p;
438:       for (j=0;j<p;j++)
439:         for (i=0;i<n;i++) S[off+i+j*n] += w[k]*R[j*n+i];
440:       w[k] *= zn[k];
441:     }
442:   }

444:   if (ds->pmode==DS_PARALLEL_DISTRIBUTED) {  /* compute final S via reduction */
445:     PetscMPIIntCast(2*mid*n*p,&len);
446:     MPIU_Allreduce(MPI_IN_PLACE,S,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)ds));
447:   }
448:   p = ctx->spls?PetscMin(ctx->spls,n):n;
449:   pp = p;
450:   do {
451:     p = pp;
452:     PetscBLASIntCast(mid*p,&colA);

454:     PetscInfo(ds,"Computing SVD of size %" PetscBLASInt_FMT "x%" PetscBLASInt_FMT "\n",rowA,colA);
455:     for (jj=0;jj<mid;jj++) {
456:       for (ii=0;ii<mid;ii++) {
457:         off = jj*p*rowA+ii*n;
458:         for (j=0;j<p;j++)
459:           for (i=0;i<n;i++) Q[off+j*rowA+i] = S[((jj+ii)*n+j)*n+i];
460:       }
461:     }
462:     PetscFPTrapPush(PETSC_FP_TRAP_OFF);
463:     PetscCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&rowA,&colA,Q,&rowA,sigma,U,&rowA,V,&colA,work,&lwork,rwork,&info));
464:     SlepcCheckLapackInfo("gesvd",info);
465:     PetscFPTrapPop();

467:     rk = colA;
468:     for (i=1;i<colA;i++) if (sigma[i]/sigma[0]<PETSC_MACHINE_EPSILON*1e4) {rk = i; break;}
469:     if (rk<colA || p==n) break;
470:     pp *= 2;
471:   } while (pp<=n);
472:   PetscInfo(ds,"Solving generalized eigenproblem of size %" PetscInt_FMT "\n",rk);
473:   for (jj=0;jj<mid;jj++) {
474:     for (ii=0;ii<mid;ii++) {
475:       off = jj*p*rowA+ii*n;
476:       for (j=0;j<p;j++)
477:         for (i=0;i<n;i++) Q[off+j*rowA+i] = S[((jj+ii+1)*n+j)*n+i];
478:     }
479:   }
480:   PetscBLASIntCast(rk,&rk_);
481:   PetscCallBLAS("BLASgemm",BLASgemm_("N","C",&rowA,&rk_,&colA,&sone,Q,&rowA,V,&colA,&szero,Z,&rowA));
482:   PetscCallBLAS("BLASgemm",BLASgemm_("C","N",&rk_,&rk_,&rowA,&sone,U,&rowA,Z,&rowA,&szero,Q,&rk_));
483:   PetscArrayzero(Z,n*mid*n*mid);
484:   for (j=0;j<rk;j++) Z[j+j*rk_] = sigma[j];
485:   PetscCallBLAS("LAPACKggev",LAPACKggev_("N","V",&rk_,Q,&rk_,Z,&rk_,alpha,beta,NULL,&ld,V,&rk_,work,&lwork,rwork,&info));
486:   for (i=0;i<rk;i++) wr[i] = (center+radius*alpha[i]/beta[i])*rgscale;
487:   PetscMalloc1(rk,&inside);
488:   RGCheckInside(ctx->rg,rk,wr,wi,inside);
489:   k=0;
490:   for (i=0;i<rk;i++)
491:     if (inside[i]==1) inside[k++] = i;
492:   /* Discard values outside region */
493:   lds = ld*mid;
494:   PetscArrayzero(Q,lds*lds);
495:   PetscArrayzero(Z,lds*lds);
496:   for (i=0;i<k;i++) Q[i+i*lds] = (center*beta[inside[i]]+radius*alpha[inside[i]])*rgscale;
497:   for (i=0;i<k;i++) Z[i+i*lds] = beta[inside[i]];
498:   for (i=0;i<k;i++) wr[i] = Q[i+i*lds]/Z[i+i*lds];
499:   for (j=0;j<k;j++) for (i=0;i<rk;i++) V[j*rk+i] = sigma[i]*V[inside[j]*rk+i];
500:   PetscBLASIntCast(k,&k_);
501:   MatDenseGetArray(ds->omat[DS_MAT_X],&X);
502:   PetscCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&k_,&rk_,&sone,U,&rowA,V,&rk_,&szero,X,&ld));
503:   /* Normalize */
504:   for (j=0;j<k;j++) {
505:     norm = BLASnrm2_(&n,X+ld*j,&one);
506:     for (i=0;i<n;i++) X[ld*j+i] /= norm;
507:   }
508:   PetscFree(inside);
509:   MatDenseRestoreArray(ds->omat[DS_MAT_X],&X);
510:   MatDenseRestoreArray(ds->omat[DS_MAT_Q],&Q);
511:   MatDenseRestoreArray(ds->omat[DS_MAT_Z],&Z);
512:   MatDenseRestoreArray(ds->omat[DS_MAT_U],&U);
513:   MatDenseRestoreArray(ds->omat[DS_MAT_V],&V);

515:   /* Newton refinement */
516:   DSNEPNewtonRefine(ds,k,wr);
517:   ds->t = k;
518:   PetscRandomDestroy(&rand);
519:   return 0;
520: }
521: #endif

523: PetscErrorCode DSSynchronize_NEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
524: {
525:   DS_NEP         *ctx = (DS_NEP*)ds->data;
526:   PetscInt       ld=ds->ld,k=0;
527:   PetscMPIInt    n,n2,rank,size,off=0;
528:   PetscScalar    *X;

530:   if (!ds->method) { /* SLP */
531:     if (ds->state>=DS_STATE_CONDENSED) k += ds->n;
532:     if (eigr) k += 1;
533:     if (eigi) k += 1;
534:     PetscMPIIntCast(1,&n);
535:     PetscMPIIntCast(ds->n,&n2);
536:   } else { /* Contour */
537:     if (ds->state>=DS_STATE_CONDENSED) k += ctx->max_mid*ds->n*ld;
538:     if (eigr) k += ctx->max_mid*ds->n;
539:     if (eigi) k += ctx->max_mid*ds->n;
540:     PetscMPIIntCast(ctx->max_mid*ds->n,&n);
541:     PetscMPIIntCast(ctx->max_mid*ds->n*ld,&n2);
542:   }
543:   DSAllocateWork_Private(ds,k,0,0);
544:   PetscMPIIntCast(k*sizeof(PetscScalar),&size);
545:   if (ds->state>=DS_STATE_CONDENSED) MatDenseGetArray(ds->omat[DS_MAT_X],&X);
546:   MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
547:   if (!rank) {
548:     if (ds->state>=DS_STATE_CONDENSED) MPI_Pack(X,n2,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
549:     if (eigr) MPI_Pack(eigr,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
550: #if !defined(PETSC_USE_COMPLEX)
551:     if (eigi) MPI_Pack(eigi,n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
552: #endif
553:   }
554:   MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
555:   if (rank) {
556:     if (ds->state>=DS_STATE_CONDENSED) MPI_Unpack(ds->work,size,&off,X,n2,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
557:     if (eigr) MPI_Unpack(ds->work,size,&off,eigr,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
558: #if !defined(PETSC_USE_COMPLEX)
559:     if (eigi) MPI_Unpack(ds->work,size,&off,eigi,n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
560: #endif
561:   }
562:   if (ds->state>=DS_STATE_CONDENSED) MatDenseRestoreArray(ds->omat[DS_MAT_X],&X);
563:   return 0;
564: }

566: static PetscErrorCode DSNEPSetFN_NEP(DS ds,PetscInt n,FN fn[])
567: {
568:   DS_NEP         *ctx = (DS_NEP*)ds->data;
569:   PetscInt       i;

573:   if (ds->ld) PetscInfo(ds,"DSNEPSetFN() called after DSAllocate()\n");
574:   for (i=0;i<n;i++) PetscObjectReference((PetscObject)fn[i]);
575:   for (i=0;i<ctx->nf;i++) FNDestroy(&ctx->f[i]);
576:   for (i=0;i<n;i++) ctx->f[i] = fn[i];
577:   ctx->nf = n;
578:   return 0;
579: }

581: /*@
582:    DSNEPSetFN - Sets a number of functions that define the nonlinear
583:    eigenproblem.

585:    Collective on ds

587:    Input Parameters:
588: +  ds - the direct solver context
589: .  n  - number of functions
590: -  fn - array of functions

592:    Notes:
593:    The nonlinear eigenproblem is defined in terms of the split nonlinear
594:    operator T(lambda) = sum_i A_i*f_i(lambda).

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

600:    Level: advanced

602: .seealso: DSNEPGetFN(), DSAllocate()
603: @*/
604: PetscErrorCode DSNEPSetFN(DS ds,PetscInt n,FN fn[])
605: {
606:   PetscInt       i;

611:   for (i=0;i<n;i++) {
614:   }
615:   PetscTryMethod(ds,"DSNEPSetFN_C",(DS,PetscInt,FN[]),(ds,n,fn));
616:   return 0;
617: }

619: static PetscErrorCode DSNEPGetFN_NEP(DS ds,PetscInt k,FN *fn)
620: {
621:   DS_NEP *ctx = (DS_NEP*)ds->data;

624:   *fn = ctx->f[k];
625:   return 0;
626: }

628: /*@
629:    DSNEPGetFN - Gets the functions associated with the nonlinear DS.

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

633:    Input Parameters:
634: +  ds - the direct solver context
635: -  k  - the index of the requested function (starting in 0)

637:    Output Parameter:
638: .  fn - the function

640:    Level: advanced

642: .seealso: DSNEPSetFN()
643: @*/
644: PetscErrorCode DSNEPGetFN(DS ds,PetscInt k,FN *fn)
645: {
648:   PetscUseMethod(ds,"DSNEPGetFN_C",(DS,PetscInt,FN*),(ds,k,fn));
649:   return 0;
650: }

652: static PetscErrorCode DSNEPGetNumFN_NEP(DS ds,PetscInt *n)
653: {
654:   DS_NEP *ctx = (DS_NEP*)ds->data;

656:   *n = ctx->nf;
657:   return 0;
658: }

660: /*@
661:    DSNEPGetNumFN - Returns the number of functions stored internally by
662:    the DS.

664:    Not collective

666:    Input Parameter:
667: .  ds - the direct solver context

669:    Output Parameters:
670: .  n - the number of functions passed in DSNEPSetFN()

672:    Level: advanced

674: .seealso: DSNEPSetFN()
675: @*/
676: PetscErrorCode DSNEPGetNumFN(DS ds,PetscInt *n)
677: {
680:   PetscUseMethod(ds,"DSNEPGetNumFN_C",(DS,PetscInt*),(ds,n));
681:   return 0;
682: }

684: static PetscErrorCode DSNEPSetMinimality_NEP(DS ds,PetscInt n)
685: {
686:   DS_NEP *ctx = (DS_NEP*)ds->data;

688:   if (n == PETSC_DECIDE || n == PETSC_DEFAULT) ctx->max_mid = 4;
689:   else {
691:     ctx->max_mid = n;
692:   }
693:   return 0;
694: }

696: /*@
697:    DSNEPSetMinimality - Sets the maximum minimality index used internally by
698:    the DSNEP.

700:    Logically Collective on ds

702:    Input Parameters:
703: +  ds - the direct solver context
704: -  n  - the maximum minimality index

706:    Options Database Key:
707: .  -ds_nep_minimality <n> - sets the maximum minimality index

709:    Notes:
710:    The maximum minimality index is used only in the contour integral method,
711:    and is related to the highest momemts used in the method. The default
712:    value is 1, an larger value might give better accuracy in some cases, but
713:    at a higher cost.

715:    Level: advanced

717: .seealso: DSNEPGetMinimality()
718: @*/
719: PetscErrorCode DSNEPSetMinimality(DS ds,PetscInt n)
720: {
723:   PetscTryMethod(ds,"DSNEPSetMinimality_C",(DS,PetscInt),(ds,n));
724:   return 0;
725: }

727: static PetscErrorCode DSNEPGetMinimality_NEP(DS ds,PetscInt *n)
728: {
729:   DS_NEP *ctx = (DS_NEP*)ds->data;

731:   *n = ctx->max_mid;
732:   return 0;
733: }

735: /*@
736:    DSNEPGetMinimality - Returns the maximum minimality index used internally by
737:    the DSNEP.

739:    Not collective

741:    Input Parameter:
742: .  ds - the direct solver context

744:    Output Parameters:
745: .  n - the maximum minimality index passed in DSNEPSetMinimality()

747:    Level: advanced

749: .seealso: DSNEPSetMinimality()
750: @*/
751: PetscErrorCode DSNEPGetMinimality(DS ds,PetscInt *n)
752: {
755:   PetscUseMethod(ds,"DSNEPGetMinimality_C",(DS,PetscInt*),(ds,n));
756:   return 0;
757: }

759: static PetscErrorCode DSNEPSetRefine_NEP(DS ds,PetscReal tol,PetscInt its)
760: {
761:   DS_NEP *ctx = (DS_NEP*)ds->data;

763:   if (tol == PETSC_DEFAULT) ctx->rtol = PETSC_MACHINE_EPSILON/PetscSqrtReal(PETSC_SQRT_MACHINE_EPSILON);
764:   else {
766:     ctx->rtol = tol;
767:   }
768:   if (its == PETSC_DECIDE || its == PETSC_DEFAULT) ctx->Nit = 3;
769:   else {
771:     ctx->Nit = its;
772:   }
773:   return 0;
774: }

776: /*@
777:    DSNEPSetRefine - Sets the tolerance and the number of iterations of Newton iterative
778:    refinement for eigenpairs.

780:    Logically Collective on ds

782:    Input Parameters:
783: +  ds  - the direct solver context
784: .  tol - the tolerance
785: -  its - the number of iterations

787:    Options Database Key:
788: +  -ds_nep_refine_tol <tol> - sets the tolerance
789: -  -ds_nep_refine_its <its> - sets the number of Newton iterations

791:    Notes:
792:    Iterative refinement of eigenpairs is currently used only in the contour
793:    integral method.

795:    Level: advanced

797: .seealso: DSNEPGetRefine()
798: @*/
799: PetscErrorCode DSNEPSetRefine(DS ds,PetscReal tol,PetscInt its)
800: {
804:   PetscTryMethod(ds,"DSNEPSetRefine_C",(DS,PetscReal,PetscInt),(ds,tol,its));
805:   return 0;
806: }

808: static PetscErrorCode DSNEPGetRefine_NEP(DS ds,PetscReal *tol,PetscInt *its)
809: {
810:   DS_NEP *ctx = (DS_NEP*)ds->data;

812:   if (tol) *tol = ctx->rtol;
813:   if (its) *its = ctx->Nit;
814:   return 0;
815: }

817: /*@
818:    DSNEPGetRefine - Returns the tolerance and the number of iterations of Newton iterative
819:    refinement for eigenpairs.

821:    Not collective

823:    Input Parameter:
824: .  ds - the direct solver context

826:    Output Parameters:
827: +  tol - the tolerance
828: -  its - the number of iterations

830:    Level: advanced

832: .seealso: DSNEPSetRefine()
833: @*/
834: PetscErrorCode DSNEPGetRefine(DS ds,PetscReal *tol,PetscInt *its)
835: {
837:   PetscUseMethod(ds,"DSNEPGetRefine_C",(DS,PetscReal*,PetscInt*),(ds,tol,its));
838:   return 0;
839: }

841: static PetscErrorCode DSNEPSetIntegrationPoints_NEP(DS ds,PetscInt ip)
842: {
843:   DS_NEP         *ctx = (DS_NEP*)ds->data;

845:   if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) ctx->nnod = 64;
846:   else {
848:     ctx->nnod = ip;
849:   }
850:   PetscLayoutDestroy(&ctx->map);  /* need to redistribute at next solve */
851:   return 0;
852: }

854: /*@
855:    DSNEPSetIntegrationPoints - Sets the number of integration points to be
856:    used in the contour integral method.

858:    Logically Collective on ds

860:    Input Parameters:
861: +  ds - the direct solver context
862: -  ip - the number of integration points

864:    Options Database Key:
865: .  -ds_nep_integration_points <ip> - sets the number of integration points

867:    Notes:
868:    This parameter is relevant only in the contour integral method.

870:    Level: advanced

872: .seealso: DSNEPGetIntegrationPoints()
873: @*/
874: PetscErrorCode DSNEPSetIntegrationPoints(DS ds,PetscInt ip)
875: {
878:   PetscTryMethod(ds,"DSNEPSetIntegrationPoints_C",(DS,PetscInt),(ds,ip));
879:   return 0;
880: }

882: static PetscErrorCode DSNEPGetIntegrationPoints_NEP(DS ds,PetscInt *ip)
883: {
884:   DS_NEP *ctx = (DS_NEP*)ds->data;

886:   *ip = ctx->nnod;
887:   return 0;
888: }

890: /*@
891:    DSNEPGetIntegrationPoints - Returns the number of integration points used
892:    in the contour integral method.

894:    Not collective

896:    Input Parameter:
897: .  ds - the direct solver context

899:    Output Parameters:
900: .  ip - the number of integration points

902:    Level: advanced

904: .seealso: DSNEPSetIntegrationPoints()
905: @*/
906: PetscErrorCode DSNEPGetIntegrationPoints(DS ds,PetscInt *ip)
907: {
910:   PetscUseMethod(ds,"DSNEPGetIntegrationPoints_C",(DS,PetscInt*),(ds,ip));
911:   return 0;
912: }

914: static PetscErrorCode DSNEPSetSamplingSize_NEP(DS ds,PetscInt p)
915: {
916:   DS_NEP *ctx = (DS_NEP*)ds->data;

918:   if (p == PETSC_DECIDE || p == PETSC_DEFAULT) ctx->spls = 0;
919:   else {
921:     ctx->spls = p;
922:   }
923:   return 0;
924: }

926: /*@
927:    DSNEPSetSamplingSize - Sets the number of sampling columns to be
928:    used in the contour integral method.

930:    Logically Collective on ds

932:    Input Parameters:
933: +  ds - the direct solver context
934: -  p - the number of columns for the sampling matrix

936:    Options Database Key:
937: .  -ds_nep_sampling_size <p> - sets the number of sampling columns

939:    Notes:
940:    This parameter is relevant only in the contour integral method.

942:    Level: advanced

944: .seealso: DSNEPGetSamplingSize()
945: @*/
946: PetscErrorCode DSNEPSetSamplingSize(DS ds,PetscInt p)
947: {
950:   PetscTryMethod(ds,"DSNEPSetSamplingSize_C",(DS,PetscInt),(ds,p));
951:   return 0;
952: }

954: static PetscErrorCode DSNEPGetSamplingSize_NEP(DS ds,PetscInt *p)
955: {
956:   DS_NEP *ctx = (DS_NEP*)ds->data;

958:   *p = ctx->spls;
959:   return 0;
960: }

962: /*@
963:    DSNEPGetSamplingSize - Returns the number of sampling columns used
964:    in the contour integral method.

966:    Not collective

968:    Input Parameter:
969: .  ds - the direct solver context

971:    Output Parameters:
972: .  p -  the number of columns for the sampling matrix

974:    Level: advanced

976: .seealso: DSNEPSetSamplingSize()
977: @*/
978: PetscErrorCode DSNEPGetSamplingSize(DS ds,PetscInt *p)
979: {
982:   PetscUseMethod(ds,"DSNEPGetSamplingSize_C",(DS,PetscInt*),(ds,p));
983:   return 0;
984: }

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

990:   dsctx->computematrix    = fun;
991:   dsctx->computematrixctx = ctx;
992:   return 0;
993: }

995: /*@C
996:    DSNEPSetComputeMatrixFunction - Sets a user-provided subroutine to compute
997:    the matrices T(lambda) or T'(lambda).

999:    Logically Collective on ds

1001:    Input Parameters:
1002: +  ds  - the direct solver context
1003: .  fun - a pointer to the user function
1004: -  ctx - a context pointer (the last parameter to the user function)

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

1009: +   ds     - the direct solver object
1010: .   lambda - point where T(lambda) or T'(lambda) must be evaluated
1011: .   deriv  - if true compute T'(lambda), otherwise compute T(lambda)
1012: .   mat    - the DS matrix where the result must be stored
1013: -   ctx    - optional context, as set by DSNEPSetComputeMatrixFunction()

1015:    Note:
1016:    The result is computed as T(lambda) = sum_i E_i*f_i(lambda), and similarly
1017:    for the derivative.

1019:    Level: developer

1021: .seealso: DSNEPGetComputeMatrixFunction()
1022: @*/
1023: PetscErrorCode DSNEPSetComputeMatrixFunction(DS ds,PetscErrorCode (*fun)(DS,PetscScalar,PetscBool,DSMatType,void*),void *ctx)
1024: {
1026:   PetscTryMethod(ds,"DSNEPSetComputeMatrixFunction_C",(DS,PetscErrorCode (*)(DS,PetscScalar,PetscBool,DSMatType,void*),void*),(ds,fun,ctx));
1027:   return 0;
1028: }

1030: static PetscErrorCode DSNEPGetComputeMatrixFunction_NEP(DS ds,PetscErrorCode (**fun)(DS,PetscScalar,PetscBool,DSMatType,void*),void **ctx)
1031: {
1032:   DS_NEP *dsctx = (DS_NEP*)ds->data;

1034:   if (fun) *fun = dsctx->computematrix;
1035:   if (ctx) *ctx = dsctx->computematrixctx;
1036:   return 0;
1037: }

1039: /*@C
1040:    DSNEPGetComputeMatrixFunction - Returns the user-provided callback function
1041:    set in DSNEPSetComputeMatrixFunction().

1043:    Not Collective

1045:    Input Parameter:
1046: .  ds  - the direct solver context

1048:    Output Parameters:
1049: +  fun - the pointer to the user function
1050: -  ctx - the context pointer

1052:    Level: developer

1054: .seealso: DSNEPSetComputeMatrixFunction()
1055: @*/
1056: PetscErrorCode DSNEPGetComputeMatrixFunction(DS ds,PetscErrorCode (**fun)(DS,PetscScalar,PetscBool,DSMatType,void*),void **ctx)
1057: {
1059:   PetscUseMethod(ds,"DSNEPGetComputeMatrixFunction_C",(DS,PetscErrorCode (**)(DS,PetscScalar,PetscBool,DSMatType,void*),void**),(ds,fun,ctx));
1060:   return 0;
1061: }

1063: static PetscErrorCode DSNEPSetRG_NEP(DS ds,RG rg)
1064: {
1065:   DS_NEP         *dsctx = (DS_NEP*)ds->data;

1067:   PetscObjectReference((PetscObject)rg);
1068:   RGDestroy(&dsctx->rg);
1069:   dsctx->rg = rg;
1070:   return 0;
1071: }

1073: /*@
1074:    DSNEPSetRG - Associates a region object to the DSNEP solver.

1076:    Logically Collective on ds

1078:    Input Parameters:
1079: +  ds  - the direct solver context
1080: -  rg  - the region context

1082:    Notes:
1083:    The region is used only in the contour integral method, and
1084:    should enclose the wanted eigenvalues.

1086:    Level: developer

1088: .seealso: DSNEPGetRG()
1089: @*/
1090: PetscErrorCode DSNEPSetRG(DS ds,RG rg)
1091: {
1093:   if (rg) {
1096:   }
1097:   PetscTryMethod(ds,"DSNEPSetRG_C",(DS,RG),(ds,rg));
1098:   return 0;
1099: }

1101: static PetscErrorCode DSNEPGetRG_NEP(DS ds,RG *rg)
1102: {
1103:   DS_NEP         *ctx = (DS_NEP*)ds->data;

1105:   if (!ctx->rg) {
1106:     RGCreate(PetscObjectComm((PetscObject)ds),&ctx->rg);
1107:     PetscObjectIncrementTabLevel((PetscObject)ctx->rg,(PetscObject)ds,1);
1108:     RGSetOptionsPrefix(ctx->rg,((PetscObject)ds)->prefix);
1109:     RGAppendOptionsPrefix(ctx->rg,"ds_nep_");
1110:     PetscObjectSetOptions((PetscObject)ctx->rg,((PetscObject)ds)->options);
1111:   }
1112:   *rg = ctx->rg;
1113:   return 0;
1114: }

1116: /*@
1117:    DSNEPGetRG - Obtain the region object associated to the DSNEP solver.

1119:    Not Collective

1121:    Input Parameter:
1122: .  ds  - the direct solver context

1124:    Output Parameter:
1125: .  rg  - the region context

1127:    Level: developer

1129: .seealso: DSNEPSetRG()
1130: @*/
1131: PetscErrorCode DSNEPGetRG(DS ds,RG *rg)
1132: {
1135:   PetscUseMethod(ds,"DSNEPGetRG_C",(DS,RG*),(ds,rg));
1136:   return 0;
1137: }

1139: PetscErrorCode DSSetFromOptions_NEP(DS ds,PetscOptionItems *PetscOptionsObject)
1140: {
1141:   PetscInt       k;
1142:   PetscBool      flg;
1143: #if defined(PETSC_USE_COMPLEX)
1144:   PetscReal      r;
1145:   PetscBool      flg1;
1146:   DS_NEP         *ctx = (DS_NEP*)ds->data;
1147: #endif

1149:   PetscOptionsHeadBegin(PetscOptionsObject,"DS NEP Options");

1151:     PetscOptionsInt("-ds_nep_minimality","Maximum minimality index","DSNEPSetMinimality",4,&k,&flg);
1152:     if (flg) DSNEPSetMinimality(ds,k);

1154:     PetscOptionsInt("-ds_nep_integration_points","Number of integration points","DSNEPSetIntegrationPoints",64,&k,&flg);
1155:     if (flg) DSNEPSetIntegrationPoints(ds,k);

1157:     PetscOptionsInt("-ds_nep_sampling_size","Number of sampling columns","DSNEPSetSamplingSize",0,&k,&flg);
1158:     if (flg) DSNEPSetSamplingSize(ds,k);

1160: #if defined(PETSC_USE_COMPLEX)
1161:     r = ctx->rtol;
1162:     PetscOptionsReal("-ds_nep_refine_tol","Refinement tolerance","DSNEPSetRefine",ctx->rtol,&r,&flg1);
1163:     k = ctx->Nit;
1164:     PetscOptionsInt("-ds_nep_refine_its","Number of iterative refinement iterations","DSNEPSetRefine",ctx->Nit,&k,&flg);
1165:     if (flg1||flg) DSNEPSetRefine(ds,r,k);

1167:     if (ds->method==1) {
1168:       if (!ctx->rg) DSNEPGetRG(ds,&ctx->rg);
1169:       RGSetFromOptions(ctx->rg);
1170:     }
1171: #endif

1173:   PetscOptionsHeadEnd();
1174:   return 0;
1175: }

1177: PetscErrorCode DSDestroy_NEP(DS ds)
1178: {
1179:   DS_NEP         *ctx = (DS_NEP*)ds->data;
1180:   PetscInt       i;

1182:   for (i=0;i<ctx->nf;i++) FNDestroy(&ctx->f[i]);
1183:   RGDestroy(&ctx->rg);
1184:   PetscLayoutDestroy(&ctx->map);
1185:   PetscFree(ds->data);
1186:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetFN_C",NULL);
1187:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetFN_C",NULL);
1188:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetNumFN_C",NULL);
1189:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetMinimality_C",NULL);
1190:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetMinimality_C",NULL);
1191:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetRefine_C",NULL);
1192:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetRefine_C",NULL);
1193:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetIntegrationPoints_C",NULL);
1194:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetIntegrationPoints_C",NULL);
1195:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetSamplingSize_C",NULL);
1196:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetSamplingSize_C",NULL);
1197:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetRG_C",NULL);
1198:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetRG_C",NULL);
1199:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetComputeMatrixFunction_C",NULL);
1200:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetComputeMatrixFunction_C",NULL);
1201:   return 0;
1202: }

1204: PetscErrorCode DSMatGetSize_NEP(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
1205: {
1206:   DS_NEP *ctx = (DS_NEP*)ds->data;

1208:   *rows = ds->n;
1209:   if (t==DS_MAT_Q || t==DS_MAT_Z || t==DS_MAT_U || t==DS_MAT_V) *rows *= ctx->max_mid;
1210:   *cols = ds->n;
1211:   if (t==DS_MAT_Q || t==DS_MAT_Z || t==DS_MAT_U || t==DS_MAT_V || t==DS_MAT_X || t==DS_MAT_Y) *cols *= ctx->max_mid;
1212:   return 0;
1213: }

1215: /*MC
1216:    DSNEP - Dense Nonlinear Eigenvalue Problem.

1218:    Level: beginner

1220:    Notes:
1221:    The problem is expressed as T(lambda)*x = 0, where T(lambda) is a
1222:    parameter-dependent matrix written as T(lambda) = sum_i E_i*f_i(lambda).
1223:    The eigenvalues lambda are the arguments returned by DSSolve()..

1225:    The coefficient matrices E_i are the extra matrices of the DS, and
1226:    the scalar functions f_i are passed via DSNEPSetFN(). Optionally, a
1227:    callback function to fill the E_i matrices can be set with
1228:    DSNEPSetComputeMatrixFunction().

1230:    Used DS matrices:
1231: +  DS_MAT_Ex - coefficient matrices of the split form of T(lambda)
1232: .  DS_MAT_A  - (workspace) T(lambda) evaluated at a given lambda (SLP only)
1233: .  DS_MAT_B  - (workspace) T'(lambda) evaluated at a given lambda (SLP only)
1234: .  DS_MAT_Q  - (workspace) left Hankel matrix (contour only)
1235: .  DS_MAT_Z  - (workspace) right Hankel matrix (contour only)
1236: .  DS_MAT_U  - (workspace) left singular vectors (contour only)
1237: .  DS_MAT_V  - (workspace) right singular vectors (contour only)
1238: -  DS_MAT_W  - (workspace) auxiliary matrix of size nxn

1240:    Implemented methods:
1241: +  0 - Successive Linear Problems (SLP), computes just one eigenpair
1242: -  1 - Contour integral, computes all eigenvalues inside a region

1244: .seealso: DSCreate(), DSSetType(), DSType, DSNEPSetFN(), DSNEPSetComputeMatrixFunction()
1245: M*/
1246: SLEPC_EXTERN PetscErrorCode DSCreate_NEP(DS ds)
1247: {
1248:   DS_NEP         *ctx;

1250:   PetscNew(&ctx);
1251:   ds->data = (void*)ctx;
1252:   ctx->max_mid = 4;
1253:   ctx->nnod    = 64;
1254:   ctx->Nit     = 3;
1255:   ctx->rtol    = PETSC_MACHINE_EPSILON/PetscSqrtReal(PETSC_SQRT_MACHINE_EPSILON);

1257:   ds->ops->allocate       = DSAllocate_NEP;
1258:   ds->ops->setfromoptions = DSSetFromOptions_NEP;
1259:   ds->ops->view           = DSView_NEP;
1260:   ds->ops->vectors        = DSVectors_NEP;
1261:   ds->ops->solve[0]       = DSSolve_NEP_SLP;
1262: #if defined(PETSC_USE_COMPLEX)
1263:   ds->ops->solve[1]       = DSSolve_NEP_Contour;
1264: #endif
1265:   ds->ops->sort           = DSSort_NEP;
1266:   ds->ops->synchronize    = DSSynchronize_NEP;
1267:   ds->ops->destroy        = DSDestroy_NEP;
1268:   ds->ops->matgetsize     = DSMatGetSize_NEP;

1270:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetFN_C",DSNEPSetFN_NEP);
1271:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetFN_C",DSNEPGetFN_NEP);
1272:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetNumFN_C",DSNEPGetNumFN_NEP);
1273:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetMinimality_C",DSNEPGetMinimality_NEP);
1274:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetMinimality_C",DSNEPSetMinimality_NEP);
1275:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetRefine_C",DSNEPGetRefine_NEP);
1276:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetRefine_C",DSNEPSetRefine_NEP);
1277:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetIntegrationPoints_C",DSNEPGetIntegrationPoints_NEP);
1278:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetIntegrationPoints_C",DSNEPSetIntegrationPoints_NEP);
1279:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetSamplingSize_C",DSNEPGetSamplingSize_NEP);
1280:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetSamplingSize_C",DSNEPSetSamplingSize_NEP);
1281:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetRG_C",DSNEPSetRG_NEP);
1282:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetRG_C",DSNEPGetRG_NEP);
1283:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetComputeMatrixFunction_C",DSNEPSetComputeMatrixFunction_NEP);
1284:   PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetComputeMatrixFunction_C",DSNEPGetComputeMatrixFunction_NEP);
1285:   return 0;
1286: }