Actual source code: ex11.c

slepc-3.11.1 2019-04-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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;

 30:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 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,"\nFiedler vector of a 2-D regular mesh, N=%D (%Dx%D grid)\n\n",N,n,m);

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

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

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

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

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

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

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

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

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

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

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

101:   EPSSolve(eps);

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

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

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

131: /*TEST

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

144: TEST*/