Actual source code: test12.c
slepc-3.7.2 2016-07-19
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Test DSNEP.\n\n";
24: #include <slepcds.h>
28: int main(int argc,char **argv)
29: {
31: DS ds;
32: FN f1,f2,f3,funs[3];
33: PetscScalar *Id,*A,*B,*wr,*wi,coeffs[2];
34: PetscReal tau=0.001,h,a=20,xi,re,im;
35: PetscInt i,n=10,ld,nev;
36: PetscViewer viewer;
37: PetscBool verbose;
39: SlepcInitialize(&argc,&argv,(char*)0,help);
40: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
41: PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
42: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NEP - dimension %D, tau=%g.\n",n,(double)tau);
43: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
45: /* Create DS object */
46: DSCreate(PETSC_COMM_WORLD,&ds);
47: DSSetType(ds,DSNEP);
48: DSSetFromOptions(ds);
50: /* Set functions (prior to DSAllocate) */
51: FNCreate(PETSC_COMM_WORLD,&f1);
52: FNSetType(f1,FNRATIONAL);
53: coeffs[0] = -1.0; coeffs[1] = 0.0;
54: FNRationalSetNumerator(f1,2,coeffs);
56: FNCreate(PETSC_COMM_WORLD,&f2);
57: FNSetType(f2,FNRATIONAL);
58: coeffs[0] = 1.0;
59: FNRationalSetNumerator(f2,1,coeffs);
61: FNCreate(PETSC_COMM_WORLD,&f3);
62: FNSetType(f3,FNEXP);
63: FNSetScale(f3,-tau,1.0);
65: funs[0] = f1;
66: funs[1] = f2;
67: funs[2] = f3;
68: DSNEPSetFN(ds,3,funs);
70: /* Set dimensions */
71: ld = n+2; /* test leading dimension larger than n */
72: DSAllocate(ds,ld);
73: DSSetDimensions(ds,n,0,0,0);
75: /* Set up viewer */
76: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
77: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
78: DSView(ds,viewer);
79: PetscViewerPopFormat(viewer);
80: if (verbose) {
81: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
82: }
84: /* Fill matrices */
85: DSGetArray(ds,DS_MAT_E0,&Id);
86: for (i=0;i<n;i++) Id[i+i*ld]=1.0;
87: DSRestoreArray(ds,DS_MAT_E0,&Id);
88: h = PETSC_PI/(PetscReal)(n+1);
89: DSGetArray(ds,DS_MAT_E1,&A);
90: for (i=0;i<n;i++) A[i+i*ld]=-2.0/(h*h)+a;
91: for (i=1;i<n;i++) {
92: A[i+(i-1)*ld]=1.0/(h*h);
93: A[(i-1)+i*ld]=1.0/(h*h);
94: }
95: DSRestoreArray(ds,DS_MAT_E1,&A);
96: DSGetArray(ds,DS_MAT_E2,&B);
97: for (i=0;i<n;i++) {
98: xi = (i+1)*h;
99: B[i+i*ld] = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
100: }
101: DSRestoreArray(ds,DS_MAT_E2,&B);
103: if (verbose) {
104: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
105: DSView(ds,viewer);
106: }
108: /* Solve */
109: PetscMalloc2(n,&wr,n,&wi);
110: DSSolve(ds,wr,wi);
111: if (verbose) {
112: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
113: DSView(ds,viewer);
114: }
116: /* Print first eigenvalue */
117: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalue =\n");
118: nev = 1;
119: for (i=0;i<nev;i++) {
120: #if defined(PETSC_USE_COMPLEX)
121: re = PetscRealPart(wr[i]);
122: im = PetscImaginaryPart(wr[i]);
123: #else
124: re = wr[i];
125: im = wi[i];
126: #endif
127: if (PetscAbs(im)<1e-10) {
128: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re);
129: } else {
130: PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im);
131: }
132: }
134: PetscFree2(wr,wi);
135: FNDestroy(&f1);
136: FNDestroy(&f2);
137: FNDestroy(&f3);
138: DSDestroy(&ds);
139: SlepcFinalize();
140: return ierr;
141: }