Actual source code: ex11.c

slepc-3.17.1 2022-04-11
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[] = "Computes the smallest nonzero eigenvalue of the Laplacian of a graph.\n\n"
 12:   "This example illustrates EPSSetDeflationSpace(). The example graph corresponds to a "
 13:   "2-D regular mesh. The command line options are:\n"
 14:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 15:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";

 17: #include <slepceps.h>

 19: int main (int argc,char **argv)
 20: {
 21:   EPS            eps;             /* eigenproblem solver context */
 22:   Mat            A;               /* operator matrix */
 23:   Vec            x;
 24:   EPSType        type;
 25:   PetscInt       N,n=10,m,i,j,II,Istart,Iend,nev;
 26:   PetscScalar    w;
 27:   PetscBool      flag,terse;

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

 31:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 32:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 33:   if (!flag) m=n;
 34:   N = n*m;
 35:   PetscPrintf(PETSC_COMM_WORLD,"\nFiedler vector of a 2-D regular mesh, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);

 37:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 38:      Compute the operator matrix that defines the eigensystem, Ax=kx
 39:      In this example, A = L(G), where L is the Laplacian of graph G, i.e.
 40:      Lii = degree of node i, Lij = -1 if edge (i,j) exists in G
 41:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

 59:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 60:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 62:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 63:                 Create the eigensolver and set various options
 64:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 66:   /*
 67:      Create eigensolver context
 68:   */
 69:   EPSCreate(PETSC_COMM_WORLD,&eps);

 71:   /*
 72:      Set operators. In this case, it is a standard eigenvalue problem
 73:   */
 74:   EPSSetOperators(eps,A,NULL);
 75:   EPSSetProblemType(eps,EPS_HEP);

 77:   /*
 78:      Select portion of spectrum
 79:   */
 80:   EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);

 82:   /*
 83:      Set solver parameters at runtime
 84:   */
 85:   EPSSetFromOptions(eps);

 87:   /*
 88:      Attach deflation space: in this case, the matrix has a constant
 89:      nullspace, [1 1 ... 1]^T is the eigenvector of the zero eigenvalue
 90:   */
 91:   MatCreateVecs(A,&x,NULL);
 92:   VecSet(x,1.0);
 93:   EPSSetDeflationSpace(eps,1,&x);
 94:   VecDestroy(&x);

 96:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97:                       Solve the eigensystem
 98:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

100:   EPSSolve(eps);

102:   /*
103:      Optional: Get some information from the solver and display it
104:   */
105:   EPSGetType(eps,&type);
106:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
107:   EPSGetDimensions(eps,&nev,NULL,NULL);
108:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);

110:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111:                     Display solution and clean up
112:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

129: /*TEST

131:    testset:
132:       args: -eps_nev 4 -terse
133:       output_file: output/ex11_1.out
134:       test:
135:          suffix: 1
136:          args: -eps_krylovschur_restart .2
137:       test:
138:          suffix: 2
139:          args: -eps_ncv 20 -eps_target 0 -st_type sinvert -st_ksp_type cg -st_pc_type jacobi

141: TEST*/