Actual source code: test2.c
slepc-3.6.3 2016-03-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, 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 DSHEP.\n\n";
24: #include <slepcds.h>
28: int main(int argc,char **argv)
29: {
31: DS ds;
32: SlepcSC sc;
33: PetscScalar *A,*eig;
34: PetscInt i,j,n=10,ld;
35: PetscViewer viewer;
36: PetscBool verbose,extrarow;
38: SlepcInitialize(&argc,&argv,(char*)0,help);
39: PetscOptionsGetInt(NULL,"-n",&n,NULL);
40: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type HEP - dimension %D.\n",n);
41: PetscOptionsHasName(NULL,"-verbose",&verbose);
42: PetscOptionsHasName(NULL,"-extrarow",&extrarow);
44: /* Create DS object */
45: DSCreate(PETSC_COMM_WORLD,&ds);
46: DSSetType(ds,DSHEP);
47: DSSetFromOptions(ds);
48: ld = n+2; /* test leading dimension larger than n */
49: DSAllocate(ds,ld);
50: DSSetDimensions(ds,n,0,0,0);
51: DSSetExtraRow(ds,extrarow);
53: /* Set up viewer */
54: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
55: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
56: DSView(ds,viewer);
57: PetscViewerPopFormat(viewer);
58: if (verbose) {
59: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
60: }
62: /* Fill with a symmetric Toeplitz matrix */
63: DSGetArray(ds,DS_MAT_A,&A);
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: if (extrarow) { A[n+(n-2)*ld]=1.0; A[n+(n-1)*ld]=1.0; }
69: DSRestoreArray(ds,DS_MAT_A,&A);
70: DSSetState(ds,DS_STATE_RAW);
71: if (verbose) {
72: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
73: DSView(ds,viewer);
74: }
76: /* Solve */
77: PetscMalloc1(n,&eig);
78: DSGetSlepcSC(ds,&sc);
79: sc->comparison = SlepcCompareLargestMagnitude;
80: sc->comparisonctx = NULL;
81: sc->map = NULL;
82: sc->mapobj = NULL;
83: DSSolve(ds,eig,NULL);
84: DSSort(ds,eig,NULL,NULL,NULL,NULL);
85: if (extrarow) { DSUpdateExtraRow(ds); }
86: if (verbose) {
87: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
88: DSView(ds,viewer);
89: }
91: /* Print eigenvalues */
92: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
93: for (i=0;i<n;i++) {
94: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)PetscRealPart(eig[i]));
95: }
97: PetscFree(eig);
98: DSDestroy(&ds);
99: SlepcFinalize();
100: return 0;
101: }