Actual source code: test1.c
slepc-3.6.1 2015-09-03
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 DSNHEP.\n\n";
24: #include <slepcds.h>
28: int main(int argc,char **argv)
29: {
31: DS ds;
32: SlepcSC sc;
33: PetscScalar *A,*wr,*wi;
34: PetscReal re,im;
35: PetscInt i,j,n=10,ld;
36: PetscViewer viewer;
37: PetscBool verbose,extrarow;
39: SlepcInitialize(&argc,&argv,(char*)0,help);
40: PetscOptionsGetInt(NULL,"-n",&n,NULL);
41: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NHEP - dimension %D.\n",n);
42: PetscOptionsHasName(NULL,"-verbose",&verbose);
43: PetscOptionsHasName(NULL,"-extrarow",&extrarow);
45: /* Create DS object */
46: DSCreate(PETSC_COMM_WORLD,&ds);
47: DSSetType(ds,DSNHEP);
48: DSSetFromOptions(ds);
49: ld = n+2; /* test leading dimension larger than n */
50: DSAllocate(ds,ld);
51: DSSetDimensions(ds,n,0,0,0);
52: DSSetExtraRow(ds,extrarow);
54: /* Set up viewer */
55: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
56: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
57: DSView(ds,viewer);
58: PetscViewerPopFormat(viewer);
59: if (verbose) {
60: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
61: }
63: /* Fill with Grcar matrix */
64: DSGetArray(ds,DS_MAT_A,&A);
65: for (i=1;i<n;i++) A[i+(i-1)*ld]=-1.0;
66: for (j=0;j<4;j++) {
67: for (i=0;i<n-j;i++) A[i+(i+j)*ld]=1.0;
68: }
69: if (extrarow) A[n+(n-1)*ld]=-1.0;
70: DSRestoreArray(ds,DS_MAT_A,&A);
71: DSSetState(ds,DS_STATE_INTERMEDIATE);
72: if (verbose) {
73: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
74: DSView(ds,viewer);
75: }
77: /* Solve */
78: PetscMalloc2(n,&wr,n,&wi);
79: DSGetSlepcSC(ds,&sc);
80: sc->comparison = SlepcCompareLargestMagnitude;
81: sc->comparisonctx = NULL;
82: sc->map = NULL;
83: sc->mapobj = NULL;
84: DSSolve(ds,wr,wi);
85: DSSort(ds,wr,wi,NULL,NULL,NULL);
86: if (extrarow) { DSUpdateExtraRow(ds); }
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: #if defined(PETSC_USE_COMPLEX)
96: re = PetscRealPart(wr[i]);
97: im = PetscImaginaryPart(wr[i]);
98: #else
99: re = wr[i];
100: im = wi[i];
101: #endif
102: if (PetscAbs(im)<1e-10) {
103: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re);
104: } else {
105: PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im);
106: }
107: }
109: PetscFree2(wr,wi);
110: DSDestroy(&ds);
111: SlepcFinalize();
112: return 0;
113: }