Actual source code: test23.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: static char help[] = "Test interface functions of DSNEP.\n\n";
13: #include <slepcds.h>
15: int main(int argc,char **argv)
16: {
17: DS ds;
18: FN f1,f2,f3,funs[3];
19: SlepcSC sc;
20: PetscScalar *Id,*A,*B,*wr,*wi,*X,*W,coeffs[2],auxr,alpha;
21: PetscReal tau=0.001,h,a=20,xi,re,im,nrm,aux;
22: PetscInt i,j,ii,jj,k,n=10,ld,nev,nfun,midx,ip,rits,meth,spls;
23: PetscViewer viewer;
24: PetscBool verbose;
25: RG rg;
26: DSMatType mat[3]={DS_MAT_E0,DS_MAT_E1,DS_MAT_E2};
27: #if defined(PETSC_USE_COMPLEX)
28: PetscBool flg;
29: #else
30: PetscScalar auxi;
31: #endif
33: SlepcInitialize(&argc,&argv,(char*)0,help);
34: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
35: PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
36: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NEP - dimension %" PetscInt_FMT ", tau=%g.\n",n,(double)tau);
37: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
39: /* Create DS object and set options */
40: DSCreate(PETSC_COMM_WORLD,&ds);
41: DSSetType(ds,DSNEP);
42: #if defined(PETSC_USE_COMPLEX)
43: DSSetMethod(ds,1); /* contour integral */
44: #endif
45: DSNEPGetRG(ds,&rg);
46: RGSetType(rg,RGELLIPSE);
47: DSNEPSetMinimality(ds,1);
48: DSNEPSetIntegrationPoints(ds,16);
49: DSNEPSetRefine(ds,PETSC_DEFAULT,2);
50: DSNEPSetSamplingSize(ds,25);
51: DSSetFromOptions(ds);
53: /* Print current options */
54: DSGetMethod(ds,&meth);
55: #if defined(PETSC_USE_COMPLEX)
57: RGIsTrivial(rg,&flg);
59: #endif
61: DSNEPGetMinimality(ds,&midx);
62: DSNEPGetIntegrationPoints(ds,&ip);
63: DSNEPGetRefine(ds,NULL,&rits);
64: DSNEPGetSamplingSize(ds,&spls);
65: if (meth==1) {
66: PetscPrintf(PETSC_COMM_WORLD,"Contour integral method with %" PetscInt_FMT " integration points, minimality index %" PetscInt_FMT ", and sampling size %" PetscInt_FMT "\n",ip,midx,spls);
67: if (rits) PetscPrintf(PETSC_COMM_WORLD,"Doing %" PetscInt_FMT " iterations of Newton refinement\n",rits);
68: }
70: /* Set functions (prior to DSAllocate) */
71: FNCreate(PETSC_COMM_WORLD,&f1);
72: FNSetType(f1,FNRATIONAL);
73: coeffs[0] = -1.0; coeffs[1] = 0.0;
74: FNRationalSetNumerator(f1,2,coeffs);
76: FNCreate(PETSC_COMM_WORLD,&f2);
77: FNSetType(f2,FNRATIONAL);
78: coeffs[0] = 1.0;
79: FNRationalSetNumerator(f2,1,coeffs);
81: FNCreate(PETSC_COMM_WORLD,&f3);
82: FNSetType(f3,FNEXP);
83: FNSetScale(f3,-tau,1.0);
85: funs[0] = f1;
86: funs[1] = f2;
87: funs[2] = f3;
88: DSNEPSetFN(ds,3,funs);
90: /* Set dimensions */
91: ld = n;
92: DSAllocate(ds,ld);
93: DSSetDimensions(ds,n,0,0);
95: /* Set up viewer */
96: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
97: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
98: PetscViewerPopFormat(viewer);
99: if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
101: /* Fill matrices */
102: DSGetArray(ds,DS_MAT_E0,&Id);
103: for (i=0;i<n;i++) Id[i+i*ld]=1.0;
104: DSRestoreArray(ds,DS_MAT_E0,&Id);
105: h = PETSC_PI/(PetscReal)(n+1);
106: DSGetArray(ds,DS_MAT_E1,&A);
107: for (i=0;i<n;i++) A[i+i*ld]=-2.0/(h*h)+a;
108: for (i=1;i<n;i++) {
109: A[i+(i-1)*ld]=1.0/(h*h);
110: A[(i-1)+i*ld]=1.0/(h*h);
111: }
112: DSRestoreArray(ds,DS_MAT_E1,&A);
113: DSGetArray(ds,DS_MAT_E2,&B);
114: for (i=0;i<n;i++) {
115: xi = (i+1)*h;
116: B[i+i*ld] = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
117: }
118: DSRestoreArray(ds,DS_MAT_E2,&B);
120: if (verbose) {
121: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
122: DSView(ds,viewer);
123: }
125: /* Solve */
126: PetscCalloc2(n,&wr,n,&wi);
127: DSGetSlepcSC(ds,&sc);
128: sc->comparison = SlepcCompareLargestMagnitude;
129: sc->comparisonctx = NULL;
130: sc->map = NULL;
131: sc->mapobj = NULL;
132: DSSolve(ds,wr,wi);
133: DSSort(ds,wr,wi,NULL,NULL,NULL);
135: if (verbose) {
136: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
137: DSView(ds,viewer);
138: }
139: DSGetDimensions(ds,NULL,NULL,NULL,&nev);
141: /* Print computed eigenvalues */
142: DSNEPGetNumFN(ds,&nfun);
143: PetscMalloc1(ld*ld,&W);
144: DSVectors(ds,DS_MAT_X,NULL,NULL);
145: DSGetArray(ds,DS_MAT_X,&X);
146: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
147: for (i=0;i<nev;i++) {
148: #if defined(PETSC_USE_COMPLEX)
149: re = PetscRealPart(wr[i]);
150: im = PetscImaginaryPart(wr[i]);
151: #else
152: re = wr[i];
153: im = wi[i];
154: #endif
155: /* Residual */
156: PetscArrayzero(W,ld*ld);
157: for (k=0;k<nfun;k++) {
158: FNEvaluateFunction(funs[k],wr[i],&alpha);
159: DSGetArray(ds,mat[k],&A);
160: for (jj=0;jj<n;jj++) for (ii=0;ii<n;ii++) W[jj*ld+ii] += alpha*A[jj*ld+ii];
161: DSRestoreArray(ds,mat[k],&A);
162: }
163: nrm = 0.0;
164: for (k=0;k<n;k++) {
165: auxr = 0.0;
166: #if !defined(PETSC_USE_COMPLEX)
167: auxi = 0.0;
168: #endif
169: for (j=0;j<n;j++) {
170: auxr += W[k+j*ld]*X[i*ld+j];
171: #if !defined(PETSC_USE_COMPLEX)
172: if (PetscAbs(wi[j])!=0.0) auxi += W[k+j*ld]*X[(i+1)*ld+j];
173: #endif
174: }
175: aux = SlepcAbsEigenvalue(auxr,auxi);
176: nrm += aux*aux;
177: }
178: nrm = PetscSqrtReal(nrm);
179: if (nrm>1000*n*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: the residual norm of the %" PetscInt_FMT "-th computed eigenpair %g\n",i,(double)nrm);
180: if (PetscAbs(im)<1e-10) PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re);
181: else PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im);
182: }
183: PetscFree(W);
184: DSRestoreArray(ds,DS_MAT_X,&X);
185: DSRestoreArray(ds,DS_MAT_W,&W);
186: PetscFree2(wr,wi);
187: FNDestroy(&f1);
188: FNDestroy(&f2);
189: FNDestroy(&f3);
190: DSDestroy(&ds);
191: SlepcFinalize();
192: return 0;
193: }
195: /*TEST
197: testset:
198: test:
199: suffix: 1
200: requires: !complex
201: test:
202: suffix: 2
203: args: -ds_nep_rg_ellipse_radius 10
204: filter: sed -e "s/[+-]0\.0*i//g" | sed -e "s/37411/37410/"
205: requires: complex
207: TEST*/