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