xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex5.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 /* ------------------------------------------------------------------- */
4360f0b76eSHong Zhang PetscErrorCode InitialConditions(DM da,Vec U)
4460f0b76eSHong Zhang {
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;
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   */
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,U,&u));
5960f0b76eSHong Zhang 
6060f0b76eSHong Zhang   /*
6160f0b76eSHong Zhang      Get local grid boundaries
6260f0b76eSHong Zhang   */
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
7266baab88SBarry Smith       if (PetscApproximateGTE(x,1.0) && PetscApproximateLTE(x,1.5) && PetscApproximateGTE(y,1.0) && PetscApproximateLTE(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;
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   */
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,U,&u));
8360f0b76eSHong Zhang   PetscFunctionReturn(0);
8460f0b76eSHong Zhang }
85c4762a1bSJed Brown 
86c4762a1bSJed Brown int main(int argc,char **argv)
87c4762a1bSJed Brown {
88c4762a1bSJed Brown   TS             ts;                  /* ODE integrator */
89c4762a1bSJed Brown   Vec            x;                   /* solution */
90c4762a1bSJed Brown   PetscErrorCode ierr;
91c4762a1bSJed Brown   DM             da;
92c4762a1bSJed Brown   AppCtx         appctx;
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95c4762a1bSJed Brown      Initialize program
96c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
98c4762a1bSJed Brown   PetscFunctionBeginUser;
99c4762a1bSJed Brown   appctx.D1    = 8.0e-5;
100c4762a1bSJed Brown   appctx.D2    = 4.0e-5;
101c4762a1bSJed Brown   appctx.gamma = .024;
102c4762a1bSJed Brown   appctx.kappa = .06;
10360f0b76eSHong Zhang   appctx.aijpc = PETSC_FALSE;
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
107c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,65,65,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da));
109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
110*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
111*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,0,"u"));
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,1,"v"));
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115c4762a1bSJed Brown      Create global vector from DMDA; this will be used to store the solution
116c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&x));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120c4762a1bSJed Brown      Create timestepping solver context
121c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSARKIMEX));
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE));
125*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts,da));
126*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
127*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,&appctx));
128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(ts,NULL,NULL,RHSJacobian,&appctx));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131c4762a1bSJed Brown      Set initial conditions
132c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(InitialConditions(da,x));
134*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,x));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137c4762a1bSJed Brown      Set solver options
138c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,2000.0));
140*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,.0001));
141*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
142*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145c4762a1bSJed Brown      Solve ODE system
146c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,x));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150c4762a1bSJed Brown      Free work space.
151c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
154*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   ierr = PetscFinalize();
157c4762a1bSJed Brown   return ierr;
158c4762a1bSJed Brown }
159c4762a1bSJed Brown 
160c4762a1bSJed Brown /*TEST
161c4762a1bSJed Brown 
16260f0b76eSHong Zhang    build:
16360f0b76eSHong Zhang      depends: reaction_diffusion.c
16460f0b76eSHong Zhang 
165c4762a1bSJed Brown    test:
166c4762a1bSJed Brown       args: -ts_view  -ts_monitor -ts_max_time 500
167c4762a1bSJed Brown       requires: double
168c4762a1bSJed Brown       timeoutfactor: 3
169c4762a1bSJed Brown 
170c4762a1bSJed Brown    test:
171c4762a1bSJed Brown       suffix: 2
172c4762a1bSJed Brown       args: -ts_view  -ts_monitor -ts_max_time 500 -ts_monitor_draw_solution
173c4762a1bSJed Brown       requires: x double
174c4762a1bSJed Brown       output_file: output/ex5_1.out
175c4762a1bSJed Brown       timeoutfactor: 3
176c4762a1bSJed Brown 
177c4762a1bSJed Brown TEST*/
178