Actual source code: test11.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 CISS solver with the problem of ex22.\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;
 38:   Mat            Id,A,B,mats[3];
 39:   FN             f1,f2,f3,funs[3];
 40:   RG             rg;
 41:   KSP            *ksp;
 42:   PC             pc;
 43:   PetscScalar    coeffs[2],b;
 44:   PetscInt       n=128,Istart,Iend,i,nsolve;
 45:   PetscReal      tau=0.001,h,a=20,xi;
 46:   PetscBool      terse;

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

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

 59:   NEPCreate(PETSC_COMM_WORLD,&nep);

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

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

 89:   /* B = diag(b(xi)) */
 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:     MatSetValues(B,1,&i,1,&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:   /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */
105:   FNCreate(PETSC_COMM_WORLD,&f1);
106:   FNSetType(f1,FNRATIONAL);
107:   coeffs[0] = -1.0; coeffs[1] = 0.0;
108:   FNRationalSetNumerator(f1,2,coeffs);

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

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

119:   /* Set the split operator */
120:   mats[0] = A;  funs[0] = f2;
121:   mats[1] = Id; funs[1] = f1;
122:   mats[2] = B;  funs[2] = f3;
123:   NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);

125:   /* Customize nonlinear solver; set runtime options */
126:   NEPSetType(nep,NEPCISS);
127:   NEPSetDimensions(nep,1,24,PETSC_DEFAULT);
128:   NEPSetTolerances(nep,1e-9,PETSC_DEFAULT);
129:   NEPGetRG(nep,&rg);
130:   RGSetType(rg,RGELLIPSE);
131:   RGEllipseSetParameters(rg,10.0,9.5,0.1);
132:   NEPCISSSetSizes(nep,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1,PETSC_DEFAULT,PETSC_TRUE);
133:   NEPCISSGetKSPs(nep,&nsolve,&ksp);
134:   for (i=0;i<nsolve;i++) {
135:     KSPSetTolerances(ksp[i],1e-12,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
136:     KSPSetType(ksp[i],KSPBCGS);
137:     KSPGetPC(ksp[i],&pc);
138:     PCSetType(pc,PCSOR);
139:   }
140:   NEPSetFromOptions(nep);

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

146:   PetscPrintf(PETSC_COMM_WORLD," Running CISS with %D KSP solvers\n",nsolve);
147:   NEPSolve(nep);

149:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150:                     Display solution and clean up
151:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

174: /*TEST

176:    build:
177:       requires: complex

179:    test:
180:       suffix: 1
181:       args: -terse
182:       requires: complex

184: TEST*/