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