Actual source code: test11.c

slepc-3.21.1 2024-04-26
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:   PetscFunctionBeginUser;
 50:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 51:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 52:   PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
 53:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau));
 54:   h = PETSC_PI/(PetscReal)(n+1);

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

 60:   PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));

 62:   /* Identity matrix */
 63:   PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id));
 64:   PetscCall(MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE));

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

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

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

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

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

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

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

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

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

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

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

167: /*TEST

169:    build:
170:       requires: complex

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

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

194: TEST*/