Actual source code: ex12.c

slepc-3.17.2 2022-08-09
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[] = "Compute all eigenvalues in an interval of a symmetric-definite problem.\n\n"
 12:   "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,k,*inertias;
 26:   PetscBool      flag,showinertia=PETSC_TRUE;
 27:   PetscReal      int0,int1,*shifts;

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

 31:   PetscOptionsGetBool(NULL,NULL,"-showinertia",&showinertia,NULL);
 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,"\nSymmetric-definite problem with two intervals, 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:   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);

 82:   /*
 83:      Set first interval and other settings for spectrum slicing
 84:   */
 85:   EPSSetType(eps,EPSKRYLOVSCHUR);
 86:   EPSSetWhichEigenpairs(eps,EPS_ALL);
 87:   EPSSetInterval(eps,1.1,1.3);
 88:   EPSGetST(eps,&st);
 89:   STSetType(st,STSINVERT);
 90:   EPSKrylovSchurGetKSP(eps,&ksp);
 91:   KSPGetPC(ksp,&pc);
 92:   KSPSetType(ksp,KSPPREONLY);
 93:   PCSetType(pc,PCCHOLESKY);
 94:   EPSSetFromOptions(eps);

 96:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97:                  Solve for first interval and display info
 98:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

100:   EPSSolve(eps);
101:   EPSGetDimensions(eps,&nev,NULL,NULL);
102:   EPSGetInterval(eps,&int0,&int1);
103:   PetscPrintf(PETSC_COMM_WORLD," Found %" PetscInt_FMT " eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1);
104:   if (showinertia) {
105:     EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias);
106:     PetscPrintf(PETSC_COMM_WORLD," Used %" PetscInt_FMT " shifts (inertia):\n",k);
107:     for (i=0;i<k;i++) PetscPrintf(PETSC_COMM_WORLD," .. %g (%" PetscInt_FMT ")\n",(double)shifts[i],inertias[i]);
108:     PetscFree(shifts);
109:     PetscFree(inertias);
110:   }

112:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113:                  Solve for second interval and display info
114:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115:   EPSSetInterval(eps,1.499,1.6);
116:   EPSSolve(eps);
117:   EPSGetDimensions(eps,&nev,NULL,NULL);
118:   EPSGetInterval(eps,&int0,&int1);
119:   PetscPrintf(PETSC_COMM_WORLD," Found %" PetscInt_FMT " eigenvalues in interval [%g,%g]\n",nev,(double)int0,(double)int1);
120:   if (showinertia) {
121:     EPSKrylovSchurGetInertias(eps,&k,&shifts,&inertias);
122:     PetscPrintf(PETSC_COMM_WORLD," Used %" PetscInt_FMT " shifts (inertia):\n",k);
123:     for (i=0;i<k;i++) PetscPrintf(PETSC_COMM_WORLD," .. %g (%" PetscInt_FMT ")\n",(double)shifts[i],inertias[i]);
124:     PetscFree(shifts);
125:     PetscFree(inertias);
126:   }

128:   EPSDestroy(&eps);
129:   MatDestroy(&A);
130:   MatDestroy(&B);
131:   SlepcFinalize();
132:   return 0;
133: }

135: /*TEST

137:    test:
138:       suffix: 1
139:       args: -showinertia 0 -eps_error_relative
140:       requires: !single

142: TEST*/