Actual source code: ex9.c
slepc-3.7.1 2016-05-27
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 a problem associated to the Brusselator wave model in chemical reactions, illustrating the use of shell matrices.\n\n"
23: "The command line options are:\n"
24: " -n <n>, where <n> = block dimension of the 2x2 block matrix.\n"
25: " -L <L>, where <L> = bifurcation parameter.\n"
26: " -alpha <alpha>, -beta <beta>, -delta1 <delta1>, -delta2 <delta2>,\n"
27: " where <alpha> <beta> <delta1> <delta2> = model parameters.\n\n";
29: #include <slepceps.h>
31: /*
32: This example computes the eigenvalues with largest real part of the
33: following matrix
35: A = [ tau1*T+(beta-1)*I alpha^2*I
36: -beta*I tau2*T-alpha^2*I ],
38: where
40: T = tridiag{1,-2,1}
41: h = 1/(n+1)
42: tau1 = delta1/(h*L)^2
43: tau2 = delta2/(h*L)^2
44: */
47: /*
48: Matrix operations
49: */
50: PetscErrorCode MatMult_Brussel(Mat,Vec,Vec);
51: PetscErrorCode MatShift_Brussel(PetscScalar*,Mat);
52: PetscErrorCode MatGetDiagonal_Brussel(Mat,Vec);
54: typedef struct {
55: Mat T;
56: Vec x1,x2,y1,y2;
57: PetscScalar alpha,beta,tau1,tau2,sigma;
58: } CTX_BRUSSEL;
62: int main(int argc,char **argv)
63: {
64: Mat A; /* eigenvalue problem matrix */
65: EPS eps; /* eigenproblem solver context */
66: EPSType type;
67: PetscScalar delta1,delta2,L,h;
68: PetscInt N=30,n,i,Istart,Iend,nev;
69: CTX_BRUSSEL *ctx;
70: PetscBool terse;
71: PetscViewer viewer;
74: SlepcInitialize(&argc,&argv,(char*)0,help);
76: PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
77: PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%D\n\n",N);
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Generate the matrix
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83: /*
84: Create shell matrix context and set default parameters
85: */
86: PetscNew(&ctx);
87: ctx->alpha = 2.0;
88: ctx->beta = 5.45;
89: delta1 = 0.008;
90: delta2 = 0.004;
91: L = 0.51302;
93: /*
94: Look the command line for user-provided parameters
95: */
96: PetscOptionsGetScalar(NULL,NULL,"-L",&L,NULL);
97: PetscOptionsGetScalar(NULL,NULL,"-alpha",&ctx->alpha,NULL);
98: PetscOptionsGetScalar(NULL,NULL,"-beta",&ctx->beta,NULL);
99: PetscOptionsGetScalar(NULL,NULL,"-delta1",&delta1,NULL);
100: PetscOptionsGetScalar(NULL,NULL,"-delta2",&delta2,NULL);
102: /*
103: Create matrix T
104: */
105: MatCreate(PETSC_COMM_WORLD,&ctx->T);
106: MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N);
107: MatSetFromOptions(ctx->T);
108: MatSetUp(ctx->T);
110: MatGetOwnershipRange(ctx->T,&Istart,&Iend);
111: for (i=Istart;i<Iend;i++) {
112: if (i>0) { MatSetValue(ctx->T,i,i-1,1.0,INSERT_VALUES); }
113: if (i<N-1) { MatSetValue(ctx->T,i,i+1,1.0,INSERT_VALUES); }
114: MatSetValue(ctx->T,i,i,-2.0,INSERT_VALUES);
115: }
116: MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY);
117: MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY);
118: MatGetLocalSize(ctx->T,&n,NULL);
120: /*
121: Fill the remaining information in the shell matrix context
122: and create auxiliary vectors
123: */
124: h = 1.0 / (PetscReal)(N+1);
125: ctx->tau1 = delta1 / ((h*L)*(h*L));
126: ctx->tau2 = delta2 / ((h*L)*(h*L));
127: ctx->sigma = 0.0;
128: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x1);
129: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x2);
130: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y1);
131: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y2);
133: /*
134: Create the shell matrix
135: */
136: MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A);
137: MatShellSetOperation(A,MATOP_MULT,(void(*)())MatMult_Brussel);
138: MatShellSetOperation(A,MATOP_SHIFT,(void(*)())MatShift_Brussel);
139: MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)())MatGetDiagonal_Brussel);
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Create the eigensolver and set various options
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145: /*
146: Create eigensolver context
147: */
148: EPSCreate(PETSC_COMM_WORLD,&eps);
150: /*
151: Set operators. In this case, it is a standard eigenvalue problem
152: */
153: EPSSetOperators(eps,A,NULL);
154: EPSSetProblemType(eps,EPS_NHEP);
156: /*
157: Ask for the rightmost eigenvalues
158: */
159: EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
161: /*
162: Set other solver options at runtime
163: */
164: EPSSetFromOptions(eps);
166: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167: Solve the eigensystem
168: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170: EPSSolve(eps);
172: /*
173: Optional: Get some information from the solver and display it
174: */
175: EPSGetType(eps,&type);
176: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
177: EPSGetDimensions(eps,&nev,NULL,NULL);
178: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
180: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181: Display solution and clean up
182: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
184: /* show detailed info unless -terse option is given by user */
185: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
186: if (terse) {
187: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
188: } else {
189: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
190: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
191: EPSReasonView(eps,viewer);
192: EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer);
193: PetscViewerPopFormat(viewer);
194: }
195: EPSDestroy(&eps);
196: MatDestroy(&A);
197: MatDestroy(&ctx->T);
198: VecDestroy(&ctx->x1);
199: VecDestroy(&ctx->x2);
200: VecDestroy(&ctx->y1);
201: VecDestroy(&ctx->y2);
202: PetscFree(ctx);
203: SlepcFinalize();
204: return ierr;
205: }
209: PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
210: {
211: PetscInt n;
212: const PetscScalar *px;
213: PetscScalar *py;
214: CTX_BRUSSEL *ctx;
215: PetscErrorCode ierr;
218: MatShellGetContext(A,(void**)&ctx);
219: MatGetLocalSize(ctx->T,&n,NULL);
220: VecGetArrayRead(x,&px);
221: VecGetArray(y,&py);
222: VecPlaceArray(ctx->x1,px);
223: VecPlaceArray(ctx->x2,px+n);
224: VecPlaceArray(ctx->y1,py);
225: VecPlaceArray(ctx->y2,py+n);
227: MatMult(ctx->T,ctx->x1,ctx->y1);
228: VecScale(ctx->y1,ctx->tau1);
229: VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1);
230: VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2);
232: MatMult(ctx->T,ctx->x2,ctx->y2);
233: VecScale(ctx->y2,ctx->tau2);
234: VecAXPY(ctx->y2,-ctx->beta,ctx->x1);
235: VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2);
237: VecRestoreArrayRead(x,&px);
238: VecRestoreArray(y,&py);
239: VecResetArray(ctx->x1);
240: VecResetArray(ctx->x2);
241: VecResetArray(ctx->y1);
242: VecResetArray(ctx->y2);
243: return(0);
244: }
248: PetscErrorCode MatShift_Brussel(PetscScalar* a,Mat Y)
249: {
250: CTX_BRUSSEL *ctx;
254: MatShellGetContext(Y,(void**)&ctx);
255: ctx->sigma += *a;
256: return(0);
257: }
261: PetscErrorCode MatGetDiagonal_Brussel(Mat A,Vec diag)
262: {
263: Vec d1,d2;
264: PetscInt n;
265: PetscScalar *pd;
266: MPI_Comm comm;
267: CTX_BRUSSEL *ctx;
271: MatShellGetContext(A,(void**)&ctx);
272: PetscObjectGetComm((PetscObject)A,&comm);
273: MatGetLocalSize(ctx->T,&n,NULL);
274: VecGetArray(diag,&pd);
275: VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd,&d1);
276: VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd+n,&d2);
278: VecSet(d1,-2.0*ctx->tau1 + ctx->beta - 1.0 + ctx->sigma);
279: VecSet(d2,-2.0*ctx->tau2 - ctx->alpha*ctx->alpha + ctx->sigma);
281: VecDestroy(&d1);
282: VecDestroy(&d2);
283: VecRestoreArray(diag,&pd);
284: return(0);
285: }