Actual source code: test8.c
slepc-3.7.2 2016-07-19
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
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 MatMult_Laplacian2D(Mat A,Vec x,Vec y);
34: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag);
38: int main(int argc,char **argv)
39: {
40: Mat A; /* operator matrix */
41: EPS eps; /* eigenproblem solver context */
42: EPSType type;
43: PetscReal tol=1000*PETSC_MACHINE_EPSILON;
44: PetscMPIInt size;
45: PetscInt N,n=10,nev;
48: SlepcInitialize(&argc,&argv,(char*)0,help);
49: MPI_Comm_size(PETSC_COMM_WORLD,&size);
50: if (size != 1) SETERRQ(PETSC_COMM_WORLD,1,"This is a uniprocessor example only");
52: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
53: N = n*n;
54: PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem (matrix-free version), N=%D (%Dx%D grid)\n\n",N,n,n);
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Compute the operator matrix that defines the eigensystem, Ax=kx
58: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60: MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,&n,&A);
61: MatSetFromOptions(A);
62: MatShellSetOperation(A,MATOP_MULT,(void(*)())MatMult_Laplacian2D);
63: MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)())MatMult_Laplacian2D);
64: MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)())MatGetDiagonal_Laplacian2D);
66: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: Create the eigensolver and set various options
68: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
70: /*
71: Create eigensolver context
72: */
73: EPSCreate(PETSC_COMM_WORLD,&eps);
75: /*
76: Set operators. In this case, it is a standard eigenvalue problem
77: */
78: EPSSetOperators(eps,A,NULL);
79: EPSSetProblemType(eps,EPS_HEP);
80: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
82: /*
83: Set solver parameters at runtime
84: */
85: EPSSetFromOptions(eps);
87: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: Solve the eigensystem
89: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91: EPSSolve(eps);
93: /*
94: Optional: Get some information from the solver and display it
95: */
96: EPSGetType(eps,&type);
97: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
98: EPSGetDimensions(eps,&nev,NULL,NULL);
99: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
101: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: Display solution and clean up
103: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
106: EPSDestroy(&eps);
107: MatDestroy(&A);
108: SlepcFinalize();
109: return ierr;
110: }
112: /*
113: Compute the matrix vector multiplication y<---T*x where T is a nx by nx
114: tridiagonal matrix with DD on the diagonal, DL on the subdiagonal, and
115: DU on the superdiagonal.
116: */
117: static void tv(int nx,const PetscScalar *x,PetscScalar *y)
118: {
119: PetscScalar dd,dl,du;
120: int j;
122: dd = 4.0;
123: dl = -1.0;
124: du = -1.0;
126: y[0] = dd*x[0] + du*x[1];
127: for (j=1;j<nx-1;j++)
128: y[j] = dl*x[j-1] + dd*x[j] + du*x[j+1];
129: y[nx-1] = dl*x[nx-2] + dd*x[nx-1];
130: }
134: /*
135: Matrix-vector product subroutine for the 2D Laplacian.
137: The matrix used is the 2 dimensional discrete Laplacian on unit square with
138: zero Dirichlet boundary condition.
140: Computes y <-- A*x, where A is the block tridiagonal matrix
142: | T -I |
143: |-I T -I |
144: A = | -I T |
145: | ... -I|
146: | -I T|
148: The subroutine TV is called to compute y<--T*x.
149: */
150: PetscErrorCode MatMult_Laplacian2D(Mat A,Vec x,Vec y)
151: {
152: void *ctx;
153: int nx,lo,i,j;
154: const PetscScalar *px;
155: PetscScalar *py;
156: PetscErrorCode ierr;
159: MatShellGetContext(A,&ctx);
160: nx = *(int*)ctx;
161: VecGetArrayRead(x,&px);
162: VecGetArray(y,&py);
164: tv(nx,&px[0],&py[0]);
165: for (i=0;i<nx;i++) py[i] -= px[nx+i];
167: for (j=2;j<nx;j++) {
168: lo = (j-1)*nx;
169: tv(nx,&px[lo],&py[lo]);
170: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i] + px[lo+nx+i];
171: }
173: lo = (nx-1)*nx;
174: tv(nx,&px[lo],&py[lo]);
175: for (i=0;i<nx;i++) py[lo+i] -= px[lo-nx+i];
177: VecRestoreArrayRead(x,&px);
178: VecRestoreArray(y,&py);
179: return(0);
180: }
184: PetscErrorCode MatGetDiagonal_Laplacian2D(Mat A,Vec diag)
185: {
189: VecSet(diag,4.0);
190: return(0);
191: }