Actual source code: test5.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[] = "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);if (ierr) return ierr;
 52:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 53:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 54:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%D, 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:   MatCreate(PETSC_COMM_WORLD,&Id);
 78:   MatSetSizes(Id,PETSC_DECIDE,PETSC_DECIDE,n,n);
 79:   MatSetFromOptions(Id);
 80:   MatSetUp(Id);
 81:   MatGetOwnershipRange(Id,&Istart,&Iend);
 82:   for (i=Istart;i<Iend;i++) {
 83:     MatSetValue(Id,i,i,1.0,INSERT_VALUES);
 84:   }
 85:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
 86:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
 87:   MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);

 89:   /* A = 1/h^2*tridiag(1,-2,1) + a*I */
 90:   MatCreate(PETSC_COMM_WORLD,&A);
 91:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 92:   MatSetFromOptions(A);
 93:   MatSetUp(A);
 94:   MatGetOwnershipRange(A,&Istart,&Iend);
 95:   for (i=Istart;i<Iend;i++) {
 96:     if (i>0) { MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES); }
 97:     if (i<n-1) { MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES); }
 98:     MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES);
 99:   }
100:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
101:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
102:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);

104:   /* B = diag(b(xi)) */
105:   MatCreate(PETSC_COMM_WORLD,&B);
106:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
107:   MatSetFromOptions(B);
108:   MatSetUp(B);
109:   MatGetOwnershipRange(B,&Istart,&Iend);
110:   for (i=Istart;i<Iend;i++) {
111:     xi = (i+1)*h;
112:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
113:     MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES);
114:   }
115:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
116:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
117:   MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);

119:   /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */
120:   FNCreate(PETSC_COMM_WORLD,&f1);
121:   FNSetType(f1,FNRATIONAL);
122:   coeffs[0] = -1.0; coeffs[1] = 0.0;
123:   FNRationalSetNumerator(f1,2,coeffs);

125:   FNCreate(PETSC_COMM_WORLD,&f2);
126:   FNSetType(f2,FNRATIONAL);
127:   coeffs[0] = 1.0;
128:   FNRationalSetNumerator(f2,1,coeffs);

130:   FNCreate(PETSC_COMM_WORLD,&f3);
131:   FNSetType(f3,FNEXP);
132:   FNSetScale(f3,-tau,1.0);

134:   /* Set the split operator */
135:   mats[0] = A;  funs[0] = f2;
136:   mats[1] = Id; funs[1] = f1;
137:   mats[2] = B;  funs[2] = f3;
138:   NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);

140:   /* Customize nonlinear solver; set runtime options */
141:   NEPSetType(nep,NEPINTERPOL);
142:   NEPSetRG(nep,rg);
143:   NEPInterpolSetPEP(nep,pep);
144:   NEPInterpolGetInterpolation(nep,&tol,&deg);
145:   NEPInterpolSetInterpolation(nep,tol,deg+2);  /* increase degree of interpolation */
146:   NEPSetFromOptions(nep);

148:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149:                       Solve the eigensystem
150:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

152:   NEPSolve(nep);
153:   NEPGetDimensions(nep,&nev,NULL,NULL);
154:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

156:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157:                     Display solution and clean up
158:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

183: /*TEST

185:    test:
186:       suffix: 1
187:       args: -nep_nev 3 -nep_target 5 -terse
188:       requires: !single

190:    test:
191:       suffix: 2_cuda
192:       args: -nep_nev 3 -nep_target 5 -mat_type aijcusparse -terse
193:       requires: cuda
194:       output_file: output/test5_1.out

196: TEST*/