1 2 /* 3 Code for manipulating distributed regular arrays in parallel. 4 */ 5 6 #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 7 #include <petscdmfield.h> 8 9 PetscErrorCode DMCreateCoordinateDM_DA(DM dm, DM *cdm) 10 { 11 PetscErrorCode ierr; 12 PetscFunctionBegin; 13 ierr = DMDACreateCompatibleDMDA(dm,dm->dim,cdm);CHKERRQ(ierr); 14 PetscFunctionReturn(0); 15 } 16 17 PetscErrorCode DMCreateCoordinateField_DA(DM dm, DMField *field) 18 { 19 PetscReal gmin[3], gmax[3]; 20 PetscScalar corners[24]; 21 PetscInt dim; 22 PetscInt i, j; 23 DM cdm; 24 PetscErrorCode ierr; 25 26 PetscFunctionBegin; 27 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 28 /* TODO: this is wrong if coordinates are not rectilinear */ 29 ierr = DMGetBoundingBox(dm,gmin,gmax);CHKERRQ(ierr); 30 for (i = 0; i < (1 << dim); i++) { 31 for (j = 0; j < dim; j++) { 32 corners[i*dim + j] = (i & (1 << j)) ? gmax[j] : gmin[j]; 33 } 34 } 35 ierr = DMClone(dm,&cdm);CHKERRQ(ierr); 36 ierr = DMFieldCreateDA(cdm,dim,corners,field);CHKERRQ(ierr); 37 ierr = DMDestroy(&cdm);CHKERRQ(ierr); 38 PetscFunctionReturn(0); 39 } 40 41 /*@C 42 DMDASetFieldName - Sets the names of individual field components in multicomponent 43 vectors associated with a DMDA. 44 45 Logically collective; name must contain a common value 46 47 Input Parameters: 48 + da - the distributed array 49 . nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the 50 number of degrees of freedom per node within the DMDA 51 - names - the name of the field (component) 52 53 Notes: 54 It must be called after having called DMSetUp(). 55 56 Level: intermediate 57 58 .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMDASetFieldNames(), DMSetUp() 59 @*/ 60 PetscErrorCode DMDASetFieldName(DM da,PetscInt nf,const char name[]) 61 { 62 PetscErrorCode ierr; 63 DM_DA *dd = (DM_DA*)da->data; 64 65 PetscFunctionBegin; 66 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 67 if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf); 68 if (!dd->fieldname) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 69 ierr = PetscFree(dd->fieldname[nf]);CHKERRQ(ierr); 70 ierr = PetscStrallocpy(name,&dd->fieldname[nf]);CHKERRQ(ierr); 71 PetscFunctionReturn(0); 72 } 73 74 /*@C 75 DMDAGetFieldNames - Gets the name of each component in the vector associated with the DMDA 76 77 Not collective; names will contain a common value 78 79 Input Parameter: 80 . dm - the DMDA object 81 82 Output Parameter: 83 . names - the names of the components, final string is NULL, will have the same number of entries as the dof used in creating the DMDA 84 85 Level: intermediate 86 87 Not supported from Fortran, use DMDAGetFieldName() 88 89 .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMDASetFieldName(), DMDASetFieldNames() 90 @*/ 91 PetscErrorCode DMDAGetFieldNames(DM da,const char * const **names) 92 { 93 DM_DA *dd = (DM_DA*)da->data; 94 95 PetscFunctionBegin; 96 *names = (const char * const *) dd->fieldname; 97 PetscFunctionReturn(0); 98 } 99 100 /*@C 101 DMDASetFieldNames - Sets the name of each component in the vector associated with the DMDA 102 103 Logically collective; names must contain a common value 104 105 Input Parameters: 106 + dm - the DMDA object 107 - names - the names of the components, final string must be NULL, must have the same number of entries as the dof used in creating the DMDA 108 109 Notes: 110 It must be called after having called DMSetUp(). 111 112 Level: intermediate 113 114 Not supported from Fortran, use DMDASetFieldName() 115 116 .seealso: DMDAGetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMDASetFieldName(), DMSetUp() 117 @*/ 118 PetscErrorCode DMDASetFieldNames(DM da,const char * const *names) 119 { 120 PetscErrorCode ierr; 121 DM_DA *dd = (DM_DA*)da->data; 122 char **fieldname; 123 PetscInt nf = 0; 124 125 PetscFunctionBegin; 126 if (!dd->fieldname) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 127 while (names[nf++]) {}; 128 if (nf != dd->w+1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid number of fields %D",nf-1); 129 ierr = PetscStrArrayallocpy(names,&fieldname);CHKERRQ(ierr); 130 ierr = PetscStrArrayDestroy(&dd->fieldname);CHKERRQ(ierr); 131 dd->fieldname = fieldname; 132 PetscFunctionReturn(0); 133 } 134 135 /*@C 136 DMDAGetFieldName - Gets the names of individual field components in multicomponent 137 vectors associated with a DMDA. 138 139 Not collective; name will contain a common value 140 141 Input Parameters: 142 + da - the distributed array 143 - nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the 144 number of degrees of freedom per node within the DMDA 145 146 Output Parameter: 147 . names - the name of the field (component) 148 149 Notes: 150 It must be called after having called DMSetUp(). 151 152 Level: intermediate 153 154 .seealso: DMDASetFieldName(), DMDASetCoordinateName(), DMDAGetCoordinateName(), DMSetUp() 155 @*/ 156 PetscErrorCode DMDAGetFieldName(DM da,PetscInt nf,const char **name) 157 { 158 DM_DA *dd = (DM_DA*)da->data; 159 160 PetscFunctionBegin; 161 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 162 PetscValidPointer(name,3); 163 if (nf < 0 || nf >= dd->w) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %D",nf); 164 if (!dd->fieldname) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 165 *name = dd->fieldname[nf]; 166 PetscFunctionReturn(0); 167 } 168 169 /*@C 170 DMDASetCoordinateName - Sets the name of the coordinate directions associated with a DMDA, for example "x" or "y" 171 172 Logically collective; name must contain a common value 173 174 Input Parameters: 175 + dm - the DM 176 . nf - coordinate number for the DMDA (0, 1, ... dim-1), 177 - name - the name of the coordinate 178 179 Notes: 180 It must be called after having called DMSetUp(). 181 182 Level: intermediate 183 184 Not supported from Fortran 185 186 .seealso: DMDAGetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName(), DMSetUp() 187 @*/ 188 PetscErrorCode DMDASetCoordinateName(DM dm,PetscInt nf,const char name[]) 189 { 190 PetscErrorCode ierr; 191 DM_DA *dd = (DM_DA*)dm->data; 192 193 PetscFunctionBegin; 194 PetscValidHeaderSpecificType(dm,DM_CLASSID,1,DMDA); 195 if (nf < 0 || nf >= dm->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf); 196 if (!dd->coordinatename) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 197 ierr = PetscFree(dd->coordinatename[nf]);CHKERRQ(ierr); 198 ierr = PetscStrallocpy(name,&dd->coordinatename[nf]);CHKERRQ(ierr); 199 PetscFunctionReturn(0); 200 } 201 202 /*@C 203 DMDAGetCoordinateName - Gets the name of a coodinate direction associated with a DMDA. 204 205 Not collective; name will contain a common value 206 207 Input Parameters: 208 + dm - the DM 209 - nf - number for the DMDA (0, 1, ... dim-1) 210 211 Output Parameter: 212 . names - the name of the coordinate direction 213 214 Notes: 215 It must be called after having called DMSetUp(). 216 217 Level: intermediate 218 219 Not supported from Fortran 220 221 .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName(), DMSetUp() 222 @*/ 223 PetscErrorCode DMDAGetCoordinateName(DM dm,PetscInt nf,const char **name) 224 { 225 DM_DA *dd = (DM_DA*)dm->data; 226 227 PetscFunctionBegin; 228 PetscValidHeaderSpecificType(dm,DM_CLASSID,1,DMDA); 229 PetscValidPointer(name,3); 230 if (nf < 0 || nf >= dm->dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %D",nf); 231 if (!dd->coordinatename) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 232 *name = dd->coordinatename[nf]; 233 PetscFunctionReturn(0); 234 } 235 236 /*@C 237 DMDAGetCorners - Returns the global (x,y,z) indices of the lower left 238 corner and size of the local region, excluding ghost points. 239 240 Not collective 241 242 Input Parameter: 243 . da - the distributed array 244 245 Output Parameters: 246 + x - the corner index for the first dimension 247 . y - the corner index for the second dimension (only used in 2D and 3D problems) 248 . z - the corner index for the third dimension (only used in 3D problems) 249 . m - the width in the first dimension 250 . n - the width in the second dimension (only used in 2D and 3D problems) 251 - p - the width in the third dimension (only used in 3D problems) 252 253 Note: 254 The corner information is independent of the number of degrees of 255 freedom per node set with the DMDACreateXX() routine. Thus the x, y, z, and 256 m, n, p can be thought of as coordinates on a logical grid, where each 257 grid point has (potentially) several degrees of freedom. 258 Any of y, z, n, and p can be passed in as NULL if not needed. 259 260 Level: beginner 261 262 .seealso: DMDAGetGhostCorners(), DMDAGetOwnershipRanges(), DMStagGetCorners() 263 @*/ 264 PetscErrorCode DMDAGetCorners(DM da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p) 265 { 266 PetscInt w; 267 DM_DA *dd = (DM_DA*)da->data; 268 269 PetscFunctionBegin; 270 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 271 /* since the xs, xe ... have all been multiplied by the number of degrees 272 of freedom per cell, w = dd->w, we divide that out before returning.*/ 273 w = dd->w; 274 if (x) *x = dd->xs/w + dd->xo; 275 /* the y and z have NOT been multiplied by w */ 276 if (y) *y = dd->ys + dd->yo; 277 if (z) *z = dd->zs + dd->zo; 278 if (m) *m = (dd->xe - dd->xs)/w; 279 if (n) *n = (dd->ye - dd->ys); 280 if (p) *p = (dd->ze - dd->zs); 281 PetscFunctionReturn(0); 282 } 283 284 PetscErrorCode DMGetLocalBoundingIndices_DMDA(DM dm, PetscReal lmin[], PetscReal lmax[]) 285 { 286 DMDALocalInfo info; 287 PetscErrorCode ierr; 288 289 PetscFunctionBegin; 290 ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr); 291 lmin[0] = info.xs; 292 lmin[1] = info.ys; 293 lmin[2] = info.zs; 294 lmax[0] = info.xs + info.xm-1; 295 lmax[1] = info.ys + info.ym-1; 296 lmax[2] = info.zs + info.zm-1; 297 PetscFunctionReturn(0); 298 } 299 300 /*@ 301 DMDAGetReducedDMDA - Deprecated; use DMDACreateCompatibleDMDA() 302 303 Level: deprecated 304 @*/ 305 PetscErrorCode DMDAGetReducedDMDA(DM da,PetscInt nfields,DM *nda) 306 { 307 PetscErrorCode ierr; 308 309 PetscFunctionBegin; 310 ierr = DMDACreateCompatibleDMDA(da,nfields,nda);CHKERRQ(ierr); 311 PetscFunctionReturn(0); 312 } 313 314 /*@ 315 DMDACreateCompatibleDMDA - Creates a DMDA with the same layout but with fewer or more fields 316 317 Collective 318 319 Input Parameters: 320 + da - the distributed array 321 - nfields - number of fields in new DMDA 322 323 Output Parameter: 324 . nda - the new DMDA 325 326 Level: intermediate 327 328 .seealso: DMDAGetGhostCorners(), DMSetCoordinates(), DMDASetUniformCoordinates(), DMGetCoordinates(), DMDAGetGhostedCoordinates(), DMStagCreateCompatibleDMStag() 329 @*/ 330 PetscErrorCode DMDACreateCompatibleDMDA(DM da,PetscInt nfields,DM *nda) 331 { 332 PetscErrorCode ierr; 333 DM_DA *dd = (DM_DA*)da->data; 334 PetscInt s,m,n,p,M,N,P,dim,Mo,No,Po; 335 const PetscInt *lx,*ly,*lz; 336 DMBoundaryType bx,by,bz; 337 DMDAStencilType stencil_type; 338 PetscInt ox,oy,oz; 339 PetscInt cl,rl; 340 341 PetscFunctionBegin; 342 dim = da->dim; 343 M = dd->M; 344 N = dd->N; 345 P = dd->P; 346 m = dd->m; 347 n = dd->n; 348 p = dd->p; 349 s = dd->s; 350 bx = dd->bx; 351 by = dd->by; 352 bz = dd->bz; 353 354 stencil_type = dd->stencil_type; 355 356 ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr); 357 if (dim == 1) { 358 ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),bx,M,nfields,s,dd->lx,nda);CHKERRQ(ierr); 359 } else if (dim == 2) { 360 ierr = DMDACreate2d(PetscObjectComm((PetscObject)da),bx,by,stencil_type,M,N,m,n,nfields,s,lx,ly,nda);CHKERRQ(ierr); 361 } else if (dim == 3) { 362 ierr = DMDACreate3d(PetscObjectComm((PetscObject)da),bx,by,bz,stencil_type,M,N,P,m,n,p,nfields,s,lx,ly,lz,nda);CHKERRQ(ierr); 363 } 364 ierr = DMSetUp(*nda);CHKERRQ(ierr); 365 if (da->coordinates) { 366 ierr = PetscObjectReference((PetscObject)da->coordinates);CHKERRQ(ierr); 367 (*nda)->coordinates = da->coordinates; 368 } 369 370 /* allow for getting a reduced DA corresponding to a domain decomposition */ 371 ierr = DMDAGetOffset(da,&ox,&oy,&oz,&Mo,&No,&Po);CHKERRQ(ierr); 372 ierr = DMDASetOffset(*nda,ox,oy,oz,Mo,No,Po);CHKERRQ(ierr); 373 374 /* allow for getting a reduced DA corresponding to a coarsened DA */ 375 ierr = DMGetCoarsenLevel(da,&cl);CHKERRQ(ierr); 376 ierr = DMGetRefineLevel(da,&rl);CHKERRQ(ierr); 377 378 (*nda)->levelup = rl; 379 (*nda)->leveldown = cl; 380 PetscFunctionReturn(0); 381 } 382 383 /*@C 384 DMDAGetCoordinateArray - Gets an array containing the coordinates of the DMDA 385 386 Not collective 387 388 Input Parameter: 389 . dm - the DM 390 391 Output Parameter: 392 . xc - the coordinates 393 394 Level: intermediate 395 396 Not supported from Fortran 397 398 .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName(), DMDARestoreCoordinateArray() 399 @*/ 400 PetscErrorCode DMDAGetCoordinateArray(DM dm,void *xc) 401 { 402 PetscErrorCode ierr; 403 DM cdm; 404 Vec x; 405 406 PetscFunctionBegin; 407 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 408 ierr = DMGetCoordinates(dm,&x);CHKERRQ(ierr); 409 ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr); 410 ierr = DMDAVecGetArray(cdm,x,xc);CHKERRQ(ierr); 411 PetscFunctionReturn(0); 412 } 413 414 /*@C 415 DMDARestoreCoordinateArray - Sets an array containing the coordinates of the DMDA 416 417 Not collective 418 419 Input Parameters: 420 + dm - the DM 421 - xc - the coordinates 422 423 Level: intermediate 424 425 Not supported from Fortran 426 427 .seealso: DMDASetCoordinateName(), DMDASetFieldName(), DMDAGetFieldName(), DMDAGetCoordinateArray() 428 @*/ 429 PetscErrorCode DMDARestoreCoordinateArray(DM dm,void *xc) 430 { 431 PetscErrorCode ierr; 432 DM cdm; 433 Vec x; 434 435 PetscFunctionBegin; 436 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 437 ierr = DMGetCoordinates(dm,&x);CHKERRQ(ierr); 438 ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr); 439 ierr = DMDAVecRestoreArray(cdm,x,xc);CHKERRQ(ierr); 440 PetscFunctionReturn(0); 441 } 442