Actual source code: planar_waveguide.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */
 21: /*
 22:    This example implements one of the problems found at
 23:        NLEVP: A Collection of Nonlinear Eigenvalue Problems,
 24:        The University of Manchester.
 25:    The details of the collection can be found at:
 26:        [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
 27:            Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.

 29:    The planar_waveguide problem is a quartic PEP with symmetric matrices,
 30:    arising from a finite element solution of the propagation constants in a
 31:    six-layer planar waveguide.
 32: */

 34: static char help[] = "NLEVP problem: planar_waveguide.\n\n"
 35:   "The command line options are:\n"
 36:   "  -n <n>, the dimension of the matrices.\n\n";

 38: #include <slepcpep.h>

 40: #define NMAT 5
 41: #define NL   6

 45: int main(int argc,char **argv)
 46: {
 47:   Mat            A[NMAT];         /* problem matrices */
 48:   PEP            pep;             /* polynomial eigenproblem solver context */
 49:   PetscInt       n=128,nlocal,k,Istart,Iend,i,j,start_ct,end_ct;
 50:   PetscReal      w=9.92918,a=0.0,b=2.0,h,deltasq;
 51:   PetscReal      nref[NL],K2[NL],q[NL],*md,*supd,*subd;
 52:   PetscScalar    v,alpha;

 55:   SlepcInitialize(&argc,&argv,(char*)0,help);

 57:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 58:   PetscPrintf(PETSC_COMM_WORLD,"\nPlanar waveguide, n=%D\n\n",n);
 59:   h = (b-a)/n;
 60:   nlocal = (n/4)-1;

 62:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 63:           Set waveguide parameters used in construction of matrices 
 64:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 66:   /* refractive indices in each layer */
 67:   nref[0] = 1.5;
 68:   nref[1] = 1.66;
 69:   nref[2] = 1.6;
 70:   nref[3] = 1.53;
 71:   nref[4] = 1.66;
 72:   nref[5] = 1.0;

 74:   for (i=0;i<NL;i++) K2[i] = w*w*nref[i]*nref[i];
 75:   deltasq = K2[0] - K2[NL-1];
 76:   for (i=0;i<NL;i++) q[i] = K2[i] - (K2[0] + K2[NL-1])/2;

 78:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 79:                      Compute the polynomial matrices 
 80:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 82:   /* initialize matrices */
 83:   for (i=0;i<NMAT;i++) {
 84:     MatCreate(PETSC_COMM_WORLD,&A[i]);
 85:     MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n+1,n+1);
 86:     MatSetFromOptions(A[i]);
 87:     MatSetUp(A[i]);
 88:   }
 89:   MatGetOwnershipRange(A[0],&Istart,&Iend);

 91:   /* A0 */
 92:   alpha = (h/6)*(deltasq*deltasq/16);
 93:   for (i=Istart;i<Iend;i++) {
 94:     v = 4.0;
 95:     if (i==0 || i==n) v = 2.0;
 96:     MatSetValue(A[0],i,i,v*alpha,INSERT_VALUES);
 97:     if (i>0) { MatSetValue(A[0],i,i-1,alpha,INSERT_VALUES); }
 98:     if (i<n) { MatSetValue(A[0],i,i+1,alpha,INSERT_VALUES); }
 99:   }

101:   /* A1 */
102:   if (Istart==0) { MatSetValue(A[1],0,0,-deltasq/4,INSERT_VALUES); }
103:   if (Iend==n+1) { MatSetValue(A[1],n,n,deltasq/4,INSERT_VALUES); }

105:   /* A2 */
106:   alpha = 1.0/h;
107:   for (i=Istart;i<Iend;i++) {
108:     v = 2.0;
109:     if (i==0 || i==n) v = 1.0;
110:     MatSetValue(A[2],i,i,v*alpha,ADD_VALUES);
111:     if (i>0) { MatSetValue(A[2],i,i-1,-alpha,ADD_VALUES); }
112:     if (i<n) { MatSetValue(A[2],i,i+1,-alpha,ADD_VALUES); }
113:   }
114:   PetscMalloc3(n+1,&md,n+1,&supd,n+1,&subd);

116:   md[0]   = 2.0*q[1];
117:   supd[1] = q[1];
118:   subd[0] = q[1];

120:   for (k=1;k<=NL-2;k++) {

122:     end_ct = k*(nlocal+1);
123:     start_ct = end_ct-nlocal;

125:     for (j=start_ct;j<end_ct;j++) {
126:       md[j] = 4*q[k];
127:       supd[j+1] = q[k];
128:       subd[j] = q[k];
129:     }

131:     if (k < 4) {  /* interface points */
132:       md[end_ct] = 4*(q[k] + q[k+1])/2.0;
133:       supd[end_ct+1] = q[k+1];
134:       subd[end_ct] = q[k+1];
135:     }

137:   }

139:   md[n] = 2*q[NL-2];
140:   supd[n] = q[NL-2];
141:   subd[n] = q[NL-2];

143:   alpha = -h/6.0;
144:   for (i=Istart;i<Iend;i++) {
145:     MatSetValue(A[2],i,i,md[i]*alpha,ADD_VALUES);
146:     if (i>0) { MatSetValue(A[2],i,i-1,subd[i]*alpha,ADD_VALUES); }
147:     if (i<n) { MatSetValue(A[2],i,i+1,supd[i+1]*alpha,ADD_VALUES); }
148:   }
149:   PetscFree3(md,supd,subd);

151:   /* A3 */
152:   if (Istart==0) { MatSetValue(A[3],0,0,1.0,INSERT_VALUES); }
153:   if (Iend==n+1) { MatSetValue(A[3],n,n,1.0,INSERT_VALUES); }

155:   /* A4 */
156:   alpha = (h/6);
157:   for (i=Istart;i<Iend;i++) {
158:     v = 4.0;
159:     if (i==0 || i==n) v = 2.0;
160:     MatSetValue(A[4],i,i,v*alpha,INSERT_VALUES);
161:     if (i>0) { MatSetValue(A[4],i,i-1,alpha,INSERT_VALUES); }
162:     if (i<n) { MatSetValue(A[4],i,i+1,alpha,INSERT_VALUES); }
163:   }

165:   /* assemble matrices */
166:   for (i=0;i<NMAT;i++) {
167:     MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY);
168:   }
169:   for (i=0;i<NMAT;i++) {
170:     MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY);
171:   }

173:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
174:                 Create the eigensolver and solve the problem
175:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

177:   PEPCreate(PETSC_COMM_WORLD,&pep);
178:   PEPSetOperators(pep,NMAT,A);
179:   PEPSetFromOptions(pep);
180:   PEPSolve(pep);

182:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183:                     Display solution and clean up
184:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
185:   
186:   PEPPrintSolution(pep,NULL);
187:   PEPDestroy(&pep);
188:   for (i=0;i<NMAT;i++) {
189:     MatDestroy(&A[i]);
190:   }
191:   SlepcFinalize();
192:   return 0;
193: }