Actual source code: test9.c
slepc-3.7.4 2017-05-17
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 DSGHEP.\n\n";
24: #include <slepcds.h>
28: int main(int argc,char **argv)
29: {
31: DS ds;
32: SlepcSC sc;
33: PetscReal re;
34: PetscScalar *A,*B,*eig;
35: PetscInt i,j,n=10,ld;
36: PetscViewer viewer;
37: PetscBool verbose;
39: SlepcInitialize(&argc,&argv,(char*)0,help);
40: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
41: PetscPrintf(PETSC_COMM_WORLD,"Solve a System of type GHEP - dimension %D.\n",n);
42: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
44: /* Create DS object */
45: DSCreate(PETSC_COMM_WORLD,&ds);
46: DSSetType(ds,DSGHEP);
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: /* Diagonal matrix */
70: for (i=0;i<n;i++) B[i+i*ld]=0.1*(i+1);
71: DSRestoreArray(ds,DS_MAT_A,&A);
72: DSRestoreArray(ds,DS_MAT_B,&B);
73: DSSetState(ds,DS_STATE_RAW);
74: if (verbose) {
75: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
76: DSView(ds,viewer);
77: }
79: /* Solve */
80: PetscMalloc1(n,&eig);
81: DSGetSlepcSC(ds,&sc);
82: sc->comparison = SlepcCompareLargestMagnitude;
83: sc->comparisonctx = NULL;
84: sc->map = NULL;
85: sc->mapobj = NULL;
86: DSSolve(ds,eig,NULL);
87: DSSort(ds,eig,NULL,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");
95: for (i=0;i<n;i++) {
96: re = PetscRealPart(eig[i]);
97: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re);
98: }
99: PetscFree(eig);
100: DSDestroy(&ds);
101: SlepcFinalize();
102: return ierr;
103: }