Actual source code: ex22.c

slepc-3.18.2 2023-01-26
Report Typos and Errors
  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: */

 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);
 49:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 50:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 51:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", 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:   MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id);
 68:   MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);

 70:   /*
 71:      A = 1/h^2*tridiag(1,-2,1) + a*I
 72:   */
 73:   MatCreate(PETSC_COMM_WORLD,&A);
 74:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 75:   MatSetFromOptions(A);
 76:   MatSetUp(A);
 77:   MatGetOwnershipRange(A,&Istart,&Iend);
 78:   for (i=Istart;i<Iend;i++) {
 79:     if (i>0) MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES);
 80:     if (i<n-1) MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES);
 81:     MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES);
 82:   }
 83:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 84:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 85:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);

 87:   /*
 88:      B = diag(b(xi))
 89:   */
 90:   MatCreate(PETSC_COMM_WORLD,&B);
 91:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 92:   MatSetFromOptions(B);
 93:   MatSetUp(B);
 94:   MatGetOwnershipRange(B,&Istart,&Iend);
 95:   for (i=Istart;i<Iend;i++) {
 96:     xi = (i+1)*h;
 97:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
 98:     MatSetValue(B,i,i,b,INSERT_VALUES);
 99:   }
100:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
101:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
102:   MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);

104:   /*
105:      Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda)
106:   */
107:   FNCreate(PETSC_COMM_WORLD,&f1);
108:   FNSetType(f1,FNRATIONAL);
109:   coeffs[0] = -1.0; coeffs[1] = 0.0;
110:   FNRationalSetNumerator(f1,2,coeffs);

112:   FNCreate(PETSC_COMM_WORLD,&f2);
113:   FNSetType(f2,FNRATIONAL);
114:   coeffs[0] = 1.0;
115:   FNRationalSetNumerator(f2,1,coeffs);

117:   FNCreate(PETSC_COMM_WORLD,&f3);
118:   FNSetType(f3,FNEXP);
119:   FNSetScale(f3,-tau,1.0);

121:   /*
122:      Set the split operator. Note that A is passed first so that
123:      SUBSET_NONZERO_PATTERN can be used
124:   */
125:   mats[0] = A;  funs[0] = f2;
126:   mats[1] = Id; funs[1] = f1;
127:   mats[2] = B;  funs[2] = f3;
128:   NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);

130:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131:              Customize nonlinear solver; set runtime options
132:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

134:   NEPSetTolerances(nep,1e-9,PETSC_DEFAULT);
135:   NEPSetDimensions(nep,1,PETSC_DEFAULT,PETSC_DEFAULT);
136:   NEPRIISetLagPreconditioner(nep,0);

138:   /*
139:      Set solver parameters at runtime
140:   */
141:   NEPSetFromOptions(nep);

143:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144:                       Solve the eigensystem
145:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

147:   NEPSolve(nep);
148:   NEPGetDimensions(nep,&nev,NULL,NULL);
149:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);

151:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152:                     Display solution and clean up
153:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

155:   /* show detailed info unless -terse option is given by user */
156:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
157:   if (terse) NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
158:   else {
159:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
160:     NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
161:     NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
162:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
163:   }
164:   NEPDestroy(&nep);
165:   MatDestroy(&Id);
166:   MatDestroy(&A);
167:   MatDestroy(&B);
168:   FNDestroy(&f1);
169:   FNDestroy(&f2);
170:   FNDestroy(&f3);
171:   SlepcFinalize();
172:   return 0;
173: }

175: /*TEST

177:    testset:
178:       suffix: 1
179:       args: -nep_type {{rii slp narnoldi}} -terse
180:       filter: sed -e "s/[+-]0\.0*i//g"
181:       requires: !single

183:    test:
184:       suffix: 1_ciss
185:       args: -nep_type ciss -nep_ciss_extraction {{ritz hankel caa}} -rg_type ellipse -rg_ellipse_center 10 -rg_ellipse_radius 9.5 -nep_ncv 24 -terse
186:       requires: complex !single

188:    test:
189:       suffix: 2
190:       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
191:       filter: sed -e "s/[+-]0\.0*i//g"
192:       requires: !single

194:    testset:
195:       args: -n 512 -nep_target 10 -nep_nev 3 -terse
196:       filter: sed -e "s/[+-]0\.0*i//g"
197:       requires: !single
198:       output_file: output/ex22_3.out
199:       test:
200:          suffix: 3
201:          args: -nep_type {{rii slp narnoldi}}
202:       test:
203:          suffix: 3_simpleu
204:          args: -nep_type {{rii slp narnoldi}} -nep_deflation_simpleu
205:       test:
206:          suffix: 3_slp_thres
207:          args: -nep_type slp -nep_slp_deflation_threshold 1e-8
208:       test:
209:          suffix: 3_rii_thres
210:          args: -nep_type rii -nep_rii_deflation_threshold 1e-8

212:    test:
213:       suffix: 4
214:       args: -nep_type interpol -rg_type interval -rg_interval_endpoints 5,20,-.1,.1 -nep_nev 3 -nep_target 5 -terse -nep_monitor draw::draw_lg
215:       filter: sed -e "s/[+-]0\.0*i//g"
216:       requires: x !single
217:       output_file: output/ex22_2.out

219: TEST*/