Lines Matching +full:- +full:m

13   DM_DA      *dd = (DM_DA *)da->data;  in DMView_DA_3d()
50 …PetscCall(PetscViewerASCIIPrintf(viewer, " Load Balance - Grid Points: Min %" PetscInt_FMT " avg… in DMView_DA_3d()
56M %" PetscInt_FMT " N %" PetscInt_FMT " P %" PetscInt_FMT " m %" PetscInt_FMT " n %" PetscInt_FMT … in DMView_DA_3d()
66 last = last - 3; in DMView_DA_3d()
76 PetscReal ymin = -1.0, ymax = (PetscReal)dd->N; in DMView_DA_3d()
77 PetscReal xmin = -1.0, xmax = (PetscReal)((dd->M + 2) * dd->P), x, y, ycoord, xcoord; in DMView_DA_3d()
94 for (k = 0; k < dd->P; k++) { in DMView_DA_3d()
96 ymax = (PetscReal)(dd->N - 1); in DMView_DA_3d()
97 …for (xmin = (PetscReal)(k * (dd->M + 1)); xmin < (PetscReal)(dd->M + (k * (dd->M + 1))); xmin++) P… in DMView_DA_3d()
98 xmin = (PetscReal)(k * (dd->M + 1)); in DMView_DA_3d()
99 xmax = xmin + (PetscReal)(dd->M - 1); in DMView_DA_3d()
100 …for (ymin = 0; ymin < (PetscReal)dd->N; ymin++) PetscCall(PetscDrawLine(draw, xmin, ymin, xmax, ym… in DMView_DA_3d()
109 for (k = 0; k < dd->P; k++) { in DMView_DA_3d()
110 if ((k >= dd->zs) && (k < dd->ze)) { in DMView_DA_3d()
112 ymin = dd->ys; in DMView_DA_3d()
113 ymax = dd->ye - 1; in DMView_DA_3d()
114 xmin = dd->xs / dd->w + (dd->M + 1) * k; in DMView_DA_3d()
115 xmax = (dd->xe - 1) / dd->w + (dd->M + 1) * k; in DMView_DA_3d()
122 xmin = dd->xs / dd->w; in DMView_DA_3d()
123 xmax = (dd->xe - 1) / dd->w; in DMView_DA_3d()
127 … PetscCall(PetscDrawString(draw, xmin + (dd->M + 1) * k + .2, ymin + .3, PETSC_DRAW_RED, node)); in DMView_DA_3d()
129 base = (dd->base + (dd->xe - dd->xs) * (dd->ye - dd->ys) * (k - dd->zs)) / dd->w; in DMView_DA_3d()
131 for (x = xmin + (dd->M + 1) * k; x <= xmax + (dd->M + 1) * k; x++) { in DMView_DA_3d()
143 for (k = 0 - dd->s; k < dd->P + dd->s; k++) { in DMView_DA_3d()
145 if ((k >= dd->Zs) && (k < dd->Ze)) { in DMView_DA_3d()
147 base = (dd->Xe - dd->Xs) * (dd->Ye - dd->Ys) * (k - dd->Zs) / dd->w; in DMView_DA_3d()
148 PetscCall(ISLocalToGlobalMappingGetBlockIndices(da->ltogmap, &idx)); in DMView_DA_3d()
151 if (k < 0) plane = dd->P + k; in DMView_DA_3d()
152 if (k >= dd->P) plane = k - dd->P; in DMView_DA_3d()
153 ymin = dd->Ys; in DMView_DA_3d()
154 ymax = dd->Ye; in DMView_DA_3d()
155 xmin = (dd->M + 1) * plane * dd->w; in DMView_DA_3d()
156 xmax = (dd->M + 1) * plane * dd->w + dd->M * dd->w; in DMView_DA_3d()
158 for (x = xmin + dd->Xs; x < xmin + dd->Xe; x += dd->w) { in DMView_DA_3d()
162 if (y < 0) ycoord = dd->N + y; in DMView_DA_3d()
163 if (y >= dd->N) ycoord = y - dd->N; in DMView_DA_3d()
165 if (x < xmin) xcoord = xmax - (xmin - x); in DMView_DA_3d()
166 if (x >= xmax) xcoord = xmin + (x - xmax); in DMView_DA_3d()
167 PetscCall(PetscDrawString(draw, xcoord / dd->w, ycoord, PETSC_DRAW_BLUE, node)); in DMView_DA_3d()
171 PetscCall(ISLocalToGlobalMappingRestoreBlockIndices(da->ltogmap, &idx)); in DMView_DA_3d()
192 DM_DA *dd = (DM_DA *)da->data; in DMSetUp_DA_3D()
193 const PetscInt M = dd->M; in DMSetUp_DA_3D() local
194 const PetscInt N = dd->N; in DMSetUp_DA_3D()
195 const PetscInt P = dd->P; in DMSetUp_DA_3D()
196 PetscMPIInt m, n, p; in DMSetUp_DA_3D() local
197 const PetscInt dof = dd->w; in DMSetUp_DA_3D()
198 const PetscInt s = dd->s; in DMSetUp_DA_3D()
199 DMBoundaryType bx = dd->bx; in DMSetUp_DA_3D()
200 DMBoundaryType by = dd->by; in DMSetUp_DA_3D()
201 DMBoundaryType bz = dd->bz; in DMSetUp_DA_3D()
202 DMDAStencilType stencil_type = dd->stencil_type; in DMSetUp_DA_3D()
203 PetscInt *lx = dd->lx; in DMSetUp_DA_3D()
204 PetscInt *ly = dd->ly; in DMSetUp_DA_3D()
205 PetscInt *lz = dd->lz; in DMSetUp_DA_3D()
226M) * ((PetscInt64)N) * ((PetscInt64)P) * ((PetscInt64)dof) <= (PetscInt64)PETSC_MPI_INT_MAX, comm,… in DMSetUp_DA_3D()
228 PetscCall(PetscMPIIntCast(dd->m, &m)); in DMSetUp_DA_3D()
229 PetscCall(PetscMPIIntCast(dd->n, &n)); in DMSetUp_DA_3D()
230 PetscCall(PetscMPIIntCast(dd->p, &p)); in DMSetUp_DA_3D()
235 if (m != PETSC_DECIDE) { in DMSetUp_DA_3D()
236 …PetscCheck(m >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors i… in DMSetUp_DA_3D()
237 …PetscCheck(m <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many processors in X directi… in DMSetUp_DA_3D()
240 …PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors i… in DMSetUp_DA_3D()
244 …PetscCheck(p >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Non-positive number of processors i… in DMSetUp_DA_3D()
247 …PetscCheck(m <= 0 || n <= 0 || p <= 0 || m * n * p == size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRA… in DMSetUp_DA_3D()
250 if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) { in DMSetUp_DA_3D()
251 m = size / (n * p); in DMSetUp_DA_3D()
252 } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { in DMSetUp_DA_3D()
253 n = size / (m * p); in DMSetUp_DA_3D()
254 } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { in DMSetUp_DA_3D()
255 p = size / (m * n); in DMSetUp_DA_3D()
256 } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) { in DMSetUp_DA_3D()
258 m = (int)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)N * p))); in DMSetUp_DA_3D()
259 if (!m) m = 1; in DMSetUp_DA_3D()
260 while (m > 0) { in DMSetUp_DA_3D()
261 n = size / (m * p); in DMSetUp_DA_3D()
262 if (m * n * p == size) break; in DMSetUp_DA_3D()
263 m--; in DMSetUp_DA_3D()
265 PetscCheck(m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad p value: p = %d", p); in DMSetUp_DA_3D()
266 if (M > N && m < n) { in DMSetUp_DA_3D()
267 PetscMPIInt _m = m; in DMSetUp_DA_3D()
268 m = n; in DMSetUp_DA_3D()
271 } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) { in DMSetUp_DA_3D()
273 m = (int)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)P * n))); in DMSetUp_DA_3D()
274 if (!m) m = 1; in DMSetUp_DA_3D()
275 while (m > 0) { in DMSetUp_DA_3D()
276 p = size / (m * n); in DMSetUp_DA_3D()
277 if (m * n * p == size) break; in DMSetUp_DA_3D()
278 m--; in DMSetUp_DA_3D()
280 PetscCheck(m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad n value: n = %d", n); in DMSetUp_DA_3D()
281 if (M > P && m < p) { in DMSetUp_DA_3D()
282 PetscMPIInt _m = m; in DMSetUp_DA_3D()
283 m = p; in DMSetUp_DA_3D()
286 } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { in DMSetUp_DA_3D()
288 n = (int)(0.5 + PetscSqrtReal(((PetscReal)N) * ((PetscReal)size) / ((PetscReal)P * m))); in DMSetUp_DA_3D()
291 p = size / (m * n); in DMSetUp_DA_3D()
292 if (m * n * p == size) break; in DMSetUp_DA_3D()
293 n--; in DMSetUp_DA_3D()
295 PetscCheck(n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "bad m value: m = %d", n); in DMSetUp_DA_3D()
301 } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) { in DMSetUp_DA_3D()
303 …n = (PetscMPIInt)(0.5 + PetscPowReal(((PetscReal)N * N) * ((PetscReal)size) / ((PetscReal)P * M), … in DMSetUp_DA_3D()
308 n--; in DMSetUp_DA_3D()
311 m = (PetscMPIInt)(0.5 + PetscSqrtReal(((PetscReal)M) * ((PetscReal)size) / ((PetscReal)P * n))); in DMSetUp_DA_3D()
312 if (!m) m = 1; in DMSetUp_DA_3D()
313 while (m > 0) { in DMSetUp_DA_3D()
314 p = size / (m * n); in DMSetUp_DA_3D()
315 if (m * n * p == size) break; in DMSetUp_DA_3D()
316 m--; in DMSetUp_DA_3D()
318 if (M > P && m < p) { in DMSetUp_DA_3D()
319 PetscMPIInt _m = m; in DMSetUp_DA_3D()
320 m = p; in DMSetUp_DA_3D()
323 …} else PetscCheck(m * n * p == size, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "… in DMSetUp_DA_3D()
325 …PetscCheck(m * n * p == size, PetscObjectComm((PetscObject)da), PETSC_ERR_PLIB, "Could not find go… in DMSetUp_DA_3D()
326 …tscCheck(M >= m, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_OUTOFRANGE, "Partition in x direc… in DMSetUp_DA_3D()
336 PetscCall(PetscMalloc1(m, &dd->lx)); in DMSetUp_DA_3D()
337 lx = dd->lx; in DMSetUp_DA_3D()
338 for (i = 0; i < m; i++) lx[i] = M / m + ((M % m) > (i % m)); in DMSetUp_DA_3D()
340 x = lx[rank % m]; in DMSetUp_DA_3D()
342 for (i = 0; i < (rank % m); i++) xs += lx[i]; in DMSetUp_DA_3D()
343 …PetscCheck(x >= s || (m <= 1 && bx != DM_BOUNDARY_PERIODIC), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFR… in DMSetUp_DA_3D()
346 PetscCall(PetscMalloc1(n, &dd->ly)); in DMSetUp_DA_3D()
347 ly = dd->ly; in DMSetUp_DA_3D()
350 y = ly[(rank % (m * n)) / m]; in DMSetUp_DA_3D()
351 …= DM_BOUNDARY_PERIODIC), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local y-width of domain y %" … in DMSetUp_DA_3D()
354 for (i = 0; i < (rank % (m * n)) / m; i++) ys += ly[i]; in DMSetUp_DA_3D()
357 PetscCall(PetscMalloc1(p, &dd->lz)); in DMSetUp_DA_3D()
358 lz = dd->lz; in DMSetUp_DA_3D()
361 z = lz[rank / (m * n)]; in DMSetUp_DA_3D()
363 … != DM_BOUNDARY_MIRROR), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local x-width of domain x %" … in DMSetUp_DA_3D()
364 … != DM_BOUNDARY_MIRROR), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local y-width of domain y %" … in DMSetUp_DA_3D()
365 … != DM_BOUNDARY_MIRROR), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local z-width of domain z %" … in DMSetUp_DA_3D()
367 /* note this is different than x- and y-, as we will handle as an important special in DMSetUp_DA_3D()
372 …= DM_BOUNDARY_PERIODIC), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local z-width of domain z %" … in DMSetUp_DA_3D()
374 for (i = 0; i < (rank / (m * n)); i++) zs += lz[i]; in DMSetUp_DA_3D()
380 if (xs - s > 0) { in DMSetUp_DA_3D()
381 Xs = xs - s; in DMSetUp_DA_3D()
382 IXs = xs - s; in DMSetUp_DA_3D()
384 if (bx) Xs = xs - s; in DMSetUp_DA_3D()
388 if (xe + s <= M) { in DMSetUp_DA_3D()
393 Xs = xs - s; in DMSetUp_DA_3D()
395 } else Xe = M; in DMSetUp_DA_3D()
396 IXe = M; in DMSetUp_DA_3D()
400 IXs = xs - s; in DMSetUp_DA_3D()
402 Xs = xs - s; in DMSetUp_DA_3D()
406 if (ys - s > 0) { in DMSetUp_DA_3D()
407 Ys = ys - s; in DMSetUp_DA_3D()
408 IYs = ys - s; in DMSetUp_DA_3D()
410 if (by) Ys = ys - s; in DMSetUp_DA_3D()
424 IYs = ys - s; in DMSetUp_DA_3D()
426 Ys = ys - s; in DMSetUp_DA_3D()
430 if (zs - s > 0) { in DMSetUp_DA_3D()
431 Zs = zs - s; in DMSetUp_DA_3D()
432 IZs = zs - s; in DMSetUp_DA_3D()
434 if (bz) Zs = zs - s; in DMSetUp_DA_3D()
448 IZs = zs - s; in DMSetUp_DA_3D()
450 Zs = zs - s; in DMSetUp_DA_3D()
464 for (i = 1; i <= size; i++) bases[i] = ldims[i - 1]; in DMSetUp_DA_3D()
465 for (i = 1; i <= size; i++) bases[i] += bases[i - 1]; in DMSetUp_DA_3D()
469 dd->Nlocal = x * y * z * dof; in DMSetUp_DA_3D()
470 PetscCall(VecCreateMPIWithArray(comm, dof, dd->Nlocal, PETSC_DECIDE, NULL, &global)); in DMSetUp_DA_3D()
471 dd->nlocal = (Xe - Xs) * (Ye - Ys) * (Ze - Zs) * dof; in DMSetUp_DA_3D()
472 PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, dof, dd->nlocal, NULL, &local)); in DMSetUp_DA_3D()
478 PetscCall(PetscMalloc1((IXe - IXs) * (IYe - IYs) * (IZe - IZs), &idx)); in DMSetUp_DA_3D()
480 left = IXs - Xs; in DMSetUp_DA_3D()
481 right = left + (IXe - IXs); in DMSetUp_DA_3D()
482 bottom = IYs - Ys; in DMSetUp_DA_3D()
483 top = bottom + (IYe - IYs); in DMSetUp_DA_3D()
484 down = IZs - Zs; in DMSetUp_DA_3D()
485 up = down + (IZe - IZs); in DMSetUp_DA_3D()
489 for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k; in DMSetUp_DA_3D()
495 left = xs - Xs; in DMSetUp_DA_3D()
497 bottom = ys - Ys; in DMSetUp_DA_3D()
499 down = zs - Zs; in DMSetUp_DA_3D()
503 for (i = (IZs - Zs); i < down; i++) { in DMSetUp_DA_3D()
505 for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k; in DMSetUp_DA_3D()
511 for (j = (IYs - Ys); j < bottom; j++) { in DMSetUp_DA_3D()
512 for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k; in DMSetUp_DA_3D()
516 for (k = IXs - Xs; k < IXe - Xs; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k; in DMSetUp_DA_3D()
519 for (j = top; j < top + IYe - ye; j++) { in DMSetUp_DA_3D()
520 for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k; in DMSetUp_DA_3D()
524 for (i = up; i < up + IZe - ze; i++) { in DMSetUp_DA_3D()
526 for (k = left; k < right; k++) idx[count++] = (i * (Ye - Ys) + j) * (Xe - Xs) + k; in DMSetUp_DA_3D()
547 n0 = rank - m * n - m - 1; in DMSetUp_DA_3D()
548 n1 = rank - m * n - m; in DMSetUp_DA_3D()
549 n2 = rank - m * n - m + 1; in DMSetUp_DA_3D()
550 n3 = rank - m * n - 1; in DMSetUp_DA_3D()
551 n4 = rank - m * n; in DMSetUp_DA_3D()
552 n5 = rank - m * n + 1; in DMSetUp_DA_3D()
553 n6 = rank - m * n + m - 1; in DMSetUp_DA_3D()
554 n7 = rank - m * n + m; in DMSetUp_DA_3D()
555 n8 = rank - m * n + m + 1; in DMSetUp_DA_3D()
557 n9 = rank - m - 1; in DMSetUp_DA_3D()
558 n10 = rank - m; in DMSetUp_DA_3D()
559 n11 = rank - m + 1; in DMSetUp_DA_3D()
560 n12 = rank - 1; in DMSetUp_DA_3D()
562 n15 = rank + m - 1; in DMSetUp_DA_3D()
563 n16 = rank + m; in DMSetUp_DA_3D()
564 n17 = rank + m + 1; in DMSetUp_DA_3D()
566 n18 = rank + m * n - m - 1; in DMSetUp_DA_3D()
567 n19 = rank + m * n - m; in DMSetUp_DA_3D()
568 n20 = rank + m * n - m + 1; in DMSetUp_DA_3D()
569 n21 = rank + m * n - 1; in DMSetUp_DA_3D()
570 n22 = rank + m * n; in DMSetUp_DA_3D()
571 n23 = rank + m * n + 1; in DMSetUp_DA_3D()
572 n24 = rank + m * n + m - 1; in DMSetUp_DA_3D()
573 n25 = rank + m * n + m; in DMSetUp_DA_3D()
574 n26 = rank + m * n + m + 1; in DMSetUp_DA_3D()
579 n0 = rank - 1 - (m * n); in DMSetUp_DA_3D()
580 n3 = rank + m - 1 - (m * n); in DMSetUp_DA_3D()
581 n6 = rank + 2 * m - 1 - (m * n); in DMSetUp_DA_3D()
582 n9 = rank - 1; in DMSetUp_DA_3D()
583 n12 = rank + m - 1; in DMSetUp_DA_3D()
584 n15 = rank + 2 * m - 1; in DMSetUp_DA_3D()
585 n18 = rank - 1 + (m * n); in DMSetUp_DA_3D()
586 n21 = rank + m - 1 + (m * n); in DMSetUp_DA_3D()
587 n24 = rank + 2 * m - 1 + (m * n); in DMSetUp_DA_3D()
590 if (xe == M) { /* First assume not corner or edge */ in DMSetUp_DA_3D()
591 n2 = rank - 2 * m + 1 - (m * n); in DMSetUp_DA_3D()
592 n5 = rank - m + 1 - (m * n); in DMSetUp_DA_3D()
593 n8 = rank + 1 - (m * n); in DMSetUp_DA_3D()
594 n11 = rank - 2 * m + 1; in DMSetUp_DA_3D()
595 n14 = rank - m + 1; in DMSetUp_DA_3D()
597 n20 = rank - 2 * m + 1 + (m * n); in DMSetUp_DA_3D()
598 n23 = rank - m + 1 + (m * n); in DMSetUp_DA_3D()
599 n26 = rank + 1 + (m * n); in DMSetUp_DA_3D()
603 n0 = rank + m * (n - 1) - 1 - (m * n); in DMSetUp_DA_3D()
604 n1 = rank + m * (n - 1) - (m * n); in DMSetUp_DA_3D()
605 n2 = rank + m * (n - 1) + 1 - (m * n); in DMSetUp_DA_3D()
606 n9 = rank + m * (n - 1) - 1; in DMSetUp_DA_3D()
607 n10 = rank + m * (n - 1); in DMSetUp_DA_3D()
608 n11 = rank + m * (n - 1) + 1; in DMSetUp_DA_3D()
609 n18 = rank + m * (n - 1) - 1 + (m * n); in DMSetUp_DA_3D()
610 n19 = rank + m * (n - 1) + (m * n); in DMSetUp_DA_3D()
611 n20 = rank + m * (n - 1) + 1 + (m * n); in DMSetUp_DA_3D()
615 n6 = rank - m * (n - 1) - 1 - (m * n); in DMSetUp_DA_3D()
616 n7 = rank - m * (n - 1) - (m * n); in DMSetUp_DA_3D()
617 n8 = rank - m * (n - 1) + 1 - (m * n); in DMSetUp_DA_3D()
618 n15 = rank - m * (n - 1) - 1; in DMSetUp_DA_3D()
619 n16 = rank - m * (n - 1); in DMSetUp_DA_3D()
620 n17 = rank - m * (n - 1) + 1; in DMSetUp_DA_3D()
621 n24 = rank - m * (n - 1) - 1 + (m * n); in DMSetUp_DA_3D()
622 n25 = rank - m * (n - 1) + (m * n); in DMSetUp_DA_3D()
623 n26 = rank - m * (n - 1) + 1 + (m * n); in DMSetUp_DA_3D()
627 n0 = size - (m * n) + rank - m - 1; in DMSetUp_DA_3D()
628 n1 = size - (m * n) + rank - m; in DMSetUp_DA_3D()
629 n2 = size - (m * n) + rank - m + 1; in DMSetUp_DA_3D()
630 n3 = size - (m * n) + rank - 1; in DMSetUp_DA_3D()
631 n4 = size - (m * n) + rank; in DMSetUp_DA_3D()
632 n5 = size - (m * n) + rank + 1; in DMSetUp_DA_3D()
633 n6 = size - (m * n) + rank + m - 1; in DMSetUp_DA_3D()
634 n7 = size - (m * n) + rank + m; in DMSetUp_DA_3D()
635 n8 = size - (m * n) + rank + m + 1; in DMSetUp_DA_3D()
639 n18 = (m * n) - (size - rank) - m - 1; in DMSetUp_DA_3D()
640 n19 = (m * n) - (size - rank) - m; in DMSetUp_DA_3D()
641 n20 = (m * n) - (size - rank) - m + 1; in DMSetUp_DA_3D()
642 n21 = (m * n) - (size - rank) - 1; in DMSetUp_DA_3D()
643 n22 = (m * n) - (size - rank); in DMSetUp_DA_3D()
644 n23 = (m * n) - (size - rank) + 1; in DMSetUp_DA_3D()
645 n24 = (m * n) - (size - rank) + m - 1; in DMSetUp_DA_3D()
646 n25 = (m * n) - (size - rank) + m; in DMSetUp_DA_3D()
647 n26 = (m * n) - (size - rank) + m + 1; in DMSetUp_DA_3D()
651 n0 = size - m * n + rank + m - 1 - m; in DMSetUp_DA_3D()
652 n3 = size - m * n + rank + m - 1; in DMSetUp_DA_3D()
653 n6 = size - m * n + rank + m - 1 + m; in DMSetUp_DA_3D()
657 n18 = m * n - (size - rank) + m - 1 - m; in DMSetUp_DA_3D()
658 n21 = m * n - (size - rank) + m - 1; in DMSetUp_DA_3D()
659 n24 = m * n - (size - rank) + m - 1 + m; in DMSetUp_DA_3D()
663 n0 = rank + m * n - 1 - m * n; in DMSetUp_DA_3D()
664 n9 = rank + m * n - 1; in DMSetUp_DA_3D()
665 n18 = rank + m * n - 1 + m * n; in DMSetUp_DA_3D()
669 n6 = rank - m * (n - 1) + m - 1 - m * n; in DMSetUp_DA_3D()
670 n15 = rank - m * (n - 1) + m - 1; in DMSetUp_DA_3D()
671 n24 = rank - m * (n - 1) + m - 1 + m * n; in DMSetUp_DA_3D()
674 if ((xe == M) && (zs == 0)) { /* Assume an edge, not corner */ in DMSetUp_DA_3D()
675 n2 = size - (m * n - rank) - (m - 1) - m; in DMSetUp_DA_3D()
676 n5 = size - (m * n - rank) - (m - 1); in DMSetUp_DA_3D()
677 n8 = size - (m * n - rank) - (m - 1) + m; in DMSetUp_DA_3D()
680 if ((xe == M) && (ze == P)) { /* Assume an edge, not corner */ in DMSetUp_DA_3D()
681 n20 = m * n - (size - rank) - (m - 1) - m; in DMSetUp_DA_3D()
682 n23 = m * n - (size - rank) - (m - 1); in DMSetUp_DA_3D()
683 n26 = m * n - (size - rank) - (m - 1) + m; in DMSetUp_DA_3D()
686 if ((xe == M) && (ys == 0)) { /* Assume an edge, not corner */ in DMSetUp_DA_3D()
687 n2 = rank + m * (n - 1) - (m - 1) - m * n; in DMSetUp_DA_3D()
688 n11 = rank + m * (n - 1) - (m - 1); in DMSetUp_DA_3D()
689 n20 = rank + m * (n - 1) - (m - 1) + m * n; in DMSetUp_DA_3D()
692 if ((xe == M) && (ye == N)) { /* Assume an edge, not corner */ in DMSetUp_DA_3D()
693 n8 = rank - m * n + 1 - m * n; in DMSetUp_DA_3D()
694 n17 = rank - m * n + 1; in DMSetUp_DA_3D()
695 n26 = rank - m * n + 1 + m * n; in DMSetUp_DA_3D()
699 n0 = size - m + rank - 1; in DMSetUp_DA_3D()
700 n1 = size - m + rank; in DMSetUp_DA_3D()
701 n2 = size - m + rank + 1; in DMSetUp_DA_3D()
705 n18 = m * n - (size - rank) + m * (n - 1) - 1; in DMSetUp_DA_3D()
706 n19 = m * n - (size - rank) + m * (n - 1); in DMSetUp_DA_3D()
707 n20 = m * n - (size - rank) + m * (n - 1) + 1; in DMSetUp_DA_3D()
711 n6 = size - (m * n - rank) - m * (n - 1) - 1; in DMSetUp_DA_3D()
712 n7 = size - (m * n - rank) - m * (n - 1); in DMSetUp_DA_3D()
713 n8 = size - (m * n - rank) - m * (n - 1) + 1; in DMSetUp_DA_3D()
717 n24 = rank - (size - m) - 1; in DMSetUp_DA_3D()
718 n25 = rank - (size - m); in DMSetUp_DA_3D()
719 n26 = rank - (size - m) + 1; in DMSetUp_DA_3D()
723 if ((xs == 0) && (ys == 0) && (zs == 0)) n0 = size - 1; in DMSetUp_DA_3D()
724 if ((xs == 0) && (ys == 0) && (ze == P)) n18 = m * n - 1; in DMSetUp_DA_3D()
725 if ((xs == 0) && (ye == N) && (zs == 0)) n6 = (size - 1) - m * (n - 1); in DMSetUp_DA_3D()
726 if ((xs == 0) && (ye == N) && (ze == P)) n24 = m - 1; in DMSetUp_DA_3D()
727 if ((xe == M) && (ys == 0) && (zs == 0)) n2 = size - m; in DMSetUp_DA_3D()
728 if ((xe == M) && (ys == 0) && (ze == P)) n20 = m * n - m; in DMSetUp_DA_3D()
729 if ((xe == M) && (ye == N) && (zs == 0)) n8 = size - m * n; in DMSetUp_DA_3D()
730 if ((xe == M) && (ye == N) && (ze == P)) n26 = 0; in DMSetUp_DA_3D()
736 if (xs == 0) n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2; in DMSetUp_DA_3D()
737 if (xe == M) n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2; in DMSetUp_DA_3D()
742 if (ys == 0) n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2; in DMSetUp_DA_3D()
743 if (ye == N) n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2; in DMSetUp_DA_3D()
748 if (zs == 0) n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2; in DMSetUp_DA_3D()
749 if (ze == P) n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2; in DMSetUp_DA_3D()
752 PetscCall(PetscMalloc1(27, &dd->neighbors)); in DMSetUp_DA_3D()
754 dd->neighbors[0] = n0; in DMSetUp_DA_3D()
755 dd->neighbors[1] = n1; in DMSetUp_DA_3D()
756 dd->neighbors[2] = n2; in DMSetUp_DA_3D()
757 dd->neighbors[3] = n3; in DMSetUp_DA_3D()
758 dd->neighbors[4] = n4; in DMSetUp_DA_3D()
759 dd->neighbors[5] = n5; in DMSetUp_DA_3D()
760 dd->neighbors[6] = n6; in DMSetUp_DA_3D()
761 dd->neighbors[7] = n7; in DMSetUp_DA_3D()
762 dd->neighbors[8] = n8; in DMSetUp_DA_3D()
763 dd->neighbors[9] = n9; in DMSetUp_DA_3D()
764 dd->neighbors[10] = n10; in DMSetUp_DA_3D()
765 dd->neighbors[11] = n11; in DMSetUp_DA_3D()
766 dd->neighbors[12] = n12; in DMSetUp_DA_3D()
767 dd->neighbors[13] = rank; in DMSetUp_DA_3D()
768 dd->neighbors[14] = n14; in DMSetUp_DA_3D()
769 dd->neighbors[15] = n15; in DMSetUp_DA_3D()
770 dd->neighbors[16] = n16; in DMSetUp_DA_3D()
771 dd->neighbors[17] = n17; in DMSetUp_DA_3D()
772 dd->neighbors[18] = n18; in DMSetUp_DA_3D()
773 dd->neighbors[19] = n19; in DMSetUp_DA_3D()
774 dd->neighbors[20] = n20; in DMSetUp_DA_3D()
775 dd->neighbors[21] = n21; in DMSetUp_DA_3D()
776 dd->neighbors[22] = n22; in DMSetUp_DA_3D()
777 dd->neighbors[23] = n23; in DMSetUp_DA_3D()
778 dd->neighbors[24] = n24; in DMSetUp_DA_3D()
779 dd->neighbors[25] = n25; in DMSetUp_DA_3D()
780 dd->neighbors[26] = n26; in DMSetUp_DA_3D()
805 …n3 = n5 = n6 = n7 = n8 = n9 = n11 = n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1; in DMSetUp_DA_3D()
808 PetscCall(PetscMalloc1((Xe - Xs) * (Ye - Ys) * (Ze - Zs), &idx)); in DMSetUp_DA_3D()
815 x_t = lx[n0 % m]; in DMSetUp_DA_3D()
816 y_t = ly[(n0 % (m * n)) / m]; in DMSetUp_DA_3D()
817 z_t = lz[n0 / (m * n)]; in DMSetUp_DA_3D()
818 s_t = bases[n0] + x_t * y_t * z_t - (s_y - i) * x_t - s_x - (s_z - k - 1) * x_t * y_t; in DMSetUp_DA_3D()
819 … if (twod && (s_t < 0)) s_t = bases[n0] + x_t * y_t * z_t - (s_y - i) * x_t - s_x; /* 2D case */ in DMSetUp_DA_3D()
824 y_t = ly[(n1 % (m * n)) / m]; in DMSetUp_DA_3D()
825 z_t = lz[n1 / (m * n)]; in DMSetUp_DA_3D()
826 s_t = bases[n1] + x_t * y_t * z_t - (s_y + 1 - i) * x_t - (s_z - k - 1) * x_t * y_t; in DMSetUp_DA_3D()
827 … if (twod && (s_t < 0)) s_t = bases[n1] + x_t * y_t * z_t - (s_y + 1 - i) * x_t; /* 2D case */ in DMSetUp_DA_3D()
831 x_t = lx[n2 % m]; in DMSetUp_DA_3D()
832 y_t = ly[(n2 % (m * n)) / m]; in DMSetUp_DA_3D()
833 z_t = lz[n2 / (m * n)]; in DMSetUp_DA_3D()
834 s_t = bases[n2] + x_t * y_t * z_t - (s_y + 1 - i) * x_t - (s_z - k - 1) * x_t * y_t; in DMSetUp_DA_3D()
835 … if (twod && (s_t < 0)) s_t = bases[n2] + x_t * y_t * z_t - (s_y + 1 - i) * x_t; /* 2D case */ in DMSetUp_DA_3D()
842 x_t = lx[n3 % m]; in DMSetUp_DA_3D()
844 z_t = lz[n3 / (m * n)]; in DMSetUp_DA_3D()
845 s_t = bases[n3] + (i + 1) * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t; in DMSetUp_DA_3D()
846 …if (twod && (s_t < 0)) s_t = bases[n3] + (i + 1) * x_t - s_x + x_t * y_t * z_t - x_t * y_t; /* 2D … in DMSetUp_DA_3D()
853 z_t = lz[n4 / (m * n)]; in DMSetUp_DA_3D()
854 s_t = bases[n4] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; in DMSetUp_DA_3D()
855 … if (twod && (s_t < 0)) s_t = bases[n4] + i * x_t + x_t * y_t * z_t - x_t * y_t; /* 2D case */ in DMSetUp_DA_3D()
858 for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + x * i + (s_z - k - 1) * x * y; in DMSetUp_DA_3D()
862 x_t = lx[n5 % m]; in DMSetUp_DA_3D()
864 z_t = lz[n5 / (m * n)]; in DMSetUp_DA_3D()
865 s_t = bases[n5] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; in DMSetUp_DA_3D()
866 … if (twod && (s_t < 0)) s_t = bases[n5] + i * x_t + x_t * y_t * z_t - x_t * y_t; /* 2D case */ in DMSetUp_DA_3D()
873 x_t = lx[n6 % m]; in DMSetUp_DA_3D()
874 y_t = ly[(n6 % (m * n)) / m]; in DMSetUp_DA_3D()
875 z_t = lz[n6 / (m * n)]; in DMSetUp_DA_3D()
876 s_t = bases[n6] + i * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t; in DMSetUp_DA_3D()
877 …if (twod && (s_t < 0)) s_t = bases[n6] + i * x_t - s_x + x_t * y_t * z_t - x_t * y_t; /* 2D case */ in DMSetUp_DA_3D()
882 y_t = ly[(n7 % (m * n)) / m]; in DMSetUp_DA_3D()
883 z_t = lz[n7 / (m * n)]; in DMSetUp_DA_3D()
884 s_t = bases[n7] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; in DMSetUp_DA_3D()
885 …if (twod && (s_t < 0)) s_t = bases[n7] + (i - 1) * x_t + x_t * y_t * z_t - x_t * y_t; /* 2D case */ in DMSetUp_DA_3D()
889 x_t = lx[n8 % m]; in DMSetUp_DA_3D()
890 y_t = ly[(n8 % (m * n)) / m]; in DMSetUp_DA_3D()
891 z_t = lz[n8 / (m * n)]; in DMSetUp_DA_3D()
892 s_t = bases[n8] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; in DMSetUp_DA_3D()
893 …if (twod && (s_t < 0)) s_t = bases[n8] + (i - 1) * x_t + x_t * y_t * z_t - x_t * y_t; /* 2D case */ in DMSetUp_DA_3D()
903 x_t = lx[n9 % m]; in DMSetUp_DA_3D()
904 y_t = ly[(n9 % (m * n)) / m]; in DMSetUp_DA_3D()
906 s_t = bases[n9] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t; in DMSetUp_DA_3D()
911 y_t = ly[(n10 % (m * n)) / m]; in DMSetUp_DA_3D()
913 s_t = bases[n10] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; in DMSetUp_DA_3D()
916 for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (s_y - i) * x; in DMSetUp_DA_3D()
919 x_t = lx[n11 % m]; in DMSetUp_DA_3D()
920 y_t = ly[(n11 % (m * n)) / m]; in DMSetUp_DA_3D()
922 s_t = bases[n11] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; in DMSetUp_DA_3D()
929 x_t = lx[n12 % m]; in DMSetUp_DA_3D()
932 s_t = bases[n12] + (i + 1) * x_t - s_x + k * x_t * y_t; in DMSetUp_DA_3D()
935 for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k * x * y + i * x; in DMSetUp_DA_3D()
943 x_t = lx[n14 % m]; in DMSetUp_DA_3D()
949 for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k * x * y + i * x; in DMSetUp_DA_3D()
955 x_t = lx[n15 % m]; in DMSetUp_DA_3D()
956 y_t = ly[(n15 % (m * n)) / m]; in DMSetUp_DA_3D()
958 s_t = bases[n15] + i * x_t - s_x + k * x_t * y_t; in DMSetUp_DA_3D()
963 y_t = ly[(n16 % (m * n)) / m]; in DMSetUp_DA_3D()
965 s_t = bases[n16] + (i - 1) * x_t + k * x_t * y_t; in DMSetUp_DA_3D()
968 for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (y - i) * x; in DMSetUp_DA_3D()
971 x_t = lx[n17 % m]; in DMSetUp_DA_3D()
972 y_t = ly[(n17 % (m * n)) / m]; in DMSetUp_DA_3D()
974 s_t = bases[n17] + (i - 1) * x_t + k * x_t * y_t; in DMSetUp_DA_3D()
984 x_t = lx[n18 % m]; in DMSetUp_DA_3D()
985 y_t = ly[(n18 % (m * n)) / m]; in DMSetUp_DA_3D()
986 /* z_t = lz[n18 / (m*n)]; */ in DMSetUp_DA_3D()
987 s_t = bases[n18] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t; in DMSetUp_DA_3D()
988 …if (twod && (s_t >= M * N * P)) s_t = bases[n18] - (s_y - i) * x_t - s_x + x_t * y_t; /* 2d case */ in DMSetUp_DA_3D()
993 y_t = ly[(n19 % (m * n)) / m]; in DMSetUp_DA_3D()
994 /* z_t = lz[n19 / (m*n)]; */ in DMSetUp_DA_3D()
995 s_t = bases[n19] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; in DMSetUp_DA_3D()
996 … if (twod && (s_t >= M * N * P)) s_t = bases[n19] - (s_y + 1 - i) * x_t + x_t * y_t; /* 2d case */ in DMSetUp_DA_3D()
1000 x_t = lx[n20 % m]; in DMSetUp_DA_3D()
1001 y_t = ly[(n20 % (m * n)) / m]; in DMSetUp_DA_3D()
1002 /* z_t = lz[n20 / (m*n)]; */ in DMSetUp_DA_3D()
1003 s_t = bases[n20] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; in DMSetUp_DA_3D()
1004 … if (twod && (s_t >= M * N * P)) s_t = bases[n20] - (s_y + 1 - i) * x_t + x_t * y_t; /* 2d case */ in DMSetUp_DA_3D()
1011 x_t = lx[n21 % m]; in DMSetUp_DA_3D()
1013 /* z_t = lz[n21 / (m*n)]; */ in DMSetUp_DA_3D()
1014 s_t = bases[n21] + (i + 1) * x_t - s_x + k * x_t * y_t; in DMSetUp_DA_3D()
1015 if (twod && (s_t >= M * N * P)) s_t = bases[n21] + (i + 1) * x_t - s_x; /* 2d case */ in DMSetUp_DA_3D()
1022 /* z_t = lz[n22 / (m*n)]; */ in DMSetUp_DA_3D()
1024 if (twod && (s_t >= M * N * P)) s_t = bases[n22] + i * x_t; /* 2d case */ in DMSetUp_DA_3D()
1027 for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + (z - k - 1) * x * y + i * x; in DMSetUp_DA_3D()
1031 x_t = lx[n23 % m]; in DMSetUp_DA_3D()
1033 /* z_t = lz[n23 / (m*n)]; */ in DMSetUp_DA_3D()
1035 if (twod && (s_t >= M * N * P)) s_t = bases[n23] + i * x_t; /* 2d case */ in DMSetUp_DA_3D()
1042 x_t = lx[n24 % m]; in DMSetUp_DA_3D()
1043 y_t = ly[(n24 % (m * n)) / m]; in DMSetUp_DA_3D()
1044 /* z_t = lz[n24 / (m*n)]; */ in DMSetUp_DA_3D()
1045 s_t = bases[n24] + i * x_t - s_x + k * x_t * y_t; in DMSetUp_DA_3D()
1046 if (twod && (s_t >= M * N * P)) s_t = bases[n24] + i * x_t - s_x; /* 2d case */ in DMSetUp_DA_3D()
1051 y_t = ly[(n25 % (m * n)) / m]; in DMSetUp_DA_3D()
1052 /* z_t = lz[n25 / (m*n)]; */ in DMSetUp_DA_3D()
1053 s_t = bases[n25] + (i - 1) * x_t + k * x_t * y_t; in DMSetUp_DA_3D()
1054 if (twod && (s_t >= M * N * P)) s_t = bases[n25] + (i - 1) * x_t; /* 2d case */ in DMSetUp_DA_3D()
1058 x_t = lx[n26 % m]; in DMSetUp_DA_3D()
1059 y_t = ly[(n26 % (m * n)) / m]; in DMSetUp_DA_3D()
1060 /* z_t = lz[n26 / (m*n)]; */ in DMSetUp_DA_3D()
1061 s_t = bases[n26] + (i - 1) * x_t + k * x_t * y_t; in DMSetUp_DA_3D()
1062 if (twod && (s_t >= M * N * P)) s_t = bases[n26] + (i - 1) * x_t; /* 2d case */ in DMSetUp_DA_3D()
1106 x_t = lx[n0 % m]; in DMSetUp_DA_3D()
1107 y_t = ly[(n0 % (m * n)) / m]; in DMSetUp_DA_3D()
1108 z_t = lz[n0 / (m * n)]; in DMSetUp_DA_3D()
1109 s_t = bases[n0] + x_t * y_t * z_t - (s_y - i) * x_t - s_x - (s_z - k - 1) * x_t * y_t; in DMSetUp_DA_3D()
1111 } else if (Xs - xs < 0 && Ys - ys < 0 && Zs - zs < 0) { in DMSetUp_DA_3D()
1112 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1116 y_t = ly[(n1 % (m * n)) / m]; in DMSetUp_DA_3D()
1117 z_t = lz[n1 / (m * n)]; in DMSetUp_DA_3D()
1118 s_t = bases[n1] + x_t * y_t * z_t - (s_y + 1 - i) * x_t - (s_z - k - 1) * x_t * y_t; in DMSetUp_DA_3D()
1120 } else if (Ys - ys < 0 && Zs - zs < 0) { in DMSetUp_DA_3D()
1121 for (j = 0; j < x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1124 x_t = lx[n2 % m]; in DMSetUp_DA_3D()
1125 y_t = ly[(n2 % (m * n)) / m]; in DMSetUp_DA_3D()
1126 z_t = lz[n2 / (m * n)]; in DMSetUp_DA_3D()
1127 s_t = bases[n2] + x_t * y_t * z_t - (s_y + 1 - i) * x_t - (s_z - k - 1) * x_t * y_t; in DMSetUp_DA_3D()
1129 } else if (xe - Xe < 0 && Ys - ys < 0 && Zs - zs < 0) { in DMSetUp_DA_3D()
1130 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1136 x_t = lx[n3 % m]; in DMSetUp_DA_3D()
1138 z_t = lz[n3 / (m * n)]; in DMSetUp_DA_3D()
1139 s_t = bases[n3] + (i + 1) * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t; in DMSetUp_DA_3D()
1141 } else if (Xs - xs < 0 && Zs - zs < 0) { in DMSetUp_DA_3D()
1142 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1148 z_t = lz[n4 / (m * n)]; in DMSetUp_DA_3D()
1149 s_t = bases[n4] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; in DMSetUp_DA_3D()
1151 } else if (Zs - zs < 0) { in DMSetUp_DA_3D()
1153 for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + x * i + (s_z - k) * x * y; in DMSetUp_DA_3D()
1155 for (j = 0; j < x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1160 x_t = lx[n5 % m]; in DMSetUp_DA_3D()
1162 z_t = lz[n5 / (m * n)]; in DMSetUp_DA_3D()
1163 s_t = bases[n5] + i * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; in DMSetUp_DA_3D()
1165 } else if (xe - Xe < 0 && Zs - zs < 0) { in DMSetUp_DA_3D()
1166 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1172 x_t = lx[n6 % m]; in DMSetUp_DA_3D()
1173 y_t = ly[(n6 % (m * n)) / m]; in DMSetUp_DA_3D()
1174 z_t = lz[n6 / (m * n)]; in DMSetUp_DA_3D()
1175 s_t = bases[n6] + i * x_t - s_x + x_t * y_t * z_t - (s_z - k) * x_t * y_t; in DMSetUp_DA_3D()
1177 } else if (Xs - xs < 0 && ye - Ye < 0 && Zs - zs < 0) { in DMSetUp_DA_3D()
1178 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1182 y_t = ly[(n7 % (m * n)) / m]; in DMSetUp_DA_3D()
1183 z_t = lz[n7 / (m * n)]; in DMSetUp_DA_3D()
1184 s_t = bases[n7] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; in DMSetUp_DA_3D()
1186 } else if (ye - Ye < 0 && Zs - zs < 0) { in DMSetUp_DA_3D()
1187 for (j = 0; j < x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1190 x_t = lx[n8 % m]; in DMSetUp_DA_3D()
1191 y_t = ly[(n8 % (m * n)) / m]; in DMSetUp_DA_3D()
1192 z_t = lz[n8 / (m * n)]; in DMSetUp_DA_3D()
1193 s_t = bases[n8] + (i - 1) * x_t + x_t * y_t * z_t - (s_z - k) * x_t * y_t; in DMSetUp_DA_3D()
1195 } else if (xe - Xe < 0 && ye - Ye < 0 && Zs - zs < 0) { in DMSetUp_DA_3D()
1196 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1205 x_t = lx[n9 % m]; in DMSetUp_DA_3D()
1206 y_t = ly[(n9 % (m * n)) / m]; in DMSetUp_DA_3D()
1208 s_t = bases[n9] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t; in DMSetUp_DA_3D()
1210 } else if (Xs - xs < 0 && Ys - ys < 0) { in DMSetUp_DA_3D()
1211 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1215 y_t = ly[(n10 % (m * n)) / m]; in DMSetUp_DA_3D()
1217 s_t = bases[n10] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; in DMSetUp_DA_3D()
1219 } else if (Ys - ys < 0) { in DMSetUp_DA_3D()
1221 for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (s_y - i + 1) * x; in DMSetUp_DA_3D()
1223 for (j = 0; j < x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1227 x_t = lx[n11 % m]; in DMSetUp_DA_3D()
1228 y_t = ly[(n11 % (m * n)) / m]; in DMSetUp_DA_3D()
1230 s_t = bases[n11] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; in DMSetUp_DA_3D()
1232 } else if (xe - Xe < 0 && Ys - ys < 0) { in DMSetUp_DA_3D()
1233 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1239 x_t = lx[n12 % m]; in DMSetUp_DA_3D()
1242 s_t = bases[n12] + (i + 1) * x_t - s_x + k * x_t * y_t; in DMSetUp_DA_3D()
1244 } else if (Xs - xs < 0) { in DMSetUp_DA_3D()
1246 for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + s_x - j - 1 + k * x * y + i * x + 1; in DMSetUp_DA_3D()
1248 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1257 x_t = lx[n14 % m]; in DMSetUp_DA_3D()
1262 } else if (xe - Xe < 0) { in DMSetUp_DA_3D()
1264 for (j = 0; j < s_x; j++) idx[nn++] = bases[rank] + x - j - 1 + k * x * y + i * x - 1; in DMSetUp_DA_3D()
1266 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1273 x_t = lx[n15 % m]; in DMSetUp_DA_3D()
1274 y_t = ly[(n15 % (m * n)) / m]; in DMSetUp_DA_3D()
1276 s_t = bases[n15] + i * x_t - s_x + k * x_t * y_t; in DMSetUp_DA_3D()
1278 } else if (Xs - xs < 0 && ye - Ye < 0) { in DMSetUp_DA_3D()
1279 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1283 y_t = ly[(n16 % (m * n)) / m]; in DMSetUp_DA_3D()
1285 s_t = bases[n16] + (i - 1) * x_t + k * x_t * y_t; in DMSetUp_DA_3D()
1287 } else if (ye - Ye < 0) { in DMSetUp_DA_3D()
1289 for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + k * x * y + (y - i - 1) * x; in DMSetUp_DA_3D()
1291 for (j = 0; j < x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1295 x_t = lx[n17 % m]; in DMSetUp_DA_3D()
1296 y_t = ly[(n17 % (m * n)) / m]; in DMSetUp_DA_3D()
1298 s_t = bases[n17] + (i - 1) * x_t + k * x_t * y_t; in DMSetUp_DA_3D()
1300 } else if (xe - Xe < 0 && ye - Ye < 0) { in DMSetUp_DA_3D()
1301 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1310 x_t = lx[n18 % m]; in DMSetUp_DA_3D()
1311 y_t = ly[(n18 % (m * n)) / m]; in DMSetUp_DA_3D()
1312 /* z_t = lz[n18 / (m*n)]; */ in DMSetUp_DA_3D()
1313 s_t = bases[n18] - (s_y - i) * x_t - s_x + (k + 1) * x_t * y_t; in DMSetUp_DA_3D()
1315 } else if (Xs - xs < 0 && Ys - ys < 0 && ze - Ze < 0) { in DMSetUp_DA_3D()
1316 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1320 y_t = ly[(n19 % (m * n)) / m]; in DMSetUp_DA_3D()
1321 /* z_t = lz[n19 / (m*n)]; */ in DMSetUp_DA_3D()
1322 s_t = bases[n19] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; in DMSetUp_DA_3D()
1324 } else if (Ys - ys < 0 && ze - Ze < 0) { in DMSetUp_DA_3D()
1325 for (j = 0; j < x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1328 x_t = lx[n20 % m]; in DMSetUp_DA_3D()
1329 y_t = ly[(n20 % (m * n)) / m]; in DMSetUp_DA_3D()
1330 /* z_t = lz[n20 / (m*n)]; */ in DMSetUp_DA_3D()
1331 s_t = bases[n20] - (s_y + 1 - i) * x_t + (k + 1) * x_t * y_t; in DMSetUp_DA_3D()
1333 } else if (xe - Xe < 0 && Ys - ys < 0 && ze - Ze < 0) { in DMSetUp_DA_3D()
1334 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1340 x_t = lx[n21 % m]; in DMSetUp_DA_3D()
1342 /* z_t = lz[n21 / (m*n)]; */ in DMSetUp_DA_3D()
1343 s_t = bases[n21] + (i + 1) * x_t - s_x + k * x_t * y_t; in DMSetUp_DA_3D()
1345 } else if (Xs - xs < 0 && ze - Ze < 0) { in DMSetUp_DA_3D()
1346 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1352 /* z_t = lz[n22 / (m*n)]; */ in DMSetUp_DA_3D()
1355 } else if (ze - Ze < 0) { in DMSetUp_DA_3D()
1357 for (j = 0; j < x; j++) idx[nn++] = bases[rank] + j + (z - k - 2) * x * y + i * x; in DMSetUp_DA_3D()
1359 for (j = 0; j < x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1364 x_t = lx[n23 % m]; in DMSetUp_DA_3D()
1366 /* z_t = lz[n23 / (m*n)]; */ in DMSetUp_DA_3D()
1369 } else if (xe - Xe < 0 && ze - Ze < 0) { in DMSetUp_DA_3D()
1370 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1376 x_t = lx[n24 % m]; in DMSetUp_DA_3D()
1377 y_t = ly[(n24 % (m * n)) / m]; in DMSetUp_DA_3D()
1378 /* z_t = lz[n24 / (m*n)]; */ in DMSetUp_DA_3D()
1379 s_t = bases[n24] + i * x_t - s_x + k * x_t * y_t; in DMSetUp_DA_3D()
1381 } else if (Xs - xs < 0 && ye - Ye < 0 && ze - Ze < 0) { in DMSetUp_DA_3D()
1382 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1386 y_t = ly[(n25 % (m * n)) / m]; in DMSetUp_DA_3D()
1387 /* z_t = lz[n25 / (m*n)]; */ in DMSetUp_DA_3D()
1388 s_t = bases[n25] + (i - 1) * x_t + k * x_t * y_t; in DMSetUp_DA_3D()
1390 } else if (ye - Ye < 0 && ze - Ze < 0) { in DMSetUp_DA_3D()
1391 for (j = 0; j < x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1394 x_t = lx[n26 % m]; in DMSetUp_DA_3D()
1395 y_t = ly[(n26 % (m * n)) / m]; in DMSetUp_DA_3D()
1396 /* z_t = lz[n26 / (m*n)]; */ in DMSetUp_DA_3D()
1397 s_t = bases[n26] + (i - 1) * x_t + k * x_t * y_t; in DMSetUp_DA_3D()
1399 } else if (xe - Xe < 0 && ye - Ye < 0 && ze - Ze < 0) { in DMSetUp_DA_3D()
1400 for (j = 0; j < s_x; j++) idx[nn++] = -1; in DMSetUp_DA_3D()
1409 PetscCall(ISLocalToGlobalMappingCreate(comm, dof, nn, idx, PETSC_OWN_POINTER, &da->ltogmap)); in DMSetUp_DA_3D()
1412 dd->m = m; in DMSetUp_DA_3D()
1413 dd->n = n; in DMSetUp_DA_3D()
1414 dd->p = p; in DMSetUp_DA_3D()
1416 dd->xs = xs * dof; in DMSetUp_DA_3D()
1417 dd->xe = xe * dof; in DMSetUp_DA_3D()
1418 dd->ys = ys; in DMSetUp_DA_3D()
1419 dd->ye = ye; in DMSetUp_DA_3D()
1420 dd->zs = zs; in DMSetUp_DA_3D()
1421 dd->ze = ze; in DMSetUp_DA_3D()
1422 dd->Xs = Xs * dof; in DMSetUp_DA_3D()
1423 dd->Xe = Xe * dof; in DMSetUp_DA_3D()
1424 dd->Ys = Ys; in DMSetUp_DA_3D()
1425 dd->Ye = Ye; in DMSetUp_DA_3D()
1426 dd->Zs = Zs; in DMSetUp_DA_3D()
1427 dd->Ze = Ze; in DMSetUp_DA_3D()
1432 dd->gtol = gtol; in DMSetUp_DA_3D()
1433 dd->base = base; in DMSetUp_DA_3D()
1434 da->ops->view = DMView_DA_3d; in DMSetUp_DA_3D()
1435 dd->ltol = NULL; in DMSetUp_DA_3D()
1436 dd->ao = NULL; in DMSetUp_DA_3D()
1441 DMDACreate3d - Creates an object that will manage the communication of three-dimensional
1447 + comm - MPI communicator
1448 . bx - type of x ghost nodes the array have.
1450 . by - type of y ghost nodes the array have.
1452 . bz - type of z ghost nodes the array have.
1454 . stencil_type - Type of stencil (`DMDA_STENCIL_STAR` or `DMDA_STENCIL_BOX`)
1455 . M - global dimension in x direction of the array
1456 . N - global dimension in y direction of the array
1457 . P - global dimension in z direction of the array
1458 . m - corresponding number of processors in x dimension (or `PETSC_DECIDE` to have calcu…
1459 . n - corresponding number of processors in y dimension (or `PETSC_DECIDE` to have calcu…
1460 . p - corresponding number of processors in z dimension (or `PETSC_DECIDE` to have calcu…
1461 . dof - number of degrees of freedom per node
1462 . s - stencil width
1463 . lx - arrays containing the number of nodes in each cell along the x coordinates, or `N…
1464 . ly - arrays containing the number of nodes in each cell along the y coordinates, or `NU…
1465 - lz - arrays containing the number of nodes in each cell along the z coordinates, or `NU…
1468 . da - the resulting distributed array object
1471 + -dm_view - Calls `DMView()` at the conclusion of `DMDACreate3d()`
1472 . -da_grid_x <nx> - number of grid points in x direction
1473 . -da_grid_y <ny> - number of grid points in y direction
1474 . -da_grid_z <nz> - number of grid points in z direction
1475 . -da_processors_x <MX> - number of processors in x direction
1476 . -da_processors_y <MY> - number of processors in y direction
1477 . -da_processors_z <MZ> - number of processors in z direction
1478 . -da_bd_x <bx> - boundary type in x direction
1479 . -da_bd_y <by> - boundary type in y direction
1480 . -da_bd_z <bz> - boundary type in x direction
1481 . -da_bd_all <bt> - boundary type in all directions
1482 . -da_refine_x <rx> - refinement ratio in x direction
1483 . -da_refine_y <ry> - refinement ratio in y direction
1484 . -da_refine_z <rz> - refinement ratio in z directio
1485 - -da_refine <n> - refine the `DMDA` n times before creating it
1490 If `lx`, `ly`, or `lz` are non-null, these must be of length as `m`, `n`, `p` and the
1491 corresponding `m`, `n`, or `p` cannot be `PETSC_DECIDE`. Sum of the `lx` entries must be `M`,
1495 standard 7-pt stencil, while `DMDA_STENCIL_BOX` with width 1 denotes
1496 the standard 27-pt stencil.
1512 …DMBoundaryType bz, DMDAStencilType stencil_type, PetscInt M, PetscInt N, PetscInt P, PetscInt m, P… in DMDACreate3d() argument
1517 PetscCall(DMDASetSizes(*da, M, N, P)); in DMDACreate3d()
1518 PetscCall(DMDASetNumProcs(*da, m, n, p)); in DMDACreate3d()