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 23c4762a1bSJed Brown int main(int argc,char **argv) 24c4762a1bSJed Brown { 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 32*327415f7SBarry Smith PetscFunctionBeginUser; 339566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,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)); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-viz_interval","Vizualization 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++) 66c4762a1bSJed Brown for (jj=-1; jj<=1; jj++) 67c4762a1bSJed Brown x[j+jj][i+ii] = GLIDER[1-jj][ii+1]; 68c4762a1bSJed Brown } 69c4762a1bSJed Brown if (has_blinker && i == blinker_loc[0] && j == blinker_loc[1]) { 70c4762a1bSJed Brown x[j-1][i] = 1; 71c4762a1bSJed Brown x[j][i] = 1; 72c4762a1bSJed Brown x[j+1][i] = 1; 73c4762a1bSJed Brown } 74c4762a1bSJed Brown } 75c4762a1bSJed Brown } 769566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xlocal, &x)); 779566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobal(da, Xlocal, ADD_VALUES, Xglobal)); 78c4762a1bSJed Brown } 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* View the initial condition */ 819566063dSJacob Faibussowitsch PetscCall(VecView(Xglobal, viewer)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown { /* Play */ 84c4762a1bSJed Brown PetscInt step; 85c4762a1bSJed Brown 86c4762a1bSJed Brown for (step=0; step<steps; step++) { 87c4762a1bSJed Brown const PetscScalar **x; 88c4762a1bSJed Brown PetscScalar **y; 89c4762a1bSJed Brown DMDALocalInfo info; 90c4762a1bSJed Brown PetscInt i,j; 91c4762a1bSJed Brown 929566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocal(da, Xglobal, INSERT_VALUES, Xlocal)); 939566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da,&info)); 9495b2e421SBarry Smith PetscCall(DMDAVecGetArrayRead(da, Xlocal, (void*)&x)); 959566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayWrite(da, Xglobal, &y)); 96c4762a1bSJed Brown for (j=info.ys; j<info.ys+info.ym; j++) { 97c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 98c4762a1bSJed Brown PetscInt live_neighbors = 0; 99c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j-1][i-1]) > 0; 100c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j-1][i]) > 0; 101c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j-1][i+1]) > 0; 102c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j][i-1]) > 0; 103c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j][i+1]) > 0; 104c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j+1][i-1]) > 0; 105c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j+1][i]) > 0; 106c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j+1][i+1]) > 0; 107c4762a1bSJed Brown if (PetscRealPart(x[j][i]) > 0) { /* Live cell */ 108c4762a1bSJed Brown switch (live_neighbors) { 109c4762a1bSJed Brown case 2: 110c4762a1bSJed Brown case 3: 111c4762a1bSJed Brown y[j][i] = 1; /* Survive */ 112c4762a1bSJed Brown break; 113c4762a1bSJed Brown default: 114c4762a1bSJed Brown y[j][i] = 0; /* Death */ 115c4762a1bSJed Brown } 116c4762a1bSJed Brown } else { /* Dead cell */ 117c4762a1bSJed Brown if (live_neighbors == 3) y[j][i] = 1; /* Birth */ 118c4762a1bSJed Brown else y[j][i] = 0; 119c4762a1bSJed Brown } 120c4762a1bSJed Brown } 121c4762a1bSJed Brown } 12295b2e421SBarry Smith PetscCall(DMDAVecRestoreArrayRead(da, Xlocal, (void*)&x)); 1239566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayWrite(da, Xglobal, &y)); 124c4762a1bSJed Brown if (step == check_step_alive || step == check_step_dead) { 125c4762a1bSJed Brown PetscScalar sum; 1269566063dSJacob Faibussowitsch PetscCall(VecSum(Xglobal, &sum)); 127c4762a1bSJed Brown if (PetscAbsScalar(sum) > 0.1) { 128c4762a1bSJed Brown if (step == check_step_dead) { 12963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Simulation alive at step %" PetscInt_FMT "\n",step)); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown } else if (step == check_step_alive) { 13263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Simulation dead at step %" PetscInt_FMT "\n",step)); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown } 135c4762a1bSJed Brown if (step % viz_interval == 0) { 1369566063dSJacob Faibussowitsch PetscCall(VecView(Xglobal, viewer)); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown } 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xglobal)); 1439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xlocal)); 1449566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1459566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 146b122ec5aSJacob Faibussowitsch return 0; 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149c4762a1bSJed Brown /*TEST 150c4762a1bSJed Brown 151c4762a1bSJed Brown test: 152c4762a1bSJed Brown requires: x 153c4762a1bSJed Brown nsize: 2 154c4762a1bSJed 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 155c4762a1bSJed Brown 156c4762a1bSJed Brown TEST*/ 157