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