Actual source code: dsnep.c
slepc-3.18.3 2023-03-24
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 <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 max_mid; /* maximum minimality index */
18: PetscInt nnod; /* number of nodes for quadrature rules */
19: PetscInt spls; /* number of sampling columns for quadrature rules */
20: PetscInt Nit; /* number of refinement iterations */
21: PetscReal rtol; /* tolerance of Newton refinement */
22: RG rg; /* region for contour integral */
23: PetscLayout map; /* used to distribute work among MPI processes */
24: void *computematrixctx;
25: PetscErrorCode (*computematrix)(DS,PetscScalar,PetscBool,DSMatType,void*);
26: } DS_NEP;
28: /*
29: DSNEPComputeMatrix - Build the matrix associated with a nonlinear operator
30: T(lambda) or its derivative T'(lambda), given the parameter lambda, where
31: T(lambda) = sum_i E_i*f_i(lambda). The result is written in mat.
32: */
33: static PetscErrorCode DSNEPComputeMatrix(DS ds,PetscScalar lambda,PetscBool deriv,DSMatType mat)
34: {
35: DS_NEP *ctx = (DS_NEP*)ds->data;
36: PetscScalar *T,alpha;
37: const PetscScalar *E;
38: PetscInt i,ld,n;
39: PetscBLASInt k,inc=1;
41: PetscLogEventBegin(DS_Other,ds,0,0,0);
42: if (ctx->computematrix) (*ctx->computematrix)(ds,lambda,deriv,mat,ctx->computematrixctx);
43: else {
44: DSGetDimensions(ds,&n,NULL,NULL,NULL);
45: DSGetLeadingDimension(ds,&ld);
46: PetscBLASIntCast(ld*n,&k);
47: MatDenseGetArray(ds->omat[mat],&T);
48: PetscArrayzero(T,k);
49: for (i=0;i<ctx->nf;i++) {
50: if (deriv) FNEvaluateDerivative(ctx->f[i],lambda,&alpha);
51: else FNEvaluateFunction(ctx->f[i],lambda,&alpha);
52: MatDenseGetArrayRead(ds->omat[DSMatExtra[i]],&E);
53: PetscCallBLAS("BLASaxpy",BLASaxpy_(&k,&alpha,E,&inc,T,&inc));
54: MatDenseRestoreArrayRead(ds->omat[DSMatExtra[i]],&E);
55: }
56: MatDenseRestoreArray(ds->omat[mat],&T);
57: }
58: PetscLogEventEnd(DS_Other,ds,0,0,0);
59: return 0;
60: }
62: PetscErrorCode DSAllocate_NEP(DS ds,PetscInt ld)
63: {
64: DS_NEP *ctx = (DS_NEP*)ds->data;
65: PetscInt i;
67: DSAllocateMat_Private(ds,DS_MAT_X);
68: for (i=0;i<ctx->nf;i++) DSAllocateMat_Private(ds,DSMatExtra[i]);
69: PetscFree(ds->perm);
70: PetscMalloc1(ld*ctx->max_mid,&ds->perm);
71: return 0;
72: }
74: PetscErrorCode DSView_NEP(DS ds,PetscViewer viewer)
75: {
76: DS_NEP *ctx = (DS_NEP*)ds->data;
77: PetscViewerFormat format;
78: PetscInt i;
79: const char *methodname[] = {
80: "Successive Linear Problems",
81: "Contour Integral"
82: };
83: const int nmeth=PETSC_STATIC_ARRAY_LENGTH(methodname);
85: PetscViewerGetFormat(viewer,&format);
86: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
87: if (ds->method<nmeth) PetscViewerASCIIPrintf(viewer,"solving the problem with: %s\n",methodname[ds->method]);
88: #if defined(PETSC_USE_COMPLEX)
89: if (ds->method==1) { /* contour integral method */
90: PetscViewerASCIIPrintf(viewer,"number of integration points: %" PetscInt_FMT "\n",ctx->nnod);
91: PetscViewerASCIIPrintf(viewer,"maximum minimality index: %" PetscInt_FMT "\n",ctx->max_mid);
92: if (ctx->spls) PetscViewerASCIIPrintf(viewer,"number of sampling columns for quadrature: %" PetscInt_FMT "\n",ctx->spls);
93: if (ctx->Nit) PetscViewerASCIIPrintf(viewer,"doing iterative refinement (%" PetscInt_FMT " its, tolerance %g)\n",ctx->Nit,(double)ctx->rtol);
94: RGView(ctx->rg,viewer);
95: }
96: #endif
97: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscViewerASCIIPrintf(viewer,"number of functions: %" PetscInt_FMT "\n",ctx->nf);
98: return 0;
99: }
100: for (i=0;i<ctx->nf;i++) {
101: FNView(ctx->f[i],viewer);
102: DSViewMat(ds,viewer,DSMatExtra[i]);
103: }
104: if (ds->state>DS_STATE_INTERMEDIATE) DSViewMat(ds,viewer,DS_MAT_X);
105: return 0;
106: }
108: PetscErrorCode DSVectors_NEP(DS ds,DSMatType mat,PetscInt *j,PetscReal *rnorm)
109: {
111: switch (mat) {
112: case DS_MAT_X:
113: break;
114: case DS_MAT_Y:
115: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_SUP,"Not implemented yet");
116: default:
117: SETERRQ(PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mat parameter");
118: }
119: return 0;
120: }
122: PetscErrorCode DSSort_NEP(DS ds,PetscScalar *wr,PetscScalar *wi,PetscScalar *rr,PetscScalar *ri,PetscInt *dummy)
123: {
124: DS_NEP *ctx = (DS_NEP*)ds->data;
125: PetscInt n,l,i,*perm,lds;
126: PetscScalar *Q;
128: if (!ds->sc) return 0;
129: if (!ds->method) return 0; /* SLP computes just one eigenvalue */
130: n = ds->n*ctx->max_mid;
131: lds = ds->ld*ctx->max_mid;
132: l = ds->l;
133: perm = ds->perm;
134: for (i=0;i<n;i++) perm[i] = i;
135: if (rr) DSSortEigenvalues_Private(ds,rr,ri,perm,PETSC_FALSE);
136: else DSSortEigenvalues_Private(ds,wr,NULL,perm,PETSC_FALSE);
137: MatDenseGetArray(ds->omat[DS_MAT_Q],&Q);
138: for (i=l;i<ds->t;i++) Q[i+i*lds] = wr[perm[i]];
139: for (i=l;i<ds->t;i++) wr[i] = Q[i+i*lds];
140: MatDenseRestoreArray(ds->omat[DS_MAT_Q],&Q);
141: /* n != ds->n */
142: DSPermuteColumns_Private(ds,0,ds->t,ds->n,DS_MAT_X,perm);
143: return 0;
144: }
146: PetscErrorCode DSSolve_NEP_SLP(DS ds,PetscScalar *wr,PetscScalar *wi)
147: {
148: PetscScalar *A,*B,*W,*X,*work,*alpha,*beta;
149: PetscScalar sigma,lambda,mu,re,re2,sone=1.0,szero=0.0;
150: PetscBLASInt info,n,ld,lrwork=0,lwork,one=1,zero=0;
151: PetscInt it,pos,j,maxit=100,result;
152: PetscReal norm,tol,done=1.0;
153: #if defined(PETSC_USE_COMPLEX)
154: PetscReal *rwork;
155: #else
156: PetscReal *alphai,im,im2;
157: #endif
159: PetscBLASIntCast(ds->n,&n);
160: PetscBLASIntCast(ds->ld,&ld);
161: #if defined(PETSC_USE_COMPLEX)
162: PetscBLASIntCast(2*ds->n+2*ds->n,&lwork);
163: PetscBLASIntCast(8*ds->n,&lrwork);
164: #else
165: PetscBLASIntCast(3*ds->n+8*ds->n,&lwork);
166: #endif
167: DSAllocateWork_Private(ds,lwork,lrwork,0);
168: alpha = ds->work;
169: beta = ds->work + ds->n;
170: #if defined(PETSC_USE_COMPLEX)
171: work = ds->work + 2*ds->n;
172: lwork -= 2*ds->n;
173: #else
174: alphai = ds->work + 2*ds->n;
175: work = ds->work + 3*ds->n;
176: lwork -= 3*ds->n;
177: #endif
178: DSAllocateMat_Private(ds,DS_MAT_A);
179: DSAllocateMat_Private(ds,DS_MAT_B);
180: DSAllocateMat_Private(ds,DS_MAT_W);
181: MatDenseGetArray(ds->omat[DS_MAT_A],&A);
182: MatDenseGetArray(ds->omat[DS_MAT_B],&B);
183: MatDenseGetArray(ds->omat[DS_MAT_W],&W);
184: MatDenseGetArray(ds->omat[DS_MAT_X],&X);
186: sigma = 0.0;
187: if (ds->sc->comparison==SlepcCompareTargetMagnitude || ds->sc->comparison==SlepcCompareTargetReal) sigma = *(PetscScalar*)ds->sc->comparisonctx;
188: lambda = sigma;
189: tol = n*PETSC_MACHINE_EPSILON/PetscSqrtReal(PETSC_SQRT_MACHINE_EPSILON);
191: for (it=0;it<maxit;it++) {
193: /* evaluate T and T' */
194: DSNEPComputeMatrix(ds,lambda,PETSC_FALSE,DS_MAT_A);
195: if (it) {
196: PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,A,&ld,X,&one,&szero,X+ld,&one));
197: norm = BLASnrm2_(&n,X+ld,&one);
198: if (norm/PetscAbsScalar(lambda)<=tol) break;
199: }
200: DSNEPComputeMatrix(ds,lambda,PETSC_TRUE,DS_MAT_B);
202: /* compute eigenvalue correction mu and eigenvector u */
203: #if defined(PETSC_USE_COMPLEX)
204: rwork = ds->rwork;
205: PetscCallBLAS("LAPACKggev",LAPACKggev_("N","V",&n,A,&ld,B,&ld,alpha,beta,NULL,&ld,W,&ld,work,&lwork,rwork,&info));
206: #else
207: PetscCallBLAS("LAPACKggev",LAPACKggev_("N","V",&n,A,&ld,B,&ld,alpha,alphai,beta,NULL,&ld,W,&ld,work,&lwork,&info));
208: #endif
209: SlepcCheckLapackInfo("ggev",info);
211: /* find smallest eigenvalue */
212: j = 0;
213: if (beta[j]==0.0) re = (PetscRealPart(alpha[j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
214: else re = alpha[j]/beta[j];
215: #if !defined(PETSC_USE_COMPLEX)
216: if (beta[j]==0.0) im = (alphai[j]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
217: else im = alphai[j]/beta[j];
218: #endif
219: pos = 0;
220: for (j=1;j<n;j++) {
221: if (beta[j]==0.0) re2 = (PetscRealPart(alpha[j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
222: else re2 = alpha[j]/beta[j];
223: #if !defined(PETSC_USE_COMPLEX)
224: if (beta[j]==0.0) im2 = (alphai[j]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
225: else im2 = alphai[j]/beta[j];
226: SlepcCompareSmallestMagnitude(re,im,re2,im2,&result,NULL);
227: #else
228: SlepcCompareSmallestMagnitude(re,0.0,re2,0.0,&result,NULL);
229: #endif
230: if (result > 0) {
231: re = re2;
232: #if !defined(PETSC_USE_COMPLEX)
233: im = im2;
234: #endif
235: pos = j;
236: }
237: }
239: #if !defined(PETSC_USE_COMPLEX)
241: #endif
242: mu = alpha[pos]/beta[pos];
243: PetscArraycpy(X,W+pos*ld,n);
244: norm = BLASnrm2_(&n,X,&one);
245: PetscCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&one,X,&n,&info));
246: SlepcCheckLapackInfo("lascl",info);
248: /* correct eigenvalue approximation */
249: lambda = lambda - mu;
250: }
251: MatDenseRestoreArray(ds->omat[DS_MAT_A],&A);
252: MatDenseRestoreArray(ds->omat[DS_MAT_B],&B);
253: MatDenseRestoreArray(ds->omat[DS_MAT_W],&W);
254: MatDenseRestoreArray(ds->omat[DS_MAT_X],&X);
257: ds->t = 1;
258: wr[0] = lambda;
259: if (wi) wi[0] = 0.0;
260: return 0;
261: }
263: #if defined(PETSC_USE_COMPLEX)
264: /*
265: Newton refinement for eigenpairs computed with contour integral.
266: k - number of eigenpairs to refine
267: wr - eigenvalues (eigenvectors are stored in DS_MAT_X)
268: */
269: static PetscErrorCode DSNEPNewtonRefine(DS ds,PetscInt k,PetscScalar *wr)
270: {
271: DS_NEP *ctx = (DS_NEP*)ds->data;
272: PetscScalar *X,*W,*U,*R,sone=1.0,szero=0.0;
273: PetscReal norm;
274: PetscInt i,j,ii,nwu=0,*p,jstart=0,jend=k;
275: const PetscInt *range;
276: PetscBLASInt n,*perm,info,ld,one=1,n1;
277: PetscMPIInt len,size,root;
278: PetscLayout map;
280: MatDenseGetArray(ds->omat[DS_MAT_X],&X);
281: PetscBLASIntCast(ds->n,&n);
282: PetscBLASIntCast(ds->ld,&ld);
283: n1 = n+1;
284: p = ds->perm;
285: PetscArrayzero(p,k);
286: DSAllocateWork_Private(ds,(n+2)*(n+1),0,n+1);
287: U = ds->work+nwu; nwu += (n+1)*(n+1);
288: R = ds->work+nwu; /*nwu += n+1;*/
289: perm = ds->iwork;
290: if (ds->pmode==DS_PARALLEL_DISTRIBUTED) {
291: PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)ds),PETSC_DECIDE,k,1,&map);
292: PetscLayoutGetRange(map,&jstart,&jend);
293: }
294: for (ii=0;ii<ctx->Nit;ii++) {
295: for (j=jstart;j<jend;j++) {
296: if (p[j]<2) {
297: DSNEPComputeMatrix(ds,wr[j],PETSC_FALSE,DS_MAT_W);
298: MatDenseGetArray(ds->omat[DS_MAT_W],&W);
299: PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,W,&ld,X+ld*j,&one,&szero,R,&one));
300: MatDenseRestoreArray(ds->omat[DS_MAT_W],&W);
301: norm = BLASnrm2_(&n,R,&one);
302: if (norm/PetscAbsScalar(wr[j]) > ctx->rtol) {
303: PetscInfo(NULL,"Refining eigenpair %" PetscInt_FMT ", residual=%g\n",j,(double)(norm/PetscAbsScalar(wr[j])));
304: p[j] = 1;
305: R[n] = 0.0;
306: MatDenseGetArray(ds->omat[DS_MAT_W],&W);
307: for (i=0;i<n;i++) {
308: PetscArraycpy(U+i*n1,W+i*ld,n);
309: U[n+i*n1] = PetscConj(X[j*ld+i]);
310: }
311: MatDenseRestoreArray(ds->omat[DS_MAT_W],&W);
312: U[n+n*n1] = 0.0;
313: DSNEPComputeMatrix(ds,wr[j],PETSC_TRUE,DS_MAT_W);
314: MatDenseGetArray(ds->omat[DS_MAT_W],&W);
315: PetscCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,W,&ld,X+ld*j,&one,&szero,U+n*(n+1),&one));
316: MatDenseRestoreArray(ds->omat[DS_MAT_W],&W);
317: /* solve system */
318: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
319: PetscCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n1,&n1,U,&n1,perm,&info));
320: SlepcCheckLapackInfo("getrf",info);
321: PetscCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&n1,&one,U,&n1,perm,R,&n1,&info));
322: SlepcCheckLapackInfo("getrs",info);
323: PetscFPTrapPop();
324: wr[j] -= R[n];
325: for (i=0;i<n;i++) X[j*ld+i] -= R[i];
326: /* normalization */
327: norm = BLASnrm2_(&n,X+ld*j,&one);
328: for (i=0;i<n;i++) X[ld*j+i] /= norm;
329: } else p[j] = 2;
330: }
331: }
332: }
333: if (ds->pmode==DS_PARALLEL_DISTRIBUTED) { /* communicate results */
334: PetscMPIIntCast(k,&len);
335: MPIU_Allreduce(MPI_IN_PLACE,p,len,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)ds));
336: MPI_Comm_size(PetscObjectComm((PetscObject)ds),&size);
337: PetscLayoutGetRanges(map,&range);
338: for (j=0;j<k;j++) {
339: if (p[j]) { /* j-th eigenpair has been refined */
340: for (root=0;root<size;root++) if (range[root+1]>j) break;
341: PetscMPIIntCast(1,&len);
342: MPI_Bcast(wr+j,len,MPIU_SCALAR,root,PetscObjectComm((PetscObject)ds));
343: PetscMPIIntCast(n,&len);
344: MPI_Bcast(X+ld*j,len,MPIU_SCALAR,root,PetscObjectComm((PetscObject)ds));
345: }
346: }
347: PetscLayoutDestroy(&map);
348: }
349: MatDenseRestoreArray(ds->omat[DS_MAT_X],&X);
350: return 0;
351: }
353: PetscErrorCode DSSolve_NEP_Contour(DS ds,PetscScalar *wr,PetscScalar *wi)
354: {
355: DS_NEP *ctx = (DS_NEP*)ds->data;
356: PetscScalar *alpha,*beta,*Q,*Z,*X,*U,*V,*W,*work,*Rc,*R,*w,*z,*zn,*S;
357: PetscScalar sone=1.0,szero=0.0,center;
358: PetscReal *rwork,norm,radius,vscale,rgscale,*sigma;
359: PetscBLASInt info,n,*perm,p,pp,ld,lwork,k_,rk_,colA,rowA,one=1;
360: PetscInt mid,lds,nnod=ctx->nnod,k,i,ii,jj,j,s,off,rk,nwu=0,nw,lrwork,*inside,kstart=0,kend=nnod;
361: PetscMPIInt len;
362: PetscBool isellipse;
363: PetscRandom rand;
366: /* Contour parameters */
367: PetscObjectTypeCompare((PetscObject)ctx->rg,RGELLIPSE,&isellipse);
369: RGEllipseGetParameters(ctx->rg,¢er,&radius,&vscale);
370: RGGetScale(ctx->rg,&rgscale);
371: if (ds->pmode==DS_PARALLEL_DISTRIBUTED) {
372: if (!ctx->map) PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)ds),PETSC_DECIDE,ctx->nnod,1,&ctx->map);
373: PetscLayoutGetRange(ctx->map,&kstart,&kend);
374: }
376: DSAllocateMat_Private(ds,DS_MAT_W); /* size n */
377: DSAllocateMat_Private(ds,DS_MAT_Q); /* size mid*n */
378: DSAllocateMat_Private(ds,DS_MAT_Z); /* size mid*n */
379: DSAllocateMat_Private(ds,DS_MAT_U); /* size mid*n */
380: DSAllocateMat_Private(ds,DS_MAT_V); /* size mid*n */
381: MatDenseGetArray(ds->omat[DS_MAT_Q],&Q);
382: MatDenseGetArray(ds->omat[DS_MAT_Z],&Z);
383: MatDenseGetArray(ds->omat[DS_MAT_U],&U);
384: MatDenseGetArray(ds->omat[DS_MAT_V],&V);
385: mid = ctx->max_mid;
386: PetscBLASIntCast(ds->n,&n);
387: p = n; /* maximum number of columns for the probing matrix */
388: PetscBLASIntCast(ds->ld,&ld);
389: PetscBLASIntCast(mid*n,&rowA);
390: PetscBLASIntCast(5*rowA,&lwork);
391: nw = n*(2*p+7*mid)+3*nnod+2*mid*n*p;
392: lrwork = mid*n*6+8*n;
393: DSAllocateWork_Private(ds,nw,lrwork,n+1);
395: sigma = ds->rwork;
396: rwork = ds->rwork+mid*n;
397: perm = ds->iwork;
398: z = ds->work+nwu; nwu += nnod; /* quadrature points */
399: zn = ds->work+nwu; nwu += nnod; /* normalized quadrature points */
400: w = ds->work+nwu; nwu += nnod; /* quadrature weights */
401: Rc = ds->work+nwu; nwu += n*p;
402: R = ds->work+nwu; nwu += n*p;
403: alpha = ds->work+nwu; nwu += mid*n;
404: beta = ds->work+nwu; nwu += mid*n;
405: S = ds->work+nwu; nwu += 2*mid*n*p;
406: work = ds->work+nwu; /*nwu += mid*n*5;*/
408: /* Compute quadrature parameters */
409: RGComputeQuadrature(ctx->rg,RG_QUADRULE_TRAPEZOIDAL,nnod,z,zn,w);
411: /* Set random matrix */
412: PetscRandomCreate(PetscObjectComm((PetscObject)ds),&rand);
413: PetscRandomSetSeed(rand,0x12345678);
414: PetscRandomSeed(rand);
415: for (j=0;j<p;j++)
416: for (i=0;i<n;i++) PetscRandomGetValue(rand,Rc+i+j*n);
417: PetscArrayzero(S,2*mid*n*p);
418: /* Loop of integration points */
419: for (k=kstart;k<kend;k++) {
420: PetscInfo(NULL,"Solving integration point %" PetscInt_FMT "\n",k);
421: PetscArraycpy(R,Rc,p*n);
422: DSNEPComputeMatrix(ds,z[k],PETSC_FALSE,DS_MAT_W);
424: /* LU factorization */
425: MatDenseGetArray(ds->omat[DS_MAT_W],&W);
426: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
427: PetscCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,W,&ld,perm,&info));
428: SlepcCheckLapackInfo("getrf",info);
429: PetscCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&n,&p,W,&ld,perm,R,&n,&info));
430: SlepcCheckLapackInfo("getrs",info);
431: PetscFPTrapPop();
432: MatDenseRestoreArray(ds->omat[DS_MAT_W],&W);
434: /* Moments computation */
435: for (s=0;s<2*ctx->max_mid;s++) {
436: off = s*n*p;
437: for (j=0;j<p;j++)
438: for (i=0;i<n;i++) S[off+i+j*n] += w[k]*R[j*n+i];
439: w[k] *= zn[k];
440: }
441: }
443: if (ds->pmode==DS_PARALLEL_DISTRIBUTED) { /* compute final S via reduction */
444: PetscMPIIntCast(2*mid*n*p,&len);
445: MPIU_Allreduce(MPI_IN_PLACE,S,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)ds));
446: }
447: p = ctx->spls?PetscMin(ctx->spls,n):n;
448: pp = p;
449: do {
450: p = pp;
451: PetscBLASIntCast(mid*p,&colA);
453: PetscInfo(ds,"Computing SVD of size %" PetscBLASInt_FMT "x%" PetscBLASInt_FMT "\n",rowA,colA);
454: for (jj=0;jj<mid;jj++) {
455: for (ii=0;ii<mid;ii++) {
456: off = jj*p*rowA+ii*n;
457: for (j=0;j<p;j++)
458: for (i=0;i<n;i++) Q[off+j*rowA+i] = S[((jj+ii)*n+j)*n+i];
459: }
460: }
461: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
462: PetscCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&rowA,&colA,Q,&rowA,sigma,U,&rowA,V,&colA,work,&lwork,rwork,&info));
463: SlepcCheckLapackInfo("gesvd",info);
464: PetscFPTrapPop();
466: rk = colA;
467: for (i=1;i<colA;i++) if (sigma[i]/sigma[0]<PETSC_MACHINE_EPSILON*1e4) {rk = i; break;}
468: if (rk<colA || p==n) break;
469: pp *= 2;
470: } while (pp<=n);
471: PetscInfo(ds,"Solving generalized eigenproblem of size %" PetscInt_FMT "\n",rk);
472: for (jj=0;jj<mid;jj++) {
473: for (ii=0;ii<mid;ii++) {
474: off = jj*p*rowA+ii*n;
475: for (j=0;j<p;j++)
476: for (i=0;i<n;i++) Q[off+j*rowA+i] = S[((jj+ii+1)*n+j)*n+i];
477: }
478: }
479: PetscBLASIntCast(rk,&rk_);
480: PetscCallBLAS("BLASgemm",BLASgemm_("N","C",&rowA,&rk_,&colA,&sone,Q,&rowA,V,&colA,&szero,Z,&rowA));
481: PetscCallBLAS("BLASgemm",BLASgemm_("C","N",&rk_,&rk_,&rowA,&sone,U,&rowA,Z,&rowA,&szero,Q,&rk_));
482: PetscArrayzero(Z,n*mid*n*mid);
483: for (j=0;j<rk;j++) Z[j+j*rk_] = sigma[j];
484: PetscCallBLAS("LAPACKggev",LAPACKggev_("N","V",&rk_,Q,&rk_,Z,&rk_,alpha,beta,NULL,&ld,V,&rk_,work,&lwork,rwork,&info));
485: for (i=0;i<rk;i++) wr[i] = (center+radius*alpha[i]/beta[i])*rgscale;
486: PetscMalloc1(rk,&inside);
487: RGCheckInside(ctx->rg,rk,wr,wi,inside);
488: k=0;
489: for (i=0;i<rk;i++)
490: if (inside[i]==1) inside[k++] = i;
491: /* Discard values outside region */
492: lds = ld*mid;
493: PetscArrayzero(Q,lds*lds);
494: PetscArrayzero(Z,lds*lds);
495: for (i=0;i<k;i++) Q[i+i*lds] = (center*beta[inside[i]]+radius*alpha[inside[i]])*rgscale;
496: for (i=0;i<k;i++) Z[i+i*lds] = beta[inside[i]];
497: for (i=0;i<k;i++) wr[i] = Q[i+i*lds]/Z[i+i*lds];
498: for (j=0;j<k;j++) for (i=0;i<rk;i++) V[j*rk+i] = sigma[i]*V[inside[j]*rk+i];
499: PetscBLASIntCast(k,&k_);
500: MatDenseGetArray(ds->omat[DS_MAT_X],&X);
501: PetscCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&k_,&rk_,&sone,U,&rowA,V,&rk_,&szero,X,&ld));
502: /* Normalize */
503: for (j=0;j<k;j++) {
504: norm = BLASnrm2_(&n,X+ld*j,&one);
505: for (i=0;i<n;i++) X[ld*j+i] /= norm;
506: }
507: PetscFree(inside);
508: MatDenseRestoreArray(ds->omat[DS_MAT_X],&X);
509: MatDenseRestoreArray(ds->omat[DS_MAT_Q],&Q);
510: MatDenseRestoreArray(ds->omat[DS_MAT_Z],&Z);
511: MatDenseRestoreArray(ds->omat[DS_MAT_U],&U);
512: MatDenseRestoreArray(ds->omat[DS_MAT_V],&V);
514: /* Newton refinement */
515: if (ctx->Nit) DSNEPNewtonRefine(ds,k,wr);
516: ds->t = k;
517: PetscRandomDestroy(&rand);
518: return 0;
519: }
520: #endif
522: #if !defined(PETSC_HAVE_MPIUNI)
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: }
565: #endif
567: static PetscErrorCode DSNEPSetFN_NEP(DS ds,PetscInt n,FN fn[])
568: {
569: DS_NEP *ctx = (DS_NEP*)ds->data;
570: PetscInt i;
574: if (ds->ld) PetscInfo(ds,"DSNEPSetFN() called after DSAllocate()\n");
575: for (i=0;i<n;i++) PetscObjectReference((PetscObject)fn[i]);
576: for (i=0;i<ctx->nf;i++) FNDestroy(&ctx->f[i]);
577: for (i=0;i<n;i++) ctx->f[i] = fn[i];
578: ctx->nf = n;
579: return 0;
580: }
582: /*@
583: DSNEPSetFN - Sets a number of functions that define the nonlinear
584: eigenproblem.
586: Collective on ds
588: Input Parameters:
589: + ds - the direct solver context
590: . n - number of functions
591: - fn - array of functions
593: Notes:
594: The nonlinear eigenproblem is defined in terms of the split nonlinear
595: operator T(lambda) = sum_i A_i*f_i(lambda).
597: This function must be called before DSAllocate(). Then DSAllocate()
598: will allocate an extra matrix A_i per each function, that can be
599: filled in the usual way.
601: Level: advanced
603: .seealso: DSNEPGetFN(), DSAllocate()
604: @*/
605: PetscErrorCode DSNEPSetFN(DS ds,PetscInt n,FN fn[])
606: {
607: PetscInt i;
612: for (i=0;i<n;i++) {
615: }
616: PetscTryMethod(ds,"DSNEPSetFN_C",(DS,PetscInt,FN[]),(ds,n,fn));
617: return 0;
618: }
620: static PetscErrorCode DSNEPGetFN_NEP(DS ds,PetscInt k,FN *fn)
621: {
622: DS_NEP *ctx = (DS_NEP*)ds->data;
625: *fn = ctx->f[k];
626: return 0;
627: }
629: /*@
630: DSNEPGetFN - Gets the functions associated with the nonlinear DS.
632: Not collective, though parallel FNs are returned if the DS is parallel
634: Input Parameters:
635: + ds - the direct solver context
636: - k - the index of the requested function (starting in 0)
638: Output Parameter:
639: . fn - the function
641: Level: advanced
643: .seealso: DSNEPSetFN()
644: @*/
645: PetscErrorCode DSNEPGetFN(DS ds,PetscInt k,FN *fn)
646: {
649: PetscUseMethod(ds,"DSNEPGetFN_C",(DS,PetscInt,FN*),(ds,k,fn));
650: return 0;
651: }
653: static PetscErrorCode DSNEPGetNumFN_NEP(DS ds,PetscInt *n)
654: {
655: DS_NEP *ctx = (DS_NEP*)ds->data;
657: *n = ctx->nf;
658: return 0;
659: }
661: /*@
662: DSNEPGetNumFN - Returns the number of functions stored internally by
663: the DS.
665: Not collective
667: Input Parameter:
668: . ds - the direct solver context
670: Output Parameters:
671: . n - the number of functions passed in DSNEPSetFN()
673: Level: advanced
675: .seealso: DSNEPSetFN()
676: @*/
677: PetscErrorCode DSNEPGetNumFN(DS ds,PetscInt *n)
678: {
681: PetscUseMethod(ds,"DSNEPGetNumFN_C",(DS,PetscInt*),(ds,n));
682: return 0;
683: }
685: static PetscErrorCode DSNEPSetMinimality_NEP(DS ds,PetscInt n)
686: {
687: DS_NEP *ctx = (DS_NEP*)ds->data;
689: if (n == PETSC_DECIDE || n == PETSC_DEFAULT) ctx->max_mid = 4;
690: else {
692: ctx->max_mid = n;
693: }
694: return 0;
695: }
697: /*@
698: DSNEPSetMinimality - Sets the maximum minimality index used internally by
699: the DSNEP.
701: Logically Collective on ds
703: Input Parameters:
704: + ds - the direct solver context
705: - n - the maximum minimality index
707: Options Database Key:
708: . -ds_nep_minimality <n> - sets the maximum minimality index
710: Notes:
711: The maximum minimality index is used only in the contour integral method,
712: and is related to the highest momemts used in the method. The default
713: value is 1, an larger value might give better accuracy in some cases, but
714: at a higher cost.
716: Level: advanced
718: .seealso: DSNEPGetMinimality()
719: @*/
720: PetscErrorCode DSNEPSetMinimality(DS ds,PetscInt n)
721: {
724: PetscTryMethod(ds,"DSNEPSetMinimality_C",(DS,PetscInt),(ds,n));
725: return 0;
726: }
728: static PetscErrorCode DSNEPGetMinimality_NEP(DS ds,PetscInt *n)
729: {
730: DS_NEP *ctx = (DS_NEP*)ds->data;
732: *n = ctx->max_mid;
733: return 0;
734: }
736: /*@
737: DSNEPGetMinimality - Returns the maximum minimality index used internally by
738: the DSNEP.
740: Not collective
742: Input Parameter:
743: . ds - the direct solver context
745: Output Parameters:
746: . n - the maximum minimality index passed in DSNEPSetMinimality()
748: Level: advanced
750: .seealso: DSNEPSetMinimality()
751: @*/
752: PetscErrorCode DSNEPGetMinimality(DS ds,PetscInt *n)
753: {
756: PetscUseMethod(ds,"DSNEPGetMinimality_C",(DS,PetscInt*),(ds,n));
757: return 0;
758: }
760: static PetscErrorCode DSNEPSetRefine_NEP(DS ds,PetscReal tol,PetscInt its)
761: {
762: DS_NEP *ctx = (DS_NEP*)ds->data;
764: if (tol == PETSC_DEFAULT) ctx->rtol = PETSC_MACHINE_EPSILON/PetscSqrtReal(PETSC_SQRT_MACHINE_EPSILON);
765: else {
767: ctx->rtol = tol;
768: }
769: if (its == PETSC_DECIDE || its == PETSC_DEFAULT) ctx->Nit = 3;
770: else {
772: ctx->Nit = its;
773: }
774: return 0;
775: }
777: /*@
778: DSNEPSetRefine - Sets the tolerance and the number of iterations of Newton iterative
779: refinement for eigenpairs.
781: Logically Collective on ds
783: Input Parameters:
784: + ds - the direct solver context
785: . tol - the tolerance
786: - its - the number of iterations
788: Options Database Key:
789: + -ds_nep_refine_tol <tol> - sets the tolerance
790: - -ds_nep_refine_its <its> - sets the number of Newton iterations
792: Notes:
793: Iterative refinement of eigenpairs is currently used only in the contour
794: integral method.
796: Level: advanced
798: .seealso: DSNEPGetRefine()
799: @*/
800: PetscErrorCode DSNEPSetRefine(DS ds,PetscReal tol,PetscInt its)
801: {
805: PetscTryMethod(ds,"DSNEPSetRefine_C",(DS,PetscReal,PetscInt),(ds,tol,its));
806: return 0;
807: }
809: static PetscErrorCode DSNEPGetRefine_NEP(DS ds,PetscReal *tol,PetscInt *its)
810: {
811: DS_NEP *ctx = (DS_NEP*)ds->data;
813: if (tol) *tol = ctx->rtol;
814: if (its) *its = ctx->Nit;
815: return 0;
816: }
818: /*@
819: DSNEPGetRefine - Returns the tolerance and the number of iterations of Newton iterative
820: refinement for eigenpairs.
822: Not collective
824: Input Parameter:
825: . ds - the direct solver context
827: Output Parameters:
828: + tol - the tolerance
829: - its - the number of iterations
831: Level: advanced
833: .seealso: DSNEPSetRefine()
834: @*/
835: PetscErrorCode DSNEPGetRefine(DS ds,PetscReal *tol,PetscInt *its)
836: {
838: PetscUseMethod(ds,"DSNEPGetRefine_C",(DS,PetscReal*,PetscInt*),(ds,tol,its));
839: return 0;
840: }
842: static PetscErrorCode DSNEPSetIntegrationPoints_NEP(DS ds,PetscInt ip)
843: {
844: DS_NEP *ctx = (DS_NEP*)ds->data;
846: if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) ctx->nnod = 64;
847: else {
849: ctx->nnod = ip;
850: }
851: PetscLayoutDestroy(&ctx->map); /* need to redistribute at next solve */
852: return 0;
853: }
855: /*@
856: DSNEPSetIntegrationPoints - Sets the number of integration points to be
857: used in the contour integral method.
859: Logically Collective on ds
861: Input Parameters:
862: + ds - the direct solver context
863: - ip - the number of integration points
865: Options Database Key:
866: . -ds_nep_integration_points <ip> - sets the number of integration points
868: Notes:
869: This parameter is relevant only in the contour integral method.
871: Level: advanced
873: .seealso: DSNEPGetIntegrationPoints()
874: @*/
875: PetscErrorCode DSNEPSetIntegrationPoints(DS ds,PetscInt ip)
876: {
879: PetscTryMethod(ds,"DSNEPSetIntegrationPoints_C",(DS,PetscInt),(ds,ip));
880: return 0;
881: }
883: static PetscErrorCode DSNEPGetIntegrationPoints_NEP(DS ds,PetscInt *ip)
884: {
885: DS_NEP *ctx = (DS_NEP*)ds->data;
887: *ip = ctx->nnod;
888: return 0;
889: }
891: /*@
892: DSNEPGetIntegrationPoints - Returns the number of integration points used
893: in the contour integral method.
895: Not collective
897: Input Parameter:
898: . ds - the direct solver context
900: Output Parameters:
901: . ip - the number of integration points
903: Level: advanced
905: .seealso: DSNEPSetIntegrationPoints()
906: @*/
907: PetscErrorCode DSNEPGetIntegrationPoints(DS ds,PetscInt *ip)
908: {
911: PetscUseMethod(ds,"DSNEPGetIntegrationPoints_C",(DS,PetscInt*),(ds,ip));
912: return 0;
913: }
915: static PetscErrorCode DSNEPSetSamplingSize_NEP(DS ds,PetscInt p)
916: {
917: DS_NEP *ctx = (DS_NEP*)ds->data;
919: if (p == PETSC_DECIDE || p == PETSC_DEFAULT) ctx->spls = 0;
920: else {
922: ctx->spls = p;
923: }
924: return 0;
925: }
927: /*@
928: DSNEPSetSamplingSize - Sets the number of sampling columns to be
929: used in the contour integral method.
931: Logically Collective on ds
933: Input Parameters:
934: + ds - the direct solver context
935: - p - the number of columns for the sampling matrix
937: Options Database Key:
938: . -ds_nep_sampling_size <p> - sets the number of sampling columns
940: Notes:
941: This parameter is relevant only in the contour integral method.
943: Level: advanced
945: .seealso: DSNEPGetSamplingSize()
946: @*/
947: PetscErrorCode DSNEPSetSamplingSize(DS ds,PetscInt p)
948: {
951: PetscTryMethod(ds,"DSNEPSetSamplingSize_C",(DS,PetscInt),(ds,p));
952: return 0;
953: }
955: static PetscErrorCode DSNEPGetSamplingSize_NEP(DS ds,PetscInt *p)
956: {
957: DS_NEP *ctx = (DS_NEP*)ds->data;
959: *p = ctx->spls;
960: return 0;
961: }
963: /*@
964: DSNEPGetSamplingSize - Returns the number of sampling columns used
965: in the contour integral method.
967: Not collective
969: Input Parameter:
970: . ds - the direct solver context
972: Output Parameters:
973: . p - the number of columns for the sampling matrix
975: Level: advanced
977: .seealso: DSNEPSetSamplingSize()
978: @*/
979: PetscErrorCode DSNEPGetSamplingSize(DS ds,PetscInt *p)
980: {
983: PetscUseMethod(ds,"DSNEPGetSamplingSize_C",(DS,PetscInt*),(ds,p));
984: return 0;
985: }
987: static PetscErrorCode DSNEPSetComputeMatrixFunction_NEP(DS ds,PetscErrorCode (*fun)(DS,PetscScalar,PetscBool,DSMatType,void*),void *ctx)
988: {
989: DS_NEP *dsctx = (DS_NEP*)ds->data;
991: dsctx->computematrix = fun;
992: dsctx->computematrixctx = ctx;
993: return 0;
994: }
996: /*@C
997: DSNEPSetComputeMatrixFunction - Sets a user-provided subroutine to compute
998: the matrices T(lambda) or T'(lambda).
1000: Logically Collective on ds
1002: Input Parameters:
1003: + ds - the direct solver context
1004: . fun - a pointer to the user function
1005: - ctx - a context pointer (the last parameter to the user function)
1007: Calling Sequence of fun:
1008: $ fun(DS ds,PetscScalar lambda,PetscBool deriv,DSMatType mat,void *ctx)
1010: + ds - the direct solver object
1011: . lambda - point where T(lambda) or T'(lambda) must be evaluated
1012: . deriv - if true compute T'(lambda), otherwise compute T(lambda)
1013: . mat - the DS matrix where the result must be stored
1014: - ctx - optional context, as set by DSNEPSetComputeMatrixFunction()
1016: Note:
1017: The result is computed as T(lambda) = sum_i E_i*f_i(lambda), and similarly
1018: for the derivative.
1020: Level: developer
1022: .seealso: DSNEPGetComputeMatrixFunction()
1023: @*/
1024: PetscErrorCode DSNEPSetComputeMatrixFunction(DS ds,PetscErrorCode (*fun)(DS,PetscScalar,PetscBool,DSMatType,void*),void *ctx)
1025: {
1027: PetscTryMethod(ds,"DSNEPSetComputeMatrixFunction_C",(DS,PetscErrorCode (*)(DS,PetscScalar,PetscBool,DSMatType,void*),void*),(ds,fun,ctx));
1028: return 0;
1029: }
1031: static PetscErrorCode DSNEPGetComputeMatrixFunction_NEP(DS ds,PetscErrorCode (**fun)(DS,PetscScalar,PetscBool,DSMatType,void*),void **ctx)
1032: {
1033: DS_NEP *dsctx = (DS_NEP*)ds->data;
1035: if (fun) *fun = dsctx->computematrix;
1036: if (ctx) *ctx = dsctx->computematrixctx;
1037: return 0;
1038: }
1040: /*@C
1041: DSNEPGetComputeMatrixFunction - Returns the user-provided callback function
1042: set in DSNEPSetComputeMatrixFunction().
1044: Not Collective
1046: Input Parameter:
1047: . ds - the direct solver context
1049: Output Parameters:
1050: + fun - the pointer to the user function
1051: - ctx - the context pointer
1053: Level: developer
1055: .seealso: DSNEPSetComputeMatrixFunction()
1056: @*/
1057: PetscErrorCode DSNEPGetComputeMatrixFunction(DS ds,PetscErrorCode (**fun)(DS,PetscScalar,PetscBool,DSMatType,void*),void **ctx)
1058: {
1060: PetscUseMethod(ds,"DSNEPGetComputeMatrixFunction_C",(DS,PetscErrorCode (**)(DS,PetscScalar,PetscBool,DSMatType,void*),void**),(ds,fun,ctx));
1061: return 0;
1062: }
1064: static PetscErrorCode DSNEPSetRG_NEP(DS ds,RG rg)
1065: {
1066: DS_NEP *dsctx = (DS_NEP*)ds->data;
1068: PetscObjectReference((PetscObject)rg);
1069: RGDestroy(&dsctx->rg);
1070: dsctx->rg = rg;
1071: return 0;
1072: }
1074: /*@
1075: DSNEPSetRG - Associates a region object to the DSNEP solver.
1077: Logically Collective on ds
1079: Input Parameters:
1080: + ds - the direct solver context
1081: - rg - the region context
1083: Notes:
1084: The region is used only in the contour integral method, and
1085: should enclose the wanted eigenvalues.
1087: Level: developer
1089: .seealso: DSNEPGetRG()
1090: @*/
1091: PetscErrorCode DSNEPSetRG(DS ds,RG rg)
1092: {
1094: if (rg) {
1097: }
1098: PetscTryMethod(ds,"DSNEPSetRG_C",(DS,RG),(ds,rg));
1099: return 0;
1100: }
1102: static PetscErrorCode DSNEPGetRG_NEP(DS ds,RG *rg)
1103: {
1104: DS_NEP *ctx = (DS_NEP*)ds->data;
1106: if (!ctx->rg) {
1107: RGCreate(PetscObjectComm((PetscObject)ds),&ctx->rg);
1108: PetscObjectIncrementTabLevel((PetscObject)ctx->rg,(PetscObject)ds,1);
1109: RGSetOptionsPrefix(ctx->rg,((PetscObject)ds)->prefix);
1110: RGAppendOptionsPrefix(ctx->rg,"ds_nep_");
1111: PetscObjectSetOptions((PetscObject)ctx->rg,((PetscObject)ds)->options);
1112: }
1113: *rg = ctx->rg;
1114: return 0;
1115: }
1117: /*@
1118: DSNEPGetRG - Obtain the region object associated to the DSNEP solver.
1120: Not Collective
1122: Input Parameter:
1123: . ds - the direct solver context
1125: Output Parameter:
1126: . rg - the region context
1128: Level: developer
1130: .seealso: DSNEPSetRG()
1131: @*/
1132: PetscErrorCode DSNEPGetRG(DS ds,RG *rg)
1133: {
1136: PetscUseMethod(ds,"DSNEPGetRG_C",(DS,RG*),(ds,rg));
1137: return 0;
1138: }
1140: PetscErrorCode DSSetFromOptions_NEP(DS ds,PetscOptionItems *PetscOptionsObject)
1141: {
1142: PetscInt k;
1143: PetscBool flg;
1144: #if defined(PETSC_USE_COMPLEX)
1145: PetscReal r;
1146: PetscBool flg1;
1147: DS_NEP *ctx = (DS_NEP*)ds->data;
1148: #endif
1150: PetscOptionsHeadBegin(PetscOptionsObject,"DS NEP Options");
1152: PetscOptionsInt("-ds_nep_minimality","Maximum minimality index","DSNEPSetMinimality",4,&k,&flg);
1153: if (flg) DSNEPSetMinimality(ds,k);
1155: PetscOptionsInt("-ds_nep_integration_points","Number of integration points","DSNEPSetIntegrationPoints",64,&k,&flg);
1156: if (flg) DSNEPSetIntegrationPoints(ds,k);
1158: PetscOptionsInt("-ds_nep_sampling_size","Number of sampling columns","DSNEPSetSamplingSize",0,&k,&flg);
1159: if (flg) DSNEPSetSamplingSize(ds,k);
1161: #if defined(PETSC_USE_COMPLEX)
1162: r = ctx->rtol;
1163: PetscOptionsReal("-ds_nep_refine_tol","Refinement tolerance","DSNEPSetRefine",ctx->rtol,&r,&flg1);
1164: k = ctx->Nit;
1165: PetscOptionsInt("-ds_nep_refine_its","Number of iterative refinement iterations","DSNEPSetRefine",ctx->Nit,&k,&flg);
1166: if (flg1||flg) DSNEPSetRefine(ds,r,k);
1168: if (ds->method==1) {
1169: if (!ctx->rg) DSNEPGetRG(ds,&ctx->rg);
1170: RGSetFromOptions(ctx->rg);
1171: }
1172: #endif
1174: PetscOptionsHeadEnd();
1175: return 0;
1176: }
1178: PetscErrorCode DSDestroy_NEP(DS ds)
1179: {
1180: DS_NEP *ctx = (DS_NEP*)ds->data;
1181: PetscInt i;
1183: for (i=0;i<ctx->nf;i++) FNDestroy(&ctx->f[i]);
1184: RGDestroy(&ctx->rg);
1185: PetscLayoutDestroy(&ctx->map);
1186: PetscFree(ds->data);
1187: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetFN_C",NULL);
1188: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetFN_C",NULL);
1189: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetNumFN_C",NULL);
1190: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetMinimality_C",NULL);
1191: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetMinimality_C",NULL);
1192: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetRefine_C",NULL);
1193: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetRefine_C",NULL);
1194: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetIntegrationPoints_C",NULL);
1195: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetIntegrationPoints_C",NULL);
1196: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetSamplingSize_C",NULL);
1197: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetSamplingSize_C",NULL);
1198: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetRG_C",NULL);
1199: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetRG_C",NULL);
1200: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetComputeMatrixFunction_C",NULL);
1201: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetComputeMatrixFunction_C",NULL);
1202: return 0;
1203: }
1205: PetscErrorCode DSMatGetSize_NEP(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
1206: {
1207: DS_NEP *ctx = (DS_NEP*)ds->data;
1209: *rows = ds->n;
1210: if (t==DS_MAT_Q || t==DS_MAT_Z || t==DS_MAT_U || t==DS_MAT_V) *rows *= ctx->max_mid;
1211: *cols = ds->n;
1212: 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;
1213: return 0;
1214: }
1216: /*MC
1217: DSNEP - Dense Nonlinear Eigenvalue Problem.
1219: Level: beginner
1221: Notes:
1222: The problem is expressed as T(lambda)*x = 0, where T(lambda) is a
1223: parameter-dependent matrix written as T(lambda) = sum_i E_i*f_i(lambda).
1224: The eigenvalues lambda are the arguments returned by DSSolve()..
1226: The coefficient matrices E_i are the extra matrices of the DS, and
1227: the scalar functions f_i are passed via DSNEPSetFN(). Optionally, a
1228: callback function to fill the E_i matrices can be set with
1229: DSNEPSetComputeMatrixFunction().
1231: Used DS matrices:
1232: + DS_MAT_Ex - coefficient matrices of the split form of T(lambda)
1233: . DS_MAT_A - (workspace) T(lambda) evaluated at a given lambda (SLP only)
1234: . DS_MAT_B - (workspace) T'(lambda) evaluated at a given lambda (SLP only)
1235: . DS_MAT_Q - (workspace) left Hankel matrix (contour only)
1236: . DS_MAT_Z - (workspace) right Hankel matrix (contour only)
1237: . DS_MAT_U - (workspace) left singular vectors (contour only)
1238: . DS_MAT_V - (workspace) right singular vectors (contour only)
1239: - DS_MAT_W - (workspace) auxiliary matrix of size nxn
1241: Implemented methods:
1242: + 0 - Successive Linear Problems (SLP), computes just one eigenpair
1243: - 1 - Contour integral, computes all eigenvalues inside a region
1245: .seealso: DSCreate(), DSSetType(), DSType, DSNEPSetFN(), DSNEPSetComputeMatrixFunction()
1246: M*/
1247: SLEPC_EXTERN PetscErrorCode DSCreate_NEP(DS ds)
1248: {
1249: DS_NEP *ctx;
1251: PetscNew(&ctx);
1252: ds->data = (void*)ctx;
1253: ctx->max_mid = 4;
1254: ctx->nnod = 64;
1255: ctx->Nit = 3;
1256: ctx->rtol = PETSC_MACHINE_EPSILON/PetscSqrtReal(PETSC_SQRT_MACHINE_EPSILON);
1258: ds->ops->allocate = DSAllocate_NEP;
1259: ds->ops->setfromoptions = DSSetFromOptions_NEP;
1260: ds->ops->view = DSView_NEP;
1261: ds->ops->vectors = DSVectors_NEP;
1262: ds->ops->solve[0] = DSSolve_NEP_SLP;
1263: #if defined(PETSC_USE_COMPLEX)
1264: ds->ops->solve[1] = DSSolve_NEP_Contour;
1265: #endif
1266: ds->ops->sort = DSSort_NEP;
1267: #if !defined(PETSC_HAVE_MPIUNI)
1268: ds->ops->synchronize = DSSynchronize_NEP;
1269: #endif
1270: ds->ops->destroy = DSDestroy_NEP;
1271: ds->ops->matgetsize = DSMatGetSize_NEP;
1273: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetFN_C",DSNEPSetFN_NEP);
1274: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetFN_C",DSNEPGetFN_NEP);
1275: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetNumFN_C",DSNEPGetNumFN_NEP);
1276: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetMinimality_C",DSNEPGetMinimality_NEP);
1277: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetMinimality_C",DSNEPSetMinimality_NEP);
1278: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetRefine_C",DSNEPGetRefine_NEP);
1279: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetRefine_C",DSNEPSetRefine_NEP);
1280: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetIntegrationPoints_C",DSNEPGetIntegrationPoints_NEP);
1281: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetIntegrationPoints_C",DSNEPSetIntegrationPoints_NEP);
1282: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetSamplingSize_C",DSNEPGetSamplingSize_NEP);
1283: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetSamplingSize_C",DSNEPSetSamplingSize_NEP);
1284: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetRG_C",DSNEPSetRG_NEP);
1285: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetRG_C",DSNEPGetRG_NEP);
1286: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetComputeMatrixFunction_C",DSNEPSetComputeMatrixFunction_NEP);
1287: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetComputeMatrixFunction_C",DSNEPGetComputeMatrixFunction_NEP);
1288: return 0;
1289: }