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