Actual source code: ex9.c
slepc-3.9.0 2018-04-12
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Solves a problem associated to the Brusselator wave model in chemical reactions, illustrating the use of shell matrices.\n\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = block dimension of the 2x2 block matrix.\n"
14: " -L <L>, where <L> = bifurcation parameter.\n"
15: " -alpha <alpha>, -beta <beta>, -delta1 <delta1>, -delta2 <delta2>,\n"
16: " where <alpha> <beta> <delta1> <delta2> = model parameters.\n\n";
18: #include <slepceps.h>
20: /*
21: This example computes the eigenvalues with largest real part of the
22: following matrix
24: A = [ tau1*T+(beta-1)*I alpha^2*I
25: -beta*I tau2*T-alpha^2*I ],
27: where
29: T = tridiag{1,-2,1}
30: h = 1/(n+1)
31: tau1 = delta1/(h*L)^2
32: tau2 = delta2/(h*L)^2
33: */
36: /*
37: Matrix operations
38: */
39: PetscErrorCode MatMult_Brussel(Mat,Vec,Vec);
40: PetscErrorCode MatGetDiagonal_Brussel(Mat,Vec);
42: typedef struct {
43: Mat T;
44: Vec x1,x2,y1,y2;
45: PetscScalar alpha,beta,tau1,tau2,sigma;
46: } CTX_BRUSSEL;
48: int main(int argc,char **argv)
49: {
50: Mat A; /* eigenvalue problem matrix */
51: EPS eps; /* eigenproblem solver context */
52: EPSType type;
53: PetscScalar delta1,delta2,L,h;
54: PetscInt N=30,n,i,Istart,Iend,nev;
55: CTX_BRUSSEL *ctx;
56: PetscBool terse;
57: PetscViewer viewer;
60: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
62: PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
63: PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%D\n\n",N);
65: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66: Generate the matrix
67: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69: /*
70: Create shell matrix context and set default parameters
71: */
72: PetscNew(&ctx);
73: ctx->alpha = 2.0;
74: ctx->beta = 5.45;
75: delta1 = 0.008;
76: delta2 = 0.004;
77: L = 0.51302;
79: /*
80: Look the command line for user-provided parameters
81: */
82: PetscOptionsGetScalar(NULL,NULL,"-L",&L,NULL);
83: PetscOptionsGetScalar(NULL,NULL,"-alpha",&ctx->alpha,NULL);
84: PetscOptionsGetScalar(NULL,NULL,"-beta",&ctx->beta,NULL);
85: PetscOptionsGetScalar(NULL,NULL,"-delta1",&delta1,NULL);
86: PetscOptionsGetScalar(NULL,NULL,"-delta2",&delta2,NULL);
88: /*
89: Create matrix T
90: */
91: MatCreate(PETSC_COMM_WORLD,&ctx->T);
92: MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N);
93: MatSetFromOptions(ctx->T);
94: MatSetUp(ctx->T);
96: MatGetOwnershipRange(ctx->T,&Istart,&Iend);
97: for (i=Istart;i<Iend;i++) {
98: if (i>0) { MatSetValue(ctx->T,i,i-1,1.0,INSERT_VALUES); }
99: if (i<N-1) { MatSetValue(ctx->T,i,i+1,1.0,INSERT_VALUES); }
100: MatSetValue(ctx->T,i,i,-2.0,INSERT_VALUES);
101: }
102: MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY);
103: MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY);
104: MatGetLocalSize(ctx->T,&n,NULL);
106: /*
107: Fill the remaining information in the shell matrix context
108: and create auxiliary vectors
109: */
110: h = 1.0 / (PetscReal)(N+1);
111: ctx->tau1 = delta1 / ((h*L)*(h*L));
112: ctx->tau2 = delta2 / ((h*L)*(h*L));
113: ctx->sigma = 0.0;
114: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x1);
115: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x2);
116: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y1);
117: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y2);
119: /*
120: Create the shell matrix
121: */
122: MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A);
123: MatShellSetOperation(A,MATOP_MULT,(void(*)())MatMult_Brussel);
124: MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)())MatGetDiagonal_Brussel);
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Create the eigensolver and set various options
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: /*
131: Create eigensolver context
132: */
133: EPSCreate(PETSC_COMM_WORLD,&eps);
135: /*
136: Set operators. In this case, it is a standard eigenvalue problem
137: */
138: EPSSetOperators(eps,A,NULL);
139: EPSSetProblemType(eps,EPS_NHEP);
141: /*
142: Ask for the rightmost eigenvalues
143: */
144: EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
146: /*
147: Set other solver options at runtime
148: */
149: EPSSetFromOptions(eps);
151: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152: Solve the eigensystem
153: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155: EPSSolve(eps);
157: /*
158: Optional: Get some information from the solver and display it
159: */
160: EPSGetType(eps,&type);
161: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
162: EPSGetDimensions(eps,&nev,NULL,NULL);
163: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
165: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166: Display solution and clean up
167: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169: /* show detailed info unless -terse option is given by user */
170: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
171: if (terse) {
172: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
173: } else {
174: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
175: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
176: EPSReasonView(eps,viewer);
177: EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer);
178: PetscViewerPopFormat(viewer);
179: }
180: EPSDestroy(&eps);
181: MatDestroy(&A);
182: MatDestroy(&ctx->T);
183: VecDestroy(&ctx->x1);
184: VecDestroy(&ctx->x2);
185: VecDestroy(&ctx->y1);
186: VecDestroy(&ctx->y2);
187: PetscFree(ctx);
188: SlepcFinalize();
189: return ierr;
190: }
192: PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
193: {
194: PetscInt n;
195: const PetscScalar *px;
196: PetscScalar *py;
197: CTX_BRUSSEL *ctx;
198: PetscErrorCode ierr;
201: MatShellGetContext(A,(void**)&ctx);
202: MatGetLocalSize(ctx->T,&n,NULL);
203: VecGetArrayRead(x,&px);
204: VecGetArray(y,&py);
205: VecPlaceArray(ctx->x1,px);
206: VecPlaceArray(ctx->x2,px+n);
207: VecPlaceArray(ctx->y1,py);
208: VecPlaceArray(ctx->y2,py+n);
210: MatMult(ctx->T,ctx->x1,ctx->y1);
211: VecScale(ctx->y1,ctx->tau1);
212: VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1);
213: VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2);
215: MatMult(ctx->T,ctx->x2,ctx->y2);
216: VecScale(ctx->y2,ctx->tau2);
217: VecAXPY(ctx->y2,-ctx->beta,ctx->x1);
218: VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2);
220: VecRestoreArrayRead(x,&px);
221: VecRestoreArray(y,&py);
222: VecResetArray(ctx->x1);
223: VecResetArray(ctx->x2);
224: VecResetArray(ctx->y1);
225: VecResetArray(ctx->y2);
226: return(0);
227: }
229: PetscErrorCode MatGetDiagonal_Brussel(Mat A,Vec diag)
230: {
231: Vec d1,d2;
232: PetscInt n;
233: PetscScalar *pd;
234: MPI_Comm comm;
235: CTX_BRUSSEL *ctx;
239: MatShellGetContext(A,(void**)&ctx);
240: PetscObjectGetComm((PetscObject)A,&comm);
241: MatGetLocalSize(ctx->T,&n,NULL);
242: VecGetArray(diag,&pd);
243: VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd,&d1);
244: VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd+n,&d2);
246: VecSet(d1,-2.0*ctx->tau1 + ctx->beta - 1.0 + ctx->sigma);
247: VecSet(d2,-2.0*ctx->tau2 - ctx->alpha*ctx->alpha + ctx->sigma);
249: VecDestroy(&d1);
250: VecDestroy(&d2);
251: VecRestoreArray(diag,&pd);
252: return(0);
253: }