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