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