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
6 dsic.upv.es!jroman 1
 
228 dsic.upv.es!jroman 2
static char help[] = "Solves a problem associated to the Brusselator wave model in chemical reactions, illustrating the use of shell matrices.\n\n"
6 dsic.upv.es!jroman 3
  "The command line options are:\n\n"
4
  "  -n <n>, where <n> = block dimension of the 2x2 block matrix.\n"
5
  "  -L <L>, where <L> = bifurcation parameter.\n"
6
  "  -alpha <alpha>, -beta <beta>, -delta1 <delta1>,  -delta2 <delta2>,\n"
7
  "       where <alpha> <beta> <delta1> <delta2> = model parameters.\n\n";
8
 
9
#include "slepceps.h"
10
 
11
/*
12
   This example computes the eigenvalues with largest real part of the
13
   following matrix
14
 
15
        A = [ tau1*T+(beta-1)*I     alpha^2*I
16
                  -beta*I        tau2*T-alpha^2*I ],
17
 
18
   where
19
 
20
        T = tridiag{1,-2,1}
21
        h = 1/(n+1)
22
        tau1 = delta1/(h*L)^2
23
        tau2 = delta2/(h*L)^2
24
 */
25
 
26
 
27
/*
28
   Matrix operations
29
*/
30
extern int MatBrussel_Mult(Mat,Vec,Vec);
31
extern int MatBrussel_Shift(PetscScalar*,Mat);
32
extern int MatBrussel_GetDiagonal(Mat,Vec);
33
 
34
typedef struct {
35
  Mat         T;
36
  Vec         x1, x2, y1, y2;
37
  PetscScalar alpha, beta, tau1, tau2, sigma;
38
} CTX_BRUSSEL;
39
 
40
#undef __FUNCT__
41
#define __FUNCT__ "main"
42
int main( int argc, char **argv )
43
{
44
  Mat         A;               /* eigenvalue problem matrix */
45
  EPS         eps;             /* eigenproblem solver context */
46
  EPSType     type;
97 dsic.upv.es!antodo 47
  PetscReal   error, tol, re, im;
48
  PetscScalar delta1, delta2, L, h, kr, ki, value[3];
982 slepc 49
  PetscInt    N=30, n, i, col[3], Istart, Iend, FirstBlock=0, LastBlock=0;
50
  int         nev, ierr, maxit, its, nconv;
6 dsic.upv.es!jroman 51
  CTX_BRUSSEL *ctx;
52
 
53
  SlepcInitialize(&argc,&argv,(char*)0,help);
54
 
55
  ierr = PetscOptionsGetInt(PETSC_NULL,"-n",&N,PETSC_NULL);CHKERRQ(ierr);
56
  ierr = PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%d\n\n",N);CHKERRQ(ierr);
57
 
58
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59
        Generate the matrix
60
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61
 
62
  /*
63
     Create shell matrix context and set default parameters
64
  */
65
  ierr = PetscNew(CTX_BRUSSEL,&ctx);CHKERRQ(ierr);
66
  ctx->alpha = 2.0;
67
  ctx->beta  = 5.45;
68
  delta1     = 0.008;
69
  delta2     = 0.004;
70
  L          = 0.51302;
71
 
72
  /*
73
     Look the command line for user-provided parameters
74
  */
75
  ierr = PetscOptionsGetScalar(PETSC_NULL,"-L",&L,PETSC_NULL);CHKERRQ(ierr);
76
  ierr = PetscOptionsGetScalar(PETSC_NULL,"-alpha",&ctx->alpha,PETSC_NULL);CHKERRQ(ierr);
77
  ierr = PetscOptionsGetScalar(PETSC_NULL,"-beta",&ctx->beta,PETSC_NULL);CHKERRQ(ierr);
78
  ierr = PetscOptionsGetScalar(PETSC_NULL,"-delta1",&delta1,PETSC_NULL);CHKERRQ(ierr);
79
  ierr = PetscOptionsGetScalar(PETSC_NULL,"-delta2",&delta2,PETSC_NULL);CHKERRQ(ierr);
80
 
81
  /*
82
     Create matrix T
83
  */
828 dsic.upv.es!antodo 84
  ierr = MatCreate(PETSC_COMM_WORLD,&ctx->T);CHKERRQ(ierr);
85
  ierr = MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr);
