Actual source code: test35.c

slepc-3.16.3 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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 BLOPEX.\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,i,j,bs;
 26:   PetscBool      flag;

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

 31:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 32:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 33:   if (!flag) m=n;
 34:   N = n*m;
 35:   PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized eigenproblem with BLOPEX, N=%D (%Dx%D grid)\n\n",N,n,m);

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

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

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

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

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

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

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

 82:   /*
 83:      Set several options
 84:   */
 85:   EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
 86:   EPSGetST(eps,&st);
 87:   STSetType(st,STPRECOND);
 88:   STGetKSP(st,&ksp);
 89:   KSPGetPC(ksp,&pc);
 90:   KSPSetType(ksp,KSPPREONLY);
 91:   PCSetType(pc,PCICC);

 93:   EPSBLOPEXSetBlockSize(eps,4);
 94:   EPSSetFromOptions(eps);

 96:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97:            Compute all eigenvalues in interval and display info
 98:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

100:   EPSSolve(eps);
101:   EPSBLOPEXGetBlockSize(eps,&bs);
102:   PetscPrintf(PETSC_COMM_WORLD," BLOPEX: using block size %D\n",bs);

104:   EPSErrorView(eps,EPS_ERROR_ABSOLUTE,NULL);

106:   EPSDestroy(&eps);
107:   MatDestroy(&A);
108:   MatDestroy(&B);
109:   SlepcFinalize();
110:   return ierr;
111: }

113: /*TEST

115:    build:
116:       requires: blopex

118:    test:
119:       suffix: 1
120:       args: -eps_nev 8 -eps_view
121:       requires: blopex
122:       filter: grep -v tolerance | sed -e "s/hermitian/symmetric/"

124: TEST*/