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