Actual source code: ex3.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[] = "Solves the same eigenproblem as in example ex2, but using a shell matrix. "
 12:   "The problem is a standard symmetric eigenproblem corresponding to the 2-D Laplacian operator.\n\n"
 13:   "The command line options are:\n"
 14:   "  -n <n>, where <n> = number of grid subdivisions in both x and y dimensions.\n\n";

 16: #include <slepceps.h>

 18: /*
 19:    User-defined routines
 20: */
 21: PetscErrorCode MatMult_Laplacian2D(Mat A,Vec x,Vec y);
 22: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag);

 24: int main(int argc,char **argv)
 25: {
 26:   Mat            A;               /* operator matrix */
 27:   EPS            eps;             /* eigenproblem solver context */
 28:   EPSType        type;
 29:   PetscMPIInt    size;
 30:   PetscInt       N,n=10,nev;
 31:   PetscBool      terse;

 33:   SlepcInitialize(&argc,&argv,(char*)0,help);
 34:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 37:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 38:   N = n*n;
 39:   PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem (matrix-free version), N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,n);

 41:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 42:        Create the operator matrix that defines the eigensystem, Ax=kx
 43:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 45:   MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,&n,&A);
 46:   MatShellSetOperation(A,MATOP_MULT,(void(*)(void))MatMult_Laplacian2D);
 47:   MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_Laplacian2D);
 48:   MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Laplacian2D);

 50:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 51:                 Create the eigensolver and set various options
 52:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 54:   /*
 55:      Create eigensolver context
 56:   */
 57:   EPSCreate(PETSC_COMM_WORLD,&eps);

 59:   /*
 60:      Set operators. In this case, it is a standard eigenvalue problem
 61:   */
 62:   EPSSetOperators(eps,A,NULL);
 63:   EPSSetProblemType(eps,EPS_HEP);

 65:   /*
 66:      Set solver parameters at runtime
 67:   */
 68:   EPSSetFromOptions(eps);

 70:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 71:                       Solve the eigensystem
 72:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 74:   EPSSolve(eps);

 76:   /*
 77:      Optional: Get some information from the solver and display it
 78:   */
 79:   EPSGetType(eps,&type);
 80:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
 81:   EPSGetDimensions(eps,&nev,NULL,NULL);
 82:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);

 84:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 85:                     Display solution and clean up
 86:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 88:   /* show detailed info unless -terse option is given by user */
 89:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
 90:   if (terse) EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
 91:   else {
 92:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
 93:     EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
 94:     EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
 95:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
 96:   }
 97:   EPSDestroy(&eps);
 98:   MatDestroy(&A);
 99:   SlepcFinalize();
100:   return 0;
101: }

103: /*
104:     Compute the matrix vector multiplication y<---T*x where T is a nx by nx
105:     tridiagonal matrix with DD on the diagonal, DL on the subdiagonal, and
106:     DU on the superdiagonal.
107:  */
108: static void tv(int nx,const PetscScalar *x,PetscScalar *y)
109: {
110:   PetscScalar dd,dl,du;
111:   int         j;

113:   dd  = 4.0;
114:   dl  = -1.0;
115:   du  = -1.0;

117:   y[0] =  dd*x[0] + du*x[1];
118:   for (j=1;j<nx-1;j++)
119:     y[j] = dl*x[j-1] + dd*x[j] + du*x[j+1];
120:   y[nx-1] = dl*x[nx-2] + dd*x[nx-1];
121: }

123: /*
124:     Matrix-vector product subroutine for the 2D Laplacian.

126:     The matrix used is the 2 dimensional discrete Laplacian on unit square with
127:     zero Dirichlet boundary condition.

129:     Computes y <-- A*x, where A is the block tridiagonal matrix

131:                  | T -I          |
132:                  |-I  T -I       |
133:              A = |   -I  T       |
134:                  |        ...  -I|
135:                  |           -I T|

137:     The subroutine TV is called to compute y<--T*x.
138:  */
139: PetscErrorCode MatMult_Laplacian2D(Mat A,Vec x,Vec y)
140: {
141:   void              *ctx;
142:   int               nx,lo,i,j;
143:   const PetscScalar *px;
144:   PetscScalar       *py;

147:   MatShellGetContext(A,&ctx);
148:   nx = *(int*)ctx;
149:   VecGetArrayRead(x,&px);
150:   VecGetArray(y,&py);

152:   tv(nx,&px[0],&py[0]);
153:   for (i=0;i<nx;i++) py[i] -= px[nx+i];

155:   for (j=2;j<nx;j++) {
156:     lo = (j-1)*nx;
157:     tv(nx,&px[lo],&py[lo]);
158:     for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i] + px[lo+nx+i];
159:   }

161:   lo = (nx-1)*nx;
162:   tv(nx,&px[lo],&py[lo]);
163:   for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i];

165:   VecRestoreArrayRead(x,&px);
166:   VecRestoreArray(y,&py);
167:   PetscFunctionReturn(0);
168: }

170: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag)
171: {
173:   VecSet(diag,4.0);
174:   PetscFunctionReturn(0);
175: }

177: /*TEST

179:    test:
180:       suffix: 1
181:       args: -n 72 -eps_nev 4 -eps_ncv 20 -terse
182:       requires: !single

184: TEST*/