Actual source code: test9.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 DSGHEP.\n\n";
24: #include slepcds.h
28: int main( int argc, char **argv )
29: {
31: DS ds;
32: PetscReal re;
33: PetscScalar *A,*B,*eig;
34: PetscInt i,j,n=10,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 System of type GHEP - dimension %D.\n",n);
41: PetscOptionsHasName(PETSC_NULL,"-verbose",&verbose);
43: /* Create DS object */
44: DSCreate(PETSC_COMM_WORLD,&ds);
45: DSSetType(ds,DSGHEP);
46: DSSetFromOptions(ds);
47: ld = n+2; /* test leading dimension larger than n */
48: DSAllocate(ds,ld);
49: DSSetDimensions(ds,n,PETSC_IGNORE,0,0);
51: /* Set up viewer */
52: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
53: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
54: DSView(ds,viewer);
55: PetscViewerPopFormat(viewer);
56: if (verbose) {
57: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
58: }
60: /* Fill with a symmetric Toeplitz matrix */
61: DSGetArray(ds,DS_MAT_A,&A);
62: DSGetArray(ds,DS_MAT_B,&B);
63: for (i=0;i<n;i++) A[i+i*ld]=2.0;
64: for (j=1;j<3;j++) {
65: for (i=0;i<n-j;i++) { A[i+(i+j)*ld]=1.0; A[(i+j)+i*ld]=1.0; }
66: }
67: for (j=1;j<3;j++) { A[0+j*ld]=-1.0*(j+2); A[j+0*ld]=-1.0*(j+2); }
68: /* Diagonal matrix */
69: for (i=0;i<n;i++) B[i+i*ld]=0.1*(i+1);
70: DSRestoreArray(ds,DS_MAT_A,&A);
71: DSRestoreArray(ds,DS_MAT_B,&B);
72: DSSetState(ds,DS_STATE_RAW);
73: if (verbose) {
74: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
75: DSView(ds,viewer);
76: }
78: /* Solve */
79: PetscMalloc(n*sizeof(PetscScalar),&eig);
80: DSSetEigenvalueComparison(ds,SlepcCompareLargestMagnitude,PETSC_NULL);
81: DSSolve(ds,eig,PETSC_NULL);
82: DSSort(ds,eig,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
83: if (verbose) {
84: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
85: DSView(ds,viewer);
86: }
87:
88: /* Print eigenvalues */
89: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n",n);
90: for (i=0;i<n;i++) {
91: re = PetscRealPart(eig[i]);
92: PetscViewerASCIIPrintf(viewer," %.5F\n",re);
93: }
94: PetscFree(eig);
95: DSDestroy(&ds);
96: SlepcFinalize();
97: return 0;
98: }