Actual source code: test34.c

slepc-3.19.0 2023-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 PRIMME.\n\n"
 12:   "This is based on ex3.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;           /* matrix */
 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;
 27:   EPSPRIMMEMethod meth;

 29:   PetscFunctionBeginUser;
 30:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

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

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

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

 47:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 48:   for (II=Istart;II<Iend;II++) {
 49:     i = II/n; j = II-i*n;
 50:     if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
 51:     if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
 52:     if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
 53:     if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
 54:     PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
 55:   }

 57:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 58:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));

 60:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 61:                 Create the eigensolver and set various options
 62:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 64:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 65:   PetscCall(EPSSetOperators(eps,A,NULL));
 66:   PetscCall(EPSSetProblemType(eps,EPS_HEP));
 67:   PetscCall(EPSSetType(eps,EPSPRIMME));

 69:   /*
 70:      Set several options
 71:   */
 72:   PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
 73:   PetscCall(EPSGetST(eps,&st));
 74:   PetscCall(STSetType(st,STPRECOND));
 75:   PetscCall(STGetKSP(st,&ksp));
 76:   PetscCall(KSPGetPC(ksp,&pc));
 77:   PetscCall(KSPSetType(ksp,KSPPREONLY));
 78:   PetscCall(PCSetType(pc,PCICC));

 80:   PetscCall(EPSPRIMMESetBlockSize(eps,4));
 81:   PetscCall(EPSPRIMMESetMethod(eps,EPS_PRIMME_GD_OLSEN_PLUSK));
 82:   PetscCall(EPSSetFromOptions(eps));

 84:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 85:                  Compute eigenvalues and display info
 86:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 88:   PetscCall(EPSSolve(eps));
 89:   PetscCall(EPSPRIMMEGetBlockSize(eps,&bs));
 90:   PetscCall(EPSPRIMMEGetMethod(eps,&meth));
 91:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," PRIMME: using block size %" PetscInt_FMT ", method %s\n",bs,EPSPRIMMEMethods[meth]));

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

 95:   PetscCall(EPSDestroy(&eps));
 96:   PetscCall(MatDestroy(&A));
 97:   PetscCall(SlepcFinalize());
 98:   return 0;
 99: }

101: /*TEST

103:    build:
104:       requires: primme

106:    testset:
107:       args: -eps_nev 4
108:       requires: primme
109:       output_file: output/test34_1.out
110:       test:
111:          suffix: 1
112:       test:
113:          suffix: 2
114:          args: -st_pc_type bjacobi -eps_target 0.01 -eps_target_real -eps_refined
115:          nsize: 2
116:       test:
117:          suffix: 3
118:          args: -eps_smallest_magnitude -eps_harmonic

120: TEST*/