Actual source code: test6.c

slepc-3.13.1 2020-04-12
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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 NArnoldi solver with a user-provided KSP.\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"
 16:   "  -initv ... set an initial vector.\n\n";

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

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

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

 26:    Discretization leads to a DDE of dimension n

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

 30:    which results in the nonlinear eigenproblem

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

 35: #include <slepcnep.h>

 37: int main(int argc,char **argv)
 38: {
 39:   NEP            nep;
 40:   KSP            ksp;
 41:   PC             pc;
 42:   Mat            Id,A,B,mats[3];
 43:   FN             f1,f2,f3,funs[3];
 44:   Vec            v0;
 45:   PetscScalar    coeffs[2],b,*pv;
 46:   PetscInt       n=128,nev,Istart,Iend,i,lag;
 47:   PetscReal      tau=0.001,h,a=20,xi;
 48:   PetscBool      terse,initv=PETSC_FALSE;
 49:   const char     *prefix;

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

 59:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 60:       Create a standalone KSP with appropriate settings
 61:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 63:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 64:   KSPSetType(ksp,KSPBCGS);
 65:   KSPGetPC(ksp,&pc);
 66:   PCSetType(pc,PCBJACOBI);
 67:   KSPSetFromOptions(ksp);

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

 73:   NEPCreate(PETSC_COMM_WORLD,&nep);

 75:   /* Identity matrix */
 76:   MatCreate(PETSC_COMM_WORLD,&Id);
 77:   MatSetSizes(Id,PETSC_DECIDE,PETSC_DECIDE,n,n);
 78:   MatSetFromOptions(Id);
 79:   MatSetUp(Id);
 80:   MatGetOwnershipRange(Id,&Istart,&Iend);
 81:   for (i=Istart;i<Iend;i++) {
 82:     MatSetValue(Id,i,i,1.0,INSERT_VALUES);
 83:   }
 84:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
 85:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
 86:   MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);

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

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

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

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

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

133:   /* Set the split operator */
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:   /* Customize nonlinear solver; set runtime options */
140:   NEPSetOptionsPrefix(nep,"check_");
141:   NEPAppendOptionsPrefix(nep,"myprefix_");
142:   NEPGetOptionsPrefix(nep,&prefix);
143:   PetscPrintf(PETSC_COMM_WORLD,"NEP prefix is currently: %s\n\n",prefix);
144:   NEPSetType(nep,NEPNARNOLDI);
145:   NEPNArnoldiSetKSP(nep,ksp);
146:   if (initv) { /* initial vector */
147:     MatCreateVecs(A,&v0,NULL);
148:     VecGetArray(v0,&pv);
149:     for (i=Istart;i<Iend;i++) pv[i-Istart] = PetscSinReal((4.0*PETSC_PI*i)/n);
150:     VecRestoreArray(v0,&pv);
151:     NEPSetInitialSpace(nep,1,&v0);
152:     VecDestroy(&v0);
153:   }
154:   NEPSetFromOptions(nep);

156:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157:                       Solve the eigensystem
158:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

160:   NEPSolve(nep);
161:   NEPGetDimensions(nep,&nev,NULL,NULL);
162:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
163:   NEPNArnoldiGetLagPreconditioner(nep,&lag);
164:   PetscPrintf(PETSC_COMM_WORLD," N-Arnoldi lag parameter: %D\n",lag);

166:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167:                     Display solution and clean up
168:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

192: /*TEST

194:    test:
195:       suffix: 1
196:       args: -check_myprefix_nep_view -check_myprefix_nep_monitor_conv -initv -terse
197:       filter: grep -v "tolerance" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g" -e "s/+0i//"
198:       requires: double

200: TEST*/