Actual source code: test13.c

slepc-3.14.1 2020-12-08
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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() operator.\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;

 49:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 50:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 51:   PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
 52:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%D, tau=%g\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:   MatCreate(PETSC_COMM_WORLD,&Id);
 63:   MatSetSizes(Id,PETSC_DECIDE,PETSC_DECIDE,n,n);
 64:   MatSetFromOptions(Id);
 65:   MatSetUp(Id);
 66:   MatGetOwnershipRange(Id,&Istart,&Iend);
 67:   for (i=Istart;i<Iend;i++) {
 68:     MatSetValue(Id,i,i,1.0,INSERT_VALUES);
 69:   }
 70:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
 71:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
 72:   MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);

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

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

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

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

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

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

127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128:                       Project the NEP
129:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

131:   NEPSetUp(nep);
132:   NEPGetBV(nep,&V);
133:   BVGetSizes(V,NULL,NULL,&nc);
134:   for (i=0;i<nc;i++) {
135:     BVGetColumn(V,i,&v);
136:     VecSetValue(v,i,1.0,INSERT_VALUES);
137:     VecAssemblyBegin(v);
138:     VecAssemblyEnd(v);
139:     BVRestoreColumn(V,i,&v);
140:   }
141:   NEPGetDS(nep,&ds);
142:   DSSetType(ds,DSNEP);
143:   DSNEPSetFN(ds,3,funs);
144:   DSAllocate(ds,nc);
145:   DSSetDimensions(ds,nc,0,0,0);
146:   NEPProjectOperator(nep,0,nc);

148:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149:                 Display projected matrices and clean up
150:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

152:   for (k=0;k<3;k++) {
153:     DSGetArray(ds,DSMatExtra[k],&M);
154:     PetscPrintf(PETSC_COMM_WORLD,"\nMatrix E%d = \n",k);
155:     for (i=0;i<nc;i++) {
156:       for (j=0;j<nc;j++) {
157:         PetscPrintf(PETSC_COMM_WORLD,"  %.5g",(double)PetscRealPart(M[i+j*nc]));
158:       }
159:       PetscPrintf(PETSC_COMM_WORLD,"\n");
160:     }
161:     DSRestoreArray(ds,DSMatExtra[k],&M);
162:   }

164:   NEPDestroy(&nep);
165:   MatDestroy(&Id);
166:   MatDestroy(&A);
167:   MatDestroy(&B);
168:   FNDestroy(&f1);
169:   FNDestroy(&f2);
170:   FNDestroy(&f3);
171:   SlepcFinalize();
172:   return ierr;
173: }

175: /*TEST

177:    test:
178:       suffix: 1
179:       args: -nep_ncv 5

181: TEST*/