Actual source code: test5.c

slepc-3.18.1 2022-11-02
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[] = "Test the INTERPOL solver with a user-provided PEP.\n\n"
 12:   "This is based on ex22.\n"
 13:   "The command line options are:\n"
 14:   "  -n <n>, where <n> = number of grid subdivisions.\n"
 15:   "  -tau <tau>, where <tau> is the delay parameter.\n\n";

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

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

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

 25:    Discretization leads to a DDE of dimension n

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

 29:    which results in the nonlinear eigenproblem

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

 34: #include <slepcnep.h>

 36: int main(int argc,char **argv)
 37: {
 38:   NEP            nep;
 39:   PEP            pep;
 40:   Mat            Id,A,B;
 41:   FN             f1,f2,f3;
 42:   RG             rg;
 43:   Mat            mats[3];
 44:   FN             funs[3];
 45:   PetscScalar    coeffs[2],b;
 46:   PetscInt       n=128,nev,Istart,Iend,i,deg;
 47:   PetscReal      tau=0.001,h,a=20,xi,tol;
 48:   PetscBool      terse;

 51:   SlepcInitialize(&argc,&argv,(char*)0,help);
 52:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 53:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 54:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau);
 55:   h = PETSC_PI/(PetscReal)(n+1);

 57:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 58:       Create a standalone PEP and RG objects with appropriate settings
 59:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 61:   PEPCreate(PETSC_COMM_WORLD,&pep);
 62:   PEPSetType(pep,PEPTOAR);
 63:   PEPSetFromOptions(pep);

 65:   RGCreate(PETSC_COMM_WORLD,&rg);
 66:   RGSetType(rg,RGINTERVAL);
 67:   RGIntervalSetEndpoints(rg,5,20,-0.1,0.1);
 68:   RGSetFromOptions(rg);

 70:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 71:      Create nonlinear eigensolver context
 72:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 74:   NEPCreate(PETSC_COMM_WORLD,&nep);

 76:   /* Identity matrix */
 77:   MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id);
 78:   MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);

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

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

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

116:   FNCreate(PETSC_COMM_WORLD,&f2);
117:   FNSetType(f2,FNRATIONAL);
118:   coeffs[0] = 1.0;
119:   FNRationalSetNumerator(f2,1,coeffs);

121:   FNCreate(PETSC_COMM_WORLD,&f3);
122:   FNSetType(f3,FNEXP);
123:   FNSetScale(f3,-tau,1.0);

125:   /* Set the split operator */
126:   mats[0] = A;  funs[0] = f2;
127:   mats[1] = Id; funs[1] = f1;
128:   mats[2] = B;  funs[2] = f3;
129:   NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);

131:   /* Customize nonlinear solver; set runtime options */
132:   NEPSetType(nep,NEPINTERPOL);
133:   NEPSetRG(nep,rg);
134:   NEPInterpolSetPEP(nep,pep);
135:   NEPInterpolGetInterpolation(nep,&tol,&deg);
136:   NEPInterpolSetInterpolation(nep,tol,deg+2);  /* increase degree of interpolation */
137:   NEPSetFromOptions(nep);

139:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140:                       Solve the eigensystem
141:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

143:   NEPSolve(nep);
144:   NEPGetDimensions(nep,&nev,NULL,NULL);
145:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);

147:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148:                     Display solution and clean up
149:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

173: /*TEST

175:    testset:
176:       args: -nep_nev 3 -nep_target 5 -terse
177:       output_file: output/test5_1.out
178:       filter: sed -e "s/[+-]0\.0*i//g"
179:       requires: !single
180:       test:
181:          suffix: 1
182:       test:
183:          suffix: 2_cuda
184:          args: -mat_type aijcusparse
185:          requires: cuda
186:       test:
187:          suffix: 3
188:          args: -nep_view_values draw
189:          requires: x

191: TEST*/