Lines Matching +full:- +full:m
13 DM_DA *dd = (DM_DA *)da->data; in DMView_DA_1d()
48 …PetscCall(PetscViewerASCIIPrintf(viewer, " Load Balance - Grid Points: Min %" PetscInt_FMT " avg… in DMView_DA_1d()
55 …r, "Processor [%d] M %" PetscInt_FMT " m %" PetscInt_FMT " w %" PetscInt_FMT " s %" PetscInt_FMT "… in DMView_DA_1d()
62 double ymin = -1, ymax = 1, xmin = -1, xmax = dd->M, x; in DMView_DA_1d()
81 …for (xmin_tmp = 0; xmin_tmp < dd->M; xmin_tmp++) PetscCall(PetscDrawLine(draw, (double)xmin_tmp, y… in DMView_DA_1d()
83 xmax = dd->M - 1; in DMView_DA_1d()
95 xmin = dd->xs / dd->w; in DMView_DA_1d()
96 xmax = (dd->xe / dd->w) - 1; in DMView_DA_1d()
102 base = dd->base / dd->w; in DMView_DA_1d()
125 DM_DA *dd = (DM_DA *)da->data; in DMSetUp_DA_1D()
126 const PetscInt M = dd->M; in DMSetUp_DA_1D() local
127 const PetscInt dof = dd->w; in DMSetUp_DA_1D()
128 const PetscInt s = dd->s; in DMSetUp_DA_1D()
130 const PetscInt *lx = dd->lx; in DMSetUp_DA_1D()
131 DMBoundaryType bx = dd->bx; in DMSetUp_DA_1D()
138 PetscInt i, *idx, nn, left, xs, xe, x, Xs, Xe, start, m, IXs, IXe; in DMSetUp_DA_1D() local
145 dd->p = 1; in DMSetUp_DA_1D()
146 dd->n = 1; in DMSetUp_DA_1D()
147 dd->m = size; in DMSetUp_DA_1D()
148 m = dd->m; in DMSetUp_DA_1D()
152 …PetscCheck(M >= m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "More processes than data points! %"… in DMSetUp_DA_1D()
153 …heck((M - 1) >= s || size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Array is too small for… in DMSetUp_DA_1D()
161 PetscCall(PetscMalloc1(m, &dd->lx)); in DMSetUp_DA_1D()
162 …PetscCall(PetscOptionsGetBool(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_partitio… in DMSetUp_DA_1D()
163 …PetscCall(PetscOptionsGetBool(((PetscObject)da)->options, ((PetscObject)da)->prefix, "-da_partitio… in DMSetUp_DA_1D()
165 xs = rank * M / m; in DMSetUp_DA_1D()
166 x = (rank + 1) * M / m - xs; in DMSetUp_DA_1D()
168 x = (M + rank) / m; in DMSetUp_DA_1D()
169 if (M / m == x) xs = rank * x; in DMSetUp_DA_1D()
170 else xs = rank * (x - 1) + (M + rank) % (x * m); in DMSetUp_DA_1D()
173 x = M / m + ((M % m) > rank); in DMSetUp_DA_1D()
174 if (rank >= (M % m)) xs = (rank * (M / m) + M % m); in DMSetUp_DA_1D()
175 else xs = rank * (M / m) + rank; in DMSetUp_DA_1D()
177 PetscCallMPI(MPI_Allgather(&xs, 1, MPIU_INT, dd->lx, 1, MPIU_INT, comm)); in DMSetUp_DA_1D()
178 for (i = 0; i < m - 1; i++) dd->lx[i] = dd->lx[i + 1] - dd->lx[i]; in DMSetUp_DA_1D()
179 dd->lx[m - 1] = M - dd->lx[m - 1]; in DMSetUp_DA_1D()
187 …k(left == M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Sum of lx across processors not equal to … in DMSetUp_DA_1D()
194 …PetscCheck((x >= s) || ((M <= 1) && (bx != DM_BOUNDARY_PERIODIC)), PETSC_COMM_SELF, PETSC_ERR_ARG_… in DMSetUp_DA_1D()
199 if (xs - sDist > 0) { in DMSetUp_DA_1D()
200 Xs = xs - sDist; in DMSetUp_DA_1D()
201 IXs = xs - sDist; in DMSetUp_DA_1D()
203 if (bx) Xs = xs - sDist; in DMSetUp_DA_1D()
207 if (xe + sDist <= M) { in DMSetUp_DA_1D()
212 else Xe = M; in DMSetUp_DA_1D()
213 IXe = M; in DMSetUp_DA_1D()
217 Xs = xs - sDist; in DMSetUp_DA_1D()
219 IXs = xs - sDist; in DMSetUp_DA_1D()
224 dd->Nlocal = dof * x; in DMSetUp_DA_1D()
225 PetscCall(VecCreateMPIWithArray(comm, dof, dd->Nlocal, PETSC_DECIDE, NULL, &global)); in DMSetUp_DA_1D()
226 dd->nlocal = dof * (Xe - Xs); in DMSetUp_DA_1D()
227 PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, dof, dd->nlocal, NULL, &local)); in DMSetUp_DA_1D()
233 PetscCall(ISCreateStride(comm, dof * (IXe - IXs), dof * (IXs - Xs), 1, &to)); in DMSetUp_DA_1D()
237 for (i = 0; i < IXs - Xs; i++) idx[i] = -1; /* prepend with -1s if needed for ghosted case*/ in DMSetUp_DA_1D()
239 nn = IXs - Xs; in DMSetUp_DA_1D()
242 if ((xs - sDist + i) >= 0) idx[nn++] = xs - sDist + i; in DMSetUp_DA_1D()
243 else idx[nn++] = M + (xs - sDist + i); in DMSetUp_DA_1D()
246 for (i = 0; i < x; i++) idx[nn++] = xs + i; /* Non-ghost points */ in DMSetUp_DA_1D()
249 if ((xe + i) < M) idx[nn++] = xe + i; in DMSetUp_DA_1D()
250 else idx[nn++] = (xe + i) - M; in DMSetUp_DA_1D()
254 if ((xs - sDist + i) >= 0) idx[nn++] = xs - sDist + i; in DMSetUp_DA_1D()
255 else idx[nn++] = sDist - i; in DMSetUp_DA_1D()
258 for (i = 0; i < x; i++) idx[nn++] = xs + i; /* Non-ghost points */ in DMSetUp_DA_1D()
261 if ((xe + i) < M) idx[nn++] = xe + i; in DMSetUp_DA_1D()
262 else idx[nn++] = M - (i + 2); in DMSetUp_DA_1D()
265 if (0 <= xs - sDist) { in DMSetUp_DA_1D()
266 for (i = 0; i < sDist; i++) idx[nn++] = xs - sDist + i; in DMSetUp_DA_1D()
273 if ((xe + sDist) <= M) { in DMSetUp_DA_1D()
276 for (i = xe; i < M; i++) idx[nn++] = i; in DMSetUp_DA_1D()
280 PetscCall(ISCreateBlock(comm, dof, nn - IXs + Xs, &idx[IXs - Xs], PETSC_USE_POINTER, &from)); in DMSetUp_DA_1D()
287 dd->xs = dof * xs; in DMSetUp_DA_1D()
288 dd->xe = dof * xe; in DMSetUp_DA_1D()
289 dd->ys = 0; in DMSetUp_DA_1D()
290 dd->ye = 1; in DMSetUp_DA_1D()
291 dd->zs = 0; in DMSetUp_DA_1D()
292 dd->ze = 1; in DMSetUp_DA_1D()
293 dd->Xs = dof * Xs; in DMSetUp_DA_1D()
294 dd->Xe = dof * Xe; in DMSetUp_DA_1D()
295 dd->Ys = 0; in DMSetUp_DA_1D()
296 dd->Ye = 1; in DMSetUp_DA_1D()
297 dd->Zs = 0; in DMSetUp_DA_1D()
298 dd->Ze = 1; in DMSetUp_DA_1D()
300 dd->gtol = gtol; in DMSetUp_DA_1D()
301 dd->base = dof * xs; in DMSetUp_DA_1D()
302 da->ops->view = DMView_DA_1d; in DMSetUp_DA_1D()
308 for (i = 0; i < Xe - IXe; i++) idx[nn++] = -1; /* pad with -1s if needed for ghosted case*/ in DMSetUp_DA_1D()
310 PetscCall(ISLocalToGlobalMappingCreate(comm, dof, nn, idx, PETSC_OWN_POINTER, &da->ltogmap)); in DMSetUp_DA_1D()
315 DMDACreate1d - Creates an object that will manage the communication of one-dimensional
321 + comm - MPI communicator
322 . bx - type of ghost cells at the boundary the array should have, if any. Use
324 . M - global dimension of the array (that is the number of grid points)
325 . dof - number of degrees of freedom per node
326 . s - stencil width
327 - lx - array containing number of nodes in the X direction on each processor,
328 or `NULL`. If non-null, must be of length as the number of processes in the MPI_Comm.
329 The sum of these entries must equal `M`
332 . da - the resulting distributed array object
335 + -dm_view - Calls `DMView()` at the conclusion of `DMDACreate1d()`
336 . -da_grid_x <nx> - number of grid points in x direction
337 . -da_refine_x <rx> - refinement factor
338 - -da_refine <n> - refine the `DMDA` n times before creating it
357 PetscErrorCode DMDACreate1d(MPI_Comm comm, DMBoundaryType bx, PetscInt M, PetscInt dof, PetscInt s,… in DMDACreate1d() argument
364 PetscCall(DMDASetSizes(*da, M, 1, 1)); in DMDACreate1d()