Actual source code: test13.c

slepc-3.20.1 2023-11-27
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 NEPProjectOperator() function.\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";

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

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

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

 25:    Discretization leads to a DDE of dimension n

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

 29:    which results in the nonlinear eigenproblem

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

 34: #include <slepcnep.h>

 36: int main(int argc,char **argv)
 37: {
 38:   NEP            nep;
 39:   Mat            Id,A,B,mats[3];
 40:   FN             f1,f2,f3,funs[3];
 41:   BV             V;
 42:   DS             ds;
 43:   Vec            v;
 44:   PetscScalar    coeffs[2],b,*M;
 45:   PetscInt       n=32,Istart,Iend,i,j,k,nc;
 46:   PetscReal      tau=0.001,h,a=20,xi;

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

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

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

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

 65:   /* A = 1/h^2*tridiag(1,-2,1) + a*I */
 66:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 67:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
 68:   PetscCall(MatSetFromOptions(A));
 69:   PetscCall(MatSetUp(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(MatSetUp(B));
 85:   PetscCall(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:     PetscCall(MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES));
 90:   }
 91:   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
 92:   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
 93:   PetscCall(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE));

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

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

106:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
107:   PetscCall(FNSetType(f3,FNEXP));
108:   PetscCall(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:   PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));
115:   PetscCall(NEPSetType(nep,NEPNARNOLDI));
116:   PetscCall(NEPSetFromOptions(nep));

118:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119:                       Project the NEP
120:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

122:   PetscCall(NEPSetUp(nep));
123:   PetscCall(NEPGetBV(nep,&V));
124:   PetscCall(BVGetSizes(V,NULL,NULL,&nc));
125:   for (i=0;i<nc;i++) {
126:     PetscCall(BVGetColumn(V,i,&v));
127:     PetscCall(VecSetValue(v,i,1.0,INSERT_VALUES));
128:     PetscCall(VecAssemblyBegin(v));
129:     PetscCall(VecAssemblyEnd(v));
130:     PetscCall(BVRestoreColumn(V,i,&v));
131:   }
132:   PetscCall(NEPGetDS(nep,&ds));
133:   PetscCall(DSSetType(ds,DSNEP));
134:   PetscCall(DSNEPSetFN(ds,3,funs));
135:   PetscCall(DSAllocate(ds,nc));
136:   PetscCall(DSSetDimensions(ds,nc,0,0));
137:   PetscCall(NEPProjectOperator(nep,0,nc));

139:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140:                 Display projected matrices and clean up
141:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

143:   for (k=0;k<3;k++) {
144:     PetscCall(DSGetArray(ds,DSMatExtra[k],&M));
145:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nMatrix E%" PetscInt_FMT " = \n",k));
146:     for (i=0;i<nc;i++) {
147:       for (j=0;j<nc;j++) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  %.5g",(double)PetscRealPart(M[i+j*nc])));
148:       PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
149:     }
150:     PetscCall(DSRestoreArray(ds,DSMatExtra[k],&M));
151:   }

153:   PetscCall(NEPDestroy(&nep));
154:   PetscCall(MatDestroy(&Id));
155:   PetscCall(MatDestroy(&A));
156:   PetscCall(MatDestroy(&B));
157:   PetscCall(FNDestroy(&f1));
158:   PetscCall(FNDestroy(&f2));
159:   PetscCall(FNDestroy(&f3));
160:   PetscCall(SlepcFinalize());
161:   return 0;
162: }

164: /*TEST

166:    test:
167:       suffix: 1
168:       args: -nep_ncv 5

170: TEST*/