Actual source code: test12.c
slepc-3.16.3 2022-04-11
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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 DSNEP.\n\n";
13: #include <slepcds.h>
15: int main(int argc,char **argv)
16: {
18: DS ds;
19: FN f1,f2,f3,funs[3],qfun;
20: SlepcSC sc;
21: PetscScalar *Id,*A,*B,*wr,*wi,*X,*W,coeffs[2],auxr,alpha;
22: PetscReal tol,tau=0.001,radius=10,h,a=20,xi,re,im,nrm,aux;
23: PetscInt i,j,ii,jj,k,n=10,ld,nev,nfun;
24: PetscViewer viewer;
25: PetscBool verbose;
26: RG rg;
27: DSMatType mat[3]={DS_MAT_E0,DS_MAT_E1,DS_MAT_E2};
28: #if !defined(PETSC_USE_COMPLEX)
29: PetscScalar auxi;
30: #endif
32: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
33: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
34: PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
35: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NEP - dimension %D, tau=%g.\n",n,(double)tau);
36: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
37: PetscOptionsGetReal(NULL,NULL,"-radius",&radius,NULL);
39: /* Create DS object */
40: DSCreate(PETSC_COMM_WORLD,&ds);
41: DSSetType(ds,DSNEP);
42: tol = 1000*n*PETSC_MACHINE_EPSILON;
43: DSNEPSetRefine(ds,tol,PETSC_DECIDE);
44: DSSetFromOptions(ds);
46: /* Set functions (prior to DSAllocate) */
47: FNCreate(PETSC_COMM_WORLD,&f1);
48: FNSetType(f1,FNRATIONAL);
49: coeffs[0] = -1.0; coeffs[1] = 0.0;
50: FNRationalSetNumerator(f1,2,coeffs);
52: FNCreate(PETSC_COMM_WORLD,&f2);
53: FNSetType(f2,FNRATIONAL);
54: coeffs[0] = 1.0;
55: FNRationalSetNumerator(f2,1,coeffs);
57: FNCreate(PETSC_COMM_WORLD,&f3);
58: FNSetType(f3,FNEXP);
59: FNSetScale(f3,-tau,1.0);
61: funs[0] = f1;
62: funs[1] = f2;
63: funs[2] = f3;
64: DSNEPSetFN(ds,3,funs);
66: /* Set dimensions */
67: ld = n+2; /* test leading dimension larger than n */
68: DSAllocate(ds,ld);
69: DSSetDimensions(ds,n,0,0);
71: /* Set region (used only in method=1) */
72: RGCreate(PETSC_COMM_WORLD,&rg);
73: RGSetType(rg,RGELLIPSE);
74: RGEllipseSetParameters(rg,0.0,radius,1.0);
75: DSNEPSetRG(ds,rg);
76: RGDestroy(&rg);
78: /* Set up viewer */
79: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
80: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
81: DSView(ds,viewer);
82: PetscViewerPopFormat(viewer);
83: if (verbose) {
84: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
85: }
87: /* Show info about functions */
88: DSNEPGetNumFN(ds,&nfun);
89: for (i=0;i<nfun;i++) {
90: PetscPrintf(PETSC_COMM_WORLD,"Function %D:\n",i);
91: DSNEPGetFN(ds,i,&qfun);
92: FNView(qfun,NULL);
93: }
95: /* Fill matrices */
96: DSGetArray(ds,DS_MAT_E0,&Id);
97: for (i=0;i<n;i++) Id[i+i*ld]=1.0;
98: DSRestoreArray(ds,DS_MAT_E0,&Id);
99: h = PETSC_PI/(PetscReal)(n+1);
100: DSGetArray(ds,DS_MAT_E1,&A);
101: for (i=0;i<n;i++) A[i+i*ld]=-2.0/(h*h)+a;
102: for (i=1;i<n;i++) {
103: A[i+(i-1)*ld]=1.0/(h*h);
104: A[(i-1)+i*ld]=1.0/(h*h);
105: }
106: DSRestoreArray(ds,DS_MAT_E1,&A);
107: DSGetArray(ds,DS_MAT_E2,&B);
108: for (i=0;i<n;i++) {
109: xi = (i+1)*h;
110: B[i+i*ld] = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
111: }
112: DSRestoreArray(ds,DS_MAT_E2,&B);
114: if (verbose) {
115: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
116: DSView(ds,viewer);
117: }
119: /* Solve */
120: PetscCalloc2(n,&wr,n,&wi);
121: DSGetSlepcSC(ds,&sc);
122: sc->comparison = SlepcCompareLargestMagnitude;
123: sc->comparisonctx = NULL;
124: sc->map = NULL;
125: sc->mapobj = NULL;
126: DSSolve(ds,wr,wi);
127: DSSort(ds,wr,wi,NULL,NULL,NULL);
129: if (verbose) {
130: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
131: DSView(ds,viewer);
132: }
133: DSGetDimensions(ds,NULL,NULL,NULL,&nev);
135: /* Print computed eigenvalues */
136: PetscMalloc1(ld*ld,&W);
137: DSVectors(ds,DS_MAT_X,NULL,NULL);
138: DSGetArray(ds,DS_MAT_X,&X);
139: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
140: for (i=0;i<nev;i++) {
141: #if defined(PETSC_USE_COMPLEX)
142: re = PetscRealPart(wr[i]);
143: im = PetscImaginaryPart(wr[i]);
144: #else
145: re = wr[i];
146: im = wi[i];
147: #endif
148: /* Residual */
149: PetscArrayzero(W,ld*ld);
150: for (k=0;k<nfun;k++) {
151: FNEvaluateFunction(funs[k],wr[i],&alpha);
152: DSGetArray(ds,mat[k],&A);
153: for (jj=0;jj<n;jj++) for (ii=0;ii<n;ii++) W[jj*ld+ii] += alpha*A[jj*ld+ii];
154: DSRestoreArray(ds,mat[k],&A);
155: }
156: nrm = 0.0;
157: for (k=0;k<n;k++) {
158: auxr = 0.0;
159: #if !defined(PETSC_USE_COMPLEX)
160: auxi = 0.0;
161: #endif
162: for (j=0;j<n;j++) {
163: auxr += W[k+j*ld]*X[i*ld+j];
164: #if !defined(PETSC_USE_COMPLEX)
165: if (PetscAbs(wi[j])!=0.0) auxi += W[k+j*ld]*X[(i+1)*ld+j];
166: #endif
167: }
168: aux = SlepcAbsEigenvalue(auxr,auxi);
169: nrm += aux*aux;
170: }
171: nrm = PetscSqrtReal(nrm);
172: if (nrm/SlepcAbsEigenvalue(wr[i],wi[i])>tol) {
173: PetscPrintf(PETSC_COMM_WORLD,"Warning: the residual norm of the %D-th computed eigenpair %g\n",i,(double)nrm);
174: }
175: if (PetscAbs(im)<1e-10) {
176: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re);
177: } else {
178: PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im);
179: }
180: }
181: PetscFree(W);
182: DSRestoreArray(ds,DS_MAT_X,&X);
183: DSRestoreArray(ds,DS_MAT_W,&W);
184: PetscFree2(wr,wi);
185: FNDestroy(&f1);
186: FNDestroy(&f2);
187: FNDestroy(&f3);
188: DSDestroy(&ds);
189: SlepcFinalize();
190: return ierr;
191: }
193: /*TEST
195: testset:
196: test:
197: filter: grep -v "solving the problem"
198: suffix: 1
199: test:
200: suffix: 2
201: args: -ds_method 1 -radius 10 -ds_nep_refine_its 1
202: filter: grep -v "solving the problem" | sed -e "s/[+-]0\.0*i//g" | sed -e "s/37411/37410/" | sed -e "s/tolerance [0-9]\.[0-9]*e[+-]\([0-9]*\)/tolerance removed/" | sed -e "s/tolerance [0-9]\.\([0-9]*\)/tolerance removed/"
203: requires: complex
205: TEST*/