6 dsic.upv.es!jroman 86
  ierr = MatSetFromOptions(ctx->T);CHKERRQ(ierr);
87
 
88
  ierr = MatGetOwnershipRange(ctx->T,&Istart,&Iend);CHKERRQ(ierr);
89
  if (Istart==0) FirstBlock=PETSC_TRUE;
90
  if (Iend==N) LastBlock=PETSC_TRUE;
91
  value[0]=1.0; value[1]=-2.0; value[2]=1.0;
92
  for( i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++ ) {
93
    col[0]=i-1; col[1]=i; col[2]=i+1;
94
    ierr = MatSetValues(ctx->T,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
95
  }
96
  if (LastBlock) {
97
    i=N-1; col[0]=N-2; col[1]=N-1;
98
    ierr = MatSetValues(ctx->T,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
99
  }
100
  if (FirstBlock) {
101
    i=0; col[0]=0; col[1]=1; value[0]=-2.0; value[1]=1.0;
102
    ierr = MatSetValues(ctx->T,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
103
  }
104
 
105
  ierr = MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
106
  ierr = MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
107
  ierr = MatGetLocalSize(ctx->T,&n,PETSC_NULL);CHKERRQ(ierr);
108
 
109
  /*
110
     Fill the remaining information in the shell matrix context
111
     and create auxiliary vectors
112
  */
113
  h = 1.0 / (double)(N+1);
114
  ctx->tau1 = delta1 / ((h*L)*(h*L));
115
  ctx->tau2 = delta2 / ((h*L)*(h*L));
116
  ctx->sigma = 0.0;
117
  ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->x1);CHKERRQ(ierr);
118
  ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->x2);CHKERRQ(ierr);
119
  ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->y1);CHKERRQ(ierr);
120
  ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n,PETSC_DECIDE,PETSC_NULL,&ctx->y2);CHKERRQ(ierr);
121
 
122
  /*
123
     Create the shell matrix
124
  */
125
  ierr = MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A);CHKERRQ(ierr);
126
  ierr = MatShellSetOperation(A,MATOP_MULT,(void(*)())MatBrussel_Mult);CHKERRQ(ierr);
127
  ierr = MatShellSetOperation(A,MATOP_SHIFT,(void(*)())MatBrussel_Shift);CHKERRQ(ierr);
128
  ierr = MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)())MatBrussel_GetDiagonal);CHKERRQ(ierr);
129
 
130
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131
                Create the eigensolver and set various options
132
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133
 
134
  /*
135
     Create eigensolver context
136
  */
137
  ierr = EPSCreate(PETSC_COMM_WORLD,&eps);CHKERRQ(ierr);
138
 
139
  /*
140
     Set operators. In this case, it is a standard eigenvalue problem
141
  */
142
  ierr = EPSSetOperators(eps,A,PETSC_NULL);CHKERRQ(ierr);
408 dsic.upv.es!antodo 143
  ierr = EPSSetProblemType(eps,EPS_NHEP);CHKERRQ(ierr);
6 dsic.upv.es!jroman 144
 
145
  /*
146
     Force to use ARPACK if it is installed and ask for the rightmost eigenvalues
147
  */
148
#if defined(SLEPC_HAVE_ARPACK)
149
  ierr = EPSSetType(eps,EPSARPACK);CHKERRQ(ierr);
147 dsic.upv.es!antodo 150
#else
151
  ierr = EPSSetType(eps, EPSLAPACK);CHKERRQ(ierr);
152
#endif
6 dsic.upv.es!jroman 153
  ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);CHKERRQ(ierr);
154
 
155
  /*
156
     Set other solver options at runtime
157
  */
158
  ierr = EPSSetFromOptions(eps);CHKERRQ(ierr);
159
 
160
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161
                      Solve the eigensystem
162
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163
 
78 dsic.upv.es!antodo 164
  ierr = EPSSolve(eps);CHKERRQ(ierr);
165
  ierr = EPSGetIterationNumber(eps, &its);CHKERRQ(ierr);
6 dsic.upv.es!jroman 166
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);CHKERRQ(ierr);
167
 
168
  /*
169
     Optional: Get some information from the solver and display it
170
  */
171
  ierr = EPSGetType(eps,&type);CHKERRQ(ierr);
172
  ierr = PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);CHKERRQ(ierr);
