Actual source code: test6.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7:
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 DSGHIEP with compact storage.\n\n";
24: #include slepcds.h
28: int main( int argc, char **argv )
29: {
31: DS ds;
32: PetscReal *T,*s,re,im;
33: PetscScalar *eigr,*eigi;
34: PetscInt i,n=10,l=2,k=5,ld;
35: PetscViewer viewer;
36: PetscBool verbose;
38: SlepcInitialize(&argc,&argv,(char*)0,help);
39: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
40: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GHIEP with compact storage - dimension %D.\n",n);
41: PetscOptionsGetInt(PETSC_NULL,"-l",&l,PETSC_NULL);
42: PetscOptionsGetInt(PETSC_NULL,"-k",&k,PETSC_NULL);
43: if (l>n || k>n || l>k) SETERRQ(PETSC_COMM_WORLD,1,"Wrong value of dimensions");
44: PetscOptionsHasName(PETSC_NULL,"-verbose",&verbose);
46: /* Create DS object */
47: DSCreate(PETSC_COMM_WORLD,&ds);
48: DSSetType(ds,DSGHIEP);
49: DSSetFromOptions(ds);
50: ld = n+2; /* test leading dimension larger than n */
51: DSAllocate(ds,ld);
52: DSSetDimensions(ds,n,PETSC_IGNORE,l,k);
53: DSSetCompact(ds,PETSC_TRUE);
55: /* Set up viewer */
56: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
57: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
58: DSView(ds,viewer);
59: PetscViewerPopFormat(viewer);
60: if (verbose) {
61: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
62: }
64: /* Fill arrow-tridiagonal matrix */
65: DSGetArrayReal(ds,DS_MAT_T,&T);
66: DSGetArrayReal(ds,DS_MAT_D,&s);
67: for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
68: for (i=k;i<n-1;i++) T[i+ld] = 1.0;
69: for (i=l;i<k;i++) T[i+2*ld] = 1.0;
70: T[2*ld+l+1] = -7; T[ld+k+1] = -7;
71: /* Signature matrix */
72: for (i=0;i<n;i++) s[i] = 1.0;
73: s[l+1] = -1.0;
74: s[k+1] = -1.0;
75: DSRestoreArrayReal(ds,DS_MAT_T,&T);
76: DSRestoreArrayReal(ds,DS_MAT_D,&s);
77: if (l==0 && k==0) {
78: DSSetState(ds,DS_STATE_INTERMEDIATE);
79: } else {
80: DSSetState(ds,DS_STATE_RAW);
81: }
82: if (verbose) {
83: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
84: DSView(ds,viewer);
85: }
87: /* Solve */
88: PetscMalloc(n*sizeof(PetscScalar),&eigr);
89: PetscMalloc(n*sizeof(PetscScalar),&eigi);
90: PetscMemzero(eigi,n*sizeof(PetscScalar));
91: DSSetEigenvalueComparison(ds,SlepcCompareLargestMagnitude,PETSC_NULL);
92: DSSolve(ds,eigr,eigi);
93: DSSort(ds,eigr,eigi,PETSC_NULL,PETSC_NULL,PETSC_NULL);
94: if (verbose) {
95: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
96: DSView(ds,viewer);
97: }
98:
99: /* Print eigenvalues */
100: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n",n);
101: for (i=0;i<n;i++) {
102: #if defined(PETSC_USE_COMPLEX)
103: re = PetscRealPart(eigr[i]);
104: im = PetscImaginaryPart(eigr[i]);
105: #else
106: re = eigr[i];
107: im = eigi[i];
108: #endif
109: if (PetscAbs(im)<1e-10) { PetscViewerASCIIPrintf(viewer," %.5F\n",re); }
110: else { PetscViewerASCIIPrintf(viewer," %.5F%+.5Fi\n",re,im); }
111: }
112: PetscFree(eigr);
113: PetscFree(eigi);
114: DSDestroy(&ds);
115: SlepcFinalize();
116: return 0;
117: }