1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests various 2-dimensional DMDA routines.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdm.h> 5c4762a1bSJed Brown #include <petscdmda.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown int main(int argc,char **argv) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown PetscMPIInt rank; 10c4762a1bSJed Brown PetscErrorCode ierr; 11c4762a1bSJed Brown PetscInt M = 10,N = 8,m = PETSC_DECIDE; 12c4762a1bSJed Brown PetscInt s =2,w=2,n = PETSC_DECIDE,nloc,l,i,j,kk; 13c4762a1bSJed Brown PetscInt Xs,Xm,Ys,Ym,iloc,*iglobal; 14c4762a1bSJed Brown const PetscInt *ltog; 15c4762a1bSJed Brown PetscInt *lx = NULL,*ly = NULL; 16c4762a1bSJed Brown PetscBool testorder = PETSC_FALSE,flg; 17c4762a1bSJed Brown DMBoundaryType bx = DM_BOUNDARY_NONE,by= DM_BOUNDARY_NONE; 18c4762a1bSJed Brown DM da; 19c4762a1bSJed Brown PetscViewer viewer; 20c4762a1bSJed Brown Vec local,global; 21c4762a1bSJed Brown PetscScalar value; 22c4762a1bSJed Brown DMDAStencilType st = DMDA_STENCIL_BOX; 23c4762a1bSJed Brown AO ao; 24c4762a1bSJed Brown 25c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",300,0,400,400,&viewer)); 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* Readoptions */ 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-NX",&M,NULL)); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-NY",&N,NULL)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-s",&s,NULL)); 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-w",&w,NULL)); 35c4762a1bSJed Brown 36c4762a1bSJed Brown flg = PETSC_FALSE; 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-xperiodic",&flg,NULL)); if (flg) bx = DM_BOUNDARY_PERIODIC; 38c4762a1bSJed Brown flg = PETSC_FALSE; 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-yperiodic",&flg,NULL)); if (flg) by = DM_BOUNDARY_PERIODIC; 40c4762a1bSJed Brown flg = PETSC_FALSE; 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-xghosted",&flg,NULL)); if (flg) bx = DM_BOUNDARY_GHOSTED; 42c4762a1bSJed Brown flg = PETSC_FALSE; 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-yghosted",&flg,NULL)); if (flg) by = DM_BOUNDARY_GHOSTED; 44c4762a1bSJed Brown flg = PETSC_FALSE; 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-star",&flg,NULL)); if (flg) st = DMDA_STENCIL_STAR; 46c4762a1bSJed Brown flg = PETSC_FALSE; 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-box",&flg,NULL)); if (flg) st = DMDA_STENCIL_BOX; 48c4762a1bSJed Brown flg = PETSC_FALSE; 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-testorder",&testorder,NULL)); 50c4762a1bSJed Brown /* 51c4762a1bSJed Brown Test putting two nodes in x and y on each processor, exact last processor 52c4762a1bSJed Brown in x and y gets the rest. 53c4762a1bSJed Brown */ 54c4762a1bSJed Brown flg = PETSC_FALSE; 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-distribute",&flg,NULL)); 56c4762a1bSJed Brown if (flg) { 572c71b3e2SJacob Faibussowitsch PetscCheckFalse(m == PETSC_DECIDE,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Must set -m option with -distribute option"); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(m,&lx)); 59c4762a1bSJed Brown for (i=0; i<m-1; i++) { lx[i] = 4;} 60c4762a1bSJed Brown lx[m-1] = M - 4*(m-1); 612c71b3e2SJacob Faibussowitsch PetscCheckFalse(n == PETSC_DECIDE,PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Must set -n option with -distribute option"); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&ly)); 63c4762a1bSJed Brown for (i=0; i<n-1; i++) { ly[i] = 2;} 64c4762a1bSJed Brown ly[n-1] = N - 2*(n-1); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* Create distributed array and get vectors */ 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,bx,by,st,M,N,m,n,w,s,lx,ly,&da)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(lx)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ly)); 73c4762a1bSJed Brown 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMView(da,viewer)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&global)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(da,&local)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Set global vector; send ghost points to local vectors */ 79c4762a1bSJed Brown value = 1; 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(global,value)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local)); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Scale local vectors according to processor rank; pass to global vector */ 85*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 86c4762a1bSJed Brown value = rank; 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(local,value)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(da,local,INSERT_VALUES,global)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(da,local,INSERT_VALUES,global)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown if (!testorder) { /* turn off printing when testing ordering mappings */ 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nGlobal Vectors:\n")); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(global,PETSC_VIEWER_STDOUT_WORLD)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n\n")); 95c4762a1bSJed Brown } 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Send ghost points to local vectors */ 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local)); 100c4762a1bSJed Brown 101c4762a1bSJed Brown flg = PETSC_FALSE; 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-local_print",&flg,NULL)); 103c4762a1bSJed Brown if (flg) { 104c4762a1bSJed Brown PetscViewer sviewer; 105c4762a1bSJed Brown 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"\nLocal Vector: processor %d\n",rank)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(local,sviewer)); 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 115a5b23f4aSJose E. Roman /* Tests mappings between application/PETSc orderings */ 116c4762a1bSJed Brown if (testorder) { 117c4762a1bSJed Brown ISLocalToGlobalMapping ltogm; 118c4762a1bSJed Brown 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalToGlobalMapping(da,<ogm)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISLocalToGlobalMappingGetSize(ltogm,&nloc)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISLocalToGlobalMappingGetIndices(ltogm,<og)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetGhostCorners(da,&Xs,&Ys,NULL,&Xm,&Ym,NULL)); 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetAO(da,&ao)); 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nloc,&iglobal)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* Set iglobal to be global indices for each processor's local and ghost nodes, 127c4762a1bSJed Brown using the DMDA ordering of grid points */ 128c4762a1bSJed Brown kk = 0; 129c4762a1bSJed Brown for (j=Ys; j<Ys+Ym; j++) { 130c4762a1bSJed Brown for (i=Xs; i<Xs+Xm; i++) { 131c4762a1bSJed Brown iloc = w*((j-Ys)*Xm + i-Xs); 132c4762a1bSJed Brown for (l=0; l<w; l++) { 133c4762a1bSJed Brown iglobal[kk++] = ltog[iloc+l]; 134c4762a1bSJed Brown } 135c4762a1bSJed Brown } 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Map this to the application ordering (which for DMDAs is just the natural ordering 139c4762a1bSJed Brown that would be used for 1 processor, numbering most rapidly by x, then y) */ 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(AOPetscToApplication(ao,nloc,iglobal)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* Then map the application ordering back to the PETSc DMDA ordering */ 143*5f80ce2aSJacob Faibussowitsch CHKERRQ(AOApplicationToPetsc(ao,nloc,iglobal)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Verify the mappings */ 146c4762a1bSJed Brown kk=0; 147c4762a1bSJed Brown for (j=Ys; j<Ys+Ym; j++) { 148c4762a1bSJed Brown for (i=Xs; i<Xs+Xm; i++) { 149c4762a1bSJed Brown iloc = w*((j-Ys)*Xm + i-Xs); 150c4762a1bSJed Brown for (l=0; l<w; l++) { 151c4762a1bSJed Brown if (iglobal[kk] != ltog[iloc+l]) { 152*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPrintf(PETSC_COMM_SELF,stdout,"[%d] Problem with mapping: j=%D, i=%D, l=%D, petsc1=%D, petsc2=%D\n",rank,j,i,l,ltog[iloc+l],iglobal[kk])); 153c4762a1bSJed Brown } 154c4762a1bSJed Brown kk++; 155c4762a1bSJed Brown } 156c4762a1bSJed Brown } 157c4762a1bSJed Brown } 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(iglobal)); 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISLocalToGlobalMappingRestoreIndices(ltogm,<og)); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* Free memory */ 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 164*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&local)); 165*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&global)); 166*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 167c4762a1bSJed Brown 168c4762a1bSJed Brown ierr = PetscFinalize(); 169c4762a1bSJed Brown return ierr; 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 172c4762a1bSJed Brown /*TEST 173c4762a1bSJed Brown 174c4762a1bSJed Brown test: 175c4762a1bSJed Brown nsize: 4 176c4762a1bSJed Brown args: -nox 177c4762a1bSJed Brown filter: grep -v -i Object 178c4762a1bSJed Brown requires: x 179c4762a1bSJed Brown 180c4762a1bSJed Brown test: 181c4762a1bSJed Brown suffix: 2 182c4762a1bSJed Brown args: -testorder -nox 183c4762a1bSJed Brown requires: x 184c4762a1bSJed Brown 185c4762a1bSJed Brown TEST*/ 186