Actual source code: test12.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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,coeffs[2];
22: PetscReal tau=0.001,h,a=20,xi,re,im,nrm,aux;
23: PetscInt i,j,n=10,ld,nev,nfun;
24: PetscViewer viewer;
25: PetscBool verbose;
27: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
28: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
29: PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
30: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NEP - dimension %D, tau=%g.\n",n,(double)tau);
31: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
33: /* Create DS object */
34: DSCreate(PETSC_COMM_WORLD,&ds);
35: DSSetType(ds,DSNEP);
36: DSSetFromOptions(ds);
38: /* Set functions (prior to DSAllocate) */
39: FNCreate(PETSC_COMM_WORLD,&f1);
40: FNSetType(f1,FNRATIONAL);
41: coeffs[0] = -1.0; coeffs[1] = 0.0;
42: FNRationalSetNumerator(f1,2,coeffs);
44: FNCreate(PETSC_COMM_WORLD,&f2);
45: FNSetType(f2,FNRATIONAL);
46: coeffs[0] = 1.0;
47: FNRationalSetNumerator(f2,1,coeffs);
49: FNCreate(PETSC_COMM_WORLD,&f3);
50: FNSetType(f3,FNEXP);
51: FNSetScale(f3,-tau,1.0);
53: funs[0] = f1;
54: funs[1] = f2;
55: funs[2] = f3;
56: DSNEPSetFN(ds,3,funs);
58: /* Set dimensions */
59: ld = n+2; /* test leading dimension larger than n */
60: DSAllocate(ds,ld);
61: DSSetDimensions(ds,n,0,0,0);
63: /* Set up viewer */
64: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
65: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
66: DSView(ds,viewer);
67: PetscViewerPopFormat(viewer);
68: if (verbose) {
69: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
70: }
72: /* Show info about functions */
73: DSNEPGetNumFN(ds,&nfun);
74: for (i=0;i<nfun;i++) {
75: PetscPrintf(PETSC_COMM_WORLD,"Function %D:\n",i);
76: DSNEPGetFN(ds,i,&qfun);
77: FNView(qfun,NULL);
78: }
80: /* Fill matrices */
81: DSGetArray(ds,DS_MAT_E0,&Id);
82: for (i=0;i<n;i++) Id[i+i*ld]=1.0;
83: DSRestoreArray(ds,DS_MAT_E0,&Id);
84: h = PETSC_PI/(PetscReal)(n+1);
85: DSGetArray(ds,DS_MAT_E1,&A);
86: for (i=0;i<n;i++) A[i+i*ld]=-2.0/(h*h)+a;
87: for (i=1;i<n;i++) {
88: A[i+(i-1)*ld]=1.0/(h*h);
89: A[(i-1)+i*ld]=1.0/(h*h);
90: }
91: DSRestoreArray(ds,DS_MAT_E1,&A);
92: DSGetArray(ds,DS_MAT_E2,&B);
93: for (i=0;i<n;i++) {
94: xi = (i+1)*h;
95: B[i+i*ld] = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
96: }
97: DSRestoreArray(ds,DS_MAT_E2,&B);
99: if (verbose) {
100: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
101: DSView(ds,viewer);
102: }
104: /* Solve */
105: PetscMalloc2(n,&wr,n,&wi);
106: DSGetSlepcSC(ds,&sc);
107: sc->comparison = SlepcCompareLargestMagnitude;
108: sc->comparisonctx = NULL;
109: sc->map = NULL;
110: sc->mapobj = NULL;
111: DSSolve(ds,wr,wi);
112: DSSort(ds,wr,wi,NULL,NULL,NULL);
114: if (verbose) {
115: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
116: DSView(ds,viewer);
117: }
119: /* Print first eigenvalue */
120: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalue =\n");
121: nev = 1;
122: for (i=0;i<nev;i++) {
123: #if defined(PETSC_USE_COMPLEX)
124: re = PetscRealPart(wr[i]);
125: im = PetscImaginaryPart(wr[i]);
126: #else
127: re = wr[i];
128: im = wi[i];
129: #endif
130: if (PetscAbs(im)<1e-10) {
131: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re);
132: } else {
133: PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im);
134: }
135: }
137: /* Eigenvectors */
138: DSVectors(ds,DS_MAT_X,NULL,NULL);
139: j = 0;
140: nrm = 0.0;
141: DSGetArray(ds,DS_MAT_X,&X);
142: for (i=0;i<n;i++) {
143: #if defined(PETSC_USE_COMPLEX)
144: aux = PetscAbsScalar(X[i+j*ld]);
145: #else
146: if (PetscAbs(wi[j])==0.0) aux = PetscAbsScalar(X[i+j*ld]);
147: else aux = SlepcAbsEigenvalue(X[i+j*ld],X[i+(j+1)*ld]);
148: #endif
149: nrm += aux*aux;
150: }
151: DSRestoreArray(ds,DS_MAT_X,&X);
152: nrm = PetscSqrtReal(nrm);
153: PetscPrintf(PETSC_COMM_WORLD,"Norm of eigenvector = %.3f\n",(double)nrm);
154: if (verbose) {
155: PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n");
156: DSView(ds,viewer);
157: }
159: PetscFree2(wr,wi);
160: FNDestroy(&f1);
161: FNDestroy(&f2);
162: FNDestroy(&f3);
163: DSDestroy(&ds);
164: SlepcFinalize();
165: return ierr;
166: }
168: /*TEST
170: test:
171: suffix: 1
173: TEST*/