1 #include <petsc-private/dmdaimpl.h> /*I "petscdmda.h" I*/ 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "FillClosureArray_Static" 5 PETSC_STATIC_INLINE PetscErrorCode FillClosureArray_Static(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, PetscInt *csize, const PetscScalar **array) 6 { 7 PetscScalar *a; 8 PetscInt pStart, pEnd, size = 0, dof, off, d, k, i; 9 PetscErrorCode ierr; 10 11 PetscFunctionBegin; 12 ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 13 for (i = 0; i < nP; ++i) { 14 const PetscInt p = points[i]; 15 16 if ((p < pStart) || (p >= pEnd)) continue; 17 ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 18 size += dof; 19 } 20 if (csize) *csize = size; 21 if (array) { 22 ierr = DMGetWorkArray(dm, size, PETSC_SCALAR, &a);CHKERRQ(ierr); 23 for (i = 0, k = 0; i < nP; ++i) { 24 const PetscInt p = points[i]; 25 26 if ((p < pStart) || (p >= pEnd)) continue; 27 ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 28 ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 29 30 for (d = 0; d < dof; ++d, ++k) a[k] = vArray[off+d]; 31 } 32 *array = a; 33 } 34 PetscFunctionReturn(0); 35 } 36 37 #undef __FUNCT__ 38 #define __FUNCT__ "FillClosureVec_Private" 39 PETSC_STATIC_INLINE PetscErrorCode FillClosureVec_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar *array, InsertMode mode) 40 { 41 PetscInt dof, off, d, k, i; 42 PetscErrorCode ierr; 43 44 PetscFunctionBegin; 45 if ((mode == INSERT_VALUES) || (mode == INSERT_ALL_VALUES)) { 46 for (i = 0, k = 0; i < nP; ++i) { 47 ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr); 48 ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr); 49 50 for (d = 0; d < dof; ++d, ++k) vArray[off+d] = array[k]; 51 } 52 } else { 53 for (i = 0, k = 0; i < nP; ++i) { 54 ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr); 55 ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr); 56 57 for (d = 0; d < dof; ++d, ++k) vArray[off+d] += array[k]; 58 } 59 } 60 PetscFunctionReturn(0); 61 } 62 63 #undef __FUNCT__ 64 #define __FUNCT__ "GetPointArray_Private" 65 PETSC_STATIC_INLINE PetscErrorCode GetPointArray_Private(DM dm,PetscInt n,PetscInt *points,PetscInt *rn,const PetscInt **rpoints) 66 { 67 PetscErrorCode ierr; 68 PetscInt *work; 69 70 PetscFunctionBegin; 71 if (rn) *rn = n; 72 if (rpoints) { 73 ierr = DMGetWorkArray(dm,n,PETSC_INT,&work);CHKERRQ(ierr); 74 ierr = PetscMemcpy(work,points,n*sizeof(PetscInt));CHKERRQ(ierr); 75 *rpoints = work; 76 } 77 PetscFunctionReturn(0); 78 } 79 80 #undef __FUNCT__ 81 #define __FUNCT__ "RestorePointArray_Private" 82 PETSC_STATIC_INLINE PetscErrorCode RestorePointArray_Private(DM dm,PetscInt *rn,const PetscInt **rpoints) 83 { 84 PetscErrorCode ierr; 85 86 PetscFunctionBegin; 87 if (rn) *rn = 0; 88 if (rpoints) { 89 ierr = DMRestoreWorkArray(dm,*rn, PETSC_INT, (void*) rpoints);CHKERRQ(ierr); 90 } 91 PetscFunctionReturn(0); 92 } 93 94 #undef __FUNCT__ 95 #define __FUNCT__ "DMDAGetTransitiveClosure" 96 PetscErrorCode DMDAGetTransitiveClosure(DM dm, PetscInt p, PetscBool useClosure, PetscInt *closureSize, PetscInt **closure) 97 { 98 DM_DA *da = (DM_DA *) dm->data; 99 PetscInt dim = da->dim; 100 PetscInt nVx, nVy, nVz, nxF, nXF, nyF, nYF, nzF, nZF; 101 PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart; 102 PetscErrorCode ierr; 103 104 PetscFunctionBegin; 105 if (!useClosure) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Star operation is not yet supported"); 106 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 107 ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 108 ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 109 ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 110 ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 111 ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, NULL);CHKERRQ(ierr); 112 ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 113 xfStart = fStart; xfEnd = xfStart+nXF; 114 yfStart = xfEnd; 115 if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd); 116 if ((p >= cStart) || (p < cEnd)) { 117 /* Cell */ 118 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 119 else if (dim == 2) { 120 /* 4 faces, 4 vertices 121 Bottom-left vertex follows same order as cells 122 Bottom y-face same order as cells 123 Left x-face follows same order as cells 124 We number the quad: 125 126 8--3--7 127 | | 128 4 0 2 129 | | 130 5--1--6 131 */ 132 PetscInt c = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1); 133 PetscInt v = cy*nVx + cx + vStart; 134 PetscInt xf = cy*nxF + cx + xfStart; 135 PetscInt yf = c + yfStart; 136 137 if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 9;} 138 if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);} 139 (*closure)[0] = p; (*closure)[1] = yf; (*closure)[2] = xf+1; (*closure)[3] = yf+nyF; (*closure)[4] = xf+0; (*closure)[5] = v+0; (*closure)[6]= v+1; (*closure)[7] = v+nVx+1; (*closure)[8] = v+nVx+0; 140 } else { 141 /* 6 faces, 12 edges, 8 vertices 142 Bottom-left vertex follows same order as cells 143 Bottom y-face same order as cells 144 Left x-face follows same order as cells 145 We number the quad: 146 147 8--3--7 148 | | 149 4 0 2 150 | | 151 5--1--6 152 */ 153 PetscInt c = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1), cz = c / ((nVx-1)*(nVy-1)); 154 PetscInt v = cy*nVx + cx + vStart; 155 PetscInt xf = cy*nxF + cx + xfStart; 156 PetscInt yf = c + yfStart; 157 158 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 159 if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 26;} 160 if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);} 161 } 162 } else if ((p >= vStart) || (p < vEnd)) { 163 /* Vertex */ 164 if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 1;} 165 if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);} 166 (*closure)[0] = p; 167 } else if ((p >= fStart) || (p < fStart + nXF)) { 168 /* X Face */ 169 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 170 else if (dim == 2) { 171 /* 2 vertices: The bottom vertex has the same numbering as the face */ 172 PetscInt f = p - xfStart; 173 174 if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 3;} 175 if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);} 176 (*closure)[0] = p; (*closure)[1] = f; (*closure)[2] = f+nVx; 177 } else if (dim == 3) { 178 /* 4 vertices */ 179 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 180 } 181 } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) { 182 /* Y Face */ 183 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 184 else if (dim == 2) { 185 /* 2 vertices: The left vertex has the same numbering as the face */ 186 PetscInt f = p - yfStart; 187 188 if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 3;} 189 if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);} 190 (*closure)[0] = p; (*closure)[1] = f; (*closure)[2]= f+1; 191 } else if (dim == 3) { 192 /* 4 vertices */ 193 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 194 } 195 } else { 196 /* Z Face */ 197 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 198 else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D"); 199 else if (dim == 3) { 200 /* 4 vertices */ 201 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 202 } 203 } 204 PetscFunctionReturn(0); 205 } 206 207 #undef __FUNCT__ 208 #define __FUNCT__ "DMDARestoreTransitiveClosure" 209 PetscErrorCode DMDARestoreTransitiveClosure(DM dm, PetscInt p, PetscBool useClosure, PetscInt *closureSize, PetscInt **closure) 210 { 211 PetscErrorCode ierr; 212 213 PetscFunctionBegin; 214 ierr = DMRestoreWorkArray(dm, 0, PETSC_INT, closure);CHKERRQ(ierr); 215 PetscFunctionReturn(0); 216 } 217 218 #undef __FUNCT__ 219 #define __FUNCT__ "DMDAGetClosure" 220 PetscErrorCode DMDAGetClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure) 221 { 222 DM_DA *da = (DM_DA*) dm->data; 223 PetscInt dim = da->dim; 224 PetscInt nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF; 225 PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart; 226 PetscErrorCode ierr; 227 228 PetscFunctionBegin; 229 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 230 if (n) PetscValidIntPointer(n,4); 231 if (closure) PetscValidPointer(closure, 5); 232 if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 233 if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection"); 234 ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 235 ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 236 ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 237 ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 238 ierr = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr); 239 ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 240 xfStart = fStart; xfEnd = xfStart+nXF; 241 yfStart = xfEnd; 242 if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd); 243 if ((p >= cStart) || (p < cEnd)) { 244 /* Cell */ 245 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 246 else if (dim == 2) { 247 /* 4 faces, 4 vertices 248 Bottom-left vertex follows same order as cells 249 Bottom y-face same order as cells 250 Left x-face follows same order as cells 251 We number the quad: 252 253 8--3--7 254 | | 255 4 0 2 256 | | 257 5--1--6 258 */ 259 PetscInt c = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1); 260 PetscInt v = cy*nVx + cx + vStart; 261 PetscInt xf = cy*nxF + cx + xfStart; 262 PetscInt yf = c + yfStart; 263 PetscInt points[9]; 264 265 /* Note: initializer list is not computable at compile time, hence we fill the array manually */ 266 points[0] = p; points[1] = yf; points[2] = xf+1; points[3] = yf+nyF; points[4] = xf+0; points[5] = v+0; points[6]= v+1; points[7] = v+nVx+1; points[8] = v+nVx+0; 267 268 ierr = GetPointArray_Private(dm,9,points,n,closure);CHKERRQ(ierr); 269 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 270 } else if ((p >= vStart) || (p < vEnd)) { 271 /* Vertex */ 272 ierr = GetPointArray_Private(dm,1,&p,n,closure);CHKERRQ(ierr); 273 } else if ((p >= fStart) || (p < fStart + nXF)) { 274 /* X Face */ 275 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 276 else if (dim == 2) { 277 /* 2 vertices: The bottom vertex has the same numbering as the face */ 278 PetscInt f = p - xfStart; 279 PetscInt points[3]; 280 281 points[0] = p; points[1] = f; points[2] = f+nVx; 282 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken"); 283 ierr = GetPointArray_Private(dm,3,points,n,closure);CHKERRQ(ierr); 284 } else if (dim == 3) { 285 /* 4 vertices */ 286 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 287 } 288 } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) { 289 /* Y Face */ 290 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 291 else if (dim == 2) { 292 /* 2 vertices: The left vertex has the same numbering as the face */ 293 PetscInt f = p - yfStart; 294 PetscInt points[3]; 295 296 points[0] = p; points[1] = f; points[2]= f+1; 297 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken"); 298 ierr = GetPointArray_Private(dm, 3, points, n, closure);CHKERRQ(ierr); 299 } else if (dim == 3) { 300 /* 4 vertices */ 301 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 302 } 303 } else { 304 /* Z Face */ 305 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 306 else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D"); 307 else if (dim == 3) { 308 /* 4 vertices */ 309 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 310 } 311 } 312 PetscFunctionReturn(0); 313 } 314 315 #undef __FUNCT__ 316 #define __FUNCT__ "DMDARestoreClosure" 317 /* Zeros n and closure. */ 318 PetscErrorCode DMDARestoreClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure) 319 { 320 PetscErrorCode ierr; 321 322 PetscFunctionBegin; 323 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 324 if (n) PetscValidIntPointer(n,4); 325 if (closure) PetscValidPointer(closure, 5); 326 ierr = RestorePointArray_Private(dm,n,closure);CHKERRQ(ierr); 327 PetscFunctionReturn(0); 328 } 329 330 #undef __FUNCT__ 331 #define __FUNCT__ "DMDAGetClosureScalar" 332 PetscErrorCode DMDAGetClosureScalar(DM dm, PetscSection section, PetscInt p, PetscScalar *vArray, PetscInt *csize, const PetscScalar **values) 333 { 334 DM_DA *da = (DM_DA*) dm->data; 335 PetscInt dim = da->dim; 336 PetscInt nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF; 337 PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart; 338 PetscErrorCode ierr; 339 340 PetscFunctionBegin; 341 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 342 PetscValidScalarPointer(vArray, 4); 343 if (values) PetscValidPointer(values, 6); 344 if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 345 if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection"); 346 ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 347 ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 348 ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 349 ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 350 ierr = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr); 351 ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 352 xfStart = fStart; xfEnd = xfStart+nXF; 353 yfStart = xfEnd; yfEnd = yfStart+nYF; 354 zfStart = yfEnd; 355 if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd); 356 if ((p >= cStart) || (p < cEnd)) { 357 /* Cell */ 358 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 359 else if (dim == 2) { 360 /* 4 faces, 4 vertices 361 Bottom-left vertex follows same order as cells 362 Bottom y-face same order as cells 363 Left x-face follows same order as cells 364 We number the quad: 365 366 8--3--7 367 | | 368 4 0 2 369 | | 370 5--1--6 371 */ 372 PetscInt c = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1); 373 PetscInt v = cy*nVx + cx + vStart; 374 PetscInt xf = cy*nxF + cx + xfStart; 375 PetscInt yf = c + yfStart; 376 PetscInt points[9]; 377 378 points[0] = p; points[1] = yf; points[2] = xf+1; points[3] = yf+nyF; points[4] = xf+0; points[5] = v+0; points[6] = v+1; points[7] = v+nVx+1; points[8] = v+nVx+0; 379 ierr = FillClosureArray_Static(dm, section, 9, points, vArray, csize, values);CHKERRQ(ierr); 380 } else { 381 /* 6 faces, 8 vertices 382 Bottom-left-back vertex follows same order as cells 383 Back z-face follows same order as cells 384 Bottom y-face follows same order as cells 385 Left x-face follows same order as cells 386 387 14-----13 388 /| /| 389 / | 2 / | 390 / 5| / | 391 10-----9 4| 392 | 11-|---12 393 |6 / | / 394 | /1 3| / 395 |/ |/ 396 7-----8 397 */ 398 PetscInt c = p - cStart; 399 PetscInt points[15]; 400 401 points[0] = p; points[1] = c+zfStart; points[2] = c+zfStart+nzF; points[3] = c+yfStart; points[4] = c+xfStart+nxF; points[5] = c+yfStart+nyF; points[6] = c+xfStart; 402 points[7] = c+vStart+0; points[8] = c+vStart+1; points[9] = c+vStart+nVx+1; points[10] = c+vStart+nVx+0; points[11] = c+vStart+nVx*nVy+0; 403 points[12] = c+vStart+nVx*nVy+1; points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0; 404 405 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken"); 406 ierr = FillClosureArray_Static(dm, section, 15, points, vArray, csize, values);CHKERRQ(ierr); 407 } 408 } else if ((p >= vStart) || (p < vEnd)) { 409 /* Vertex */ 410 ierr = FillClosureArray_Static(dm, section, 1, &p, vArray, csize, values);CHKERRQ(ierr); 411 } else if ((p >= fStart) || (p < fStart + nXF)) { 412 /* X Face */ 413 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 414 else if (dim == 2) { 415 /* 2 vertices: The bottom vertex has the same numbering as the face */ 416 PetscInt f = p - xfStart; 417 PetscInt points[3]; 418 419 points[0] = p; points[1] = f; points[2] = f+nVx; 420 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken"); 421 ierr = FillClosureArray_Static(dm, section, 3, points, vArray, csize, values);CHKERRQ(ierr); 422 } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 423 } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) { 424 /* Y Face */ 425 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 426 else if (dim == 2) { 427 /* 2 vertices: The left vertex has the same numbering as the face */ 428 PetscInt f = p - yfStart; 429 PetscInt points[3]; 430 431 points[0] = p; points[1] = f; points[2] = f+1; 432 SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken"); 433 ierr = FillClosureArray_Static(dm, section, 3, points, vArray, csize, values);CHKERRQ(ierr); 434 } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 435 } else { 436 /* Z Face */ 437 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 438 else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D"); 439 else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 440 } 441 PetscFunctionReturn(0); 442 } 443 444 #undef __FUNCT__ 445 #define __FUNCT__ "DMDAVecGetClosure" 446 PetscErrorCode DMDAVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt p, PetscInt *csize, const PetscScalar **values) 447 { 448 PetscScalar *vArray; 449 PetscErrorCode ierr; 450 451 PetscFunctionBegin; 452 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 453 PetscValidHeaderSpecific(v, VEC_CLASSID, 3); 454 if (values) PetscValidPointer(values, 6); 455 ierr = VecGetArray(v, &vArray);CHKERRQ(ierr); 456 ierr = DMDAGetClosureScalar(dm, section, p, vArray, csize, values);CHKERRQ(ierr); 457 ierr = VecRestoreArray(v, &vArray);CHKERRQ(ierr); 458 PetscFunctionReturn(0); 459 } 460 461 #undef __FUNCT__ 462 #define __FUNCT__ "DMDARestoreClosureScalar" 463 PetscErrorCode DMDARestoreClosureScalar(DM dm, PetscSection section, PetscInt p, PetscScalar *vArray, PetscInt *csize, const PetscScalar **values) 464 { 465 PetscErrorCode ierr; 466 467 PetscFunctionBegin; 468 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 469 PetscValidPointer(values, 6); 470 ierr = DMRestoreWorkArray(dm, csize ? *csize : 0, PETSC_SCALAR, (void*) values);CHKERRQ(ierr); 471 PetscFunctionReturn(0); 472 } 473 474 #undef __FUNCT__ 475 #define __FUNCT__ "DMDAVecRestoreClosure" 476 PetscErrorCode DMDAVecRestoreClosure(DM dm, PetscSection section, Vec v, PetscInt p, PetscInt *csize, const PetscScalar **values) 477 { 478 PetscErrorCode ierr; 479 480 PetscFunctionBegin; 481 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 482 PetscValidHeaderSpecific(v, VEC_CLASSID, 3); 483 PetscValidPointer(values, 5); 484 ierr = DMDARestoreClosureScalar(dm, section, p, NULL, csize, values);CHKERRQ(ierr); 485 PetscFunctionReturn(0); 486 } 487 488 #undef __FUNCT__ 489 #define __FUNCT__ "DMDASetClosureScalar" 490 PetscErrorCode DMDASetClosureScalar(DM dm, PetscSection section, PetscInt p,PetscScalar *vArray, const PetscScalar *values, InsertMode mode) 491 { 492 DM_DA *da = (DM_DA*) dm->data; 493 PetscInt dim = da->dim; 494 PetscInt nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF, nCx, nCy; 495 PetscInt pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart; 496 PetscErrorCode ierr; 497 498 PetscFunctionBegin; 499 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 500 PetscValidScalarPointer(values, 4); 501 PetscValidPointer(values, 5); 502 if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 503 if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection"); 504 ierr = DMDAGetHeightStratum(dm, -1, &pStart, &pEnd);CHKERRQ(ierr); 505 ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 506 ierr = DMDAGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 507 ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr); 508 ierr = DMDAGetNumCells(dm, &nCx, &nCy, NULL, NULL);CHKERRQ(ierr); 509 ierr = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr); 510 ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr); 511 xfStart = fStart; xfEnd = xfStart+nXF; 512 yfStart = xfEnd; yfEnd = yfStart+nYF; 513 zfStart = yfEnd; 514 if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd); 515 if ((p >= cStart) || (p < cEnd)) { 516 /* Cell */ 517 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 518 else if (dim == 2) { 519 /* 4 faces, 4 vertices 520 Bottom-left vertex follows same order as cells 521 Bottom y-face same order as cells 522 Left x-face follows same order as cells 523 We number the quad: 524 525 8--3--7 526 | | 527 4 0 2 528 | | 529 5--1--6 530 */ 531 PetscInt c = p - cStart, l = c/nCx; 532 PetscInt points[9]; 533 534 points[0] = p; 535 points[1] = yfStart+c+0; points[2] = xfStart+l+c+1; points[3] = yfStart+nyF+c+0; points[4] = xfStart+l+c+0; 536 points[5] = vStart+l+c+0; points[6] = vStart+l+c+1; points[7] = vStart+nVx+l+c+1; points[8] = vStart+nVx+l+c+0; 537 ierr = FillClosureVec_Private(dm, section, 9, points, vArray, values, mode);CHKERRQ(ierr); 538 } else { 539 /* 6 faces, 8 vertices 540 Bottom-left-back vertex follows same order as cells 541 Back z-face follows same order as cells 542 Bottom y-face follows same order as cells 543 Left x-face follows same order as cells 544 545 14-----13 546 /| /| 547 / | 2 / | 548 / 5| / | 549 10-----9 4| 550 | 11-|---12 551 |6 / | / 552 | /1 3| / 553 |/ |/ 554 7-----8 555 */ 556 PetscInt c = p - cStart; 557 PetscInt points[15]; 558 559 points[0] = p; points[1] = c+zfStart; points[2] = c+zfStart+nzF; points[3] = c+yfStart; points[4] = c+xfStart+nxF; points[5] = c+yfStart+nyF; points[6] = c+xfStart; 560 points[7] = c+vStart+0; points[8] = c+vStart+1; points[9] = c+vStart+nVx+1; points[10] = c+vStart+nVx+0; points[11] = c+vStart+nVx*nVy+0; points[12] = c+vStart+nVx*nVy+1; 561 points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0; 562 ierr = FillClosureVec_Private(dm, section, 15, points, vArray, values, mode);CHKERRQ(ierr); 563 } 564 } else if ((p >= vStart) || (p < vEnd)) { 565 /* Vertex */ 566 ierr = FillClosureVec_Private(dm, section, 1, &p, vArray, values, mode);CHKERRQ(ierr); 567 } else if ((p >= fStart) || (p < fStart + nXF)) { 568 /* X Face */ 569 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 570 else if (dim == 2) { 571 /* 2 vertices: The bottom vertex has the same numbering as the face */ 572 PetscInt f = p - xfStart; 573 PetscInt points[3]; 574 575 points[0] = p; points[1] = f; points[2] = f+nVx; 576 ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr); 577 } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 578 } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) { 579 /* Y Face */ 580 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 581 else if (dim == 2) { 582 /* 2 vertices: The left vertex has the same numbering as the face */ 583 PetscInt f = p - yfStart; 584 PetscInt points[3]; 585 586 points[0] = p; points[1] = f; points[2] = f+1; 587 ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr); 588 } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 589 } else { 590 /* Z Face */ 591 if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D"); 592 else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D"); 593 else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented"); 594 } 595 PetscFunctionReturn(0); 596 } 597 598 #undef __FUNCT__ 599 #define __FUNCT__ "DMDAVecSetClosure" 600 PetscErrorCode DMDAVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar *values, InsertMode mode) 601 { 602 PetscScalar *vArray; 603 PetscErrorCode ierr; 604 605 PetscFunctionBegin; 606 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 607 PetscValidHeaderSpecific(v, VEC_CLASSID, 3); 608 PetscValidPointer(values, 5); 609 ierr = VecGetArray(v,&vArray);CHKERRQ(ierr); 610 ierr = DMDASetClosureScalar(dm,section,p,vArray,values,mode);CHKERRQ(ierr); 611 ierr = VecRestoreArray(v,&vArray);CHKERRQ(ierr); 612 PetscFunctionReturn(0); 613 } 614 615 #undef __FUNCT__ 616 #define __FUNCT__ "DMDACnvertToCell" 617 /*@ 618 DMDAConvertToCell - Convert (i,j,k) to local cell number 619 620 Not Collective 621 622 Input Parameter: 623 + da - the distributed array 624 = s - A MatStencil giving (i,j,k) 625 626 Output Parameter: 627 . cell - the local cell number 628 629 Level: developer 630 631 .seealso: DMDAVecGetClosure() 632 @*/ 633 PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell) 634 { 635 DM_DA *da = (DM_DA*) dm->data; 636 const PetscInt dim = da->dim; 637 const PetscInt mx = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys /*, mz = da->Ze - da->Zs*/; 638 const PetscInt il = s.i - da->Xs/da->w, jl = dim > 1 ? s.j - da->Ys : 0, kl = dim > 2 ? s.k - da->Zs : 0; 639 640 PetscFunctionBegin; 641 *cell = -1; 642 if ((s.i < da->Xs/da->w) || (s.i >= da->Xe/da->w)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil i %d should be in [%d, %d)", s.i, da->Xs/da->w, da->Xe/da->w); 643 if ((dim > 1) && ((s.j < da->Ys) || (s.j >= da->Ye))) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil j %d should be in [%d, %d)", s.j, da->Ys, da->Ye); 644 if ((dim > 2) && ((s.k < da->Zs) || (s.k >= da->Ze))) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil k %d should be in [%d, %d)", s.k, da->Zs, da->Ze); 645 *cell = (kl*my + jl)*mx + il; 646 PetscFunctionReturn(0); 647 } 648 649 #undef __FUNCT__ 650 #define __FUNCT__ "DMDAComputeCellGeometry_2D" 651 PetscErrorCode DMDAComputeCellGeometry_2D(DM dm, const PetscScalar vertices[], const PetscReal refPoint[], PetscReal J[], PetscReal invJ[], PetscReal *detJ) 652 { 653 const PetscScalar x0 = vertices[0]; 654 const PetscScalar y0 = vertices[1]; 655 const PetscScalar x1 = vertices[2]; 656 const PetscScalar y1 = vertices[3]; 657 const PetscScalar x2 = vertices[4]; 658 const PetscScalar y2 = vertices[5]; 659 const PetscScalar x3 = vertices[6]; 660 const PetscScalar y3 = vertices[7]; 661 const PetscScalar f_01 = x2 - x1 - x3 + x0; 662 const PetscScalar g_01 = y2 - y1 - y3 + y0; 663 const PetscScalar x = refPoint[0]; 664 const PetscScalar y = refPoint[1]; 665 PetscReal invDet; 666 PetscErrorCode ierr; 667 668 PetscFunctionBegin; 669 #if 0 670 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell (%g,%g)--(%g,%g)--(%g,%g)--(%g,%g)\n", 671 PetscRealPart(x0),PetscRealPart(y0),PetscRealPart(x1),PetscRealPart(y1),PetscRealPart(x2),PetscRealPart(y2),PetscRealPart(x3),PetscRealPart(y3));CHKERRQ(ierr); 672 ierr = PetscPrintf(PETSC_COMM_SELF, "Ref Point (%g,%g)\n", PetscRealPart(x), PetscRealPart(y));CHKERRQ(ierr); 673 #endif 674 J[0] = PetscRealPart(x1 - x0 + f_01*y) * 0.5; J[1] = PetscRealPart(x3 - x0 + f_01*x) * 0.5; 675 J[2] = PetscRealPart(y1 - y0 + g_01*y) * 0.5; J[3] = PetscRealPart(y3 - y0 + g_01*x) * 0.5; 676 *detJ = J[0]*J[3] - J[1]*J[2]; 677 invDet = 1.0/(*detJ); 678 if (invJ) { 679 invJ[0] = invDet*J[3]; invJ[1] = -invDet*J[1]; 680 invJ[2] = -invDet*J[2]; invJ[3] = invDet*J[0]; 681 } 682 ierr = PetscLogFlops(30);CHKERRQ(ierr); 683 PetscFunctionReturn(0); 684 } 685 686 #undef __FUNCT__ 687 #define __FUNCT__ "DMDAComputeCellGeometry" 688 PetscErrorCode DMDAComputeCellGeometry(DM dm, PetscInt cell, PetscQuadrature *quad, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal detJ[]) 689 { 690 DM cdm; 691 Vec coordinates; 692 const PetscScalar *vertices = NULL; 693 PetscInt csize, dim, d, q; 694 PetscErrorCode ierr; 695 696 PetscFunctionBegin; 697 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 698 ierr = DMDAGetInfo(dm, &dim, 0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 699 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 700 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 701 ierr = DMDAVecGetClosure(cdm, NULL, coordinates, cell, &csize, &vertices);CHKERRQ(ierr); 702 for (d = 0; d < dim; ++d) v0[d] = PetscRealPart(vertices[d]); 703 switch (dim) { 704 case 2: 705 for (q = 0; q < quad->numPoints; ++q) { 706 ierr = DMDAComputeCellGeometry_2D(dm, vertices, &quad->points[q*dim], J, invJ, detJ);CHKERRQ(ierr); 707 } 708 break; 709 default: 710 SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Dimension %d not supported", dim); 711 } 712 ierr = DMDAVecRestoreClosure(cdm, NULL, coordinates, cell, &csize, &vertices);CHKERRQ(ierr); 713 PetscFunctionReturn(0); 714 } 715