Actual source code: test33.c

slepc-3.14.0 2020-09-30
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 interface to external package BLZPACK.\n\n"
 12:   "This is based on ex12.c. The command line options are:\n"
 13:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 14:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";

 16: #include <slepceps.h>

 18: int main(int argc,char **argv)
 19: {
 20:   Mat            A,B;         /* matrices */
 21:   EPS            eps;         /* eigenproblem solver context */
 22:   ST             st;          /* spectral transformation context */
 23:   KSP            ksp;
 24:   PC             pc;
 25:   PetscInt       N,n=35,m,Istart,Iend,II,nev,i,j,bs,nsteps;
 26:   PetscBool      flag;
 27:   PetscReal      int0,int1;

 30:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 32:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 33:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 34:   if (!flag) m=n;
 35:   N = n*m;
 36:   PetscPrintf(PETSC_COMM_WORLD,"\nSpectrum-slicing test, N=%D (%Dx%D grid)\n\n",N,n,m);

 38:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 39:      Compute the matrices that define the eigensystem, Ax=kBx
 40:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 42:   MatCreate(PETSC_COMM_WORLD,&A);
 43:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 44:   MatSetFromOptions(A);
 45:   MatSetUp(A);

 47:   MatCreate(PETSC_COMM_WORLD,&B);
 48:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
 49:   MatSetFromOptions(B);
 50:   MatSetUp(B);

 52:   MatGetOwnershipRange(A,&Istart,&Iend);
 53:   for (II=Istart;II<Iend;II++) {
 54:     i = II/n; j = II-i*n;
 55:     if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
 56:     if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
 57:     if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
 58:     if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
 59:     MatSetValue(A,II,II,4.0,INSERT_VALUES);
 60:     MatSetValue(B,II,II,2.0,INSERT_VALUES);
 61:   }
 62:   if (Istart==0) {
 63:     MatSetValue(B,0,0,6.0,INSERT_VALUES);
 64:     MatSetValue(B,0,1,-1.0,INSERT_VALUES);
 65:     MatSetValue(B,1,0,-1.0,INSERT_VALUES);
 66:     MatSetValue(B,1,1,1.0,INSERT_VALUES);
 67:   }

 69:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 70:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 71:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 72:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 74:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 75:                 Create the eigensolver and set various options
 76:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 78:   EPSCreate(PETSC_COMM_WORLD,&eps);
 79:   EPSSetOperators(eps,A,B);
 80:   EPSSetProblemType(eps,EPS_GHEP);
 81:   EPSSetType(eps,EPSBLZPACK);

 83:   /*
 84:      Set interval and other settings for spectrum slicing
 85:   */
 86:   EPSSetWhichEigenpairs(eps,EPS_ALL);
 87:   int0 = 1.1; int1 = 1.2;
 88:   EPSSetInterval(eps,int0,int1);
 89:   EPSSetDimensions(eps,100,PETSC_DEFAULT,PETSC_DEFAULT);
 90:   EPSGetST(eps,&st);
 91:   STSetType(st,STSINVERT);
 92:   STGetKSP(st,&ksp);
 93:   KSPGetPC(ksp,&pc);
 94:   KSPSetType(ksp,KSPPREONLY);
 95:   PCSetType(pc,PCCHOLESKY);

 97:   EPSBlzpackSetBlockSize(eps,4);
 98:   EPSBlzpackSetNSteps(eps,20);
 99:   EPSSetFromOptions(eps);

101:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102:            Compute all eigenvalues in interval and display info
103:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

105:   EPSSolve(eps);
106:   EPSBlzpackGetBlockSize(eps,&bs);
107:   EPSBlzpackGetNSteps(eps,&nsteps);
108:   PetscPrintf(PETSC_COMM_WORLD," Using block size %D, maximum number of steps %D\n",bs,nsteps);
109:   EPSGetDimensions(eps,&nev,NULL,NULL);
110:   EPSGetInterval(eps,&int0,&int1);
111:   PetscPrintf(PETSC_COMM_WORLD," Found %D eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1);

113:   EPSErrorView(eps,EPS_ERROR_ABSOLUTE,NULL);

115:   EPSDestroy(&eps);
116:   MatDestroy(&A);
117:   MatDestroy(&B);
118:   SlepcFinalize();
119:   return ierr;
120: }

122: /*TEST

124:    build:
125:       requires: blzpack

127:    test:
128:       suffix: 1
129:       requires: blzpack

131: TEST*/