173
  ierr = EPSGetDimensions(eps,&nev,PETSC_NULL);CHKERRQ(ierr);
174
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);CHKERRQ(ierr);
175
  ierr = EPSGetTolerances(eps,&tol,&maxit);CHKERRQ(ierr);
176
  ierr = PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);CHKERRQ(ierr);
177
 
178
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179
                    Display solution and clean up
180
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181
 
182
  /*
183
     Get number of converged eigenpairs
184
  */
132 dsic.upv.es!antodo 185
  ierr = EPSGetConverged(eps,&nconv);CHKERRQ(ierr);
186
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);CHKERRQ(ierr);
6 dsic.upv.es!jroman 187
 
132 dsic.upv.es!antodo 188
  if (nconv>0) {
6 dsic.upv.es!jroman 189
    /*
190
       Display eigenvalues and relative errors
191
    */
192
    ierr = PetscPrintf(PETSC_COMM_WORLD,
97 dsic.upv.es!antodo 193
         "           k             ||Ax-kx||/||kx||\n"
194
         "  --------------------- ------------------\n" );CHKERRQ(ierr);
132 dsic.upv.es!antodo 195
    for( i=0; i<nconv; i++ ) {
97 dsic.upv.es!antodo 196
      /*
197
         Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
198
         ki (imaginary part)
199
      */
200
      ierr = EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
201
 
202
      /*
203
         Compute the relative error associated to each eigenpair
204
      */
205
      ierr = EPSComputeRelativeError(eps,i,&error);CHKERRQ(ierr);
206
 
6 dsic.upv.es!jroman 207
#if defined(PETSC_USE_COMPLEX)
97 dsic.upv.es!antodo 208
      re = PetscRealPart(kr);
209
      im = PetscImaginaryPart(kr);
6 dsic.upv.es!jroman 210
#else
97 dsic.upv.es!antodo 211
      re = kr;
212
      im = ki;
6 dsic.upv.es!jroman 213
#endif
97 dsic.upv.es!antodo 214
      if( im != 0.0 ) {
215
        ierr = PetscPrintf(PETSC_COMM_WORLD," % 6f %+6f i",re,im);CHKERRQ(ierr);
216
      } else {
132 dsic.upv.es!antodo 217
        ierr = PetscPrintf(PETSC_COMM_WORLD,"       % 6f      ",re); CHKERRQ(ierr);
97 dsic.upv.es!antodo 218
      }
219
      ierr = PetscPrintf(PETSC_COMM_WORLD," % 12f\n",error);CHKERRQ(ierr);
6 dsic.upv.es!jroman 220
    }
221
    ierr = PetscPrintf(PETSC_COMM_WORLD,"\n" );CHKERRQ(ierr);
222
  }
223
 
224
  /*
225
     Free work space
226
  */
227
  ierr = EPSDestroy(eps);CHKERRQ(ierr);
228
  ierr = MatDestroy(A);CHKERRQ(ierr);
229
  ierr = MatDestroy(ctx->T);CHKERRQ(ierr);
230
  ierr = VecDestroy(ctx->x1);CHKERRQ(ierr);
231
  ierr = VecDestroy(ctx->x2);CHKERRQ(ierr);
232
  ierr = VecDestroy(ctx->y1);CHKERRQ(ierr);
233
  ierr = VecDestroy(ctx->y2);CHKERRQ(ierr);
234
  ierr = PetscFree(ctx);CHKERRQ(ierr);
235
  ierr = SlepcFinalize();CHKERRQ(ierr);
236
  return 0;
237
}
238
 
