Actual source code: test6.c

slepc-3.17.2 2022-08-09
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 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;

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

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

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

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

 72:   NEPCreate(PETSC_COMM_WORLD,&nep);

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

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

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

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

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

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

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

129:   /* Customize nonlinear solver; set runtime options */
130:   NEPSetOptionsPrefix(nep,"check_");
131:   NEPAppendOptionsPrefix(nep,"myprefix_");
132:   NEPGetOptionsPrefix(nep,&prefix);
133:   PetscPrintf(PETSC_COMM_WORLD,"NEP prefix is currently: %s\n\n",prefix);
134:   NEPSetType(nep,NEPNARNOLDI);
135:   NEPNArnoldiSetKSP(nep,ksp);
136:   if (initv) { /* initial vector */
137:     MatCreateVecs(A,&v0,NULL);
138:     VecGetArray(v0,&pv);
139:     for (i=Istart;i<Iend;i++) pv[i-Istart] = PetscSinReal((4.0*PETSC_PI*i)/n);
140:     VecRestoreArray(v0,&pv);
141:     NEPSetInitialSpace(nep,1,&v0);
142:     VecDestroy(&v0);
143:   }
144:   NEPSetFromOptions(nep);

146:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147:                       Solve the eigensystem
148:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

150:   NEPSolve(nep);
151:   NEPGetDimensions(nep,&nev,NULL,NULL);
152:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
153:   NEPNArnoldiGetLagPreconditioner(nep,&lag);
154:   PetscPrintf(PETSC_COMM_WORLD," N-Arnoldi lag parameter: %" PetscInt_FMT "\n",lag);

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) NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
163:   else {
164:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
165:     NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
166:     NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
167:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
168:   }
169:   NEPDestroy(&nep);
170:   KSPDestroy(&ksp);
171:   MatDestroy(&Id);
172:   MatDestroy(&A);
173:   MatDestroy(&B);
174:   FNDestroy(&f1);
175:   FNDestroy(&f2);
176:   FNDestroy(&f3);
177:   SlepcFinalize();
178:   return 0;
179: }

181: /*TEST

183:    test:
184:       suffix: 1
185:       args: -check_myprefix_nep_view -check_myprefix_nep_monitor_conv -initv -terse
186:       filter: grep -v "tolerance" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g" -e "s/+0i//g"
187:       requires: double

189: TEST*/