Actual source code: test5.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.\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 re,im;
34: PetscScalar *A,*B,*eigr,*eigi;
35: PetscInt i,j,n=10,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 - dimension %D.\n",n);
42: PetscOptionsHasName(NULL,"-verbose",&verbose);
44: /* Create DS object */
45: DSCreate(PETSC_COMM_WORLD,&ds);
46: DSSetType(ds,DSGHIEP);
47: DSSetFromOptions(ds);
48: ld = n+2; /* test leading dimension larger than n */
49: DSAllocate(ds,ld);
50: DSSetDimensions(ds,n,0,0,0);
52: /* Set up viewer */
53: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
54: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
55: DSView(ds,viewer);
56: PetscViewerPopFormat(viewer);
57: if (verbose) {
58: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
59: }
61: /* Fill with a symmetric Toeplitz matrix */
62: DSGetArray(ds,DS_MAT_A,&A);
63: DSGetArray(ds,DS_MAT_B,&B);
64: for (i=0;i<n;i++) A[i+i*ld]=2.0;
65: for (j=1;j<3;j++) {
66: for (i=0;i<n-j;i++) { A[i+(i+j)*ld]=1.0; A[(i+j)+i*ld]=1.0; }
67: }
68: for (j=1;j<3;j++) { A[0+j*ld]=-1.0*(j+2); A[j+0*ld]=-1.0*(j+2); }
69: /* Signature matrix */
70: for (i=0;i<n;i++) B[i+i*ld]=1.0;
71: B[0] = -1.0;
72: B[n-1+(n-1)*ld] = -1.0;
73: DSRestoreArray(ds,DS_MAT_A,&A);
74: DSRestoreArray(ds,DS_MAT_B,&B);
75: DSSetState(ds,DS_STATE_RAW);
76: if (verbose) {
77: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
78: DSView(ds,viewer);
79: }
81: /* Solve */
82: PetscMalloc(n*sizeof(PetscScalar),&eigr);
83: PetscMalloc(n*sizeof(PetscScalar),&eigi);
84: PetscMemzero(eigi,n*sizeof(PetscScalar));
85: DSSetEigenvalueComparison(ds,SlepcCompareLargestMagnitude,NULL);
86: DSSolve(ds,eigr,eigi);
87: DSSort(ds,eigr,eigi,NULL,NULL,NULL);
88: if (verbose) {
89: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
90: DSView(ds,viewer);
91: }
93: /* Print eigenvalues */
94: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n",n);
95: for (i=0;i<n;i++) {
96: #if defined(PETSC_USE_COMPLEX)
97: re = PetscRealPart(eigr[i]);
98: im = PetscImaginaryPart(eigr[i]);
99: #else
100: re = eigr[i];
101: im = eigi[i];
102: #endif
103: if (PetscAbs(im)<1e-10) {
104: PetscViewerASCIIPrintf(viewer," %.5F\n",re);
105: } else {
106: PetscViewerASCIIPrintf(viewer," %.5F%+.5Fi\n",re,im);
107: }
108: }
109: PetscFree(eigr);
110: PetscFree(eigi);
111: DSDestroy(&ds);
112: SlepcFinalize();
113: return 0;
114: }