1c4762a1bSJed Brown static const char help[] = "Test DMDAGetOwnershipRanges()\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdm.h> 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 PetscViewer vw; 11c4762a1bSJed Brown PetscInt dim = 2,m,n,p; 12c4762a1bSJed Brown const PetscInt *lx,*ly,*lz; 13c4762a1bSJed Brown PetscMPIInt rank; 14c4762a1bSJed Brown 15c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 16*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,0,"-dim",&dim,0)); 17c4762a1bSJed Brown switch (dim) { 18c4762a1bSJed Brown case 2: 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR, 3,5,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da)); 20c4762a1bSJed Brown break; 21c4762a1bSJed Brown case 3: 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR, 3,5,7,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,NULL,&da)); 23c4762a1bSJed Brown break; 2498921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"No support for %D dimensions",dim); 25c4762a1bSJed Brown } 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da, 0, 0,0,0, &m,&n,&p, 0,0, 0,0,0,0)); 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetOwnershipRanges(da,&lx,&ly,&lz)); 30*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 31c4762a1bSJed Brown 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&vw)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(vw,"[%d] lx ly%s\n",rank,dim>2 ? " lz" : "")); 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscIntView(m,lx,vw)); 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscIntView(n,ly,vw)); 36*5f80ce2aSJacob Faibussowitsch if (dim > 2) CHKERRQ(PetscIntView(n,lz,vw)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&vw)); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 39c4762a1bSJed Brown 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 41c4762a1bSJed Brown ierr = PetscFinalize(); 42c4762a1bSJed Brown return ierr; 43c4762a1bSJed Brown } 44c4762a1bSJed Brown 45c4762a1bSJed Brown /*TEST 46c4762a1bSJed Brown 47c4762a1bSJed Brown test: 48c4762a1bSJed Brown nsize: 12 49c4762a1bSJed Brown args: -dm_view -dim 3 -da_grid_x 11 -da_grid_y 5 -da_grid_z 7 50c4762a1bSJed Brown 51c4762a1bSJed Brown TEST*/ 52