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 | This example implements one of the problems found at | ||
12 | NLEVP: A Collection of Nonlinear Eigenvalue Problems, | ||
13 | The University of Manchester. | ||
14 | The details of the collection can be found at: | ||
15 | [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue | ||
16 | Problems", ACM Trans. Math. Software 39(2), Article 7, 2013. | ||
17 | |||
18 | The planar_waveguide problem is a quartic PEP with symmetric matrices, | ||
19 | arising from a finite element solution of the propagation constants in a | ||
20 | six-layer planar waveguide. | ||
21 | */ | ||
22 | |||
23 | static char help[] = "FEM solution of the propagation constants in a six-layer planar waveguide.\n\n" | ||
24 | "The command line options are:\n" | ||
25 | " -n <n>, the dimension of the matrices.\n\n"; | ||
26 | |||
27 | #include <slepcpep.h> | ||
28 | |||
29 | #define NMAT 5 | ||
30 | #define NL 6 | ||
31 | |||
32 | 20 | int main(int argc,char **argv) | |
33 | { | ||
34 | 20 | Mat A[NMAT]; /* problem matrices */ | |
35 | 20 | PEP pep; /* polynomial eigenproblem solver context */ | |
36 | 20 | PetscInt n=128,nlocal,k,Istart,Iend,i,j,start_ct,end_ct; | |
37 | 20 | PetscReal w=9.92918,a=0.0,b=2.0,h,deltasq; | |
38 | 20 | PetscReal nref[NL],K2[NL],q[NL],*md,*supd,*subd; | |
39 | 20 | PetscScalar v,alpha; | |
40 | 20 | PetscBool terse; | |
41 | |||
42 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBeginUser; |
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.
|
20 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
44 | |||
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(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
46 | 20 | n = (n/4)*4; | |
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(PetscPrintf(PETSC_COMM_WORLD,"\nPlanar waveguide, n=%" PetscInt_FMT "\n\n",n+1)); |
48 | 20 | h = (b-a)/n; | |
49 | 20 | nlocal = (n/4)-1; | |
50 | |||
51 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
52 | Set waveguide parameters used in construction of matrices | ||
53 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
54 | |||
55 | /* refractive indices in each layer */ | ||
56 | 20 | nref[0] = 1.5; | |
57 | 20 | nref[1] = 1.66; | |
58 | 20 | nref[2] = 1.6; | |
59 | 20 | nref[3] = 1.53; | |
60 | 20 | nref[4] = 1.66; | |
61 | 20 | nref[5] = 1.0; | |
62 | |||
63 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
140 | for (i=0;i<NL;i++) K2[i] = w*w*nref[i]*nref[i]; |
64 | 20 | deltasq = K2[0] - K2[NL-1]; | |
65 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
140 | for (i=0;i<NL;i++) q[i] = K2[i] - (K2[0] + K2[NL-1])/2; |
66 | |||
67 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
68 | Compute the polynomial matrices | ||
69 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
70 | |||
71 | /* initialize matrices */ | ||
72 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
120 | for (i=0;i<NMAT;i++) { |
73 |
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.
|
100 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A[i])); |
74 |
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.
|
100 | PetscCall(MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n+1,n+1)); |
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.
|
100 | PetscCall(MatSetFromOptions(A[i])); |
76 | } | ||
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(MatGetOwnershipRange(A[0],&Istart,&Iend)); |
78 | |||
79 | /* A0 */ | ||
80 | 20 | alpha = (h/6)*(deltasq*deltasq/16); | |
81 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2600 | for (i=Istart;i<Iend;i++) { |
82 | 2580 | v = 4.0; | |
83 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
2580 | if (i==0 || i==n) v = 2.0; |
84 |
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.
|
2580 | PetscCall(MatSetValue(A[0],i,i,v*alpha,INSERT_VALUES)); |
85 |
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.
|
2580 | if (i>0) PetscCall(MatSetValue(A[0],i,i-1,alpha,INSERT_VALUES)); |
86 |
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.
|
2580 | if (i<n) PetscCall(MatSetValue(A[0],i,i+1,alpha,INSERT_VALUES)); |
87 | } | ||
88 | |||
89 | /* A1 */ | ||
90 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ 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 (Istart==0) PetscCall(MatSetValue(A[1],0,0,-deltasq/4,INSERT_VALUES)); |
91 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ 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 (Iend==n+1) PetscCall(MatSetValue(A[1],n,n,deltasq/4,INSERT_VALUES)); |
92 | |||
93 | /* A2 */ | ||
94 | 20 | alpha = 1.0/h; | |
95 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2600 | for (i=Istart;i<Iend;i++) { |
96 | 2580 | v = 2.0; | |
97 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
2580 | if (i==0 || i==n) v = 1.0; |
98 |
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.
|
2580 | PetscCall(MatSetValue(A[2],i,i,v*alpha,ADD_VALUES)); |
99 |
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.
|
2580 | if (i>0) PetscCall(MatSetValue(A[2],i,i-1,-alpha,ADD_VALUES)); |
100 |
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.
|
2580 | if (i<n) PetscCall(MatSetValue(A[2],i,i+1,-alpha,ADD_VALUES)); |
101 | } | ||
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.
|
20 | PetscCall(PetscMalloc3(n+1,&md,n+1,&supd,n+1,&subd)); |
103 | |||
104 | 20 | md[0] = 2.0*q[1]; | |
105 | 20 | supd[1] = q[1]; | |
106 | 20 | subd[0] = q[1]; | |
107 | |||
108 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
100 | for (k=1;k<=NL-2;k++) { |
109 | |||
110 | 80 | end_ct = k*(nlocal+1); | |
111 | 80 | start_ct = end_ct-nlocal; | |
112 | |||
113 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2560 | for (j=start_ct;j<end_ct;j++) { |
114 | 2480 | md[j] = 4*q[k]; | |
115 | 2480 | supd[j+1] = q[k]; | |
116 | 2480 | subd[j] = q[k]; | |
117 | } | ||
118 | |||
119 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
80 | if (k < 4) { /* interface points */ |
120 | 60 | md[end_ct] = 4*(q[k] + q[k+1])/2.0; | |
121 | 60 | supd[end_ct+1] = q[k+1]; | |
122 | 60 | subd[end_ct] = q[k+1]; | |
123 | } | ||
124 | |||
125 | } | ||
126 | |||
127 | 20 | md[n] = 2*q[NL-2]; | |
128 | 20 | supd[n] = q[NL-2]; | |
129 | 20 | subd[n] = q[NL-2]; | |
130 | |||
131 | 20 | alpha = -h/6.0; | |
132 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2600 | for (i=Istart;i<Iend;i++) { |
133 |
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.
|
2580 | PetscCall(MatSetValue(A[2],i,i,md[i]*alpha,ADD_VALUES)); |
134 |
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.
|
2580 | if (i>0) PetscCall(MatSetValue(A[2],i,i-1,subd[i-1]*alpha,ADD_VALUES)); |
135 |
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.
|
2580 | if (i<n) PetscCall(MatSetValue(A[2],i,i+1,supd[i+1]*alpha,ADD_VALUES)); |
136 | } | ||
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(PetscFree3(md,supd,subd)); |
138 | |||
139 | /* A3 */ | ||
140 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ 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 (Istart==0) PetscCall(MatSetValue(A[3],0,0,1.0,INSERT_VALUES)); |
141 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ 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 (Iend==n+1) PetscCall(MatSetValue(A[3],n,n,1.0,INSERT_VALUES)); |
142 | |||
143 | /* A4 */ | ||
144 | 20 | alpha = (h/6); | |
145 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2600 | for (i=Istart;i<Iend;i++) { |
146 | 2580 | v = 4.0; | |
147 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
2580 | if (i==0 || i==n) v = 2.0; |
148 |
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.
|
2580 | PetscCall(MatSetValue(A[4],i,i,v*alpha,INSERT_VALUES)); |
149 |
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.
|
2580 | if (i>0) PetscCall(MatSetValue(A[4],i,i-1,alpha,INSERT_VALUES)); |
150 |
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.
|
2580 | if (i<n) PetscCall(MatSetValue(A[4],i,i+1,alpha,INSERT_VALUES)); |
151 | } | ||
152 | |||
153 | /* assemble matrices */ | ||
154 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
120 | for (i=0;i<NMAT;i++) PetscCall(MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY)); |
155 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
120 | for (i=0;i<NMAT;i++) PetscCall(MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY)); |
156 | |||
157 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
158 | Create the eigensolver and solve the problem | ||
159 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
160 | |||
161 |
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(PEPCreate(PETSC_COMM_WORLD,&pep)); |
162 |
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(PEPSetOperators(pep,NMAT,A)); |
163 |
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(PEPSetFromOptions(pep)); |
164 |
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(PEPSolve(pep)); |
165 | |||
166 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
167 | Display solution and clean up | ||
168 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
169 | |||
170 | /* show detailed info unless -terse option is given by user */ | ||
171 |
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,"-terse",&terse)); |
172 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ 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 (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL)); |
173 | else { | ||
174 | ✗ | PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); | |
175 | ✗ | PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD)); | |
176 | ✗ | PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD)); | |
177 | ✗ | PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); | |
178 | } | ||
179 |
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(PEPDestroy(&pep)); |
180 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
120 | for (i=0;i<NMAT;i++) PetscCall(MatDestroy(&A[i])); |
181 |
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()); |
182 | return 0; | ||
183 | } | ||
184 | |||
185 | /*TEST | ||
186 | |||
187 | test: | ||
188 | suffix: 1 | ||
189 | args: -pep_type {{toar linear}} -pep_nev 4 -st_type sinvert -terse | ||
190 | requires: !single | ||
191 | |||
192 | TEST*/ | ||
193 |