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