Line | Branch | Exec | Source |
---|---|---|---|
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.\n\n"; | ||
12 | |||
13 | #include <slepcds.h> | ||
14 | |||
15 | 60 | int main(int argc,char **argv) | |
16 | { | ||
17 | 60 | DS ds; | |
18 | 60 | SlepcSC sc; | |
19 | 60 | PetscReal re,im; | |
20 | 60 | PetscScalar *A,*B,*Q,*eigr,*eigi,d; | |
21 | 60 | PetscInt i,j,n=10,ld; | |
22 | 60 | PetscViewer viewer; | |
23 | 60 | PetscBool verbose,extrarow; | |
24 | |||
25 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
60 | PetscFunctionBeginUser; |
26 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
27 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
28 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GHIEP - dimension %" PetscInt_FMT ".\n",n)); |
29 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose)); |
30 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow)); |
31 | |||
32 | /* Create DS object */ | ||
33 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSCreate(PETSC_COMM_WORLD,&ds)); |
34 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSSetType(ds,DSGHIEP)); |
35 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSSetFromOptions(ds)); |
36 | 60 | ld = n+2; /* test leading dimension larger than n */ | |
37 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSAllocate(ds,ld)); |
38 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSSetDimensions(ds,n,0,0)); |
39 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSSetExtraRow(ds,extrarow)); |
40 | |||
41 | /* Set up viewer */ | ||
42 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer)); |
43 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL)); |
44 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSView(ds,viewer)); |
45 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(PetscViewerPopFormat(viewer)); |
46 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
60 | if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB)); |
47 | |||
48 | /* Fill with a symmetric Toeplitz matrix */ | ||
49 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSGetArray(ds,DS_MAT_A,&A)); |
50 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSGetArray(ds,DS_MAT_B,&B)); |
51 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
660 | for (i=0;i<n;i++) A[i+i*ld]=2.0; |
52 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
180 | for (j=1;j<3;j++) { |
53 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1140 | for (i=0;i<n-j;i++) { A[i+(i+j)*ld]=1.0; A[(i+j)+i*ld]=1.0; } |
54 | } | ||
55 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
180 | for (j=1;j<3;j++) { A[0+j*ld]=-1.0*(j+2); A[j+0*ld]=-1.0*(j+2); } |
56 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
60 | if (extrarow) A[n+(n-1)*ld]=-1.0; |
57 | /* Signature matrix */ | ||
58 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
660 | for (i=0;i<n;i++) B[i+i*ld]=1.0; |
59 | 60 | B[0] = -1.0; | |
60 | 60 | B[n-1+(n-1)*ld] = -1.0; | |
61 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSRestoreArray(ds,DS_MAT_A,&A)); |
62 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSRestoreArray(ds,DS_MAT_B,&B)); |
63 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSSetState(ds,DS_STATE_RAW)); |
64 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
60 | if (verbose) { |
65 | ✗ | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n")); | |
66 | ✗ | PetscCall(DSView(ds,viewer)); | |
67 | } | ||
68 | |||
69 | /* Solve */ | ||
70 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(PetscCalloc2(n,&eigr,n,&eigi)); |
71 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSGetSlepcSC(ds,&sc)); |
72 | 60 | sc->comparison = SlepcCompareLargestMagnitude; | |
73 | 60 | sc->comparisonctx = NULL; | |
74 | 60 | sc->map = NULL; | |
75 | 60 | sc->mapobj = NULL; | |
76 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSSolve(ds,eigr,eigi)); |
77 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSSort(ds,eigr,eigi,NULL,NULL,NULL)); |
78 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
60 | if (extrarow) PetscCall(DSUpdateExtraRow(ds)); |
79 | |||
80 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
60 | if (verbose) { |
81 | ✗ | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n")); | |
82 | ✗ | PetscCall(DSView(ds,viewer)); | |
83 | } | ||
84 | |||
85 | /* Print eigenvalues */ | ||
86 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n")); |
87 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
660 | for (i=0;i<n;i++) { |
88 | #if defined(PETSC_USE_COMPLEX) | ||
89 | 300 | re = PetscRealPart(eigr[i]); | |
90 | 300 | im = PetscImaginaryPart(eigr[i]); | |
91 | #else | ||
92 | 300 | re = eigr[i]; | |
93 | 300 | im = eigi[i]; | |
94 | #endif | ||
95 |
8/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
600 | if (PetscAbs(im)<1e-10) PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re)); |
96 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
600 | else PetscCall(PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im)); |
97 | } | ||
98 | |||
99 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
60 | if (extrarow) { |
100 | /* Check that extra row is correct */ | ||
101 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(DSGetArray(ds,DS_MAT_A,&A)); |
102 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(DSGetArray(ds,DS_MAT_Q,&Q)); |
103 | d = 0.0; | ||
104 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
330 | for (i=0;i<n;i++) d += A[n+i*ld]+Q[n-1+i*ld]; |
105 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
30 | if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in the extra row of %g\n",(double)PetscAbsScalar(d))); |
106 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(DSRestoreArray(ds,DS_MAT_A,&A)); |
107 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(DSRestoreArray(ds,DS_MAT_Q,&Q)); |
108 | } | ||
109 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(PetscFree2(eigr,eigi)); |
110 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(DSDestroy(&ds)); |
111 |
3/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
60 | PetscCall(SlepcFinalize()); |
112 | return 0; | ||
113 | } | ||
114 | |||
115 | /*TEST | ||
116 | |||
117 | testset: | ||
118 | args: -ds_method {{0 1 2}} | ||
119 | filter: grep -v "solving the problem" | sed -e "s/extrarow//" | ||
120 | output_file: output/test5_1.out | ||
121 | requires: !single | ||
122 | test: | ||
123 | suffix: 1 | ||
124 | test: | ||
125 | suffix: 2 | ||
126 | args: -extrarow | ||
127 | |||
128 | TEST*/ | ||
129 |