Actual source code: test10.c

slepc-3.19.2 2023-09-05
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:   PetscInt       N,n=10,m,i,j,II,Istart,Iend,nev;
 25:   PetscScalar    w;
 26:   PetscBool      flag;

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

 31:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 32:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
 33:   if (!flag) m=n;
 34:   N = n*m;
 35:   PetscCall(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:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 44:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
 45:   PetscCall(MatSetFromOptions(A));
 46:   PetscCall(MatSetUp(A));

 48:   PetscCall(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) { PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES)); w=w+1.0; }
 53:     if (i<m-1) { PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES)); w=w+1.0; }
 54:     if (j>0) { PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES)); w=w+1.0; }
 55:     if (j<n-1) { PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES)); w=w+1.0; }
 56:     PetscCall(MatSetValue(A,II,II,w,INSERT_VALUES));
 57:   }

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

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

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

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

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

 82:   /*
 83:      Set solver parameters at runtime
 84:   */
 85:   PetscCall(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:   PetscCall(MatCreateVecs(A,&x,NULL));
 92:   PetscCall(VecSet(x,1.0));
 93:   PetscCall(EPSSetDeflationSpace(eps,1,&x));
 94:   PetscCall(VecDestroy(&x));

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

100:   PetscCall(EPSSolve(eps));
101:   PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
102:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));

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

108:   PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
109:   PetscCall(EPSDestroy(&eps));
110:   PetscCall(MatDestroy(&A));
111:   PetscCall(SlepcFinalize());
112:   return 0;
113: }

115: /*TEST

117:    testset:
118:       args: -eps_nev 4 -m 11 -eps_max_it 500
119:       output_file: output/test10_1.out
120:       test:
121:          suffix: 1
122:          args: -eps_type {{krylovschur arnoldi gd jd rqcg}}
123:       test:
124:          suffix: 1_lobpcg
125:          args: -eps_type lobpcg -eps_lobpcg_blocksize 6
126:       test:
127:          suffix: 1_lanczos
128:          args: -eps_type lanczos -eps_lanczos_reorthog local
129:          requires: !single
130:       test:
131:          suffix: 1_gd2
132:          args: -eps_type gd -eps_gd_double_expansion

134: TEST*/