Actual source code: test4.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2013, 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 DSGNHEP.\n\n";
24: #include <slepcds.h>
25: #include <slepc-private/dsimpl.h> /* for the definition of SlepcCompare* */
29: int main(int argc,char **argv)
30: {
32: DS ds;
33: PetscScalar *A,*B,*wr,*wi;
34: PetscReal re,im;
35: PetscInt i,j,n=10,ld;
36: PetscViewer viewer;
37: PetscBool verbose;
39: SlepcInitialize(&argc,&argv,(char*)0,help);
40: PetscOptionsGetInt(NULL,"-n",&n,NULL);
41: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GNHEP - dimension %D.\n",n);
42: PetscOptionsHasName(NULL,"-verbose",&verbose);
44: /* Create DS object */
45: DSCreate(PETSC_COMM_WORLD,&ds);
46: DSSetType(ds,DSGNHEP);
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 A with Grcar matrix */
62: DSGetArray(ds,DS_MAT_A,&A);
63: PetscMemzero(A,sizeof(PetscScalar)*ld*n);
64: for (i=1;i<n;i++) A[i+(i-1)*ld]=-1.0;
65: for (j=0;j<4;j++) {
66: for (i=0;i<n-j;i++) A[i+(i+j)*ld]=1.0;
67: }
68: DSRestoreArray(ds,DS_MAT_A,&A);
69: /* Fill B with an identity matrix */
70: DSGetArray(ds,DS_MAT_B,&B);
71: PetscMemzero(B,sizeof(PetscScalar)*ld*n);
72: for (i=0;i<n;i++) B[i+i*ld]=1.0;
73: DSRestoreArray(ds,DS_MAT_B,&B);
75: if (verbose) {
76: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
77: DSView(ds,viewer);
78: }
80: /* Solve */
81: PetscMalloc(n*sizeof(PetscScalar),&wr);
82: PetscMalloc(n*sizeof(PetscScalar),&wi);
83: DSSetEigenvalueComparison(ds,SlepcCompareLargestMagnitude,NULL);
84: DSSolve(ds,wr,wi);
85: DSSort(ds,wr,wi,NULL,NULL,NULL);
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",n);
93: for (i=0;i<n;i++) {
94: #if defined(PETSC_USE_COMPLEX)
95: re = PetscRealPart(wr[i]);
96: im = PetscImaginaryPart(wr[i]);
97: #else
98: re = wr[i];
99: im = wi[i];
100: #endif
101: if (PetscAbs(im)<1e-10) {
102: PetscViewerASCIIPrintf(viewer," %.5F\n",re);
103: } else {
104: PetscViewerASCIIPrintf(viewer," %.5F%+.5Fi\n",re,im);
105: }
106: }
108: PetscFree(wr);
109: PetscFree(wi);
110: DSDestroy(&ds);
111: SlepcFinalize();
112: return 0;
113: }