Subversion Repositories slepc-dev

Rev

Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1376 slepc 1
/*
2
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 3
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2116 eromero 4
   Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
6 dsic.upv.es!jroman 5
 
1672 slepc 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.
11
 
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.
16
 
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/>.
1376 slepc 19
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20
*/
21
 
228 dsic.upv.es!jroman 22
static char help[] = "Solves the same eigenproblem as in example ex2, but using a shell matrix. "
161 dsic.upv.es!antodo 23
  "The problem is a standard symmetric eigenproblem corresponding to the 2-D Laplacian operator.\n\n"
1156 slepc 24
  "The command line options are:\n"
6 dsic.upv.es!jroman 25
  "  -n <n>, where <n> = number of grid subdivisions in both x and y dimensions.\n\n";
26
 
27
#include "slepceps.h"
28
#include "petscblaslapack.h"
29
 
30
/*
31
   User-defined routines
32
*/
983 slepc 33
PetscErrorCode MatLaplacian2D_Mult( Mat A, Vec x, Vec y );
1732 antodo 34
PetscErrorCode MatLaplacian2D_GetDiagonal( Mat A, Vec diag );
6 dsic.upv.es!jroman 35
 
36
#undef __FUNCT__
37
#define __FUNCT__ "main"
38
int main( int argc, char **argv )
39
{
1505 slepc 40
  Mat            A;               /* operator matrix */
41
  EPS            eps;             /* eigenproblem solver context */
1507 slepc 42
  const EPSType  type;
1505 slepc 43
  PetscReal      error, tol, re, im;
44
  PetscScalar    kr, ki;
45
  PetscMPIInt    size;
983 slepc 46
  PetscErrorCode ierr;
1505 slepc 47
  PetscInt       N, n=10, nev, maxit, i, its, nconv;
6 dsic.upv.es!jroman 48
 
49
  SlepcInitialize(&argc,&argv,(char*)0,help);
50
  ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
51
  if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
52
 
53
  ierr = PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);CHKERRQ(ierr);
54
  N = n*n;
55
  ierr = PetscPrintf(PETSC_COMM_WORLD,"\n2-D Laplacian Eigenproblem (matrix-free version), N=%d (%dx%d grid)\n\n",N,n,n);CHKERRQ(ierr);
56
 
57
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58
     Compute the operator matrix that defines the eigensystem, Ax=kx
59
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60
 
61
  ierr = MatCreateShell(PETSC_COMM_WORLD,N,N,N,N,&n,&A);CHKERRQ(ierr);
62
  ierr = MatSetFromOptions(A);CHKERRQ(ierr);
63
  ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)())MatLaplacian2D_Mult);CHKERRQ(ierr);
58 dsic.upv.es!jroman 64
  ierr = MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)())MatLaplacian2D_Mult);CHKERRQ(ierr);
1732 antodo 65
  ierr = MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)())MatLaplacian2D_GetDiagonal);CHKERRQ(ierr);
6 dsic.upv.es!jroman 66
 
67
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68
                Create the eigensolver and set various options
69
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
70
 
71
  /*
72
     Create eigensolver context
73
  */
74
  ierr = EPSCreate(PETSC_COMM_WORLD,&eps);CHKERRQ(ierr);
75
 
76
  /*
77
     Set operators. In this case, it is a standard eigenvalue problem
78
  */
79
  ierr = EPSSetOperators(eps,A,PETSC_NULL);CHKERRQ(ierr);
408 dsic.upv.es!antodo 80
  ierr = EPSSetProblemType(eps,EPS_HEP);CHKERRQ(ierr);
6 dsic.upv.es!jroman 81
 
82
  /*
83
     Set solver parameters at runtime
84
  */
85
  ierr = EPSSetFromOptions(eps);CHKERRQ(ierr);
86
 
87
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88
                      Solve the eigensystem
89
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90
 
78 dsic.upv.es!antodo 91
  ierr = EPSSolve(eps);CHKERRQ(ierr);
92
  ierr = EPSGetIterationNumber(eps, &its);CHKERRQ(ierr);
6 dsic.upv.es!jroman 93
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);CHKERRQ(ierr);
94
 
95
  /*
96
     Optional: Get some information from the solver and display it
97
  */
98
  ierr = EPSGetType(eps,&type);CHKERRQ(ierr);
99
  ierr = PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);CHKERRQ(ierr);
1575 slepc 100
  ierr = EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
6 dsic.upv.es!jroman 101
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);CHKERRQ(ierr);
102
  ierr = EPSGetTolerances(eps,&tol,&maxit);CHKERRQ(ierr);
103
  ierr = PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);CHKERRQ(ierr);
104
 
105
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106
                    Display solution and clean up
107
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108
 
109
  /*
97 dsic.upv.es!antodo 110
     Get number of converged approximate eigenpairs
6 dsic.upv.es!jroman 111
  */
132 dsic.upv.es!antodo 112
  ierr = EPSGetConverged(eps,&nconv);CHKERRQ(ierr);
113
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %d\n\n",nconv);
97 dsic.upv.es!antodo 114
         CHKERRQ(ierr);
6 dsic.upv.es!jroman 115
 
