Actual source code: test11.c

slepc-3.16.0 2021-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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:   NEPCISSExtraction ext;
 44:   PetscScalar       coeffs[2],b;
 45:   PetscInt          n=128,Istart,Iend,i,nsolve;
 46:   PetscReal         tau=0.001,h,a=20,xi;
 47:   PetscBool         flg,terse;
 48:   PetscErrorCode    ierr;

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

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

 60:   NEPCreate(PETSC_COMM_WORLD,&nep);

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

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

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

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

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

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

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

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

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

147:   PetscObjectTypeCompare((PetscObject)nep,NEPCISS,&flg);
148:   if (flg) {
149:     NEPCISSGetExtraction(nep,&ext);
150:     PetscPrintf(PETSC_COMM_WORLD," Running CISS with %D KSP solvers (%s extraction)\n",nsolve,NEPCISSExtractions[ext]);
151:   }
152:   NEPSolve(nep);

154:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155:                     Display solution and clean up
156:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

179: /*TEST

181:    build:
182:       requires: complex

184:    testset:
185:       args: -nep_ciss_extraction {{ritz hankel caa}} -terse
186:       requires: complex !single
187:       output_file: output/test11_1.out
188:       filter: sed -e "s/([A-Z]* extraction)//"
189:       test:
190:          suffix: 1
191:          args: -nep_ciss_delta 1e-10
192:       test:
193:          suffix: 2
194:          nsize: 2
195:          args: -nep_ciss_partitions 2
196:       test:
197:          suffix: 3
198:          args: -nep_ciss_moments 4 -nep_ciss_blocksize 5 -nep_ciss_refine_inner 1 -nep_ciss_refine_blocksize 2

200:    test:
201:       suffix: 4
202:       args: -terse -nep_view
203:       requires: complex !single
204:       filter: grep -v tolerance

206: TEST*/