xref: /petsc/src/dm/tutorials/ex2.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1c4762a1bSJed Brown static char help[] = "Demonstrates Conway's Game of Life using a 2d DMDA.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown  At each step in time, the following transitions occur:
5c4762a1bSJed Brown 
6c4762a1bSJed Brown     Any live cell with fewer than two live neighbours dies, as if by underpopulation.
7c4762a1bSJed Brown     Any live cell with two or three live neighbours lives on to the next generation.
8c4762a1bSJed Brown     Any live cell with more than three live neighbours dies, as if by overpopulation.
9c4762a1bSJed Brown     Any dead cell with exactly three live neighbours becomes a live cell, as if by reproduction.
10c4762a1bSJed Brown 
11c4762a1bSJed Brown  https://en.wikipedia.org/wiki/Conway%27s_Game_of_Life
12c4762a1bSJed Brown */
13c4762a1bSJed Brown 
14c4762a1bSJed Brown #include <petscdm.h>
15c4762a1bSJed Brown #include <petscdmda.h>
16c4762a1bSJed Brown 
17c4762a1bSJed Brown static const int GLIDER[3][3] = {
18c4762a1bSJed Brown   {0, 1, 0},
19c4762a1bSJed Brown   {0, 1, 1},
20c4762a1bSJed Brown   {1, 0, 1}
21c4762a1bSJed Brown };
22c4762a1bSJed Brown 
main(int argc,char ** argv)23d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
24d71ae5a4SJacob Faibussowitsch {
25c4762a1bSJed Brown   DM          da;
26c4762a1bSJed Brown   PetscViewer viewer;
27c4762a1bSJed Brown   Vec         Xlocal, Xglobal;
28c4762a1bSJed Brown   PetscInt    glider_loc[2] = {10, 20}, blinker_loc[2] = {20, 10}, two, steps = 100, viz_interval = 1;
29c4762a1bSJed Brown   PetscInt    check_step_alive = -1, check_step_dead = -1;
30c4762a1bSJed Brown   PetscBool   has_glider, has_blinker;
31c4762a1bSJed Brown 
32327415f7SBarry Smith   PetscFunctionBeginUser;
33c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
34d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Conway's Game of Life", "");
35c4762a1bSJed Brown   {
369566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-glider", "Coordinate at which to center a glider", NULL, glider_loc, (two = 2, &two), &has_glider));
379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-blinker", "Coordinate at which to center a blinker", NULL, blinker_loc, (two = 2, &two), &has_blinker));
389566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-steps", "Number of steps to take", NULL, steps, &steps, NULL));
39d5b43468SJose E. Roman     PetscCall(PetscOptionsInt("-viz_interval", "Visualization interval", NULL, viz_interval, &viz_interval, NULL));
409566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-check_step_alive", "Step on which to check that the simulation is alive", NULL, check_step_alive, &check_step_alive, NULL));
419566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-check_step_dead", "Step on which to check that the simulation is dead", NULL, check_step_dead, &check_step_dead, NULL));
42c4762a1bSJed Brown   }
43d0609cedSBarry Smith   PetscOptionsEnd();
44c4762a1bSJed Brown 
459566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD, NULL, "Life", PETSC_DECIDE, PETSC_DECIDE, 1000, 1000, &viewer));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Create distributed array and get vectors */
489566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, DM_BOUNDARY_PERIODIC, DMDA_STENCIL_BOX, 30, 30, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da));
499566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
509566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
519566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(da, &Xlocal));
529566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &Xglobal));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   { /* Initialize */
55c4762a1bSJed Brown     DMDALocalInfo info;
56c4762a1bSJed Brown     PetscScalar **x;
57c4762a1bSJed Brown     PetscInt      i, j;
58c4762a1bSJed Brown 
599566063dSJacob Faibussowitsch     PetscCall(DMDAGetLocalInfo(da, &info));
609566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da, Xlocal, &x));
61c4762a1bSJed Brown     for (j = info.ys; j < info.ys + info.ym; j++) {
62c4762a1bSJed Brown       for (i = info.xs; i < info.xs + info.xm; i++) {
63c4762a1bSJed Brown         if (has_glider && i == glider_loc[0] && j == glider_loc[1]) {
64c4762a1bSJed Brown           PetscInt ii, jj;
65c4762a1bSJed Brown           for (ii = -1; ii <= 1; ii++)
669371c9d4SSatish Balay             for (jj = -1; jj <= 1; jj++) x[j + jj][i + ii] = GLIDER[1 - jj][ii + 1];
67c4762a1bSJed Brown         }
68c4762a1bSJed Brown         if (has_blinker && i == blinker_loc[0] && j == blinker_loc[1]) {
69c4762a1bSJed Brown           x[j - 1][i] = 1;
70c4762a1bSJed Brown           x[j][i]     = 1;
71c4762a1bSJed Brown           x[j + 1][i] = 1;
72c4762a1bSJed Brown         }
73c4762a1bSJed Brown       }
74c4762a1bSJed Brown     }
759566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da, Xlocal, &x));
769566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobal(da, Xlocal, ADD_VALUES, Xglobal));
77c4762a1bSJed Brown   }
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* View the initial condition */
809566063dSJacob Faibussowitsch   PetscCall(VecView(Xglobal, viewer));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   { /* Play */
83c4762a1bSJed Brown     PetscInt step;
84c4762a1bSJed Brown 
85c4762a1bSJed Brown     for (step = 0; step < steps; step++) {
86c4762a1bSJed Brown       const PetscScalar **x;
87c4762a1bSJed Brown       PetscScalar       **y;
88c4762a1bSJed Brown       DMDALocalInfo       info;
89c4762a1bSJed Brown       PetscInt            i, j;
90c4762a1bSJed Brown 
919566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocal(da, Xglobal, INSERT_VALUES, Xlocal));
929566063dSJacob Faibussowitsch       PetscCall(DMDAGetLocalInfo(da, &info));
9395b2e421SBarry Smith       PetscCall(DMDAVecGetArrayRead(da, Xlocal, (void *)&x));
949566063dSJacob Faibussowitsch       PetscCall(DMDAVecGetArrayWrite(da, Xglobal, &y));
95c4762a1bSJed Brown       for (j = info.ys; j < info.ys + info.ym; j++) {
96c4762a1bSJed Brown         for (i = info.xs; i < info.xs + info.xm; i++) {
97c4762a1bSJed Brown           PetscInt live_neighbors = 0;
98c4762a1bSJed Brown           live_neighbors += PetscRealPart(x[j - 1][i - 1]) > 0;
99c4762a1bSJed Brown           live_neighbors += PetscRealPart(x[j - 1][i]) > 0;
100c4762a1bSJed Brown           live_neighbors += PetscRealPart(x[j - 1][i + 1]) > 0;
101c4762a1bSJed Brown           live_neighbors += PetscRealPart(x[j][i - 1]) > 0;
102c4762a1bSJed Brown           live_neighbors += PetscRealPart(x[j][i + 1]) > 0;
103c4762a1bSJed Brown           live_neighbors += PetscRealPart(x[j + 1][i - 1]) > 0;
104c4762a1bSJed Brown           live_neighbors += PetscRealPart(x[j + 1][i]) > 0;
105c4762a1bSJed Brown           live_neighbors += PetscRealPart(x[j + 1][i + 1]) > 0;
106c4762a1bSJed Brown           if (PetscRealPart(x[j][i]) > 0) { /* Live cell */
107c4762a1bSJed Brown             switch (live_neighbors) {
108c4762a1bSJed Brown             case 2:
109c4762a1bSJed Brown             case 3:
110c4762a1bSJed Brown               y[j][i] = 1; /* Survive */
111c4762a1bSJed Brown               break;
112d71ae5a4SJacob Faibussowitsch             default:
113d71ae5a4SJacob Faibussowitsch               y[j][i] = 0; /* Death */
114c4762a1bSJed Brown             }
115c4762a1bSJed Brown           } else {                                /* Dead cell */
116c4762a1bSJed Brown             if (live_neighbors == 3) y[j][i] = 1; /* Birth */
117c4762a1bSJed Brown             else y[j][i] = 0;
118c4762a1bSJed Brown           }
119c4762a1bSJed Brown         }
120c4762a1bSJed Brown       }
12195b2e421SBarry Smith       PetscCall(DMDAVecRestoreArrayRead(da, Xlocal, (void *)&x));
1229566063dSJacob Faibussowitsch       PetscCall(DMDAVecRestoreArrayWrite(da, Xglobal, &y));
123c4762a1bSJed Brown       if (step == check_step_alive || step == check_step_dead) {
124c4762a1bSJed Brown         PetscScalar sum;
1259566063dSJacob Faibussowitsch         PetscCall(VecSum(Xglobal, &sum));
126c4762a1bSJed Brown         if (PetscAbsScalar(sum) > 0.1) {
12748a46eb9SPierre Jolivet           if (step == check_step_dead) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Simulation alive at step %" PetscInt_FMT "\n", step));
128c4762a1bSJed Brown         } else if (step == check_step_alive) {
12963a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Simulation dead at step %" PetscInt_FMT "\n", step));
130c4762a1bSJed Brown         }
131c4762a1bSJed Brown       }
13248a46eb9SPierre Jolivet       if (step % viz_interval == 0) PetscCall(VecView(Xglobal, viewer));
133c4762a1bSJed Brown     }
134c4762a1bSJed Brown   }
135c4762a1bSJed Brown 
1369566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
1379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Xglobal));
1389566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Xlocal));
1399566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1409566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
141b122ec5aSJacob Faibussowitsch   return 0;
142c4762a1bSJed Brown }
143c4762a1bSJed Brown 
144c4762a1bSJed Brown /*TEST
145c4762a1bSJed Brown 
146c4762a1bSJed Brown    test:
147c4762a1bSJed Brown       requires: x
148c4762a1bSJed Brown       nsize: 2
149c4762a1bSJed Brown       args: -glider 5,6 -blinker 12,12 -steps 35 -check_step_alive 31 -check_step_dead 32 -da_grid_x 20 -da_grid_y 20 -nox
150*3886731fSPierre Jolivet       output_file: output/empty.out
151c4762a1bSJed Brown 
152c4762a1bSJed Brown TEST*/
153