Actual source code: ex3.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[] = "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;

 34:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 35:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 36:   if (size != 1) SETERRQ(PETSC_COMM_WORLD,1,"This is a uniprocessor example only");

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

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

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

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

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

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

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

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

 75:   EPSSolve(eps);

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

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

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

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

115:   dd  = 4.0;
116:   dl  = -1.0;
117:   du  = -1.0;

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

125: /*
126:     Matrix-vector product subroutine for the 2D Laplacian.

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

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

133:                  | T -I          |
134:                  |-I  T -I       |
135:              A = |   -I  T       |
136:                  |        ...  -I|
137:                  |           -I T|

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

150:   MatShellGetContext(A,&ctx);
151:   nx = *(int*)ctx;
152:   VecGetArrayRead(x,&px);
153:   VecGetArray(y,&py);

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

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

164:   lo = (nx-1)*nx;
165:   tv(nx,&px[lo],&py[lo]);
166:   for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i];

168:   VecRestoreArrayRead(x,&px);
169:   VecRestoreArray(y,&py);
170:   return(0);
171: }

173: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag)
174: {

178:   VecSet(diag,4.0);
179:   return(0);
180: }

182: /*TEST

184:    test:
185:       suffix: 1
186:       args: -n 72 -eps_nev 4 -eps_ncv 20 -terse
187:       requires: !single

189: TEST*/