Actual source code: test3.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Test DSHEP with compact storage.\n\n";
13: #include <slepcds.h>
15: int main(int argc,char **argv)
16: {
18: DS ds;
19: SlepcSC sc;
20: PetscReal *T;
21: PetscScalar *eig;
22: PetscInt i,n=10,l=2,k=5,ld;
23: PetscViewer viewer;
24: PetscBool verbose,extrarow;
26: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
27: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
28: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type HEP with compact storage - dimension %D.\n",n);
29: PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
30: PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
31: if (l>n || k>n || l>k) SETERRQ(PETSC_COMM_WORLD,1,"Wrong value of dimensions");
32: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
33: PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow);
35: /* Create DS object */
36: DSCreate(PETSC_COMM_WORLD,&ds);
37: DSSetType(ds,DSHEP);
38: DSSetFromOptions(ds);
39: ld = n+2; /* test leading dimension larger than n */
40: DSAllocate(ds,ld);
41: DSSetDimensions(ds,n,0,l,k);
42: DSSetCompact(ds,PETSC_TRUE);
43: DSSetExtraRow(ds,extrarow);
45: /* Set up viewer */
46: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
47: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
48: DSView(ds,viewer);
49: PetscViewerPopFormat(viewer);
50: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
52: /* Fill arrow-tridiagonal matrix */
53: DSGetArrayReal(ds,DS_MAT_T,&T);
54: for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
55: for (i=l;i<n-1;i++) T[i+ld] = 1.0;
56: if (extrarow) T[n-1+ld] = 1.0;
57: DSRestoreArrayReal(ds,DS_MAT_T,&T);
58: if (l==0 && k==0) {
59: DSSetState(ds,DS_STATE_INTERMEDIATE);
60: } else {
61: DSSetState(ds,DS_STATE_RAW);
62: }
63: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
64: DSView(ds,viewer);
66: /* Solve */
67: PetscMalloc1(n,&eig);
68: DSGetSlepcSC(ds,&sc);
69: sc->comparison = SlepcCompareLargestMagnitude;
70: sc->comparisonctx = NULL;
71: sc->map = NULL;
72: sc->mapobj = NULL;
73: DSSolve(ds,eig,NULL);
74: DSSort(ds,eig,NULL,NULL,NULL,NULL);
75: if (extrarow) { DSUpdateExtraRow(ds); }
76: if (verbose) {
77: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
78: DSView(ds,viewer);
79: }
81: /* Print eigenvalues */
82: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
83: for (i=0;i<n;i++) {
84: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)PetscRealPart(eig[i]));
85: }
87: PetscFree(eig);
88: DSDestroy(&ds);
89: SlepcFinalize();
90: return ierr;
91: }
93: /*TEST
95: test:
96: suffix: 1
97: requires: !single
98: args: -n 9 -ds_method {{0 1 2}}
99: filter: grep -v "solving the problem"
101: TEST*/