Actual source code: planar_waveguide.c
slepc-main 2024-12-17
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: */
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.
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: */
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";
27: #include <slepcpep.h>
29: #define NMAT 5
30: #define NL 6
32: int main(int argc,char **argv)
33: {
34: Mat A[NMAT]; /* problem matrices */
35: PEP pep; /* polynomial eigenproblem solver context */
36: PetscInt n=128,nlocal,k,Istart,Iend,i,j,start_ct,end_ct;
37: PetscReal w=9.92918,a=0.0,b=2.0,h,deltasq;
38: PetscReal nref[NL],K2[NL],q[NL],*md,*supd,*subd;
39: PetscScalar v,alpha;
40: PetscBool terse;
42: PetscFunctionBeginUser;
43: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
45: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
46: n = (n/4)*4;
47: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nPlanar waveguide, n=%" PetscInt_FMT "\n\n",n+1));
48: h = (b-a)/n;
49: nlocal = (n/4)-1;
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Set waveguide parameters used in construction of matrices
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: /* refractive indices in each layer */
56: nref[0] = 1.5;
57: nref[1] = 1.66;
58: nref[2] = 1.6;
59: nref[3] = 1.53;
60: nref[4] = 1.66;
61: nref[5] = 1.0;
63: for (i=0;i<NL;i++) K2[i] = w*w*nref[i]*nref[i];
64: deltasq = K2[0] - K2[NL-1];
65: for (i=0;i<NL;i++) q[i] = K2[i] - (K2[0] + K2[NL-1])/2;
67: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68: Compute the polynomial matrices
69: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71: /* initialize matrices */
72: for (i=0;i<NMAT;i++) {
73: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[i]));
74: PetscCall(MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n+1,n+1));
75: PetscCall(MatSetFromOptions(A[i]));
76: }
77: PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
79: /* A0 */
80: alpha = (h/6)*(deltasq*deltasq/16);
81: for (i=Istart;i<Iend;i++) {
82: v = 4.0;
83: if (i==0 || i==n) v = 2.0;
84: PetscCall(MatSetValue(A[0],i,i,v*alpha,INSERT_VALUES));
85: if (i>0) PetscCall(MatSetValue(A[0],i,i-1,alpha,INSERT_VALUES));
86: if (i<n) PetscCall(MatSetValue(A[0],i,i+1,alpha,INSERT_VALUES));
87: }
89: /* A1 */
90: if (Istart==0) PetscCall(MatSetValue(A[1],0,0,-deltasq/4,INSERT_VALUES));
91: if (Iend==n+1) PetscCall(MatSetValue(A[1],n,n,deltasq/4,INSERT_VALUES));
93: /* A2 */
94: alpha = 1.0/h;
95: for (i=Istart;i<Iend;i++) {
96: v = 2.0;
97: if (i==0 || i==n) v = 1.0;
98: PetscCall(MatSetValue(A[2],i,i,v*alpha,ADD_VALUES));
99: if (i>0) PetscCall(MatSetValue(A[2],i,i-1,-alpha,ADD_VALUES));
100: if (i<n) PetscCall(MatSetValue(A[2],i,i+1,-alpha,ADD_VALUES));
101: }
102: PetscCall(PetscMalloc3(n+1,&md,n+1,&supd,n+1,&subd));
104: md[0] = 2.0*q[1];
105: supd[1] = q[1];
106: subd[0] = q[1];
108: for (k=1;k<=NL-2;k++) {
110: end_ct = k*(nlocal+1);
111: start_ct = end_ct-nlocal;
113: for (j=start_ct;j<end_ct;j++) {
114: md[j] = 4*q[k];
115: supd[j+1] = q[k];
116: subd[j] = q[k];
117: }
119: if (k < 4) { /* interface points */
120: md[end_ct] = 4*(q[k] + q[k+1])/2.0;
121: supd[end_ct+1] = q[k+1];
122: subd[end_ct] = q[k+1];
123: }
125: }
127: md[n] = 2*q[NL-2];
128: supd[n] = q[NL-2];
129: subd[n] = q[NL-2];
131: alpha = -h/6.0;
132: for (i=Istart;i<Iend;i++) {
133: PetscCall(MatSetValue(A[2],i,i,md[i]*alpha,ADD_VALUES));
134: if (i>0) PetscCall(MatSetValue(A[2],i,i-1,subd[i-1]*alpha,ADD_VALUES));
135: if (i<n) PetscCall(MatSetValue(A[2],i,i+1,supd[i+1]*alpha,ADD_VALUES));
136: }
137: PetscCall(PetscFree3(md,supd,subd));
139: /* A3 */
140: if (Istart==0) PetscCall(MatSetValue(A[3],0,0,1.0,INSERT_VALUES));
141: if (Iend==n+1) PetscCall(MatSetValue(A[3],n,n,1.0,INSERT_VALUES));
143: /* A4 */
144: alpha = (h/6);
145: for (i=Istart;i<Iend;i++) {
146: v = 4.0;
147: if (i==0 || i==n) v = 2.0;
148: PetscCall(MatSetValue(A[4],i,i,v*alpha,INSERT_VALUES));
149: if (i>0) PetscCall(MatSetValue(A[4],i,i-1,alpha,INSERT_VALUES));
150: if (i<n) PetscCall(MatSetValue(A[4],i,i+1,alpha,INSERT_VALUES));
151: }
153: /* assemble matrices */
154: for (i=0;i<NMAT;i++) PetscCall(MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY));
155: for (i=0;i<NMAT;i++) PetscCall(MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY));
157: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158: Create the eigensolver and solve the problem
159: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
162: PetscCall(PEPSetOperators(pep,NMAT,A));
163: PetscCall(PEPSetFromOptions(pep));
164: PetscCall(PEPSolve(pep));
166: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167: Display solution and clean up
168: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170: /* show detailed info unless -terse option is given by user */
171: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
172: 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: PetscCall(PEPDestroy(&pep));
180: for (i=0;i<NMAT;i++) PetscCall(MatDestroy(&A[i]));
181: PetscCall(SlepcFinalize());
182: return 0;
183: }
185: /*TEST
187: test:
188: suffix: 1
189: args: -pep_type {{toar linear}} -pep_nev 4 -st_type sinvert -terse
190: requires: !single
192: TEST*/