132 dsic.upv.es!antodo 116
  if (nconv>0) {
6 dsic.upv.es!jroman 117
    /*
118
       Display eigenvalues and relative errors
119
    */
120
    ierr = PetscPrintf(PETSC_COMM_WORLD,
97 dsic.upv.es!antodo 121
         "           k          ||Ax-kx||/||kx||\n"
122
         "   ----------------- ------------------\n" );CHKERRQ(ierr);
123
 
132 dsic.upv.es!antodo 124
    for( i=0; i<nconv; i++ ) {
97 dsic.upv.es!antodo 125
      /*
126
        Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
127
        ki (imaginary part)
128
      */
132 dsic.upv.es!antodo 129
      ierr = EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
97 dsic.upv.es!antodo 130
      /*
131
         Compute the relative error associated to each eigenpair
132
      */
133
      ierr = EPSComputeRelativeError(eps,i,&error);CHKERRQ(ierr);
134
 
135
#ifdef PETSC_USE_COMPLEX
132 dsic.upv.es!antodo 136
      re = PetscRealPart(kr);
137
      im = PetscImaginaryPart(kr);
138
#else
139
      re = kr;
140
      im = ki;
97 dsic.upv.es!antodo 141
#endif
132 dsic.upv.es!antodo 142
      if (im!=0.0) {
1162 slepc 143
        ierr = PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g\n",re,im,error);CHKERRQ(ierr);
97 dsic.upv.es!antodo 144
      } else {
1162 slepc 145
        ierr = PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g\n",re,error);CHKERRQ(ierr);
97 dsic.upv.es!antodo 146
      }
6 dsic.upv.es!jroman 147
    }
148
    ierr = PetscPrintf(PETSC_COMM_WORLD,"\n" );CHKERRQ(ierr);
149
  }
150
 
151
  /*
152
     Free work space
153
  */
154
  ierr = EPSDestroy(eps);CHKERRQ(ierr);
155
  ierr = MatDestroy(A);CHKERRQ(ierr);
156
  ierr = SlepcFinalize();CHKERRQ(ierr);
157
  return 0;
158
}
159
 
160
/*
161
    Compute the matrix vector multiplication y<---T*x where T is a nx by nx
162
    tridiagonal matrix with DD on the diagonal, DL on the subdiagonal, and
163
    DU on the superdiagonal.
164
 */  
165
static void tv( int nx, PetscScalar *x, PetscScalar *y )
166
{
167
  PetscScalar dd, dl, du;
168
  int         j;
169
 
170
  dd  = 4.0;
171
  dl  = -1.0;
172
  du  = -1.0;
173
 
174
  y[0] =  dd*x[0] + du*x[1];
175
  for( j=1; j<nx-1; j++ )
176
    y[j] = dl*x[j-1] + dd*x[j] + du*x[j+1];
177
  y[nx-1] = dl*x[nx-2] + dd*x[nx-1];
178
}
179
 
180
#undef __FUNCT__
181
#define __FUNCT__ "MatLaplacian2D_Mult"
182
/*
183
    Matrix-vector product subroutine for the 2D Laplacian.
184
 
185
    The matrix used is the 2 dimensional discrete Laplacian on unit square with
186
    zero Dirichlet boundary condition.
187
 
188
    Computes y <-- A*x, where A is the block tridiagonal matrix
189
 
190
                 | T -I          |
191
                 |-I  T -I       |
192
             A = |   -I  T       |
193
                 |        ...  -I|
194
                 |           -I T|
195
 
196
    The subroutine TV is called to compute y<--T*x.
197
 */
983 slepc 198
PetscErrorCode MatLaplacian2D_Mult( Mat A, Vec x, Vec y )
6 dsic.upv.es!jroman 199
{
983 slepc 200
  void           *ctx;
201
  PetscErrorCode ierr;
202
  int            nx, lo, j, one=1;
203
  PetscScalar    *px, *py, dmone=-1.0;
6 dsic.upv.es!jroman 204
 
983 slepc 205
  PetscFunctionBegin;
6 dsic.upv.es!jroman 206
  ierr = MatShellGetContext( A, &ctx ); CHKERRQ(ierr);
207
  nx = *(int *)ctx;
208
  ierr = VecGetArray( x, &px ); CHKERRQ(ierr);
209
  ierr = VecGetArray( y, &py ); CHKERRQ(ierr);
210
 
211
  tv( nx, &px[0], &py[0] );
801 dsic.upv.es!antodo 212
  BLASaxpy_( &nx, &dmone, &px[nx], &one, &py[0], &one );
6 dsic.upv.es!jroman 213
 
214
  for( j=2; j<nx; j++ ) {
215
    lo = (j-1)*nx;
216
    tv( nx, &px[lo], &py[lo]);
801 dsic.upv.es!antodo 217
    BLASaxpy_( &nx, &dmone, &px[lo-nx], &one, &py[lo], &one );
218
    BLASaxpy_( &nx, &dmone, &px[lo+nx], &one, &py[lo], &one );
6 dsic.upv.es!jroman 219
  }
220
 
221
  lo = (nx-1)*nx;
222
  tv( nx, &px[lo], &py[lo]);
801 dsic.upv.es!antodo 223
  BLASaxpy_( &nx, &dmone, &px[lo-nx], &one, &py[lo], &one );
6 dsic.upv.es!jroman 224
 
225
  ierr = VecRestoreArray( x, &px ); CHKERRQ(ierr);
226
  ierr = VecRestoreArray( y, &py ); CHKERRQ(ierr);
983 slepc 227
  PetscFunctionReturn(0);
6 dsic.upv.es!jroman 228
}
229
 
1732 antodo 230
#undef __FUNCT__
231
#define __FUNCT__ "MatLaplacian2D_GetDiagonal"
232
PetscErrorCode MatLaplacian2D_GetDiagonal( Mat A, Vec diag )
233
{
234
  PetscErrorCode ierr;
235
 
236
  PetscFunctionBegin;
237
  ierr = VecSet(diag,4.0);CHKERRQ(ierr);
238
  PetscFunctionReturn(0);
239
}
240
 
241