239
#undef __FUNC__
240
#define __FUNC__ "MatBrussel_Mult"
241
int MatBrussel_Mult(Mat A,Vec x,Vec y)
242
{
982 slepc 243
  int         ierr;
244
  PetscInt    n;
828 dsic.upv.es!antodo 245
  PetscScalar *px, *py;
6 dsic.upv.es!jroman 246
  CTX_BRUSSEL *ctx;
247
 
248
  ierr = MatShellGetContext(A,(void**)&ctx);CHKERRQ(ierr);
249
  ierr = MatGetLocalSize(ctx->T,&n,PETSC_NULL);CHKERRQ(ierr);
250
  ierr = VecGetArray(x,&px);CHKERRQ(ierr);
251
  ierr = VecGetArray(y,&py);CHKERRQ(ierr);
252
  ierr = VecPlaceArray(ctx->x1,px);CHKERRQ(ierr);
253
  ierr = VecPlaceArray(ctx->x2,px+n);CHKERRQ(ierr);
254
  ierr = VecPlaceArray(ctx->y1,py);CHKERRQ(ierr);
255
  ierr = VecPlaceArray(ctx->y2,py+n);CHKERRQ(ierr);
256
 
257
  ierr = MatMult(ctx->T,ctx->x1,ctx->y1);CHKERRQ(ierr);
828 dsic.upv.es!antodo 258
  ierr = VecScale(ctx->y1,ctx->tau1);CHKERRQ(ierr);
259
  ierr = VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1);CHKERRQ(ierr);
260
  ierr = VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2);CHKERRQ(ierr);
6 dsic.upv.es!jroman 261
 
262
  ierr = MatMult(ctx->T,ctx->x2,ctx->y2);CHKERRQ(ierr);
828 dsic.upv.es!antodo 263
  ierr = VecScale(ctx->y2,ctx->tau2);CHKERRQ(ierr);
264
  ierr = VecAXPY(ctx->y2,-ctx->beta,ctx->x1);CHKERRQ(ierr);
265
  ierr = VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2);CHKERRQ(ierr);
6 dsic.upv.es!jroman 266
 
267
  ierr = VecRestoreArray(x,&px);CHKERRQ(ierr);
268
  ierr = VecRestoreArray(y,&py);CHKERRQ(ierr);
864 dsic.upv.es!jroman 269
  ierr = VecResetArray(ctx->x1);CHKERRQ(ierr);
270
  ierr = VecResetArray(ctx->x2);CHKERRQ(ierr);
271
  ierr = VecResetArray(ctx->y1);CHKERRQ(ierr);
272
  ierr = VecResetArray(ctx->y2);CHKERRQ(ierr);
6 dsic.upv.es!jroman 273
 
274
  return 0;
275
}
276
 
277
#undef __FUNCT__
278
#define __FUNCT__ "MatBrussel_Shift"
279
int MatBrussel_Shift( PetscScalar* a, Mat Y )
280
{
281
  CTX_BRUSSEL *ctx;
282
  int         ierr;
283
 
284
  ierr = MatShellGetContext( Y, (void**)&ctx ); CHKERRQ(ierr);
285
  ctx->sigma += *a;
286
  PetscFunctionReturn(0);
287
}
288
 
289
#undef __FUNC__
290
#define __FUNC__ "MatBrussel_GetDiagonal"
291
int MatBrussel_GetDiagonal(Mat A,Vec diag)
292
{
293
  Vec         d1, d2;
982 slepc 294
  int         ierr;
295
  PetscInt    n;
828 dsic.upv.es!antodo 296
  PetscScalar *pd;
6 dsic.upv.es!jroman 297
  MPI_Comm    comm;
298
  CTX_BRUSSEL *ctx;
299
 
300
  ierr = MatShellGetContext(A,(void**)&ctx);CHKERRQ(ierr);
301
  ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
302
  ierr = MatGetLocalSize(ctx->T,&n,PETSC_NULL);CHKERRQ(ierr);
303
  ierr = VecGetArray(diag,&pd);CHKERRQ(ierr);
304
  ierr = VecCreateMPIWithArray(comm,n,PETSC_DECIDE,pd,&d1);CHKERRQ(ierr);
305
  ierr = VecCreateMPIWithArray(comm,n,PETSC_DECIDE,pd+n,&d2);CHKERRQ(ierr);
306
 
828 dsic.upv.es!antodo 307
  ierr = VecSet(d1,-2.0*ctx->tau1 + ctx->beta - 1.0 + ctx->sigma);CHKERRQ(ierr);
308
  ierr = VecSet(d2,-2.0*ctx->tau2 - ctx->alpha*ctx->alpha + ctx->sigma);CHKERRQ(ierr);
6 dsic.upv.es!jroman 309
 
310
  ierr = VecDestroy(d1);CHKERRQ(ierr);
311
  ierr = VecDestroy(d2);CHKERRQ(ierr);
312
  ierr = VecRestoreArray(diag,&pd);CHKERRQ(ierr);
313
 
314
  return 0;
315
}
316
 
317