1*c4762a1bSJed Brown static char help[] = "Demonstrates Conway's Game of Life using a 2d DMDA.\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown /* 4*c4762a1bSJed Brown At each step in time, the following transitions occur: 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown Any live cell with fewer than two live neighbours dies, as if by underpopulation. 7*c4762a1bSJed Brown Any live cell with two or three live neighbours lives on to the next generation. 8*c4762a1bSJed Brown Any live cell with more than three live neighbours dies, as if by overpopulation. 9*c4762a1bSJed Brown Any dead cell with exactly three live neighbours becomes a live cell, as if by reproduction. 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown https://en.wikipedia.org/wiki/Conway%27s_Game_of_Life 12*c4762a1bSJed Brown */ 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown #include <petscdm.h> 15*c4762a1bSJed Brown #include <petscdmda.h> 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown static const int GLIDER[3][3] = { 18*c4762a1bSJed Brown {0, 1, 0}, 19*c4762a1bSJed Brown {0, 1, 1}, 20*c4762a1bSJed Brown {1, 0, 1} 21*c4762a1bSJed Brown }; 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown int main(int argc,char **argv) 24*c4762a1bSJed Brown { 25*c4762a1bSJed Brown PetscErrorCode ierr; 26*c4762a1bSJed Brown DM da; 27*c4762a1bSJed Brown PetscViewer viewer; 28*c4762a1bSJed Brown Vec Xlocal, Xglobal; 29*c4762a1bSJed Brown PetscInt glider_loc[2] = {10, 20}, blinker_loc[2] = {20, 10}, two, steps = 100, viz_interval = 1; 30*c4762a1bSJed Brown PetscInt check_step_alive = -1, check_step_dead = -1; 31*c4762a1bSJed Brown PetscBool has_glider, has_blinker; 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 34*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Conway's Game of Life","");CHKERRQ(ierr); 35*c4762a1bSJed Brown { 36*c4762a1bSJed Brown ierr = PetscOptionsIntArray("-glider","Coordinate at which to center a glider",NULL,glider_loc,(two=2,&two),&has_glider);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = PetscOptionsIntArray("-blinker","Coordinate at which to center a blinker",NULL,blinker_loc,(two=2,&two),&has_blinker);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = PetscOptionsInt("-steps","Number of steps to take",NULL,steps,&steps,NULL);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = PetscOptionsInt("-viz_interval","Vizualization interval",NULL,viz_interval,&viz_interval,NULL);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = PetscOptionsInt("-check_step_alive","Step on which to check that the simulation is alive",NULL,check_step_alive,&check_step_alive,NULL);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = PetscOptionsInt("-check_step_dead","Step on which to check that the simulation is dead",NULL,check_step_dead,&check_step_dead,NULL);CHKERRQ(ierr); 42*c4762a1bSJed Brown } 43*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,NULL,"Life",PETSC_DECIDE,PETSC_DECIDE,1000,1000,&viewer);CHKERRQ(ierr); 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown /* Create distributed array and get vectors */ 48*c4762a1bSJed Brown ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_BOX,30,30,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 51*c4762a1bSJed Brown ierr = DMCreateLocalVector(da,&Xlocal);CHKERRQ(ierr); 52*c4762a1bSJed Brown ierr = DMCreateGlobalVector(da,&Xglobal);CHKERRQ(ierr); 53*c4762a1bSJed Brown 54*c4762a1bSJed Brown { /* Initialize */ 55*c4762a1bSJed Brown DMDALocalInfo info; 56*c4762a1bSJed Brown PetscScalar **x; 57*c4762a1bSJed Brown PetscInt i,j; 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = DMDAVecGetArray(da, Xlocal, &x);CHKERRQ(ierr); 61*c4762a1bSJed Brown for (j=info.ys; j<info.ys+info.ym; j++) { 62*c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 63*c4762a1bSJed Brown if (has_glider && i == glider_loc[0] && j == glider_loc[1]) { 64*c4762a1bSJed Brown PetscInt ii,jj; 65*c4762a1bSJed Brown for (ii=-1; ii<=1; ii++) 66*c4762a1bSJed Brown for (jj=-1; jj<=1; jj++) 67*c4762a1bSJed Brown x[j+jj][i+ii] = GLIDER[1-jj][ii+1]; 68*c4762a1bSJed Brown } 69*c4762a1bSJed Brown if (has_blinker && i == blinker_loc[0] && j == blinker_loc[1]) { 70*c4762a1bSJed Brown x[j-1][i] = 1; 71*c4762a1bSJed Brown x[j][i] = 1; 72*c4762a1bSJed Brown x[j+1][i] = 1; 73*c4762a1bSJed Brown } 74*c4762a1bSJed Brown } 75*c4762a1bSJed Brown } 76*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(da, Xlocal, &x);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = DMLocalToGlobal(da, Xlocal, ADD_VALUES, Xglobal);CHKERRQ(ierr); 78*c4762a1bSJed Brown } 79*c4762a1bSJed Brown 80*c4762a1bSJed Brown /* View the initial condition */ 81*c4762a1bSJed Brown ierr = VecView(Xglobal, viewer);CHKERRQ(ierr); 82*c4762a1bSJed Brown 83*c4762a1bSJed Brown { /* Play */ 84*c4762a1bSJed Brown PetscInt step; 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown for (step=0; step<steps; step++) { 87*c4762a1bSJed Brown const PetscScalar **x; 88*c4762a1bSJed Brown PetscScalar **y; 89*c4762a1bSJed Brown DMDALocalInfo info; 90*c4762a1bSJed Brown PetscInt i,j; 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown ierr = DMGlobalToLocal(da, Xglobal, INSERT_VALUES, Xlocal);CHKERRQ(ierr); 93*c4762a1bSJed Brown ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 94*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(da, Xlocal, &x);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = DMDAVecGetArrayWrite(da, Xglobal, &y);CHKERRQ(ierr); 96*c4762a1bSJed Brown for (j=info.ys; j<info.ys+info.ym; j++) { 97*c4762a1bSJed Brown for (i=info.xs; i<info.xs+info.xm; i++) { 98*c4762a1bSJed Brown PetscInt live_neighbors = 0; 99*c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j-1][i-1]) > 0; 100*c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j-1][i]) > 0; 101*c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j-1][i+1]) > 0; 102*c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j][i-1]) > 0; 103*c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j][i+1]) > 0; 104*c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j+1][i-1]) > 0; 105*c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j+1][i]) > 0; 106*c4762a1bSJed Brown live_neighbors += PetscRealPart(x[j+1][i+1]) > 0; 107*c4762a1bSJed Brown if (PetscRealPart(x[j][i]) > 0) { /* Live cell */ 108*c4762a1bSJed Brown switch (live_neighbors) { 109*c4762a1bSJed Brown case 2: 110*c4762a1bSJed Brown case 3: 111*c4762a1bSJed Brown y[j][i] = 1; /* Survive */ 112*c4762a1bSJed Brown break; 113*c4762a1bSJed Brown default: 114*c4762a1bSJed Brown y[j][i] = 0; /* Death */ 115*c4762a1bSJed Brown } 116*c4762a1bSJed Brown } else { /* Dead cell */ 117*c4762a1bSJed Brown if (live_neighbors == 3) y[j][i] = 1; /* Birth */ 118*c4762a1bSJed Brown else y[j][i] = 0; 119*c4762a1bSJed Brown } 120*c4762a1bSJed Brown } 121*c4762a1bSJed Brown } 122*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(da, Xlocal, &x);CHKERRQ(ierr); 123*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayWrite(da, Xglobal, &y);CHKERRQ(ierr); 124*c4762a1bSJed Brown if (step == check_step_alive || step == check_step_dead) { 125*c4762a1bSJed Brown PetscScalar sum; 126*c4762a1bSJed Brown ierr = VecSum(Xglobal, &sum);CHKERRQ(ierr); 127*c4762a1bSJed Brown if (PetscAbsScalar(sum) > 0.1) { 128*c4762a1bSJed Brown if (step == check_step_dead) { 129*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Simulation alive at step %D\n",step);CHKERRQ(ierr); 130*c4762a1bSJed Brown } 131*c4762a1bSJed Brown } else if (step == check_step_alive) { 132*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Simulation dead at step %D\n",step);CHKERRQ(ierr); 133*c4762a1bSJed Brown } 134*c4762a1bSJed Brown } 135*c4762a1bSJed Brown if (step % viz_interval == 0) { 136*c4762a1bSJed Brown ierr = VecView(Xglobal, viewer);CHKERRQ(ierr); 137*c4762a1bSJed Brown } 138*c4762a1bSJed Brown } 139*c4762a1bSJed Brown } 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 142*c4762a1bSJed Brown ierr = VecDestroy(&Xglobal);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = VecDestroy(&Xlocal);CHKERRQ(ierr); 144*c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 145*c4762a1bSJed Brown ierr = PetscFinalize(); 146*c4762a1bSJed Brown return ierr; 147*c4762a1bSJed Brown } 148*c4762a1bSJed Brown 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown /*TEST 151*c4762a1bSJed Brown 152*c4762a1bSJed Brown test: 153*c4762a1bSJed Brown requires: x 154*c4762a1bSJed Brown nsize: 2 155*c4762a1bSJed 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 156*c4762a1bSJed Brown 157*c4762a1bSJed Brown TEST*/ 158