147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 347c6ae99SBarry Smith Code for manipulating distributed regular arrays in parallel. 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 64035e84dSBarry Smith #include <petsc-private/dmdaimpl.h> /*I "petscdmda.h" I*/ 747c6ae99SBarry Smith 847c6ae99SBarry Smith #undef __FUNCT__ 96636e97aSMatthew G Knepley #define __FUNCT__ "DMCreateCoordinateDM_DA" 106636e97aSMatthew G Knepley PetscErrorCode DMCreateCoordinateDM_DA(DM dm, DM *cdm) 1147c6ae99SBarry Smith { 1247c6ae99SBarry Smith PetscErrorCode ierr; 13e1d31a9fSPeter Brune DM_DA *da = (DM_DA*) dm->data; 14*3f2d3b52SPeter Brune PetscFunctionBegin; 15e1d31a9fSPeter Brune ierr = DMDAGetReducedDMDA(dm,da->dim,cdm);CHKERRQ(ierr); 1647c6ae99SBarry Smith PetscFunctionReturn(0); 1747c6ae99SBarry Smith } 1847c6ae99SBarry Smith 1947c6ae99SBarry Smith #undef __FUNCT__ 20aa219208SBarry Smith #define __FUNCT__ "DMDASetFieldName" 2147c6ae99SBarry Smith /*@C 22aa219208SBarry Smith DMDASetFieldName - Sets the names of individual field components in multicomponent 23aa219208SBarry Smith vectors associated with a DMDA. 2447c6ae99SBarry Smith 2547c6ae99SBarry Smith Not Collective 2647c6ae99SBarry Smith 2747c6ae99SBarry Smith Input Parameters: 2847c6ae99SBarry Smith + da - the distributed array 29aa219208SBarry Smith . nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the 30aa219208SBarry Smith number of degrees of freedom per node within the DMDA 3147c6ae99SBarry Smith - names - the name of the field (component) 3247c6ae99SBarry Smith 3347c6ae99SBarry Smith Level: intermediate 3447c6ae99SBarry Smith 3547c6ae99SBarry Smith .keywords: distributed array, get, component name 3647c6ae99SBarry Smith 37109c9344SBarry Smith .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName() 3847c6ae99SBarry Smith @*/ 397087cfbeSBarry Smith PetscErrorCode DMDASetFieldName(DM da,PetscInt nf,const char name[]) 4047c6ae99SBarry Smith { 4147c6ae99SBarry Smith PetscErrorCode ierr; 4247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 4347c6ae99SBarry Smith 4447c6ae99SBarry Smith PetscFunctionBegin; 4547c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 4647c6ae99SBarry Smith if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf); 47c31cb41cSBarry Smith ierr = PetscFree(dd->fieldname[nf]);CHKERRQ(ierr); 4847c6ae99SBarry Smith ierr = PetscStrallocpy(name,&dd->fieldname[nf]);CHKERRQ(ierr); 4947c6ae99SBarry Smith PetscFunctionReturn(0); 5047c6ae99SBarry Smith } 5147c6ae99SBarry Smith 5247c6ae99SBarry Smith #undef __FUNCT__ 53aa219208SBarry Smith #define __FUNCT__ "DMDAGetFieldName" 5447c6ae99SBarry Smith /*@C 55aa219208SBarry Smith DMDAGetFieldName - Gets the names of individual field components in multicomponent 56aa219208SBarry Smith vectors associated with a DMDA. 5747c6ae99SBarry Smith 5847c6ae99SBarry Smith Not Collective 5947c6ae99SBarry Smith 6047c6ae99SBarry Smith Input Parameter: 6147c6ae99SBarry Smith + da - the distributed array 62aa219208SBarry Smith - nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the 63aa219208SBarry Smith number of degrees of freedom per node within the DMDA 6447c6ae99SBarry Smith 6547c6ae99SBarry Smith Output Parameter: 6647c6ae99SBarry Smith . names - the name of the field (component) 6747c6ae99SBarry Smith 6847c6ae99SBarry Smith Level: intermediate 6947c6ae99SBarry Smith 7047c6ae99SBarry Smith .keywords: distributed array, get, component name 7147c6ae99SBarry Smith 72109c9344SBarry Smith .seealso: DMDASetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName() 7347c6ae99SBarry Smith @*/ 747087cfbeSBarry Smith PetscErrorCode DMDAGetFieldName(DM da,PetscInt nf,const char **name) 7547c6ae99SBarry Smith { 7647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 7747c6ae99SBarry Smith 7847c6ae99SBarry Smith PetscFunctionBegin; 7947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 8047c6ae99SBarry Smith PetscValidPointer(name,3); 8147c6ae99SBarry Smith if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf); 8247c6ae99SBarry Smith *name = dd->fieldname[nf]; 8347c6ae99SBarry Smith PetscFunctionReturn(0); 8447c6ae99SBarry Smith } 8547c6ae99SBarry Smith 8647c6ae99SBarry Smith #undef __FUNCT__ 87109c9344SBarry Smith #define __FUNCT__ "DMDASetCoordinateName" 88109c9344SBarry Smith /*@C 89109c9344SBarry Smith DMDASetCoordinateName - Sets the name of the coordinate directions associated with a DMDA, for example "x" or "y" 90109c9344SBarry Smith 91109c9344SBarry Smith Not Collective 92109c9344SBarry Smith 93109c9344SBarry Smith Input Parameters: 94109c9344SBarry Smith + da - the distributed array 95109c9344SBarry Smith . nf - coordinate number for the DMDA (0, 1, ... dim-1), 96109c9344SBarry Smith - name - the name of the coordinate 97109c9344SBarry Smith 98109c9344SBarry Smith Level: intermediate 99109c9344SBarry Smith 100109c9344SBarry Smith .keywords: distributed array, get, component name 101109c9344SBarry Smith 102109c9344SBarry Smith .seealso: DMDAGetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName() 103109c9344SBarry Smith @*/ 104109c9344SBarry Smith PetscErrorCode DMDASetCoordinateName(DM da,PetscInt nf,const char name[]) 105109c9344SBarry Smith { 106109c9344SBarry Smith PetscErrorCode ierr; 107109c9344SBarry Smith DM_DA *dd = (DM_DA*)da->data; 108109c9344SBarry Smith 109109c9344SBarry Smith PetscFunctionBegin; 110109c9344SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 111109c9344SBarry Smith if (nf < 0 || nf >= dd->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf); 112109c9344SBarry Smith ierr = PetscFree(dd->coordinatename[nf]);CHKERRQ(ierr); 113109c9344SBarry Smith ierr = PetscStrallocpy(name,&dd->coordinatename[nf]);CHKERRQ(ierr); 114109c9344SBarry Smith PetscFunctionReturn(0); 115109c9344SBarry Smith } 116109c9344SBarry Smith 117109c9344SBarry Smith #undef __FUNCT__ 118109c9344SBarry Smith #define __FUNCT__ "DMDAGetCoordinateName" 119109c9344SBarry Smith /*@C 120109c9344SBarry Smith DMDAGetCoordinateName - Gets the name of a coodinate direction associated with a DMDA. 121109c9344SBarry Smith 122109c9344SBarry Smith Not Collective 123109c9344SBarry Smith 124109c9344SBarry Smith Input Parameter: 125109c9344SBarry Smith + da - the distributed array 126109c9344SBarry Smith - nf - number for the DMDA (0, 1, ... dim-1) 127109c9344SBarry Smith 128109c9344SBarry Smith Output Parameter: 129109c9344SBarry Smith . names - the name of the coordinate direction 130109c9344SBarry Smith 131109c9344SBarry Smith Level: intermediate 132109c9344SBarry Smith 133109c9344SBarry Smith .keywords: distributed array, get, component name 134109c9344SBarry Smith 135109c9344SBarry Smith .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName() 136109c9344SBarry Smith @*/ 137109c9344SBarry Smith PetscErrorCode DMDAGetCoordinateName(DM da,PetscInt nf,const char **name) 138109c9344SBarry Smith { 139109c9344SBarry Smith DM_DA *dd = (DM_DA*)da->data; 140109c9344SBarry Smith 141109c9344SBarry Smith PetscFunctionBegin; 142109c9344SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 143109c9344SBarry Smith PetscValidPointer(name,3); 144109c9344SBarry Smith if (nf < 0 || nf >= dd->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf); 145109c9344SBarry Smith *name = dd->coordinatename[nf]; 146109c9344SBarry Smith PetscFunctionReturn(0); 147109c9344SBarry Smith } 148109c9344SBarry Smith 149109c9344SBarry Smith #undef __FUNCT__ 150aa219208SBarry Smith #define __FUNCT__ "DMDAGetCorners" 15147c6ae99SBarry Smith /*@ 152aa219208SBarry Smith DMDAGetCorners - Returns the global (x,y,z) indices of the lower left 15347c6ae99SBarry Smith corner of the local region, excluding ghost points. 15447c6ae99SBarry Smith 15547c6ae99SBarry Smith Not Collective 15647c6ae99SBarry Smith 15747c6ae99SBarry Smith Input Parameter: 15847c6ae99SBarry Smith . da - the distributed array 15947c6ae99SBarry Smith 16047c6ae99SBarry Smith Output Parameters: 16147c6ae99SBarry Smith + x,y,z - the corner indices (where y and z are optional; these are used 16247c6ae99SBarry Smith for 2D and 3D problems) 16347c6ae99SBarry Smith - m,n,p - widths in the corresponding directions (where n and p are optional; 16447c6ae99SBarry Smith these are used for 2D and 3D problems) 16547c6ae99SBarry Smith 16647c6ae99SBarry Smith Note: 16747c6ae99SBarry Smith The corner information is independent of the number of degrees of 168aa219208SBarry Smith freedom per node set with the DMDACreateXX() routine. Thus the x, y, z, and 16947c6ae99SBarry Smith m, n, p can be thought of as coordinates on a logical grid, where each 17047c6ae99SBarry Smith grid point has (potentially) several degrees of freedom. 1710298fd71SBarry Smith Any of y, z, n, and p can be passed in as NULL if not needed. 17247c6ae99SBarry Smith 17347c6ae99SBarry Smith Level: beginner 17447c6ae99SBarry Smith 17547c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices 17647c6ae99SBarry Smith 177aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetOwnershipRanges() 17847c6ae99SBarry Smith @*/ 1797087cfbeSBarry Smith PetscErrorCode DMDAGetCorners(DM da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p) 18047c6ae99SBarry Smith { 18147c6ae99SBarry Smith PetscInt w; 18247c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 18347c6ae99SBarry Smith 18447c6ae99SBarry Smith PetscFunctionBegin; 18547c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 18647c6ae99SBarry Smith /* since the xs, xe ... have all been multiplied by the number of degrees 18747c6ae99SBarry Smith of freedom per cell, w = dd->w, we divide that out before returning.*/ 18847c6ae99SBarry Smith w = dd->w; 189d886c4f4SPeter Brune if (x) *x = dd->xs/w + dd->xo; if (m) *m = (dd->xe - dd->xs)/w; 19047c6ae99SBarry Smith /* the y and z have NOT been multiplied by w */ 191d886c4f4SPeter Brune if (y) *y = dd->ys + dd->yo; if (n) *n = (dd->ye - dd->ys); 192d886c4f4SPeter Brune if (z) *z = dd->zs + dd->zo; if (p) *p = (dd->ze - dd->zs); 19347c6ae99SBarry Smith PetscFunctionReturn(0); 19447c6ae99SBarry Smith } 19547c6ae99SBarry Smith 19647c6ae99SBarry Smith #undef __FUNCT__ 197aa219208SBarry Smith #define __FUNCT__ "DMDAGetLocalBoundingBox" 19847c6ae99SBarry Smith /*@ 199aa219208SBarry Smith DMDAGetLocalBoundingBox - Returns the local bounding box for the DMDA. 20047c6ae99SBarry Smith 20147c6ae99SBarry Smith Not Collective 20247c6ae99SBarry Smith 20347c6ae99SBarry Smith Input Parameter: 20447c6ae99SBarry Smith . da - the distributed array 20547c6ae99SBarry Smith 20647c6ae99SBarry Smith Output Parameters: 20747c6ae99SBarry Smith + lmin - local minimum coordinates (length dim, optional) 20847c6ae99SBarry Smith - lmax - local maximim coordinates (length dim, optional) 20947c6ae99SBarry Smith 21047c6ae99SBarry Smith Level: beginner 21147c6ae99SBarry Smith 21247c6ae99SBarry Smith .keywords: distributed array, get, coordinates 21347c6ae99SBarry Smith 2142150357eSBarry Smith .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetBoundingBox() 21547c6ae99SBarry Smith @*/ 2167087cfbeSBarry Smith PetscErrorCode DMDAGetLocalBoundingBox(DM da,PetscReal lmin[],PetscReal lmax[]) 21747c6ae99SBarry Smith { 21847c6ae99SBarry Smith PetscErrorCode ierr; 2190298fd71SBarry Smith Vec coords = NULL; 22047c6ae99SBarry Smith PetscInt dim,i,j; 22147c6ae99SBarry Smith const PetscScalar *local_coords; 222ea345e14SBarry Smith PetscReal min[3]={PETSC_MAX_REAL,PETSC_MAX_REAL,PETSC_MAX_REAL},max[3]={PETSC_MIN_REAL,PETSC_MIN_REAL,PETSC_MIN_REAL}; 22347c6ae99SBarry Smith PetscInt N,Ni; 22447c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 22547c6ae99SBarry Smith 22647c6ae99SBarry Smith PetscFunctionBegin; 22747c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 22847c6ae99SBarry Smith dim = dd->dim; 2296636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&coords);CHKERRQ(ierr); 2307324c66bSJed Brown if (coords) { 23147c6ae99SBarry Smith ierr = VecGetArrayRead(coords,&local_coords);CHKERRQ(ierr); 23247c6ae99SBarry Smith ierr = VecGetLocalSize(coords,&N);CHKERRQ(ierr); 23347c6ae99SBarry Smith Ni = N/dim; 23447c6ae99SBarry Smith for (i=0; i<Ni; i++) { 2357324c66bSJed Brown for (j=0; j<3; j++) { 2367324c66bSJed Brown min[j] = j < dim ? PetscMin(min[j],PetscRealPart(local_coords[i*dim+j])) : 0; 2372197688cSJed Brown max[j] = j < dim ? PetscMax(max[j],PetscRealPart(local_coords[i*dim+j])) : 0; 23847c6ae99SBarry Smith } 23947c6ae99SBarry Smith } 24047c6ae99SBarry Smith ierr = VecRestoreArrayRead(coords,&local_coords);CHKERRQ(ierr); 2417324c66bSJed Brown } else { /* Just use grid indices */ 2427324c66bSJed Brown DMDALocalInfo info; 2437324c66bSJed Brown ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 2447324c66bSJed Brown min[0] = info.xs; 2457324c66bSJed Brown min[1] = info.ys; 2467324c66bSJed Brown min[2] = info.zs; 2477324c66bSJed Brown max[0] = info.xs + info.xm-1; 2487324c66bSJed Brown max[1] = info.ys + info.ym-1; 2497324c66bSJed Brown max[2] = info.zs + info.zm-1; 2507324c66bSJed Brown } 25147c6ae99SBarry Smith if (lmin) {ierr = PetscMemcpy(lmin,min,dim*sizeof(PetscReal));CHKERRQ(ierr);} 25247c6ae99SBarry Smith if (lmax) {ierr = PetscMemcpy(lmax,max,dim*sizeof(PetscReal));CHKERRQ(ierr);} 25347c6ae99SBarry Smith PetscFunctionReturn(0); 25447c6ae99SBarry Smith } 25547c6ae99SBarry Smith 25647c6ae99SBarry Smith #undef __FUNCT__ 257aa219208SBarry Smith #define __FUNCT__ "DMDAGetBoundingBox" 25847c6ae99SBarry Smith /*@ 259aa219208SBarry Smith DMDAGetBoundingBox - Returns the global bounding box for the DMDA. 26047c6ae99SBarry Smith 261aa219208SBarry Smith Collective on DMDA 26247c6ae99SBarry Smith 26347c6ae99SBarry Smith Input Parameter: 26447c6ae99SBarry Smith . da - the distributed array 26547c6ae99SBarry Smith 26647c6ae99SBarry Smith Output Parameters: 26747c6ae99SBarry Smith + gmin - global minimum coordinates (length dim, optional) 26847c6ae99SBarry Smith - gmax - global maximim coordinates (length dim, optional) 26947c6ae99SBarry Smith 27047c6ae99SBarry Smith Level: beginner 27147c6ae99SBarry Smith 27247c6ae99SBarry Smith .keywords: distributed array, get, coordinates 27347c6ae99SBarry Smith 2742150357eSBarry Smith .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetLocalBoundingBox() 27547c6ae99SBarry Smith @*/ 2767087cfbeSBarry Smith PetscErrorCode DMDAGetBoundingBox(DM da,PetscReal gmin[],PetscReal gmax[]) 27747c6ae99SBarry Smith { 27847c6ae99SBarry Smith PetscErrorCode ierr; 27947c6ae99SBarry Smith PetscMPIInt count; 28047c6ae99SBarry Smith PetscReal lmin[3],lmax[3]; 28147c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 28247c6ae99SBarry Smith 28347c6ae99SBarry Smith PetscFunctionBegin; 28447c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 2854dc2109aSBarry Smith ierr = PetscMPIIntCast(dd->dim,&count);CHKERRQ(ierr); 286aa219208SBarry Smith ierr = DMDAGetLocalBoundingBox(da,lmin,lmax);CHKERRQ(ierr); 287ce94432eSBarry Smith if (gmin) {ierr = MPI_Allreduce(lmin,gmin,count,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);} 288ce94432eSBarry Smith if (gmax) {ierr = MPI_Allreduce(lmax,gmax,count,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)da));CHKERRQ(ierr);} 28947c6ae99SBarry Smith PetscFunctionReturn(0); 29047c6ae99SBarry Smith } 291bc2bf880SBarry Smith 292bc2bf880SBarry Smith #undef __FUNCT__ 293db05f41bSBarry Smith #define __FUNCT__ "DMDAGetReducedDMDA" 294bc2bf880SBarry Smith /*@ 295db05f41bSBarry Smith DMDAGetReducedDMDA - Gets the DMDA with the same layout but with fewer or more fields 296bc2bf880SBarry Smith 297bc2bf880SBarry Smith Collective on DMDA 298bc2bf880SBarry Smith 299bc2bf880SBarry Smith Input Parameter: 300bc2bf880SBarry Smith + da - the distributed array 301bc2bf880SBarry Smith . nfields - number of fields in new DMDA 302bc2bf880SBarry Smith 303bc2bf880SBarry Smith Output Parameter: 304bc2bf880SBarry Smith . nda - the new DMDA 305bc2bf880SBarry Smith 306bc2bf880SBarry Smith Level: intermediate 307bc2bf880SBarry Smith 308bc2bf880SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates 309bc2bf880SBarry Smith 3102150357eSBarry Smith .seealso: DMDAGetGhostCorners(), DMSetCoordinates(), DMDASetUniformCoordinates(), DMGetCoordinates(), DMDAGetGhostedCoordinates() 311bc2bf880SBarry Smith @*/ 312db05f41bSBarry Smith PetscErrorCode DMDAGetReducedDMDA(DM da,PetscInt nfields,DM *nda) 313bc2bf880SBarry Smith { 314bc2bf880SBarry Smith PetscErrorCode ierr; 315bc2bf880SBarry Smith DM_DA *dd = (DM_DA*)da->data; 31695c13181SPeter Brune PetscInt s,m,n,p,M,N,P,dim,Mo,No,Po; 317320964c4SBlaise Bourdin const PetscInt *lx,*ly,*lz; 318bc2bf880SBarry Smith DMDABoundaryType bx,by,bz; 319320964c4SBlaise Bourdin DMDAStencilType stencil_type; 32095c13181SPeter Brune PetscInt ox,oy,oz; 32195c13181SPeter Brune PetscInt cl,rl; 322320964c4SBlaise Bourdin 323320964c4SBlaise Bourdin PetscFunctionBegin; 32495c13181SPeter Brune dim = dd->dim; 32595c13181SPeter Brune M = dd->M; 32695c13181SPeter Brune N = dd->N; 32795c13181SPeter Brune P = dd->P; 32895c13181SPeter Brune m = dd->m; 32995c13181SPeter Brune n = dd->n; 33095c13181SPeter Brune p = dd->p; 33195c13181SPeter Brune s = dd->s; 33295c13181SPeter Brune bx = dd->bx; 33395c13181SPeter Brune by = dd->by; 33495c13181SPeter Brune bz = dd->bz; 3358865f1eaSKarl Rupp 33695c13181SPeter Brune stencil_type = dd->stencil_type; 3378865f1eaSKarl Rupp 338320964c4SBlaise Bourdin ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr); 339320964c4SBlaise Bourdin if (dim == 1) { 340ce94432eSBarry Smith ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),bx,M,nfields,s,dd->lx,nda);CHKERRQ(ierr); 341320964c4SBlaise Bourdin } else if (dim == 2) { 342ce94432eSBarry Smith ierr = DMDACreate2d(PetscObjectComm((PetscObject)da),bx,by,stencil_type,M,N,m,n,nfields,s,lx,ly,nda);CHKERRQ(ierr); 343320964c4SBlaise Bourdin } else if (dim == 3) { 344ce94432eSBarry Smith ierr = DMDACreate3d(PetscObjectComm((PetscObject)da),bx,by,bz,stencil_type,M,N,P,m,n,p,nfields,s,lx,ly,lz,nda);CHKERRQ(ierr); 345bc2bf880SBarry Smith } 3466636e97aSMatthew G Knepley if (da->coordinates) { 3476636e97aSMatthew G Knepley ierr = PetscObjectReference((PetscObject)da->coordinates);CHKERRQ(ierr); 3488865f1eaSKarl Rupp 3496636e97aSMatthew G Knepley (*nda)->coordinates = da->coordinates; 350bc2bf880SBarry Smith } 35195c13181SPeter Brune 35295c13181SPeter Brune /* allow for getting a reduced DA corresponding to a domain decomposition */ 35395c13181SPeter Brune ierr = DMDAGetOffset(da,&ox,&oy,&oz,&Mo,&No,&Po);CHKERRQ(ierr); 35495c13181SPeter Brune ierr = DMDASetOffset(*nda,ox,oy,oz,Mo,No,Po);CHKERRQ(ierr); 35595c13181SPeter Brune 35695c13181SPeter Brune /* allow for getting a reduced DA corresponding to a coarsened DA */ 35795c13181SPeter Brune ierr = DMGetCoarsenLevel(da,&cl);CHKERRQ(ierr); 35895c13181SPeter Brune ierr = DMGetRefineLevel(da,&rl);CHKERRQ(ierr); 3598865f1eaSKarl Rupp 36095c13181SPeter Brune (*nda)->levelup = rl; 36195c13181SPeter Brune (*nda)->leveldown = cl; 362bc2bf880SBarry Smith PetscFunctionReturn(0); 363bc2bf880SBarry Smith } 364bc2bf880SBarry Smith 365