Actual source code: test23.c
slepc-3.18.2 2023-01-26
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
34: SlepcInitialize(&argc,&argv,(char*)0,help);
35: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
36: PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
37: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NEP - dimension %" PetscInt_FMT ", tau=%g.\n",n,(double)tau);
38: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
40: /* Create DS object and set options */
41: DSCreate(PETSC_COMM_WORLD,&ds);
42: DSSetType(ds,DSNEP);
43: #if defined(PETSC_USE_COMPLEX)
44: DSSetMethod(ds,1); /* contour integral */
45: #endif
46: DSNEPGetRG(ds,&rg);
47: RGSetType(rg,RGELLIPSE);
48: DSNEPSetMinimality(ds,1);
49: DSNEPSetIntegrationPoints(ds,16);
50: DSNEPSetRefine(ds,PETSC_DEFAULT,2);
51: DSNEPSetSamplingSize(ds,25);
52: DSSetFromOptions(ds);
54: /* Print current options */
55: DSGetMethod(ds,&meth);
56: #if defined(PETSC_USE_COMPLEX)
58: RGIsTrivial(rg,&flg);
60: #endif
62: DSNEPGetMinimality(ds,&midx);
63: DSNEPGetIntegrationPoints(ds,&ip);
64: DSNEPGetRefine(ds,NULL,&rits);
65: DSNEPGetSamplingSize(ds,&spls);
66: if (meth==1) {
67: 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);
68: if (rits) PetscPrintf(PETSC_COMM_WORLD,"Doing %" PetscInt_FMT " iterations of Newton refinement\n",rits);
69: }
71: /* Set functions (prior to DSAllocate) */
72: FNCreate(PETSC_COMM_WORLD,&f1);
73: FNSetType(f1,FNRATIONAL);
74: coeffs[0] = -1.0; coeffs[1] = 0.0;
75: FNRationalSetNumerator(f1,2,coeffs);
77: FNCreate(PETSC_COMM_WORLD,&f2);
78: FNSetType(f2,FNRATIONAL);
79: coeffs[0] = 1.0;
80: FNRationalSetNumerator(f2,1,coeffs);
82: FNCreate(PETSC_COMM_WORLD,&f3);
83: FNSetType(f3,FNEXP);
84: FNSetScale(f3,-tau,1.0);
86: funs[0] = f1;
87: funs[1] = f2;
88: funs[2] = f3;
89: DSNEPSetFN(ds,3,funs);
91: /* Set dimensions */
92: ld = n;
93: DSAllocate(ds,ld);
94: DSSetDimensions(ds,n,0,0);
96: /* Set up viewer */
97: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
98: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
99: PetscViewerPopFormat(viewer);
100: if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
102: /* Fill matrices */
103: DSGetArray(ds,DS_MAT_E0,&Id);
104: for (i=0;i<n;i++) Id[i+i*ld]=1.0;
105: DSRestoreArray(ds,DS_MAT_E0,&Id);
106: h = PETSC_PI/(PetscReal)(n+1);
107: DSGetArray(ds,DS_MAT_E1,&A);
108: for (i=0;i<n;i++) A[i+i*ld]=-2.0/(h*h)+a;
109: for (i=1;i<n;i++) {
110: A[i+(i-1)*ld]=1.0/(h*h);
111: A[(i-1)+i*ld]=1.0/(h*h);
112: }
113: DSRestoreArray(ds,DS_MAT_E1,&A);
114: DSGetArray(ds,DS_MAT_E2,&B);
115: for (i=0;i<n;i++) {
116: xi = (i+1)*h;
117: B[i+i*ld] = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
118: }
119: DSRestoreArray(ds,DS_MAT_E2,&B);
121: if (verbose) {
122: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
123: DSView(ds,viewer);
124: }
126: /* Solve */
127: PetscCalloc2(n,&wr,n,&wi);
128: DSGetSlepcSC(ds,&sc);
129: sc->comparison = SlepcCompareLargestMagnitude;
130: sc->comparisonctx = NULL;
131: sc->map = NULL;
132: sc->mapobj = NULL;
133: DSSolve(ds,wr,wi);
134: DSSort(ds,wr,wi,NULL,NULL,NULL);
136: if (verbose) {
137: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
138: DSView(ds,viewer);
139: }
140: DSGetDimensions(ds,NULL,NULL,NULL,&nev);
142: /* Print computed eigenvalues */
143: DSNEPGetNumFN(ds,&nfun);
144: PetscMalloc1(ld*ld,&W);
145: DSVectors(ds,DS_MAT_X,NULL,NULL);
146: DSGetArray(ds,DS_MAT_X,&X);
147: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
148: for (i=0;i<nev;i++) {
149: #if defined(PETSC_USE_COMPLEX)
150: re = PetscRealPart(wr[i]);
151: im = PetscImaginaryPart(wr[i]);
152: #else
153: re = wr[i];
154: im = wi[i];
155: #endif
156: /* Residual */
157: PetscArrayzero(W,ld*ld);
158: for (k=0;k<nfun;k++) {
159: FNEvaluateFunction(funs[k],wr[i],&alpha);
160: DSGetArray(ds,mat[k],&A);
161: for (jj=0;jj<n;jj++) for (ii=0;ii<n;ii++) W[jj*ld+ii] += alpha*A[jj*ld+ii];
162: DSRestoreArray(ds,mat[k],&A);
163: }
164: nrm = 0.0;
165: for (k=0;k<n;k++) {
166: auxr = 0.0;
167: #if !defined(PETSC_USE_COMPLEX)
168: auxi = 0.0;
169: #endif
170: for (j=0;j<n;j++) {
171: auxr += W[k+j*ld]*X[i*ld+j];
172: #if !defined(PETSC_USE_COMPLEX)
173: if (PetscAbs(wi[j])!=0.0) auxi += W[k+j*ld]*X[(i+1)*ld+j];
174: #endif
175: }
176: aux = SlepcAbsEigenvalue(auxr,auxi);
177: nrm += aux*aux;
178: }
179: nrm = PetscSqrtReal(nrm);
180: 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);
181: if (PetscAbs(im)<1e-10) PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re);
182: else PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im);
183: }
184: PetscFree(W);
185: DSRestoreArray(ds,DS_MAT_X,&X);
186: DSRestoreArray(ds,DS_MAT_W,&W);
187: PetscFree2(wr,wi);
188: FNDestroy(&f1);
189: FNDestroy(&f2);
190: FNDestroy(&f3);
191: DSDestroy(&ds);
192: SlepcFinalize();
193: return 0;
194: }
196: /*TEST
198: testset:
199: test:
200: suffix: 1
201: requires: !complex
202: test:
203: suffix: 2
204: args: -ds_nep_rg_ellipse_radius 10
205: filter: sed -e "s/[+-]0\.0*i//g" | sed -e "s/37411/37410/"
206: requires: complex
208: TEST*/