Actual source code: dsnep.c
slepc-3.17.0 2022-03-31
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,*E,alpha;
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: DSGetArray(ds,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: E = ds->mat[DSMatExtra[i]];
53: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&k,&alpha,E,&inc,T,&inc));
54: }
55: DSRestoreArray(ds,mat,&T);
56: }
57: PetscLogEventEnd(DS_Other,ds,0,0,0);
58: PetscFunctionReturn(0);
59: }
61: PetscErrorCode DSAllocate_NEP(DS ds,PetscInt ld)
62: {
63: DS_NEP *ctx = (DS_NEP*)ds->data;
64: PetscInt i;
66: DSAllocateMat_Private(ds,DS_MAT_X);
67: for (i=0;i<ctx->nf;i++) DSAllocateMat_Private(ds,DSMatExtra[i]);
68: PetscFree(ds->perm);
69: PetscMalloc1(ld*ctx->max_mid,&ds->perm);
70: PetscLogObjectMemory((PetscObject)ds,ld*ctx->max_mid*sizeof(PetscInt));
71: PetscFunctionReturn(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=sizeof(methodname)/sizeof(methodname[0]);
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: PetscFunctionReturn(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: PetscFunctionReturn(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: PetscFunctionReturn(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 *A;
128: if (!ds->sc) PetscFunctionReturn(0);
129: n = ds->n*ctx->max_mid;
130: lds = ds->ld*ctx->max_mid;
131: l = ds->l;
132: A = ds->mat[DS_MAT_A];
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: for (i=l;i<ds->t;i++) A[i+i*lds] = wr[perm[i]];
138: for (i=l;i<ds->t;i++) wr[i] = A[i+i*lds];
139: /* n != ds->n */
140: DSPermuteColumns_Private(ds,0,ds->t,ds->n,DS_MAT_X,perm);
141: PetscFunctionReturn(0);
142: }
144: PetscErrorCode DSSolve_NEP_SLP(DS ds,PetscScalar *wr,PetscScalar *wi)
145: {
146: PetscScalar *A,*B,*W,*X,*work,*alpha,*beta;
147: PetscScalar sigma,lambda,mu,re,re2,sone=1.0,szero=0.0;
148: PetscBLASInt info,n,ld,lrwork=0,lwork,one=1,zero=0;
149: PetscInt it,pos,j,maxit=100,result;
150: PetscReal norm,tol,done=1.0;
151: #if defined(PETSC_USE_COMPLEX)
152: PetscReal *rwork;
153: #else
154: PetscReal *alphai,im,im2;
155: #endif
157: if (!ds->mat[DS_MAT_A]) DSAllocateMat_Private(ds,DS_MAT_A);
158: if (!ds->mat[DS_MAT_B]) DSAllocateMat_Private(ds,DS_MAT_B);
159: if (!ds->mat[DS_MAT_W]) DSAllocateMat_Private(ds,DS_MAT_W);
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: A = ds->mat[DS_MAT_A];
180: B = ds->mat[DS_MAT_B];
181: W = ds->mat[DS_MAT_W];
182: X = ds->mat[DS_MAT_X];
184: sigma = 0.0;
185: if (ds->sc->comparison==SlepcCompareTargetMagnitude || ds->sc->comparison==SlepcCompareTargetReal) sigma = *(PetscScalar*)ds->sc->comparisonctx;
186: lambda = sigma;
187: tol = n*PETSC_MACHINE_EPSILON/PetscSqrtReal(PETSC_SQRT_MACHINE_EPSILON);
189: for (it=0;it<maxit;it++) {
191: /* evaluate T and T' */
192: DSNEPComputeMatrix(ds,lambda,PETSC_FALSE,DS_MAT_A);
193: if (it) {
194: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,A,&ld,X,&one,&szero,X+ld,&one));
195: norm = BLASnrm2_(&n,X+ld,&one);
196: if (norm/PetscAbsScalar(lambda)<=tol) break;
197: }
198: DSNEPComputeMatrix(ds,lambda,PETSC_TRUE,DS_MAT_B);
200: /* compute eigenvalue correction mu and eigenvector u */
201: #if defined(PETSC_USE_COMPLEX)
202: rwork = ds->rwork;
203: PetscStackCallBLAS("LAPACKggev",LAPACKggev_("N","V",&n,A,&ld,B,&ld,alpha,beta,NULL,&ld,W,&ld,work,&lwork,rwork,&info));
204: #else
205: PetscStackCallBLAS("LAPACKggev",LAPACKggev_("N","V",&n,A,&ld,B,&ld,alpha,alphai,beta,NULL,&ld,W,&ld,work,&lwork,&info));
206: #endif
207: SlepcCheckLapackInfo("ggev",info);
209: /* find smallest eigenvalue */
210: j = 0;
211: if (beta[j]==0.0) re = (PetscRealPart(alpha[j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
212: else re = alpha[j]/beta[j];
213: #if !defined(PETSC_USE_COMPLEX)
214: if (beta[j]==0.0) im = (alphai[j]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
215: else im = alphai[j]/beta[j];
216: #endif
217: pos = 0;
218: for (j=1;j<n;j++) {
219: if (beta[j]==0.0) re2 = (PetscRealPart(alpha[j])>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
220: else re2 = alpha[j]/beta[j];
221: #if !defined(PETSC_USE_COMPLEX)
222: if (beta[j]==0.0) im2 = (alphai[j]>0.0)? PETSC_MAX_REAL: PETSC_MIN_REAL;
223: else im2 = alphai[j]/beta[j];
224: SlepcCompareSmallestMagnitude(re,im,re2,im2,&result,NULL);
225: #else
226: SlepcCompareSmallestMagnitude(re,0.0,re2,0.0,&result,NULL);
227: #endif
228: if (result > 0) {
229: re = re2;
230: #if !defined(PETSC_USE_COMPLEX)
231: im = im2;
232: #endif
233: pos = j;
234: }
235: }
237: #if !defined(PETSC_USE_COMPLEX)
239: #endif
240: mu = alpha[pos]/beta[pos];
241: PetscArraycpy(X,W+pos*ld,n);
242: norm = BLASnrm2_(&n,X,&one);
243: PetscStackCallBLAS("LAPACKlascl",LAPACKlascl_("G",&zero,&zero,&norm,&done,&n,&one,X,&n,&info));
244: SlepcCheckLapackInfo("lascl",info);
246: /* correct eigenvalue approximation */
247: lambda = lambda - mu;
248: }
251: ds->t = 1;
252: wr[0] = lambda;
253: if (wi) wi[0] = 0.0;
254: PetscFunctionReturn(0);
255: }
257: #if defined(PETSC_USE_COMPLEX)
258: /*
259: Newton refinement for eigenpairs computed with contour integral.
260: k - number of eigenpairs to refine
261: wr - eigenvalues (eigenvectors are stored in DS_MAT_X)
262: */
263: static PetscErrorCode DSNEPNewtonRefine(DS ds,PetscInt k,PetscScalar *wr)
264: {
265: DS_NEP *ctx = (DS_NEP*)ds->data;
266: PetscScalar *X,*W,*U,*R,sone=1.0,szero=0.0;
267: PetscReal norm;
268: PetscInt i,j,ii,nwu=0,*p,jstart=0,jend=k;
269: const PetscInt *range;
270: PetscBLASInt n,*perm,info,ld,one=1,n1;
271: PetscMPIInt len,size,root;
272: PetscLayout map;
274: X = ds->mat[DS_MAT_X];
275: W = ds->mat[DS_MAT_W];
276: PetscBLASIntCast(ds->n,&n);
277: PetscBLASIntCast(ds->ld,&ld);
278: n1 = n+1;
279: p = ds->perm;
280: PetscArrayzero(p,k);
281: DSAllocateWork_Private(ds,(n+2)*(n+1),0,n+1);
282: U = ds->work+nwu; nwu += (n+1)*(n+1);
283: R = ds->work+nwu; /*nwu += n+1;*/
284: perm = ds->iwork;
285: if (ds->pmode==DS_PARALLEL_DISTRIBUTED) {
286: PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)ds),PETSC_DECIDE,k,1,&map);
287: PetscLayoutGetRange(map,&jstart,&jend);
288: }
289: for (ii=0;ii<ctx->Nit;ii++) {
290: for (j=jstart;j<jend;j++) {
291: if (p[j]<2) {
292: DSNEPComputeMatrix(ds,wr[j],PETSC_FALSE,DS_MAT_W);
293: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,W,&ld,X+ld*j,&one,&szero,R,&one));
294: norm = BLASnrm2_(&n,R,&one);
295: if (norm/PetscAbsScalar(wr[j]) > ctx->rtol) {
296: PetscInfo(NULL,"Refining eigenpair %" PetscInt_FMT ", residual=%g\n",j,(double)norm/PetscAbsScalar(wr[j]));
297: p[j] = 1;
298: R[n] = 0.0;
299: for (i=0;i<n;i++) {
300: PetscArraycpy(U+i*n1,W+i*ld,n);
301: U[n+i*n1] = PetscConj(X[j*ld+i]);
302: }
303: U[n+n*n1] = 0.0;
304: DSNEPComputeMatrix(ds,wr[j],PETSC_TRUE,DS_MAT_W);
305: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,W,&ld,X+ld*j,&one,&szero,U+n*(n+1),&one));
306: /* solve system */
307: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
308: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n1,&n1,U,&n1,perm,&info));
309: SlepcCheckLapackInfo("getrf",info);
310: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&n1,&one,U,&n1,perm,R,&n1,&info));
311: SlepcCheckLapackInfo("getrs",info);
312: PetscFPTrapPop();
313: wr[j] -= R[n];
314: for (i=0;i<n;i++) X[j*ld+i] -= R[i];
315: /* normalization */
316: norm = BLASnrm2_(&n,X+ld*j,&one);
317: for (i=0;i<n;i++) X[ld*j+i] /= norm;
318: } else p[j] = 2;
319: }
320: }
321: }
322: if (ds->pmode==DS_PARALLEL_DISTRIBUTED) { /* communicate results */
323: PetscMPIIntCast(k,&len);
324: MPIU_Allreduce(MPI_IN_PLACE,p,len,MPIU_INT,MPIU_SUM,PetscObjectComm((PetscObject)ds));
325: MPI_Comm_size(PetscObjectComm((PetscObject)ds),&size);
326: PetscLayoutGetRanges(map,&range);
327: for (j=0;j<k;j++) {
328: if (p[j]) { /* j-th eigenpair has been refined */
329: for (root=0;root<size;root++) if (range[root+1]>j) break;
330: PetscMPIIntCast(1,&len);
331: MPI_Bcast(wr+j,len,MPIU_SCALAR,root,PetscObjectComm((PetscObject)ds));
332: PetscMPIIntCast(n,&len);
333: MPI_Bcast(X+ld*j,len,MPIU_SCALAR,root,PetscObjectComm((PetscObject)ds));
334: }
335: }
336: PetscLayoutDestroy(&map);
337: }
338: PetscFunctionReturn(0);
339: }
341: PetscErrorCode DSSolve_NEP_Contour(DS ds,PetscScalar *wr,PetscScalar *wi)
342: {
343: DS_NEP *ctx = (DS_NEP*)ds->data;
344: PetscScalar *alpha,*beta,*A,*B,*X,*W,*work,*Rc,*R,*w,*z,*zn,*S,*U,*V;
345: PetscScalar sone=1.0,szero=0.0,center;
346: PetscReal *rwork,norm,radius,vscale,rgscale,*sigma;
347: PetscBLASInt info,n,*perm,p,pp,ld,lwork,k_,rk_,colA,rowA,one=1;
348: PetscInt mid,lds,nnod=ctx->nnod,k,i,ii,jj,j,s,off,rk,nwu=0,nw,lrwork,*inside,kstart=0,kend=nnod;
349: PetscMPIInt len;
350: PetscBool isellipse;
351: PetscRandom rand;
354: /* Contour parameters */
355: PetscObjectTypeCompare((PetscObject)ctx->rg,RGELLIPSE,&isellipse);
357: RGEllipseGetParameters(ctx->rg,¢er,&radius,&vscale);
358: RGGetScale(ctx->rg,&rgscale);
359: if (ds->pmode==DS_PARALLEL_DISTRIBUTED) {
360: if (!ctx->map) PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)ds),PETSC_DECIDE,ctx->nnod,1,&ctx->map);
361: PetscLayoutGetRange(ctx->map,&kstart,&kend);
362: }
364: if (!ds->mat[DS_MAT_A]) DSAllocateMat_Private(ds,DS_MAT_A); /* size mid*n */
365: if (!ds->mat[DS_MAT_B]) DSAllocateMat_Private(ds,DS_MAT_B); /* size mid*n */
366: if (!ds->mat[DS_MAT_W]) DSAllocateMat_Private(ds,DS_MAT_W); /* size mid*n */
367: if (!ds->mat[DS_MAT_U]) DSAllocateMat_Private(ds,DS_MAT_U); /* size mid*n */
368: if (!ds->mat[DS_MAT_V]) DSAllocateMat_Private(ds,DS_MAT_V); /* size n */
369: A = ds->mat[DS_MAT_A];
370: B = ds->mat[DS_MAT_B];
371: W = ds->mat[DS_MAT_W];
372: U = ds->mat[DS_MAT_U];
373: V = ds->mat[DS_MAT_V];
374: X = ds->mat[DS_MAT_X];
375: mid = ctx->max_mid;
376: PetscBLASIntCast(ds->n,&n);
377: p = n; /* maximum number of columns for the probing matrix */
378: PetscBLASIntCast(ds->ld,&ld);
379: PetscBLASIntCast(mid*n,&rowA);
380: PetscBLASIntCast(5*rowA,&lwork);
381: nw = n*(2*p+7*mid)+3*nnod+2*mid*n*p;
382: lrwork = mid*n*6+8*n;
383: DSAllocateWork_Private(ds,nw,lrwork,n+1);
385: sigma = ds->rwork;
386: rwork = ds->rwork+mid*n;
387: perm = ds->iwork;
388: z = ds->work+nwu; nwu += nnod; /* quadrature points */
389: zn = ds->work+nwu; nwu += nnod; /* normalized quadrature points */
390: w = ds->work+nwu; nwu += nnod; /* quadrature weights */
391: Rc = ds->work+nwu; nwu += n*p;
392: R = ds->work+nwu; nwu += n*p;
393: alpha = ds->work+nwu; nwu += mid*n;
394: beta = ds->work+nwu; nwu += mid*n;
395: S = ds->work+nwu; nwu += 2*mid*n*p;
396: work = ds->work+nwu; /*nwu += mid*n*5;*/
398: /* Compute quadrature parameters */
399: RGComputeQuadrature(ctx->rg,RG_QUADRULE_TRAPEZOIDAL,nnod,z,zn,w);
401: /* Set random matrix */
402: PetscRandomCreate(PetscObjectComm((PetscObject)ds),&rand);
403: PetscRandomSetSeed(rand,0x12345678);
404: PetscRandomSeed(rand);
405: for (j=0;j<p;j++)
406: for (i=0;i<n;i++) PetscRandomGetValue(rand,Rc+i+j*n);
407: PetscArrayzero(S,2*mid*n*p);
408: /* Loop of integration points */
409: for (k=kstart;k<kend;k++) {
410: PetscInfo(NULL,"Solving integration point %" PetscInt_FMT "\n",k);
411: PetscArraycpy(R,Rc,p*n);
412: DSNEPComputeMatrix(ds,z[k],PETSC_FALSE,DS_MAT_V);
414: /* LU factorization */
415: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
416: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n,&n,V,&ld,perm,&info));
417: SlepcCheckLapackInfo("getrf",info);
418: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&n,&p,V,&ld,perm,R,&n,&info));
419: SlepcCheckLapackInfo("getrs",info);
420: PetscFPTrapPop();
422: /* Moments computation */
423: for (s=0;s<2*ctx->max_mid;s++) {
424: off = s*n*p;
425: for (j=0;j<p;j++)
426: for (i=0;i<n;i++) S[off+i+j*n] += w[k]*R[j*n+i];
427: w[k] *= zn[k];
428: }
429: }
431: if (ds->pmode==DS_PARALLEL_DISTRIBUTED) { /* compute final S via reduction */
432: PetscMPIIntCast(2*mid*n*p,&len);
433: MPIU_Allreduce(MPI_IN_PLACE,S,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)ds));
434: }
435: p = ctx->spls?PetscMin(ctx->spls,n):n;
436: pp = p;
437: do {
438: p = pp;
439: PetscBLASIntCast(mid*p,&colA);
441: PetscInfo(ds,"Computing SVD of size %" PetscBLASInt_FMT "x%" PetscBLASInt_FMT "\n",rowA,colA);
442: for (jj=0;jj<mid;jj++) {
443: for (ii=0;ii<mid;ii++) {
444: off = jj*p*rowA+ii*n;
445: for (j=0;j<p;j++)
446: for (i=0;i<n;i++) A[off+j*rowA+i] = S[((jj+ii)*n+j)*n+i];
447: }
448: }
449: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
450: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&rowA,&colA,A,&rowA,sigma,U,&rowA,W,&colA,work,&lwork,rwork,&info));
451: SlepcCheckLapackInfo("gesvd",info);
452: PetscFPTrapPop();
454: rk = colA;
455: for (i=1;i<colA;i++) if (sigma[i]/sigma[0]<PETSC_MACHINE_EPSILON*1e4) {rk = i; break;}
456: if (rk<colA || p==n) break;
457: pp *= 2;
458: } while (pp<=n);
459: PetscInfo(ds,"Solving generalized eigenproblem of size %" PetscInt_FMT "\n",rk);
460: for (jj=0;jj<mid;jj++) {
461: for (ii=0;ii<mid;ii++) {
462: off = jj*p*rowA+ii*n;
463: for (j=0;j<p;j++)
464: for (i=0;i<n;i++) A[off+j*rowA+i] = S[((jj+ii+1)*n+j)*n+i];
465: }
466: }
467: PetscBLASIntCast(rk,&rk_);
468: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","C",&rowA,&rk_,&colA,&sone,A,&rowA,W,&colA,&szero,B,&rowA));
469: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&rk_,&rk_,&rowA,&sone,U,&rowA,B,&rowA,&szero,A,&rk_));
470: PetscArrayzero(B,n*mid*n*mid);
471: for (j=0;j<rk;j++) B[j+j*rk_] = sigma[j];
472: PetscStackCallBLAS("LAPACKggev",LAPACKggev_("N","V",&rk_,A,&rk_,B,&rk_,alpha,beta,NULL,&ld,W,&rk_,work,&lwork,rwork,&info));
473: for (i=0;i<rk;i++) wr[i] = (center+radius*alpha[i]/beta[i])*rgscale;
474: PetscMalloc1(rk,&inside);
475: RGCheckInside(ctx->rg,rk,wr,wi,inside);
476: k=0;
477: for (i=0;i<rk;i++)
478: if (inside[i]==1) inside[k++] = i;
479: /* Discard values outside region */
480: lds = ld*mid;
481: PetscArrayzero(A,lds*lds);
482: PetscArrayzero(B,lds*lds);
483: for (i=0;i<k;i++) A[i+i*lds] = (center*beta[inside[i]]+radius*alpha[inside[i]])*rgscale;
484: for (i=0;i<k;i++) B[i+i*lds] = beta[inside[i]];
485: for (i=0;i<k;i++) wr[i] = A[i+i*lds]/B[i+i*lds];
486: for (j=0;j<k;j++) for (i=0;i<rk;i++) W[j*rk+i] = sigma[i]*W[inside[j]*rk+i];
487: PetscBLASIntCast(k,&k_);
488: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n,&k_,&rk_,&sone,U,&rowA,W,&rk_,&szero,X,&ld));
490: /* Normalize */
491: for (j=0;j<k;j++) {
492: norm = BLASnrm2_(&n,X+ld*j,&one);
493: for (i=0;i<n;i++) X[ld*j+i] /= norm;
494: }
495: PetscFree(inside);
496: /* Newton refinement */
497: DSNEPNewtonRefine(ds,k,wr);
498: ds->t = k;
499: PetscRandomDestroy(&rand);
500: PetscFunctionReturn(0);
501: }
502: #endif
504: PetscErrorCode DSSynchronize_NEP(DS ds,PetscScalar eigr[],PetscScalar eigi[])
505: {
506: PetscInt k=0;
507: PetscMPIInt n,rank,size,off=0;
509: if (ds->state>=DS_STATE_CONDENSED) k += ds->n;
510: if (eigr) k += 1;
511: if (eigi) k += 1;
512: DSAllocateWork_Private(ds,k,0,0);
513: PetscMPIIntCast(k*sizeof(PetscScalar),&size);
514: PetscMPIIntCast(ds->n,&n);
515: MPI_Comm_rank(PetscObjectComm((PetscObject)ds),&rank);
516: if (!rank) {
517: if (ds->state>=DS_STATE_CONDENSED) MPI_Pack(ds->mat[DS_MAT_X],n,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
518: if (eigr) MPI_Pack(eigr,1,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
519: #if !defined(PETSC_USE_COMPLEX)
520: if (eigi) MPI_Pack(eigi,1,MPIU_SCALAR,ds->work,size,&off,PetscObjectComm((PetscObject)ds));
521: #endif
522: }
523: MPI_Bcast(ds->work,size,MPI_BYTE,0,PetscObjectComm((PetscObject)ds));
524: if (rank) {
525: if (ds->state>=DS_STATE_CONDENSED) MPI_Unpack(ds->work,size,&off,ds->mat[DS_MAT_X],n,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
526: if (eigr) MPI_Unpack(ds->work,size,&off,eigr,1,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
527: #if !defined(PETSC_USE_COMPLEX)
528: if (eigi) MPI_Unpack(ds->work,size,&off,eigi,1,MPIU_SCALAR,PetscObjectComm((PetscObject)ds));
529: #endif
530: }
531: PetscFunctionReturn(0);
532: }
534: static PetscErrorCode DSNEPSetFN_NEP(DS ds,PetscInt n,FN fn[])
535: {
536: DS_NEP *ctx = (DS_NEP*)ds->data;
537: PetscInt i;
541: if (ds->ld) PetscInfo(ds,"DSNEPSetFN() called after DSAllocate()\n");
542: for (i=0;i<n;i++) PetscObjectReference((PetscObject)fn[i]);
543: for (i=0;i<ctx->nf;i++) FNDestroy(&ctx->f[i]);
544: for (i=0;i<n;i++) ctx->f[i] = fn[i];
545: ctx->nf = n;
546: PetscFunctionReturn(0);
547: }
549: /*@
550: DSNEPSetFN - Sets a number of functions that define the nonlinear
551: eigenproblem.
553: Collective on ds
555: Input Parameters:
556: + ds - the direct solver context
557: . n - number of functions
558: - fn - array of functions
560: Notes:
561: The nonlinear eigenproblem is defined in terms of the split nonlinear
562: operator T(lambda) = sum_i A_i*f_i(lambda).
564: This function must be called before DSAllocate(). Then DSAllocate()
565: will allocate an extra matrix A_i per each function, that can be
566: filled in the usual way.
568: Level: advanced
570: .seealso: DSNEPGetFN(), DSAllocate()
571: @*/
572: PetscErrorCode DSNEPSetFN(DS ds,PetscInt n,FN fn[])
573: {
574: PetscInt i;
579: for (i=0;i<n;i++) {
582: }
583: PetscTryMethod(ds,"DSNEPSetFN_C",(DS,PetscInt,FN[]),(ds,n,fn));
584: PetscFunctionReturn(0);
585: }
587: static PetscErrorCode DSNEPGetFN_NEP(DS ds,PetscInt k,FN *fn)
588: {
589: DS_NEP *ctx = (DS_NEP*)ds->data;
592: *fn = ctx->f[k];
593: PetscFunctionReturn(0);
594: }
596: /*@
597: DSNEPGetFN - Gets the functions associated with the nonlinear DS.
599: Not collective, though parallel FNs are returned if the DS is parallel
601: Input Parameters:
602: + ds - the direct solver context
603: - k - the index of the requested function (starting in 0)
605: Output Parameter:
606: . fn - the function
608: Level: advanced
610: .seealso: DSNEPSetFN()
611: @*/
612: PetscErrorCode DSNEPGetFN(DS ds,PetscInt k,FN *fn)
613: {
616: PetscUseMethod(ds,"DSNEPGetFN_C",(DS,PetscInt,FN*),(ds,k,fn));
617: PetscFunctionReturn(0);
618: }
620: static PetscErrorCode DSNEPGetNumFN_NEP(DS ds,PetscInt *n)
621: {
622: DS_NEP *ctx = (DS_NEP*)ds->data;
624: *n = ctx->nf;
625: PetscFunctionReturn(0);
626: }
628: /*@
629: DSNEPGetNumFN - Returns the number of functions stored internally by
630: the DS.
632: Not collective
634: Input Parameter:
635: . ds - the direct solver context
637: Output Parameters:
638: . n - the number of functions passed in DSNEPSetFN()
640: Level: advanced
642: .seealso: DSNEPSetFN()
643: @*/
644: PetscErrorCode DSNEPGetNumFN(DS ds,PetscInt *n)
645: {
648: PetscUseMethod(ds,"DSNEPGetNumFN_C",(DS,PetscInt*),(ds,n));
649: PetscFunctionReturn(0);
650: }
652: static PetscErrorCode DSNEPSetMinimality_NEP(DS ds,PetscInt n)
653: {
654: DS_NEP *ctx = (DS_NEP*)ds->data;
656: if (n == PETSC_DECIDE || n == PETSC_DEFAULT) ctx->max_mid = 4;
657: else {
659: ctx->max_mid = n;
660: }
661: PetscFunctionReturn(0);
662: }
664: /*@
665: DSNEPSetMinimality - Sets the maximum minimality index used internally by
666: the DSNEP.
668: Logically Collective on ds
670: Input Parameters:
671: + ds - the direct solver context
672: - n - the maximum minimality index
674: Options Database Key:
675: . -ds_nep_minimality <n> - sets the maximum minimality index
677: Notes:
678: The maximum minimality index is used only in the contour integral method,
679: and is related to the highest momemts used in the method. The default
680: value is 1, an larger value might give better accuracy in some cases, but
681: at a higher cost.
683: Level: advanced
685: .seealso: DSNEPGetMinimality()
686: @*/
687: PetscErrorCode DSNEPSetMinimality(DS ds,PetscInt n)
688: {
691: PetscTryMethod(ds,"DSNEPSetMinimality_C",(DS,PetscInt),(ds,n));
692: PetscFunctionReturn(0);
693: }
695: static PetscErrorCode DSNEPGetMinimality_NEP(DS ds,PetscInt *n)
696: {
697: DS_NEP *ctx = (DS_NEP*)ds->data;
699: *n = ctx->max_mid;
700: PetscFunctionReturn(0);
701: }
703: /*@
704: DSNEPGetMinimality - Returns the maximum minimality index used internally by
705: the DSNEP.
707: Not collective
709: Input Parameter:
710: . ds - the direct solver context
712: Output Parameters:
713: . n - the maximum minimality index passed in DSNEPSetMinimality()
715: Level: advanced
717: .seealso: DSNEPSetMinimality()
718: @*/
719: PetscErrorCode DSNEPGetMinimality(DS ds,PetscInt *n)
720: {
723: PetscUseMethod(ds,"DSNEPGetMinimality_C",(DS,PetscInt*),(ds,n));
724: PetscFunctionReturn(0);
725: }
727: static PetscErrorCode DSNEPSetRefine_NEP(DS ds,PetscReal tol,PetscInt its)
728: {
729: DS_NEP *ctx = (DS_NEP*)ds->data;
731: if (tol == PETSC_DEFAULT) ctx->rtol = PETSC_MACHINE_EPSILON/PetscSqrtReal(PETSC_SQRT_MACHINE_EPSILON);
732: else {
734: ctx->rtol = tol;
735: }
736: if (its == PETSC_DECIDE || its == PETSC_DEFAULT) ctx->Nit = 3;
737: else {
739: ctx->Nit = its;
740: }
741: PetscFunctionReturn(0);
742: }
744: /*@
745: DSNEPSetRefine - Sets the tolerance and the number of iterations of Newton iterative
746: refinement for eigenpairs.
748: Logically Collective on ds
750: Input Parameters:
751: + ds - the direct solver context
752: . tol - the tolerance
753: - its - the number of iterations
755: Options Database Key:
756: + -ds_nep_refine_tol <tol> - sets the tolerance
757: - -ds_nep_refine_its <its> - sets the number of Newton iterations
759: Notes:
760: Iterative refinement of eigenpairs is currently used only in the contour
761: integral method.
763: Level: advanced
765: .seealso: DSNEPGetRefine()
766: @*/
767: PetscErrorCode DSNEPSetRefine(DS ds,PetscReal tol,PetscInt its)
768: {
772: PetscTryMethod(ds,"DSNEPSetRefine_C",(DS,PetscReal,PetscInt),(ds,tol,its));
773: PetscFunctionReturn(0);
774: }
776: static PetscErrorCode DSNEPGetRefine_NEP(DS ds,PetscReal *tol,PetscInt *its)
777: {
778: DS_NEP *ctx = (DS_NEP*)ds->data;
780: if (tol) *tol = ctx->rtol;
781: if (its) *its = ctx->Nit;
782: PetscFunctionReturn(0);
783: }
785: /*@
786: DSNEPGetRefine - Returns the tolerance and the number of iterations of Newton iterative
787: refinement for eigenpairs.
789: Not collective
791: Input Parameter:
792: . ds - the direct solver context
794: Output Parameters:
795: + tol - the tolerance
796: - its - the number of iterations
798: Level: advanced
800: .seealso: DSNEPSetRefine()
801: @*/
802: PetscErrorCode DSNEPGetRefine(DS ds,PetscReal *tol,PetscInt *its)
803: {
805: PetscUseMethod(ds,"DSNEPGetRefine_C",(DS,PetscReal*,PetscInt*),(ds,tol,its));
806: PetscFunctionReturn(0);
807: }
809: static PetscErrorCode DSNEPSetIntegrationPoints_NEP(DS ds,PetscInt ip)
810: {
811: DS_NEP *ctx = (DS_NEP*)ds->data;
813: if (ip == PETSC_DECIDE || ip == PETSC_DEFAULT) ctx->nnod = 64;
814: else {
816: ctx->nnod = ip;
817: }
818: PetscLayoutDestroy(&ctx->map); /* need to redistribute at next solve */
819: PetscFunctionReturn(0);
820: }
822: /*@
823: DSNEPSetIntegrationPoints - Sets the number of integration points to be
824: used in the contour integral method.
826: Logically Collective on ds
828: Input Parameters:
829: + ds - the direct solver context
830: - ip - the number of integration points
832: Options Database Key:
833: . -ds_nep_integration_points <ip> - sets the number of integration points
835: Notes:
836: This parameter is relevant only in the contour integral method.
838: Level: advanced
840: .seealso: DSNEPGetIntegrationPoints()
841: @*/
842: PetscErrorCode DSNEPSetIntegrationPoints(DS ds,PetscInt ip)
843: {
846: PetscTryMethod(ds,"DSNEPSetIntegrationPoints_C",(DS,PetscInt),(ds,ip));
847: PetscFunctionReturn(0);
848: }
850: static PetscErrorCode DSNEPGetIntegrationPoints_NEP(DS ds,PetscInt *ip)
851: {
852: DS_NEP *ctx = (DS_NEP*)ds->data;
854: *ip = ctx->nnod;
855: PetscFunctionReturn(0);
856: }
858: /*@
859: DSNEPGetIntegrationPoints - Returns the number of integration points used
860: in the contour integral method.
862: Not collective
864: Input Parameter:
865: . ds - the direct solver context
867: Output Parameters:
868: . ip - the number of integration points
870: Level: advanced
872: .seealso: DSNEPSetIntegrationPoints()
873: @*/
874: PetscErrorCode DSNEPGetIntegrationPoints(DS ds,PetscInt *ip)
875: {
878: PetscUseMethod(ds,"DSNEPGetIntegrationPoints_C",(DS,PetscInt*),(ds,ip));
879: PetscFunctionReturn(0);
880: }
882: static PetscErrorCode DSNEPSetSamplingSize_NEP(DS ds,PetscInt p)
883: {
884: DS_NEP *ctx = (DS_NEP*)ds->data;
886: if (p == PETSC_DECIDE || p == PETSC_DEFAULT) ctx->spls = 0;
887: else {
889: ctx->spls = p;
890: }
891: PetscFunctionReturn(0);
892: }
894: /*@
895: DSNEPSetSamplingSize - Sets the number of sampling columns to be
896: used in the contour integral method.
898: Logically Collective on ds
900: Input Parameters:
901: + ds - the direct solver context
902: - p - the number of columns for the sampling matrix
904: Options Database Key:
905: . -ds_nep_sampling_size <p> - sets the number of sampling columns
907: Notes:
908: This parameter is relevant only in the contour integral method.
910: Level: advanced
912: .seealso: DSNEPGetSamplingSize()
913: @*/
914: PetscErrorCode DSNEPSetSamplingSize(DS ds,PetscInt p)
915: {
918: PetscTryMethod(ds,"DSNEPSetSamplingSize_C",(DS,PetscInt),(ds,p));
919: PetscFunctionReturn(0);
920: }
922: static PetscErrorCode DSNEPGetSamplingSize_NEP(DS ds,PetscInt *p)
923: {
924: DS_NEP *ctx = (DS_NEP*)ds->data;
926: *p = ctx->spls;
927: PetscFunctionReturn(0);
928: }
930: /*@
931: DSNEPGetSamplingSize - Returns the number of sampling columns used
932: in the contour integral method.
934: Not collective
936: Input Parameter:
937: . ds - the direct solver context
939: Output Parameters:
940: . p - the number of columns for the sampling matrix
942: Level: advanced
944: .seealso: DSNEPSetSamplingSize()
945: @*/
946: PetscErrorCode DSNEPGetSamplingSize(DS ds,PetscInt *p)
947: {
950: PetscUseMethod(ds,"DSNEPGetSamplingSize_C",(DS,PetscInt*),(ds,p));
951: PetscFunctionReturn(0);
952: }
954: static PetscErrorCode DSNEPSetComputeMatrixFunction_NEP(DS ds,PetscErrorCode (*fun)(DS,PetscScalar,PetscBool,DSMatType,void*),void *ctx)
955: {
956: DS_NEP *dsctx = (DS_NEP*)ds->data;
958: dsctx->computematrix = fun;
959: dsctx->computematrixctx = ctx;
960: PetscFunctionReturn(0);
961: }
963: /*@C
964: DSNEPSetComputeMatrixFunction - Sets a user-provided subroutine to compute
965: the matrices T(lambda) or T'(lambda).
967: Logically Collective on ds
969: Input Parameters:
970: + ds - the direct solver context
971: . fun - a pointer to the user function
972: - ctx - a context pointer (the last parameter to the user function)
974: Calling Sequence of fun:
975: $ fun(DS ds,PetscScalar lambda,PetscBool deriv,DSMatType mat,void *ctx)
977: + ds - the direct solver object
978: . lambda - point where T(lambda) or T'(lambda) must be evaluated
979: . deriv - if true compute T'(lambda), otherwise compute T(lambda)
980: . mat - the DS matrix where the result must be stored
981: - ctx - optional context, as set by DSNEPSetComputeMatrixFunction()
983: Note:
984: The result is computed as T(lambda) = sum_i E_i*f_i(lambda), and similarly
985: for the derivative.
987: Level: developer
989: .seealso: DSNEPGetComputeMatrixFunction()
990: @*/
991: PetscErrorCode DSNEPSetComputeMatrixFunction(DS ds,PetscErrorCode (*fun)(DS,PetscScalar,PetscBool,DSMatType,void*),void *ctx)
992: {
994: PetscTryMethod(ds,"DSNEPSetComputeMatrixFunction_C",(DS,PetscErrorCode (*)(DS,PetscScalar,PetscBool,DSMatType,void*),void*),(ds,fun,ctx));
995: PetscFunctionReturn(0);
996: }
998: static PetscErrorCode DSNEPGetComputeMatrixFunction_NEP(DS ds,PetscErrorCode (**fun)(DS,PetscScalar,PetscBool,DSMatType,void*),void **ctx)
999: {
1000: DS_NEP *dsctx = (DS_NEP*)ds->data;
1002: if (fun) *fun = dsctx->computematrix;
1003: if (ctx) *ctx = dsctx->computematrixctx;
1004: PetscFunctionReturn(0);
1005: }
1007: /*@C
1008: DSNEPGetComputeMatrixFunction - Returns the user-provided callback function
1009: set in DSNEPSetComputeMatrixFunction().
1011: Not Collective
1013: Input Parameter:
1014: . ds - the direct solver context
1016: Output Parameters:
1017: + fun - the pointer to the user function
1018: - ctx - the context pointer
1020: Level: developer
1022: .seealso: DSNEPSetComputeMatrixFunction()
1023: @*/
1024: PetscErrorCode DSNEPGetComputeMatrixFunction(DS ds,PetscErrorCode (**fun)(DS,PetscScalar,PetscBool,DSMatType,void*),void **ctx)
1025: {
1027: PetscUseMethod(ds,"DSNEPGetComputeMatrixFunction_C",(DS,PetscErrorCode (**)(DS,PetscScalar,PetscBool,DSMatType,void*),void**),(ds,fun,ctx));
1028: PetscFunctionReturn(0);
1029: }
1031: static PetscErrorCode DSNEPSetRG_NEP(DS ds,RG rg)
1032: {
1033: DS_NEP *dsctx = (DS_NEP*)ds->data;
1035: PetscObjectReference((PetscObject)rg);
1036: RGDestroy(&dsctx->rg);
1037: dsctx->rg = rg;
1038: PetscLogObjectParent((PetscObject)ds,(PetscObject)dsctx->rg);
1039: PetscFunctionReturn(0);
1040: }
1042: /*@
1043: DSNEPSetRG - Associates a region object to the DSNEP solver.
1045: Logically Collective on ds
1047: Input Parameters:
1048: + ds - the direct solver context
1049: - rg - the region context
1051: Notes:
1052: The region is used only in the contour integral method, and
1053: should enclose the wanted eigenvalues.
1055: Level: developer
1057: .seealso: DSNEPGetRG()
1058: @*/
1059: PetscErrorCode DSNEPSetRG(DS ds,RG rg)
1060: {
1062: if (rg) {
1065: }
1066: PetscTryMethod(ds,"DSNEPSetRG_C",(DS,RG),(ds,rg));
1067: PetscFunctionReturn(0);
1068: }
1070: static PetscErrorCode DSNEPGetRG_NEP(DS ds,RG *rg)
1071: {
1072: DS_NEP *ctx = (DS_NEP*)ds->data;
1074: if (!ctx->rg) {
1075: RGCreate(PetscObjectComm((PetscObject)ds),&ctx->rg);
1076: PetscObjectIncrementTabLevel((PetscObject)ctx->rg,(PetscObject)ds,1);
1077: RGSetOptionsPrefix(ctx->rg,((PetscObject)ds)->prefix);
1078: RGAppendOptionsPrefix(ctx->rg,"ds_nep_");
1079: PetscLogObjectParent((PetscObject)ds,(PetscObject)ctx->rg);
1080: PetscObjectSetOptions((PetscObject)ctx->rg,((PetscObject)ds)->options);
1081: }
1082: *rg = ctx->rg;
1083: PetscFunctionReturn(0);
1084: }
1086: /*@
1087: DSNEPGetRG - Obtain the region object associated to the DSNEP solver.
1089: Not Collective
1091: Input Parameter:
1092: . ds - the direct solver context
1094: Output Parameter:
1095: . rg - the region context
1097: Level: developer
1099: .seealso: DSNEPSetRG()
1100: @*/
1101: PetscErrorCode DSNEPGetRG(DS ds,RG *rg)
1102: {
1105: PetscUseMethod(ds,"DSNEPGetRG_C",(DS,RG*),(ds,rg));
1106: PetscFunctionReturn(0);
1107: }
1109: PetscErrorCode DSSetFromOptions_NEP(PetscOptionItems *PetscOptionsObject,DS ds)
1110: {
1111: PetscInt k;
1112: PetscBool flg;
1113: #if defined(PETSC_USE_COMPLEX)
1114: PetscReal r;
1115: PetscBool flg1;
1116: DS_NEP *ctx = (DS_NEP*)ds->data;
1117: #endif
1119: PetscOptionsHead(PetscOptionsObject,"DS NEP Options");
1121: PetscOptionsInt("-ds_nep_minimality","Maximum minimality index","DSNEPSetMinimality",4,&k,&flg);
1122: if (flg) DSNEPSetMinimality(ds,k);
1124: PetscOptionsInt("-ds_nep_integration_points","Number of integration points","DSNEPSetIntegrationPoints",64,&k,&flg);
1125: if (flg) DSNEPSetIntegrationPoints(ds,k);
1127: PetscOptionsInt("-ds_nep_sampling_size","Number of sampling columns","DSNEPSetSamplingSize",0,&k,&flg);
1128: if (flg) DSNEPSetSamplingSize(ds,k);
1130: #if defined(PETSC_USE_COMPLEX)
1131: r = ctx->rtol;
1132: PetscOptionsReal("-ds_nep_refine_tol","Refinement tolerance","DSNEPSetRefine",ctx->rtol,&r,&flg1);
1133: k = ctx->Nit;
1134: PetscOptionsInt("-ds_nep_refine_its","Number of iterative refinement iterations","DSNEPSetRefine",ctx->Nit,&k,&flg);
1135: if (flg1||flg) DSNEPSetRefine(ds,r,k);
1137: if (ds->method==1) {
1138: if (!ctx->rg) DSNEPGetRG(ds,&ctx->rg);
1139: RGSetFromOptions(ctx->rg);
1140: }
1141: #endif
1143: PetscOptionsTail();
1144: PetscFunctionReturn(0);
1145: }
1147: PetscErrorCode DSDestroy_NEP(DS ds)
1148: {
1149: DS_NEP *ctx = (DS_NEP*)ds->data;
1150: PetscInt i;
1152: for (i=0;i<ctx->nf;i++) FNDestroy(&ctx->f[i]);
1153: RGDestroy(&ctx->rg);
1154: PetscLayoutDestroy(&ctx->map);
1155: PetscFree(ds->data);
1156: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetFN_C",NULL);
1157: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetFN_C",NULL);
1158: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetNumFN_C",NULL);
1159: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetMinimality_C",NULL);
1160: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetMinimality_C",NULL);
1161: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetRefine_C",NULL);
1162: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetRefine_C",NULL);
1163: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetIntegrationPoints_C",NULL);
1164: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetIntegrationPoints_C",NULL);
1165: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetSamplingSize_C",NULL);
1166: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetSamplingSize_C",NULL);
1167: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetRG_C",NULL);
1168: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetRG_C",NULL);
1169: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetComputeMatrixFunction_C",NULL);
1170: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetComputeMatrixFunction_C",NULL);
1171: PetscFunctionReturn(0);
1172: }
1174: PetscErrorCode DSMatGetSize_NEP(DS ds,DSMatType t,PetscInt *rows,PetscInt *cols)
1175: {
1176: DS_NEP *ctx = (DS_NEP*)ds->data;
1178: *rows = ds->n;
1179: *cols = ds->n;
1180: if (t==DS_MAT_X || t==DS_MAT_Y) *cols *= ctx->max_mid;
1181: PetscFunctionReturn(0);
1182: }
1184: /*MC
1185: DSNEP - Dense Nonlinear Eigenvalue Problem.
1187: Level: beginner
1189: Notes:
1190: The problem is expressed as T(lambda)*x = 0, where T(lambda) is a
1191: parameter-dependent matrix written as T(lambda) = sum_i E_i*f_i(lambda).
1192: The eigenvalues lambda are the arguments returned by DSSolve()..
1194: The coefficient matrices E_i are the extra matrices of the DS, and
1195: the scalar functions f_i are passed via DSNEPSetFN(). Optionally, a
1196: callback function to fill the E_i matrices can be set with
1197: DSNEPSetComputeMatrixFunction().
1199: Used DS matrices:
1200: + DS_MAT_Ex - coefficient matrices of the split form of T(lambda)
1201: . DS_MAT_A - (workspace) T(lambda) evaluated at a given lambda
1202: . DS_MAT_B - (workspace) T'(lambda) evaluated at a given lambda
1203: - DS_MAT_W - (workspace) eigenvectors of linearization in SLP
1205: Implemented methods:
1206: + 0 - Successive Linear Problems (SLP), computes just one eigenpair
1207: - 1 - Contour integral, computes all eigenvalues inside a region
1209: .seealso: DSCreate(), DSSetType(), DSType, DSNEPSetFN(), DSNEPSetComputeMatrixFunction()
1210: M*/
1211: SLEPC_EXTERN PetscErrorCode DSCreate_NEP(DS ds)
1212: {
1213: DS_NEP *ctx;
1215: PetscNewLog(ds,&ctx);
1216: ds->data = (void*)ctx;
1217: ctx->max_mid = 4;
1218: ctx->nnod = 64;
1219: ctx->Nit = 3;
1220: ctx->rtol = PETSC_MACHINE_EPSILON/PetscSqrtReal(PETSC_SQRT_MACHINE_EPSILON);
1222: ds->ops->allocate = DSAllocate_NEP;
1223: ds->ops->setfromoptions = DSSetFromOptions_NEP;
1224: ds->ops->view = DSView_NEP;
1225: ds->ops->vectors = DSVectors_NEP;
1226: ds->ops->solve[0] = DSSolve_NEP_SLP;
1227: #if defined(PETSC_USE_COMPLEX)
1228: ds->ops->solve[1] = DSSolve_NEP_Contour;
1229: #endif
1230: ds->ops->sort = DSSort_NEP;
1231: ds->ops->synchronize = DSSynchronize_NEP;
1232: ds->ops->destroy = DSDestroy_NEP;
1233: ds->ops->matgetsize = DSMatGetSize_NEP;
1235: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetFN_C",DSNEPSetFN_NEP);
1236: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetFN_C",DSNEPGetFN_NEP);
1237: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetNumFN_C",DSNEPGetNumFN_NEP);
1238: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetMinimality_C",DSNEPGetMinimality_NEP);
1239: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetMinimality_C",DSNEPSetMinimality_NEP);
1240: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetRefine_C",DSNEPGetRefine_NEP);
1241: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetRefine_C",DSNEPSetRefine_NEP);
1242: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetIntegrationPoints_C",DSNEPGetIntegrationPoints_NEP);
1243: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetIntegrationPoints_C",DSNEPSetIntegrationPoints_NEP);
1244: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetSamplingSize_C",DSNEPGetSamplingSize_NEP);
1245: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetSamplingSize_C",DSNEPSetSamplingSize_NEP);
1246: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetRG_C",DSNEPSetRG_NEP);
1247: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetRG_C",DSNEPGetRG_NEP);
1248: PetscObjectComposeFunction((PetscObject)ds,"DSNEPSetComputeMatrixFunction_C",DSNEPSetComputeMatrixFunction_NEP);
1249: PetscObjectComposeFunction((PetscObject)ds,"DSNEPGetComputeMatrixFunction_C",DSNEPGetComputeMatrixFunction_NEP);
1250: PetscFunctionReturn(0);
1251: }