Actual source code: test3.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 DSHEP with compact storage.\n\n";
24: #include slepcds.h
28: int main( int argc, char **argv )
29: {
31: DS ds;
32: PetscReal *T;
33: PetscScalar *eig;
34: PetscInt i,n=10,l=2,k=5,ld;
35: PetscViewer viewer;
36: PetscBool verbose,extrarow;
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 HEP 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);
45: PetscOptionsHasName(PETSC_NULL,"-extrarow",&extrarow);
47: /* Create DS object */
48: DSCreate(PETSC_COMM_WORLD,&ds);
49: DSSetType(ds,DSHEP);
50: DSSetFromOptions(ds);
51: ld = n+2; /* test leading dimension larger than n */
52: DSAllocate(ds,ld);
53: DSSetDimensions(ds,n,PETSC_IGNORE,l,k);
54: DSSetCompact(ds,PETSC_TRUE);
55: DSSetExtraRow(ds,extrarow);
57: /* Set up viewer */
58: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
59: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
60: DSView(ds,viewer);
61: PetscViewerPopFormat(viewer);
62: if (verbose) {
63: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
64: }
66: /* Fill arrow-tridiagonal matrix */
67: DSGetArrayReal(ds,DS_MAT_T,&T);
68: for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
69: for (i=l;i<n-1;i++) T[i+ld] = 1.0;
70: if (extrarow) T[n-1+ld] = 1.0;
71: DSRestoreArrayReal(ds,DS_MAT_T,&T);
72: if (l==0 && k==0) {
73: DSSetState(ds,DS_STATE_INTERMEDIATE);
74: } else {
75: DSSetState(ds,DS_STATE_RAW);
76: }
77: if (verbose) {
78: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
79: DSView(ds,viewer);
80: }
82: /* Solve */
83: PetscMalloc(n*sizeof(PetscScalar),&eig);
84: DSSetEigenvalueComparison(ds,SlepcCompareLargestMagnitude,PETSC_NULL);
85: DSSolve(ds,eig,PETSC_NULL);
86: DSSort(ds,eig,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
87: if (extrarow) { DSUpdateExtraRow(ds); }
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: PetscViewerASCIIPrintf(viewer," %.5F\n",PetscRealPart(eig[i]));
97: }
99: PetscFree(eig);
100: DSDestroy(&ds);
101: SlepcFinalize();
102: return 0;
103: }