Actual source code: ex3.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:       
  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: static char help[] = "Solves the same eigenproblem as in example ex2, but using a shell matrix. "
 23:   "The problem is a standard symmetric eigenproblem corresponding to the 2-D Laplacian operator.\n\n"
 24:   "The command line options are:\n"
 25:   "  -n <n>, where <n> = number of grid subdivisions in both x and y dimensions.\n\n";

 27: #include <slepceps.h>
 28: #include <petscblaslapack.h>

 30: /* 
 31:    User-defined routines
 32: */
 33: PetscErrorCode MatLaplacian2D_Mult(Mat A,Vec x,Vec y);
 34: PetscErrorCode MatLaplacian2D_GetDiagonal(Mat A,Vec diag);

 38: int main(int argc,char **argv)
 39: {
 40:   Mat            A;               /* operator matrix */
 41:   EPS            eps;             /* eigenproblem solver context */
 42:   const EPSType  type;
 43:   PetscMPIInt    size;
 44:   PetscInt       N,n=10,nev;

 47:   SlepcInitialize(&argc,&argv,(char*)0,help);
 48:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 49:   if (size != 1) SETERRQ(PETSC_COMM_WORLD,1,"This is a uniprocessor example only");

 51:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 52:   N = n*n;
 53:   PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem (matrix-free version), N=%D (%Dx%D grid)\n\n",N,n,n);

 55:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 56:      Compute the operator matrix that defines the eigensystem, Ax=kx
 57:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 59:   MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,&n,&A);
 60:   MatSetFromOptions(A);
 61:   MatShellSetOperation(A,MATOP_MULT,(void(*)())MatLaplacian2D_Mult);
 62:   MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)())MatLaplacian2D_Mult);
 63:   MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)())MatLaplacian2D_GetDiagonal);

 65:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 66:                 Create the eigensolver and set various options
 67:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

 80:   /*
 81:      Set solver parameters at runtime
 82:   */
 83:   EPSSetFromOptions(eps);

 85:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 86:                       Solve the eigensystem
 87:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 89:   EPSSolve(eps);

 91:   /*
 92:      Optional: Get some information from the solver and display it
 93:   */
 94:   EPSGetType(eps,&type);
 95:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
 96:   EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
 97:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

 99:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
100:                     Display solution and clean up
101:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

103:   EPSPrintSolution(eps,PETSC_NULL);
104:   EPSDestroy(&eps);
105:   MatDestroy(&A);
106:   SlepcFinalize();
107:   return 0;
108: }

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

120:   dd  = 4.0;
121:   dl  = -1.0;
122:   du  = -1.0;

124:   y[0] =  dd*x[0] + du*x[1];
125:   for (j=1;j<nx-1;j++)
126:     y[j] = dl*x[j-1] + dd*x[j] + du*x[j+1];
127:   y[nx-1] = dl*x[nx-2] + dd*x[nx-1];
128: }

132: /*
133:     Matrix-vector product subroutine for the 2D Laplacian.

135:     The matrix used is the 2 dimensional discrete Laplacian on unit square with
136:     zero Dirichlet boundary condition.
137:  
138:     Computes y <-- A*x, where A is the block tridiagonal matrix
139:  
140:                  | T -I          | 
141:                  |-I  T -I       |
142:              A = |   -I  T       |
143:                  |        ...  -I|
144:                  |           -I T|
145:  
146:     The subroutine TV is called to compute y<--T*x.
147:  */
148: PetscErrorCode MatLaplacian2D_Mult(Mat A,Vec x,Vec y)
149: {
150:   void              *ctx;
151:   int               nx,lo,j,one=1;
152:   const PetscScalar *px;
153:   PetscScalar       *py,dmone=-1.0;
154:   PetscErrorCode    ierr;
155: 
157:   MatShellGetContext(A,&ctx);
158:   nx = *(int*)ctx;
159:   VecGetArrayRead(x,&px);
160:   VecGetArray(y,&py);

162:   tv(nx,&px[0],&py[0]);
163:   BLASaxpy_(&nx,&dmone,&px[nx],&one,&py[0],&one);

165:   for (j=2;j<nx;j++) {
166:     lo = (j-1)*nx;
167:     tv(nx,&px[lo],&py[lo]);
168:     BLASaxpy_(&nx,&dmone,&px[lo-nx],&one,&py[lo],&one);
169:     BLASaxpy_(&nx,&dmone,&px[lo+nx],&one,&py[lo],&one);
170:   }

172:   lo = (nx-1)*nx;
173:   tv(nx,&px[lo],&py[lo]);
174:   BLASaxpy_(&nx,&dmone,&px[lo-nx],&one,&py[lo],&one);

176:   VecRestoreArrayRead(x,&px);
177:   VecRestoreArray(y,&py);
178:   return(0);
179: }

183: PetscErrorCode MatLaplacian2D_GetDiagonal(Mat A,Vec diag)
184: {

188:   VecSet(diag,4.0);
189:   return(0);
190: }