Actual source code: test35.c

slepc-3.22.1 2024-10-28
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:   PetscFunctionBeginUser;
 29:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));

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

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

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

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

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

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

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

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

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

 91:   PetscCall(EPSBLOPEXSetBlockSize(eps,4));
 92:   PetscCall(EPSSetFromOptions(eps));

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

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

102:   PetscCall(EPSErrorView(eps,EPS_ERROR_ABSOLUTE,NULL));

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

111: /*TEST

113:    build:
114:       requires: blopex

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

122: TEST*/