1c4762a1bSJed Brown static char help[] = "Demonstrates Pattern Formation with Reaction-Diffusion Equations.\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown /*F
4c4762a1bSJed Brown This example is taken from the book, Numerical Solution of Time-Dependent Advection-Diffusion-Reaction Equations by
5c4762a1bSJed Brown W. Hundsdorf and J.G. Verwer, Page 21, Pattern Formation with Reaction-Diffusion Equations
6c4762a1bSJed Brown \begin{eqnarray*}
7c4762a1bSJed Brown u_t = D_1 (u_{xx} + u_{yy}) - u*v^2 + \gamma(1 -u) \\
8c4762a1bSJed Brown v_t = D_2 (v_{xx} + v_{yy}) + u*v^2 - (\gamma + \kappa)v
9c4762a1bSJed Brown \end{eqnarray*}
10c4762a1bSJed Brown Unlike in the book this uses periodic boundary conditions instead of Neumann
11c4762a1bSJed Brown (since they are easier for finite differences).
12c4762a1bSJed Brown F*/
13c4762a1bSJed Brown
14c4762a1bSJed Brown /*
15c4762a1bSJed Brown Helpful runtime monitor options:
16c4762a1bSJed Brown -ts_monitor_draw_solution
17c4762a1bSJed Brown -draw_save -draw_save_movie
18c4762a1bSJed Brown
19c4762a1bSJed Brown Helpful runtime linear solver options:
20c4762a1bSJed Brown -pc_type mg -pc_mg_galerkin pmat -da_refine 1 -snes_monitor -ksp_monitor -ts_view (note that these Jacobians are so well-conditioned multigrid may not be the best solver)
21c4762a1bSJed Brown
22c4762a1bSJed Brown Point your browser to localhost:8080 to monitor the simulation
23c4762a1bSJed Brown ./ex5 -ts_view_pre saws -stack_view saws -draw_save -draw_save_single_file -x_virtual -ts_monitor_draw_solution -saws_root .
24c4762a1bSJed Brown
25c4762a1bSJed Brown */
26c4762a1bSJed Brown
27c4762a1bSJed Brown /*
28c4762a1bSJed Brown
29c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
30c4762a1bSJed Brown Include "petscts.h" so that we can use SNES numerical (ODE) integrators. Note that this
31c4762a1bSJed Brown file automatically includes:
32c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors
33c4762a1bSJed Brown petscmat.h - matrices petscis.h - index sets
34c4762a1bSJed Brown petscksp.h - Krylov subspace methods petscpc.h - preconditioners
35c4762a1bSJed Brown petscviewer.h - viewers petscsnes.h - nonlinear solvers
36c4762a1bSJed Brown */
3760f0b76eSHong Zhang #include "reaction_diffusion.h"
38c4762a1bSJed Brown #include <petscdm.h>
39c4762a1bSJed Brown #include <petscdmda.h>
40c4762a1bSJed Brown
4160f0b76eSHong Zhang /* ------------------------------------------------------------------- */
InitialConditions(DM da,Vec U)42d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(DM da, Vec U)
43d71ae5a4SJacob Faibussowitsch {
4460f0b76eSHong Zhang PetscInt i, j, xs, ys, xm, ym, Mx, My;
4560f0b76eSHong Zhang Field **u;
4660f0b76eSHong Zhang PetscReal hx, hy, x, y;
47c4762a1bSJed Brown
4860f0b76eSHong Zhang PetscFunctionBegin;
499566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
5060f0b76eSHong Zhang
51*57508eceSPierre Jolivet hx = 2.5 / (PetscReal)Mx;
52*57508eceSPierre Jolivet hy = 2.5 / (PetscReal)My;
53c4762a1bSJed Brown
54c4762a1bSJed Brown /*
5560f0b76eSHong Zhang Get pointers to actual vector data
56c4762a1bSJed Brown */
579566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u));
5860f0b76eSHong Zhang
5960f0b76eSHong Zhang /*
6060f0b76eSHong Zhang Get local grid boundaries
6160f0b76eSHong Zhang */
629566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
6360f0b76eSHong Zhang
6460f0b76eSHong Zhang /*
6560f0b76eSHong Zhang Compute function over the locally owned part of the grid
6660f0b76eSHong Zhang */
6760f0b76eSHong Zhang for (j = ys; j < ys + ym; j++) {
6860f0b76eSHong Zhang y = j * hy;
6960f0b76eSHong Zhang for (i = xs; i < xs + xm; i++) {
7060f0b76eSHong Zhang x = i * hx;
719371c9d4SSatish Balay if (PetscApproximateGTE(x, 1.0) && PetscApproximateLTE(x, 1.5) && PetscApproximateGTE(y, 1.0) && PetscApproximateLTE(y, 1.5))
729371c9d4SSatish Balay u[j][i].v = PetscPowReal(PetscSinReal(4.0 * PETSC_PI * x), 2.0) * PetscPowReal(PetscSinReal(4.0 * PETSC_PI * y), 2.0) / 4.0;
7360f0b76eSHong Zhang else u[j][i].v = 0.0;
7460f0b76eSHong Zhang
7560f0b76eSHong Zhang u[j][i].u = 1.0 - 2.0 * u[j][i].v;
7660f0b76eSHong Zhang }
7760f0b76eSHong Zhang }
7860f0b76eSHong Zhang
7960f0b76eSHong Zhang /*
8060f0b76eSHong Zhang Restore access to vector
8160f0b76eSHong Zhang */
829566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u));
833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8460f0b76eSHong Zhang }
85c4762a1bSJed Brown
main(int argc,char ** argv)86d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
87d71ae5a4SJacob Faibussowitsch {
88c4762a1bSJed Brown TS ts; /* ODE integrator */
89c4762a1bSJed Brown Vec x; /* solution */
90c4762a1bSJed Brown DM da;
91c4762a1bSJed Brown AppCtx appctx;
92c4762a1bSJed Brown
93c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94c4762a1bSJed Brown Initialize program
95c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96327415f7SBarry Smith PetscFunctionBeginUser;
97c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
98c4762a1bSJed Brown appctx.D1 = 8.0e-5;
99c4762a1bSJed Brown appctx.D2 = 4.0e-5;
100c4762a1bSJed Brown appctx.gamma = .024;
101c4762a1bSJed Brown appctx.kappa = .06;
10260f0b76eSHong Zhang appctx.aijpc = PETSC_FALSE;
103c4762a1bSJed Brown
104c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors
106c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1079566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, DM_BOUNDARY_PERIODIC, DMDA_STENCIL_STAR, 65, 65, PETSC_DECIDE, PETSC_DECIDE, 2, 1, NULL, NULL, &da));
1089566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
1099566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
1109566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "u"));
1119566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "v"));
112c4762a1bSJed Brown
113c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114c4762a1bSJed Brown Create global vector from DMDA; this will be used to store the solution
115c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1169566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x));
117c4762a1bSJed Brown
118c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119c4762a1bSJed Brown Create timestepping solver context
120c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1219566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1229566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSARKIMEX));
1239566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetFullyImplicit(ts, PETSC_TRUE));
1249566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da));
1259566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1269566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &appctx));
1279566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSJacobian, &appctx));
128c4762a1bSJed Brown
129c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130c4762a1bSJed Brown Set initial conditions
131c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1329566063dSJacob Faibussowitsch PetscCall(InitialConditions(da, x));
1339566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x));
134c4762a1bSJed Brown
135c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136c4762a1bSJed Brown Set solver options
137c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1389566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 2000.0));
1399566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .0001));
1409566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1419566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
142c4762a1bSJed Brown
143c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144c4762a1bSJed Brown Solve ODE system
145c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1469566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x));
147c4762a1bSJed Brown
148c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149c4762a1bSJed Brown Free work space.
150c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1529566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
1539566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
154c4762a1bSJed Brown
1559566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
156b122ec5aSJacob Faibussowitsch return 0;
157c4762a1bSJed Brown }
158c4762a1bSJed Brown
159c4762a1bSJed Brown /*TEST
160c4762a1bSJed Brown
16160f0b76eSHong Zhang build:
16260f0b76eSHong Zhang depends: reaction_diffusion.c
16360f0b76eSHong Zhang
164c4762a1bSJed Brown test:
165c4762a1bSJed Brown args: -ts_view -ts_monitor -ts_max_time 500
166c4762a1bSJed Brown requires: double
167c4762a1bSJed Brown timeoutfactor: 3
168c4762a1bSJed Brown
169c4762a1bSJed Brown test:
170c4762a1bSJed Brown suffix: 2
171c4762a1bSJed Brown args: -ts_view -ts_monitor -ts_max_time 500 -ts_monitor_draw_solution
172c4762a1bSJed Brown requires: x double
173c4762a1bSJed Brown output_file: output/ex5_1.out
174c4762a1bSJed Brown timeoutfactor: 3
175c4762a1bSJed Brown
176c4762a1bSJed Brown TEST*/
177