Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
11 : static char help[] = "Test DSGHIEP with compact storage.\n\n";
12 :
13 : #include <slepcds.h>
14 :
15 6 : int main(int argc,char **argv)
16 : {
17 6 : DS ds;
18 6 : SlepcSC sc;
19 6 : PetscReal *T,*s,re,im;
20 6 : PetscScalar *eigr,*eigi;
21 6 : PetscInt i,n=10,l=2,k=5,ld;
22 6 : PetscViewer viewer;
23 6 : PetscBool verbose;
24 :
25 6 : PetscFunctionBeginUser;
26 6 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
27 6 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
28 6 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GHIEP with compact storage - dimension %" PetscInt_FMT ".\n",n));
29 6 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL));
30 6 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
31 6 : PetscCheck(l<=n && k<=n && l<=k,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Wrong value of dimensions");
32 6 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
33 :
34 : /* Create DS object */
35 6 : PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
36 6 : PetscCall(DSSetType(ds,DSGHIEP));
37 6 : PetscCall(DSSetFromOptions(ds));
38 6 : ld = n+2; /* test leading dimension larger than n */
39 6 : PetscCall(DSAllocate(ds,ld));
40 6 : PetscCall(DSSetDimensions(ds,n,l,k));
41 6 : PetscCall(DSSetCompact(ds,PETSC_TRUE));
42 :
43 : /* Set up viewer */
44 6 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
45 6 : PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
46 6 : PetscCall(DSView(ds,viewer));
47 6 : PetscCall(PetscViewerPopFormat(viewer));
48 6 : PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
49 :
50 : /* Fill arrow-tridiagonal matrix */
51 6 : PetscCall(DSGetArrayReal(ds,DS_MAT_T,&T));
52 6 : PetscCall(DSGetArrayReal(ds,DS_MAT_D,&s));
53 51 : for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
54 18 : for (i=k;i<n-1;i++) T[i+ld] = 1.0;
55 15 : for (i=l;i<k;i++) T[i+2*ld] = 1.0;
56 6 : T[2*ld+l+1] = -7; T[ld+k+1] = -7;
57 : /* Signature matrix */
58 51 : for (i=0;i<n;i++) s[i] = 1.0;
59 6 : s[l+1] = -1.0;
60 6 : s[k+1] = -1.0;
61 6 : PetscCall(DSRestoreArrayReal(ds,DS_MAT_T,&T));
62 6 : PetscCall(DSRestoreArrayReal(ds,DS_MAT_D,&s));
63 6 : if (l==0 && k==0) PetscCall(DSSetState(ds,DS_STATE_INTERMEDIATE));
64 6 : else PetscCall(DSSetState(ds,DS_STATE_RAW));
65 6 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
66 6 : PetscCall(DSView(ds,viewer));
67 :
68 : /* Solve */
69 6 : PetscCall(PetscCalloc2(n,&eigr,n,&eigi));
70 6 : PetscCall(DSGetSlepcSC(ds,&sc));
71 6 : sc->comparison = SlepcCompareLargestMagnitude;
72 6 : sc->comparisonctx = NULL;
73 6 : sc->map = NULL;
74 6 : sc->mapobj = NULL;
75 6 : PetscCall(DSSolve(ds,eigr,eigi));
76 6 : PetscCall(DSSort(ds,eigr,eigi,NULL,NULL,NULL));
77 6 : if (verbose) {
78 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
79 0 : PetscCall(DSView(ds,viewer));
80 : }
81 :
82 : /* Print eigenvalues */
83 6 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n"));
84 51 : for (i=0;i<n;i++) {
85 : #if defined(PETSC_USE_COMPLEX)
86 45 : re = PetscRealPart(eigr[i]);
87 45 : im = PetscImaginaryPart(eigr[i]);
88 : #else
89 : re = eigr[i];
90 : im = eigi[i];
91 : #endif
92 45 : if (PetscAbs(im)<1e-10) PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re));
93 45 : else PetscCall(PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im));
94 : }
95 6 : PetscCall(PetscFree2(eigr,eigi));
96 6 : PetscCall(DSDestroy(&ds));
97 6 : PetscCall(SlepcFinalize());
98 : return 0;
99 : }
100 :
101 : /*TEST
102 :
103 : test:
104 : suffix: 1
105 : requires: !single
106 : args: -ds_method {{0 1 2}}
107 : filter: grep -v "solving the problem"
108 :
109 : test:
110 : suffix: 2
111 : args: -n 5 -k 4 -l 4 -ds_method {{0 1 2}}
112 : filter: grep -v "solving the problem"
113 :
114 : TEST*/
|