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 DSNHEP.\n\n"; | ||
12 | |||
13 | #include <slepcds.h> | ||
14 | |||
15 | 20 | int main(int argc,char **argv) | |
16 | { | ||
17 | 20 | DS ds; | |
18 | 20 | SlepcSC sc; | |
19 | 20 | DSType type; | |
20 | 20 | DSStateType state; | |
21 | 20 | PetscScalar *A,*X,*Q,*wr,*wi,d; | |
22 | 20 | PetscReal re,im,rnorm,aux; | |
23 | 20 | PetscInt i,j,n=10,ld,method; | |
24 | 20 | PetscViewer viewer; | |
25 | 20 | PetscBool verbose,extrarow; | |
26 | |||
27 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBeginUser; |
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.
|
20 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
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.
|
20 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
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.
|
20 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NHEP - dimension %" PetscInt_FMT ".\n",n)); |
31 |
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.
|
20 | PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose)); |
32 |
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.
|
20 | PetscCall(PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow)); |
33 | |||
34 | /* Create DS object */ | ||
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.
|
20 | PetscCall(DSCreate(PETSC_COMM_WORLD,&ds)); |
36 |
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.
|
20 | PetscCall(DSSetType(ds,DSNHEP)); |
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.
|
20 | PetscCall(DSSetFromOptions(ds)); |
38 | 20 | ld = n+2; /* test leading dimension larger than n */ | |
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.
|
20 | PetscCall(DSAllocate(ds,ld)); |
40 |
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.
|
20 | PetscCall(DSSetDimensions(ds,n,0,0)); |
41 |
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.
|
20 | PetscCall(DSSetExtraRow(ds,extrarow)); |
42 | |||
43 | /* Set up viewer */ | ||
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.
|
20 | PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&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.
|
20 | PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL)); |
46 |
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.
|
20 | PetscCall(DSView(ds,viewer)); |
47 |
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.
|
20 | PetscCall(PetscViewerPopFormat(viewer)); |
48 |
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.
|
20 | if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB)); |
49 | |||
50 | /* Fill with Grcar matrix */ | ||
51 |
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.
|
20 | PetscCall(DSGetArray(ds,DS_MAT_A,&A)); |
52 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
200 | for (i=1;i<n;i++) A[i+(i-1)*ld]=-1.0; |
53 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
100 | for (j=0;j<4;j++) { |
54 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
760 | for (i=0;i<n-j;i++) A[i+(i+j)*ld]=1.0; |
55 | } | ||
56 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
20 | if (extrarow) A[n+(n-1)*ld]=-1.0; |
57 |
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.
|
20 | PetscCall(DSRestoreArray(ds,DS_MAT_A,&A)); |
58 |
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.
|
20 | PetscCall(DSSetState(ds,DS_STATE_INTERMEDIATE)); |
59 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
20 | if (verbose) { |
60 | ✗ | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n")); | |
61 | ✗ | PetscCall(DSView(ds,viewer)); | |
62 | } | ||
63 | |||
64 | /* Solve */ | ||
65 |
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.
|
20 | PetscCall(PetscMalloc2(n,&wr,n,&wi)); |
66 |
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.
|
20 | PetscCall(DSGetSlepcSC(ds,&sc)); |
67 | 20 | sc->comparison = SlepcCompareLargestMagnitude; | |
68 | 20 | sc->comparisonctx = NULL; | |
69 | 20 | sc->map = NULL; | |
70 | 20 | sc->mapobj = NULL; | |
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.
|
20 | PetscCall(DSSolve(ds,wr,wi)); |
72 |
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.
|
20 | PetscCall(DSSort(ds,wr,wi,NULL,NULL,NULL)); |
73 |
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.
|
20 | if (extrarow) PetscCall(DSUpdateExtraRow(ds)); |
74 | |||
75 |
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.
|
20 | PetscCall(DSGetType(ds,&type)); |
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.
|
20 | PetscCall(DSGetMethod(ds,&method)); |
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.
|
20 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"DS of type %s, method used=%" PetscInt_FMT "\n",type,method)); |
78 |
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.
|
20 | PetscCall(DSGetState(ds,&state)); |
79 |
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.
|
20 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"State after solve: %s\n",DSStateTypes[state])); |
80 | |||
81 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
20 | if (verbose) { |
82 | ✗ | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n")); | |
83 | ✗ | PetscCall(DSView(ds,viewer)); | |
84 | } | ||
85 | |||
86 | /* Print eigenvalues */ | ||
87 |
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.
|
20 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n")); |
88 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
220 | for (i=0;i<n;i++) { |
89 | #if defined(PETSC_USE_COMPLEX) | ||
90 | 100 | re = PetscRealPart(wr[i]); | |
91 | 100 | im = PetscImaginaryPart(wr[i]); | |
92 | #else | ||
93 | 100 | re = wr[i]; | |
94 | 100 | im = wi[i]; | |
95 | #endif | ||
96 |
3/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
200 | if (PetscAbs(im)<1e-10) PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re)); |
97 |
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.
|
200 | else PetscCall(PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im)); |
98 | } | ||
99 | |||
100 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
20 | if (extrarow) { |
101 | /* Check that extra row is correct */ | ||
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.
|
10 | PetscCall(DSGetArray(ds,DS_MAT_A,&A)); |
103 |
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.
|
10 | PetscCall(DSGetArray(ds,DS_MAT_Q,&Q)); |
104 | d = 0.0; | ||
105 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
110 | for (i=0;i<n;i++) d += A[n+i*ld]+Q[n-1+i*ld]; |
106 |
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.
|
10 | 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))); |
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.
|
10 | PetscCall(DSRestoreArray(ds,DS_MAT_A,&A)); |
108 |
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.
|
10 | PetscCall(DSRestoreArray(ds,DS_MAT_Q,&Q)); |
109 | } | ||
110 | |||
111 | /* Eigenvectors */ | ||
112 | 20 | j = 2; | |
113 |
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.
|
20 | PetscCall(DSVectors(ds,DS_MAT_X,&j,&rnorm)); /* third eigenvector */ |
114 |
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.
|
20 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Value of rnorm for 3rd vector = %.3f\n",(double)rnorm)); |
115 |
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.
|
20 | PetscCall(DSVectors(ds,DS_MAT_X,NULL,NULL)); /* all eigenvectors */ |
116 | 20 | j = 0; | |
117 | 20 | rnorm = 0.0; | |
118 |
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.
|
20 | PetscCall(DSGetArray(ds,DS_MAT_X,&X)); |
119 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
220 | for (i=0;i<n;i++) { |
120 | #if defined(PETSC_USE_COMPLEX) | ||
121 | 100 | aux = PetscAbsScalar(X[i+j*ld]); | |
122 | #else | ||
123 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
|
100 | if (PetscAbs(wi[j])==0.0) aux = PetscAbsScalar(X[i+j*ld]); |
124 | 100 | else aux = SlepcAbsEigenvalue(X[i+j*ld],X[i+(j+1)*ld]); | |
125 | #endif | ||
126 | 200 | rnorm += aux*aux; | |
127 | } | ||
128 |
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.
|
20 | PetscCall(DSRestoreArray(ds,DS_MAT_X,&X)); |
129 | 20 | rnorm = PetscSqrtReal(rnorm); | |
130 |
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.
|
20 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st vector = %.3f\n",(double)rnorm)); |
131 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
20 | if (verbose) { |
132 | ✗ | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n")); | |
133 | ✗ | PetscCall(DSView(ds,viewer)); | |
134 | } | ||
135 | |||
136 |
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.
|
20 | PetscCall(PetscFree2(wr,wi)); |
137 |
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.
|
20 | PetscCall(DSDestroy(&ds)); |
138 |
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.
|
20 | PetscCall(SlepcFinalize()); |
139 | return 0; | ||
140 | } | ||
141 | |||
142 | /*TEST | ||
143 | |||
144 | testset: | ||
145 | filter: sed -e "s/[+-]\([0-9]\.[0-9]*i\)/+-\\1/" | sed -e "s/extrarow//" | ||
146 | output_file: output/test1_1.out | ||
147 | requires: !single | ||
148 | test: | ||
149 | suffix: 1 | ||
150 | test: | ||
151 | suffix: 2 | ||
152 | args: -extrarow | ||
153 | |||
154 | TEST*/ | ||
155 |