Actual source code: ex13.c

slepc-3.18.0 2022-10-01
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[] = "Generalized Symmetric eigenproblem.\n\n"
 12:   "The problem is Ax = lambda Bx, with:\n"
 13:   "   A = Laplacian operator in 2-D\n"
 14:   "   B = diagonal matrix with all values equal to 4 except nulldim zeros\n\n"
 15:   "The command line options are:\n"
 16:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 17:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n"
 18:   "  -nulldim <k>, where <k> = dimension of the nullspace of B.\n\n";

 20: #include <slepceps.h>

 22: int main(int argc,char **argv)
 23: {
 24:   Mat            A,B;         /* matrices */
 25:   EPS            eps;         /* eigenproblem solver context */
 26:   EPSType        type;
 27:   PetscInt       N,n=10,m,Istart,Iend,II,nev,i,j,nulldim=0;
 28:   PetscBool      flag,terse;

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

 33:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 34:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 35:   if (!flag) m=n;
 36:   N = n*m;
 37:   PetscOptionsGetInt(NULL,NULL,"-nulldim",&nulldim,NULL);
 38:   PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized Symmetric Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid), null(B)=%" PetscInt_FMT "\n\n",N,n,m,nulldim);

 40:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 41:      Compute the matrices that define the eigensystem, Ax=kBx
 42:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

 49:   MatCreate(PETSC_COMM_WORLD,&B);
 50:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
 51:   MatSetFromOptions(B);
 52:   MatSetUp(B);

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

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

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

 74:   /*
 75:      Create eigensolver context
 76:   */
 77:   EPSCreate(PETSC_COMM_WORLD,&eps);

 79:   /*
 80:      Set operators. In this case, it is a generalized eigenvalue problem
 81:   */
 82:   EPSSetOperators(eps,A,B);
 83:   EPSSetProblemType(eps,EPS_GHEP);

 85:   /*
 86:      Set solver parameters at runtime
 87:   */
 88:   EPSSetFromOptions(eps);

 90:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 91:                       Solve the eigensystem
 92:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 94:   EPSSolve(eps);

 96:   /*
 97:      Optional: Get some information from the solver and display it
 98:   */
 99:   EPSGetType(eps,&type);
100:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
101:   EPSGetDimensions(eps,&nev,NULL,NULL);
102:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);

104:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105:                     Display solution and clean up
106:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

108:   /* show detailed info unless -terse option is given by user */
109:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
110:   if (terse) EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
111:   else {
112:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
113:     EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
114:     EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
115:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
116:   }
117:   EPSDestroy(&eps);
118:   MatDestroy(&A);
119:   MatDestroy(&B);
120:   SlepcFinalize();
121:   return 0;
122: }

124: /*TEST

126:    test:
127:       suffix: 1
128:       args: -eps_nev 4 -eps_ncv 22 -eps_tol 1e-5 -st_type sinvert -terse
129:       filter: grep -v Solution

131:    test:
132:       suffix: 2
133:       args: -n 110 -nulldim 6 -eps_nev 4 -eps_ncv 18 -eps_tol 1e-5 -eps_purify 1 -st_type sinvert -st_matstructure {{different subset}} -terse
134:       requires: !single

136:    test:
137:       suffix: 3
138:       args: -eps_nev 3 -eps_tol 1e-5 -mat_type sbaij -st_type sinvert -terse

140:    test:
141:       suffix: 4
142:       args: -eps_nev 4 -eps_tol 1e-4 -eps_smallest_real -eps_type {{gd lobpcg rqcg}} -terse
143:       output_file: output/ex13_1.out
144:       filter: grep -v Solution

146:    test:
147:       suffix: 5_primme
148:       args: -n 10 -m 12 -eps_nev 4 -eps_target 0.9 -eps_max_it 15000 -eps_type primme -st_pc_type jacobi -terse
149:       requires: primme defined(SLEPC_HAVE_PRIMME3) !single

151:    test:
152:       suffix: 6
153:       nsize: 2
154:       args: -eps_type ciss -rg_type ellipse -rg_ellipse_center 1.4 -rg_ellipse_radius 0.1 -eps_ciss_partitions 2 -terse
155:       requires: !single

157: TEST*/