Actual source code: test13.c
slepc-3.7.3 2016-09-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, 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 with block size larger than one.\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,ld,bs,maxbw=3,nblks=8;
35: PetscViewer viewer;
36: PetscBool verbose;
38: SlepcInitialize(&argc,&argv,(char*)0,help);
39: PetscOptionsGetInt(NULL,NULL,"-maxbw",&maxbw,NULL);
40: PetscOptionsGetInt(NULL,NULL,"-nblks",&nblks,NULL);
41: n = maxbw*nblks;
42: bs = maxbw;
43: PetscPrintf(PETSC_COMM_WORLD,"Solve a block HEP Dense System - dimension %D (bandwidth=%D, blocks=%D).\n",n,maxbw,nblks);
44: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
46: /* Create DS object */
47: DSCreate(PETSC_COMM_WORLD,&ds);
48: DSSetType(ds,DSHEP);
49: DSSetMethod(ds,3); /* Select block divide-and-conquer */
50: DSSetBlockSize(ds,bs);
51: DSSetFromOptions(ds);
52: ld = n;
53: DSAllocate(ds,ld);
54: DSSetDimensions(ds,n,0,0,0);
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 with a symmetric band Toeplitz matrix */
66: DSGetArray(ds,DS_MAT_A,&A);
67: for (i=0;i<n;i++) A[i+i*ld]=2.0;
68: for (j=1;j<=bs;j++) {
69: for (i=0;i<n-j;i++) { A[i+(i+j)*ld]=1.0; A[(i+j)+i*ld]=1.0; }
70: }
71: DSRestoreArray(ds,DS_MAT_A,&A);
72: DSSetState(ds,DS_STATE_RAW);
73: if (verbose) {
74: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
75: DSView(ds,viewer);
76: }
78: /* Solve */
79: PetscMalloc1(n,&eig);
80: DSGetSlepcSC(ds,&sc);
81: sc->comparison = SlepcCompareSmallestReal;
82: sc->comparisonctx = NULL;
83: sc->map = NULL;
84: sc->mapobj = NULL;
85: DSSolve(ds,eig,NULL);
86: DSSort(ds,eig,NULL,NULL,NULL,NULL);
87: if (verbose) {
88: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
89: DSView(ds,viewer);
90: }
92: /* Print eigenvalues */
93: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
94: for (i=0;i<n;i++) {
95: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)PetscRealPart(eig[i]));
96: }
98: PetscFree(eig);
99: DSDestroy(&ds);
100: SlepcFinalize();
101: return ierr;
102: }