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 <petscbt.h> 8 #include <petscsf.h> 9 #include <petscds.h> 10 #include <petscfe.h> 11 12 /* 13 This allows the DMDA vectors to properly tell MATLAB their dimensions 14 */ 15 #if defined(PETSC_HAVE_MATLAB_ENGINE) 16 #include <engine.h> /* MATLAB include file */ 17 #include <mex.h> /* MATLAB include file */ 18 static PetscErrorCode VecMatlabEnginePut_DA2d(PetscObject obj, void *mengine) { 19 PetscInt n, m; 20 Vec vec = (Vec)obj; 21 PetscScalar *array; 22 mxArray *mat; 23 DM da; 24 25 PetscFunctionBegin; 26 PetscCall(VecGetDM(vec, &da)); 27 PetscCheck(da, PetscObjectComm((PetscObject)vec), PETSC_ERR_ARG_WRONGSTATE, "Vector not associated with a DMDA"); 28 PetscCall(DMDAGetGhostCorners(da, 0, 0, 0, &m, &n, 0)); 29 30 PetscCall(VecGetArray(vec, &array)); 31 #if !defined(PETSC_USE_COMPLEX) 32 mat = mxCreateDoubleMatrix(m, n, mxREAL); 33 #else 34 mat = mxCreateDoubleMatrix(m, n, mxCOMPLEX); 35 #endif 36 PetscCall(PetscArraycpy(mxGetPr(mat), array, n * m)); 37 PetscCall(PetscObjectName(obj)); 38 engPutVariable((Engine *)mengine, obj->name, mat); 39 40 PetscCall(VecRestoreArray(vec, &array)); 41 PetscFunctionReturn(0); 42 } 43 #endif 44 45 PetscErrorCode DMCreateLocalVector_DA(DM da, Vec *g) { 46 DM_DA *dd = (DM_DA *)da->data; 47 48 PetscFunctionBegin; 49 PetscValidHeaderSpecific(da, DM_CLASSID, 1); 50 PetscValidPointer(g, 2); 51 PetscCall(VecCreate(PETSC_COMM_SELF, g)); 52 PetscCall(VecSetSizes(*g, dd->nlocal, PETSC_DETERMINE)); 53 PetscCall(VecSetBlockSize(*g, dd->w)); 54 PetscCall(VecSetType(*g, da->vectype)); 55 if (dd->nlocal < da->bind_below) { 56 PetscCall(VecSetBindingPropagates(*g, PETSC_TRUE)); 57 PetscCall(VecBindToCPU(*g, PETSC_TRUE)); 58 } 59 PetscCall(VecSetDM(*g, da)); 60 #if defined(PETSC_HAVE_MATLAB_ENGINE) 61 if (dd->w == 1 && da->dim == 2) { PetscCall(PetscObjectComposeFunction((PetscObject)*g, "PetscMatlabEnginePut_C", VecMatlabEnginePut_DA2d)); } 62 #endif 63 PetscFunctionReturn(0); 64 } 65 66 /*@ 67 DMDAGetNumCells - Get the number of cells in the local piece of the DMDA. This includes ghost cells. 68 69 Input Parameter: 70 . dm - The DM object 71 72 Output Parameters: 73 + numCellsX - The number of local cells in the x-direction 74 . numCellsY - The number of local cells in the y-direction 75 . numCellsZ - The number of local cells in the z-direction 76 - numCells - The number of local cells 77 78 Level: developer 79 80 .seealso: `DMDAGetCellPoint()` 81 @*/ 82 PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells) { 83 DM_DA *da = (DM_DA *)dm->data; 84 const PetscInt dim = dm->dim; 85 const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 86 const PetscInt nC = (mx) * (dim > 1 ? (my) * (dim > 2 ? (mz) : 1) : 1); 87 88 PetscFunctionBegin; 89 PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 90 if (numCellsX) { 91 PetscValidIntPointer(numCellsX, 2); 92 *numCellsX = mx; 93 } 94 if (numCellsY) { 95 PetscValidIntPointer(numCellsY, 3); 96 *numCellsY = my; 97 } 98 if (numCellsZ) { 99 PetscValidIntPointer(numCellsZ, 4); 100 *numCellsZ = mz; 101 } 102 if (numCells) { 103 PetscValidIntPointer(numCells, 5); 104 *numCells = nC; 105 } 106 PetscFunctionReturn(0); 107 } 108 109 /*@ 110 DMDAGetCellPoint - Get the DM point corresponding to the tuple (i, j, k) in the DMDA 111 112 Input Parameters: 113 + dm - The DM object 114 - i,j,k - The global indices for the cell 115 116 Output Parameters: 117 . point - The local DM point 118 119 Level: developer 120 121 .seealso: `DMDAGetNumCells()` 122 @*/ 123 PetscErrorCode DMDAGetCellPoint(DM dm, PetscInt i, PetscInt j, PetscInt k, PetscInt *point) { 124 const PetscInt dim = dm->dim; 125 DMDALocalInfo info; 126 127 PetscFunctionBegin; 128 PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 129 PetscValidIntPointer(point, 5); 130 PetscCall(DMDAGetLocalInfo(dm, &info)); 131 if (dim > 0) PetscCheck(!(i < info.gxs) && !(i >= info.gxs + info.gxm), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "X index %" PetscInt_FMT " not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", i, info.gxs, info.gxs + info.gxm); 132 if (dim > 1) PetscCheck(!(j < info.gys) && !(j >= info.gys + info.gym), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Y index %" PetscInt_FMT " not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", j, info.gys, info.gys + info.gym); 133 if (dim > 2) PetscCheck(!(k < info.gzs) && !(k >= info.gzs + info.gzm), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Z index %" PetscInt_FMT PetscInt_FMT " not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", k, info.gzs, info.gzs + info.gzm); 134 *point = i + (dim > 1 ? (j + (dim > 2 ? k * info.gym : 0)) * info.gxm : 0); 135 PetscFunctionReturn(0); 136 } 137 138 PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices) { 139 DM_DA *da = (DM_DA *)dm->data; 140 const PetscInt dim = dm->dim; 141 const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 142 const PetscInt nVx = mx + 1; 143 const PetscInt nVy = dim > 1 ? (my + 1) : 1; 144 const PetscInt nVz = dim > 2 ? (mz + 1) : 1; 145 const PetscInt nV = nVx * nVy * nVz; 146 147 PetscFunctionBegin; 148 if (numVerticesX) { 149 PetscValidIntPointer(numVerticesX, 2); 150 *numVerticesX = nVx; 151 } 152 if (numVerticesY) { 153 PetscValidIntPointer(numVerticesY, 3); 154 *numVerticesY = nVy; 155 } 156 if (numVerticesZ) { 157 PetscValidIntPointer(numVerticesZ, 4); 158 *numVerticesZ = nVz; 159 } 160 if (numVertices) { 161 PetscValidIntPointer(numVertices, 5); 162 *numVertices = nV; 163 } 164 PetscFunctionReturn(0); 165 } 166 167 PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces) { 168 DM_DA *da = (DM_DA *)dm->data; 169 const PetscInt dim = dm->dim; 170 const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 171 const PetscInt nxF = (dim > 1 ? (my) * (dim > 2 ? (mz) : 1) : 1); 172 const PetscInt nXF = (mx + 1) * nxF; 173 const PetscInt nyF = mx * (dim > 2 ? mz : 1); 174 const PetscInt nYF = dim > 1 ? (my + 1) * nyF : 0; 175 const PetscInt nzF = mx * (dim > 1 ? my : 0); 176 const PetscInt nZF = dim > 2 ? (mz + 1) * nzF : 0; 177 178 PetscFunctionBegin; 179 if (numXFacesX) { 180 PetscValidIntPointer(numXFacesX, 2); 181 *numXFacesX = nxF; 182 } 183 if (numXFaces) { 184 PetscValidIntPointer(numXFaces, 3); 185 *numXFaces = nXF; 186 } 187 if (numYFacesY) { 188 PetscValidIntPointer(numYFacesY, 4); 189 *numYFacesY = nyF; 190 } 191 if (numYFaces) { 192 PetscValidIntPointer(numYFaces, 5); 193 *numYFaces = nYF; 194 } 195 if (numZFacesZ) { 196 PetscValidIntPointer(numZFacesZ, 6); 197 *numZFacesZ = nzF; 198 } 199 if (numZFaces) { 200 PetscValidIntPointer(numZFaces, 7); 201 *numZFaces = nZF; 202 } 203 PetscFunctionReturn(0); 204 } 205 206 PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd) { 207 const PetscInt dim = dm->dim; 208 PetscInt nC, nV, nXF, nYF, nZF; 209 210 PetscFunctionBegin; 211 if (pStart) PetscValidIntPointer(pStart, 3); 212 if (pEnd) PetscValidIntPointer(pEnd, 4); 213 PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC)); 214 PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV)); 215 PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF)); 216 if (height == 0) { 217 /* Cells */ 218 if (pStart) *pStart = 0; 219 if (pEnd) *pEnd = nC; 220 } else if (height == 1) { 221 /* Faces */ 222 if (pStart) *pStart = nC + nV; 223 if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF; 224 } else if (height == dim) { 225 /* Vertices */ 226 if (pStart) *pStart = nC; 227 if (pEnd) *pEnd = nC + nV; 228 } else if (height < 0) { 229 /* All points */ 230 if (pStart) *pStart = 0; 231 if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF; 232 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %" PetscInt_FMT " in the DA", height); 233 PetscFunctionReturn(0); 234 } 235 236 PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd) { 237 const PetscInt dim = dm->dim; 238 PetscInt nC, nV, nXF, nYF, nZF; 239 240 PetscFunctionBegin; 241 if (pStart) PetscValidIntPointer(pStart, 3); 242 if (pEnd) PetscValidIntPointer(pEnd, 4); 243 PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC)); 244 PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV)); 245 PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF)); 246 if (depth == dim) { 247 /* Cells */ 248 if (pStart) *pStart = 0; 249 if (pEnd) *pEnd = nC; 250 } else if (depth == dim - 1) { 251 /* Faces */ 252 if (pStart) *pStart = nC + nV; 253 if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF; 254 } else if (depth == 0) { 255 /* Vertices */ 256 if (pStart) *pStart = nC; 257 if (pEnd) *pEnd = nC + nV; 258 } else if (depth < 0) { 259 /* All points */ 260 if (pStart) *pStart = 0; 261 if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF; 262 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %" PetscInt_FMT " in the DA", depth); 263 PetscFunctionReturn(0); 264 } 265 266 PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize) { 267 const PetscInt dim = dm->dim; 268 PetscInt nC, nV, nXF, nYF, nZF; 269 270 PetscFunctionBegin; 271 *coneSize = 0; 272 PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC)); 273 PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV)); 274 PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF)); 275 switch (dim) { 276 case 2: 277 if (p >= 0) { 278 if (p < nC) { 279 *coneSize = 4; 280 } else if (p < nC + nV) { 281 *coneSize = 0; 282 } else if (p < nC + nV + nXF + nYF + nZF) { 283 *coneSize = 2; 284 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", p, nC + nV + nXF + nYF + nZF); 285 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %" PetscInt_FMT " is invalid", p); 286 break; 287 case 3: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 288 } 289 PetscFunctionReturn(0); 290 } 291 292 PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[]) { 293 const PetscInt dim = dm->dim; 294 PetscInt nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF; 295 296 PetscFunctionBegin; 297 if (!cone) PetscCall(DMGetWorkArray(dm, 6, MPIU_INT, cone)); 298 PetscCall(DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC)); 299 PetscCall(DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV)); 300 PetscCall(DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF)); 301 switch (dim) { 302 case 2: 303 if (p >= 0) { 304 if (p < nC) { 305 const PetscInt cy = p / nCx; 306 const PetscInt cx = p % nCx; 307 308 (*cone)[0] = cy * nxF + cx + nC + nV; 309 (*cone)[1] = cx * nyF + cy + nyF + nC + nV + nXF; 310 (*cone)[2] = cy * nxF + cx + nxF + nC + nV; 311 (*cone)[3] = cx * nyF + cy + nC + nV + nXF; 312 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones"); 313 } else if (p < nC + nV) { 314 } else if (p < nC + nV + nXF) { 315 const PetscInt fy = (p - nC + nV) / nxF; 316 const PetscInt fx = (p - nC + nV) % nxF; 317 318 (*cone)[0] = fy * nVx + fx + nC; 319 (*cone)[1] = fy * nVx + fx + 1 + nC; 320 } else if (p < nC + nV + nXF + nYF) { 321 const PetscInt fx = (p - nC + nV + nXF) / nyF; 322 const PetscInt fy = (p - nC + nV + nXF) % nyF; 323 324 (*cone)[0] = fy * nVx + fx + nC; 325 (*cone)[1] = fy * nVx + fx + nVx + nC; 326 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", p, nC + nV + nXF + nYF + nZF); 327 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %" PetscInt_FMT " is invalid", p); 328 break; 329 case 3: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 330 } 331 PetscFunctionReturn(0); 332 } 333 334 PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[]) { 335 PetscFunctionBegin; 336 PetscCall(DMGetWorkArray(dm, 6, MPIU_INT, cone)); 337 PetscFunctionReturn(0); 338 } 339 340 PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu) { 341 DM_DA *da = (DM_DA *)dm->data; 342 Vec coordinates; 343 PetscSection section; 344 PetscScalar *coords; 345 PetscReal h[3]; 346 PetscInt dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k; 347 348 PetscFunctionBegin; 349 PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 350 PetscCall(DMDAGetInfo(dm, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 351 PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "The following code only works for dim <= 3"); 352 h[0] = (xu - xl) / M; 353 h[1] = (yu - yl) / N; 354 h[2] = (zu - zl) / P; 355 PetscCall(DMDAGetDepthStratum(dm, 0, &vStart, &vEnd)); 356 PetscCall(DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV)); 357 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion)); 358 PetscCall(PetscSectionSetNumFields(section, 1)); 359 PetscCall(PetscSectionSetFieldComponents(section, 0, dim)); 360 PetscCall(PetscSectionSetChart(section, vStart, vEnd)); 361 for (v = vStart; v < vEnd; ++v) { PetscCall(PetscSectionSetDof(section, v, dim)); } 362 PetscCall(PetscSectionSetUp(section)); 363 PetscCall(PetscSectionGetStorageSize(section, &size)); 364 PetscCall(VecCreateSeq(PETSC_COMM_SELF, size, &coordinates)); 365 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 366 PetscCall(VecGetArray(coordinates, &coords)); 367 for (k = 0; k < nVz; ++k) { 368 PetscInt ind[3], d, off; 369 370 ind[0] = 0; 371 ind[1] = 0; 372 ind[2] = k + da->zs; 373 for (j = 0; j < nVy; ++j) { 374 ind[1] = j + da->ys; 375 for (i = 0; i < nVx; ++i) { 376 const PetscInt vertex = (k * nVy + j) * nVx + i + vStart; 377 378 PetscCall(PetscSectionGetOffset(section, vertex, &off)); 379 ind[0] = i + da->xs; 380 for (d = 0; d < dim; ++d) { coords[off + d] = h[d] * ind[d]; } 381 } 382 } 383 } 384 PetscCall(VecRestoreArray(coordinates, &coords)); 385 PetscCall(DMSetCoordinateSection(dm, PETSC_DETERMINE, section)); 386 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 387 PetscCall(PetscSectionDestroy(§ion)); 388 PetscCall(VecDestroy(&coordinates)); 389 PetscFunctionReturn(0); 390 } 391 392 /* ------------------------------------------------------------------- */ 393 394 /*@C 395 DMDAGetArray - Gets a work array for a DMDA 396 397 Input Parameters: 398 + da - information about my local patch 399 - ghosted - do you want arrays for the ghosted or nonghosted patch 400 401 Output Parameters: 402 . vptr - array data structured 403 404 Note: The vector values are NOT initialized and may have garbage in them, so you may need 405 to zero them. 406 407 Level: advanced 408 409 .seealso: `DMDARestoreArray()` 410 411 @*/ 412 PetscErrorCode DMDAGetArray(DM da, PetscBool ghosted, void *vptr) { 413 PetscInt j, i, xs, ys, xm, ym, zs, zm; 414 char *iarray_start; 415 void **iptr = (void **)vptr; 416 DM_DA *dd = (DM_DA *)da->data; 417 418 PetscFunctionBegin; 419 PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 420 if (ghosted) { 421 for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 422 if (dd->arrayghostedin[i]) { 423 *iptr = dd->arrayghostedin[i]; 424 iarray_start = (char *)dd->startghostedin[i]; 425 dd->arrayghostedin[i] = NULL; 426 dd->startghostedin[i] = NULL; 427 428 goto done; 429 } 430 } 431 xs = dd->Xs; 432 ys = dd->Ys; 433 zs = dd->Zs; 434 xm = dd->Xe - dd->Xs; 435 ym = dd->Ye - dd->Ys; 436 zm = dd->Ze - dd->Zs; 437 } else { 438 for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 439 if (dd->arrayin[i]) { 440 *iptr = dd->arrayin[i]; 441 iarray_start = (char *)dd->startin[i]; 442 dd->arrayin[i] = NULL; 443 dd->startin[i] = NULL; 444 445 goto done; 446 } 447 } 448 xs = dd->xs; 449 ys = dd->ys; 450 zs = dd->zs; 451 xm = dd->xe - dd->xs; 452 ym = dd->ye - dd->ys; 453 zm = dd->ze - dd->zs; 454 } 455 456 switch (da->dim) { 457 case 1: { 458 void *ptr; 459 460 PetscCall(PetscMalloc(xm * sizeof(PetscScalar), &iarray_start)); 461 462 ptr = (void *)(iarray_start - xs * sizeof(PetscScalar)); 463 *iptr = (void *)ptr; 464 break; 465 } 466 case 2: { 467 void **ptr; 468 469 PetscCall(PetscMalloc((ym + 1) * sizeof(void *) + xm * ym * sizeof(PetscScalar), &iarray_start)); 470 471 ptr = (void **)(iarray_start + xm * ym * sizeof(PetscScalar) - ys * sizeof(void *)); 472 for (j = ys; j < ys + ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar) * (xm * (j - ys) - xs); 473 *iptr = (void *)ptr; 474 break; 475 } 476 case 3: { 477 void ***ptr, **bptr; 478 479 PetscCall(PetscMalloc((zm + 1) * sizeof(void **) + (ym * zm + 1) * sizeof(void *) + xm * ym * zm * sizeof(PetscScalar), &iarray_start)); 480 481 ptr = (void ***)(iarray_start + xm * ym * zm * sizeof(PetscScalar) - zs * sizeof(void *)); 482 bptr = (void **)(iarray_start + xm * ym * zm * sizeof(PetscScalar) + zm * sizeof(void **)); 483 for (i = zs; i < zs + zm; i++) ptr[i] = bptr + ((i - zs) * ym - ys); 484 for (i = zs; i < zs + zm; i++) { 485 for (j = ys; j < ys + ym; j++) { ptr[i][j] = iarray_start + sizeof(PetscScalar) * (xm * ym * (i - zs) + xm * (j - ys) - xs); } 486 } 487 *iptr = (void *)ptr; 488 break; 489 } 490 default: SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", da->dim); 491 } 492 493 done: 494 /* add arrays to the checked out list */ 495 if (ghosted) { 496 for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 497 if (!dd->arrayghostedout[i]) { 498 dd->arrayghostedout[i] = *iptr; 499 dd->startghostedout[i] = iarray_start; 500 break; 501 } 502 } 503 } else { 504 for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 505 if (!dd->arrayout[i]) { 506 dd->arrayout[i] = *iptr; 507 dd->startout[i] = iarray_start; 508 break; 509 } 510 } 511 } 512 PetscFunctionReturn(0); 513 } 514 515 /*@C 516 DMDARestoreArray - Restores an array of derivative types for a DMDA 517 518 Input Parameters: 519 + da - information about my local patch 520 . ghosted - do you want arrays for the ghosted or nonghosted patch 521 - vptr - array data structured 522 523 Level: advanced 524 525 .seealso: `DMDAGetArray()` 526 527 @*/ 528 PetscErrorCode DMDARestoreArray(DM da, PetscBool ghosted, void *vptr) { 529 PetscInt i; 530 void **iptr = (void **)vptr, *iarray_start = NULL; 531 DM_DA *dd = (DM_DA *)da->data; 532 533 PetscFunctionBegin; 534 PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 535 if (ghosted) { 536 for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 537 if (dd->arrayghostedout[i] == *iptr) { 538 iarray_start = dd->startghostedout[i]; 539 dd->arrayghostedout[i] = NULL; 540 dd->startghostedout[i] = NULL; 541 break; 542 } 543 } 544 for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 545 if (!dd->arrayghostedin[i]) { 546 dd->arrayghostedin[i] = *iptr; 547 dd->startghostedin[i] = iarray_start; 548 break; 549 } 550 } 551 } else { 552 for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 553 if (dd->arrayout[i] == *iptr) { 554 iarray_start = dd->startout[i]; 555 dd->arrayout[i] = NULL; 556 dd->startout[i] = NULL; 557 break; 558 } 559 } 560 for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 561 if (!dd->arrayin[i]) { 562 dd->arrayin[i] = *iptr; 563 dd->startin[i] = iarray_start; 564 break; 565 } 566 } 567 } 568 PetscFunctionReturn(0); 569 } 570