1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests periodic boundary conditions for DMDA1d with periodic boundary conditions.\n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscdmda.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown int main(int argc,char **argv) 7*c4762a1bSJed Brown { 8*c4762a1bSJed Brown PetscErrorCode ierr; 9*c4762a1bSJed Brown DM da; 10*c4762a1bSJed Brown Mat A; 11*c4762a1bSJed Brown const PetscInt dfill[4] = {0,1,0,1},ofill[4] = {0,1,1,0}; 12*c4762a1bSJed Brown 13*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 14*c4762a1bSJed Brown ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,7,2,1,NULL,&da);CHKERRQ(ierr); 15*c4762a1bSJed Brown ierr = DMSetFromOptions(da);CHKERRQ(ierr); 16*c4762a1bSJed Brown ierr = DMDASetBlockFills(da,dfill,ofill);CHKERRQ(ierr); 17*c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 18*c4762a1bSJed Brown ierr = DMCreateMatrix(da,&A);CHKERRQ(ierr); 19*c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 20*c4762a1bSJed Brown ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 21*c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 22*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 23*c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = PetscFinalize(); 25*c4762a1bSJed Brown return ierr; 26*c4762a1bSJed Brown } 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown /*TEST 30*c4762a1bSJed Brown 31*c4762a1bSJed Brown test: 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown test: 34*c4762a1bSJed Brown suffix: 2 35*c4762a1bSJed Brown nsize: 2 36*c4762a1bSJed Brown 37*c4762a1bSJed Brown TEST*/ 38