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 */ 38*60f0b76eSHong Zhang #include "reaction_diffusion.h" 39c4762a1bSJed Brown #include <petscdm.h> 40c4762a1bSJed Brown #include <petscdmda.h> 41c4762a1bSJed Brown 42*60f0b76eSHong Zhang /* ------------------------------------------------------------------- */ 43*60f0b76eSHong Zhang PetscErrorCode InitialConditions(DM da,Vec U) 44*60f0b76eSHong Zhang { 45*60f0b76eSHong Zhang PetscErrorCode ierr; 46*60f0b76eSHong Zhang PetscInt i,j,xs,ys,xm,ym,Mx,My; 47*60f0b76eSHong Zhang Field **u; 48*60f0b76eSHong Zhang PetscReal hx,hy,x,y; 49c4762a1bSJed Brown 50*60f0b76eSHong Zhang PetscFunctionBegin; 51*60f0b76eSHong Zhang ierr = 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);CHKERRQ(ierr); 52*60f0b76eSHong Zhang 53*60f0b76eSHong Zhang hx = 2.5/(PetscReal)(Mx); 54*60f0b76eSHong Zhang hy = 2.5/(PetscReal)(My); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* 57*60f0b76eSHong Zhang Get pointers to actual vector data 58c4762a1bSJed Brown */ 59*60f0b76eSHong Zhang ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr); 60*60f0b76eSHong Zhang 61*60f0b76eSHong Zhang /* 62*60f0b76eSHong Zhang Get local grid boundaries 63*60f0b76eSHong Zhang */ 64*60f0b76eSHong Zhang ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr); 65*60f0b76eSHong Zhang 66*60f0b76eSHong Zhang /* 67*60f0b76eSHong Zhang Compute function over the locally owned part of the grid 68*60f0b76eSHong Zhang */ 69*60f0b76eSHong Zhang for (j=ys; j<ys+ym; j++) { 70*60f0b76eSHong Zhang y = j*hy; 71*60f0b76eSHong Zhang for (i=xs; i<xs+xm; i++) { 72*60f0b76eSHong Zhang x = i*hx; 73*60f0b76eSHong Zhang if (PetscGTE(x,1.0) && PetscLTE(x,1.5) && PetscGTE(y,1.0) && PetscLTE(y,1.5)) u[j][i].v = PetscPowReal(PetscSinReal(4.0*PETSC_PI*x),2.0)*PetscPowReal(PetscSinReal(4.0*PETSC_PI*y),2.0)/4.0; 74*60f0b76eSHong Zhang else u[j][i].v = 0.0; 75*60f0b76eSHong Zhang 76*60f0b76eSHong Zhang u[j][i].u = 1.0 - 2.0*u[j][i].v; 77*60f0b76eSHong Zhang } 78*60f0b76eSHong Zhang } 79*60f0b76eSHong Zhang 80*60f0b76eSHong Zhang /* 81*60f0b76eSHong Zhang Restore access to vector 82*60f0b76eSHong Zhang */ 83*60f0b76eSHong Zhang ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr); 84*60f0b76eSHong Zhang PetscFunctionReturn(0); 85*60f0b76eSHong Zhang } 86c4762a1bSJed Brown 87c4762a1bSJed Brown int main(int argc,char **argv) 88c4762a1bSJed Brown { 89c4762a1bSJed Brown TS ts; /* ODE integrator */ 90c4762a1bSJed Brown Vec x; /* solution */ 91c4762a1bSJed Brown PetscErrorCode ierr; 92c4762a1bSJed Brown DM da; 93c4762a1bSJed Brown AppCtx appctx; 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 96c4762a1bSJed Brown Initialize program 97c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 98c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 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; 104*60f0b76eSHong Zhang appctx.aijpc = PETSC_FALSE; 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 107c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 108c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 109c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,65,65,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = DMDASetFieldName(da,0,"u");CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = DMDASetFieldName(da,1,"v");CHKERRQ(ierr); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 116c4762a1bSJed Brown Create global vector from DMDA; this will be used to store the solution 117c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 118c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 121c4762a1bSJed Brown Create timestepping solver context 122c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 123c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = TSSetDM(ts,da);CHKERRQ(ierr); 127c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&appctx);CHKERRQ(ierr); 129c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,NULL,NULL,RHSJacobian,&appctx);CHKERRQ(ierr); 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 132c4762a1bSJed Brown Set initial conditions 133c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 134c4762a1bSJed Brown ierr = InitialConditions(da,x);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 138c4762a1bSJed Brown Set solver options 139c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 140c4762a1bSJed Brown ierr = TSSetMaxTime(ts,2000.0);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = TSSetTimeStep(ts,.0001);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 146c4762a1bSJed Brown Solve ODE system 147c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 148c4762a1bSJed Brown ierr = TSSolve(ts,x);CHKERRQ(ierr); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 151c4762a1bSJed Brown Free work space. 152c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 153c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 156c4762a1bSJed Brown 157c4762a1bSJed Brown ierr = PetscFinalize(); 158c4762a1bSJed Brown return ierr; 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown /*TEST 162c4762a1bSJed Brown 163*60f0b76eSHong Zhang build: 164*60f0b76eSHong Zhang depends: reaction_diffusion.c 165*60f0b76eSHong 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