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