1*47c6ae99SBarry Smith #define PETSCDM_DLL 2*47c6ae99SBarry Smith 3*47c6ae99SBarry Smith /* 4*47c6ae99SBarry Smith Code for manipulating distributed regular arrays in parallel. 5*47c6ae99SBarry Smith */ 6*47c6ae99SBarry Smith 7*47c6ae99SBarry Smith #include "private/daimpl.h" /*I "petscda.h" I*/ 8*47c6ae99SBarry Smith 9*47c6ae99SBarry Smith #undef __FUNCT__ 10*47c6ae99SBarry Smith #define __FUNCT__ "DASetCoordinates" 11*47c6ae99SBarry Smith /*@ 12*47c6ae99SBarry Smith DASetCoordinates - Sets into the DA a vector that indicates the 13*47c6ae99SBarry Smith coordinates of the local nodes (NOT including ghost nodes). 14*47c6ae99SBarry Smith 15*47c6ae99SBarry Smith Collective on DA 16*47c6ae99SBarry Smith 17*47c6ae99SBarry Smith Input Parameter: 18*47c6ae99SBarry Smith + da - the distributed array 19*47c6ae99SBarry Smith - c - coordinate vector 20*47c6ae99SBarry Smith 21*47c6ae99SBarry Smith Note: 22*47c6ae99SBarry Smith The coordinates should NOT include those for all ghost points 23*47c6ae99SBarry Smith 24*47c6ae99SBarry Smith Level: intermediate 25*47c6ae99SBarry Smith 26*47c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates 27*47c6ae99SBarry Smith 28*47c6ae99SBarry Smith .seealso: DAGetGhostCorners(), DAGetCoordinates(), DASetUniformCoordinates(). DAGetGhostedCoordinates(), DAGetCoordinateDA() 29*47c6ae99SBarry Smith @*/ 30*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetCoordinates(DA da,Vec c) 31*47c6ae99SBarry Smith { 32*47c6ae99SBarry Smith PetscErrorCode ierr; 33*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 34*47c6ae99SBarry Smith 35*47c6ae99SBarry Smith PetscFunctionBegin; 36*47c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 37*47c6ae99SBarry Smith PetscValidHeaderSpecific(c,VEC_CLASSID,2); 38*47c6ae99SBarry Smith ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 39*47c6ae99SBarry Smith if (dd->coordinates) {ierr = VecDestroy(dd->coordinates);CHKERRQ(ierr);} 40*47c6ae99SBarry Smith dd->coordinates = c; 41*47c6ae99SBarry Smith ierr = VecSetBlockSize(c,dd->dim);CHKERRQ(ierr); 42*47c6ae99SBarry Smith if (dd->ghosted_coordinates) { /* The ghosted coordinates are no longer valid */ 43*47c6ae99SBarry Smith ierr = VecDestroy(dd->ghosted_coordinates);CHKERRQ(ierr); 44*47c6ae99SBarry Smith dd->ghosted_coordinates = PETSC_NULL; 45*47c6ae99SBarry Smith } 46*47c6ae99SBarry Smith PetscFunctionReturn(0); 47*47c6ae99SBarry Smith } 48*47c6ae99SBarry Smith 49*47c6ae99SBarry Smith #undef __FUNCT__ 50*47c6ae99SBarry Smith #define __FUNCT__ "DAGetCoordinates" 51*47c6ae99SBarry Smith /*@ 52*47c6ae99SBarry Smith DAGetCoordinates - Gets the node coordinates associated with a DA. 53*47c6ae99SBarry Smith 54*47c6ae99SBarry Smith Not Collective 55*47c6ae99SBarry Smith 56*47c6ae99SBarry Smith Input Parameter: 57*47c6ae99SBarry Smith . da - the distributed array 58*47c6ae99SBarry Smith 59*47c6ae99SBarry Smith Output Parameter: 60*47c6ae99SBarry Smith . c - coordinate vector 61*47c6ae99SBarry Smith 62*47c6ae99SBarry Smith Note: 63*47c6ae99SBarry Smith Each process has only the coordinates for its local nodes (does NOT have the 64*47c6ae99SBarry Smith coordinates for the ghost nodes). 65*47c6ae99SBarry Smith 66*47c6ae99SBarry Smith For two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 67*47c6ae99SBarry Smith and (x_0,y_0,z_0,x_1,y_1,z_1...) 68*47c6ae99SBarry Smith 69*47c6ae99SBarry Smith Level: intermediate 70*47c6ae99SBarry Smith 71*47c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates 72*47c6ae99SBarry Smith 73*47c6ae99SBarry Smith .seealso: DAGetGhostCorners(), DASetCoordinates(), DASetUniformCoordinates(), DAGetGhostedCoordinates(), DAGetCoordinateDA() 74*47c6ae99SBarry Smith @*/ 75*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetCoordinates(DA da,Vec *c) 76*47c6ae99SBarry Smith { 77*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 78*47c6ae99SBarry Smith PetscFunctionBegin; 79*47c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 80*47c6ae99SBarry Smith PetscValidPointer(c,2); 81*47c6ae99SBarry Smith *c = dd->coordinates; 82*47c6ae99SBarry Smith PetscFunctionReturn(0); 83*47c6ae99SBarry Smith } 84*47c6ae99SBarry Smith 85*47c6ae99SBarry Smith #undef __FUNCT__ 86*47c6ae99SBarry Smith #define __FUNCT__ "DAGetCoordinateDA" 87*47c6ae99SBarry Smith /*@ 88*47c6ae99SBarry Smith DAGetCoordinateDA - Gets the DA that scatters between global and local DA coordinates 89*47c6ae99SBarry Smith 90*47c6ae99SBarry Smith Collective on DA 91*47c6ae99SBarry Smith 92*47c6ae99SBarry Smith Input Parameter: 93*47c6ae99SBarry Smith . da - the distributed array 94*47c6ae99SBarry Smith 95*47c6ae99SBarry Smith Output Parameter: 96*47c6ae99SBarry Smith . dac - coordinate DA 97*47c6ae99SBarry Smith 98*47c6ae99SBarry Smith Level: intermediate 99*47c6ae99SBarry Smith 100*47c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates 101*47c6ae99SBarry Smith 102*47c6ae99SBarry Smith .seealso: DAGetGhostCorners(), DASetCoordinates(), DASetUniformCoordinates(), DAGetCoordinates(), DAGetGhostedCoordinates() 103*47c6ae99SBarry Smith @*/ 104*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetCoordinateDA(DA da,DA *cda) 105*47c6ae99SBarry Smith { 106*47c6ae99SBarry Smith PetscMPIInt size; 107*47c6ae99SBarry Smith PetscErrorCode ierr; 108*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 109*47c6ae99SBarry Smith 110*47c6ae99SBarry Smith PetscFunctionBegin; 111*47c6ae99SBarry Smith if (!dd->da_coordinates) { 112*47c6ae99SBarry Smith ierr = MPI_Comm_size(((PetscObject)da)->comm,&size);CHKERRQ(ierr); 113*47c6ae99SBarry Smith if (dd->dim == 1) { 114*47c6ae99SBarry Smith PetscInt s,m,*lc,l; 115*47c6ae99SBarry Smith DAPeriodicType pt; 116*47c6ae99SBarry Smith ierr = DAGetInfo(da,0,&m,0,0,0,0,0,0,&s,&pt,0);CHKERRQ(ierr); 117*47c6ae99SBarry Smith ierr = DAGetCorners(da,0,0,0,&l,0,0);CHKERRQ(ierr); 118*47c6ae99SBarry Smith ierr = PetscMalloc(size*sizeof(PetscInt),&lc);CHKERRQ(ierr); 119*47c6ae99SBarry Smith ierr = MPI_Allgather(&l,1,MPIU_INT,lc,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr); 120*47c6ae99SBarry Smith ierr = DACreate1d(((PetscObject)da)->comm,pt,m,1,s,lc,&dd->da_coordinates);CHKERRQ(ierr); 121*47c6ae99SBarry Smith ierr = PetscFree(lc);CHKERRQ(ierr); 122*47c6ae99SBarry Smith } else if (dd->dim == 2) { 123*47c6ae99SBarry Smith PetscInt i,s,m,*lc,*ld,l,k,n,M,N; 124*47c6ae99SBarry Smith DAPeriodicType pt; 125*47c6ae99SBarry Smith ierr = DAGetInfo(da,0,&m,&n,0,&M,&N,0,0,&s,&pt,0);CHKERRQ(ierr); 126*47c6ae99SBarry Smith ierr = DAGetCorners(da,0,0,0,&l,&k,0);CHKERRQ(ierr); 127*47c6ae99SBarry Smith ierr = PetscMalloc2(size,PetscInt,&lc,size,PetscInt,&ld);CHKERRQ(ierr); 128*47c6ae99SBarry Smith /* only first M values in lc matter */ 129*47c6ae99SBarry Smith ierr = MPI_Allgather(&l,1,MPIU_INT,lc,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr); 130*47c6ae99SBarry Smith /* every Mth value in ld matters */ 131*47c6ae99SBarry Smith ierr = MPI_Allgather(&k,1,MPIU_INT,ld,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr); 132*47c6ae99SBarry Smith for ( i=0; i<N; i++) { 133*47c6ae99SBarry Smith ld[i] = ld[M*i]; 134*47c6ae99SBarry Smith } 135*47c6ae99SBarry Smith ierr = DACreate2d(((PetscObject)da)->comm,pt,DA_STENCIL_BOX,m,n,M,N,2,s,lc,ld,&dd->da_coordinates);CHKERRQ(ierr); 136*47c6ae99SBarry Smith ierr = PetscFree2(lc,ld);CHKERRQ(ierr); 137*47c6ae99SBarry Smith } else if (dd->dim == 3) { 138*47c6ae99SBarry Smith PetscInt i,s,m,*lc,*ld,*le,l,k,q,n,M,N,P,p; 139*47c6ae99SBarry Smith DAPeriodicType pt; 140*47c6ae99SBarry Smith ierr = DAGetInfo(da,0,&m,&n,&p,&M,&N,&P,0,&s,&pt,0);CHKERRQ(ierr); 141*47c6ae99SBarry Smith ierr = DAGetCorners(da,0,0,0,&l,&k,&q);CHKERRQ(ierr); 142*47c6ae99SBarry Smith ierr = PetscMalloc3(size,PetscInt,&lc,size,PetscInt,&ld,size,PetscInt,&le);CHKERRQ(ierr); 143*47c6ae99SBarry Smith /* only first M values in lc matter */ 144*47c6ae99SBarry Smith ierr = MPI_Allgather(&l,1,MPIU_INT,lc,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr); 145*47c6ae99SBarry Smith /* every Mth value in ld matters */ 146*47c6ae99SBarry Smith ierr = MPI_Allgather(&k,1,MPIU_INT,ld,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr); 147*47c6ae99SBarry Smith for ( i=0; i<N; i++) { 148*47c6ae99SBarry Smith ld[i] = ld[M*i]; 149*47c6ae99SBarry Smith } 150*47c6ae99SBarry Smith ierr = MPI_Allgather(&q,1,MPIU_INT,le,1,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr); 151*47c6ae99SBarry Smith for ( i=0; i<P; i++) { 152*47c6ae99SBarry Smith le[i] = le[M*N*i]; 153*47c6ae99SBarry Smith } 154*47c6ae99SBarry Smith ierr = DACreate3d(((PetscObject)da)->comm,pt,DA_STENCIL_BOX,m,n,p,M,N,P,3,s,lc,ld,le,&dd->da_coordinates);CHKERRQ(ierr); 155*47c6ae99SBarry Smith ierr = PetscFree3(lc,ld,le);CHKERRQ(ierr); 156*47c6ae99SBarry Smith } 157*47c6ae99SBarry Smith } 158*47c6ae99SBarry Smith *cda = dd->da_coordinates; 159*47c6ae99SBarry Smith PetscFunctionReturn(0); 160*47c6ae99SBarry Smith } 161*47c6ae99SBarry Smith 162*47c6ae99SBarry Smith 163*47c6ae99SBarry Smith #undef __FUNCT__ 164*47c6ae99SBarry Smith #define __FUNCT__ "DAGetGhostedCoordinates" 165*47c6ae99SBarry Smith /*@ 166*47c6ae99SBarry Smith DAGetGhostedCoordinates - Gets the node coordinates associated with a DA. 167*47c6ae99SBarry Smith 168*47c6ae99SBarry Smith Collective on DA 169*47c6ae99SBarry Smith 170*47c6ae99SBarry Smith Input Parameter: 171*47c6ae99SBarry Smith . da - the distributed array 172*47c6ae99SBarry Smith 173*47c6ae99SBarry Smith Output Parameter: 174*47c6ae99SBarry Smith . c - coordinate vector 175*47c6ae99SBarry Smith 176*47c6ae99SBarry Smith Note: 177*47c6ae99SBarry Smith Each process has only the coordinates for its local AND ghost nodes 178*47c6ae99SBarry Smith 179*47c6ae99SBarry Smith For two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 180*47c6ae99SBarry Smith and (x_0,y_0,z_0,x_1,y_1,z_1...) 181*47c6ae99SBarry Smith 182*47c6ae99SBarry Smith Level: intermediate 183*47c6ae99SBarry Smith 184*47c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates 185*47c6ae99SBarry Smith 186*47c6ae99SBarry Smith .seealso: DAGetGhostCorners(), DASetCoordinates(), DASetUniformCoordinates(), DAGetCoordinates(), DAGetCoordinateDA() 187*47c6ae99SBarry Smith @*/ 188*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetGhostedCoordinates(DA da,Vec *c) 189*47c6ae99SBarry Smith { 190*47c6ae99SBarry Smith PetscErrorCode ierr; 191*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 192*47c6ae99SBarry Smith 193*47c6ae99SBarry Smith PetscFunctionBegin; 194*47c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 195*47c6ae99SBarry Smith PetscValidPointer(c,2); 196*47c6ae99SBarry Smith if (!dd->coordinates) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"You must call DASetCoordinates() before this call"); 197*47c6ae99SBarry Smith if (!dd->ghosted_coordinates) { 198*47c6ae99SBarry Smith DA dac; 199*47c6ae99SBarry Smith ierr = DAGetCoordinateDA(da,&dac);CHKERRQ(ierr); 200*47c6ae99SBarry Smith ierr = DACreateLocalVector(dac,&dd->ghosted_coordinates);CHKERRQ(ierr); 201*47c6ae99SBarry Smith ierr = DAGlobalToLocalBegin(dac,dd->coordinates,INSERT_VALUES,dd->ghosted_coordinates);CHKERRQ(ierr); 202*47c6ae99SBarry Smith ierr = DAGlobalToLocalEnd(dac,dd->coordinates,INSERT_VALUES,dd->ghosted_coordinates);CHKERRQ(ierr); 203*47c6ae99SBarry Smith } 204*47c6ae99SBarry Smith *c = dd->ghosted_coordinates; 205*47c6ae99SBarry Smith PetscFunctionReturn(0); 206*47c6ae99SBarry Smith } 207*47c6ae99SBarry Smith 208*47c6ae99SBarry Smith #undef __FUNCT__ 209*47c6ae99SBarry Smith #define __FUNCT__ "DASetFieldName" 210*47c6ae99SBarry Smith /*@C 211*47c6ae99SBarry Smith DASetFieldName - Sets the names of individual field components in multicomponent 212*47c6ae99SBarry Smith vectors associated with a DA. 213*47c6ae99SBarry Smith 214*47c6ae99SBarry Smith Not Collective 215*47c6ae99SBarry Smith 216*47c6ae99SBarry Smith Input Parameters: 217*47c6ae99SBarry Smith + da - the distributed array 218*47c6ae99SBarry Smith . nf - field number for the DA (0, 1, ... dof-1), where dof indicates the 219*47c6ae99SBarry Smith number of degrees of freedom per node within the DA 220*47c6ae99SBarry Smith - names - the name of the field (component) 221*47c6ae99SBarry Smith 222*47c6ae99SBarry Smith Level: intermediate 223*47c6ae99SBarry Smith 224*47c6ae99SBarry Smith .keywords: distributed array, get, component name 225*47c6ae99SBarry Smith 226*47c6ae99SBarry Smith .seealso: DAGetFieldName() 227*47c6ae99SBarry Smith @*/ 228*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetFieldName(DA da,PetscInt nf,const char name[]) 229*47c6ae99SBarry Smith { 230*47c6ae99SBarry Smith PetscErrorCode ierr; 231*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 232*47c6ae99SBarry Smith 233*47c6ae99SBarry Smith PetscFunctionBegin; 234*47c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 235*47c6ae99SBarry Smith if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf); 236*47c6ae99SBarry Smith if (dd->fieldname[nf]) {ierr = PetscFree(dd->fieldname[nf]);CHKERRQ(ierr);} 237*47c6ae99SBarry Smith ierr = PetscStrallocpy(name,&dd->fieldname[nf]);CHKERRQ(ierr); 238*47c6ae99SBarry Smith PetscFunctionReturn(0); 239*47c6ae99SBarry Smith } 240*47c6ae99SBarry Smith 241*47c6ae99SBarry Smith #undef __FUNCT__ 242*47c6ae99SBarry Smith #define __FUNCT__ "DAGetFieldName" 243*47c6ae99SBarry Smith /*@C 244*47c6ae99SBarry Smith DAGetFieldName - Gets the names of individual field components in multicomponent 245*47c6ae99SBarry Smith vectors associated with a DA. 246*47c6ae99SBarry Smith 247*47c6ae99SBarry Smith Not Collective 248*47c6ae99SBarry Smith 249*47c6ae99SBarry Smith Input Parameter: 250*47c6ae99SBarry Smith + da - the distributed array 251*47c6ae99SBarry Smith - nf - field number for the DA (0, 1, ... dof-1), where dof indicates the 252*47c6ae99SBarry Smith number of degrees of freedom per node within the DA 253*47c6ae99SBarry Smith 254*47c6ae99SBarry Smith Output Parameter: 255*47c6ae99SBarry Smith . names - the name of the field (component) 256*47c6ae99SBarry Smith 257*47c6ae99SBarry Smith Level: intermediate 258*47c6ae99SBarry Smith 259*47c6ae99SBarry Smith .keywords: distributed array, get, component name 260*47c6ae99SBarry Smith 261*47c6ae99SBarry Smith .seealso: DASetFieldName() 262*47c6ae99SBarry Smith @*/ 263*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetFieldName(DA da,PetscInt nf,const char **name) 264*47c6ae99SBarry Smith { 265*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 266*47c6ae99SBarry Smith 267*47c6ae99SBarry Smith PetscFunctionBegin; 268*47c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 269*47c6ae99SBarry Smith PetscValidPointer(name,3); 270*47c6ae99SBarry Smith if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf); 271*47c6ae99SBarry Smith *name = dd->fieldname[nf]; 272*47c6ae99SBarry Smith PetscFunctionReturn(0); 273*47c6ae99SBarry Smith } 274*47c6ae99SBarry Smith 275*47c6ae99SBarry Smith #undef __FUNCT__ 276*47c6ae99SBarry Smith #define __FUNCT__ "DAGetCorners" 277*47c6ae99SBarry Smith /*@ 278*47c6ae99SBarry Smith DAGetCorners - Returns the global (x,y,z) indices of the lower left 279*47c6ae99SBarry Smith corner of the local region, excluding ghost points. 280*47c6ae99SBarry Smith 281*47c6ae99SBarry Smith Not Collective 282*47c6ae99SBarry Smith 283*47c6ae99SBarry Smith Input Parameter: 284*47c6ae99SBarry Smith . da - the distributed array 285*47c6ae99SBarry Smith 286*47c6ae99SBarry Smith Output Parameters: 287*47c6ae99SBarry Smith + x,y,z - the corner indices (where y and z are optional; these are used 288*47c6ae99SBarry Smith for 2D and 3D problems) 289*47c6ae99SBarry Smith - m,n,p - widths in the corresponding directions (where n and p are optional; 290*47c6ae99SBarry Smith these are used for 2D and 3D problems) 291*47c6ae99SBarry Smith 292*47c6ae99SBarry Smith Note: 293*47c6ae99SBarry Smith The corner information is independent of the number of degrees of 294*47c6ae99SBarry Smith freedom per node set with the DACreateXX() routine. Thus the x, y, z, and 295*47c6ae99SBarry Smith m, n, p can be thought of as coordinates on a logical grid, where each 296*47c6ae99SBarry Smith grid point has (potentially) several degrees of freedom. 297*47c6ae99SBarry Smith Any of y, z, n, and p can be passed in as PETSC_NULL if not needed. 298*47c6ae99SBarry Smith 299*47c6ae99SBarry Smith Level: beginner 300*47c6ae99SBarry Smith 301*47c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices 302*47c6ae99SBarry Smith 303*47c6ae99SBarry Smith .seealso: DAGetGhostCorners(), DAGetOwnershipRanges() 304*47c6ae99SBarry Smith @*/ 305*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetCorners(DA da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p) 306*47c6ae99SBarry Smith { 307*47c6ae99SBarry Smith PetscInt w; 308*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 309*47c6ae99SBarry Smith 310*47c6ae99SBarry Smith PetscFunctionBegin; 311*47c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 312*47c6ae99SBarry Smith /* since the xs, xe ... have all been multiplied by the number of degrees 313*47c6ae99SBarry Smith of freedom per cell, w = dd->w, we divide that out before returning.*/ 314*47c6ae99SBarry Smith w = dd->w; 315*47c6ae99SBarry Smith if (x) *x = dd->xs/w; if(m) *m = (dd->xe - dd->xs)/w; 316*47c6ae99SBarry Smith /* the y and z have NOT been multiplied by w */ 317*47c6ae99SBarry Smith if (y) *y = dd->ys; if (n) *n = (dd->ye - dd->ys); 318*47c6ae99SBarry Smith if (z) *z = dd->zs; if (p) *p = (dd->ze - dd->zs); 319*47c6ae99SBarry Smith PetscFunctionReturn(0); 320*47c6ae99SBarry Smith } 321*47c6ae99SBarry Smith 322*47c6ae99SBarry Smith #undef __FUNCT__ 323*47c6ae99SBarry Smith #define __FUNCT__ "DAGetLocalBoundingBox" 324*47c6ae99SBarry Smith /*@ 325*47c6ae99SBarry Smith DAGetLocalBoundingBox - Returns the local bounding box for the DA. 326*47c6ae99SBarry Smith 327*47c6ae99SBarry Smith Not Collective 328*47c6ae99SBarry Smith 329*47c6ae99SBarry Smith Input Parameter: 330*47c6ae99SBarry Smith . da - the distributed array 331*47c6ae99SBarry Smith 332*47c6ae99SBarry Smith Output Parameters: 333*47c6ae99SBarry Smith + lmin - local minimum coordinates (length dim, optional) 334*47c6ae99SBarry Smith - lmax - local maximim coordinates (length dim, optional) 335*47c6ae99SBarry Smith 336*47c6ae99SBarry Smith Level: beginner 337*47c6ae99SBarry Smith 338*47c6ae99SBarry Smith .keywords: distributed array, get, coordinates 339*47c6ae99SBarry Smith 340*47c6ae99SBarry Smith .seealso: DAGetCoordinateDA(), DAGetCoordinates(), DAGetBoundingBox() 341*47c6ae99SBarry Smith @*/ 342*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetLocalBoundingBox(DA da,PetscReal lmin[],PetscReal lmax[]) 343*47c6ae99SBarry Smith { 344*47c6ae99SBarry Smith PetscErrorCode ierr; 345*47c6ae99SBarry Smith Vec coords = PETSC_NULL; 346*47c6ae99SBarry Smith PetscInt dim,i,j; 347*47c6ae99SBarry Smith const PetscScalar *local_coords; 348*47c6ae99SBarry Smith PetscReal min[3]={PETSC_MAX,PETSC_MAX,PETSC_MAX},max[3]={PETSC_MIN,PETSC_MIN,PETSC_MIN}; 349*47c6ae99SBarry Smith PetscInt N,Ni; 350*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 351*47c6ae99SBarry Smith 352*47c6ae99SBarry Smith PetscFunctionBegin; 353*47c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 354*47c6ae99SBarry Smith dim = dd->dim; 355*47c6ae99SBarry Smith ierr = DAGetCoordinates(da,&coords);CHKERRQ(ierr); 356*47c6ae99SBarry Smith ierr = VecGetArrayRead(coords,&local_coords);CHKERRQ(ierr); 357*47c6ae99SBarry Smith ierr = VecGetLocalSize(coords,&N);CHKERRQ(ierr); 358*47c6ae99SBarry Smith Ni = N/dim; 359*47c6ae99SBarry Smith for (i=0; i<Ni; i++) { 360*47c6ae99SBarry Smith for (j=0; j<dim; j++) { 361*47c6ae99SBarry Smith min[j] = PetscMin(min[j],PetscRealPart(local_coords[i*dim+j]));CHKERRQ(ierr); 362*47c6ae99SBarry Smith max[j] = PetscMax(min[j],PetscRealPart(local_coords[i*dim+j]));CHKERRQ(ierr); 363*47c6ae99SBarry Smith } 364*47c6ae99SBarry Smith } 365*47c6ae99SBarry Smith ierr = VecRestoreArrayRead(coords,&local_coords);CHKERRQ(ierr); 366*47c6ae99SBarry Smith if (lmin) {ierr = PetscMemcpy(lmin,min,dim*sizeof(PetscReal));CHKERRQ(ierr);} 367*47c6ae99SBarry Smith if (lmax) {ierr = PetscMemcpy(lmax,max,dim*sizeof(PetscReal));CHKERRQ(ierr);} 368*47c6ae99SBarry Smith PetscFunctionReturn(0); 369*47c6ae99SBarry Smith } 370*47c6ae99SBarry Smith 371*47c6ae99SBarry Smith #undef __FUNCT__ 372*47c6ae99SBarry Smith #define __FUNCT__ "DAGetBoundingBox" 373*47c6ae99SBarry Smith /*@ 374*47c6ae99SBarry Smith DAGetBoundingBox - Returns the global bounding box for the DA. 375*47c6ae99SBarry Smith 376*47c6ae99SBarry Smith Collective on DA 377*47c6ae99SBarry Smith 378*47c6ae99SBarry Smith Input Parameter: 379*47c6ae99SBarry Smith . da - the distributed array 380*47c6ae99SBarry Smith 381*47c6ae99SBarry Smith Output Parameters: 382*47c6ae99SBarry Smith + gmin - global minimum coordinates (length dim, optional) 383*47c6ae99SBarry Smith - gmax - global maximim coordinates (length dim, optional) 384*47c6ae99SBarry Smith 385*47c6ae99SBarry Smith Level: beginner 386*47c6ae99SBarry Smith 387*47c6ae99SBarry Smith .keywords: distributed array, get, coordinates 388*47c6ae99SBarry Smith 389*47c6ae99SBarry Smith .seealso: DAGetCoordinateDA(), DAGetCoordinates(), DAGetLocalBoundingBox() 390*47c6ae99SBarry Smith @*/ 391*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetBoundingBox(DA da,PetscReal gmin[],PetscReal gmax[]) 392*47c6ae99SBarry Smith { 393*47c6ae99SBarry Smith PetscErrorCode ierr; 394*47c6ae99SBarry Smith PetscMPIInt count; 395*47c6ae99SBarry Smith PetscReal lmin[3],lmax[3]; 396*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 397*47c6ae99SBarry Smith 398*47c6ae99SBarry Smith PetscFunctionBegin; 399*47c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 400*47c6ae99SBarry Smith count = PetscMPIIntCast(dd->dim); 401*47c6ae99SBarry Smith ierr = DAGetLocalBoundingBox(da,lmin,lmax);CHKERRQ(ierr); 402*47c6ae99SBarry Smith if (gmin) {ierr = MPI_Allreduce(lmin,gmin,count,MPIU_REAL,MPI_MIN,((PetscObject)da)->comm);CHKERRQ(ierr);} 403*47c6ae99SBarry Smith if (gmax) {ierr = MPI_Allreduce(lmax,gmax,count,MPIU_REAL,MPI_MAX,((PetscObject)da)->comm);CHKERRQ(ierr);} 404*47c6ae99SBarry Smith PetscFunctionReturn(0); 405*47c6ae99SBarry Smith } 406