Actual source code: ex22.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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: */

 11: static char help[] = "Delay differential equation.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions.\n"
 14:   "  -tau <tau>, where <tau> is the delay parameter.\n\n";

 16: /*
 17:    Solve parabolic partial differential equation with time delay tau

 19:             u_t = u_xx + a*u(t) + b*u(t-tau)
 20:             u(0,t) = u(pi,t) = 0

 22:    with a = 20 and b(x) = -4.1+x*(1-exp(x-pi)).

 24:    Discretization leads to a DDE of dimension n

 26:             -u' = A*u(t) + B*u(t-tau)

 28:    which results in the nonlinear eigenproblem

 30:             (-lambda*I + A + exp(-tau*lambda)*B)*u = 0
 31: */

 33: #include <slepcnep.h>

 35: int main(int argc,char **argv)
 36: {
 37:   NEP            nep;             /* nonlinear eigensolver context */
 38:   Mat            Id,A,B;          /* problem matrices */
 39:   FN             f1,f2,f3;        /* functions to define the nonlinear operator */
 40:   Mat            mats[3];
 41:   FN             funs[3];
 42:   PetscScalar    coeffs[2],b;
 43:   PetscInt       n=128,nev,Istart,Iend,i;
 44:   PetscReal      tau=0.001,h,a=20,xi;
 45:   PetscBool      terse;

 48:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 49:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 50:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 51:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%D, tau=%g\n\n",n,(double)tau);
 52:   h = PETSC_PI/(PetscReal)(n+1);

 54:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 55:      Create nonlinear eigensolver context
 56:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 58:   NEPCreate(PETSC_COMM_WORLD,&nep);

 60:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 61:      Create problem matrices and coefficient functions. Pass them to NEP
 62:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 64:   /*
 65:      Identity matrix
 66:   */
 67:   MatCreate(PETSC_COMM_WORLD,&Id);
 68:   MatSetSizes(Id,PETSC_DECIDE,PETSC_DECIDE,n,n);
 69:   MatSetFromOptions(Id);
 70:   MatSetUp(Id);
 71:   MatGetOwnershipRange(Id,&Istart,&Iend);
 72:   for (i=Istart;i<Iend;i++) {
 73:     MatSetValue(Id,i,i,1.0,INSERT_VALUES);
 74:   }
 75:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
 76:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
 77:   MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);

 79:   /*
 80:      A = 1/h^2*tridiag(1,-2,1) + a*I
 81:   */
 82:   MatCreate(PETSC_COMM_WORLD,&A);
 83:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 84:   MatSetFromOptions(A);
 85:   MatSetUp(A);
 86:   MatGetOwnershipRange(A,&Istart,&Iend);
 87:   for (i=Istart;i<Iend;i++) {
 88:     if (i>0) { MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES); }
 89:     if (i<n-1) { MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES); }
 90:     MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES);
 91:   }
 92:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 93:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 94:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);

 96:   /*
 97:      B = diag(b(xi))
 98:   */
 99:   MatCreate(PETSC_COMM_WORLD,&B);
100:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
101:   MatSetFromOptions(B);
102:   MatSetUp(B);
103:   MatGetOwnershipRange(B,&Istart,&Iend);
104:   for (i=Istart;i<Iend;i++) {
105:     xi = (i+1)*h;
106:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
107:     MatSetValue(B,i,i,b,INSERT_VALUES);
108:   }
109:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
110:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
111:   MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);

113:   /*
114:      Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda)
115:   */
116:   FNCreate(PETSC_COMM_WORLD,&f1);
117:   FNSetType(f1,FNRATIONAL);
118:   coeffs[0] = -1.0; coeffs[1] = 0.0;
119:   FNRationalSetNumerator(f1,2,coeffs);

121:   FNCreate(PETSC_COMM_WORLD,&f2);
122:   FNSetType(f2,FNRATIONAL);
123:   coeffs[0] = 1.0;
124:   FNRationalSetNumerator(f2,1,coeffs);

126:   FNCreate(PETSC_COMM_WORLD,&f3);
127:   FNSetType(f3,FNEXP);
128:   FNSetScale(f3,-tau,1.0);

130:   /*
131:      Set the split operator. Note that A is passed first so that
132:      SUBSET_NONZERO_PATTERN can be used
133:   */
134:   mats[0] = A;  funs[0] = f2;
135:   mats[1] = Id; funs[1] = f1;
136:   mats[2] = B;  funs[2] = f3;
137:   NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);

139:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140:              Customize nonlinear solver; set runtime options
141:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

143:   NEPSetTolerances(nep,1e-9,PETSC_DEFAULT);
144:   NEPSetDimensions(nep,1,PETSC_DEFAULT,PETSC_DEFAULT);
145:   NEPRIISetLagPreconditioner(nep,0);

147:   /*
148:      Set solver parameters at runtime
149:   */
150:   NEPSetFromOptions(nep);

152:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153:                       Solve the eigensystem
154:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

156:   NEPSolve(nep);
157:   NEPGetDimensions(nep,&nev,NULL,NULL);
158:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

160:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161:                     Display solution and clean up
162:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

164:   /* show detailed info unless -terse option is given by user */
165:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
166:   if (terse) {
167:     NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
168:   } else {
169:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
170:     NEPReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
171:     NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
172:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
173:   }
174:   NEPDestroy(&nep);
175:   MatDestroy(&Id);
176:   MatDestroy(&A);
177:   MatDestroy(&B);
178:   FNDestroy(&f1);
179:   FNDestroy(&f2);
180:   FNDestroy(&f3);
181:   SlepcFinalize();
182:   return ierr;
183: }

185: /*TEST

187:    testset:
188:       args: -terse
189:       requires: !single
190:       output_file: output/ex22_1.out
191:       test:
192:          suffix: 1
193:          args: -nep_type {{rii slp narnoldi}}

195:    test:
196:       suffix: 1_ciss
197:       args: -nep_type ciss -rg_type ellipse -rg_ellipse_center 10 -rg_ellipse_radius 9.5 -nep_ncv 24 -terse
198:       requires: complex

200:    test:
201:       suffix: 2
202:       args: -nep_type interpol -nep_interpol_pep_extract {{none norm residual}} -rg_type interval -rg_interval_endpoints 5,20,-.1,.1 -nep_nev 3 -nep_target 5 -terse
203:       requires: !single

205:    test:
206:       suffix: 3
207:       args: -n 512 -nep_target 10 -nep_nev 3 -nep_type {{rii slp narnoldi}} -terse
208:       requires: !single

210: TEST*/