1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Demonstrates Pattern Formation with Reaction-Diffusion Equations.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*F 5c4762a1bSJed Brown This example is taken from the book, Numerical Solution of Time-Dependent Advection-Diffusion-Reaction Equations by 6c4762a1bSJed Brown W. Hundsdorf and J.G. Verwer, Page 21, Pattern Formation with Reaction-Diffusion Equations 7c4762a1bSJed Brown \begin{eqnarray*} 8c4762a1bSJed Brown u_t = D_1 (u_{xx} + u_{yy}) - u*v^2 + \gamma(1 -u) \\ 9c4762a1bSJed Brown v_t = D_2 (v_{xx} + v_{yy}) + u*v^2 - (\gamma + \kappa)v 10c4762a1bSJed Brown \end{eqnarray*} 11c4762a1bSJed Brown Unlike in the book this uses periodic boundary conditions instead of Neumann 12c4762a1bSJed Brown (since they are easier for finite differences). 13c4762a1bSJed Brown F*/ 14c4762a1bSJed Brown 15c4762a1bSJed Brown /* 16c4762a1bSJed Brown Helpful runtime monitor options: 17c4762a1bSJed Brown -ts_monitor_draw_solution 18c4762a1bSJed Brown -draw_save -draw_save_movie 19c4762a1bSJed Brown 20c4762a1bSJed Brown Helpful runtime linear solver options: 21c4762a1bSJed 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) 22c4762a1bSJed Brown 23c4762a1bSJed Brown Point your browser to localhost:8080 to monitor the simulation 24c4762a1bSJed Brown ./ex5 -ts_view_pre saws -stack_view saws -draw_save -draw_save_single_file -x_virtual -ts_monitor_draw_solution -saws_root . 25c4762a1bSJed Brown 26c4762a1bSJed Brown */ 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* 29c4762a1bSJed Brown 30c4762a1bSJed Brown Include "petscdmda.h" so that we can use distributed arrays (DMDAs). 31c4762a1bSJed Brown Include "petscts.h" so that we can use SNES numerical (ODE) integrators. Note that this 32c4762a1bSJed Brown file automatically includes: 33c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 34c4762a1bSJed Brown petscmat.h - matrices petscis.h - index sets 35c4762a1bSJed Brown petscksp.h - Krylov subspace methods petscpc.h - preconditioners 36c4762a1bSJed Brown petscviewer.h - viewers petscsnes.h - nonlinear solvers 37c4762a1bSJed Brown */ 3860f0b76eSHong Zhang #include "reaction_diffusion.h" 39c4762a1bSJed Brown #include <petscdm.h> 40c4762a1bSJed Brown #include <petscdmda.h> 41c4762a1bSJed Brown 4260f0b76eSHong Zhang /* ------------------------------------------------------------------- */ 43d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(DM da, Vec U) 44d71ae5a4SJacob Faibussowitsch { 4560f0b76eSHong Zhang PetscInt i, j, xs, ys, xm, ym, Mx, My; 4660f0b76eSHong Zhang Field **u; 4760f0b76eSHong Zhang PetscReal hx, hy, x, y; 48c4762a1bSJed Brown 4960f0b76eSHong Zhang PetscFunctionBegin; 509566063dSJacob 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)); 5160f0b76eSHong Zhang 5260f0b76eSHong Zhang hx = 2.5 / (PetscReal)(Mx); 5360f0b76eSHong Zhang hy = 2.5 / (PetscReal)(My); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* 5660f0b76eSHong Zhang Get pointers to actual vector data 57c4762a1bSJed Brown */ 589566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 5960f0b76eSHong Zhang 6060f0b76eSHong Zhang /* 6160f0b76eSHong Zhang Get local grid boundaries 6260f0b76eSHong Zhang */ 639566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 6460f0b76eSHong Zhang 6560f0b76eSHong Zhang /* 6660f0b76eSHong Zhang Compute function over the locally owned part of the grid 6760f0b76eSHong Zhang */ 6860f0b76eSHong Zhang for (j = ys; j < ys + ym; j++) { 6960f0b76eSHong Zhang y = j * hy; 7060f0b76eSHong Zhang for (i = xs; i < xs + xm; i++) { 7160f0b76eSHong Zhang x = i * hx; 729371c9d4SSatish Balay if (PetscApproximateGTE(x, 1.0) && PetscApproximateLTE(x, 1.5) && PetscApproximateGTE(y, 1.0) && PetscApproximateLTE(y, 1.5)) 739371c9d4SSatish 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; 7460f0b76eSHong Zhang else u[j][i].v = 0.0; 7560f0b76eSHong Zhang 7660f0b76eSHong Zhang u[j][i].u = 1.0 - 2.0 * u[j][i].v; 7760f0b76eSHong Zhang } 7860f0b76eSHong Zhang } 7960f0b76eSHong Zhang 8060f0b76eSHong Zhang /* 8160f0b76eSHong Zhang Restore access to vector 8260f0b76eSHong Zhang */ 839566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 84*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8560f0b76eSHong Zhang } 86c4762a1bSJed Brown 87d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 88d71ae5a4SJacob Faibussowitsch { 89c4762a1bSJed Brown TS ts; /* ODE integrator */ 90c4762a1bSJed Brown Vec x; /* solution */ 91c4762a1bSJed Brown DM da; 92c4762a1bSJed Brown AppCtx appctx; 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 95c4762a1bSJed Brown Initialize program 96c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 97327415f7SBarry Smith PetscFunctionBeginUser; 989566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 99c4762a1bSJed Brown PetscFunctionBeginUser; 100c4762a1bSJed Brown appctx.D1 = 8.0e-5; 101c4762a1bSJed Brown appctx.D2 = 4.0e-5; 102c4762a1bSJed Brown appctx.gamma = .024; 103c4762a1bSJed Brown appctx.kappa = .06; 10460f0b76eSHong Zhang appctx.aijpc = PETSC_FALSE; 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 107c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 108c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1099566063dSJacob 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)); 1109566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 1119566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1129566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "u")); 1139566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "v")); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 116c4762a1bSJed Brown Create global vector from DMDA; this will be used to store the solution 117c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1189566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 121c4762a1bSJed Brown Create timestepping solver context 122c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1239566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1249566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSARKIMEX)); 1259566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetFullyImplicit(ts, PETSC_TRUE)); 1269566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 1279566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1289566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &appctx)); 1299566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSJacobian, &appctx)); 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 132c4762a1bSJed Brown Set initial conditions 133c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1349566063dSJacob Faibussowitsch PetscCall(InitialConditions(da, x)); 1359566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 138c4762a1bSJed Brown Set solver options 139c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1409566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 2000.0)); 1419566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .0001)); 1429566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 1439566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 146c4762a1bSJed Brown Solve ODE system 147c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1489566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 151c4762a1bSJed Brown Free work space. 152c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1549566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1559566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 156c4762a1bSJed Brown 1579566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 158b122ec5aSJacob Faibussowitsch return 0; 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown /*TEST 162c4762a1bSJed Brown 16360f0b76eSHong Zhang build: 16460f0b76eSHong Zhang depends: reaction_diffusion.c 16560f0b76eSHong Zhang 166c4762a1bSJed Brown test: 167c4762a1bSJed Brown args: -ts_view -ts_monitor -ts_max_time 500 168c4762a1bSJed Brown requires: double 169c4762a1bSJed Brown timeoutfactor: 3 170c4762a1bSJed Brown 171c4762a1bSJed Brown test: 172c4762a1bSJed Brown suffix: 2 173c4762a1bSJed Brown args: -ts_view -ts_monitor -ts_max_time 500 -ts_monitor_draw_solution 174c4762a1bSJed Brown requires: x double 175c4762a1bSJed Brown output_file: output/ex5_1.out 176c4762a1bSJed Brown timeoutfactor: 3 177c4762a1bSJed Brown 178c4762a1bSJed Brown TEST*/ 179