Actual source code: test35.c

slepc-3.17.0 2022-03-31
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 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;

 28:   SlepcInitialize(&argc,&argv,(char*)0,help);

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

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

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

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

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

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

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

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

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

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

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

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

103:   EPSErrorView(eps,EPS_ERROR_ABSOLUTE,NULL);

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

112: /*TEST

114:    build:
115:       requires: blopex

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

123: TEST*/