Actual source code: ex35.cxx
1: static const char help[] = "Time-dependent Brusselator reaction-diffusion PDE in 1d. Demonstrates IMEX methods and uses MOAB.\n";
2: /*
3: u_t - alpha u_xx = A + u^2 v - (B+1) u
4: v_t - alpha v_xx = B u - u^2 v
5: 0 < x < 1;
6: A = 1, B = 3, alpha = 1/50
8: Initial conditions:
9: u(x,0) = 1 + sin(2 pi x)
10: v(x,0) = 3
12: Boundary conditions:
13: u(0,t) = u(1,t) = 1
14: v(0,t) = v(1,t) = 3
15: */
17: // PETSc includes:
18: #include <petscts.h>
19: #include <petscdmmoab.h>
21: typedef struct {
22: PetscScalar u,v;
23: } Field;
25: struct pUserCtx {
26: PetscReal A,B; /* Reaction coefficients */
27: PetscReal alpha; /* Diffusion coefficient */
28: Field leftbc; /* Dirichlet boundary conditions at left boundary */
29: Field rightbc; /* Dirichlet boundary conditions at right boundary */
30: PetscInt n,npts; /* Number of mesh points */
31: PetscInt ntsteps; /* Number of time steps */
32: PetscInt nvars; /* Number of variables in the equation system */
33: PetscBool io;
34: };
35: typedef pUserCtx* UserCtx;
37: PetscErrorCode Initialize_AppContext(UserCtx *puser)
38: {
39: UserCtx user;
40: PetscErrorCode ierr;
42: PetscNew(&user);
44: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","ex35.cxx");
45: {
46: user->nvars = 2;
47: user->A = 1;
48: user->B = 3;
49: user->alpha = 0.02;
50: user->leftbc.u = 1;
51: user->rightbc.u = 1;
52: user->leftbc.v = 3;
53: user->rightbc.v = 3;
54: user->n = 10;
55: user->ntsteps = 10000;
56: user->io = PETSC_FALSE;
57: PetscOptionsReal("-A","Reaction rate","ex35.cxx",user->A,&user->A,NULL);
58: PetscOptionsReal("-B","Reaction rate","ex35.cxx",user->B,&user->B,NULL);
59: PetscOptionsReal("-alpha","Diffusion coefficient","ex35.cxx",user->alpha,&user->alpha,NULL);
60: PetscOptionsScalar("-uleft","Dirichlet boundary condition","ex35.cxx",user->leftbc.u,&user->leftbc.u,NULL);
61: PetscOptionsScalar("-uright","Dirichlet boundary condition","ex35.cxx",user->rightbc.u,&user->rightbc.u,NULL);
62: PetscOptionsScalar("-vleft","Dirichlet boundary condition","ex35.cxx",user->leftbc.v,&user->leftbc.v,NULL);
63: PetscOptionsScalar("-vright","Dirichlet boundary condition","ex35.cxx",user->rightbc.v,&user->rightbc.v,NULL);
64: PetscOptionsInt("-n","Number of 1-D elements","ex35.cxx",user->n,&user->n,NULL);
65: PetscOptionsInt("-ndt","Number of time steps","ex35.cxx",user->ntsteps,&user->ntsteps,NULL);
66: PetscOptionsBool("-io","Write the mesh and solution output to a file.","ex35.cxx",user->io,&user->io,NULL);
67: user->npts = user->n+1;
68: }
69: PetscOptionsEnd();
71: *puser = user;
72: return 0;
73: }
75: PetscErrorCode Destroy_AppContext(UserCtx *user)
76: {
77: PetscFree(*user);
78: return 0;
79: }
81: static PetscErrorCode FormInitialSolution(TS,Vec,void*);
82: static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
83: static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
84: static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
86: /****************
87: * *
88: * MAIN *
89: * *
90: ****************/
91: int main(int argc,char **argv)
92: {
93: TS ts; /* nonlinear solver */
94: Vec X; /* solution, residual vectors */
95: Mat J; /* Jacobian matrix */
96: PetscInt steps,mx;
97: PetscReal hx,dt,ftime;
98: UserCtx user; /* user-defined work context */
99: TSConvergedReason reason;
100: DM dm;
101: const char *fields[2] = {"U","V"};
103: PetscInitialize(&argc,&argv,(char *)0,help);
105: /* Initialize the user context struct */
106: Initialize_AppContext(&user);
108: /* Fill in the user defined work context: */
109: DMMoabCreateBoxMesh(PETSC_COMM_WORLD, 1, PETSC_FALSE, NULL, user->n, 1, &dm);
110: DMMoabSetFieldNames(dm, user->nvars, fields);
111: DMMoabSetBlockSize(dm, user->nvars);
112: DMSetFromOptions(dm);
114: /* SetUp the data structures for DMMOAB */
115: DMSetUp(dm);
117: /* Create timestepping solver context */
118: TSCreate(PETSC_COMM_WORLD,&ts);
119: TSSetDM(ts, dm);
120: TSSetType(ts,TSARKIMEX);
121: TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1);
122: DMSetMatType(dm,MATBAIJ);
123: DMCreateMatrix(dm,&J);
125: TSSetRHSFunction(ts,NULL,FormRHSFunction,user);
126: TSSetIFunction(ts,NULL,FormIFunction,user);
127: TSSetIJacobian(ts,J,J,FormIJacobian,user);
129: ftime = 10.0;
130: TSSetMaxSteps(ts,user->ntsteps);
131: TSSetMaxTime(ts,ftime);
132: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
134: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135: Create the solution vector and set the initial conditions
136: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137: DMCreateGlobalVector(dm, &X);
139: FormInitialSolution(ts,X,user);
140: TSSetSolution(ts,X);
141: VecGetSize(X,&mx);
142: hx = 1.0/(PetscReal)(mx/2-1);
143: dt = 0.4 * PetscSqr(hx) / user->alpha; /* Diffusive stability limit */
144: TSSetTimeStep(ts,dt);
146: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147: Set runtime options
148: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149: TSSetFromOptions(ts);
151: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152: Solve nonlinear system
153: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154: TSSolve(ts,X);
155: TSGetSolveTime(ts,&ftime);
156: TSGetStepNumber(ts,&steps);
157: TSGetConvergedReason(ts,&reason);
158: PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],ftime,steps);
160: if (user->io) {
161: /* Print the numerical solution to screen and then dump to file */
162: VecView(X,PETSC_VIEWER_STDOUT_WORLD);
164: /* Write out the solution along with the mesh */
165: DMMoabSetGlobalFieldVector(dm, X);
166: #ifdef MOAB_HAVE_HDF5
167: DMMoabOutput(dm, "ex35.h5m", "");
168: #else
169: /* MOAB does not support true parallel writers that aren't HDF5 based
170: And so if you are using VTK as the output format in parallel,
171: the data could be jumbled due to the order in which the processors
172: write out their parts of the mesh and solution tags
173: */
174: DMMoabOutput(dm, "ex35.vtk", "");
175: #endif
176: }
178: /* Free work space.
179: Free all PETSc related resources: */
180: MatDestroy(&J);
181: VecDestroy(&X);
182: TSDestroy(&ts);
183: DMDestroy(&dm);
185: /* Free all MOAB related resources: */
186: Destroy_AppContext(&user);
187: PetscFinalize();
188: return 0;
189: }
191: /*
192: IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
193: */
194: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr)
195: {
196: UserCtx user = (UserCtx)ptr;
197: PetscInt dof;
198: PetscReal hx;
199: DM dm;
200: const moab::Range *vlocal;
201: PetscBool vonboundary;
203: TSGetDM(ts, &dm);
205: /* get the essential MOAB mesh related quantities needed for FEM assembly */
206: DMMoabGetLocalVertices(dm, &vlocal, NULL);
208: /* compute local element sizes - structured grid */
209: hx = 1.0/user->n;
211: /* Compute function over the locally owned part of the grid
212: Assemble the operator by looping over edges and computing
213: contribution for each vertex dof */
214: for (moab::Range::iterator iter = vlocal->begin(); iter != vlocal->end(); iter++) {
215: const moab::EntityHandle vhandle = *iter;
217: DMMoabGetDofsBlocked(dm, 1, &vhandle, &dof);
219: /* check if vertex is on the boundary */
220: DMMoabIsEntityOnBoundary(dm,vhandle,&vonboundary);
222: if (vonboundary) {
223: const PetscScalar bcvals[2][2] = {{hx,0},{0,hx}};
224: MatSetValuesBlocked(Jpre,1,&dof,1,&dof,&bcvals[0][0],INSERT_VALUES);
225: }
226: else {
227: const PetscInt row = dof,col[] = {dof-1,dof,dof+1};
228: const PetscScalar dxxL = -user->alpha/hx,dxx0 = 2.*user->alpha/hx,dxxR = -user->alpha/hx;
229: const PetscScalar vals[2][3][2] = {{{dxxL,0},{a *hx+dxx0,0},{dxxR,0}},
230: {{0,dxxL},{0,a*hx+dxx0},{0,dxxR}}};
231: MatSetValuesBlocked(Jpre,1,&row,3,col,&vals[0][0][0],INSERT_VALUES);
232: }
233: }
235: MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
236: MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
237: if (J != Jpre) {
238: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
239: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
240: }
241: return 0;
242: }
244: static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
245: {
246: UserCtx user = (UserCtx)ptr;
247: DM dm;
248: PetscReal hx;
249: const Field *x;
250: Field *f;
251: PetscInt dof;
252: const moab::Range *ownedvtx;
254: hx = 1.0/user->n;
255: TSGetDM(ts,&dm);
257: /* Get pointers to vector data */
258: VecSet(F,0.0);
260: DMMoabVecGetArrayRead(dm, X, &x);
261: DMMoabVecGetArray(dm, F, &f);
263: DMMoabGetLocalVertices(dm, &ownedvtx, NULL);
265: /* Compute function over the locally owned part of the grid */
266: for (moab::Range::iterator iter = ownedvtx->begin(); iter != ownedvtx->end(); iter++) {
267: const moab::EntityHandle vhandle = *iter;
268: DMMoabGetDofsBlockedLocal(dm, 1, &vhandle, &dof);
270: PetscScalar u = x[dof].u,v = x[dof].v;
271: f[dof].u = hx*(user->A + u*u*v - (user->B+1)*u);
272: f[dof].v = hx*(user->B*u - u*u*v);
273: }
275: /* Restore vectors */
276: DMMoabVecRestoreArrayRead(dm, X, &x);
277: DMMoabVecRestoreArray(dm, F, &f);
278: return 0;
279: }
281: static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
282: {
283: UserCtx user = (UserCtx)ctx;
284: DM dm;
285: Field *x,*xdot,*f;
286: PetscReal hx;
287: Vec Xloc;
288: PetscInt i,bcindx;
289: PetscBool elem_on_boundary;
290: const moab::Range *vlocal;
292: hx = 1.0/user->n;
293: TSGetDM(ts, &dm);
295: /* get the essential MOAB mesh related quantities needed for FEM assembly */
296: DMMoabGetLocalVertices(dm, &vlocal, NULL);
298: /* reset the residual vector */
299: VecSet(F,0.0);
301: DMGetLocalVector(dm,&Xloc);
302: DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);
303: DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);
305: /* get the local representation of the arrays from Vectors */
306: DMMoabVecGetArrayRead(dm, Xloc, &x);
307: DMMoabVecGetArrayRead(dm, Xdot, &xdot);
308: DMMoabVecGetArray(dm, F, &f);
310: /* loop over local elements */
311: for (moab::Range::iterator iter = vlocal->begin(); iter != vlocal->end(); iter++) {
312: const moab::EntityHandle vhandle = *iter;
314: DMMoabGetDofsBlockedLocal(dm,1,&vhandle,&i);
316: /* check if vertex is on the boundary */
317: DMMoabIsEntityOnBoundary(dm,vhandle,&elem_on_boundary);
319: if (elem_on_boundary) {
320: DMMoabGetDofsBlocked(dm, 1, &vhandle, &bcindx);
321: if (bcindx == 0) { /* Apply left BC */
322: f[i].u = hx * (x[i].u - user->leftbc.u);
323: f[i].v = hx * (x[i].v - user->leftbc.v);
324: } else { /* Apply right BC */
325: f[i].u = hx * (x[i].u - user->rightbc.u);
326: f[i].v = hx * (x[i].v - user->rightbc.v);
327: }
328: }
329: else {
330: f[i].u = hx * xdot[i].u - user->alpha * (x[i-1].u - 2.*x[i].u + x[i+1].u) / hx;
331: f[i].v = hx * xdot[i].v - user->alpha * (x[i-1].v - 2.*x[i].v + x[i+1].v) / hx;
332: }
333: }
335: /* Restore data */
336: DMMoabVecRestoreArrayRead(dm, Xloc, &x);
337: DMMoabVecRestoreArrayRead(dm, Xdot, &xdot);
338: DMMoabVecRestoreArray(dm, F, &f);
339: DMRestoreLocalVector(dm, &Xloc);
340: return 0;
341: }
343: PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
344: {
345: UserCtx user = (UserCtx)ctx;
346: PetscReal vpos[3];
347: DM dm;
348: Field *x;
349: const moab::Range *vowned;
350: PetscInt dof;
351: moab::Range::iterator iter;
353: TSGetDM(ts, &dm);
355: /* get the essential MOAB mesh related quantities needed for FEM assembly */
356: DMMoabGetLocalVertices(dm, &vowned, NULL);
358: VecSet(X, 0.0);
360: /* Get pointers to vector data */
361: DMMoabVecGetArray(dm, X, &x);
363: /* Compute function over the locally owned part of the grid */
364: for (moab::Range::iterator iter = vowned->begin(); iter != vowned->end(); iter++) {
365: const moab::EntityHandle vhandle = *iter;
366: DMMoabGetDofsBlockedLocal(dm, 1, &vhandle, &dof);
368: /* compute the mid-point of the element and use a 1-point lumped quadrature */
369: DMMoabGetVertexCoordinates(dm,1,&vhandle,vpos);
371: PetscReal xi = vpos[0];
372: x[dof].u = user->leftbc.u*(1.-xi) + user->rightbc.u*xi + PetscSinReal(2.*PETSC_PI*xi);
373: x[dof].v = user->leftbc.v*(1.-xi) + user->rightbc.v*xi;
374: }
376: /* Restore vectors */
377: DMMoabVecRestoreArray(dm, X, &x);
378: return 0;
379: }
381: /*TEST
383: build:
384: requires: moab
386: test:
387: args: -n 20 -ts_type rosw -ts_rosw_type 2p -ts_dt 5e-2 -ts_adapt_type none
389: test:
390: suffix: 2
391: nsize: 2
392: args: -n 50 -ts_type glee -ts_adapt_type none -ts_dt 0.1 -io
393: TODO:
395: TEST*/