Actual source code: test11.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 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;

 49:   SlepcInitialize(&argc,&argv,(char*)0,help);
 50:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 51:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 52:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", 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:   MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id);
 63:   MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);

 65:   /* A = 1/h^2*tridiag(1,-2,1) + a*I */
 66:   MatCreate(PETSC_COMM_WORLD,&A);
 67:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 68:   MatSetFromOptions(A);
 69:   MatSetUp(A);
 70:   MatGetOwnershipRange(A,&Istart,&Iend);
 71:   for (i=Istart;i<Iend;i++) {
 72:     if (i>0) MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES);
 73:     if (i<n-1) MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES);
 74:     MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES);
 75:   }
 76:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 77:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 78:   MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);

 80:   /* B = diag(b(xi)) */
 81:   MatCreate(PETSC_COMM_WORLD,&B);
 82:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 83:   MatSetFromOptions(B);
 84:   MatSetUp(B);
 85:   MatGetOwnershipRange(B,&Istart,&Iend);
 86:   for (i=Istart;i<Iend;i++) {
 87:     xi = (i+1)*h;
 88:     b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
 89:     MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES);
 90:   }
 91:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 92:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 93:   MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);

 95:   /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */
 96:   FNCreate(PETSC_COMM_WORLD,&f1);
 97:   FNSetType(f1,FNRATIONAL);
 98:   coeffs[0] = -1.0; coeffs[1] = 0.0;
 99:   FNRationalSetNumerator(f1,2,coeffs);

101:   FNCreate(PETSC_COMM_WORLD,&f2);
102:   FNSetType(f2,FNRATIONAL);
103:   coeffs[0] = 1.0;
104:   FNRationalSetNumerator(f2,1,coeffs);

106:   FNCreate(PETSC_COMM_WORLD,&f3);
107:   FNSetType(f3,FNEXP);
108:   FNSetScale(f3,-tau,1.0);

110:   /* Set the split operator */
111:   mats[0] = A;  funs[0] = f2;
112:   mats[1] = Id; funs[1] = f1;
113:   mats[2] = B;  funs[2] = f3;
114:   NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);

116:   /* Customize nonlinear solver; set runtime options */
117:   NEPSetType(nep,NEPCISS);
118:   NEPSetDimensions(nep,1,24,PETSC_DEFAULT);
119:   NEPSetTolerances(nep,1e-9,PETSC_DEFAULT);
120:   NEPGetRG(nep,&rg);
121:   RGSetType(rg,RGELLIPSE);
122:   RGEllipseSetParameters(rg,10.0,9.5,0.1);
123:   NEPCISSSetSizes(nep,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1,PETSC_DEFAULT,PETSC_TRUE);
124:   NEPCISSGetKSPs(nep,&nsolve,&ksp);
125:   for (i=0;i<nsolve;i++) {
126:     KSPSetTolerances(ksp[i],1e-12,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
127:     KSPSetType(ksp[i],KSPBCGS);
128:     KSPGetPC(ksp[i],&pc);
129:     PCSetType(pc,PCSOR);
130:   }
131:   NEPSetFromOptions(nep);

133:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134:                       Solve the eigensystem
135:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

137:   PetscObjectTypeCompare((PetscObject)nep,NEPCISS,&flg);
138:   if (flg) {
139:     NEPCISSGetExtraction(nep,&ext);
140:     PetscPrintf(PETSC_COMM_WORLD," Running CISS with %" PetscInt_FMT " KSP solvers (%s extraction)\n",nsolve,NEPCISSExtractions[ext]);
141:   }
142:   NEPSolve(nep);

144:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145:                     Display solution and clean up
146:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

148:   /* show detailed info unless -terse option is given by user */
149:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
150:   if (terse) NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
151:   else {
152:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
153:     NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
154:     NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
155:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
156:   }
157:   NEPDestroy(&nep);
158:   MatDestroy(&Id);
159:   MatDestroy(&A);
160:   MatDestroy(&B);
161:   FNDestroy(&f1);
162:   FNDestroy(&f2);
163:   FNDestroy(&f3);
164:   SlepcFinalize();
165:   return 0;
166: }

168: /*TEST

170:    build:
171:       requires: complex

173:    testset:
174:       args: -nep_ciss_extraction {{ritz hankel caa}} -terse
175:       requires: complex !single
176:       output_file: output/test11_1.out
177:       filter: sed -e "s/([A-Z]* extraction)//"
178:       test:
179:          suffix: 1
180:          args: -nep_ciss_delta 1e-10
181:       test:
182:          suffix: 2
183:          nsize: 2
184:          args: -nep_ciss_partitions 2
185:       test:
186:          suffix: 3
187:          args: -nep_ciss_moments 4 -nep_ciss_blocksize 5 -nep_ciss_refine_inner 1 -nep_ciss_refine_blocksize 2

189:    test:
190:       suffix: 4
191:       args: -terse -nep_view
192:       requires: complex !single
193:       filter: grep -v tolerance | grep -v threshold

195: TEST*/