Actual source code: test12.c
slepc-3.5.2 2014-10-10
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2014, 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,"-n",&n,NULL);
41: PetscOptionsGetReal(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,"-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: FNSetParameters(f1,2,coeffs,0,NULL);
56: FNCreate(PETSC_COMM_WORLD,&f2);
57: FNSetType(f2,FNRATIONAL);
58: coeffs[0] = 1.0;
59: FNSetParameters(f2,1,coeffs,0,NULL);
61: FNCreate(PETSC_COMM_WORLD,&f3);
62: FNSetType(f3,FNEXP);
63: coeffs[0] = -tau;
64: FNSetParameters(f3,1,coeffs,0,NULL);
66: funs[0] = f1;
67: funs[1] = f2;
68: funs[2] = f3;
69: DSSetFN(ds,3,funs);
71: /* Set dimensions */
72: ld = n+2; /* test leading dimension larger than n */
73: DSAllocate(ds,ld);
74: DSSetDimensions(ds,n,0,0,0);
76: /* Set up viewer */
77: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
78: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
79: DSView(ds,viewer);
80: PetscViewerPopFormat(viewer);
81: if (verbose) {
82: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
83: }
85: /* Fill matrices */
86: DSGetArray(ds,DS_MAT_E0,&Id);
87: for (i=0;i<n;i++) Id[i+i*ld]=1.0;
88: DSRestoreArray(ds,DS_MAT_E0,&Id);
89: h = PETSC_PI/(PetscReal)(n+1);
90: DSGetArray(ds,DS_MAT_E1,&A);
91: for (i=0;i<n;i++) A[i+i*ld]=-2.0/(h*h)+a;
92: for (i=1;i<n;i++) {
93: A[i+(i-1)*ld]=1.0/(h*h);
94: A[(i-1)+i*ld]=1.0/(h*h);
95: }
96: DSRestoreArray(ds,DS_MAT_E1,&A);
97: DSGetArray(ds,DS_MAT_E2,&B);
98: for (i=0;i<n;i++) {
99: xi = (i+1)*h;
100: B[i+i*ld] = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
101: }
102: DSRestoreArray(ds,DS_MAT_E2,&B);
104: if (verbose) {
105: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
106: DSView(ds,viewer);
107: }
109: /* Solve */
110: PetscMalloc2(n,&wr,n,&wi);
111: DSSolve(ds,wr,wi);
112: if (verbose) {
113: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
114: DSView(ds,viewer);
115: }
117: /* Print first eigenvalue */
118: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalue =\n",n);
119: nev = 1;
120: for (i=0;i<nev;i++) {
121: #if defined(PETSC_USE_COMPLEX)
122: re = PetscRealPart(wr[i]);
123: im = PetscImaginaryPart(wr[i]);
124: #else
125: re = wr[i];
126: im = wi[i];
127: #endif
128: if (PetscAbs(im)<1e-10) {
129: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re);
130: } else {
131: PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im);
132: }
133: }
135: PetscFree2(wr,wi);
136: FNDestroy(&f1);
137: FNDestroy(&f2);
138: FNDestroy(&f3);
139: DSDestroy(&ds);
140: SlepcFinalize();
141: return 0;
142: }