Actual source code: plex.c
1: #include <petsc/private/dmpleximpl.h>
2: #include <petsc/private/isimpl.h>
3: #include <petsc/private/vecimpl.h>
4: #include <petsc/private/glvisvecimpl.h>
5: #include <petscsf.h>
6: #include <petscds.h>
7: #include <petscdraw.h>
8: #include <petscdmfield.h>
9: #include <petscdmplextransform.h>
11: /* Logging support */
12: PetscLogEvent DMPLEX_Interpolate, DMPLEX_Partition, DMPLEX_Distribute, DMPLEX_DistributeCones, DMPLEX_DistributeLabels, DMPLEX_DistributeSF, DMPLEX_DistributeOverlap, DMPLEX_DistributeField, DMPLEX_DistributeData, DMPLEX_Migrate, DMPLEX_InterpolateSF, DMPLEX_GlobalToNaturalBegin, DMPLEX_GlobalToNaturalEnd, DMPLEX_NaturalToGlobalBegin, DMPLEX_NaturalToGlobalEnd, DMPLEX_Stratify, DMPLEX_Symmetrize, DMPLEX_Preallocate, DMPLEX_ResidualFEM, DMPLEX_JacobianFEM, DMPLEX_InterpolatorFEM, DMPLEX_InjectorFEM, DMPLEX_IntegralFEM, DMPLEX_CreateGmsh, DMPLEX_RebalanceSharedPoints, DMPLEX_PartSelf, DMPLEX_PartLabelInvert, DMPLEX_PartLabelCreateSF, DMPLEX_PartStratSF, DMPLEX_CreatePointSF,DMPLEX_LocatePoints;
14: PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
16: /*@
17: DMPlexIsSimplex - Is the first cell in this mesh a simplex?
19: Input Parameter:
20: . dm - The DMPlex object
22: Output Parameter:
23: . simplex - Flag checking for a simplex
25: Note: This just gives the first range of cells found. If the mesh has several cell types, it will only give the first.
26: If the mesh has no cells, this returns PETSC_FALSE.
28: Level: intermediate
30: .seealso DMPlexGetSimplexOrBoxCells(), DMPlexGetCellType(), DMPlexGetHeightStratum(), DMPolytopeTypeGetNumVertices()
31: @*/
32: PetscErrorCode DMPlexIsSimplex(DM dm, PetscBool *simplex)
33: {
34: DMPolytopeType ct;
35: PetscInt cStart, cEnd;
39: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
40: if (cEnd <= cStart) {*simplex = PETSC_FALSE; return(0);}
41: DMPlexGetCellType(dm, cStart, &ct);
42: *simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
43: return(0);
44: }
46: /*@
47: DMPlexGetSimplexOrBoxCells - Get the range of cells which are neither prisms nor ghost FV cells
49: Input Parameters:
50: + dm - The DMPlex object
51: - height - The cell height in the Plex, 0 is the default
53: Output Parameters:
54: + cStart - The first "normal" cell
55: - cEnd - The upper bound on "normal"" cells
57: Note: This just gives the first range of cells found. If the mesh has several cell types, it will only give the first.
59: Level: developer
61: .seealso DMPlexConstructGhostCells(), DMPlexGetGhostCellStratum()
62: @*/
63: PetscErrorCode DMPlexGetSimplexOrBoxCells(DM dm, PetscInt height, PetscInt *cStart, PetscInt *cEnd)
64: {
65: DMPolytopeType ct = DM_POLYTOPE_UNKNOWN;
66: PetscInt cS, cE, c;
70: DMPlexGetHeightStratum(dm, PetscMax(height, 0), &cS, &cE);
71: for (c = cS; c < cE; ++c) {
72: DMPolytopeType cct;
74: DMPlexGetCellType(dm, c, &cct);
75: if ((PetscInt) cct < 0) break;
76: switch (cct) {
77: case DM_POLYTOPE_POINT:
78: case DM_POLYTOPE_SEGMENT:
79: case DM_POLYTOPE_TRIANGLE:
80: case DM_POLYTOPE_QUADRILATERAL:
81: case DM_POLYTOPE_TETRAHEDRON:
82: case DM_POLYTOPE_HEXAHEDRON:
83: ct = cct;
84: break;
85: default: break;
86: }
87: if (ct != DM_POLYTOPE_UNKNOWN) break;
88: }
89: if (ct != DM_POLYTOPE_UNKNOWN) {
90: DMLabel ctLabel;
92: DMPlexGetCellTypeLabel(dm, &ctLabel);
93: DMLabelGetStratumBounds(ctLabel, ct, &cS, &cE);
94: }
95: if (cStart) *cStart = cS;
96: if (cEnd) *cEnd = cE;
97: return(0);
98: }
100: PetscErrorCode DMPlexGetFieldType_Internal(DM dm, PetscSection section, PetscInt field, PetscInt *sStart, PetscInt *sEnd, PetscViewerVTKFieldType *ft)
101: {
102: PetscInt cdim, pStart, pEnd, vStart, vEnd, cStart, cEnd;
103: PetscInt vcdof[2] = {0,0}, globalvcdof[2];
107: *ft = PETSC_VTK_INVALID;
108: DMGetCoordinateDim(dm, &cdim);
109: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
110: DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);
111: PetscSectionGetChart(section, &pStart, &pEnd);
112: if (field >= 0) {
113: if ((vStart >= pStart) && (vStart < pEnd)) {PetscSectionGetFieldDof(section, vStart, field, &vcdof[0]);}
114: if ((cStart >= pStart) && (cStart < pEnd)) {PetscSectionGetFieldDof(section, cStart, field, &vcdof[1]);}
115: } else {
116: if ((vStart >= pStart) && (vStart < pEnd)) {PetscSectionGetDof(section, vStart, &vcdof[0]);}
117: if ((cStart >= pStart) && (cStart < pEnd)) {PetscSectionGetDof(section, cStart, &vcdof[1]);}
118: }
119: MPI_Allreduce(vcdof, globalvcdof, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));
120: if (globalvcdof[0]) {
121: *sStart = vStart;
122: *sEnd = vEnd;
123: if (globalvcdof[0] == cdim) *ft = PETSC_VTK_POINT_VECTOR_FIELD;
124: else *ft = PETSC_VTK_POINT_FIELD;
125: } else if (globalvcdof[1]) {
126: *sStart = cStart;
127: *sEnd = cEnd;
128: if (globalvcdof[1] == cdim) *ft = PETSC_VTK_CELL_VECTOR_FIELD;
129: else *ft = PETSC_VTK_CELL_FIELD;
130: } else {
131: if (field >= 0) {
132: const char *fieldname;
134: PetscSectionGetFieldName(section, field, &fieldname);
135: PetscInfo2((PetscObject) dm, "Could not classify VTK output type of section field %D \"%s\"\n", field, fieldname);
136: } else {
137: PetscInfo((PetscObject) dm, "Could not classify VTK output type of section\"%s\"\n");
138: }
139: }
140: return(0);
141: }
143: static PetscErrorCode VecView_Plex_Local_Draw(Vec v, PetscViewer viewer)
144: {
145: DM dm;
146: PetscSection s;
147: PetscDraw draw, popup;
148: DM cdm;
149: PetscSection coordSection;
150: Vec coordinates;
151: const PetscScalar *coords, *array;
152: PetscReal bound[4] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
153: PetscReal vbound[2], time;
154: PetscBool isnull, flg;
155: PetscInt dim, Nf, f, Nc, comp, vStart, vEnd, cStart, cEnd, c, N, level, step, w = 0;
156: const char *name;
157: char title[PETSC_MAX_PATH_LEN];
158: PetscErrorCode ierr;
161: PetscViewerDrawGetDraw(viewer, 0, &draw);
162: PetscDrawIsNull(draw, &isnull);
163: if (isnull) return(0);
165: VecGetDM(v, &dm);
166: DMGetCoordinateDim(dm, &dim);
167: if (dim != 2) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Cannot draw meshes of dimension %D. Use PETSCVIEWERGLVIS", dim);
168: DMGetLocalSection(dm, &s);
169: PetscSectionGetNumFields(s, &Nf);
170: DMGetCoarsenLevel(dm, &level);
171: DMGetCoordinateDM(dm, &cdm);
172: DMGetLocalSection(cdm, &coordSection);
173: DMGetCoordinatesLocal(dm, &coordinates);
174: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
175: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
177: PetscObjectGetName((PetscObject) v, &name);
178: DMGetOutputSequenceNumber(dm, &step, &time);
180: VecGetLocalSize(coordinates, &N);
181: VecGetArrayRead(coordinates, &coords);
182: for (c = 0; c < N; c += dim) {
183: bound[0] = PetscMin(bound[0], PetscRealPart(coords[c])); bound[2] = PetscMax(bound[2], PetscRealPart(coords[c]));
184: bound[1] = PetscMin(bound[1], PetscRealPart(coords[c+1])); bound[3] = PetscMax(bound[3], PetscRealPart(coords[c+1]));
185: }
186: VecRestoreArrayRead(coordinates, &coords);
187: PetscDrawClear(draw);
189: /* Could implement something like DMDASelectFields() */
190: for (f = 0; f < Nf; ++f) {
191: DM fdm = dm;
192: Vec fv = v;
193: IS fis;
194: char prefix[PETSC_MAX_PATH_LEN];
195: const char *fname;
197: PetscSectionGetFieldComponents(s, f, &Nc);
198: PetscSectionGetFieldName(s, f, &fname);
200: if (v->hdr.prefix) {PetscStrncpy(prefix, v->hdr.prefix,sizeof(prefix));}
201: else {prefix[0] = '\0';}
202: if (Nf > 1) {
203: DMCreateSubDM(dm, 1, &f, &fis, &fdm);
204: VecGetSubVector(v, fis, &fv);
205: PetscStrlcat(prefix, fname,sizeof(prefix));
206: PetscStrlcat(prefix, "_",sizeof(prefix));
207: }
208: for (comp = 0; comp < Nc; ++comp, ++w) {
209: PetscInt nmax = 2;
211: PetscViewerDrawGetDraw(viewer, w, &draw);
212: if (Nc > 1) {PetscSNPrintf(title, sizeof(title), "%s:%s_%D Step: %D Time: %.4g", name, fname, comp, step, time);}
213: else {PetscSNPrintf(title, sizeof(title), "%s:%s Step: %D Time: %.4g", name, fname, step, time);}
214: PetscDrawSetTitle(draw, title);
216: /* TODO Get max and min only for this component */
217: PetscOptionsGetRealArray(NULL, prefix, "-vec_view_bounds", vbound, &nmax, &flg);
218: if (!flg) {
219: VecMin(fv, NULL, &vbound[0]);
220: VecMax(fv, NULL, &vbound[1]);
221: if (vbound[1] <= vbound[0]) vbound[1] = vbound[0] + 1.0;
222: }
223: PetscDrawGetPopup(draw, &popup);
224: PetscDrawScalePopup(popup, vbound[0], vbound[1]);
225: PetscDrawSetCoordinates(draw, bound[0], bound[1], bound[2], bound[3]);
227: VecGetArrayRead(fv, &array);
228: for (c = cStart; c < cEnd; ++c) {
229: PetscScalar *coords = NULL, *a = NULL;
230: PetscInt numCoords, color[4] = {-1,-1,-1,-1};
232: DMPlexPointLocalRead(fdm, c, array, &a);
233: if (a) {
234: color[0] = PetscDrawRealToColor(PetscRealPart(a[comp]), vbound[0], vbound[1]);
235: color[1] = color[2] = color[3] = color[0];
236: } else {
237: PetscScalar *vals = NULL;
238: PetscInt numVals, va;
240: DMPlexVecGetClosure(fdm, NULL, fv, c, &numVals, &vals);
241: if (numVals % Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The number of components %D does not divide the number of values in the closure %D", Nc, numVals);
242: switch (numVals/Nc) {
243: case 3: /* P1 Triangle */
244: case 4: /* P1 Quadrangle */
245: for (va = 0; va < numVals/Nc; ++va) color[va] = PetscDrawRealToColor(PetscRealPart(vals[va*Nc+comp]), vbound[0], vbound[1]);
246: break;
247: case 6: /* P2 Triangle */
248: case 8: /* P2 Quadrangle */
249: for (va = 0; va < numVals/(Nc*2); ++va) color[va] = PetscDrawRealToColor(PetscRealPart(vals[va*Nc+comp + numVals/(Nc*2)]), vbound[0], vbound[1]);
250: break;
251: default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of values for cell closure %D cannot be handled", numVals/Nc);
252: }
253: DMPlexVecRestoreClosure(fdm, NULL, fv, c, &numVals, &vals);
254: }
255: DMPlexVecGetClosure(dm, coordSection, coordinates, c, &numCoords, &coords);
256: switch (numCoords) {
257: case 6:
258: case 12: /* Localized triangle */
259: PetscDrawTriangle(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), color[0], color[1], color[2]);
260: break;
261: case 8:
262: case 16: /* Localized quadrilateral */
263: PetscDrawTriangle(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), color[0], color[1], color[2]);
264: PetscDrawTriangle(draw, PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[6]), PetscRealPart(coords[7]), PetscRealPart(coords[0]), PetscRealPart(coords[1]), color[2], color[3], color[0]);
265: break;
266: default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot draw cells with %D coordinates", numCoords);
267: }
268: DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, &numCoords, &coords);
269: }
270: VecRestoreArrayRead(fv, &array);
271: PetscDrawFlush(draw);
272: PetscDrawPause(draw);
273: PetscDrawSave(draw);
274: }
275: if (Nf > 1) {
276: VecRestoreSubVector(v, fis, &fv);
277: ISDestroy(&fis);
278: DMDestroy(&fdm);
279: }
280: }
281: return(0);
282: }
284: static PetscErrorCode VecView_Plex_Local_VTK(Vec v, PetscViewer viewer)
285: {
286: DM dm;
287: Vec locv;
288: const char *name;
289: PetscSection section;
290: PetscInt pStart, pEnd;
291: PetscInt numFields;
292: PetscViewerVTKFieldType ft;
293: PetscErrorCode ierr;
296: VecGetDM(v, &dm);
297: DMCreateLocalVector(dm, &locv); /* VTK viewer requires exclusive ownership of the vector */
298: PetscObjectGetName((PetscObject) v, &name);
299: PetscObjectSetName((PetscObject) locv, name);
300: VecCopy(v, locv);
301: DMGetLocalSection(dm, §ion);
302: PetscSectionGetNumFields(section, &numFields);
303: if (!numFields) {
304: DMPlexGetFieldType_Internal(dm, section, PETSC_DETERMINE, &pStart, &pEnd, &ft);
305: PetscViewerVTKAddField(viewer, (PetscObject) dm, DMPlexVTKWriteAll, PETSC_DEFAULT, ft, PETSC_TRUE,(PetscObject) locv);
306: } else {
307: PetscInt f;
309: for (f = 0; f < numFields; f++) {
310: DMPlexGetFieldType_Internal(dm, section, f, &pStart, &pEnd, &ft);
311: if (ft == PETSC_VTK_INVALID) continue;
312: PetscObjectReference((PetscObject)locv);
313: PetscViewerVTKAddField(viewer, (PetscObject) dm, DMPlexVTKWriteAll, f, ft, PETSC_TRUE,(PetscObject) locv);
314: }
315: VecDestroy(&locv);
316: }
317: return(0);
318: }
320: PetscErrorCode VecView_Plex_Local(Vec v, PetscViewer viewer)
321: {
322: DM dm;
323: PetscBool isvtk, ishdf5, isdraw, isglvis;
327: VecGetDM(v, &dm);
328: if (!dm) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
329: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);
330: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
331: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW, &isdraw);
332: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERGLVIS, &isglvis);
333: if (isvtk || ishdf5 || isdraw || isglvis) {
334: PetscInt i,numFields;
335: PetscObject fe;
336: PetscBool fem = PETSC_FALSE;
337: Vec locv = v;
338: const char *name;
339: PetscInt step;
340: PetscReal time;
342: DMGetNumFields(dm, &numFields);
343: for (i=0; i<numFields; i++) {
344: DMGetField(dm, i, NULL, &fe);
345: if (fe->classid == PETSCFE_CLASSID) { fem = PETSC_TRUE; break; }
346: }
347: if (fem) {
348: PetscObject isZero;
350: DMGetLocalVector(dm, &locv);
351: PetscObjectGetName((PetscObject) v, &name);
352: PetscObjectSetName((PetscObject) locv, name);
353: PetscObjectQuery((PetscObject) v, "__Vec_bc_zero__", &isZero);
354: PetscObjectCompose((PetscObject) locv, "__Vec_bc_zero__", isZero);
355: VecCopy(v, locv);
356: DMGetOutputSequenceNumber(dm, NULL, &time);
357: DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locv, time, NULL, NULL, NULL);
358: }
359: if (isvtk) {
360: VecView_Plex_Local_VTK(locv, viewer);
361: } else if (ishdf5) {
362: #if defined(PETSC_HAVE_HDF5)
363: VecView_Plex_Local_HDF5_Internal(locv, viewer);
364: #else
365: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
366: #endif
367: } else if (isdraw) {
368: VecView_Plex_Local_Draw(locv, viewer);
369: } else if (isglvis) {
370: DMGetOutputSequenceNumber(dm, &step, NULL);
371: PetscViewerGLVisSetSnapId(viewer, step);
372: VecView_GLVis(locv, viewer);
373: }
374: if (fem) {
375: PetscObjectCompose((PetscObject) locv, "__Vec_bc_zero__", NULL);
376: DMRestoreLocalVector(dm, &locv);
377: }
378: } else {
379: PetscBool isseq;
381: PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);
382: if (isseq) {VecView_Seq(v, viewer);}
383: else {VecView_MPI(v, viewer);}
384: }
385: return(0);
386: }
388: PetscErrorCode VecView_Plex(Vec v, PetscViewer viewer)
389: {
390: DM dm;
391: PetscBool isvtk, ishdf5, isdraw, isglvis, isexodusii;
395: VecGetDM(v, &dm);
396: if (!dm) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
397: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);
398: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
399: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW, &isdraw);
400: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERGLVIS, &isglvis);
401: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWEREXODUSII, &isexodusii);
402: if (isvtk || isdraw || isglvis) {
403: Vec locv;
404: PetscObject isZero;
405: const char *name;
407: DMGetLocalVector(dm, &locv);
408: PetscObjectGetName((PetscObject) v, &name);
409: PetscObjectSetName((PetscObject) locv, name);
410: DMGlobalToLocalBegin(dm, v, INSERT_VALUES, locv);
411: DMGlobalToLocalEnd(dm, v, INSERT_VALUES, locv);
412: PetscObjectQuery((PetscObject) v, "__Vec_bc_zero__", &isZero);
413: PetscObjectCompose((PetscObject) locv, "__Vec_bc_zero__", isZero);
414: VecView_Plex_Local(locv, viewer);
415: PetscObjectCompose((PetscObject) locv, "__Vec_bc_zero__", NULL);
416: DMRestoreLocalVector(dm, &locv);
417: } else if (ishdf5) {
418: #if defined(PETSC_HAVE_HDF5)
419: VecView_Plex_HDF5_Internal(v, viewer);
420: #else
421: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
422: #endif
423: } else if (isexodusii) {
424: #if defined(PETSC_HAVE_EXODUSII)
425: VecView_PlexExodusII_Internal(v, viewer);
426: #else
427: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "ExodusII not supported in this build.\nPlease reconfigure using --download-exodusii");
428: #endif
429: } else {
430: PetscBool isseq;
432: PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);
433: if (isseq) {VecView_Seq(v, viewer);}
434: else {VecView_MPI(v, viewer);}
435: }
436: return(0);
437: }
439: PetscErrorCode VecView_Plex_Native(Vec originalv, PetscViewer viewer)
440: {
441: DM dm;
442: MPI_Comm comm;
443: PetscViewerFormat format;
444: Vec v;
445: PetscBool isvtk, ishdf5;
446: PetscErrorCode ierr;
449: VecGetDM(originalv, &dm);
450: PetscObjectGetComm((PetscObject) originalv, &comm);
451: if (!dm) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
452: PetscViewerGetFormat(viewer, &format);
453: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
454: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);
455: if (format == PETSC_VIEWER_NATIVE) {
456: /* Natural ordering is the common case for DMDA, NATIVE means plain vector, for PLEX is the opposite */
457: /* this need a better fix */
458: if (dm->useNatural) {
459: if (dm->sfNatural) {
460: const char *vecname;
461: PetscInt n, nroots;
463: VecGetLocalSize(originalv, &n);
464: PetscSFGetGraph(dm->sfNatural, &nroots, NULL, NULL, NULL);
465: if (n == nroots) {
466: DMGetGlobalVector(dm, &v);
467: DMPlexGlobalToNaturalBegin(dm, originalv, v);
468: DMPlexGlobalToNaturalEnd(dm, originalv, v);
469: PetscObjectGetName((PetscObject) originalv, &vecname);
470: PetscObjectSetName((PetscObject) v, vecname);
471: } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "DM global to natural SF only handles global vectors");
472: } else SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created");
473: } else v = originalv;
474: } else v = originalv;
476: if (ishdf5) {
477: #if defined(PETSC_HAVE_HDF5)
478: VecView_Plex_HDF5_Native_Internal(v, viewer);
479: #else
480: SETERRQ(comm, PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
481: #endif
482: } else if (isvtk) {
483: SETERRQ(comm, PETSC_ERR_SUP, "VTK format does not support viewing in natural order. Please switch to HDF5.");
484: } else {
485: PetscBool isseq;
487: PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);
488: if (isseq) {VecView_Seq(v, viewer);}
489: else {VecView_MPI(v, viewer);}
490: }
491: if (v != originalv) {DMRestoreGlobalVector(dm, &v);}
492: return(0);
493: }
495: PetscErrorCode VecLoad_Plex_Local(Vec v, PetscViewer viewer)
496: {
497: DM dm;
498: PetscBool ishdf5;
502: VecGetDM(v, &dm);
503: if (!dm) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
504: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
505: if (ishdf5) {
506: DM dmBC;
507: Vec gv;
508: const char *name;
510: DMGetOutputDM(dm, &dmBC);
511: DMGetGlobalVector(dmBC, &gv);
512: PetscObjectGetName((PetscObject) v, &name);
513: PetscObjectSetName((PetscObject) gv, name);
514: VecLoad_Default(gv, viewer);
515: DMGlobalToLocalBegin(dmBC, gv, INSERT_VALUES, v);
516: DMGlobalToLocalEnd(dmBC, gv, INSERT_VALUES, v);
517: DMRestoreGlobalVector(dmBC, &gv);
518: } else {
519: VecLoad_Default(v, viewer);
520: }
521: return(0);
522: }
524: PetscErrorCode VecLoad_Plex(Vec v, PetscViewer viewer)
525: {
526: DM dm;
527: PetscBool ishdf5,isexodusii;
531: VecGetDM(v, &dm);
532: if (!dm) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
533: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
534: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWEREXODUSII, &isexodusii);
535: if (ishdf5) {
536: #if defined(PETSC_HAVE_HDF5)
537: VecLoad_Plex_HDF5_Internal(v, viewer);
538: #else
539: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
540: #endif
541: } else if (isexodusii) {
542: #if defined(PETSC_HAVE_EXODUSII)
543: VecLoad_PlexExodusII_Internal(v, viewer);
544: #else
545: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "ExodusII not supported in this build.\nPlease reconfigure using --download-exodusii");
546: #endif
547: } else {
548: VecLoad_Default(v, viewer);
549: }
550: return(0);
551: }
553: PetscErrorCode VecLoad_Plex_Native(Vec originalv, PetscViewer viewer)
554: {
555: DM dm;
556: PetscViewerFormat format;
557: PetscBool ishdf5;
558: PetscErrorCode ierr;
561: VecGetDM(originalv, &dm);
562: if (!dm) SETERRQ(PetscObjectComm((PetscObject) originalv), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
563: PetscViewerGetFormat(viewer, &format);
564: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
565: if (format == PETSC_VIEWER_NATIVE) {
566: if (dm->useNatural) {
567: if (dm->sfNatural) {
568: if (ishdf5) {
569: #if defined(PETSC_HAVE_HDF5)
570: Vec v;
571: const char *vecname;
573: DMGetGlobalVector(dm, &v);
574: PetscObjectGetName((PetscObject) originalv, &vecname);
575: PetscObjectSetName((PetscObject) v, vecname);
576: VecLoad_Plex_HDF5_Native_Internal(v, viewer);
577: DMPlexNaturalToGlobalBegin(dm, v, originalv);
578: DMPlexNaturalToGlobalEnd(dm, v, originalv);
579: DMRestoreGlobalVector(dm, &v);
580: #else
581: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
582: #endif
583: } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Reading in natural order is not supported for anything but HDF5.");
584: }
585: } else {
586: VecLoad_Default(originalv, viewer);
587: }
588: }
589: return(0);
590: }
592: PETSC_UNUSED static PetscErrorCode DMPlexView_Ascii_Geometry(DM dm, PetscViewer viewer)
593: {
594: PetscSection coordSection;
595: Vec coordinates;
596: DMLabel depthLabel, celltypeLabel;
597: const char *name[4];
598: const PetscScalar *a;
599: PetscInt dim, pStart, pEnd, cStart, cEnd, c;
600: PetscErrorCode ierr;
603: DMGetDimension(dm, &dim);
604: DMGetCoordinatesLocal(dm, &coordinates);
605: DMGetCoordinateSection(dm, &coordSection);
606: DMPlexGetDepthLabel(dm, &depthLabel);
607: DMPlexGetCellTypeLabel(dm, &celltypeLabel);
608: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
609: PetscSectionGetChart(coordSection, &pStart, &pEnd);
610: VecGetArrayRead(coordinates, &a);
611: name[0] = "vertex";
612: name[1] = "edge";
613: name[dim-1] = "face";
614: name[dim] = "cell";
615: for (c = cStart; c < cEnd; ++c) {
616: PetscInt *closure = NULL;
617: PetscInt closureSize, cl, ct;
619: DMLabelGetValue(celltypeLabel, c, &ct);
620: PetscViewerASCIIPrintf(viewer, "Geometry for cell %D polytope type %s:\n", c, DMPolytopeTypes[ct]);
621: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
622: PetscViewerASCIIPushTab(viewer);
623: for (cl = 0; cl < closureSize*2; cl += 2) {
624: PetscInt point = closure[cl], depth, dof, off, d, p;
626: if ((point < pStart) || (point >= pEnd)) continue;
627: PetscSectionGetDof(coordSection, point, &dof);
628: if (!dof) continue;
629: DMLabelGetValue(depthLabel, point, &depth);
630: PetscSectionGetOffset(coordSection, point, &off);
631: PetscViewerASCIIPrintf(viewer, "%s %D coords:", name[depth], point);
632: for (p = 0; p < dof/dim; ++p) {
633: PetscViewerASCIIPrintf(viewer, " (");
634: for (d = 0; d < dim; ++d) {
635: if (d > 0) {PetscViewerASCIIPrintf(viewer, ", ");}
636: PetscViewerASCIIPrintf(viewer, "%g", (double) PetscRealPart(a[off+p*dim+d]));
637: }
638: PetscViewerASCIIPrintf(viewer, ")");
639: }
640: PetscViewerASCIIPrintf(viewer, "\n");
641: }
642: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
643: PetscViewerASCIIPopTab(viewer);
644: }
645: VecRestoreArrayRead(coordinates, &a);
646: return(0);
647: }
649: static PetscErrorCode DMPlexView_Ascii(DM dm, PetscViewer viewer)
650: {
651: DM_Plex *mesh = (DM_Plex*) dm->data;
652: DM cdm;
653: PetscSection coordSection;
654: Vec coordinates;
655: PetscViewerFormat format;
656: PetscErrorCode ierr;
659: DMGetCoordinateDM(dm, &cdm);
660: DMGetLocalSection(cdm, &coordSection);
661: DMGetCoordinatesLocal(dm, &coordinates);
662: PetscViewerGetFormat(viewer, &format);
663: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
664: const char *name;
665: PetscInt dim, cellHeight, maxConeSize, maxSupportSize;
666: PetscInt pStart, pEnd, p, numLabels, l;
667: PetscMPIInt rank, size;
669: MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
670: MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
671: PetscObjectGetName((PetscObject) dm, &name);
672: DMPlexGetChart(dm, &pStart, &pEnd);
673: DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);
674: DMGetDimension(dm, &dim);
675: DMPlexGetVTKCellHeight(dm, &cellHeight);
676: if (name) {PetscViewerASCIIPrintf(viewer, "%s in %D dimension%s:\n", name, dim, dim == 1 ? "" : "s");}
677: else {PetscViewerASCIIPrintf(viewer, "Mesh in %D dimension%s:\n", dim, dim == 1 ? "" : "s");}
678: if (cellHeight) {PetscViewerASCIIPrintf(viewer, " Cells are at height %D\n", cellHeight);}
679: PetscViewerASCIIPrintf(viewer, "Supports:\n", name);
680: PetscViewerASCIIPushSynchronized(viewer);
681: PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Max support size: %D\n", rank, maxSupportSize);
682: for (p = pStart; p < pEnd; ++p) {
683: PetscInt dof, off, s;
685: PetscSectionGetDof(mesh->supportSection, p, &dof);
686: PetscSectionGetOffset(mesh->supportSection, p, &off);
687: for (s = off; s < off+dof; ++s) {
688: PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D ----> %D\n", rank, p, mesh->supports[s]);
689: }
690: }
691: PetscViewerFlush(viewer);
692: PetscViewerASCIIPrintf(viewer, "Cones:\n", name);
693: PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Max cone size: %D\n", rank, maxConeSize);
694: for (p = pStart; p < pEnd; ++p) {
695: PetscInt dof, off, c;
697: PetscSectionGetDof(mesh->coneSection, p, &dof);
698: PetscSectionGetOffset(mesh->coneSection, p, &off);
699: for (c = off; c < off+dof; ++c) {
700: PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D <---- %D (%D)\n", rank, p, mesh->cones[c], mesh->coneOrientations[c]);
701: }
702: }
703: PetscViewerFlush(viewer);
704: PetscViewerASCIIPopSynchronized(viewer);
705: if (coordSection && coordinates) {
706: PetscSectionVecView(coordSection, coordinates, viewer);
707: }
708: DMGetNumLabels(dm, &numLabels);
709: if (numLabels) {PetscViewerASCIIPrintf(viewer, "Labels:\n");}
710: for (l = 0; l < numLabels; ++l) {
711: DMLabel label;
712: PetscBool isdepth;
713: const char *name;
715: DMGetLabelName(dm, l, &name);
716: PetscStrcmp(name, "depth", &isdepth);
717: if (isdepth) continue;
718: DMGetLabel(dm, name, &label);
719: DMLabelView(label, viewer);
720: }
721: if (size > 1) {
722: PetscSF sf;
724: DMGetPointSF(dm, &sf);
725: PetscSFView(sf, viewer);
726: }
727: PetscViewerFlush(viewer);
728: } else if (format == PETSC_VIEWER_ASCII_LATEX) {
729: const char *name, *color;
730: const char *defcolors[3] = {"gray", "orange", "green"};
731: const char *deflcolors[4] = {"blue", "cyan", "red", "magenta"};
732: char lname[PETSC_MAX_PATH_LEN];
733: PetscReal scale = 2.0;
734: PetscReal tikzscale = 1.0;
735: PetscBool useNumbers = PETSC_TRUE, drawNumbers[4], drawColors[4], useLabels, useColors, plotEdges, drawHasse = PETSC_FALSE;
736: double tcoords[3];
737: PetscScalar *coords;
738: PetscInt numLabels, l, numColors, numLColors, dim, d, depth, cStart, cEnd, c, vStart, vEnd, v, eStart = 0, eEnd = 0, e, p, n;
739: PetscMPIInt rank, size;
740: char **names, **colors, **lcolors;
741: PetscBool flg, lflg;
742: PetscBT wp = NULL;
743: PetscInt pEnd, pStart;
745: DMGetDimension(dm, &dim);
746: DMPlexGetDepth(dm, &depth);
747: DMGetNumLabels(dm, &numLabels);
748: numLabels = PetscMax(numLabels, 10);
749: numColors = 10;
750: numLColors = 10;
751: PetscCalloc3(numLabels, &names, numColors, &colors, numLColors, &lcolors);
752: PetscOptionsGetReal(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_scale", &scale, NULL);
753: PetscOptionsGetReal(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_tikzscale", &tikzscale, NULL);
754: PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_numbers", &useNumbers, NULL);
755: for (d = 0; d < 4; ++d) drawNumbers[d] = useNumbers;
756: for (d = 0; d < 4; ++d) drawColors[d] = PETSC_TRUE;
757: n = 4;
758: PetscOptionsGetBoolArray(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_numbers_depth", drawNumbers, &n, &flg);
759: if (flg && n != dim+1) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Number of flags %D != %D dim+1", n, dim+1);
760: PetscOptionsGetBoolArray(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_colors_depth", drawColors, &n, &flg);
761: if (flg && n != dim+1) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Number of flags %D != %D dim+1", n, dim+1);
762: PetscOptionsGetStringArray(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_labels", names, &numLabels, &useLabels);
763: if (!useLabels) numLabels = 0;
764: PetscOptionsGetStringArray(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_colors", colors, &numColors, &useColors);
765: if (!useColors) {
766: numColors = 3;
767: for (c = 0; c < numColors; ++c) {PetscStrallocpy(defcolors[c], &colors[c]);}
768: }
769: PetscOptionsGetStringArray(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_lcolors", lcolors, &numLColors, &useColors);
770: if (!useColors) {
771: numLColors = 4;
772: for (c = 0; c < numLColors; ++c) {PetscStrallocpy(deflcolors[c], &lcolors[c]);}
773: }
774: PetscOptionsGetString(((PetscObject) viewer)->options, ((PetscObject) viewer)->prefix, "-dm_plex_view_label_filter", lname, sizeof(lname), &lflg);
775: plotEdges = (PetscBool)(depth > 1 && drawNumbers[1] && dim < 3);
776: PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_edges", &plotEdges, &flg);
777: if (flg && plotEdges && depth < dim) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Mesh must be interpolated");
778: if (depth < dim) plotEdges = PETSC_FALSE;
779: PetscOptionsGetBool(((PetscObject) viewer)->options, ((PetscObject) viewer)->prefix, "-dm_plex_view_hasse", &drawHasse, NULL);
781: /* filter points with labelvalue != labeldefaultvalue */
782: DMPlexGetChart(dm, &pStart, &pEnd);
783: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
784: DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);
785: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
786: if (lflg) {
787: DMLabel lbl;
789: DMGetLabel(dm, lname, &lbl);
790: if (lbl) {
791: PetscInt val, defval;
793: DMLabelGetDefaultValue(lbl, &defval);
794: PetscBTCreate(pEnd-pStart, &wp);
795: for (c = pStart; c < pEnd; c++) {
796: PetscInt *closure = NULL;
797: PetscInt closureSize;
799: DMLabelGetValue(lbl, c, &val);
800: if (val == defval) continue;
802: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
803: for (p = 0; p < closureSize*2; p += 2) {
804: PetscBTSet(wp, closure[p] - pStart);
805: }
806: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
807: }
808: }
809: }
811: MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
812: MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
813: PetscObjectGetName((PetscObject) dm, &name);
814: PetscViewerASCIIPrintf(viewer, "\
815: \\documentclass[tikz]{standalone}\n\n\
816: \\usepackage{pgflibraryshapes}\n\
817: \\usetikzlibrary{backgrounds}\n\
818: \\usetikzlibrary{arrows}\n\
819: \\begin{document}\n");
820: if (size > 1) {
821: PetscViewerASCIIPrintf(viewer, "%s for process ", name);
822: for (p = 0; p < size; ++p) {
823: if (p > 0 && p == size-1) {
824: PetscViewerASCIIPrintf(viewer, ", and ", colors[p%numColors], p);
825: } else if (p > 0) {
826: PetscViewerASCIIPrintf(viewer, ", ", colors[p%numColors], p);
827: }
828: PetscViewerASCIIPrintf(viewer, "{\\textcolor{%s}%D}", colors[p%numColors], p);
829: }
830: PetscViewerASCIIPrintf(viewer, ".\n\n\n");
831: }
832: if (drawHasse) {
833: PetscInt maxStratum = PetscMax(vEnd-vStart, PetscMax(eEnd-eStart, cEnd-cStart));
835: PetscViewerASCIIPrintf(viewer, "\\newcommand{\\vStart}{%D}\n", vStart);
836: PetscViewerASCIIPrintf(viewer, "\\newcommand{\\vEnd}{%D}\n", vEnd-1);
837: PetscViewerASCIIPrintf(viewer, "\\newcommand{\\numVertices}{%D}\n", vEnd-vStart);
838: PetscViewerASCIIPrintf(viewer, "\\newcommand{\\vShift}{%.2f}\n", 3 + (maxStratum-(vEnd-vStart))/2.);
839: PetscViewerASCIIPrintf(viewer, "\\newcommand{\\eStart}{%D}\n", eStart);
840: PetscViewerASCIIPrintf(viewer, "\\newcommand{\\eEnd}{%D}\n", eEnd-1);
841: PetscViewerASCIIPrintf(viewer, "\\newcommand{\\eShift}{%.2f}\n", 3 + (maxStratum-(eEnd-eStart))/2.);
842: PetscViewerASCIIPrintf(viewer, "\\newcommand{\\numEdges}{%D}\n", eEnd-eStart);
843: PetscViewerASCIIPrintf(viewer, "\\newcommand{\\cStart}{%D}\n", cStart);
844: PetscViewerASCIIPrintf(viewer, "\\newcommand{\\cEnd}{%D}\n", cEnd-1);
845: PetscViewerASCIIPrintf(viewer, "\\newcommand{\\numCells}{%D}\n", cEnd-cStart);
846: PetscViewerASCIIPrintf(viewer, "\\newcommand{\\cShift}{%.2f}\n", 3 + (maxStratum-(cEnd-cStart))/2.);
847: }
848: PetscViewerASCIIPrintf(viewer, "\\begin{tikzpicture}[scale = %g,font=\\fontsize{8}{8}\\selectfont]\n", (double) tikzscale);
850: /* Plot vertices */
851: VecGetArray(coordinates, &coords);
852: PetscViewerASCIIPushSynchronized(viewer);
853: for (v = vStart; v < vEnd; ++v) {
854: PetscInt off, dof, d;
855: PetscBool isLabeled = PETSC_FALSE;
857: if (wp && !PetscBTLookup(wp,v - pStart)) continue;
858: PetscSectionGetDof(coordSection, v, &dof);
859: PetscSectionGetOffset(coordSection, v, &off);
860: PetscViewerASCIISynchronizedPrintf(viewer, "\\path (");
861: if (PetscUnlikely(dof > 3)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"coordSection vertex %D has dof %D > 3",v,dof);
862: for (d = 0; d < dof; ++d) {
863: tcoords[d] = (double) (scale*PetscRealPart(coords[off+d]));
864: tcoords[d] = PetscAbs(tcoords[d]) < 1e-10 ? 0.0 : tcoords[d];
865: }
866: /* Rotate coordinates since PGF makes z point out of the page instead of up */
867: if (dim == 3) {PetscReal tmp = tcoords[1]; tcoords[1] = tcoords[2]; tcoords[2] = -tmp;}
868: for (d = 0; d < dof; ++d) {
869: if (d > 0) {PetscViewerASCIISynchronizedPrintf(viewer, ",");}
870: PetscViewerASCIISynchronizedPrintf(viewer, "%g", (double) tcoords[d]);
871: }
872: if (drawHasse) color = colors[0%numColors];
873: else color = colors[rank%numColors];
874: for (l = 0; l < numLabels; ++l) {
875: PetscInt val;
876: DMGetLabelValue(dm, names[l], v, &val);
877: if (val >= 0) {color = lcolors[l%numLColors]; isLabeled = PETSC_TRUE; break;}
878: }
879: if (drawNumbers[0]) {
880: PetscViewerASCIISynchronizedPrintf(viewer, ") node(%D_%d) [draw,shape=circle,color=%s] {%D};\n", v, rank, color, v);
881: } else if (drawColors[0]) {
882: PetscViewerASCIISynchronizedPrintf(viewer, ") node(%D_%d) [fill,inner sep=%dpt,shape=circle,color=%s] {};\n", v, rank, !isLabeled ? 1 : 2, color);
883: } else {
884: PetscViewerASCIISynchronizedPrintf(viewer, ") node(%D_%d) [] {};\n", v, rank);
885: }
886: }
887: VecRestoreArray(coordinates, &coords);
888: PetscViewerFlush(viewer);
889: /* Plot edges */
890: if (plotEdges) {
891: VecGetArray(coordinates, &coords);
892: PetscViewerASCIIPrintf(viewer, "\\path\n");
893: for (e = eStart; e < eEnd; ++e) {
894: const PetscInt *cone;
895: PetscInt coneSize, offA, offB, dof, d;
897: if (wp && !PetscBTLookup(wp,e - pStart)) continue;
898: DMPlexGetConeSize(dm, e, &coneSize);
899: if (coneSize != 2) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Edge %D cone should have two vertices, not %D", e, coneSize);
900: DMPlexGetCone(dm, e, &cone);
901: PetscSectionGetDof(coordSection, cone[0], &dof);
902: PetscSectionGetOffset(coordSection, cone[0], &offA);
903: PetscSectionGetOffset(coordSection, cone[1], &offB);
904: PetscViewerASCIISynchronizedPrintf(viewer, "(");
905: for (d = 0; d < dof; ++d) {
906: tcoords[d] = (double) (0.5*scale*PetscRealPart(coords[offA+d]+coords[offB+d]));
907: tcoords[d] = PetscAbs(tcoords[d]) < 1e-10 ? 0.0 : tcoords[d];
908: }
909: /* Rotate coordinates since PGF makes z point out of the page instead of up */
910: if (dim == 3) {PetscReal tmp = tcoords[1]; tcoords[1] = tcoords[2]; tcoords[2] = -tmp;}
911: for (d = 0; d < dof; ++d) {
912: if (d > 0) {PetscViewerASCIISynchronizedPrintf(viewer, ",");}
913: PetscViewerASCIISynchronizedPrintf(viewer, "%g", (double)tcoords[d]);
914: }
915: if (drawHasse) color = colors[1%numColors];
916: else color = colors[rank%numColors];
917: for (l = 0; l < numLabels; ++l) {
918: PetscInt val;
919: DMGetLabelValue(dm, names[l], v, &val);
920: if (val >= 0) {color = lcolors[l%numLColors]; break;}
921: }
922: PetscViewerASCIISynchronizedPrintf(viewer, ") node(%D_%d) [draw,shape=circle,color=%s] {%D} --\n", e, rank, color, e);
923: }
924: VecRestoreArray(coordinates, &coords);
925: PetscViewerFlush(viewer);
926: PetscViewerASCIIPrintf(viewer, "(0,0);\n");
927: }
928: /* Plot cells */
929: if (dim == 3 || !drawNumbers[1]) {
930: for (e = eStart; e < eEnd; ++e) {
931: const PetscInt *cone;
933: if (wp && !PetscBTLookup(wp,e - pStart)) continue;
934: color = colors[rank%numColors];
935: for (l = 0; l < numLabels; ++l) {
936: PetscInt val;
937: DMGetLabelValue(dm, names[l], e, &val);
938: if (val >= 0) {color = lcolors[l%numLColors]; break;}
939: }
940: DMPlexGetCone(dm, e, &cone);
941: PetscViewerASCIISynchronizedPrintf(viewer, "\\draw[color=%s] (%D_%d) -- (%D_%d);\n", color, cone[0], rank, cone[1], rank);
942: }
943: } else {
944: DMPolytopeType ct;
946: /* Drawing a 2D polygon */
947: for (c = cStart; c < cEnd; ++c) {
948: if (wp && !PetscBTLookup(wp, c - pStart)) continue;
949: DMPlexGetCellType(dm, c, &ct);
950: if (ct == DM_POLYTOPE_SEG_PRISM_TENSOR ||
951: ct == DM_POLYTOPE_TRI_PRISM_TENSOR ||
952: ct == DM_POLYTOPE_QUAD_PRISM_TENSOR) {
953: const PetscInt *cone;
954: PetscInt coneSize, e;
956: DMPlexGetCone(dm, c, &cone);
957: DMPlexGetConeSize(dm, c, &coneSize);
958: for (e = 0; e < coneSize; ++e) {
959: const PetscInt *econe;
961: DMPlexGetCone(dm, cone[e], &econe);
962: PetscViewerASCIISynchronizedPrintf(viewer, "\\draw[color=%s] (%D_%d) -- (%D_%d) -- (%D_%d);\n", colors[rank%numColors], econe[0], rank, cone[e], rank, econe[1], rank);
963: }
964: } else {
965: PetscInt *closure = NULL;
966: PetscInt closureSize, Nv = 0, v;
968: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
969: for (p = 0; p < closureSize*2; p += 2) {
970: const PetscInt point = closure[p];
972: if ((point >= vStart) && (point < vEnd)) closure[Nv++] = point;
973: }
974: PetscViewerASCIISynchronizedPrintf(viewer, "\\draw[color=%s] ", colors[rank%numColors]);
975: for (v = 0; v <= Nv; ++v) {
976: const PetscInt vertex = closure[v%Nv];
978: if (v > 0) {
979: if (plotEdges) {
980: const PetscInt *edge;
981: PetscInt endpoints[2], ne;
983: endpoints[0] = closure[v-1]; endpoints[1] = vertex;
984: DMPlexGetJoin(dm, 2, endpoints, &ne, &edge);
985: if (ne != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find edge for vertices %D, %D", endpoints[0], endpoints[1]);
986: PetscViewerASCIISynchronizedPrintf(viewer, " -- (%D_%d) -- ", edge[0], rank);
987: DMPlexRestoreJoin(dm, 2, endpoints, &ne, &edge);
988: } else {
989: PetscViewerASCIISynchronizedPrintf(viewer, " -- ");
990: }
991: }
992: PetscViewerASCIISynchronizedPrintf(viewer, "(%D_%d)", vertex, rank);
993: }
994: PetscViewerASCIISynchronizedPrintf(viewer, ";\n");
995: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
996: }
997: }
998: }
999: VecGetArray(coordinates, &coords);
1000: for (c = cStart; c < cEnd; ++c) {
1001: double ccoords[3] = {0.0, 0.0, 0.0};
1002: PetscBool isLabeled = PETSC_FALSE;
1003: PetscInt *closure = NULL;
1004: PetscInt closureSize, dof, d, n = 0;
1006: if (wp && !PetscBTLookup(wp,c - pStart)) continue;
1007: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
1008: PetscViewerASCIISynchronizedPrintf(viewer, "\\path (");
1009: for (p = 0; p < closureSize*2; p += 2) {
1010: const PetscInt point = closure[p];
1011: PetscInt off;
1013: if ((point < vStart) || (point >= vEnd)) continue;
1014: PetscSectionGetDof(coordSection, point, &dof);
1015: PetscSectionGetOffset(coordSection, point, &off);
1016: for (d = 0; d < dof; ++d) {
1017: tcoords[d] = (double) (scale*PetscRealPart(coords[off+d]));
1018: tcoords[d] = PetscAbs(tcoords[d]) < 1e-10 ? 0.0 : tcoords[d];
1019: }
1020: /* Rotate coordinates since PGF makes z point out of the page instead of up */
1021: if (dof == 3) {PetscReal tmp = tcoords[1]; tcoords[1] = tcoords[2]; tcoords[2] = -tmp;}
1022: for (d = 0; d < dof; ++d) {ccoords[d] += tcoords[d];}
1023: ++n;
1024: }
1025: for (d = 0; d < dof; ++d) {ccoords[d] /= n;}
1026: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
1027: for (d = 0; d < dof; ++d) {
1028: if (d > 0) {PetscViewerASCIISynchronizedPrintf(viewer, ",");}
1029: PetscViewerASCIISynchronizedPrintf(viewer, "%g", (double) ccoords[d]);
1030: }
1031: if (drawHasse) color = colors[depth%numColors];
1032: else color = colors[rank%numColors];
1033: for (l = 0; l < numLabels; ++l) {
1034: PetscInt val;
1035: DMGetLabelValue(dm, names[l], c, &val);
1036: if (val >= 0) {color = lcolors[l%numLColors]; isLabeled = PETSC_TRUE; break;}
1037: }
1038: if (drawNumbers[dim]) {
1039: PetscViewerASCIISynchronizedPrintf(viewer, ") node(%D_%d) [draw,shape=circle,color=%s] {%D};\n", c, rank, color, c);
1040: } else if (drawColors[dim]) {
1041: PetscViewerASCIISynchronizedPrintf(viewer, ") node(%D_%d) [fill,inner sep=%dpt,shape=circle,color=%s] {};\n", c, rank, !isLabeled ? 1 : 2, color);
1042: } else {
1043: PetscViewerASCIISynchronizedPrintf(viewer, ") node(%D_%d) [] {};\n", c, rank);
1044: }
1045: }
1046: VecRestoreArray(coordinates, &coords);
1047: if (drawHasse) {
1048: color = colors[depth%numColors];
1049: PetscViewerASCIIPrintf(viewer, "%% Cells\n");
1050: PetscViewerASCIIPrintf(viewer, "\\foreach \\c in {\\cStart,...,\\cEnd}\n");
1051: PetscViewerASCIIPrintf(viewer, "{\n");
1052: PetscViewerASCIIPrintf(viewer, " \\node(\\c_%d) [draw,shape=circle,color=%s,minimum size = 6mm] at (\\cShift+\\c-\\cStart,0) {\\c};\n", rank, color);
1053: PetscViewerASCIIPrintf(viewer, "}\n");
1055: color = colors[1%numColors];
1056: PetscViewerASCIIPrintf(viewer, "%% Edges\n");
1057: PetscViewerASCIIPrintf(viewer, "\\foreach \\e in {\\eStart,...,\\eEnd}\n");
1058: PetscViewerASCIIPrintf(viewer, "{\n");
1059: PetscViewerASCIIPrintf(viewer, " \\node(\\e_%d) [draw,shape=circle,color=%s,minimum size = 6mm] at (\\eShift+\\e-\\eStart,1) {\\e};\n", rank, color);
1060: PetscViewerASCIIPrintf(viewer, "}\n");
1062: color = colors[0%numColors];
1063: PetscViewerASCIIPrintf(viewer, "%% Vertices\n");
1064: PetscViewerASCIIPrintf(viewer, "\\foreach \\v in {\\vStart,...,\\vEnd}\n");
1065: PetscViewerASCIIPrintf(viewer, "{\n");
1066: PetscViewerASCIIPrintf(viewer, " \\node(\\v_%d) [draw,shape=circle,color=%s,minimum size = 6mm] at (\\vShift+\\v-\\vStart,2) {\\v};\n", rank, color);
1067: PetscViewerASCIIPrintf(viewer, "}\n");
1069: for (p = pStart; p < pEnd; ++p) {
1070: const PetscInt *cone;
1071: PetscInt coneSize, cp;
1073: DMPlexGetCone(dm, p, &cone);
1074: DMPlexGetConeSize(dm, p, &coneSize);
1075: for (cp = 0; cp < coneSize; ++cp) {
1076: PetscViewerASCIIPrintf(viewer, "\\draw[->, shorten >=1pt] (%D_%d) -- (%D_%d);\n", cone[cp], rank, p, rank);
1077: }
1078: }
1079: }
1080: PetscViewerFlush(viewer);
1081: PetscViewerASCIIPopSynchronized(viewer);
1082: PetscViewerASCIIPrintf(viewer, "\\end{tikzpicture}\n");
1083: PetscViewerASCIIPrintf(viewer, "\\end{document}\n", name);
1084: for (l = 0; l < numLabels; ++l) {PetscFree(names[l]);}
1085: for (c = 0; c < numColors; ++c) {PetscFree(colors[c]);}
1086: for (c = 0; c < numLColors; ++c) {PetscFree(lcolors[c]);}
1087: PetscFree3(names, colors, lcolors);
1088: PetscBTDestroy(&wp);
1089: } else if (format == PETSC_VIEWER_LOAD_BALANCE) {
1090: Vec cown,acown;
1091: VecScatter sct;
1092: ISLocalToGlobalMapping g2l;
1093: IS gid,acis;
1094: MPI_Comm comm,ncomm = MPI_COMM_NULL;
1095: MPI_Group ggroup,ngroup;
1096: PetscScalar *array,nid;
1097: const PetscInt *idxs;
1098: PetscInt *idxs2,*start,*adjacency,*work;
1099: PetscInt64 lm[3],gm[3];
1100: PetscInt i,c,cStart,cEnd,cum,numVertices,ect,ectn,cellHeight;
1101: PetscMPIInt d1,d2,rank;
1103: PetscObjectGetComm((PetscObject)dm,&comm);
1104: MPI_Comm_rank(comm,&rank);
1105: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
1106: MPI_Comm_split_type(comm,MPI_COMM_TYPE_SHARED,rank,MPI_INFO_NULL,&ncomm);
1107: #endif
1108: if (ncomm != MPI_COMM_NULL) {
1109: MPI_Comm_group(comm,&ggroup);
1110: MPI_Comm_group(ncomm,&ngroup);
1111: d1 = 0;
1112: MPI_Group_translate_ranks(ngroup,1,&d1,ggroup,&d2);
1113: nid = d2;
1114: MPI_Group_free(&ggroup);
1115: MPI_Group_free(&ngroup);
1116: MPI_Comm_free(&ncomm);
1117: } else nid = 0.0;
1119: /* Get connectivity */
1120: DMPlexGetVTKCellHeight(dm,&cellHeight);
1121: DMPlexCreatePartitionerGraph(dm,cellHeight,&numVertices,&start,&adjacency,&gid);
1123: /* filter overlapped local cells */
1124: DMPlexGetHeightStratum(dm,cellHeight,&cStart,&cEnd);
1125: ISGetIndices(gid,&idxs);
1126: ISGetLocalSize(gid,&cum);
1127: PetscMalloc1(cum,&idxs2);
1128: for (c = cStart, cum = 0; c < cEnd; c++) {
1129: if (idxs[c-cStart] < 0) continue;
1130: idxs2[cum++] = idxs[c-cStart];
1131: }
1132: ISRestoreIndices(gid,&idxs);
1133: if (numVertices != cum) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D != %D",numVertices,cum);
1134: ISDestroy(&gid);
1135: ISCreateGeneral(comm,numVertices,idxs2,PETSC_OWN_POINTER,&gid);
1137: /* support for node-aware cell locality */
1138: ISCreateGeneral(comm,start[numVertices],adjacency,PETSC_USE_POINTER,&acis);
1139: VecCreateSeq(PETSC_COMM_SELF,start[numVertices],&acown);
1140: VecCreateMPI(comm,numVertices,PETSC_DECIDE,&cown);
1141: VecGetArray(cown,&array);
1142: for (c = 0; c < numVertices; c++) array[c] = nid;
1143: VecRestoreArray(cown,&array);
1144: VecScatterCreate(cown,acis,acown,NULL,&sct);
1145: VecScatterBegin(sct,cown,acown,INSERT_VALUES,SCATTER_FORWARD);
1146: VecScatterEnd(sct,cown,acown,INSERT_VALUES,SCATTER_FORWARD);
1147: ISDestroy(&acis);
1148: VecScatterDestroy(&sct);
1149: VecDestroy(&cown);
1151: /* compute edgeCut */
1152: for (c = 0, cum = 0; c < numVertices; c++) cum = PetscMax(cum,start[c+1]-start[c]);
1153: PetscMalloc1(cum,&work);
1154: ISLocalToGlobalMappingCreateIS(gid,&g2l);
1155: ISLocalToGlobalMappingSetType(g2l,ISLOCALTOGLOBALMAPPINGHASH);
1156: ISDestroy(&gid);
1157: VecGetArray(acown,&array);
1158: for (c = 0, ect = 0, ectn = 0; c < numVertices; c++) {
1159: PetscInt totl;
1161: totl = start[c+1]-start[c];
1162: ISGlobalToLocalMappingApply(g2l,IS_GTOLM_MASK,totl,adjacency+start[c],NULL,work);
1163: for (i = 0; i < totl; i++) {
1164: if (work[i] < 0) {
1165: ect += 1;
1166: ectn += (array[i + start[c]] != nid) ? 0 : 1;
1167: }
1168: }
1169: }
1170: PetscFree(work);
1171: VecRestoreArray(acown,&array);
1172: lm[0] = numVertices > 0 ? numVertices : PETSC_MAX_INT;
1173: lm[1] = -numVertices;
1174: MPIU_Allreduce(lm,gm,2,MPIU_INT64,MPI_MIN,comm);
1175: PetscViewerASCIIPrintf(viewer," Cell balance: %.2f (max %D, min %D",-((double)gm[1])/((double)gm[0]),-(PetscInt)gm[1],(PetscInt)gm[0]);
1176: lm[0] = ect; /* edgeCut */
1177: lm[1] = ectn; /* node-aware edgeCut */
1178: lm[2] = numVertices > 0 ? 0 : 1; /* empty processes */
1179: MPIU_Allreduce(lm,gm,3,MPIU_INT64,MPI_SUM,comm);
1180: PetscViewerASCIIPrintf(viewer,", empty %D)\n",(PetscInt)gm[2]);
1181: #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
1182: PetscViewerASCIIPrintf(viewer," Edge Cut: %D (on node %.3f)\n",(PetscInt)(gm[0]/2),gm[0] ? ((double)(gm[1]))/((double)gm[0]) : 1.);
1183: #else
1184: PetscViewerASCIIPrintf(viewer," Edge Cut: %D (on node %.3f)\n",(PetscInt)(gm[0]/2),0.0);
1185: #endif
1186: ISLocalToGlobalMappingDestroy(&g2l);
1187: PetscFree(start);
1188: PetscFree(adjacency);
1189: VecDestroy(&acown);
1190: } else {
1191: const char *name;
1192: PetscInt *sizes, *hybsizes, *ghostsizes;
1193: PetscInt locDepth, depth, cellHeight, dim, d;
1194: PetscInt pStart, pEnd, p, gcStart, gcEnd, gcNum;
1195: PetscInt numLabels, l;
1196: DMPolytopeType ct0 = DM_POLYTOPE_UNKNOWN;
1197: MPI_Comm comm;
1198: PetscMPIInt size, rank;
1200: PetscObjectGetComm((PetscObject) dm, &comm);
1201: MPI_Comm_size(comm, &size);
1202: MPI_Comm_rank(comm, &rank);
1203: DMGetDimension(dm, &dim);
1204: DMPlexGetVTKCellHeight(dm, &cellHeight);
1205: PetscObjectGetName((PetscObject) dm, &name);
1206: if (name) {PetscViewerASCIIPrintf(viewer, "%s in %D dimension%s:\n", name, dim, dim == 1 ? "" : "s");}
1207: else {PetscViewerASCIIPrintf(viewer, "Mesh in %D dimension%s:\n", dim, dim == 1 ? "" : "s");}
1208: if (cellHeight) {PetscViewerASCIIPrintf(viewer, " Cells are at height %D\n", cellHeight);}
1209: DMPlexGetDepth(dm, &locDepth);
1210: MPIU_Allreduce(&locDepth, &depth, 1, MPIU_INT, MPI_MAX, comm);
1211: DMPlexGetGhostCellStratum(dm, &gcStart, &gcEnd);
1212: gcNum = gcEnd - gcStart;
1213: PetscCalloc3(size,&sizes,size,&hybsizes,size,&ghostsizes);
1214: for (d = 0; d <= depth; d++) {
1215: PetscInt Nc[2] = {0, 0}, ict;
1217: DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
1218: if (pStart < pEnd) {DMPlexGetCellType(dm, pStart, &ct0);}
1219: ict = ct0;
1220: MPI_Bcast(&ict, 1, MPIU_INT, 0, comm);
1221: ct0 = (DMPolytopeType) ict;
1222: for (p = pStart; p < pEnd; ++p) {
1223: DMPolytopeType ct;
1225: DMPlexGetCellType(dm, p, &ct);
1226: if (ct == ct0) ++Nc[0];
1227: else ++Nc[1];
1228: }
1229: MPI_Gather(&Nc[0], 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm);
1230: MPI_Gather(&Nc[1], 1, MPIU_INT, hybsizes, 1, MPIU_INT, 0, comm);
1231: if (d == depth) {MPI_Gather(&gcNum, 1, MPIU_INT, ghostsizes, 1, MPIU_INT, 0, comm);}
1232: PetscViewerASCIIPrintf(viewer, " %D-cells:", (depth == 1) && d ? dim : d);
1233: for (p = 0; p < size; ++p) {
1234: if (rank == 0) {
1235: PetscViewerASCIIPrintf(viewer, " %D", sizes[p]+hybsizes[p]);
1236: if (hybsizes[p] > 0) {PetscViewerASCIIPrintf(viewer, " (%D)", hybsizes[p]);}
1237: if (ghostsizes[p] > 0) {PetscViewerASCIIPrintf(viewer, " [%D]", ghostsizes[p]);}
1238: }
1239: }
1240: PetscViewerASCIIPrintf(viewer, "\n");
1241: }
1242: PetscFree3(sizes,hybsizes,ghostsizes);
1243: {
1244: const PetscReal *maxCell;
1245: const PetscReal *L;
1246: const DMBoundaryType *bd;
1247: PetscBool per, localized;
1249: DMGetPeriodicity(dm, &per, &maxCell, &L, &bd);
1250: DMGetCoordinatesLocalized(dm, &localized);
1251: if (per) {
1252: PetscViewerASCIIPrintf(viewer, "Periodic mesh (");
1253: PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);
1254: for (d = 0; d < dim; ++d) {
1255: if (bd && d > 0) {PetscViewerASCIIPrintf(viewer, ", ");}
1256: if (bd) {PetscViewerASCIIPrintf(viewer, "%s", DMBoundaryTypes[bd[d]]);}
1257: }
1258: PetscViewerASCIIPrintf(viewer, ") coordinates %s\n", localized ? "localized" : "not localized");
1259: PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);
1260: }
1261: }
1262: DMGetNumLabels(dm, &numLabels);
1263: if (numLabels) {PetscViewerASCIIPrintf(viewer, "Labels:\n");}
1264: for (l = 0; l < numLabels; ++l) {
1265: DMLabel label;
1266: const char *name;
1267: IS valueIS;
1268: const PetscInt *values;
1269: PetscInt numValues, v;
1271: DMGetLabelName(dm, l, &name);
1272: DMGetLabel(dm, name, &label);
1273: DMLabelGetNumValues(label, &numValues);
1274: PetscViewerASCIIPrintf(viewer, " %s: %D strata with value/size (", name, numValues);
1275: DMLabelGetValueIS(label, &valueIS);
1276: ISGetIndices(valueIS, &values);
1277: PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);
1278: for (v = 0; v < numValues; ++v) {
1279: PetscInt size;
1281: DMLabelGetStratumSize(label, values[v], &size);
1282: if (v > 0) {PetscViewerASCIIPrintf(viewer, ", ");}
1283: PetscViewerASCIIPrintf(viewer, "%D (%D)", values[v], size);
1284: }
1285: PetscViewerASCIIPrintf(viewer, ")\n");
1286: PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);
1287: ISRestoreIndices(valueIS, &values);
1288: ISDestroy(&valueIS);
1289: }
1290: {
1291: char **labelNames;
1292: PetscInt Nl = numLabels;
1293: PetscBool flg;
1295: PetscMalloc1(Nl, &labelNames);
1296: PetscOptionsGetStringArray(((PetscObject) dm)->options, ((PetscObject) dm)->prefix, "-dm_plex_view_labels", labelNames, &Nl, &flg);
1297: for (l = 0; l < Nl; ++l) {
1298: DMLabel label;
1300: DMHasLabel(dm, labelNames[l], &flg);
1301: if (flg) {
1302: DMGetLabel(dm, labelNames[l], &label);
1303: DMLabelView(label, viewer);
1304: }
1305: PetscFree(labelNames[l]);
1306: }
1307: PetscFree(labelNames);
1308: }
1309: /* If no fields are specified, people do not want to see adjacency */
1310: if (dm->Nf) {
1311: PetscInt f;
1313: for (f = 0; f < dm->Nf; ++f) {
1314: const char *name;
1316: PetscObjectGetName(dm->fields[f].disc, &name);
1317: if (numLabels) {PetscViewerASCIIPrintf(viewer, "Field %s:\n", name);}
1318: PetscViewerASCIIPushTab(viewer);
1319: if (dm->fields[f].label) {DMLabelView(dm->fields[f].label, viewer);}
1320: if (dm->fields[f].adjacency[0]) {
1321: if (dm->fields[f].adjacency[1]) {PetscViewerASCIIPrintf(viewer, "adjacency FVM++\n");}
1322: else {PetscViewerASCIIPrintf(viewer, "adjacency FVM\n");}
1323: } else {
1324: if (dm->fields[f].adjacency[1]) {PetscViewerASCIIPrintf(viewer, "adjacency FEM\n");}
1325: else {PetscViewerASCIIPrintf(viewer, "adjacency FUNKY\n");}
1326: }
1327: PetscViewerASCIIPopTab(viewer);
1328: }
1329: }
1330: DMGetCoarseDM(dm, &cdm);
1331: if (cdm) {
1332: PetscViewerASCIIPushTab(viewer);
1333: DMPlexView_Ascii(cdm, viewer);
1334: PetscViewerASCIIPopTab(viewer);
1335: }
1336: }
1337: return(0);
1338: }
1340: static PetscErrorCode DMPlexDrawCell(DM dm, PetscDraw draw, PetscInt cell, const PetscScalar coords[])
1341: {
1342: DMPolytopeType ct;
1343: PetscMPIInt rank;
1347: MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
1348: DMPlexGetCellType(dm, cell, &ct);
1349: switch (ct) {
1350: case DM_POLYTOPE_TRIANGLE:
1351: PetscDrawTriangle(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]),
1352: PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS-2) + 2,
1353: PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS-2) + 2,
1354: PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS-2) + 2);
1355: PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PETSC_DRAW_BLACK);
1356: PetscDrawLine(draw, PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_BLACK);
1357: PetscDrawLine(draw, PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[0]), PetscRealPart(coords[1]), PETSC_DRAW_BLACK);
1358: break;
1359: case DM_POLYTOPE_QUADRILATERAL:
1360: PetscDrawTriangle(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]),
1361: PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS-2) + 2,
1362: PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS-2) + 2,
1363: PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS-2) + 2);
1364: PetscDrawTriangle(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[6]), PetscRealPart(coords[7]),
1365: PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS-2) + 2,
1366: PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS-2) + 2,
1367: PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS-2) + 2);
1368: PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PETSC_DRAW_BLACK);
1369: PetscDrawLine(draw, PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_BLACK);
1370: PetscDrawLine(draw, PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[6]), PetscRealPart(coords[7]), PETSC_DRAW_BLACK);
1371: PetscDrawLine(draw, PetscRealPart(coords[6]), PetscRealPart(coords[7]), PetscRealPart(coords[0]), PetscRealPart(coords[1]), PETSC_DRAW_BLACK);
1372: break;
1373: default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot draw cells of type %s", DMPolytopeTypes[ct]);
1374: }
1375: return(0);
1376: }
1378: static PetscErrorCode DMPlexDrawCellHighOrder(DM dm, PetscDraw draw, PetscInt cell, const PetscScalar coords[], PetscInt edgeDiv, PetscReal refCoords[], PetscReal edgeCoords[])
1379: {
1380: DMPolytopeType ct;
1381: PetscReal centroid[2] = {0., 0.};
1382: PetscMPIInt rank;
1383: PetscInt fillColor, v, e, d;
1387: MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
1388: DMPlexGetCellType(dm, cell, &ct);
1389: fillColor = PETSC_DRAW_WHITE + rank % (PETSC_DRAW_BASIC_COLORS-2) + 2;
1390: switch (ct) {
1391: case DM_POLYTOPE_TRIANGLE:
1392: {
1393: PetscReal refVertices[6] = {-1., -1., 1., -1., -1., 1.};
1395: for (v = 0; v < 3; ++v) {centroid[0] += PetscRealPart(coords[v*2+0])/3.;centroid[1] += PetscRealPart(coords[v*2+1])/3.;}
1396: for (e = 0; e < 3; ++e) {
1397: refCoords[0] = refVertices[e*2+0];
1398: refCoords[1] = refVertices[e*2+1];
1399: for (d = 1; d <= edgeDiv; ++d) {
1400: refCoords[d*2+0] = refCoords[0] + (refVertices[(e+1)%3 * 2 + 0] - refCoords[0])*d/edgeDiv;
1401: refCoords[d*2+1] = refCoords[1] + (refVertices[(e+1)%3 * 2 + 1] - refCoords[1])*d/edgeDiv;
1402: }
1403: DMPlexReferenceToCoordinates(dm, cell, edgeDiv+1, refCoords, edgeCoords);
1404: for (d = 0; d < edgeDiv; ++d) {
1405: PetscDrawTriangle(draw, centroid[0], centroid[1], edgeCoords[d*2+0], edgeCoords[d*2+1], edgeCoords[(d+1)*2+0], edgeCoords[(d+1)*2+1], fillColor, fillColor, fillColor);
1406: PetscDrawLine(draw, edgeCoords[d*2+0], edgeCoords[d*2+1], edgeCoords[(d+1)*2+0], edgeCoords[(d+1)*2+1], PETSC_DRAW_BLACK);
1407: }
1408: }
1409: }
1410: break;
1411: default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot draw cells of type %s", DMPolytopeTypes[ct]);
1412: }
1413: return(0);
1414: }
1416: static PetscErrorCode DMPlexView_Draw(DM dm, PetscViewer viewer)
1417: {
1418: PetscDraw draw;
1419: DM cdm;
1420: PetscSection coordSection;
1421: Vec coordinates;
1422: const PetscScalar *coords;
1423: PetscReal xyl[2],xyr[2],bound[4] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
1424: PetscReal *refCoords, *edgeCoords;
1425: PetscBool isnull, drawAffine = PETSC_TRUE;
1426: PetscInt dim, vStart, vEnd, cStart, cEnd, c, N, edgeDiv = 4;
1427: PetscErrorCode ierr;
1430: DMGetCoordinateDim(dm, &dim);
1431: if (dim != 2) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Cannot draw meshes of dimension %D", dim);
1432: PetscOptionsGetBool(((PetscObject) dm)->options, ((PetscObject) dm)->prefix, "-dm_view_draw_affine", &drawAffine, NULL);
1433: if (!drawAffine) {PetscMalloc2((edgeDiv+1)*dim, &refCoords, (edgeDiv+1)*dim, &edgeCoords);}
1434: DMGetCoordinateDM(dm, &cdm);
1435: DMGetLocalSection(cdm, &coordSection);
1436: DMGetCoordinatesLocal(dm, &coordinates);
1437: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1438: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
1440: PetscViewerDrawGetDraw(viewer, 0, &draw);
1441: PetscDrawIsNull(draw, &isnull);
1442: if (isnull) return(0);
1443: PetscDrawSetTitle(draw, "Mesh");
1445: VecGetLocalSize(coordinates, &N);
1446: VecGetArrayRead(coordinates, &coords);
1447: for (c = 0; c < N; c += dim) {
1448: bound[0] = PetscMin(bound[0], PetscRealPart(coords[c])); bound[2] = PetscMax(bound[2], PetscRealPart(coords[c]));
1449: bound[1] = PetscMin(bound[1], PetscRealPart(coords[c+1])); bound[3] = PetscMax(bound[3], PetscRealPart(coords[c+1]));
1450: }
1451: VecRestoreArrayRead(coordinates, &coords);
1452: MPIU_Allreduce(&bound[0],xyl,2,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)dm));
1453: MPIU_Allreduce(&bound[2],xyr,2,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)dm));
1454: PetscDrawSetCoordinates(draw, xyl[0], xyl[1], xyr[0], xyr[1]);
1455: PetscDrawClear(draw);
1457: for (c = cStart; c < cEnd; ++c) {
1458: PetscScalar *coords = NULL;
1459: PetscInt numCoords;
1461: DMPlexVecGetClosureAtDepth_Internal(dm, coordSection, coordinates, c, 0, &numCoords, &coords);
1462: if (drawAffine) {
1463: DMPlexDrawCell(dm, draw, c, coords);
1464: } else {
1465: DMPlexDrawCellHighOrder(dm, draw, c, coords, edgeDiv, refCoords, edgeCoords);
1466: }
1467: DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, &numCoords, &coords);
1468: }
1469: if (!drawAffine) {PetscFree2(refCoords, edgeCoords);}
1470: PetscDrawFlush(draw);
1471: PetscDrawPause(draw);
1472: PetscDrawSave(draw);
1473: return(0);
1474: }
1476: #if defined(PETSC_HAVE_EXODUSII)
1477: #include <exodusII.h>
1478: #include <petscviewerexodusii.h>
1479: #endif
1481: PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer)
1482: {
1483: PetscBool iascii, ishdf5, isvtk, isdraw, flg, isglvis, isexodus;
1484: char name[PETSC_MAX_PATH_LEN];
1490: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
1491: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);
1492: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
1493: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW, &isdraw);
1494: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERGLVIS, &isglvis);
1495: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWEREXODUSII, &isexodus);
1496: if (iascii) {
1497: PetscViewerFormat format;
1498: PetscViewerGetFormat(viewer, &format);
1499: if (format == PETSC_VIEWER_ASCII_GLVIS) {
1500: DMPlexView_GLVis(dm, viewer);
1501: } else {
1502: DMPlexView_Ascii(dm, viewer);
1503: }
1504: } else if (ishdf5) {
1505: #if defined(PETSC_HAVE_HDF5)
1506: DMPlexView_HDF5_Internal(dm, viewer);
1507: #else
1508: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1509: #endif
1510: } else if (isvtk) {
1511: DMPlexVTKWriteAll((PetscObject) dm,viewer);
1512: } else if (isdraw) {
1513: DMPlexView_Draw(dm, viewer);
1514: } else if (isglvis) {
1515: DMPlexView_GLVis(dm, viewer);
1516: #if defined(PETSC_HAVE_EXODUSII)
1517: } else if (isexodus) {
1518: /*
1519: exodusII requires that all sets be part of exactly one cell set.
1520: If the dm does not have a "Cell Sets" label defined, we create one
1521: with ID 1, containig all cells.
1522: Note that if the Cell Sets label is defined but does not cover all cells,
1523: we may still have a problem. This should probably be checked here or in the viewer;
1524: */
1525: PetscInt numCS;
1526: DMGetLabelSize(dm,"Cell Sets",&numCS);
1527: if (!numCS) {
1528: PetscInt cStart, cEnd, c;
1529: DMCreateLabel(dm, "Cell Sets");
1530: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
1531: for (c = cStart; c < cEnd; ++c) {DMSetLabelValue(dm, "Cell Sets", c, 1);}
1532: }
1533: DMView_PlexExodusII(dm, viewer);
1534: #endif
1535: } else {
1536: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMPlex writing", ((PetscObject)viewer)->type_name);
1537: }
1538: /* Optionally view the partition */
1539: PetscOptionsHasName(((PetscObject) dm)->options, ((PetscObject) dm)->prefix, "-dm_partition_view", &flg);
1540: if (flg) {
1541: Vec ranks;
1542: DMPlexCreateRankField(dm, &ranks);
1543: VecView(ranks, viewer);
1544: VecDestroy(&ranks);
1545: }
1546: /* Optionally view a label */
1547: PetscOptionsGetString(((PetscObject) dm)->options, ((PetscObject) dm)->prefix, "-dm_label_view", name, sizeof(name), &flg);
1548: if (flg) {
1549: DMLabel label;
1550: Vec val;
1552: DMGetLabel(dm, name, &label);
1553: if (!label) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Label %s provided to -dm_label_view does not exist in this DM", name);
1554: DMPlexCreateLabelField(dm, label, &val);
1555: VecView(val, viewer);
1556: VecDestroy(&val);
1557: }
1558: return(0);
1559: }
1561: /*@
1562: DMPlexTopologyView - Saves a DMPlex topology into a file
1564: Collective on DM
1566: Input Parameters:
1567: + dm - The DM whose topology is to be saved
1568: - viewer - The PetscViewer for saving
1570: Level: advanced
1572: .seealso: DMView(), DMPlexCoordinatesView(), DMPlexLabelsView(), DMPlexTopologyLoad()
1573: @*/
1574: PetscErrorCode DMPlexTopologyView(DM dm, PetscViewer viewer)
1575: {
1576: PetscBool ishdf5;
1582: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
1583: if (ishdf5) {
1584: #if defined(PETSC_HAVE_HDF5)
1585: PetscViewerFormat format;
1586: PetscViewerGetFormat(viewer, &format);
1587: if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
1588: IS globalPointNumbering;
1590: DMPlexCreatePointNumbering(dm, &globalPointNumbering);
1591: DMPlexTopologyView_HDF5_Internal(dm, globalPointNumbering, viewer);
1592: ISDestroy(&globalPointNumbering);
1593: } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 output.", PetscViewerFormats[format]);
1594: #else
1595: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1596: #endif
1597: }
1598: return(0);
1599: }
1601: /*@
1602: DMPlexCoordinatesView - Saves DMPlex coordinates into a file
1604: Collective on DM
1606: Input Parameters:
1607: + dm - The DM whose coordinates are to be saved
1608: - viewer - The PetscViewer for saving
1610: Level: advanced
1612: .seealso: DMView(), DMPlexTopologyView(), DMPlexLabelsView(), DMPlexCoordinatesLoad()
1613: @*/
1614: PetscErrorCode DMPlexCoordinatesView(DM dm, PetscViewer viewer)
1615: {
1616: PetscBool ishdf5;
1622: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
1623: if (ishdf5) {
1624: #if defined(PETSC_HAVE_HDF5)
1625: PetscViewerFormat format;
1626: PetscViewerGetFormat(viewer, &format);
1627: if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
1628: DMPlexCoordinatesView_HDF5_Internal(dm, viewer);
1629: } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
1630: #else
1631: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1632: #endif
1633: }
1634: return(0);
1635: }
1637: /*@
1638: DMPlexLabelsView - Saves DMPlex labels into a file
1640: Collective on DM
1642: Input Parameters:
1643: + dm - The DM whose labels are to be saved
1644: - viewer - The PetscViewer for saving
1646: Level: advanced
1648: .seealso: DMView(), DMPlexTopologyView(), DMPlexCoordinatesView(), DMPlexLabelsLoad()
1649: @*/
1650: PetscErrorCode DMPlexLabelsView(DM dm, PetscViewer viewer)
1651: {
1652: PetscBool ishdf5;
1658: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
1659: if (ishdf5) {
1660: #if defined(PETSC_HAVE_HDF5)
1661: IS globalPointNumbering;
1662: PetscViewerFormat format;
1664: PetscViewerGetFormat(viewer, &format);
1665: if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
1666: DMPlexCreatePointNumbering(dm, &globalPointNumbering);
1667: DMPlexLabelsView_HDF5_Internal(dm, globalPointNumbering, viewer);
1668: ISDestroy(&globalPointNumbering);
1669: } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
1670: #else
1671: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1672: #endif
1673: }
1674: return(0);
1675: }
1677: /*@
1678: DMPlexSectionView - Saves a section associated with a DMPlex
1680: Collective on DM
1682: Input Parameters:
1683: + dm - The DM that contains the topology on which the section to be saved is defined
1684: . viewer - The PetscViewer for saving
1685: - sectiondm - The DM that contains the section to be saved
1687: Level: advanced
1689: Notes:
1690: This function is a wrapper around PetscSectionView(); in addition to the raw section, it saves information that associates the section points to the topology (dm) points. When the topology (dm) and the section are later loaded with DMPlexTopologyLoad() and DMPlexSectionLoad(), respectively, this information is used to match section points with topology points.
1692: In general dm and sectiondm are two different objects, the former carrying the topology and the latter carrying the section, and have been given a topology name and a section name, respectively, with PetscObjectSetName(). In practice, however, they can be the same object if it carries both topology and section; in that case the name of the object is used as both the topology name and the section name.
1694: .seealso: DMView(), DMPlexTopologyView(), DMPlexCoordinatesView(), DMPlexLabelsView(), DMPlexGlobalVectorView(), DMPlexLocalVectorView(), PetscSectionView(), DMPlexSectionLoad()
1695: @*/
1696: PetscErrorCode DMPlexSectionView(DM dm, PetscViewer viewer, DM sectiondm)
1697: {
1698: PetscBool ishdf5;
1705: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
1706: if (ishdf5) {
1707: #if defined(PETSC_HAVE_HDF5)
1708: DMPlexSectionView_HDF5_Internal(dm, viewer, sectiondm);
1709: #else
1710: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1711: #endif
1712: }
1713: return(0);
1714: }
1716: /*@
1717: DMPlexGlobalVectorView - Saves a global vector
1719: Collective on DM
1721: Input Parameters:
1722: + dm - The DM that represents the topology
1723: . viewer - The PetscViewer to save data with
1724: . sectiondm - The DM that contains the global section on which vec is defined
1725: - vec - The global vector to be saved
1727: Level: advanced
1729: Notes:
1730: In general dm and sectiondm are two different objects, the former carrying the topology and the latter carrying the section, and have been given a topology name and a section name, respectively, with PetscObjectSetName(). In practice, however, they can be the same object if it carries both topology and section; in that case the name of the object is used as both the topology name and the section name.
1732: Typical calling sequence
1733: $ DMCreate(PETSC_COMM_WORLD, &dm);
1734: $ DMSetType(dm, DMPLEX);
1735: $ PetscObjectSetName((PetscObject)dm, "topologydm_name");
1736: $ DMClone(dm, §iondm);
1737: $ PetscObjectSetName((PetscObject)sectiondm, "sectiondm_name");
1738: $ PetscSectionCreate(PETSC_COMM_WORLD, §ion);
1739: $ DMPlexGetChart(sectiondm, &pStart, &pEnd);
1740: $ PetscSectionSetChart(section, pStart, pEnd);
1741: $ PetscSectionSetUp(section);
1742: $ DMSetLocalSection(sectiondm, section);
1743: $ PetscSectionDestroy(§ion);
1744: $ DMGetGlobalVector(sectiondm, &vec);
1745: $ PetscObjectSetName((PetscObject)vec, "vec_name");
1746: $ DMPlexTopologyView(dm, viewer);
1747: $ DMPlexSectionView(dm, viewer, sectiondm);
1748: $ DMPlexGlobalVectorView(dm, viewer, sectiondm, vec);
1749: $ DMRestoreGlobalVector(sectiondm, &vec);
1750: $ DMDestroy(§iondm);
1751: $ DMDestroy(&dm);
1753: .seealso: DMPlexTopologyView(), DMPlexSectionView(), DMPlexLocalVectorView(), DMPlexGlobalVectorLoad(), DMPlexLocalVectorLoad()
1754: @*/
1755: PetscErrorCode DMPlexGlobalVectorView(DM dm, PetscViewer viewer, DM sectiondm, Vec vec)
1756: {
1757: PetscBool ishdf5;
1758: PetscErrorCode ierr;
1765: /* Check consistency */
1766: {
1767: PetscSection section;
1768: PetscBool includesConstraints;
1769: PetscInt m, m1;
1771: VecGetLocalSize(vec, &m1);
1772: DMGetGlobalSection(sectiondm, §ion);
1773: PetscSectionGetIncludesConstraints(section, &includesConstraints);
1774: if (includesConstraints) {PetscSectionGetStorageSize(section, &m);}
1775: else {PetscSectionGetConstrainedStorageSize(section, &m);}
1776: if (m1 != m) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Global vector size (%D) != global section storage size (%D)", m1, m);
1777: }
1778: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5);
1779: if (ishdf5) {
1780: #if defined(PETSC_HAVE_HDF5)
1781: DMPlexGlobalVectorView_HDF5_Internal(dm, viewer, sectiondm, vec);
1782: #else
1783: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1784: #endif
1785: }
1786: return(0);
1787: }
1789: /*@
1790: DMPlexLocalVectorView - Saves a local vector
1792: Collective on DM
1794: Input Parameters:
1795: + dm - The DM that represents the topology
1796: . viewer - The PetscViewer to save data with
1797: . sectiondm - The DM that contains the local section on which vec is defined; may be the same as dm
1798: - vec - The local vector to be saved
1800: Level: advanced
1802: Notes:
1803: In general dm and sectiondm are two different objects, the former carrying the topology and the latter carrying the section, and have been given a topology name and a section name, respectively, with PetscObjectSetName(). In practice, however, they can be the same object if it carries both topology and section; in that case the name of the object is used as both the topology name and the section name.
1805: Typical calling sequence
1806: $ DMCreate(PETSC_COMM_WORLD, &dm);
1807: $ DMSetType(dm, DMPLEX);
1808: $ PetscObjectSetName((PetscObject)dm, "topologydm_name");
1809: $ DMClone(dm, §iondm);
1810: $ PetscObjectSetName((PetscObject)sectiondm, "sectiondm_name");
1811: $ PetscSectionCreate(PETSC_COMM_WORLD, §ion);
1812: $ DMPlexGetChart(sectiondm, &pStart, &pEnd);
1813: $ PetscSectionSetChart(section, pStart, pEnd);
1814: $ PetscSectionSetUp(section);
1815: $ DMSetLocalSection(sectiondm, section);
1816: $ DMGetLocalVector(sectiondm, &vec);
1817: $ PetscObjectSetName((PetscObject)vec, "vec_name");
1818: $ DMPlexTopologyView(dm, viewer);
1819: $ DMPlexSectionView(dm, viewer, sectiondm);
1820: $ DMPlexLocalVectorView(dm, viewer, sectiondm, vec);
1821: $ DMRestoreLocalVector(sectiondm, &vec);
1822: $ DMDestroy(§iondm);
1823: $ DMDestroy(&dm);
1825: .seealso: DMPlexTopologyView(), DMPlexSectionView(), DMPlexGlobalVectorView(), DMPlexGlobalVectorLoad(), DMPlexLocalVectorLoad()
1826: @*/
1827: PetscErrorCode DMPlexLocalVectorView(DM dm, PetscViewer viewer, DM sectiondm, Vec vec)
1828: {
1829: PetscBool ishdf5;
1830: PetscErrorCode ierr;
1837: /* Check consistency */
1838: {
1839: PetscSection section;
1840: PetscBool includesConstraints;
1841: PetscInt m, m1;
1843: VecGetLocalSize(vec, &m1);
1844: DMGetLocalSection(sectiondm, §ion);
1845: PetscSectionGetIncludesConstraints(section, &includesConstraints);
1846: if (includesConstraints) {PetscSectionGetStorageSize(section, &m);}
1847: else {PetscSectionGetConstrainedStorageSize(section, &m);}
1848: if (m1 != m) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Local vector size (%D) != local section storage size (%D)", m1, m);
1849: }
1850: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5);
1851: if (ishdf5) {
1852: #if defined(PETSC_HAVE_HDF5)
1853: DMPlexLocalVectorView_HDF5_Internal(dm, viewer, sectiondm, vec);
1854: #else
1855: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1856: #endif
1857: }
1858: return(0);
1859: }
1861: PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer)
1862: {
1863: PetscBool ishdf5;
1869: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
1870: if (ishdf5) {
1871: #if defined(PETSC_HAVE_HDF5)
1872: PetscViewerFormat format;
1873: PetscViewerGetFormat(viewer, &format);
1874: if (format == PETSC_VIEWER_HDF5_XDMF || format == PETSC_VIEWER_HDF5_VIZ) {
1875: DMPlexLoad_HDF5_Xdmf_Internal(dm, viewer);
1876: } else if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
1877: DMPlexLoad_HDF5_Internal(dm, viewer);
1878: } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
1879: return(0);
1880: #else
1881: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1882: #endif
1883: } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMPlex loading", ((PetscObject)viewer)->type_name);
1884: }
1886: /*@
1887: DMPlexTopologyLoad - Loads a topology into a DMPlex
1889: Collective on DM
1891: Input Parameters:
1892: + dm - The DM into which the topology is loaded
1893: - viewer - The PetscViewer for the saved topology
1895: Output Parameters:
1896: . globalToLocalPointSF - The PetscSF that pushes points in [0, N) to the associated points in the loaded plex, where N is the global number of points; NULL if unneeded
1898: Level: advanced
1900: .seealso: DMLoad(), DMPlexCoordinatesLoad(), DMPlexLabelsLoad(), DMView(), PetscViewerHDF5Open(), PetscViewerPushFormat()
1901: @*/
1902: PetscErrorCode DMPlexTopologyLoad(DM dm, PetscViewer viewer, PetscSF *globalToLocalPointSF)
1903: {
1904: PetscBool ishdf5;
1911: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
1912: if (ishdf5) {
1913: #if defined(PETSC_HAVE_HDF5)
1914: PetscViewerFormat format;
1915: PetscViewerGetFormat(viewer, &format);
1916: if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
1917: DMPlexTopologyLoad_HDF5_Internal(dm, viewer, globalToLocalPointSF);
1918: } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
1919: #else
1920: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1921: #endif
1922: }
1923: return(0);
1924: }
1926: /*@
1927: DMPlexCoordinatesLoad - Loads coordinates into a DMPlex
1929: Collective on DM
1931: Input Parameters:
1932: + dm - The DM into which the coordinates are loaded
1933: - viewer - The PetscViewer for the saved coordinates
1935: Level: advanced
1937: .seealso: DMLoad(), DMPlexTopologyLoad(), DMPlexLabelsLoad(), DMView(), PetscViewerHDF5Open(), PetscViewerPushFormat()
1938: @*/
1939: PetscErrorCode DMPlexCoordinatesLoad(DM dm, PetscViewer viewer)
1940: {
1941: PetscBool ishdf5;
1947: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
1948: if (ishdf5) {
1949: #if defined(PETSC_HAVE_HDF5)
1950: PetscViewerFormat format;
1951: PetscViewerGetFormat(viewer, &format);
1952: if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
1953: DMPlexCoordinatesLoad_HDF5_Internal(dm, viewer);
1954: } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
1955: #else
1956: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1957: #endif
1958: }
1959: return(0);
1960: }
1962: /*@
1963: DMPlexLabelsLoad - Loads labels into a DMPlex
1965: Collective on DM
1967: Input Parameters:
1968: + dm - The DM into which the labels are loaded
1969: - viewer - The PetscViewer for the saved labels
1971: Level: advanced
1973: .seealso: DMLoad(), DMPlexTopologyLoad(), DMPlexCoordinatesLoad(), DMView(), PetscViewerHDF5Open(), PetscViewerPushFormat()
1974: @*/
1975: PetscErrorCode DMPlexLabelsLoad(DM dm, PetscViewer viewer)
1976: {
1977: PetscBool ishdf5;
1983: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
1984: if (ishdf5) {
1985: #if defined(PETSC_HAVE_HDF5)
1986: PetscViewerFormat format;
1988: PetscViewerGetFormat(viewer, &format);
1989: if (format == PETSC_VIEWER_HDF5_PETSC || format == PETSC_VIEWER_DEFAULT || format == PETSC_VIEWER_NATIVE) {
1990: DMPlexLabelsLoad_HDF5_Internal(dm, viewer);
1991: } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "PetscViewerFormat %s not supported for HDF5 input.", PetscViewerFormats[format]);
1992: #else
1993: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1994: #endif
1995: }
1996: return(0);
1997: }
1999: /*@
2000: DMPlexSectionLoad - Loads section into a DMPlex
2002: Collective on DM
2004: Input Parameters:
2005: + dm - The DM that represents the topology
2006: . viewer - The PetscViewer that represents the on-disk section (sectionA)
2007: . sectiondm - The DM into which the on-disk section (sectionA) is migrated
2008: - globalToLocalPointSF - The SF returned by DMPlexTopologyLoad() when loading dm from viewer
2010: Output Parameters
2011: + globalDofSF - The SF that migrates any on-disk Vec data associated with sectionA into a global Vec associated with the sectiondm's global section (NULL if not needed)
2012: - localDofSF - The SF that migrates any on-disk Vec data associated with sectionA into a local Vec associated with the sectiondm's local section (NULL if not needed)
2014: Level: advanced
2016: Notes:
2017: This function is a wrapper around PetscSectionLoad(); it loads, in addition to the raw section, a list of global point numbers that associates each on-disk section point with a global point number in [0, NX), where NX is the number of topology points in dm. Noting that globalToLocalPointSF associates each topology point in dm with a global number in [0, NX), one can readily establish an association of the on-disk section points with the topology points.
2019: In general dm and sectiondm are two different objects, the former carrying the topology and the latter carrying the section, and have been given a topology name and a section name, respectively, with PetscObjectSetName(). In practice, however, they can be the same object if it carries both topology and section; in that case the name of the object is used as both the topology name and the section name.
2021: The output parameter, globalDofSF (localDofSF), can later be used with DMPlexGlobalVectorLoad() (DMPlexLocalVectorLoad()) to load on-disk vectors into global (local) vectors associated with sectiondm's global (local) section.
2023: Example using 2 processes:
2024: $ NX (number of points on dm): 4
2025: $ sectionA : the on-disk section
2026: $ vecA : a vector associated with sectionA
2027: $ sectionB : sectiondm's local section constructed in this function
2028: $ vecB (local) : a vector associated with sectiondm's local section
2029: $ vecB (global) : a vector associated with sectiondm's global section
2030: $
2031: $ rank 0 rank 1
2032: $ vecA (global) : [.0 .4 .1 | .2 .3] <- to be loaded in DMPlexGlobalVectorLoad() or DMPlexLocalVectorLoad()
2033: $ sectionA->atlasOff : 0 2 | 1 <- loaded in PetscSectionLoad()
2034: $ sectionA->atlasDof : 1 3 | 1 <- loaded in PetscSectionLoad()
2035: $ sectionA's global point numbers: 0 2 | 3 <- loaded in DMPlexSectionLoad()
2036: $ [0, NX) : 0 1 | 2 3 <- conceptual partition used in globalToLocalPointSF
2037: $ sectionB's global point numbers: 0 1 3 | 3 2 <- associated with [0, NX) by globalToLocalPointSF
2038: $ sectionB->atlasDof : 1 0 1 | 1 3
2039: $ sectionB->atlasOff (no perm) : 0 1 1 | 0 1
2040: $ vecB (local) : [.0 .4] | [.4 .1 .2 .3] <- to be constructed by calling DMPlexLocalVectorLoad() with localDofSF
2041: $ vecB (global) : [.0 .4 | .1 .2 .3] <- to be constructed by calling DMPlexGlobalVectorLoad() with globalDofSF
2042: $
2043: $ where "|" represents a partition of loaded data, and global point 3 is assumed to be owned by rank 0.
2045: .seealso: DMLoad(), DMPlexTopologyLoad(), DMPlexCoordinatesLoad(), DMPlexLabelsLoad(), DMPlexGlobalVectorLoad(), DMPlexLocalVectorLoad(), PetscSectionLoad(), DMPlexSectionView()
2046: @*/
2047: PetscErrorCode DMPlexSectionLoad(DM dm, PetscViewer viewer, DM sectiondm, PetscSF globalToLocalPointSF, PetscSF *globalDofSF, PetscSF *localDofSF)
2048: {
2049: PetscBool ishdf5;
2059: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
2060: if (ishdf5) {
2061: #if defined(PETSC_HAVE_HDF5)
2062: DMPlexSectionLoad_HDF5_Internal(dm, viewer, sectiondm, globalToLocalPointSF, globalDofSF, localDofSF);
2063: #else
2064: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2065: #endif
2066: }
2067: return(0);
2068: }
2070: /*@
2071: DMPlexGlobalVectorLoad - Loads on-disk vector data into a global vector
2073: Collective on DM
2075: Input Parameters:
2076: + dm - The DM that represents the topology
2077: . viewer - The PetscViewer that represents the on-disk vector data
2078: . sectiondm - The DM that contains the global section on which vec is defined
2079: . sf - The SF that migrates the on-disk vector data into vec
2080: - vec - The global vector to set values of
2082: Level: advanced
2084: Notes:
2085: In general dm and sectiondm are two different objects, the former carrying the topology and the latter carrying the section, and have been given a topology name and a section name, respectively, with PetscObjectSetName(). In practice, however, they can be the same object if it carries both topology and section; in that case the name of the object is used as both the topology name and the section name.
2087: Typical calling sequence
2088: $ DMCreate(PETSC_COMM_WORLD, &dm);
2089: $ DMSetType(dm, DMPLEX);
2090: $ PetscObjectSetName((PetscObject)dm, "topologydm_name");
2091: $ DMPlexTopologyLoad(dm, viewer, &sfX);
2092: $ DMClone(dm, §iondm);
2093: $ PetscObjectSetName((PetscObject)sectiondm, "sectiondm_name");
2094: $ DMPlexSectionLoad(dm, viewer, sectiondm, sfX, &gsf, NULL);
2095: $ DMGetGlobalVector(sectiondm, &vec);
2096: $ PetscObjectSetName((PetscObject)vec, "vec_name");
2097: $ DMPlexGlobalVectorLoad(dm, viewer, sectiondm, gsf, vec);
2098: $ DMRestoreGlobalVector(sectiondm, &vec);
2099: $ PetscSFDestroy(&gsf);
2100: $ PetscSFDestroy(&sfX);
2101: $ DMDestroy(§iondm);
2102: $ DMDestroy(&dm);
2104: .seealso: DMPlexTopologyLoad(), DMPlexSectionLoad(), DMPlexLocalVectorLoad(), DMPlexGlobalVectorView(), DMPlexLocalVectorView()
2105: @*/
2106: PetscErrorCode DMPlexGlobalVectorLoad(DM dm, PetscViewer viewer, DM sectiondm, PetscSF sf, Vec vec)
2107: {
2108: PetscBool ishdf5;
2109: PetscErrorCode ierr;
2117: /* Check consistency */
2118: {
2119: PetscSection section;
2120: PetscBool includesConstraints;
2121: PetscInt m, m1;
2123: VecGetLocalSize(vec, &m1);
2124: DMGetGlobalSection(sectiondm, §ion);
2125: PetscSectionGetIncludesConstraints(section, &includesConstraints);
2126: if (includesConstraints) {PetscSectionGetStorageSize(section, &m);}
2127: else {PetscSectionGetConstrainedStorageSize(section, &m);}
2128: if (m1 != m) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Global vector size (%D) != global section storage size (%D)", m1, m);
2129: }
2130: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
2131: if (ishdf5) {
2132: #if defined(PETSC_HAVE_HDF5)
2133: DMPlexVecLoad_HDF5_Internal(dm, viewer, sectiondm, sf, vec);
2134: #else
2135: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2136: #endif
2137: }
2138: return(0);
2139: }
2141: /*@
2142: DMPlexLocalVectorLoad - Loads on-disk vector data into a local vector
2144: Collective on DM
2146: Input Parameters:
2147: + dm - The DM that represents the topology
2148: . viewer - The PetscViewer that represents the on-disk vector data
2149: . sectiondm - The DM that contains the local section on which vec is defined
2150: . sf - The SF that migrates the on-disk vector data into vec
2151: - vec - The local vector to set values of
2153: Level: advanced
2155: Notes:
2156: In general dm and sectiondm are two different objects, the former carrying the topology and the latter carrying the section, and have been given a topology name and a section name, respectively, with PetscObjectSetName(). In practice, however, they can be the same object if it carries both topology and section; in that case the name of the object is used as both the topology name and the section name.
2158: Typical calling sequence
2159: $ DMCreate(PETSC_COMM_WORLD, &dm);
2160: $ DMSetType(dm, DMPLEX);
2161: $ PetscObjectSetName((PetscObject)dm, "topologydm_name");
2162: $ DMPlexTopologyLoad(dm, viewer, &sfX);
2163: $ DMClone(dm, §iondm);
2164: $ PetscObjectSetName((PetscObject)sectiondm, "sectiondm_name");
2165: $ DMPlexSectionLoad(dm, viewer, sectiondm, sfX, NULL, &lsf);
2166: $ DMGetLocalVector(sectiondm, &vec);
2167: $ PetscObjectSetName((PetscObject)vec, "vec_name");
2168: $ DMPlexLocalVectorLoad(dm, viewer, sectiondm, lsf, vec);
2169: $ DMRestoreLocalVector(sectiondm, &vec);
2170: $ PetscSFDestroy(&lsf);
2171: $ PetscSFDestroy(&sfX);
2172: $ DMDestroy(§iondm);
2173: $ DMDestroy(&dm);
2175: .seealso: DMPlexTopologyLoad(), DMPlexSectionLoad(), DMPlexGlobalVectorLoad(), DMPlexGlobalVectorView(), DMPlexLocalVectorView()
2176: @*/
2177: PetscErrorCode DMPlexLocalVectorLoad(DM dm, PetscViewer viewer, DM sectiondm, PetscSF sf, Vec vec)
2178: {
2179: PetscBool ishdf5;
2180: PetscErrorCode ierr;
2188: /* Check consistency */
2189: {
2190: PetscSection section;
2191: PetscBool includesConstraints;
2192: PetscInt m, m1;
2194: VecGetLocalSize(vec, &m1);
2195: DMGetLocalSection(sectiondm, §ion);
2196: PetscSectionGetIncludesConstraints(section, &includesConstraints);
2197: if (includesConstraints) {PetscSectionGetStorageSize(section, &m);}
2198: else {PetscSectionGetConstrainedStorageSize(section, &m);}
2199: if (m1 != m) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Local vector size (%D) != local section storage size (%D)", m1, m);
2200: }
2201: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
2202: if (ishdf5) {
2203: #if defined(PETSC_HAVE_HDF5)
2204: DMPlexVecLoad_HDF5_Internal(dm, viewer, sectiondm, sf, vec);
2205: #else
2206: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
2207: #endif
2208: }
2209: return(0);
2210: }
2212: PetscErrorCode DMDestroy_Plex(DM dm)
2213: {
2214: DM_Plex *mesh = (DM_Plex*) dm->data;
2218: PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",NULL);
2219: PetscObjectComposeFunction((PetscObject)dm,"DMPlexInsertBoundaryValues_C", NULL);
2220: PetscObjectComposeFunction((PetscObject)dm,"DMCreateNeumannOverlap_C", NULL);
2221: PetscObjectComposeFunction((PetscObject)dm,"DMInterpolateSolution_C", NULL);
2222: if (--mesh->refct > 0) return(0);
2223: PetscSectionDestroy(&mesh->coneSection);
2224: PetscFree(mesh->cones);
2225: PetscFree(mesh->coneOrientations);
2226: PetscSectionDestroy(&mesh->supportSection);
2227: PetscSectionDestroy(&mesh->subdomainSection);
2228: PetscFree(mesh->supports);
2229: PetscFree(mesh->facesTmp);
2230: PetscFree(mesh->tetgenOpts);
2231: PetscFree(mesh->triangleOpts);
2232: PetscFree(mesh->transformType);
2233: PetscPartitionerDestroy(&mesh->partitioner);
2234: DMLabelDestroy(&mesh->subpointMap);
2235: ISDestroy(&mesh->subpointIS);
2236: ISDestroy(&mesh->globalVertexNumbers);
2237: ISDestroy(&mesh->globalCellNumbers);
2238: PetscSectionDestroy(&mesh->anchorSection);
2239: ISDestroy(&mesh->anchorIS);
2240: PetscSectionDestroy(&mesh->parentSection);
2241: PetscFree(mesh->parents);
2242: PetscFree(mesh->childIDs);
2243: PetscSectionDestroy(&mesh->childSection);
2244: PetscFree(mesh->children);
2245: DMDestroy(&mesh->referenceTree);
2246: PetscGridHashDestroy(&mesh->lbox);
2247: PetscFree(mesh->neighbors);
2248: /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
2249: PetscFree(mesh);
2250: return(0);
2251: }
2253: PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J)
2254: {
2255: PetscSection sectionGlobal;
2256: PetscInt bs = -1, mbs;
2257: PetscInt localSize;
2258: PetscBool isShell, isBlock, isSeqBlock, isMPIBlock, isSymBlock, isSymSeqBlock, isSymMPIBlock, isMatIS;
2259: PetscErrorCode ierr;
2260: MatType mtype;
2261: ISLocalToGlobalMapping ltog;
2264: MatInitializePackage();
2265: mtype = dm->mattype;
2266: DMGetGlobalSection(dm, §ionGlobal);
2267: /* PetscSectionGetStorageSize(sectionGlobal, &localSize); */
2268: PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize);
2269: MatCreate(PetscObjectComm((PetscObject)dm), J);
2270: MatSetSizes(*J, localSize, localSize, PETSC_DETERMINE, PETSC_DETERMINE);
2271: MatSetType(*J, mtype);
2272: MatSetFromOptions(*J);
2273: MatGetBlockSize(*J, &mbs);
2274: if (mbs > 1) bs = mbs;
2275: PetscStrcmp(mtype, MATSHELL, &isShell);
2276: PetscStrcmp(mtype, MATBAIJ, &isBlock);
2277: PetscStrcmp(mtype, MATSEQBAIJ, &isSeqBlock);
2278: PetscStrcmp(mtype, MATMPIBAIJ, &isMPIBlock);
2279: PetscStrcmp(mtype, MATSBAIJ, &isSymBlock);
2280: PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock);
2281: PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock);
2282: PetscStrcmp(mtype, MATIS, &isMatIS);
2283: if (!isShell) {
2284: PetscSection subSection;
2285: PetscBool fillMatrix = (PetscBool)(!dm->prealloc_only && !isMatIS);
2286: PetscInt *dnz, *onz, *dnzu, *onzu, bsLocal[2], bsMinMax[2], *ltogidx, lsize;
2287: PetscInt pStart, pEnd, p, dof, cdof;
2289: /* Set localtoglobalmapping on the matrix for MatSetValuesLocal() to work (it also creates the local matrices in case of MATIS) */
2290: if (isMatIS) { /* need a different l2g map than the one computed by DMGetLocalToGlobalMapping */
2291: PetscSection section;
2292: PetscInt size;
2294: DMGetLocalSection(dm, §ion);
2295: PetscSectionGetStorageSize(section, &size);
2296: PetscMalloc1(size,<ogidx);
2297: DMPlexGetSubdomainSection(dm, &subSection);
2298: } else {
2299: DMGetLocalToGlobalMapping(dm,<og);
2300: }
2301: PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);
2302: for (p = pStart, lsize = 0; p < pEnd; ++p) {
2303: PetscInt bdof;
2305: PetscSectionGetDof(sectionGlobal, p, &dof);
2306: PetscSectionGetConstraintDof(sectionGlobal, p, &cdof);
2307: dof = dof < 0 ? -(dof+1) : dof;
2308: bdof = cdof && (dof-cdof) ? 1 : dof;
2309: if (dof) {
2310: if (bs < 0) {bs = bdof;}
2311: else if (bs != bdof) {bs = 1; if (!isMatIS) break;}
2312: }
2313: if (isMatIS) {
2314: PetscInt loff,c,off;
2315: PetscSectionGetOffset(subSection, p, &loff);
2316: PetscSectionGetOffset(sectionGlobal, p, &off);
2317: for (c = 0; c < dof-cdof; ++c, ++lsize) ltogidx[loff+c] = off > -1 ? off+c : -(off+1)+c;
2318: }
2319: }
2320: /* Must have same blocksize on all procs (some might have no points) */
2321: bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; bsLocal[1] = bs;
2322: PetscGlobalMinMaxInt(PetscObjectComm((PetscObject) dm), bsLocal, bsMinMax);
2323: if (bsMinMax[0] != bsMinMax[1]) {bs = 1;}
2324: else {bs = bsMinMax[0];}
2325: bs = PetscMax(1,bs);
2326: if (isMatIS) { /* Must reduce indices by blocksize */
2327: PetscInt l;
2329: lsize = lsize/bs;
2330: if (bs > 1) for (l = 0; l < lsize; ++l) ltogidx[l] = ltogidx[l*bs]/bs;
2331: ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), bs, lsize, ltogidx, PETSC_OWN_POINTER, <og);
2332: }
2333: MatSetLocalToGlobalMapping(*J,ltog,ltog);
2334: if (isMatIS) {
2335: ISLocalToGlobalMappingDestroy(<og);
2336: }
2337: PetscCalloc4(localSize/bs, &dnz, localSize/bs, &onz, localSize/bs, &dnzu, localSize/bs, &onzu);
2338: DMPlexPreallocateOperator(dm, bs, dnz, onz, dnzu, onzu, *J, fillMatrix);
2339: PetscFree4(dnz, onz, dnzu, onzu);
2340: }
2341: MatSetDM(*J, dm);
2342: return(0);
2343: }
2345: /*@
2346: DMPlexGetSubdomainSection - Returns the section associated with the subdomain
2348: Not collective
2350: Input Parameter:
2351: . mesh - The DMPlex
2353: Output Parameters:
2354: . subsection - The subdomain section
2356: Level: developer
2358: .seealso:
2359: @*/
2360: PetscErrorCode DMPlexGetSubdomainSection(DM dm, PetscSection *subsection)
2361: {
2362: DM_Plex *mesh = (DM_Plex*) dm->data;
2367: if (!mesh->subdomainSection) {
2368: PetscSection section;
2369: PetscSF sf;
2371: PetscSFCreate(PETSC_COMM_SELF,&sf);
2372: DMGetLocalSection(dm,§ion);
2373: PetscSectionCreateGlobalSection(section,sf,PETSC_FALSE,PETSC_TRUE,&mesh->subdomainSection);
2374: PetscSFDestroy(&sf);
2375: }
2376: *subsection = mesh->subdomainSection;
2377: return(0);
2378: }
2380: /*@
2381: DMPlexGetChart - Return the interval for all mesh points [pStart, pEnd)
2383: Not collective
2385: Input Parameter:
2386: . mesh - The DMPlex
2388: Output Parameters:
2389: + pStart - The first mesh point
2390: - pEnd - The upper bound for mesh points
2392: Level: beginner
2394: .seealso: DMPlexCreate(), DMPlexSetChart()
2395: @*/
2396: PetscErrorCode DMPlexGetChart(DM dm, PetscInt *pStart, PetscInt *pEnd)
2397: {
2398: DM_Plex *mesh = (DM_Plex*) dm->data;
2403: PetscSectionGetChart(mesh->coneSection, pStart, pEnd);
2404: return(0);
2405: }
2407: /*@
2408: DMPlexSetChart - Set the interval for all mesh points [pStart, pEnd)
2410: Not collective
2412: Input Parameters:
2413: + mesh - The DMPlex
2414: . pStart - The first mesh point
2415: - pEnd - The upper bound for mesh points
2417: Output Parameters:
2419: Level: beginner
2421: .seealso: DMPlexCreate(), DMPlexGetChart()
2422: @*/
2423: PetscErrorCode DMPlexSetChart(DM dm, PetscInt pStart, PetscInt pEnd)
2424: {
2425: DM_Plex *mesh = (DM_Plex*) dm->data;
2430: PetscSectionSetChart(mesh->coneSection, pStart, pEnd);
2431: PetscSectionSetChart(mesh->supportSection, pStart, pEnd);
2432: return(0);
2433: }
2435: /*@
2436: DMPlexGetConeSize - Return the number of in-edges for this point in the DAG
2438: Not collective
2440: Input Parameters:
2441: + mesh - The DMPlex
2442: - p - The point, which must lie in the chart set with DMPlexSetChart()
2444: Output Parameter:
2445: . size - The cone size for point p
2447: Level: beginner
2449: .seealso: DMPlexCreate(), DMPlexSetConeSize(), DMPlexSetChart()
2450: @*/
2451: PetscErrorCode DMPlexGetConeSize(DM dm, PetscInt p, PetscInt *size)
2452: {
2453: DM_Plex *mesh = (DM_Plex*) dm->data;
2459: PetscSectionGetDof(mesh->coneSection, p, size);
2460: return(0);
2461: }
2463: /*@
2464: DMPlexSetConeSize - Set the number of in-edges for this point in the DAG
2466: Not collective
2468: Input Parameters:
2469: + mesh - The DMPlex
2470: . p - The point, which must lie in the chart set with DMPlexSetChart()
2471: - size - The cone size for point p
2473: Output Parameter:
2475: Note:
2476: This should be called after DMPlexSetChart().
2478: Level: beginner
2480: .seealso: DMPlexCreate(), DMPlexGetConeSize(), DMPlexSetChart()
2481: @*/
2482: PetscErrorCode DMPlexSetConeSize(DM dm, PetscInt p, PetscInt size)
2483: {
2484: DM_Plex *mesh = (DM_Plex*) dm->data;
2489: PetscSectionSetDof(mesh->coneSection, p, size);
2491: mesh->maxConeSize = PetscMax(mesh->maxConeSize, size);
2492: return(0);
2493: }
2495: /*@
2496: DMPlexAddConeSize - Add the given number of in-edges to this point in the DAG
2498: Not collective
2500: Input Parameters:
2501: + mesh - The DMPlex
2502: . p - The point, which must lie in the chart set with DMPlexSetChart()
2503: - size - The additional cone size for point p
2505: Output Parameter:
2507: Note:
2508: This should be called after DMPlexSetChart().
2510: Level: beginner
2512: .seealso: DMPlexCreate(), DMPlexSetConeSize(), DMPlexGetConeSize(), DMPlexSetChart()
2513: @*/
2514: PetscErrorCode DMPlexAddConeSize(DM dm, PetscInt p, PetscInt size)
2515: {
2516: DM_Plex *mesh = (DM_Plex*) dm->data;
2517: PetscInt csize;
2522: PetscSectionAddDof(mesh->coneSection, p, size);
2523: PetscSectionGetDof(mesh->coneSection, p, &csize);
2525: mesh->maxConeSize = PetscMax(mesh->maxConeSize, csize);
2526: return(0);
2527: }
2529: /*@C
2530: DMPlexGetCone - Return the points on the in-edges for this point in the DAG
2532: Not collective
2534: Input Parameters:
2535: + dm - The DMPlex
2536: - p - The point, which must lie in the chart set with DMPlexSetChart()
2538: Output Parameter:
2539: . cone - An array of points which are on the in-edges for point p
2541: Level: beginner
2543: Fortran Notes:
2544: Since it returns an array, this routine is only available in Fortran 90, and you must
2545: include petsc.h90 in your code.
2546: You must also call DMPlexRestoreCone() after you finish using the returned array.
2547: DMPlexRestoreCone() is not needed/available in C.
2549: .seealso: DMPlexGetConeSize(), DMPlexSetCone(), DMPlexGetConeTuple(), DMPlexSetChart()
2550: @*/
2551: PetscErrorCode DMPlexGetCone(DM dm, PetscInt p, const PetscInt *cone[])
2552: {
2553: DM_Plex *mesh = (DM_Plex*) dm->data;
2554: PetscInt off;
2560: PetscSectionGetOffset(mesh->coneSection, p, &off);
2561: *cone = &mesh->cones[off];
2562: return(0);
2563: }
2565: /*@C
2566: DMPlexGetConeTuple - Return the points on the in-edges of several points in the DAG
2568: Not collective
2570: Input Parameters:
2571: + dm - The DMPlex
2572: - p - The IS of points, which must lie in the chart set with DMPlexSetChart()
2574: Output Parameters:
2575: + pConesSection - PetscSection describing the layout of pCones
2576: - pCones - An array of points which are on the in-edges for the point set p
2578: Level: intermediate
2580: .seealso: DMPlexCreate(), DMPlexGetCone(), DMPlexGetConeRecursive(), DMPlexSetChart()
2581: @*/
2582: PetscErrorCode DMPlexGetConeTuple(DM dm, IS p, PetscSection *pConesSection, IS *pCones)
2583: {
2584: PetscSection cs, newcs;
2585: PetscInt *cones;
2586: PetscInt *newarr=NULL;
2587: PetscInt n;
2588: PetscErrorCode ierr;
2591: DMPlexGetCones(dm, &cones);
2592: DMPlexGetConeSection(dm, &cs);
2593: PetscSectionExtractDofsFromArray(cs, MPIU_INT, cones, p, &newcs, pCones ? ((void**)&newarr) : NULL);
2594: if (pConesSection) *pConesSection = newcs;
2595: if (pCones) {
2596: PetscSectionGetStorageSize(newcs, &n);
2597: ISCreateGeneral(PetscObjectComm((PetscObject)p), n, newarr, PETSC_OWN_POINTER, pCones);
2598: }
2599: return(0);
2600: }
2602: /*@
2603: DMPlexGetConeRecursiveVertices - Expand each given point into its cone points and do that recursively until we end up just with vertices.
2605: Not collective
2607: Input Parameters:
2608: + dm - The DMPlex
2609: - points - The IS of points, which must lie in the chart set with DMPlexSetChart()
2611: Output Parameter:
2612: . expandedPoints - An array of vertices recursively expanded from input points
2614: Level: advanced
2616: Notes:
2617: Like DMPlexGetConeRecursive but returns only the 0-depth IS (i.e. vertices only) and no sections.
2618: There is no corresponding Restore function, just call ISDestroy() on the returned IS to deallocate.
2620: .seealso: DMPlexCreate(), DMPlexGetCone(), DMPlexGetConeTuple(), DMPlexGetConeRecursive(), DMPlexRestoreConeRecursive(), DMPlexGetDepth()
2621: @*/
2622: PetscErrorCode DMPlexGetConeRecursiveVertices(DM dm, IS points, IS *expandedPoints)
2623: {
2624: IS *expandedPointsAll;
2625: PetscInt depth;
2626: PetscErrorCode ierr;
2632: DMPlexGetConeRecursive(dm, points, &depth, &expandedPointsAll, NULL);
2633: *expandedPoints = expandedPointsAll[0];
2634: PetscObjectReference((PetscObject)expandedPointsAll[0]);
2635: DMPlexRestoreConeRecursive(dm, points, &depth, &expandedPointsAll, NULL);
2636: return(0);
2637: }
2639: /*@
2640: DMPlexGetConeRecursive - Expand each given point into its cone points and do that recursively until we end up just with vertices (DAG points of depth 0, i.e. without cones).
2642: Not collective
2644: Input Parameters:
2645: + dm - The DMPlex
2646: - points - The IS of points, which must lie in the chart set with DMPlexSetChart()
2648: Output Parameters:
2649: + depth - (optional) Size of the output arrays, equal to DMPlex depth, returned by DMPlexGetDepth()
2650: . expandedPoints - (optional) An array of index sets with recursively expanded cones
2651: - sections - (optional) An array of sections which describe mappings from points to their cone points
2653: Level: advanced
2655: Notes:
2656: Like DMPlexGetConeTuple() but recursive.
2658: Array expandedPoints has size equal to depth. Each expandedPoints[d] contains DAG points with maximum depth d, recursively cone-wise expanded from the input points.
2659: For example, for d=0 it contains only vertices, for d=1 it can contain vertices and edges, etc.
2661: Array section has size equal to depth. Each PetscSection sections[d] realizes mapping from expandedPoints[d+1] (section points) to expandedPoints[d] (section dofs) as follows:
2662: (1) DAG points in expandedPoints[d+1] with depth d+1 to their cone points in expandedPoints[d];
2663: (2) DAG points in expandedPoints[d+1] with depth in [0,d] to the same points in expandedPoints[d].
2665: .seealso: DMPlexCreate(), DMPlexGetCone(), DMPlexGetConeTuple(), DMPlexRestoreConeRecursive(), DMPlexGetConeRecursiveVertices(), DMPlexGetDepth()
2666: @*/
2667: PetscErrorCode DMPlexGetConeRecursive(DM dm, IS points, PetscInt *depth, IS *expandedPoints[], PetscSection *sections[])
2668: {
2669: const PetscInt *arr0=NULL, *cone=NULL;
2670: PetscInt *arr=NULL, *newarr=NULL;
2671: PetscInt d, depth_, i, n, newn, cn, co, start, end;
2672: IS *expandedPoints_;
2673: PetscSection *sections_;
2674: PetscErrorCode ierr;
2682: ISGetLocalSize(points, &n);
2683: ISGetIndices(points, &arr0);
2684: DMPlexGetDepth(dm, &depth_);
2685: PetscCalloc1(depth_, &expandedPoints_);
2686: PetscCalloc1(depth_, §ions_);
2687: arr = (PetscInt*) arr0; /* this is ok because first generation of arr is not modified */
2688: for (d=depth_-1; d>=0; d--) {
2689: PetscSectionCreate(PETSC_COMM_SELF, §ions_[d]);
2690: PetscSectionSetChart(sections_[d], 0, n);
2691: for (i=0; i<n; i++) {
2692: DMPlexGetDepthStratum(dm, d+1, &start, &end);
2693: if (arr[i] >= start && arr[i] < end) {
2694: DMPlexGetConeSize(dm, arr[i], &cn);
2695: PetscSectionSetDof(sections_[d], i, cn);
2696: } else {
2697: PetscSectionSetDof(sections_[d], i, 1);
2698: }
2699: }
2700: PetscSectionSetUp(sections_[d]);
2701: PetscSectionGetStorageSize(sections_[d], &newn);
2702: PetscMalloc1(newn, &newarr);
2703: for (i=0; i<n; i++) {
2704: PetscSectionGetDof(sections_[d], i, &cn);
2705: PetscSectionGetOffset(sections_[d], i, &co);
2706: if (cn > 1) {
2707: DMPlexGetCone(dm, arr[i], &cone);
2708: PetscMemcpy(&newarr[co], cone, cn*sizeof(PetscInt));
2709: } else {
2710: newarr[co] = arr[i];
2711: }
2712: }
2713: ISCreateGeneral(PETSC_COMM_SELF, newn, newarr, PETSC_OWN_POINTER, &expandedPoints_[d]);
2714: arr = newarr;
2715: n = newn;
2716: }
2717: ISRestoreIndices(points, &arr0);
2718: *depth = depth_;
2719: if (expandedPoints) *expandedPoints = expandedPoints_;
2720: else {
2721: for (d=0; d<depth_; d++) {ISDestroy(&expandedPoints_[d]);}
2722: PetscFree(expandedPoints_);
2723: }
2724: if (sections) *sections = sections_;
2725: else {
2726: for (d=0; d<depth_; d++) {PetscSectionDestroy(§ions_[d]);}
2727: PetscFree(sections_);
2728: }
2729: return(0);
2730: }
2732: /*@
2733: DMPlexRestoreConeRecursive - Deallocates arrays created by DMPlexGetConeRecursive
2735: Not collective
2737: Input Parameters:
2738: + dm - The DMPlex
2739: - points - The IS of points, which must lie in the chart set with DMPlexSetChart()
2741: Output Parameters:
2742: + depth - (optional) Size of the output arrays, equal to DMPlex depth, returned by DMPlexGetDepth()
2743: . expandedPoints - (optional) An array of recursively expanded cones
2744: - sections - (optional) An array of sections which describe mappings from points to their cone points
2746: Level: advanced
2748: Notes:
2749: See DMPlexGetConeRecursive() for details.
2751: .seealso: DMPlexCreate(), DMPlexGetCone(), DMPlexGetConeTuple(), DMPlexGetConeRecursive(), DMPlexGetConeRecursiveVertices(), DMPlexGetDepth()
2752: @*/
2753: PetscErrorCode DMPlexRestoreConeRecursive(DM dm, IS points, PetscInt *depth, IS *expandedPoints[], PetscSection *sections[])
2754: {
2755: PetscInt d, depth_;
2756: PetscErrorCode ierr;
2759: DMPlexGetDepth(dm, &depth_);
2760: if (depth && *depth != depth_) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "depth changed since last call to DMPlexGetConeRecursive");
2761: if (depth) *depth = 0;
2762: if (expandedPoints) {
2763: for (d=0; d<depth_; d++) {ISDestroy(&((*expandedPoints)[d]));}
2764: PetscFree(*expandedPoints);
2765: }
2766: if (sections) {
2767: for (d=0; d<depth_; d++) {PetscSectionDestroy(&((*sections)[d]));}
2768: PetscFree(*sections);
2769: }
2770: return(0);
2771: }
2773: /*@
2774: DMPlexSetCone - Set the points on the in-edges for this point in the DAG; that is these are the points that cover the specific point
2776: Not collective
2778: Input Parameters:
2779: + mesh - The DMPlex
2780: . p - The point, which must lie in the chart set with DMPlexSetChart()
2781: - cone - An array of points which are on the in-edges for point p
2783: Output Parameter:
2785: Note:
2786: This should be called after all calls to DMPlexSetConeSize() and DMSetUp().
2788: Level: beginner
2790: .seealso: DMPlexCreate(), DMPlexGetCone(), DMPlexSetChart(), DMPlexSetConeSize(), DMSetUp(), DMPlexSetSupport(), DMPlexSetSupportSize()
2791: @*/
2792: PetscErrorCode DMPlexSetCone(DM dm, PetscInt p, const PetscInt cone[])
2793: {
2794: DM_Plex *mesh = (DM_Plex*) dm->data;
2795: PetscInt pStart, pEnd;
2796: PetscInt dof, off, c;
2801: PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd);
2802: PetscSectionGetDof(mesh->coneSection, p, &dof);
2804: PetscSectionGetOffset(mesh->coneSection, p, &off);
2805: if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D is not in the valid range [%D, %D)", p, pStart, pEnd);
2806: for (c = 0; c < dof; ++c) {
2807: if ((cone[c] < pStart) || (cone[c] >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone point %D is not in the valid range [%D, %D)", cone[c], pStart, pEnd);
2808: mesh->cones[off+c] = cone[c];
2809: }
2810: return(0);
2811: }
2813: /*@C
2814: DMPlexGetConeOrientation - Return the orientations on the in-edges for this point in the DAG
2816: Not collective
2818: Input Parameters:
2819: + mesh - The DMPlex
2820: - p - The point, which must lie in the chart set with DMPlexSetChart()
2822: Output Parameter:
2823: . coneOrientation - An array of orientations which are on the in-edges for point p. An orientation is an
2824: integer giving the prescription for cone traversal.
2826: Level: beginner
2828: Notes:
2829: The number indexes the symmetry transformations for the cell type (see manual). Orientation 0 is always
2830: the identity transformation. Negative orientation indicates reflection so that -(o+1) is the reflection
2831: of o, however it is not necessarily the inverse. To get the inverse, use DMPolytopeTypeComposeOrientationInv()
2832: with the identity.
2834: Fortran Notes:
2835: Since it returns an array, this routine is only available in Fortran 90, and you must
2836: include petsc.h90 in your code.
2837: You must also call DMPlexRestoreConeOrientation() after you finish using the returned array.
2838: DMPlexRestoreConeOrientation() is not needed/available in C.
2840: .seealso: DMPolytopeTypeComposeOrientation(), DMPolytopeTypeComposeOrientationInv(), DMPlexCreate(), DMPlexGetCone(), DMPlexSetCone(), DMPlexSetChart()
2841: @*/
2842: PetscErrorCode DMPlexGetConeOrientation(DM dm, PetscInt p, const PetscInt *coneOrientation[])
2843: {
2844: DM_Plex *mesh = (DM_Plex*) dm->data;
2845: PetscInt off;
2850: if (PetscDefined(USE_DEBUG)) {
2851: PetscInt dof;
2852: PetscSectionGetDof(mesh->coneSection, p, &dof);
2854: }
2855: PetscSectionGetOffset(mesh->coneSection, p, &off);
2857: *coneOrientation = &mesh->coneOrientations[off];
2858: return(0);
2859: }
2861: /*@
2862: DMPlexSetConeOrientation - Set the orientations on the in-edges for this point in the DAG
2864: Not collective
2866: Input Parameters:
2867: + mesh - The DMPlex
2868: . p - The point, which must lie in the chart set with DMPlexSetChart()
2869: - coneOrientation - An array of orientations
2870: Output Parameter:
2872: Notes:
2873: This should be called after all calls to DMPlexSetConeSize() and DMSetUp().
2875: The meaning of coneOrientation is detailed in DMPlexGetConeOrientation().
2877: Level: beginner
2879: .seealso: DMPlexCreate(), DMPlexGetConeOrientation(), DMPlexSetCone(), DMPlexSetChart(), DMPlexSetConeSize(), DMSetUp()
2880: @*/
2881: PetscErrorCode DMPlexSetConeOrientation(DM dm, PetscInt p, const PetscInt coneOrientation[])
2882: {
2883: DM_Plex *mesh = (DM_Plex*) dm->data;
2884: PetscInt pStart, pEnd;
2885: PetscInt dof, off, c;
2890: PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd);
2891: PetscSectionGetDof(mesh->coneSection, p, &dof);
2893: PetscSectionGetOffset(mesh->coneSection, p, &off);
2894: if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D is not in the valid range [%D, %D)", p, pStart, pEnd);
2895: for (c = 0; c < dof; ++c) {
2896: PetscInt cdof, o = coneOrientation[c];
2898: PetscSectionGetDof(mesh->coneSection, mesh->cones[off+c], &cdof);
2899: if (o && ((o < -(cdof+1)) || (o >= cdof))) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone orientation %D is not in the valid range [%D. %D)", o, -(cdof+1), cdof);
2900: mesh->coneOrientations[off+c] = o;
2901: }
2902: return(0);
2903: }
2905: /*@
2906: DMPlexInsertCone - Insert a point into the in-edges for the point p in the DAG
2908: Not collective
2910: Input Parameters:
2911: + mesh - The DMPlex
2912: . p - The point, which must lie in the chart set with DMPlexSetChart()
2913: . conePos - The local index in the cone where the point should be put
2914: - conePoint - The mesh point to insert
2916: Level: beginner
2918: .seealso: DMPlexCreate(), DMPlexGetCone(), DMPlexSetChart(), DMPlexSetConeSize(), DMSetUp()
2919: @*/
2920: PetscErrorCode DMPlexInsertCone(DM dm, PetscInt p, PetscInt conePos, PetscInt conePoint)
2921: {
2922: DM_Plex *mesh = (DM_Plex*) dm->data;
2923: PetscInt pStart, pEnd;
2924: PetscInt dof, off;
2929: PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd);
2930: if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D is not in the valid range [%D, %D)", p, pStart, pEnd);
2931: if ((conePoint < pStart) || (conePoint >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone point %D is not in the valid range [%D, %D)", conePoint, pStart, pEnd);
2932: PetscSectionGetDof(mesh->coneSection, p, &dof);
2933: PetscSectionGetOffset(mesh->coneSection, p, &off);
2934: if ((conePos < 0) || (conePos >= dof)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone position %D of point %D is not in the valid range [0, %D)", conePos, p, dof);
2935: mesh->cones[off+conePos] = conePoint;
2936: return(0);
2937: }
2939: /*@
2940: DMPlexInsertConeOrientation - Insert a point orientation for the in-edge for the point p in the DAG
2942: Not collective
2944: Input Parameters:
2945: + mesh - The DMPlex
2946: . p - The point, which must lie in the chart set with DMPlexSetChart()
2947: . conePos - The local index in the cone where the point should be put
2948: - coneOrientation - The point orientation to insert
2950: Level: beginner
2952: Notes:
2953: The meaning of coneOrientation values is detailed in DMPlexGetConeOrientation().
2955: .seealso: DMPlexCreate(), DMPlexGetCone(), DMPlexSetChart(), DMPlexSetConeSize(), DMSetUp()
2956: @*/
2957: PetscErrorCode DMPlexInsertConeOrientation(DM dm, PetscInt p, PetscInt conePos, PetscInt coneOrientation)
2958: {
2959: DM_Plex *mesh = (DM_Plex*) dm->data;
2960: PetscInt pStart, pEnd;
2961: PetscInt dof, off;
2966: PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd);
2967: if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D is not in the valid range [%D, %D)", p, pStart, pEnd);
2968: PetscSectionGetDof(mesh->coneSection, p, &dof);
2969: PetscSectionGetOffset(mesh->coneSection, p, &off);
2970: if ((conePos < 0) || (conePos >= dof)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cone position %D of point %D is not in the valid range [0, %D)", conePos, p, dof);
2971: mesh->coneOrientations[off+conePos] = coneOrientation;
2972: return(0);
2973: }
2975: /*@
2976: DMPlexGetSupportSize - Return the number of out-edges for this point in the DAG
2978: Not collective
2980: Input Parameters:
2981: + mesh - The DMPlex
2982: - p - The point, which must lie in the chart set with DMPlexSetChart()
2984: Output Parameter:
2985: . size - The support size for point p
2987: Level: beginner
2989: .seealso: DMPlexCreate(), DMPlexSetConeSize(), DMPlexSetChart(), DMPlexGetConeSize()
2990: @*/
2991: PetscErrorCode DMPlexGetSupportSize(DM dm, PetscInt p, PetscInt *size)
2992: {
2993: DM_Plex *mesh = (DM_Plex*) dm->data;
2999: PetscSectionGetDof(mesh->supportSection, p, size);
3000: return(0);
3001: }
3003: /*@
3004: DMPlexSetSupportSize - Set the number of out-edges for this point in the DAG
3006: Not collective
3008: Input Parameters:
3009: + mesh - The DMPlex
3010: . p - The point, which must lie in the chart set with DMPlexSetChart()
3011: - size - The support size for point p
3013: Output Parameter:
3015: Note:
3016: This should be called after DMPlexSetChart().
3018: Level: beginner
3020: .seealso: DMPlexCreate(), DMPlexGetSupportSize(), DMPlexSetChart()
3021: @*/
3022: PetscErrorCode DMPlexSetSupportSize(DM dm, PetscInt p, PetscInt size)
3023: {
3024: DM_Plex *mesh = (DM_Plex*) dm->data;
3029: PetscSectionSetDof(mesh->supportSection, p, size);
3031: mesh->maxSupportSize = PetscMax(mesh->maxSupportSize, size);
3032: return(0);
3033: }
3035: /*@C
3036: DMPlexGetSupport - Return the points on the out-edges for this point in the DAG
3038: Not collective
3040: Input Parameters:
3041: + mesh - The DMPlex
3042: - p - The point, which must lie in the chart set with DMPlexSetChart()
3044: Output Parameter:
3045: . support - An array of points which are on the out-edges for point p
3047: Level: beginner
3049: Fortran Notes:
3050: Since it returns an array, this routine is only available in Fortran 90, and you must
3051: include petsc.h90 in your code.
3052: You must also call DMPlexRestoreSupport() after you finish using the returned array.
3053: DMPlexRestoreSupport() is not needed/available in C.
3055: .seealso: DMPlexGetSupportSize(), DMPlexSetSupport(), DMPlexGetCone(), DMPlexSetChart()
3056: @*/
3057: PetscErrorCode DMPlexGetSupport(DM dm, PetscInt p, const PetscInt *support[])
3058: {
3059: DM_Plex *mesh = (DM_Plex*) dm->data;
3060: PetscInt off;
3066: PetscSectionGetOffset(mesh->supportSection, p, &off);
3067: *support = &mesh->supports[off];
3068: return(0);
3069: }
3071: /*@
3072: DMPlexSetSupport - Set the points on the out-edges for this point in the DAG, that is the list of points that this point covers
3074: Not collective
3076: Input Parameters:
3077: + mesh - The DMPlex
3078: . p - The point, which must lie in the chart set with DMPlexSetChart()
3079: - support - An array of points which are on the out-edges for point p
3081: Output Parameter:
3083: Note:
3084: This should be called after all calls to DMPlexSetSupportSize() and DMSetUp().
3086: Level: beginner
3088: .seealso: DMPlexSetCone(), DMPlexSetConeSize(), DMPlexCreate(), DMPlexGetSupport(), DMPlexSetChart(), DMPlexSetSupportSize(), DMSetUp()
3089: @*/
3090: PetscErrorCode DMPlexSetSupport(DM dm, PetscInt p, const PetscInt support[])
3091: {
3092: DM_Plex *mesh = (DM_Plex*) dm->data;
3093: PetscInt pStart, pEnd;
3094: PetscInt dof, off, c;
3099: PetscSectionGetChart(mesh->supportSection, &pStart, &pEnd);
3100: PetscSectionGetDof(mesh->supportSection, p, &dof);
3102: PetscSectionGetOffset(mesh->supportSection, p, &off);
3103: if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D is not in the valid range [%D, %D)", p, pStart, pEnd);
3104: for (c = 0; c < dof; ++c) {
3105: if ((support[c] < pStart) || (support[c] >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Support point %D is not in the valid range [%D, %D)", support[c], pStart, pEnd);
3106: mesh->supports[off+c] = support[c];
3107: }
3108: return(0);
3109: }
3111: /*@
3112: DMPlexInsertSupport - Insert a point into the out-edges for the point p in the DAG
3114: Not collective
3116: Input Parameters:
3117: + mesh - The DMPlex
3118: . p - The point, which must lie in the chart set with DMPlexSetChart()
3119: . supportPos - The local index in the cone where the point should be put
3120: - supportPoint - The mesh point to insert
3122: Level: beginner
3124: .seealso: DMPlexCreate(), DMPlexGetCone(), DMPlexSetChart(), DMPlexSetConeSize(), DMSetUp()
3125: @*/
3126: PetscErrorCode DMPlexInsertSupport(DM dm, PetscInt p, PetscInt supportPos, PetscInt supportPoint)
3127: {
3128: DM_Plex *mesh = (DM_Plex*) dm->data;
3129: PetscInt pStart, pEnd;
3130: PetscInt dof, off;
3135: PetscSectionGetChart(mesh->supportSection, &pStart, &pEnd);
3136: PetscSectionGetDof(mesh->supportSection, p, &dof);
3137: PetscSectionGetOffset(mesh->supportSection, p, &off);
3138: if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D is not in the valid range [%D, %D)", p, pStart, pEnd);
3139: if ((supportPoint < pStart) || (supportPoint >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Support point %D is not in the valid range [%D, %D)", supportPoint, pStart, pEnd);
3140: if (supportPos >= dof) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Support position %D of point %D is not in the valid range [0, %D)", supportPos, p, dof);
3141: mesh->supports[off+supportPos] = supportPoint;
3142: return(0);
3143: }
3145: /* Converts an orientation o in the current numbering to the previous scheme used in Plex */
3146: PetscInt DMPolytopeConvertNewOrientation_Internal(DMPolytopeType ct, PetscInt o)
3147: {
3148: switch (ct) {
3149: case DM_POLYTOPE_SEGMENT:
3150: if (o == -1) return -2;
3151: break;
3152: case DM_POLYTOPE_TRIANGLE:
3153: if (o == -3) return -1;
3154: if (o == -2) return -3;
3155: if (o == -1) return -2;
3156: break;
3157: case DM_POLYTOPE_QUADRILATERAL:
3158: if (o == -4) return -2;
3159: if (o == -3) return -1;
3160: if (o == -2) return -4;
3161: if (o == -1) return -3;
3162: break;
3163: default: return o;
3164: }
3165: return o;
3166: }
3168: /* Converts an orientation o in the previous scheme used in Plex to the current numbering */
3169: PetscInt DMPolytopeConvertOldOrientation_Internal(DMPolytopeType ct, PetscInt o)
3170: {
3171: switch (ct) {
3172: case DM_POLYTOPE_SEGMENT:
3173: if ((o == -2) || (o == 1)) return -1;
3174: if (o == -1) return 0;
3175: break;
3176: case DM_POLYTOPE_TRIANGLE:
3177: if (o == -3) return -2;
3178: if (o == -2) return -1;
3179: if (o == -1) return -3;
3180: break;
3181: case DM_POLYTOPE_QUADRILATERAL:
3182: if (o == -4) return -2;
3183: if (o == -3) return -1;
3184: if (o == -2) return -4;
3185: if (o == -1) return -3;
3186: break;
3187: default: return o;
3188: }
3189: return o;
3190: }
3192: /* Takes in a mesh whose orientations are in the previous scheme and converts them all to the current numbering */
3193: PetscErrorCode DMPlexConvertOldOrientations_Internal(DM dm)
3194: {
3195: PetscInt pStart, pEnd, p;
3199: DMPlexGetChart(dm, &pStart, &pEnd);
3200: for (p = pStart; p < pEnd; ++p) {
3201: const PetscInt *cone, *ornt;
3202: PetscInt coneSize, c;
3204: DMPlexGetConeSize(dm, p, &coneSize);
3205: DMPlexGetCone(dm, p, &cone);
3206: DMPlexGetConeOrientation(dm, p, &ornt);
3207: for (c = 0; c < coneSize; ++c) {
3208: DMPolytopeType ct;
3209: const PetscInt o = ornt[c];
3211: DMPlexGetCellType(dm, cone[c], &ct);
3212: switch (ct) {
3213: case DM_POLYTOPE_SEGMENT:
3214: if ((o == -2) || (o == 1)) {DMPlexInsertConeOrientation(dm, p, c, -1);}
3215: if (o == -1) {DMPlexInsertConeOrientation(dm, p, c, 0);}
3216: break;
3217: case DM_POLYTOPE_TRIANGLE:
3218: if (o == -3) {DMPlexInsertConeOrientation(dm, p, c, -2);}
3219: if (o == -2) {DMPlexInsertConeOrientation(dm, p, c, -1);}
3220: if (o == -1) {DMPlexInsertConeOrientation(dm, p, c, -3);}
3221: break;
3222: case DM_POLYTOPE_QUADRILATERAL:
3223: if (o == -4) {DMPlexInsertConeOrientation(dm, p, c, -2);}
3224: if (o == -3) {DMPlexInsertConeOrientation(dm, p, c, -1);}
3225: if (o == -2) {DMPlexInsertConeOrientation(dm, p, c, -4);}
3226: if (o == -1) {DMPlexInsertConeOrientation(dm, p, c, -3);}
3227: break;
3228: default: break;
3229: }
3230: }
3231: }
3232: return(0);
3233: }
3235: static PetscErrorCode DMPlexGetTransitiveClosure_Depth1_Private(DM dm, PetscInt p, PetscInt ornt, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
3236: {
3237: DMPolytopeType ct = DM_POLYTOPE_UNKNOWN;
3238: PetscInt *closure;
3239: const PetscInt *tmp = NULL, *tmpO = NULL;
3240: PetscInt off = 0, tmpSize, t;
3241: PetscErrorCode ierr;
3244: if (ornt) {
3245: DMPlexGetCellType(dm, p, &ct);
3246: if (ct == DM_POLYTOPE_FV_GHOST || ct == DM_POLYTOPE_INTERIOR_GHOST || ct == DM_POLYTOPE_UNKNOWN) ct = DM_POLYTOPE_UNKNOWN;
3247: }
3248: if (*points) {
3249: closure = *points;
3250: } else {
3251: PetscInt maxConeSize, maxSupportSize;
3252: DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);
3253: DMGetWorkArray(dm, 2*(PetscMax(maxConeSize, maxSupportSize)+1), MPIU_INT, &closure);
3254: }
3255: if (useCone) {
3256: DMPlexGetConeSize(dm, p, &tmpSize);
3257: DMPlexGetCone(dm, p, &tmp);
3258: DMPlexGetConeOrientation(dm, p, &tmpO);
3259: } else {
3260: DMPlexGetSupportSize(dm, p, &tmpSize);
3261: DMPlexGetSupport(dm, p, &tmp);
3262: }
3263: if (ct == DM_POLYTOPE_UNKNOWN) {
3264: closure[off++] = p;
3265: closure[off++] = 0;
3266: for (t = 0; t < tmpSize; ++t) {
3267: closure[off++] = tmp[t];
3268: closure[off++] = tmpO ? tmpO[t] : 0;
3269: }
3270: } else {
3271: const PetscInt *arr = DMPolytopeTypeGetArrangment(ct, ornt);
3273: /* We assume that cells with a valid type have faces with a valid type */
3274: closure[off++] = p;
3275: closure[off++] = ornt;
3276: for (t = 0; t < tmpSize; ++t) {
3277: DMPolytopeType ft;
3279: DMPlexGetCellType(dm, tmp[t], &ft);
3280: closure[off++] = tmp[arr[t]];
3281: closure[off++] = tmpO ? DMPolytopeTypeComposeOrientation(ft, ornt, tmpO[t]) : 0;
3282: }
3283: }
3284: if (numPoints) *numPoints = tmpSize+1;
3285: if (points) *points = closure;
3286: return(0);
3287: }
3289: /* We need a special tensor verison becasue we want to allow duplicate points in the endcaps for hybrid cells */
3290: static PetscErrorCode DMPlexTransitiveClosure_Tensor_Internal(DM dm, PetscInt point, DMPolytopeType ct, PetscInt o, PetscBool useCone, PetscInt *numPoints, PetscInt **points)
3291: {
3292: const PetscInt *arr = DMPolytopeTypeGetArrangment(ct, o);
3293: const PetscInt *cone, *ornt;
3294: PetscInt *pts, *closure = NULL;
3295: DMPolytopeType ft;
3296: PetscInt maxConeSize, maxSupportSize, coneSeries, supportSeries, maxSize;
3297: PetscInt dim, coneSize, c, d, clSize, cl;
3298: PetscErrorCode ierr;
3301: DMGetDimension(dm, &dim);
3302: DMPlexGetConeSize(dm, point, &coneSize);
3303: DMPlexGetCone(dm, point, &cone);
3304: DMPlexGetConeOrientation(dm, point, &ornt);
3305: DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);
3306: coneSeries = (maxConeSize > 1) ? ((PetscPowInt(maxConeSize, dim+1)-1)/(maxConeSize-1)) : dim+1;
3307: supportSeries = (maxSupportSize > 1) ? ((PetscPowInt(maxSupportSize, dim+1)-1)/(maxSupportSize-1)) : dim+1;
3308: maxSize = PetscMax(coneSeries, supportSeries);
3309: if (*points) {pts = *points;}
3310: else {DMGetWorkArray(dm, 2*maxSize, MPIU_INT, &pts);}
3311: c = 0;
3312: pts[c++] = point;
3313: pts[c++] = o;
3314: DMPlexGetCellType(dm, cone[arr[0*2+0]], &ft);
3315: DMPlexGetTransitiveClosure_Internal(dm, cone[arr[0*2+0]], DMPolytopeTypeComposeOrientation(ft, arr[0*2+1], ornt[0]), useCone, &clSize, &closure);
3316: for (cl = 0; cl < clSize*2; cl += 2) {pts[c++] = closure[cl]; pts[c++] = closure[cl+1];}
3317: DMPlexGetTransitiveClosure_Internal(dm, cone[arr[1*2+0]], DMPolytopeTypeComposeOrientation(ft, arr[1*2+1], ornt[1]), useCone, &clSize, &closure);
3318: for (cl = 0; cl < clSize*2; cl += 2) {pts[c++] = closure[cl]; pts[c++] = closure[cl+1];}
3319: DMPlexRestoreTransitiveClosure(dm, cone[0], useCone, &clSize, &closure);
3320: for (d = 2; d < coneSize; ++d) {
3321: DMPlexGetCellType(dm, cone[arr[d*2+0]], &ft);
3322: pts[c++] = cone[arr[d*2+0]];
3323: pts[c++] = DMPolytopeTypeComposeOrientation(ft, arr[d*2+1], ornt[d]);
3324: }
3325: if (dim >= 3) {
3326: for (d = 2; d < coneSize; ++d) {
3327: const PetscInt fpoint = cone[arr[d*2+0]];
3328: const PetscInt *fcone, *fornt;
3329: PetscInt fconeSize, fc, i;
3331: DMPlexGetCellType(dm, fpoint, &ft);
3332: const PetscInt *farr = DMPolytopeTypeGetArrangment(ft, DMPolytopeTypeComposeOrientation(ft, arr[d*2+1], ornt[d]));
3333: DMPlexGetConeSize(dm, fpoint, &fconeSize);
3334: DMPlexGetCone(dm, fpoint, &fcone);
3335: DMPlexGetConeOrientation(dm, fpoint, &fornt);
3336: for (fc = 0; fc < fconeSize; ++fc) {
3337: const PetscInt cp = fcone[farr[fc*2+0]];
3338: const PetscInt co = farr[fc*2+1];
3340: for (i = 0; i < c; i += 2) if (pts[i] == cp) break;
3341: if (i == c) {
3342: DMPlexGetCellType(dm, cp, &ft);
3343: pts[c++] = cp;
3344: pts[c++] = DMPolytopeTypeComposeOrientation(ft, co, fornt[farr[fc*2+0]]);
3345: }
3346: }
3347: }
3348: }
3349: *numPoints = c/2;
3350: *points = pts;
3351: return(0);
3352: }
3354: PetscErrorCode DMPlexGetTransitiveClosure_Internal(DM dm, PetscInt p, PetscInt ornt, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
3355: {
3356: DMPolytopeType ct;
3357: PetscInt *closure, *fifo;
3358: PetscInt closureSize = 0, fifoStart = 0, fifoSize = 0;
3359: PetscInt maxConeSize, maxSupportSize, coneSeries, supportSeries;
3360: PetscInt depth, maxSize;
3364: DMPlexGetDepth(dm, &depth);
3365: if (depth == 1) {
3366: DMPlexGetTransitiveClosure_Depth1_Private(dm, p, ornt, useCone, numPoints, points);
3367: return(0);
3368: }
3369: DMPlexGetCellType(dm, p, &ct);
3370: if (ct == DM_POLYTOPE_FV_GHOST || ct == DM_POLYTOPE_INTERIOR_GHOST || ct == DM_POLYTOPE_UNKNOWN) ct = DM_POLYTOPE_UNKNOWN;
3371: if (ct == DM_POLYTOPE_SEG_PRISM_TENSOR || ct == DM_POLYTOPE_TRI_PRISM_TENSOR || ct == DM_POLYTOPE_QUAD_PRISM_TENSOR) {
3372: DMPlexTransitiveClosure_Tensor_Internal(dm, p, ct, ornt, useCone, numPoints, points);
3373: return(0);
3374: }
3375: DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);
3376: coneSeries = (maxConeSize > 1) ? ((PetscPowInt(maxConeSize, depth+1)-1)/(maxConeSize-1)) : depth+1;
3377: supportSeries = (maxSupportSize > 1) ? ((PetscPowInt(maxSupportSize, depth+1)-1)/(maxSupportSize-1)) : depth+1;
3378: maxSize = PetscMax(coneSeries, supportSeries);
3379: DMGetWorkArray(dm, 3*maxSize, MPIU_INT, &fifo);
3380: if (*points) {closure = *points;}
3381: else {DMGetWorkArray(dm, 2*maxSize, MPIU_INT, &closure);}
3382: closure[closureSize++] = p;
3383: closure[closureSize++] = ornt;
3384: fifo[fifoSize++] = p;
3385: fifo[fifoSize++] = ornt;
3386: fifo[fifoSize++] = ct;
3387: /* Should kick out early when depth is reached, rather than checking all vertices for empty cones */
3388: while (fifoSize - fifoStart) {
3389: const PetscInt q = fifo[fifoStart++];
3390: const PetscInt o = fifo[fifoStart++];
3391: const DMPolytopeType qt = (DMPolytopeType) fifo[fifoStart++];
3392: const PetscInt *qarr = DMPolytopeTypeGetArrangment(qt, o);
3393: const PetscInt *tmp, *tmpO;
3394: PetscInt tmpSize, t;
3396: if (PetscDefined(USE_DEBUG)) {
3397: PetscInt nO = DMPolytopeTypeGetNumArrangments(qt)/2;
3398: if (o && (o >= nO || o < -nO)) SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid orientation %D not in [%D,%D) for %s %D", o, -nO, nO, DMPolytopeTypes[qt], q);
3399: }
3400: if (useCone) {
3401: DMPlexGetConeSize(dm, q, &tmpSize);
3402: DMPlexGetCone(dm, q, &tmp);
3403: DMPlexGetConeOrientation(dm, q, &tmpO);
3404: } else {
3405: DMPlexGetSupportSize(dm, q, &tmpSize);
3406: DMPlexGetSupport(dm, q, &tmp);
3407: tmpO = NULL;
3408: }
3409: for (t = 0; t < tmpSize; ++t) {
3410: const PetscInt ip = useCone && qarr ? qarr[t*2] : t;
3411: const PetscInt io = useCone && qarr ? qarr[t*2+1] : 0;
3412: const PetscInt cp = tmp[ip];
3413: DMPlexGetCellType(dm, cp, &ct);
3414: const PetscInt co = tmpO ? DMPolytopeTypeComposeOrientation(ct, io, tmpO[ip]) : 0;
3415: PetscInt c;
3417: /* Check for duplicate */
3418: for (c = 0; c < closureSize; c += 2) {
3419: if (closure[c] == cp) break;
3420: }
3421: if (c == closureSize) {
3422: closure[closureSize++] = cp;
3423: closure[closureSize++] = co;
3424: fifo[fifoSize++] = cp;
3425: fifo[fifoSize++] = co;
3426: fifo[fifoSize++] = ct;
3427: }
3428: }
3429: }
3430: DMRestoreWorkArray(dm, 3*maxSize, MPIU_INT, &fifo);
3431: if (numPoints) *numPoints = closureSize/2;
3432: if (points) *points = closure;
3433: return(0);
3434: }
3436: /*@C
3437: DMPlexGetTransitiveClosure - Return the points on the transitive closure of the in-edges or out-edges for this point in the DAG
3439: Not collective
3441: Input Parameters:
3442: + dm - The DMPlex
3443: . p - The mesh point
3444: - useCone - PETSC_TRUE for the closure, otherwise return the star
3446: Input/Output Parameter:
3447: . points - The points and point orientations, interleaved as pairs [p0, o0, p1, o1, ...];
3448: if NULL on input, internal storage will be returned, otherwise the provided array is used
3450: Output Parameter:
3451: . numPoints - The number of points in the closure, so points[] is of size 2*numPoints
3453: Note:
3454: If using internal storage (points is NULL on input), each call overwrites the last output.
3456: Fortran Notes:
3457: Since it returns an array, this routine is only available in Fortran 90, and you must include petsc.h90 in your code.
3459: The numPoints argument is not present in the Fortran 90 binding since it is internal to the array.
3461: Level: beginner
3463: .seealso: DMPlexRestoreTransitiveClosure(), DMPlexCreate(), DMPlexSetCone(), DMPlexSetChart(), DMPlexGetCone()
3464: @*/
3465: PetscErrorCode DMPlexGetTransitiveClosure(DM dm, PetscInt p, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
3466: {
3473: DMPlexGetTransitiveClosure_Internal(dm, p, 0, useCone, numPoints, points);
3474: return(0);
3475: }
3477: /*@C
3478: DMPlexRestoreTransitiveClosure - Restore the array of points on the transitive closure of the in-edges or out-edges for this point in the DAG
3480: Not collective
3482: Input Parameters:
3483: + dm - The DMPlex
3484: . p - The mesh point
3485: . useCone - PETSC_TRUE for the closure, otherwise return the star
3486: . numPoints - The number of points in the closure, so points[] is of size 2*numPoints
3487: - points - The points and point orientations, interleaved as pairs [p0, o0, p1, o1, ...]
3489: Note:
3490: If not using internal storage (points is not NULL on input), this call is unnecessary
3492: Fortran Notes:
3493: Since it returns an array, this routine is only available in Fortran 90, and you must include petsc.h90 in your code.
3495: The numPoints argument is not present in the Fortran 90 binding since it is internal to the array.
3497: Level: beginner
3499: .seealso: DMPlexGetTransitiveClosure(), DMPlexCreate(), DMPlexSetCone(), DMPlexSetChart(), DMPlexGetCone()
3500: @*/
3501: PetscErrorCode DMPlexRestoreTransitiveClosure(DM dm, PetscInt p, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
3502: {
3507: if (numPoints) *numPoints = 0;
3508: DMRestoreWorkArray(dm, 0, MPIU_INT, points);
3509: return(0);
3510: }
3512: /*@
3513: DMPlexGetMaxSizes - Return the maximum number of in-edges (cone) and out-edges (support) for any point in the DAG
3515: Not collective
3517: Input Parameter:
3518: . mesh - The DMPlex
3520: Output Parameters:
3521: + maxConeSize - The maximum number of in-edges
3522: - maxSupportSize - The maximum number of out-edges
3524: Level: beginner
3526: .seealso: DMPlexCreate(), DMPlexSetConeSize(), DMPlexSetChart()
3527: @*/
3528: PetscErrorCode DMPlexGetMaxSizes(DM dm, PetscInt *maxConeSize, PetscInt *maxSupportSize)
3529: {
3530: DM_Plex *mesh = (DM_Plex*) dm->data;
3534: if (maxConeSize) *maxConeSize = mesh->maxConeSize;
3535: if (maxSupportSize) *maxSupportSize = mesh->maxSupportSize;
3536: return(0);
3537: }
3539: PetscErrorCode DMSetUp_Plex(DM dm)
3540: {
3541: DM_Plex *mesh = (DM_Plex*) dm->data;
3542: PetscInt size;
3547: PetscSectionSetUp(mesh->coneSection);
3548: PetscSectionGetStorageSize(mesh->coneSection, &size);
3549: PetscMalloc1(size, &mesh->cones);
3550: PetscCalloc1(size, &mesh->coneOrientations);
3551: PetscLogObjectMemory((PetscObject) dm, size*2*sizeof(PetscInt));
3552: if (mesh->maxSupportSize) {
3553: PetscSectionSetUp(mesh->supportSection);
3554: PetscSectionGetStorageSize(mesh->supportSection, &size);
3555: PetscMalloc1(size, &mesh->supports);
3556: PetscLogObjectMemory((PetscObject) dm, size*sizeof(PetscInt));
3557: }
3558: return(0);
3559: }
3561: PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
3562: {
3566: if (subdm) {DMClone(dm, subdm);}
3567: DMCreateSectionSubDM(dm, numFields, fields, is, subdm);
3568: if (subdm) {(*subdm)->useNatural = dm->useNatural;}
3569: if (dm->useNatural && dm->sfMigration) {
3570: PetscSF sfMigrationInv,sfNatural;
3571: PetscSection section, sectionSeq;
3573: (*subdm)->sfMigration = dm->sfMigration;
3574: PetscObjectReference((PetscObject) dm->sfMigration);
3575: DMGetLocalSection((*subdm), §ion);
3576: PetscSFCreateInverseSF((*subdm)->sfMigration, &sfMigrationInv);
3577: PetscSectionCreate(PetscObjectComm((PetscObject) (*subdm)), §ionSeq);
3578: PetscSFDistributeSection(sfMigrationInv, section, NULL, sectionSeq);
3580: DMPlexCreateGlobalToNaturalSF(*subdm, sectionSeq, (*subdm)->sfMigration, &sfNatural);
3581: (*subdm)->sfNatural = sfNatural;
3582: PetscSectionDestroy(§ionSeq);
3583: PetscSFDestroy(&sfMigrationInv);
3584: }
3585: return(0);
3586: }
3588: PetscErrorCode DMCreateSuperDM_Plex(DM dms[], PetscInt len, IS **is, DM *superdm)
3589: {
3591: PetscInt i = 0;
3594: DMClone(dms[0], superdm);
3595: DMCreateSectionSuperDM(dms, len, is, superdm);
3596: (*superdm)->useNatural = PETSC_FALSE;
3597: for (i = 0; i < len; i++) {
3598: if (dms[i]->useNatural && dms[i]->sfMigration) {
3599: PetscSF sfMigrationInv,sfNatural;
3600: PetscSection section, sectionSeq;
3602: (*superdm)->sfMigration = dms[i]->sfMigration;
3603: PetscObjectReference((PetscObject) dms[i]->sfMigration);
3604: (*superdm)->useNatural = PETSC_TRUE;
3605: DMGetLocalSection((*superdm), §ion);
3606: PetscSFCreateInverseSF((*superdm)->sfMigration, &sfMigrationInv);
3607: PetscSectionCreate(PetscObjectComm((PetscObject) (*superdm)), §ionSeq);
3608: PetscSFDistributeSection(sfMigrationInv, section, NULL, sectionSeq);
3610: DMPlexCreateGlobalToNaturalSF(*superdm, sectionSeq, (*superdm)->sfMigration, &sfNatural);
3611: (*superdm)->sfNatural = sfNatural;
3612: PetscSectionDestroy(§ionSeq);
3613: PetscSFDestroy(&sfMigrationInv);
3614: break;
3615: }
3616: }
3617: return(0);
3618: }
3620: /*@
3621: DMPlexSymmetrize - Create support (out-edge) information from cone (in-edge) information
3623: Not collective
3625: Input Parameter:
3626: . mesh - The DMPlex
3628: Output Parameter:
3630: Note:
3631: This should be called after all calls to DMPlexSetCone()
3633: Level: beginner
3635: .seealso: DMPlexCreate(), DMPlexSetChart(), DMPlexSetConeSize(), DMPlexSetCone()
3636: @*/
3637: PetscErrorCode DMPlexSymmetrize(DM dm)
3638: {
3639: DM_Plex *mesh = (DM_Plex*) dm->data;
3640: PetscInt *offsets;
3641: PetscInt supportSize;
3642: PetscInt pStart, pEnd, p;
3647: if (mesh->supports) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Supports were already setup in this DMPlex");
3648: PetscLogEventBegin(DMPLEX_Symmetrize,dm,0,0,0);
3649: /* Calculate support sizes */
3650: DMPlexGetChart(dm, &pStart, &pEnd);
3651: for (p = pStart; p < pEnd; ++p) {
3652: PetscInt dof, off, c;
3654: PetscSectionGetDof(mesh->coneSection, p, &dof);
3655: PetscSectionGetOffset(mesh->coneSection, p, &off);
3656: for (c = off; c < off+dof; ++c) {
3657: PetscSectionAddDof(mesh->supportSection, mesh->cones[c], 1);
3658: }
3659: }
3660: for (p = pStart; p < pEnd; ++p) {
3661: PetscInt dof;
3663: PetscSectionGetDof(mesh->supportSection, p, &dof);
3665: mesh->maxSupportSize = PetscMax(mesh->maxSupportSize, dof);
3666: }
3667: PetscSectionSetUp(mesh->supportSection);
3668: /* Calculate supports */
3669: PetscSectionGetStorageSize(mesh->supportSection, &supportSize);
3670: PetscMalloc1(supportSize, &mesh->supports);
3671: PetscCalloc1(pEnd - pStart, &offsets);
3672: for (p = pStart; p < pEnd; ++p) {
3673: PetscInt dof, off, c;
3675: PetscSectionGetDof(mesh->coneSection, p, &dof);
3676: PetscSectionGetOffset(mesh->coneSection, p, &off);
3677: for (c = off; c < off+dof; ++c) {
3678: const PetscInt q = mesh->cones[c];
3679: PetscInt offS;
3681: PetscSectionGetOffset(mesh->supportSection, q, &offS);
3683: mesh->supports[offS+offsets[q]] = p;
3684: ++offsets[q];
3685: }
3686: }
3687: PetscFree(offsets);
3688: PetscLogEventEnd(DMPLEX_Symmetrize,dm,0,0,0);
3689: return(0);
3690: }
3692: static PetscErrorCode DMPlexCreateDepthStratum(DM dm, DMLabel label, PetscInt depth, PetscInt pStart, PetscInt pEnd)
3693: {
3694: IS stratumIS;
3698: if (pStart >= pEnd) return(0);
3699: if (PetscDefined(USE_DEBUG)) {
3700: PetscInt qStart, qEnd, numLevels, level;
3701: PetscBool overlap = PETSC_FALSE;
3702: DMLabelGetNumValues(label, &numLevels);
3703: for (level = 0; level < numLevels; level++) {
3704: DMLabelGetStratumBounds(label, level, &qStart, &qEnd);
3705: if ((pStart >= qStart && pStart < qEnd) || (pEnd > qStart && pEnd <= qEnd)) {overlap = PETSC_TRUE; break;}
3706: }
3707: if (overlap) SETERRQ6(PETSC_COMM_SELF, PETSC_ERR_PLIB, "New depth %D range [%D,%D) overlaps with depth %D range [%D,%D)", depth, pStart, pEnd, level, qStart, qEnd);
3708: }
3709: ISCreateStride(PETSC_COMM_SELF, pEnd-pStart, pStart, 1, &stratumIS);
3710: DMLabelSetStratumIS(label, depth, stratumIS);
3711: ISDestroy(&stratumIS);
3712: return(0);
3713: }
3715: /*@
3716: DMPlexStratify - The DAG for most topologies is a graded poset (https://en.wikipedia.org/wiki/Graded_poset), and
3717: can be illustrated by a Hasse Diagram (https://en.wikipedia.org/wiki/Hasse_diagram). The strata group all points of the
3718: same grade, and this function calculates the strata. This grade can be seen as the height (or depth) of the point in
3719: the DAG.
3721: Collective on dm
3723: Input Parameter:
3724: . mesh - The DMPlex
3726: Output Parameter:
3728: Notes:
3729: Concretely, DMPlexStratify() creates a new label named "depth" containing the depth in the DAG of each point. For cell-vertex
3730: meshes, vertices are depth 0 and cells are depth 1. For fully interpolated meshes, depth 0 for vertices, 1 for edges, and so on
3731: until cells have depth equal to the dimension of the mesh. The depth label can be accessed through DMPlexGetDepthLabel() or DMPlexGetDepthStratum(), or
3732: manually via DMGetLabel(). The height is defined implicitly by height = maxDimension - depth, and can be accessed
3733: via DMPlexGetHeightStratum(). For example, cells have height 0 and faces have height 1.
3735: The depth of a point is calculated by executing a breadth-first search (BFS) on the DAG. This could produce surprising results
3736: if run on a partially interpolated mesh, meaning one that had some edges and faces, but not others. For example, suppose that
3737: we had a mesh consisting of one triangle (c0) and three vertices (v0, v1, v2), and only one edge is on the boundary so we choose
3738: to interpolate only that one (e0), so that
3739: $ cone(c0) = {e0, v2}
3740: $ cone(e0) = {v0, v1}
3741: If DMPlexStratify() is run on this mesh, it will give depths
3742: $ depth 0 = {v0, v1, v2}
3743: $ depth 1 = {e0, c0}
3744: where the triangle has been given depth 1, instead of 2, because it is reachable from vertex v2.
3746: DMPlexStratify() should be called after all calls to DMPlexSymmetrize()
3748: Level: beginner
3750: .seealso: DMPlexCreate(), DMPlexSymmetrize(), DMPlexComputeCellTypes()
3751: @*/
3752: PetscErrorCode DMPlexStratify(DM dm)
3753: {
3754: DM_Plex *mesh = (DM_Plex*) dm->data;
3755: DMLabel label;
3756: PetscInt pStart, pEnd, p;
3757: PetscInt numRoots = 0, numLeaves = 0;
3762: PetscLogEventBegin(DMPLEX_Stratify,dm,0,0,0);
3764: /* Create depth label */
3765: DMPlexGetChart(dm, &pStart, &pEnd);
3766: DMCreateLabel(dm, "depth");
3767: DMPlexGetDepthLabel(dm, &label);
3769: {
3770: /* Initialize roots and count leaves */
3771: PetscInt sMin = PETSC_MAX_INT;
3772: PetscInt sMax = PETSC_MIN_INT;
3773: PetscInt coneSize, supportSize;
3775: for (p = pStart; p < pEnd; ++p) {
3776: DMPlexGetConeSize(dm, p, &coneSize);
3777: DMPlexGetSupportSize(dm, p, &supportSize);
3778: if (!coneSize && supportSize) {
3779: sMin = PetscMin(p, sMin);
3780: sMax = PetscMax(p, sMax);
3781: ++numRoots;
3782: } else if (!supportSize && coneSize) {
3783: ++numLeaves;
3784: } else if (!supportSize && !coneSize) {
3785: /* Isolated points */
3786: sMin = PetscMin(p, sMin);
3787: sMax = PetscMax(p, sMax);
3788: }
3789: }
3790: DMPlexCreateDepthStratum(dm, label, 0, sMin, sMax+1);
3791: }
3793: if (numRoots + numLeaves == (pEnd - pStart)) {
3794: PetscInt sMin = PETSC_MAX_INT;
3795: PetscInt sMax = PETSC_MIN_INT;
3796: PetscInt coneSize, supportSize;
3798: for (p = pStart; p < pEnd; ++p) {
3799: DMPlexGetConeSize(dm, p, &coneSize);
3800: DMPlexGetSupportSize(dm, p, &supportSize);
3801: if (!supportSize && coneSize) {
3802: sMin = PetscMin(p, sMin);
3803: sMax = PetscMax(p, sMax);
3804: }
3805: }
3806: DMPlexCreateDepthStratum(dm, label, 1, sMin, sMax+1);
3807: } else {
3808: PetscInt level = 0;
3809: PetscInt qStart, qEnd, q;
3811: DMLabelGetStratumBounds(label, level, &qStart, &qEnd);
3812: while (qEnd > qStart) {
3813: PetscInt sMin = PETSC_MAX_INT;
3814: PetscInt sMax = PETSC_MIN_INT;
3816: for (q = qStart; q < qEnd; ++q) {
3817: const PetscInt *support;
3818: PetscInt supportSize, s;
3820: DMPlexGetSupportSize(dm, q, &supportSize);
3821: DMPlexGetSupport(dm, q, &support);
3822: for (s = 0; s < supportSize; ++s) {
3823: sMin = PetscMin(support[s], sMin);
3824: sMax = PetscMax(support[s], sMax);
3825: }
3826: }
3827: DMLabelGetNumValues(label, &level);
3828: DMPlexCreateDepthStratum(dm, label, level, sMin, sMax+1);
3829: DMLabelGetStratumBounds(label, level, &qStart, &qEnd);
3830: }
3831: }
3832: { /* just in case there is an empty process */
3833: PetscInt numValues, maxValues = 0, v;
3835: DMLabelGetNumValues(label, &numValues);
3836: MPI_Allreduce(&numValues,&maxValues,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)dm));
3837: for (v = numValues; v < maxValues; v++) {
3838: DMLabelAddStratum(label, v);
3839: }
3840: }
3841: PetscObjectStateGet((PetscObject) label, &mesh->depthState);
3842: PetscLogEventEnd(DMPLEX_Stratify,dm,0,0,0);
3843: return(0);
3844: }
3846: PetscErrorCode DMPlexComputeCellType_Internal(DM dm, PetscInt p, PetscInt pdepth, DMPolytopeType *pt)
3847: {
3848: DMPolytopeType ct = DM_POLYTOPE_UNKNOWN;
3849: PetscInt dim, depth, pheight, coneSize;
3853: DMGetDimension(dm, &dim);
3854: DMPlexGetDepth(dm, &depth);
3855: DMPlexGetConeSize(dm, p, &coneSize);
3856: pheight = depth - pdepth;
3857: if (depth <= 1) {
3858: switch (pdepth) {
3859: case 0: ct = DM_POLYTOPE_POINT;break;
3860: case 1:
3861: switch (coneSize) {
3862: case 2: ct = DM_POLYTOPE_SEGMENT;break;
3863: case 3: ct = DM_POLYTOPE_TRIANGLE;break;
3864: case 4:
3865: switch (dim) {
3866: case 2: ct = DM_POLYTOPE_QUADRILATERAL;break;
3867: case 3: ct = DM_POLYTOPE_TETRAHEDRON;break;
3868: default: break;
3869: }
3870: break;
3871: case 5: ct = DM_POLYTOPE_PYRAMID;break;
3872: case 6: ct = DM_POLYTOPE_TRI_PRISM_TENSOR;break;
3873: case 8: ct = DM_POLYTOPE_HEXAHEDRON;break;
3874: default: break;
3875: }
3876: }
3877: } else {
3878: if (pdepth == 0) {
3879: ct = DM_POLYTOPE_POINT;
3880: } else if (pheight == 0) {
3881: switch (dim) {
3882: case 1:
3883: switch (coneSize) {
3884: case 2: ct = DM_POLYTOPE_SEGMENT;break;
3885: default: break;
3886: }
3887: break;
3888: case 2:
3889: switch (coneSize) {
3890: case 3: ct = DM_POLYTOPE_TRIANGLE;break;
3891: case 4: ct = DM_POLYTOPE_QUADRILATERAL;break;
3892: default: break;
3893: }
3894: break;
3895: case 3:
3896: switch (coneSize) {
3897: case 4: ct = DM_POLYTOPE_TETRAHEDRON;break;
3898: case 5:
3899: {
3900: const PetscInt *cone;
3901: PetscInt faceConeSize;
3903: DMPlexGetCone(dm, p, &cone);
3904: DMPlexGetConeSize(dm, cone[0], &faceConeSize);
3905: switch (faceConeSize) {
3906: case 3: ct = DM_POLYTOPE_TRI_PRISM_TENSOR;break;
3907: case 4: ct = DM_POLYTOPE_PYRAMID;break;
3908: }
3909: }
3910: break;
3911: case 6: ct = DM_POLYTOPE_HEXAHEDRON;break;
3912: default: break;
3913: }
3914: break;
3915: default: break;
3916: }
3917: } else if (pheight > 0) {
3918: switch (coneSize) {
3919: case 2: ct = DM_POLYTOPE_SEGMENT;break;
3920: case 3: ct = DM_POLYTOPE_TRIANGLE;break;
3921: case 4: ct = DM_POLYTOPE_QUADRILATERAL;break;
3922: default: break;
3923: }
3924: }
3925: }
3926: *pt = ct;
3927: return(0);
3928: }
3930: /*@
3931: DMPlexComputeCellTypes - Infer the polytope type of every cell using its dimension and cone size.
3933: Collective on dm
3935: Input Parameter:
3936: . mesh - The DMPlex
3938: DMPlexComputeCellTypes() should be called after all calls to DMPlexSymmetrize() and DMPlexStratify()
3940: Level: developer
3942: Note: This function is normally called automatically by Plex when a cell type is requested. It creates an
3943: internal DMLabel named "celltype" which can be directly accessed using DMGetLabel(). A user may disable
3944: automatic creation by creating the label manually, using DMCreateLabel(dm, "celltype").
3946: .seealso: DMPlexCreate(), DMPlexSymmetrize(), DMPlexStratify(), DMGetLabel(), DMCreateLabel()
3947: @*/
3948: PetscErrorCode DMPlexComputeCellTypes(DM dm)
3949: {
3950: DM_Plex *mesh;
3951: DMLabel ctLabel;
3952: PetscInt pStart, pEnd, p;
3957: mesh = (DM_Plex *) dm->data;
3958: DMCreateLabel(dm, "celltype");
3959: DMPlexGetCellTypeLabel(dm, &ctLabel);
3960: DMPlexGetChart(dm, &pStart, &pEnd);
3961: for (p = pStart; p < pEnd; ++p) {
3962: DMPolytopeType ct = DM_POLYTOPE_UNKNOWN;
3963: PetscInt pdepth;
3965: DMPlexGetPointDepth(dm, p, &pdepth);
3966: DMPlexComputeCellType_Internal(dm, p, pdepth, &ct);
3967: if (ct == DM_POLYTOPE_UNKNOWN) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Point %D is screwed up", p);
3968: DMLabelSetValue(ctLabel, p, ct);
3969: }
3970: PetscObjectStateGet((PetscObject) ctLabel, &mesh->celltypeState);
3971: PetscObjectViewFromOptions((PetscObject) ctLabel, NULL, "-dm_plex_celltypes_view");
3972: return(0);
3973: }
3975: /*@C
3976: DMPlexGetJoin - Get an array for the join of the set of points
3978: Not Collective
3980: Input Parameters:
3981: + dm - The DMPlex object
3982: . numPoints - The number of input points for the join
3983: - points - The input points
3985: Output Parameters:
3986: + numCoveredPoints - The number of points in the join
3987: - coveredPoints - The points in the join
3989: Level: intermediate
3991: Note: Currently, this is restricted to a single level join
3993: Fortran Notes:
3994: Since it returns an array, this routine is only available in Fortran 90, and you must
3995: include petsc.h90 in your code.
3997: The numCoveredPoints argument is not present in the Fortran 90 binding since it is internal to the array.
3999: .seealso: DMPlexRestoreJoin(), DMPlexGetMeet()
4000: @*/
4001: PetscErrorCode DMPlexGetJoin(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt **coveredPoints)
4002: {
4003: DM_Plex *mesh = (DM_Plex*) dm->data;
4004: PetscInt *join[2];
4005: PetscInt joinSize, i = 0;
4006: PetscInt dof, off, p, c, m;
4014: DMGetWorkArray(dm, mesh->maxSupportSize, MPIU_INT, &join[0]);
4015: DMGetWorkArray(dm, mesh->maxSupportSize, MPIU_INT, &join[1]);
4016: /* Copy in support of first point */
4017: PetscSectionGetDof(mesh->supportSection, points[0], &dof);
4018: PetscSectionGetOffset(mesh->supportSection, points[0], &off);
4019: for (joinSize = 0; joinSize < dof; ++joinSize) {
4020: join[i][joinSize] = mesh->supports[off+joinSize];
4021: }
4022: /* Check each successive support */
4023: for (p = 1; p < numPoints; ++p) {
4024: PetscInt newJoinSize = 0;
4026: PetscSectionGetDof(mesh->supportSection, points[p], &dof);
4027: PetscSectionGetOffset(mesh->supportSection, points[p], &off);
4028: for (c = 0; c < dof; ++c) {
4029: const PetscInt point = mesh->supports[off+c];
4031: for (m = 0; m < joinSize; ++m) {
4032: if (point == join[i][m]) {
4033: join[1-i][newJoinSize++] = point;
4034: break;
4035: }
4036: }
4037: }
4038: joinSize = newJoinSize;
4039: i = 1-i;
4040: }
4041: *numCoveredPoints = joinSize;
4042: *coveredPoints = join[i];
4043: DMRestoreWorkArray(dm, mesh->maxSupportSize, MPIU_INT, &join[1-i]);
4044: return(0);
4045: }
4047: /*@C
4048: DMPlexRestoreJoin - Restore an array for the join of the set of points
4050: Not Collective
4052: Input Parameters:
4053: + dm - The DMPlex object
4054: . numPoints - The number of input points for the join
4055: - points - The input points
4057: Output Parameters:
4058: + numCoveredPoints - The number of points in the join
4059: - coveredPoints - The points in the join
4061: Fortran Notes:
4062: Since it returns an array, this routine is only available in Fortran 90, and you must
4063: include petsc.h90 in your code.
4065: The numCoveredPoints argument is not present in the Fortran 90 binding since it is internal to the array.
4067: Level: intermediate
4069: .seealso: DMPlexGetJoin(), DMPlexGetFullJoin(), DMPlexGetMeet()
4070: @*/
4071: PetscErrorCode DMPlexRestoreJoin(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt **coveredPoints)
4072: {
4080: DMRestoreWorkArray(dm, 0, MPIU_INT, (void*) coveredPoints);
4081: if (numCoveredPoints) *numCoveredPoints = 0;
4082: return(0);
4083: }
4085: /*@C
4086: DMPlexGetFullJoin - Get an array for the join of the set of points
4088: Not Collective
4090: Input Parameters:
4091: + dm - The DMPlex object
4092: . numPoints - The number of input points for the join
4093: - points - The input points
4095: Output Parameters:
4096: + numCoveredPoints - The number of points in the join
4097: - coveredPoints - The points in the join
4099: Fortran Notes:
4100: Since it returns an array, this routine is only available in Fortran 90, and you must
4101: include petsc.h90 in your code.
4103: The numCoveredPoints argument is not present in the Fortran 90 binding since it is internal to the array.
4105: Level: intermediate
4107: .seealso: DMPlexGetJoin(), DMPlexRestoreJoin(), DMPlexGetMeet()
4108: @*/
4109: PetscErrorCode DMPlexGetFullJoin(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt **coveredPoints)
4110: {
4111: DM_Plex *mesh = (DM_Plex*) dm->data;
4112: PetscInt *offsets, **closures;
4113: PetscInt *join[2];
4114: PetscInt depth = 0, maxSize, joinSize = 0, i = 0;
4115: PetscInt p, d, c, m, ms;
4124: DMPlexGetDepth(dm, &depth);
4125: PetscCalloc1(numPoints, &closures);
4126: DMGetWorkArray(dm, numPoints*(depth+2), MPIU_INT, &offsets);
4127: ms = mesh->maxSupportSize;
4128: maxSize = (ms > 1) ? ((PetscPowInt(ms,depth+1)-1)/(ms-1)) : depth + 1;
4129: DMGetWorkArray(dm, maxSize, MPIU_INT, &join[0]);
4130: DMGetWorkArray(dm, maxSize, MPIU_INT, &join[1]);
4132: for (p = 0; p < numPoints; ++p) {
4133: PetscInt closureSize;
4135: DMPlexGetTransitiveClosure(dm, points[p], PETSC_FALSE, &closureSize, &closures[p]);
4137: offsets[p*(depth+2)+0] = 0;
4138: for (d = 0; d < depth+1; ++d) {
4139: PetscInt pStart, pEnd, i;
4141: DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
4142: for (i = offsets[p*(depth+2)+d]; i < closureSize; ++i) {
4143: if ((pStart > closures[p][i*2]) || (pEnd <= closures[p][i*2])) {
4144: offsets[p*(depth+2)+d+1] = i;
4145: break;
4146: }
4147: }
4148: if (i == closureSize) offsets[p*(depth+2)+d+1] = i;
4149: }
4150: if (offsets[p*(depth+2)+depth+1] != closureSize) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Total size of closure %D should be %D", offsets[p*(depth+2)+depth+1], closureSize);
4151: }
4152: for (d = 0; d < depth+1; ++d) {
4153: PetscInt dof;
4155: /* Copy in support of first point */
4156: dof = offsets[d+1] - offsets[d];
4157: for (joinSize = 0; joinSize < dof; ++joinSize) {
4158: join[i][joinSize] = closures[0][(offsets[d]+joinSize)*2];
4159: }
4160: /* Check each successive cone */
4161: for (p = 1; p < numPoints && joinSize; ++p) {
4162: PetscInt newJoinSize = 0;
4164: dof = offsets[p*(depth+2)+d+1] - offsets[p*(depth+2)+d];
4165: for (c = 0; c < dof; ++c) {
4166: const PetscInt point = closures[p][(offsets[p*(depth+2)+d]+c)*2];
4168: for (m = 0; m < joinSize; ++m) {
4169: if (point == join[i][m]) {
4170: join[1-i][newJoinSize++] = point;
4171: break;
4172: }
4173: }
4174: }
4175: joinSize = newJoinSize;
4176: i = 1-i;
4177: }
4178: if (joinSize) break;
4179: }
4180: *numCoveredPoints = joinSize;
4181: *coveredPoints = join[i];
4182: for (p = 0; p < numPoints; ++p) {
4183: DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_FALSE, NULL, &closures[p]);
4184: }
4185: PetscFree(closures);
4186: DMRestoreWorkArray(dm, numPoints*(depth+2), MPIU_INT, &offsets);
4187: DMRestoreWorkArray(dm, mesh->maxSupportSize, MPIU_INT, &join[1-i]);
4188: return(0);
4189: }
4191: /*@C
4192: DMPlexGetMeet - Get an array for the meet of the set of points
4194: Not Collective
4196: Input Parameters:
4197: + dm - The DMPlex object
4198: . numPoints - The number of input points for the meet
4199: - points - The input points
4201: Output Parameters:
4202: + numCoveredPoints - The number of points in the meet
4203: - coveredPoints - The points in the meet
4205: Level: intermediate
4207: Note: Currently, this is restricted to a single level meet
4209: Fortran Notes:
4210: Since it returns an array, this routine is only available in Fortran 90, and you must
4211: include petsc.h90 in your code.
4213: The numCoveredPoints argument is not present in the Fortran 90 binding since it is internal to the array.
4215: .seealso: DMPlexRestoreMeet(), DMPlexGetJoin()
4216: @*/
4217: PetscErrorCode DMPlexGetMeet(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveringPoints, const PetscInt **coveringPoints)
4218: {
4219: DM_Plex *mesh = (DM_Plex*) dm->data;
4220: PetscInt *meet[2];
4221: PetscInt meetSize, i = 0;
4222: PetscInt dof, off, p, c, m;
4230: DMGetWorkArray(dm, mesh->maxConeSize, MPIU_INT, &meet[0]);
4231: DMGetWorkArray(dm, mesh->maxConeSize, MPIU_INT, &meet[1]);
4232: /* Copy in cone of first point */
4233: PetscSectionGetDof(mesh->coneSection, points[0], &dof);
4234: PetscSectionGetOffset(mesh->coneSection, points[0], &off);
4235: for (meetSize = 0; meetSize < dof; ++meetSize) {
4236: meet[i][meetSize] = mesh->cones[off+meetSize];
4237: }
4238: /* Check each successive cone */
4239: for (p = 1; p < numPoints; ++p) {
4240: PetscInt newMeetSize = 0;
4242: PetscSectionGetDof(mesh->coneSection, points[p], &dof);
4243: PetscSectionGetOffset(mesh->coneSection, points[p], &off);
4244: for (c = 0; c < dof; ++c) {
4245: const PetscInt point = mesh->cones[off+c];
4247: for (m = 0; m < meetSize; ++m) {
4248: if (point == meet[i][m]) {
4249: meet[1-i][newMeetSize++] = point;
4250: break;
4251: }
4252: }
4253: }
4254: meetSize = newMeetSize;
4255: i = 1-i;
4256: }
4257: *numCoveringPoints = meetSize;
4258: *coveringPoints = meet[i];
4259: DMRestoreWorkArray(dm, mesh->maxConeSize, MPIU_INT, &meet[1-i]);
4260: return(0);
4261: }
4263: /*@C
4264: DMPlexRestoreMeet - Restore an array for the meet of the set of points
4266: Not Collective
4268: Input Parameters:
4269: + dm - The DMPlex object
4270: . numPoints - The number of input points for the meet
4271: - points - The input points
4273: Output Parameters:
4274: + numCoveredPoints - The number of points in the meet
4275: - coveredPoints - The points in the meet
4277: Level: intermediate
4279: Fortran Notes:
4280: Since it returns an array, this routine is only available in Fortran 90, and you must
4281: include petsc.h90 in your code.
4283: The numCoveredPoints argument is not present in the Fortran 90 binding since it is internal to the array.
4285: .seealso: DMPlexGetMeet(), DMPlexGetFullMeet(), DMPlexGetJoin()
4286: @*/
4287: PetscErrorCode DMPlexRestoreMeet(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt **coveredPoints)
4288: {
4296: DMRestoreWorkArray(dm, 0, MPIU_INT, (void*) coveredPoints);
4297: if (numCoveredPoints) *numCoveredPoints = 0;
4298: return(0);
4299: }
4301: /*@C
4302: DMPlexGetFullMeet - Get an array for the meet of the set of points
4304: Not Collective
4306: Input Parameters:
4307: + dm - The DMPlex object
4308: . numPoints - The number of input points for the meet
4309: - points - The input points
4311: Output Parameters:
4312: + numCoveredPoints - The number of points in the meet
4313: - coveredPoints - The points in the meet
4315: Level: intermediate
4317: Fortran Notes:
4318: Since it returns an array, this routine is only available in Fortran 90, and you must
4319: include petsc.h90 in your code.
4321: The numCoveredPoints argument is not present in the Fortran 90 binding since it is internal to the array.
4323: .seealso: DMPlexGetMeet(), DMPlexRestoreMeet(), DMPlexGetJoin()
4324: @*/
4325: PetscErrorCode DMPlexGetFullMeet(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt **coveredPoints)
4326: {
4327: DM_Plex *mesh = (DM_Plex*) dm->data;
4328: PetscInt *offsets, **closures;
4329: PetscInt *meet[2];
4330: PetscInt height = 0, maxSize, meetSize = 0, i = 0;
4331: PetscInt p, h, c, m, mc;
4340: DMPlexGetDepth(dm, &height);
4341: PetscMalloc1(numPoints, &closures);
4342: DMGetWorkArray(dm, numPoints*(height+2), MPIU_INT, &offsets);
4343: mc = mesh->maxConeSize;
4344: maxSize = (mc > 1) ? ((PetscPowInt(mc,height+1)-1)/(mc-1)) : height + 1;
4345: DMGetWorkArray(dm, maxSize, MPIU_INT, &meet[0]);
4346: DMGetWorkArray(dm, maxSize, MPIU_INT, &meet[1]);
4348: for (p = 0; p < numPoints; ++p) {
4349: PetscInt closureSize;
4351: DMPlexGetTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closures[p]);
4353: offsets[p*(height+2)+0] = 0;
4354: for (h = 0; h < height+1; ++h) {
4355: PetscInt pStart, pEnd, i;
4357: DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);
4358: for (i = offsets[p*(height+2)+h]; i < closureSize; ++i) {
4359: if ((pStart > closures[p][i*2]) || (pEnd <= closures[p][i*2])) {
4360: offsets[p*(height+2)+h+1] = i;
4361: break;
4362: }
4363: }
4364: if (i == closureSize) offsets[p*(height+2)+h+1] = i;
4365: }
4366: if (offsets[p*(height+2)+height+1] != closureSize) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Total size of closure %D should be %D", offsets[p*(height+2)+height+1], closureSize);
4367: }
4368: for (h = 0; h < height+1; ++h) {
4369: PetscInt dof;
4371: /* Copy in cone of first point */
4372: dof = offsets[h+1] - offsets[h];
4373: for (meetSize = 0; meetSize < dof; ++meetSize) {
4374: meet[i][meetSize] = closures[0][(offsets[h]+meetSize)*2];
4375: }
4376: /* Check each successive cone */
4377: for (p = 1; p < numPoints && meetSize; ++p) {
4378: PetscInt newMeetSize = 0;
4380: dof = offsets[p*(height+2)+h+1] - offsets[p*(height+2)+h];
4381: for (c = 0; c < dof; ++c) {
4382: const PetscInt point = closures[p][(offsets[p*(height+2)+h]+c)*2];
4384: for (m = 0; m < meetSize; ++m) {
4385: if (point == meet[i][m]) {
4386: meet[1-i][newMeetSize++] = point;
4387: break;
4388: }
4389: }
4390: }
4391: meetSize = newMeetSize;
4392: i = 1-i;
4393: }
4394: if (meetSize) break;
4395: }
4396: *numCoveredPoints = meetSize;
4397: *coveredPoints = meet[i];
4398: for (p = 0; p < numPoints; ++p) {
4399: DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_TRUE, NULL, &closures[p]);
4400: }
4401: PetscFree(closures);
4402: DMRestoreWorkArray(dm, numPoints*(height+2), MPIU_INT, &offsets);
4403: DMRestoreWorkArray(dm, mesh->maxConeSize, MPIU_INT, &meet[1-i]);
4404: return(0);
4405: }
4407: /*@C
4408: DMPlexEqual - Determine if two DMs have the same topology
4410: Not Collective
4412: Input Parameters:
4413: + dmA - A DMPlex object
4414: - dmB - A DMPlex object
4416: Output Parameters:
4417: . equal - PETSC_TRUE if the topologies are identical
4419: Level: intermediate
4421: Notes:
4422: We are not solving graph isomorphism, so we do not permutation.
4424: .seealso: DMPlexGetCone()
4425: @*/
4426: PetscErrorCode DMPlexEqual(DM dmA, DM dmB, PetscBool *equal)
4427: {
4428: PetscInt depth, depthB, pStart, pEnd, pStartB, pEndB, p;
4436: *equal = PETSC_FALSE;
4437: DMPlexGetDepth(dmA, &depth);
4438: DMPlexGetDepth(dmB, &depthB);
4439: if (depth != depthB) return(0);
4440: DMPlexGetChart(dmA, &pStart, &pEnd);
4441: DMPlexGetChart(dmB, &pStartB, &pEndB);
4442: if ((pStart != pStartB) || (pEnd != pEndB)) return(0);
4443: for (p = pStart; p < pEnd; ++p) {
4444: const PetscInt *cone, *coneB, *ornt, *orntB, *support, *supportB;
4445: PetscInt coneSize, coneSizeB, c, supportSize, supportSizeB, s;
4447: DMPlexGetConeSize(dmA, p, &coneSize);
4448: DMPlexGetCone(dmA, p, &cone);
4449: DMPlexGetConeOrientation(dmA, p, &ornt);
4450: DMPlexGetConeSize(dmB, p, &coneSizeB);
4451: DMPlexGetCone(dmB, p, &coneB);
4452: DMPlexGetConeOrientation(dmB, p, &orntB);
4453: if (coneSize != coneSizeB) return(0);
4454: for (c = 0; c < coneSize; ++c) {
4455: if (cone[c] != coneB[c]) return(0);
4456: if (ornt[c] != orntB[c]) return(0);
4457: }
4458: DMPlexGetSupportSize(dmA, p, &supportSize);
4459: DMPlexGetSupport(dmA, p, &support);
4460: DMPlexGetSupportSize(dmB, p, &supportSizeB);
4461: DMPlexGetSupport(dmB, p, &supportB);
4462: if (supportSize != supportSizeB) return(0);
4463: for (s = 0; s < supportSize; ++s) {
4464: if (support[s] != supportB[s]) return(0);
4465: }
4466: }
4467: *equal = PETSC_TRUE;
4468: return(0);
4469: }
4471: /*@C
4472: DMPlexGetNumFaceVertices - Returns the number of vertices on a face
4474: Not Collective
4476: Input Parameters:
4477: + dm - The DMPlex
4478: . cellDim - The cell dimension
4479: - numCorners - The number of vertices on a cell
4481: Output Parameters:
4482: . numFaceVertices - The number of vertices on a face
4484: Level: developer
4486: Notes:
4487: Of course this can only work for a restricted set of symmetric shapes
4489: .seealso: DMPlexGetCone()
4490: @*/
4491: PetscErrorCode DMPlexGetNumFaceVertices(DM dm, PetscInt cellDim, PetscInt numCorners, PetscInt *numFaceVertices)
4492: {
4493: MPI_Comm comm;
4497: PetscObjectGetComm((PetscObject)dm,&comm);
4499: switch (cellDim) {
4500: case 0:
4501: *numFaceVertices = 0;
4502: break;
4503: case 1:
4504: *numFaceVertices = 1;
4505: break;
4506: case 2:
4507: switch (numCorners) {
4508: case 3: /* triangle */
4509: *numFaceVertices = 2; /* Edge has 2 vertices */
4510: break;
4511: case 4: /* quadrilateral */
4512: *numFaceVertices = 2; /* Edge has 2 vertices */
4513: break;
4514: case 6: /* quadratic triangle, tri and quad cohesive Lagrange cells */
4515: *numFaceVertices = 3; /* Edge has 3 vertices */
4516: break;
4517: case 9: /* quadratic quadrilateral, quadratic quad cohesive Lagrange cells */
4518: *numFaceVertices = 3; /* Edge has 3 vertices */
4519: break;
4520: default:
4521: SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of face corners %D for dimension %D", numCorners, cellDim);
4522: }
4523: break;
4524: case 3:
4525: switch (numCorners) {
4526: case 4: /* tetradehdron */
4527: *numFaceVertices = 3; /* Face has 3 vertices */
4528: break;
4529: case 6: /* tet cohesive cells */
4530: *numFaceVertices = 4; /* Face has 4 vertices */
4531: break;
4532: case 8: /* hexahedron */
4533: *numFaceVertices = 4; /* Face has 4 vertices */
4534: break;
4535: case 9: /* tet cohesive Lagrange cells */
4536: *numFaceVertices = 6; /* Face has 6 vertices */
4537: break;
4538: case 10: /* quadratic tetrahedron */
4539: *numFaceVertices = 6; /* Face has 6 vertices */
4540: break;
4541: case 12: /* hex cohesive Lagrange cells */
4542: *numFaceVertices = 6; /* Face has 6 vertices */
4543: break;
4544: case 18: /* quadratic tet cohesive Lagrange cells */
4545: *numFaceVertices = 6; /* Face has 6 vertices */
4546: break;
4547: case 27: /* quadratic hexahedron, quadratic hex cohesive Lagrange cells */
4548: *numFaceVertices = 9; /* Face has 9 vertices */
4549: break;
4550: default:
4551: SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of face corners %D for dimension %D", numCorners, cellDim);
4552: }
4553: break;
4554: default:
4555: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid cell dimension %D", cellDim);
4556: }
4557: return(0);
4558: }
4560: /*@
4561: DMPlexGetDepthLabel - Get the DMLabel recording the depth of each point
4563: Not Collective
4565: Input Parameter:
4566: . dm - The DMPlex object
4568: Output Parameter:
4569: . depthLabel - The DMLabel recording point depth
4571: Level: developer
4573: .seealso: DMPlexGetDepth(), DMPlexGetHeightStratum(), DMPlexGetDepthStratum(), DMPlexGetPointDepth(),
4574: @*/
4575: PetscErrorCode DMPlexGetDepthLabel(DM dm, DMLabel *depthLabel)
4576: {
4580: *depthLabel = dm->depthLabel;
4581: return(0);
4582: }
4584: /*@
4585: DMPlexGetDepth - Get the depth of the DAG representing this mesh
4587: Not Collective
4589: Input Parameter:
4590: . dm - The DMPlex object
4592: Output Parameter:
4593: . depth - The number of strata (breadth first levels) in the DAG
4595: Level: developer
4597: Notes:
4598: This returns maximum of point depths over all points, i.e. maximum value of the label returned by DMPlexGetDepthLabel().
4599: The point depth is described more in detail in DMPlexGetDepthStratum().
4600: An empty mesh gives -1.
4602: .seealso: DMPlexGetDepthLabel(), DMPlexGetDepthStratum(), DMPlexGetPointDepth(), DMPlexSymmetrize()
4603: @*/
4604: PetscErrorCode DMPlexGetDepth(DM dm, PetscInt *depth)
4605: {
4606: DMLabel label;
4607: PetscInt d = 0;
4613: DMPlexGetDepthLabel(dm, &label);
4614: if (label) {DMLabelGetNumValues(label, &d);}
4615: *depth = d-1;
4616: return(0);
4617: }
4619: /*@
4620: DMPlexGetDepthStratum - Get the bounds [start, end) for all points at a certain depth.
4622: Not Collective
4624: Input Parameters:
4625: + dm - The DMPlex object
4626: - stratumValue - The requested depth
4628: Output Parameters:
4629: + start - The first point at this depth
4630: - end - One beyond the last point at this depth
4632: Notes:
4633: Depth indexing is related to topological dimension. Depth stratum 0 contains the lowest topological dimension points,
4634: often "vertices". If the mesh is "interpolated" (see DMPlexInterpolate()), then depth stratum 1 contains the next
4635: higher dimension, e.g., "edges".
4637: Level: developer
4639: .seealso: DMPlexGetHeightStratum(), DMPlexGetDepth(), DMPlexGetDepthLabel(), DMPlexGetPointDepth(), DMPlexSymmetrize(), DMPlexInterpolate()
4640: @*/
4641: PetscErrorCode DMPlexGetDepthStratum(DM dm, PetscInt stratumValue, PetscInt *start, PetscInt *end)
4642: {
4643: DMLabel label;
4644: PetscInt pStart, pEnd;
4651: DMPlexGetChart(dm, &pStart, &pEnd);
4652: if (pStart == pEnd) return(0);
4653: if (stratumValue < 0) {
4654: if (start) *start = pStart;
4655: if (end) *end = pEnd;
4656: return(0);
4657: }
4658: DMPlexGetDepthLabel(dm, &label);
4659: if (!label) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "No label named depth was found");
4660: DMLabelGetStratumBounds(label, stratumValue, start, end);
4661: return(0);
4662: }
4664: /*@
4665: DMPlexGetHeightStratum - Get the bounds [start, end) for all points at a certain height.
4667: Not Collective
4669: Input Parameters:
4670: + dm - The DMPlex object
4671: - stratumValue - The requested height
4673: Output Parameters:
4674: + start - The first point at this height
4675: - end - One beyond the last point at this height
4677: Notes:
4678: Height indexing is related to topological codimension. Height stratum 0 contains the highest topological dimension
4679: points, often called "cells" or "elements". If the mesh is "interpolated" (see DMPlexInterpolate()), then height
4680: stratum 1 contains the boundary of these "cells", often called "faces" or "facets".
4682: Level: developer
4684: .seealso: DMPlexGetDepthStratum(), DMPlexGetDepth(), DMPlexGetPointHeight()
4685: @*/
4686: PetscErrorCode DMPlexGetHeightStratum(DM dm, PetscInt stratumValue, PetscInt *start, PetscInt *end)
4687: {
4688: DMLabel label;
4689: PetscInt depth, pStart, pEnd;
4696: DMPlexGetChart(dm, &pStart, &pEnd);
4697: if (pStart == pEnd) return(0);
4698: if (stratumValue < 0) {
4699: if (start) *start = pStart;
4700: if (end) *end = pEnd;
4701: return(0);
4702: }
4703: DMPlexGetDepthLabel(dm, &label);
4704: if (!label) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "No label named depth was found");
4705: DMLabelGetNumValues(label, &depth);
4706: DMLabelGetStratumBounds(label, depth-1-stratumValue, start, end);
4707: return(0);
4708: }
4710: /*@
4711: DMPlexGetPointDepth - Get the depth of a given point
4713: Not Collective
4715: Input Parameters:
4716: + dm - The DMPlex object
4717: - point - The point
4719: Output Parameter:
4720: . depth - The depth of the point
4722: Level: intermediate
4724: .seealso: DMPlexGetCellType(), DMPlexGetDepthLabel(), DMPlexGetDepth(), DMPlexGetPointHeight()
4725: @*/
4726: PetscErrorCode DMPlexGetPointDepth(DM dm, PetscInt point, PetscInt *depth)
4727: {
4733: DMLabelGetValue(dm->depthLabel, point, depth);
4734: return(0);
4735: }
4737: /*@
4738: DMPlexGetPointHeight - Get the height of a given point
4740: Not Collective
4742: Input Parameters:
4743: + dm - The DMPlex object
4744: - point - The point
4746: Output Parameter:
4747: . height - The height of the point
4749: Level: intermediate
4751: .seealso: DMPlexGetCellType(), DMPlexGetDepthLabel(), DMPlexGetDepth(), DMPlexGetPointDepth()
4752: @*/
4753: PetscErrorCode DMPlexGetPointHeight(DM dm, PetscInt point, PetscInt *height)
4754: {
4755: PetscInt n, pDepth;
4761: DMLabelGetNumValues(dm->depthLabel, &n);
4762: DMLabelGetValue(dm->depthLabel, point, &pDepth);
4763: *height = n - 1 - pDepth; /* DAG depth is n-1 */
4764: return(0);
4765: }
4767: /*@
4768: DMPlexGetCellTypeLabel - Get the DMLabel recording the polytope type of each cell
4770: Not Collective
4772: Input Parameter:
4773: . dm - The DMPlex object
4775: Output Parameter:
4776: . celltypeLabel - The DMLabel recording cell polytope type
4778: Note: This function will trigger automatica computation of cell types. This can be disabled by calling
4779: DMCreateLabel(dm, "celltype") beforehand.
4781: Level: developer
4783: .seealso: DMPlexGetCellType(), DMPlexGetDepthLabel(), DMCreateLabel()
4784: @*/
4785: PetscErrorCode DMPlexGetCellTypeLabel(DM dm, DMLabel *celltypeLabel)
4786: {
4792: if (!dm->celltypeLabel) {DMPlexComputeCellTypes(dm);}
4793: *celltypeLabel = dm->celltypeLabel;
4794: return(0);
4795: }
4797: /*@
4798: DMPlexGetCellType - Get the polytope type of a given cell
4800: Not Collective
4802: Input Parameters:
4803: + dm - The DMPlex object
4804: - cell - The cell
4806: Output Parameter:
4807: . celltype - The polytope type of the cell
4809: Level: intermediate
4811: .seealso: DMPlexGetCellTypeLabel(), DMPlexGetDepthLabel(), DMPlexGetDepth()
4812: @*/
4813: PetscErrorCode DMPlexGetCellType(DM dm, PetscInt cell, DMPolytopeType *celltype)
4814: {
4815: DMLabel label;
4816: PetscInt ct;
4822: DMPlexGetCellTypeLabel(dm, &label);
4823: DMLabelGetValue(label, cell, &ct);
4824: if (ct < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cell %D has not been assigned a cell type", cell);
4825: *celltype = (DMPolytopeType) ct;
4826: return(0);
4827: }
4829: /*@
4830: DMPlexSetCellType - Set the polytope type of a given cell
4832: Not Collective
4834: Input Parameters:
4835: + dm - The DMPlex object
4836: . cell - The cell
4837: - celltype - The polytope type of the cell
4839: Note: By default, cell types will be automatically computed using DMPlexComputeCellTypes() before this function
4840: is executed. This function will override the computed type. However, if automatic classification will not succeed
4841: and a user wants to manually specify all types, the classification must be disabled by calling
4842: DMCreaateLabel(dm, "celltype") before getting or setting any cell types.
4844: Level: advanced
4846: .seealso: DMPlexGetCellTypeLabel(), DMPlexGetDepthLabel(), DMPlexGetDepth(), DMPlexComputeCellTypes(), DMCreateLabel()
4847: @*/
4848: PetscErrorCode DMPlexSetCellType(DM dm, PetscInt cell, DMPolytopeType celltype)
4849: {
4850: DMLabel label;
4855: DMPlexGetCellTypeLabel(dm, &label);
4856: DMLabelSetValue(label, cell, celltype);
4857: return(0);
4858: }
4860: PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm)
4861: {
4862: PetscSection section, s;
4863: Mat m;
4864: PetscInt maxHeight;
4868: DMClone(dm, cdm);
4869: DMPlexGetMaxProjectionHeight(dm, &maxHeight);
4870: DMPlexSetMaxProjectionHeight(*cdm, maxHeight);
4871: PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion);
4872: DMSetLocalSection(*cdm, section);
4873: PetscSectionDestroy(§ion);
4874: PetscSectionCreate(PETSC_COMM_SELF, &s);
4875: MatCreate(PETSC_COMM_SELF, &m);
4876: DMSetDefaultConstraints(*cdm, s, m);
4877: PetscSectionDestroy(&s);
4878: MatDestroy(&m);
4880: DMSetNumFields(*cdm, 1);
4881: DMCreateDS(*cdm);
4882: return(0);
4883: }
4885: PetscErrorCode DMCreateCoordinateField_Plex(DM dm, DMField *field)
4886: {
4887: Vec coordsLocal;
4888: DM coordsDM;
4892: *field = NULL;
4893: DMGetCoordinatesLocal(dm,&coordsLocal);
4894: DMGetCoordinateDM(dm,&coordsDM);
4895: if (coordsLocal && coordsDM) {
4896: DMFieldCreateDS(coordsDM, 0, coordsLocal, field);
4897: }
4898: return(0);
4899: }
4901: /*@C
4902: DMPlexGetConeSection - Return a section which describes the layout of cone data
4904: Not Collective
4906: Input Parameters:
4907: . dm - The DMPlex object
4909: Output Parameter:
4910: . section - The PetscSection object
4912: Level: developer
4914: .seealso: DMPlexGetSupportSection(), DMPlexGetCones(), DMPlexGetConeOrientations()
4915: @*/
4916: PetscErrorCode DMPlexGetConeSection(DM dm, PetscSection *section)
4917: {
4918: DM_Plex *mesh = (DM_Plex*) dm->data;
4922: if (section) *section = mesh->coneSection;
4923: return(0);
4924: }
4926: /*@C
4927: DMPlexGetSupportSection - Return a section which describes the layout of support data
4929: Not Collective
4931: Input Parameters:
4932: . dm - The DMPlex object
4934: Output Parameter:
4935: . section - The PetscSection object
4937: Level: developer
4939: .seealso: DMPlexGetConeSection()
4940: @*/
4941: PetscErrorCode DMPlexGetSupportSection(DM dm, PetscSection *section)
4942: {
4943: DM_Plex *mesh = (DM_Plex*) dm->data;
4947: if (section) *section = mesh->supportSection;
4948: return(0);
4949: }
4951: /*@C
4952: DMPlexGetCones - Return cone data
4954: Not Collective
4956: Input Parameters:
4957: . dm - The DMPlex object
4959: Output Parameter:
4960: . cones - The cone for each point
4962: Level: developer
4964: .seealso: DMPlexGetConeSection()
4965: @*/
4966: PetscErrorCode DMPlexGetCones(DM dm, PetscInt *cones[])
4967: {
4968: DM_Plex *mesh = (DM_Plex*) dm->data;
4972: if (cones) *cones = mesh->cones;
4973: return(0);
4974: }
4976: /*@C
4977: DMPlexGetConeOrientations - Return cone orientation data
4979: Not Collective
4981: Input Parameters:
4982: . dm - The DMPlex object
4984: Output Parameter:
4985: . coneOrientations - The array of cone orientations for all points
4987: Level: developer
4989: Notes:
4990: The PetscSection returned by DMPlexGetConeSection() partitions coneOrientations into cone orientations of particular points as returned by DMPlexGetConeOrientation().
4992: The meaning of coneOrientations values is detailed in DMPlexGetConeOrientation().
4994: .seealso: DMPlexGetConeSection(), DMPlexGetConeOrientation()
4995: @*/
4996: PetscErrorCode DMPlexGetConeOrientations(DM dm, PetscInt *coneOrientations[])
4997: {
4998: DM_Plex *mesh = (DM_Plex*) dm->data;
5002: if (coneOrientations) *coneOrientations = mesh->coneOrientations;
5003: return(0);
5004: }
5006: /******************************** FEM Support **********************************/
5008: /*
5009: Returns number of components and tensor degree for the field. For interpolated meshes, line should be a point
5010: representing a line in the section.
5011: */
5012: static PetscErrorCode PetscSectionFieldGetTensorDegree_Private(PetscSection section,PetscInt field,PetscInt line,PetscBool vertexchart,PetscInt *Nc,PetscInt *k)
5013: {
5017: PetscSectionGetFieldComponents(section, field, Nc);
5018: if (line < 0) {
5019: *k = 0;
5020: *Nc = 0;
5021: } else if (vertexchart) { /* If we only have a vertex chart, we must have degree k=1 */
5022: *k = 1;
5023: } else { /* Assume the full interpolated mesh is in the chart; lines in particular */
5024: /* An order k SEM disc has k-1 dofs on an edge */
5025: PetscSectionGetFieldDof(section, line, field, k);
5026: *k = *k / *Nc + 1;
5027: }
5028: return(0);
5029: }
5031: /*@
5033: DMPlexSetClosurePermutationTensor - Create a permutation from the default (BFS) point ordering in the closure, to a
5034: lexicographic ordering over the tensor product cell (i.e., line, quad, hex, etc.), and set this permutation in the
5035: section provided (or the section of the DM).
5037: Input Parameters:
5038: + dm - The DM
5039: . point - Either a cell (highest dim point) or an edge (dim 1 point), or PETSC_DETERMINE
5040: - section - The PetscSection to reorder, or NULL for the default section
5042: Note: The point is used to determine the number of dofs/field on an edge. For SEM, this is related to the polynomial
5043: degree of the basis.
5045: Example:
5046: A typical interpolated single-quad mesh might order points as
5047: .vb
5048: [c0, v1, v2, v3, v4, e5, e6, e7, e8]
5050: v4 -- e6 -- v3
5051: | |
5052: e7 c0 e8
5053: | |
5054: v1 -- e5 -- v2
5055: .ve
5057: (There is no significance to the ordering described here.) The default section for a Q3 quad might typically assign
5058: dofs in the order of points, e.g.,
5059: .vb
5060: c0 -> [0,1,2,3]
5061: v1 -> [4]
5062: ...
5063: e5 -> [8, 9]
5064: .ve
5066: which corresponds to the dofs
5067: .vb
5068: 6 10 11 7
5069: 13 2 3 15
5070: 12 0 1 14
5071: 4 8 9 5
5072: .ve
5074: The closure in BFS ordering works through height strata (cells, edges, vertices) to produce the ordering
5075: .vb
5076: 0 1 2 3 8 9 14 15 11 10 13 12 4 5 7 6
5077: .ve
5079: After calling DMPlexSetClosurePermutationTensor(), the closure will be ordered lexicographically,
5080: .vb
5081: 4 8 9 5 12 0 1 14 13 2 3 15 6 10 11 7
5082: .ve
5084: Level: developer
5086: .seealso: DMGetLocalSection(), PetscSectionSetClosurePermutation(), DMSetGlobalSection()
5087: @*/
5088: PetscErrorCode DMPlexSetClosurePermutationTensor(DM dm, PetscInt point, PetscSection section)
5089: {
5090: DMLabel label;
5091: PetscInt dim, depth = -1, eStart = -1, Nf;
5092: PetscBool vertexchart;
5096: DMGetDimension(dm, &dim);
5097: if (dim < 1) return(0);
5098: if (point < 0) {
5099: PetscInt sStart,sEnd;
5101: DMPlexGetDepthStratum(dm, 1, &sStart, &sEnd);
5102: point = sEnd-sStart ? sStart : point;
5103: }
5104: DMPlexGetDepthLabel(dm, &label);
5105: if (point >= 0) { DMLabelGetValue(label, point, &depth); }
5106: if (!section) {DMGetLocalSection(dm, §ion);}
5107: if (depth == 1) {eStart = point;}
5108: else if (depth == dim) {
5109: const PetscInt *cone;
5111: DMPlexGetCone(dm, point, &cone);
5112: if (dim == 2) eStart = cone[0];
5113: else if (dim == 3) {
5114: const PetscInt *cone2;
5115: DMPlexGetCone(dm, cone[0], &cone2);
5116: eStart = cone2[0];
5117: } else SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %D of depth %D cannot be used to bootstrap spectral ordering for dim %D", point, depth, dim);
5118: } else if (depth >= 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %D of depth %D cannot be used to bootstrap spectral ordering for dim %D", point, depth, dim);
5119: { /* Determine whether the chart covers all points or just vertices. */
5120: PetscInt pStart,pEnd,cStart,cEnd;
5121: DMPlexGetDepthStratum(dm,0,&pStart,&pEnd);
5122: PetscSectionGetChart(section,&cStart,&cEnd);
5123: if (pStart == cStart && pEnd == cEnd) vertexchart = PETSC_TRUE; /* Just vertices */
5124: else vertexchart = PETSC_FALSE; /* Assume all interpolated points are in chart */
5125: }
5126: PetscSectionGetNumFields(section, &Nf);
5127: for (PetscInt d=1; d<=dim; d++) {
5128: PetscInt k, f, Nc, c, i, j, size = 0, offset = 0, foffset = 0;
5129: PetscInt *perm;
5131: for (f = 0; f < Nf; ++f) {
5132: PetscSectionFieldGetTensorDegree_Private(section,f,eStart,vertexchart,&Nc,&k);
5133: size += PetscPowInt(k+1, d)*Nc;
5134: }
5135: PetscMalloc1(size, &perm);
5136: for (f = 0; f < Nf; ++f) {
5137: switch (d) {
5138: case 1:
5139: PetscSectionFieldGetTensorDegree_Private(section,f,eStart,vertexchart,&Nc,&k);
5140: /*
5141: Original ordering is [ edge of length k-1; vtx0; vtx1 ]
5142: We want [ vtx0; edge of length k-1; vtx1 ]
5143: */
5144: for (c=0; c<Nc; c++,offset++) perm[offset] = (k-1)*Nc + c + foffset;
5145: for (i=0; i<k-1; i++) for (c=0; c<Nc; c++,offset++) perm[offset] = i*Nc + c + foffset;
5146: for (c=0; c<Nc; c++,offset++) perm[offset] = k*Nc + c + foffset;
5147: foffset = offset;
5148: break;
5149: case 2:
5150: /* The original quad closure is oriented clockwise, {f, e_b, e_r, e_t, e_l, v_lb, v_rb, v_tr, v_tl} */
5151: PetscSectionFieldGetTensorDegree_Private(section,f,eStart,vertexchart,&Nc,&k);
5152: /* The SEM order is
5154: v_lb, {e_b}, v_rb,
5155: e^{(k-1)-i}_l, {f^{i*(k-1)}}, e^i_r,
5156: v_lt, reverse {e_t}, v_rt
5157: */
5158: {
5159: const PetscInt of = 0;
5160: const PetscInt oeb = of + PetscSqr(k-1);
5161: const PetscInt oer = oeb + (k-1);
5162: const PetscInt oet = oer + (k-1);
5163: const PetscInt oel = oet + (k-1);
5164: const PetscInt ovlb = oel + (k-1);
5165: const PetscInt ovrb = ovlb + 1;
5166: const PetscInt ovrt = ovrb + 1;
5167: const PetscInt ovlt = ovrt + 1;
5168: PetscInt o;
5170: /* bottom */
5171: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovlb*Nc + c + foffset;
5172: for (o = oeb; o < oer; ++o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
5173: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovrb*Nc + c + foffset;
5174: /* middle */
5175: for (i = 0; i < k-1; ++i) {
5176: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oel+(k-2)-i)*Nc + c + foffset;
5177: for (o = of+(k-1)*i; o < of+(k-1)*(i+1); ++o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
5178: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oer+i)*Nc + c + foffset;
5179: }
5180: /* top */
5181: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovlt*Nc + c + foffset;
5182: for (o = oel-1; o >= oet; --o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
5183: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovrt*Nc + c + foffset;
5184: foffset = offset;
5185: }
5186: break;
5187: case 3:
5188: /* The original hex closure is
5190: {c,
5191: f_b, f_t, f_f, f_b, f_r, f_l,
5192: e_bl, e_bb, e_br, e_bf, e_tf, e_tr, e_tb, e_tl, e_rf, e_lf, e_lb, e_rb,
5193: v_blf, v_blb, v_brb, v_brf, v_tlf, v_trf, v_trb, v_tlb}
5194: */
5195: PetscSectionFieldGetTensorDegree_Private(section,f,eStart,vertexchart,&Nc,&k);
5196: /* The SEM order is
5197: Bottom Slice
5198: v_blf, {e^{(k-1)-n}_bf}, v_brf,
5199: e^{i}_bl, f^{n*(k-1)+(k-1)-i}_b, e^{(k-1)-i}_br,
5200: v_blb, {e_bb}, v_brb,
5202: Middle Slice (j)
5203: {e^{(k-1)-j}_lf}, {f^{j*(k-1)+n}_f}, e^j_rf,
5204: f^{i*(k-1)+j}_l, {c^{(j*(k-1) + i)*(k-1)+n}_t}, f^{j*(k-1)+i}_r,
5205: e^j_lb, {f^{j*(k-1)+(k-1)-n}_b}, e^{(k-1)-j}_rb,
5207: Top Slice
5208: v_tlf, {e_tf}, v_trf,
5209: e^{(k-1)-i}_tl, {f^{i*(k-1)}_t}, e^{i}_tr,
5210: v_tlb, {e^{(k-1)-n}_tb}, v_trb,
5211: */
5212: {
5213: const PetscInt oc = 0;
5214: const PetscInt ofb = oc + PetscSqr(k-1)*(k-1);
5215: const PetscInt oft = ofb + PetscSqr(k-1);
5216: const PetscInt off = oft + PetscSqr(k-1);
5217: const PetscInt ofk = off + PetscSqr(k-1);
5218: const PetscInt ofr = ofk + PetscSqr(k-1);
5219: const PetscInt ofl = ofr + PetscSqr(k-1);
5220: const PetscInt oebl = ofl + PetscSqr(k-1);
5221: const PetscInt oebb = oebl + (k-1);
5222: const PetscInt oebr = oebb + (k-1);
5223: const PetscInt oebf = oebr + (k-1);
5224: const PetscInt oetf = oebf + (k-1);
5225: const PetscInt oetr = oetf + (k-1);
5226: const PetscInt oetb = oetr + (k-1);
5227: const PetscInt oetl = oetb + (k-1);
5228: const PetscInt oerf = oetl + (k-1);
5229: const PetscInt oelf = oerf + (k-1);
5230: const PetscInt oelb = oelf + (k-1);
5231: const PetscInt oerb = oelb + (k-1);
5232: const PetscInt ovblf = oerb + (k-1);
5233: const PetscInt ovblb = ovblf + 1;
5234: const PetscInt ovbrb = ovblb + 1;
5235: const PetscInt ovbrf = ovbrb + 1;
5236: const PetscInt ovtlf = ovbrf + 1;
5237: const PetscInt ovtrf = ovtlf + 1;
5238: const PetscInt ovtrb = ovtrf + 1;
5239: const PetscInt ovtlb = ovtrb + 1;
5240: PetscInt o, n;
5242: /* Bottom Slice */
5243: /* bottom */
5244: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovblf*Nc + c + foffset;
5245: for (o = oetf-1; o >= oebf; --o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
5246: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovbrf*Nc + c + foffset;
5247: /* middle */
5248: for (i = 0; i < k-1; ++i) {
5249: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oebl+i)*Nc + c + foffset;
5250: for (n = 0; n < k-1; ++n) {o = ofb+n*(k-1)+i; for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;}
5251: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oebr+(k-2)-i)*Nc + c + foffset;
5252: }
5253: /* top */
5254: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovblb*Nc + c + foffset;
5255: for (o = oebb; o < oebr; ++o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
5256: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovbrb*Nc + c + foffset;
5258: /* Middle Slice */
5259: for (j = 0; j < k-1; ++j) {
5260: /* bottom */
5261: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oelf+(k-2)-j)*Nc + c + foffset;
5262: for (o = off+j*(k-1); o < off+(j+1)*(k-1); ++o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
5263: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oerf+j)*Nc + c + foffset;
5264: /* middle */
5265: for (i = 0; i < k-1; ++i) {
5266: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (ofl+i*(k-1)+j)*Nc + c + foffset;
5267: for (n = 0; n < k-1; ++n) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oc+(j*(k-1)+i)*(k-1)+n)*Nc + c + foffset;
5268: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (ofr+j*(k-1)+i)*Nc + c + foffset;
5269: }
5270: /* top */
5271: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oelb+j)*Nc + c + foffset;
5272: for (o = ofk+j*(k-1)+(k-2); o >= ofk+j*(k-1); --o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
5273: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oerb+(k-2)-j)*Nc + c + foffset;
5274: }
5276: /* Top Slice */
5277: /* bottom */
5278: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtlf*Nc + c + foffset;
5279: for (o = oetf; o < oetr; ++o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
5280: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtrf*Nc + c + foffset;
5281: /* middle */
5282: for (i = 0; i < k-1; ++i) {
5283: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oetl+(k-2)-i)*Nc + c + foffset;
5284: for (n = 0; n < k-1; ++n) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oft+i*(k-1)+n)*Nc + c + foffset;
5285: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oetr+i)*Nc + c + foffset;
5286: }
5287: /* top */
5288: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtlb*Nc + c + foffset;
5289: for (o = oetl-1; o >= oetb; --o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
5290: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtrb*Nc + c + foffset;
5292: foffset = offset;
5293: }
5294: break;
5295: default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No spectral ordering for dimension %D", d);
5296: }
5297: }
5298: if (offset != size) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Number of permutation entries %D != %D", offset, size);
5299: /* Check permutation */
5300: {
5301: PetscInt *check;
5303: PetscMalloc1(size, &check);
5304: for (i = 0; i < size; ++i) {check[i] = -1; if (perm[i] < 0 || perm[i] >= size) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Invalid permutation index p[%D] = %D", i, perm[i]);}
5305: for (i = 0; i < size; ++i) check[perm[i]] = i;
5306: for (i = 0; i < size; ++i) {if (check[i] < 0) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Missing permutation index %D", i);}
5307: PetscFree(check);
5308: }
5309: PetscSectionSetClosurePermutation_Internal(section, (PetscObject) dm, d, size, PETSC_OWN_POINTER, perm);
5310: }
5311: return(0);
5312: }
5314: PetscErrorCode DMPlexGetPointDualSpaceFEM(DM dm, PetscInt point, PetscInt field, PetscDualSpace *dspace)
5315: {
5316: PetscDS prob;
5317: PetscInt depth, Nf, h;
5318: DMLabel label;
5322: DMGetDS(dm, &prob);
5323: Nf = prob->Nf;
5324: label = dm->depthLabel;
5325: *dspace = NULL;
5326: if (field < Nf) {
5327: PetscObject disc = prob->disc[field];
5329: if (disc->classid == PETSCFE_CLASSID) {
5330: PetscDualSpace dsp;
5332: PetscFEGetDualSpace((PetscFE)disc,&dsp);
5333: DMLabelGetNumValues(label,&depth);
5334: DMLabelGetValue(label,point,&h);
5335: h = depth - 1 - h;
5336: if (h) {
5337: PetscDualSpaceGetHeightSubspace(dsp,h,dspace);
5338: } else {
5339: *dspace = dsp;
5340: }
5341: }
5342: }
5343: return(0);
5344: }
5346: PETSC_STATIC_INLINE PetscErrorCode DMPlexVecGetClosure_Depth1_Static(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt *csize, PetscScalar *values[])
5347: {
5348: PetscScalar *array, *vArray;
5349: const PetscInt *cone, *coneO;
5350: PetscInt pStart, pEnd, p, numPoints, size = 0, offset = 0;
5351: PetscErrorCode ierr;
5354: PetscSectionGetChart(section, &pStart, &pEnd);
5355: DMPlexGetConeSize(dm, point, &numPoints);
5356: DMPlexGetCone(dm, point, &cone);
5357: DMPlexGetConeOrientation(dm, point, &coneO);
5358: if (!values || !*values) {
5359: if ((point >= pStart) && (point < pEnd)) {
5360: PetscInt dof;
5362: PetscSectionGetDof(section, point, &dof);
5363: size += dof;
5364: }
5365: for (p = 0; p < numPoints; ++p) {
5366: const PetscInt cp = cone[p];
5367: PetscInt dof;
5369: if ((cp < pStart) || (cp >= pEnd)) continue;
5370: PetscSectionGetDof(section, cp, &dof);
5371: size += dof;
5372: }
5373: if (!values) {
5374: if (csize) *csize = size;
5375: return(0);
5376: }
5377: DMGetWorkArray(dm, size, MPIU_SCALAR, &array);
5378: } else {
5379: array = *values;
5380: }
5381: size = 0;
5382: VecGetArray(v, &vArray);
5383: if ((point >= pStart) && (point < pEnd)) {
5384: PetscInt dof, off, d;
5385: PetscScalar *varr;
5387: PetscSectionGetDof(section, point, &dof);
5388: PetscSectionGetOffset(section, point, &off);
5389: varr = &vArray[off];
5390: for (d = 0; d < dof; ++d, ++offset) {
5391: array[offset] = varr[d];
5392: }
5393: size += dof;
5394: }
5395: for (p = 0; p < numPoints; ++p) {
5396: const PetscInt cp = cone[p];
5397: PetscInt o = coneO[p];
5398: PetscInt dof, off, d;
5399: PetscScalar *varr;
5401: if ((cp < pStart) || (cp >= pEnd)) continue;
5402: PetscSectionGetDof(section, cp, &dof);
5403: PetscSectionGetOffset(section, cp, &off);
5404: varr = &vArray[off];
5405: if (o >= 0) {
5406: for (d = 0; d < dof; ++d, ++offset) {
5407: array[offset] = varr[d];
5408: }
5409: } else {
5410: for (d = dof-1; d >= 0; --d, ++offset) {
5411: array[offset] = varr[d];
5412: }
5413: }
5414: size += dof;
5415: }
5416: VecRestoreArray(v, &vArray);
5417: if (!*values) {
5418: if (csize) *csize = size;
5419: *values = array;
5420: } else {
5421: if (size > *csize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Size of input array %D < actual size %D", *csize, size);
5422: *csize = size;
5423: }
5424: return(0);
5425: }
5427: /* Compress out points not in the section */
5428: PETSC_STATIC_INLINE PetscErrorCode CompressPoints_Private(PetscSection section, PetscInt *numPoints, PetscInt points[])
5429: {
5430: const PetscInt np = *numPoints;
5431: PetscInt pStart, pEnd, p, q;
5434: PetscSectionGetChart(section, &pStart, &pEnd);
5435: for (p = 0, q = 0; p < np; ++p) {
5436: const PetscInt r = points[p*2];
5437: if ((r >= pStart) && (r < pEnd)) {
5438: points[q*2] = r;
5439: points[q*2+1] = points[p*2+1];
5440: ++q;
5441: }
5442: }
5443: *numPoints = q;
5444: return 0;
5445: }
5447: /* Compressed closure does not apply closure permutation */
5448: PetscErrorCode DMPlexGetCompressedClosure(DM dm, PetscSection section, PetscInt point, PetscInt *numPoints, PetscInt **points, PetscSection *clSec, IS *clPoints, const PetscInt **clp)
5449: {
5450: const PetscInt *cla = NULL;
5451: PetscInt np, *pts = NULL;
5455: PetscSectionGetClosureIndex(section, (PetscObject) dm, clSec, clPoints);
5456: if (*clPoints) {
5457: PetscInt dof, off;
5459: PetscSectionGetDof(*clSec, point, &dof);
5460: PetscSectionGetOffset(*clSec, point, &off);
5461: ISGetIndices(*clPoints, &cla);
5462: np = dof/2;
5463: pts = (PetscInt *) &cla[off];
5464: } else {
5465: DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &np, &pts);
5466: CompressPoints_Private(section, &np, pts);
5467: }
5468: *numPoints = np;
5469: *points = pts;
5470: *clp = cla;
5471: return(0);
5472: }
5474: PetscErrorCode DMPlexRestoreCompressedClosure(DM dm, PetscSection section, PetscInt point, PetscInt *numPoints, PetscInt **points, PetscSection *clSec, IS *clPoints, const PetscInt **clp)
5475: {
5479: if (!*clPoints) {
5480: DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, numPoints, points);
5481: } else {
5482: ISRestoreIndices(*clPoints, clp);
5483: }
5484: *numPoints = 0;
5485: *points = NULL;
5486: *clSec = NULL;
5487: *clPoints = NULL;
5488: *clp = NULL;
5489: return(0);
5490: }
5492: PETSC_STATIC_INLINE PetscErrorCode DMPlexVecGetClosure_Static(DM dm, PetscSection section, PetscInt numPoints, const PetscInt points[], const PetscInt clperm[], const PetscScalar vArray[], PetscInt *size, PetscScalar array[])
5493: {
5494: PetscInt offset = 0, p;
5495: const PetscInt **perms = NULL;
5496: const PetscScalar **flips = NULL;
5497: PetscErrorCode ierr;
5500: *size = 0;
5501: PetscSectionGetPointSyms(section,numPoints,points,&perms,&flips);
5502: for (p = 0; p < numPoints; p++) {
5503: const PetscInt point = points[2*p];
5504: const PetscInt *perm = perms ? perms[p] : NULL;
5505: const PetscScalar *flip = flips ? flips[p] : NULL;
5506: PetscInt dof, off, d;
5507: const PetscScalar *varr;
5509: PetscSectionGetDof(section, point, &dof);
5510: PetscSectionGetOffset(section, point, &off);
5511: varr = &vArray[off];
5512: if (clperm) {
5513: if (perm) {
5514: for (d = 0; d < dof; d++) array[clperm[offset + perm[d]]] = varr[d];
5515: } else {
5516: for (d = 0; d < dof; d++) array[clperm[offset + d ]] = varr[d];
5517: }
5518: if (flip) {
5519: for (d = 0; d < dof; d++) array[clperm[offset + d ]] *= flip[d];
5520: }
5521: } else {
5522: if (perm) {
5523: for (d = 0; d < dof; d++) array[offset + perm[d]] = varr[d];
5524: } else {
5525: for (d = 0; d < dof; d++) array[offset + d ] = varr[d];
5526: }
5527: if (flip) {
5528: for (d = 0; d < dof; d++) array[offset + d ] *= flip[d];
5529: }
5530: }
5531: offset += dof;
5532: }
5533: PetscSectionRestorePointSyms(section,numPoints,points,&perms,&flips);
5534: *size = offset;
5535: return(0);
5536: }
5538: PETSC_STATIC_INLINE PetscErrorCode DMPlexVecGetClosure_Fields_Static(DM dm, PetscSection section, PetscInt numPoints, const PetscInt points[], PetscInt numFields, const PetscInt clperm[], const PetscScalar vArray[], PetscInt *size, PetscScalar array[])
5539: {
5540: PetscInt offset = 0, f;
5541: PetscErrorCode ierr;
5544: *size = 0;
5545: for (f = 0; f < numFields; ++f) {
5546: PetscInt p;
5547: const PetscInt **perms = NULL;
5548: const PetscScalar **flips = NULL;
5550: PetscSectionGetFieldPointSyms(section,f,numPoints,points,&perms,&flips);
5551: for (p = 0; p < numPoints; p++) {
5552: const PetscInt point = points[2*p];
5553: PetscInt fdof, foff, b;
5554: const PetscScalar *varr;
5555: const PetscInt *perm = perms ? perms[p] : NULL;
5556: const PetscScalar *flip = flips ? flips[p] : NULL;
5558: PetscSectionGetFieldDof(section, point, f, &fdof);
5559: PetscSectionGetFieldOffset(section, point, f, &foff);
5560: varr = &vArray[foff];
5561: if (clperm) {
5562: if (perm) {for (b = 0; b < fdof; b++) {array[clperm[offset + perm[b]]] = varr[b];}}
5563: else {for (b = 0; b < fdof; b++) {array[clperm[offset + b ]] = varr[b];}}
5564: if (flip) {for (b = 0; b < fdof; b++) {array[clperm[offset + b ]] *= flip[b];}}
5565: } else {
5566: if (perm) {for (b = 0; b < fdof; b++) {array[offset + perm[b]] = varr[b];}}
5567: else {for (b = 0; b < fdof; b++) {array[offset + b ] = varr[b];}}
5568: if (flip) {for (b = 0; b < fdof; b++) {array[offset + b ] *= flip[b];}}
5569: }
5570: offset += fdof;
5571: }
5572: PetscSectionRestoreFieldPointSyms(section,f,numPoints,points,&perms,&flips);
5573: }
5574: *size = offset;
5575: return(0);
5576: }
5578: /*@C
5579: DMPlexVecGetClosure - Get an array of the values on the closure of 'point'
5581: Not collective
5583: Input Parameters:
5584: + dm - The DM
5585: . section - The section describing the layout in v, or NULL to use the default section
5586: . v - The local vector
5587: - point - The point in the DM
5589: Input/Output Parameters:
5590: + csize - The size of the input values array, or NULL; on output the number of values in the closure
5591: - values - An array to use for the values, or NULL to have it allocated automatically;
5592: if the user provided NULL, it is a borrowed array and should not be freed
5594: $ Note that DMPlexVecGetClosure/DMPlexVecRestoreClosure only allocates the values array if it set to NULL in the
5595: $ calling function. This is because DMPlexVecGetClosure() is typically called in the inner loop of a Vec or Mat
5596: $ assembly function, and a user may already have allocated storage for this operation.
5597: $
5598: $ A typical use could be
5599: $
5600: $ values = NULL;
5601: $ DMPlexVecGetClosure(dm, NULL, v, p, &clSize, &values);
5602: $ for (cl = 0; cl < clSize; ++cl) {
5603: $ <Compute on closure>
5604: $ }
5605: $ DMPlexVecRestoreClosure(dm, NULL, v, p, &clSize, &values);
5606: $
5607: $ or
5608: $
5609: $ PetscMalloc1(clMaxSize, &values);
5610: $ for (p = pStart; p < pEnd; ++p) {
5611: $ clSize = clMaxSize;
5612: $ DMPlexVecGetClosure(dm, NULL, v, p, &clSize, &values);
5613: $ for (cl = 0; cl < clSize; ++cl) {
5614: $ <Compute on closure>
5615: $ }
5616: $ }
5617: $ PetscFree(values);
5619: Fortran Notes:
5620: Since it returns an array, this routine is only available in Fortran 90, and you must
5621: include petsc.h90 in your code.
5623: The csize argument is not present in the Fortran 90 binding since it is internal to the array.
5625: Level: intermediate
5627: .seealso DMPlexVecRestoreClosure(), DMPlexVecSetClosure(), DMPlexMatSetClosure()
5628: @*/
5629: PetscErrorCode DMPlexVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt *csize, PetscScalar *values[])
5630: {
5631: PetscSection clSection;
5632: IS clPoints;
5633: PetscInt *points = NULL;
5634: const PetscInt *clp, *perm;
5635: PetscInt depth, numFields, numPoints, asize;
5636: PetscErrorCode ierr;
5640: if (!section) {DMGetLocalSection(dm, §ion);}
5643: DMPlexGetDepth(dm, &depth);
5644: PetscSectionGetNumFields(section, &numFields);
5645: if (depth == 1 && numFields < 2) {
5646: DMPlexVecGetClosure_Depth1_Static(dm, section, v, point, csize, values);
5647: return(0);
5648: }
5649: /* Get points */
5650: DMPlexGetCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
5651: /* Get sizes */
5652: asize = 0;
5653: for (PetscInt p = 0; p < numPoints*2; p += 2) {
5654: PetscInt dof;
5655: PetscSectionGetDof(section, points[p], &dof);
5656: asize += dof;
5657: }
5658: if (values) {
5659: const PetscScalar *vArray;
5660: PetscInt size;
5662: if (*values) {
5663: if (PetscUnlikely(*csize < asize)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Provided array size %D not sufficient to hold closure size %D", *csize, asize);
5664: } else {DMGetWorkArray(dm, asize, MPIU_SCALAR, values);}
5665: PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject) dm, depth, asize, &perm);
5666: VecGetArrayRead(v, &vArray);
5667: /* Get values */
5668: if (numFields > 0) {DMPlexVecGetClosure_Fields_Static(dm, section, numPoints, points, numFields, perm, vArray, &size, *values);}
5669: else {DMPlexVecGetClosure_Static(dm, section, numPoints, points, perm, vArray, &size, *values);}
5670: if (PetscUnlikely(asize != size)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Section size %D does not match Vec closure size %D", asize, size);
5671: /* Cleanup array */
5672: VecRestoreArrayRead(v, &vArray);
5673: }
5674: if (csize) *csize = asize;
5675: /* Cleanup points */
5676: DMPlexRestoreCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
5677: return(0);
5678: }
5680: PetscErrorCode DMPlexVecGetClosureAtDepth_Internal(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt depth, PetscInt *csize, PetscScalar *values[])
5681: {
5682: DMLabel depthLabel;
5683: PetscSection clSection;
5684: IS clPoints;
5685: PetscScalar *array;
5686: const PetscScalar *vArray;
5687: PetscInt *points = NULL;
5688: const PetscInt *clp, *perm = NULL;
5689: PetscInt mdepth, numFields, numPoints, Np = 0, p, clsize, size;
5690: PetscErrorCode ierr;
5694: if (!section) {DMGetLocalSection(dm, §ion);}
5697: DMPlexGetDepth(dm, &mdepth);
5698: DMPlexGetDepthLabel(dm, &depthLabel);
5699: PetscSectionGetNumFields(section, &numFields);
5700: if (mdepth == 1 && numFields < 2) {
5701: DMPlexVecGetClosure_Depth1_Static(dm, section, v, point, csize, values);
5702: return(0);
5703: }
5704: /* Get points */
5705: DMPlexGetCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
5706: for (clsize=0,p=0; p<Np; p++) {
5707: PetscInt dof;
5708: PetscSectionGetDof(section, points[2*p], &dof);
5709: clsize += dof;
5710: }
5711: PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject) dm, depth, clsize, &perm);
5712: /* Filter points */
5713: for (p = 0; p < numPoints*2; p += 2) {
5714: PetscInt dep;
5716: DMLabelGetValue(depthLabel, points[p], &dep);
5717: if (dep != depth) continue;
5718: points[Np*2+0] = points[p];
5719: points[Np*2+1] = points[p+1];
5720: ++Np;
5721: }
5722: /* Get array */
5723: if (!values || !*values) {
5724: PetscInt asize = 0, dof;
5726: for (p = 0; p < Np*2; p += 2) {
5727: PetscSectionGetDof(section, points[p], &dof);
5728: asize += dof;
5729: }
5730: if (!values) {
5731: DMPlexRestoreCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
5732: if (csize) *csize = asize;
5733: return(0);
5734: }
5735: DMGetWorkArray(dm, asize, MPIU_SCALAR, &array);
5736: } else {
5737: array = *values;
5738: }
5739: VecGetArrayRead(v, &vArray);
5740: /* Get values */
5741: if (numFields > 0) {DMPlexVecGetClosure_Fields_Static(dm, section, Np, points, numFields, perm, vArray, &size, array);}
5742: else {DMPlexVecGetClosure_Static(dm, section, Np, points, perm, vArray, &size, array);}
5743: /* Cleanup points */
5744: DMPlexRestoreCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
5745: /* Cleanup array */
5746: VecRestoreArrayRead(v, &vArray);
5747: if (!*values) {
5748: if (csize) *csize = size;
5749: *values = array;
5750: } else {
5751: if (size > *csize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Size of input array %D < actual size %D", *csize, size);
5752: *csize = size;
5753: }
5754: return(0);
5755: }
5757: /*@C
5758: DMPlexVecRestoreClosure - Restore the array of the values on the closure of 'point'
5760: Not collective
5762: Input Parameters:
5763: + dm - The DM
5764: . section - The section describing the layout in v, or NULL to use the default section
5765: . v - The local vector
5766: . point - The point in the DM
5767: . csize - The number of values in the closure, or NULL
5768: - values - The array of values, which is a borrowed array and should not be freed
5770: Note that the array values are discarded and not copied back into v. In order to copy values back to v, use DMPlexVecSetClosure()
5772: Fortran Notes:
5773: Since it returns an array, this routine is only available in Fortran 90, and you must
5774: include petsc.h90 in your code.
5776: The csize argument is not present in the Fortran 90 binding since it is internal to the array.
5778: Level: intermediate
5780: .seealso DMPlexVecGetClosure(), DMPlexVecSetClosure(), DMPlexMatSetClosure()
5781: @*/
5782: PetscErrorCode DMPlexVecRestoreClosure(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt *csize, PetscScalar *values[])
5783: {
5784: PetscInt size = 0;
5788: /* Should work without recalculating size */
5789: DMRestoreWorkArray(dm, size, MPIU_SCALAR, (void*) values);
5790: *values = NULL;
5791: return(0);
5792: }
5794: PETSC_STATIC_INLINE void add (PetscScalar *x, PetscScalar y) {*x += y;}
5795: PETSC_STATIC_INLINE void insert(PetscScalar *x, PetscScalar y) {*x = y;}
5797: PETSC_STATIC_INLINE PetscErrorCode updatePoint_private(PetscSection section, PetscInt point, PetscInt dof, void (*fuse)(PetscScalar*, PetscScalar), PetscBool setBC, const PetscInt perm[], const PetscScalar flip[], const PetscInt clperm[], const PetscScalar values[], PetscInt offset, PetscScalar array[])
5798: {
5799: PetscInt cdof; /* The number of constraints on this point */
5800: const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
5801: PetscScalar *a;
5802: PetscInt off, cind = 0, k;
5803: PetscErrorCode ierr;
5806: PetscSectionGetConstraintDof(section, point, &cdof);
5807: PetscSectionGetOffset(section, point, &off);
5808: a = &array[off];
5809: if (!cdof || setBC) {
5810: if (clperm) {
5811: if (perm) {for (k = 0; k < dof; ++k) {fuse(&a[k], values[clperm[offset+perm[k]]] * (flip ? flip[perm[k]] : 1.));}}
5812: else {for (k = 0; k < dof; ++k) {fuse(&a[k], values[clperm[offset+ k ]] * (flip ? flip[ k ] : 1.));}}
5813: } else {
5814: if (perm) {for (k = 0; k < dof; ++k) {fuse(&a[k], values[offset+perm[k]] * (flip ? flip[perm[k]] : 1.));}}
5815: else {for (k = 0; k < dof; ++k) {fuse(&a[k], values[offset+ k ] * (flip ? flip[ k ] : 1.));}}
5816: }
5817: } else {
5818: PetscSectionGetConstraintIndices(section, point, &cdofs);
5819: if (clperm) {
5820: if (perm) {for (k = 0; k < dof; ++k) {
5821: if ((cind < cdof) && (k == cdofs[cind])) {++cind; continue;}
5822: fuse(&a[k], values[clperm[offset+perm[k]]] * (flip ? flip[perm[k]] : 1.));
5823: }
5824: } else {
5825: for (k = 0; k < dof; ++k) {
5826: if ((cind < cdof) && (k == cdofs[cind])) {++cind; continue;}
5827: fuse(&a[k], values[clperm[offset+ k ]] * (flip ? flip[ k ] : 1.));
5828: }
5829: }
5830: } else {
5831: if (perm) {
5832: for (k = 0; k < dof; ++k) {
5833: if ((cind < cdof) && (k == cdofs[cind])) {++cind; continue;}
5834: fuse(&a[k], values[offset+perm[k]] * (flip ? flip[perm[k]] : 1.));
5835: }
5836: } else {
5837: for (k = 0; k < dof; ++k) {
5838: if ((cind < cdof) && (k == cdofs[cind])) {++cind; continue;}
5839: fuse(&a[k], values[offset+ k ] * (flip ? flip[ k ] : 1.));
5840: }
5841: }
5842: }
5843: }
5844: return(0);
5845: }
5847: PETSC_STATIC_INLINE PetscErrorCode updatePointBC_private(PetscSection section, PetscInt point, PetscInt dof, void (*fuse)(PetscScalar*, PetscScalar), const PetscInt perm[], const PetscScalar flip[], const PetscInt clperm[], const PetscScalar values[], PetscInt offset, PetscScalar array[])
5848: {
5849: PetscInt cdof; /* The number of constraints on this point */
5850: const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
5851: PetscScalar *a;
5852: PetscInt off, cind = 0, k;
5853: PetscErrorCode ierr;
5856: PetscSectionGetConstraintDof(section, point, &cdof);
5857: PetscSectionGetOffset(section, point, &off);
5858: a = &array[off];
5859: if (cdof) {
5860: PetscSectionGetConstraintIndices(section, point, &cdofs);
5861: if (clperm) {
5862: if (perm) {
5863: for (k = 0; k < dof; ++k) {
5864: if ((cind < cdof) && (k == cdofs[cind])) {
5865: fuse(&a[k], values[clperm[offset+perm[k]]] * (flip ? flip[perm[k]] : 1.));
5866: cind++;
5867: }
5868: }
5869: } else {
5870: for (k = 0; k < dof; ++k) {
5871: if ((cind < cdof) && (k == cdofs[cind])) {
5872: fuse(&a[k], values[clperm[offset+ k ]] * (flip ? flip[ k ] : 1.));
5873: cind++;
5874: }
5875: }
5876: }
5877: } else {
5878: if (perm) {
5879: for (k = 0; k < dof; ++k) {
5880: if ((cind < cdof) && (k == cdofs[cind])) {
5881: fuse(&a[k], values[offset+perm[k]] * (flip ? flip[perm[k]] : 1.));
5882: cind++;
5883: }
5884: }
5885: } else {
5886: for (k = 0; k < dof; ++k) {
5887: if ((cind < cdof) && (k == cdofs[cind])) {
5888: fuse(&a[k], values[offset+ k ] * (flip ? flip[ k ] : 1.));
5889: cind++;
5890: }
5891: }
5892: }
5893: }
5894: }
5895: return(0);
5896: }
5898: PETSC_STATIC_INLINE PetscErrorCode updatePointFields_private(PetscSection section, PetscInt point, const PetscInt *perm, const PetscScalar *flip, PetscInt f, void (*fuse)(PetscScalar*, PetscScalar), PetscBool setBC, const PetscInt clperm[], const PetscScalar values[], PetscInt *offset, PetscScalar array[])
5899: {
5900: PetscScalar *a;
5901: PetscInt fdof, foff, fcdof, foffset = *offset;
5902: const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
5903: PetscInt cind = 0, b;
5904: PetscErrorCode ierr;
5907: PetscSectionGetFieldDof(section, point, f, &fdof);
5908: PetscSectionGetFieldConstraintDof(section, point, f, &fcdof);
5909: PetscSectionGetFieldOffset(section, point, f, &foff);
5910: a = &array[foff];
5911: if (!fcdof || setBC) {
5912: if (clperm) {
5913: if (perm) {for (b = 0; b < fdof; b++) {fuse(&a[b], values[clperm[foffset+perm[b]]] * (flip ? flip[perm[b]] : 1.));}}
5914: else {for (b = 0; b < fdof; b++) {fuse(&a[b], values[clperm[foffset+ b ]] * (flip ? flip[ b ] : 1.));}}
5915: } else {
5916: if (perm) {for (b = 0; b < fdof; b++) {fuse(&a[b], values[foffset+perm[b]] * (flip ? flip[perm[b]] : 1.));}}
5917: else {for (b = 0; b < fdof; b++) {fuse(&a[b], values[foffset+ b ] * (flip ? flip[ b ] : 1.));}}
5918: }
5919: } else {
5920: PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs);
5921: if (clperm) {
5922: if (perm) {
5923: for (b = 0; b < fdof; b++) {
5924: if ((cind < fcdof) && (b == fcdofs[cind])) {++cind; continue;}
5925: fuse(&a[b], values[clperm[foffset+perm[b]]] * (flip ? flip[perm[b]] : 1.));
5926: }
5927: } else {
5928: for (b = 0; b < fdof; b++) {
5929: if ((cind < fcdof) && (b == fcdofs[cind])) {++cind; continue;}
5930: fuse(&a[b], values[clperm[foffset+ b ]] * (flip ? flip[ b ] : 1.));
5931: }
5932: }
5933: } else {
5934: if (perm) {
5935: for (b = 0; b < fdof; b++) {
5936: if ((cind < fcdof) && (b == fcdofs[cind])) {++cind; continue;}
5937: fuse(&a[b], values[foffset+perm[b]] * (flip ? flip[perm[b]] : 1.));
5938: }
5939: } else {
5940: for (b = 0; b < fdof; b++) {
5941: if ((cind < fcdof) && (b == fcdofs[cind])) {++cind; continue;}
5942: fuse(&a[b], values[foffset+ b ] * (flip ? flip[ b ] : 1.));
5943: }
5944: }
5945: }
5946: }
5947: *offset += fdof;
5948: return(0);
5949: }
5951: PETSC_STATIC_INLINE PetscErrorCode updatePointFieldsBC_private(PetscSection section, PetscInt point, const PetscInt perm[], const PetscScalar flip[], PetscInt f, PetscInt Ncc, const PetscInt comps[], void (*fuse)(PetscScalar*, PetscScalar), const PetscInt clperm[], const PetscScalar values[], PetscInt *offset, PetscScalar array[])
5952: {
5953: PetscScalar *a;
5954: PetscInt fdof, foff, fcdof, foffset = *offset;
5955: const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
5956: PetscInt Nc, cind = 0, ncind = 0, b;
5957: PetscBool ncSet, fcSet;
5958: PetscErrorCode ierr;
5961: PetscSectionGetFieldComponents(section, f, &Nc);
5962: PetscSectionGetFieldDof(section, point, f, &fdof);
5963: PetscSectionGetFieldConstraintDof(section, point, f, &fcdof);
5964: PetscSectionGetFieldOffset(section, point, f, &foff);
5965: a = &array[foff];
5966: if (fcdof) {
5967: /* We just override fcdof and fcdofs with Ncc and comps */
5968: PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs);
5969: if (clperm) {
5970: if (perm) {
5971: if (comps) {
5972: for (b = 0; b < fdof; b++) {
5973: ncSet = fcSet = PETSC_FALSE;
5974: if (b%Nc == comps[ncind]) {ncind = (ncind+1)%Ncc; ncSet = PETSC_TRUE;}
5975: if ((cind < fcdof) && (b == fcdofs[cind])) {++cind; fcSet = PETSC_TRUE;}
5976: if (ncSet && fcSet) {fuse(&a[b], values[clperm[foffset+perm[b]]] * (flip ? flip[perm[b]] : 1.));}
5977: }
5978: } else {
5979: for (b = 0; b < fdof; b++) {
5980: if ((cind < fcdof) && (b == fcdofs[cind])) {
5981: fuse(&a[b], values[clperm[foffset+perm[b]]] * (flip ? flip[perm[b]] : 1.));
5982: ++cind;
5983: }
5984: }
5985: }
5986: } else {
5987: if (comps) {
5988: for (b = 0; b < fdof; b++) {
5989: ncSet = fcSet = PETSC_FALSE;
5990: if (b%Nc == comps[ncind]) {ncind = (ncind+1)%Ncc; ncSet = PETSC_TRUE;}
5991: if ((cind < fcdof) && (b == fcdofs[cind])) {++cind; fcSet = PETSC_TRUE;}
5992: if (ncSet && fcSet) {fuse(&a[b], values[clperm[foffset+ b ]] * (flip ? flip[ b ] : 1.));}
5993: }
5994: } else {
5995: for (b = 0; b < fdof; b++) {
5996: if ((cind < fcdof) && (b == fcdofs[cind])) {
5997: fuse(&a[b], values[clperm[foffset+ b ]] * (flip ? flip[ b ] : 1.));
5998: ++cind;
5999: }
6000: }
6001: }
6002: }
6003: } else {
6004: if (perm) {
6005: if (comps) {
6006: for (b = 0; b < fdof; b++) {
6007: ncSet = fcSet = PETSC_FALSE;
6008: if (b%Nc == comps[ncind]) {ncind = (ncind+1)%Ncc; ncSet = PETSC_TRUE;}
6009: if ((cind < fcdof) && (b == fcdofs[cind])) {++cind; fcSet = PETSC_TRUE;}
6010: if (ncSet && fcSet) {fuse(&a[b], values[foffset+perm[b]] * (flip ? flip[perm[b]] : 1.));}
6011: }
6012: } else {
6013: for (b = 0; b < fdof; b++) {
6014: if ((cind < fcdof) && (b == fcdofs[cind])) {
6015: fuse(&a[b], values[foffset+perm[b]] * (flip ? flip[perm[b]] : 1.));
6016: ++cind;
6017: }
6018: }
6019: }
6020: } else {
6021: if (comps) {
6022: for (b = 0; b < fdof; b++) {
6023: ncSet = fcSet = PETSC_FALSE;
6024: if (b%Nc == comps[ncind]) {ncind = (ncind+1)%Ncc; ncSet = PETSC_TRUE;}
6025: if ((cind < fcdof) && (b == fcdofs[cind])) {++cind; fcSet = PETSC_TRUE;}
6026: if (ncSet && fcSet) {fuse(&a[b], values[foffset+ b ] * (flip ? flip[ b ] : 1.));}
6027: }
6028: } else {
6029: for (b = 0; b < fdof; b++) {
6030: if ((cind < fcdof) && (b == fcdofs[cind])) {
6031: fuse(&a[b], values[foffset+ b ] * (flip ? flip[ b ] : 1.));
6032: ++cind;
6033: }
6034: }
6035: }
6036: }
6037: }
6038: }
6039: *offset += fdof;
6040: return(0);
6041: }
6043: PETSC_STATIC_INLINE PetscErrorCode DMPlexVecSetClosure_Depth1_Static(DM dm, PetscSection section, Vec v, PetscInt point, const PetscScalar values[], InsertMode mode)
6044: {
6045: PetscScalar *array;
6046: const PetscInt *cone, *coneO;
6047: PetscInt pStart, pEnd, p, numPoints, off, dof;
6048: PetscErrorCode ierr;
6051: PetscSectionGetChart(section, &pStart, &pEnd);
6052: DMPlexGetConeSize(dm, point, &numPoints);
6053: DMPlexGetCone(dm, point, &cone);
6054: DMPlexGetConeOrientation(dm, point, &coneO);
6055: VecGetArray(v, &array);
6056: for (p = 0, off = 0; p <= numPoints; ++p, off += dof) {
6057: const PetscInt cp = !p ? point : cone[p-1];
6058: const PetscInt o = !p ? 0 : coneO[p-1];
6060: if ((cp < pStart) || (cp >= pEnd)) {dof = 0; continue;}
6061: PetscSectionGetDof(section, cp, &dof);
6062: /* ADD_VALUES */
6063: {
6064: const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
6065: PetscScalar *a;
6066: PetscInt cdof, coff, cind = 0, k;
6068: PetscSectionGetConstraintDof(section, cp, &cdof);
6069: PetscSectionGetOffset(section, cp, &coff);
6070: a = &array[coff];
6071: if (!cdof) {
6072: if (o >= 0) {
6073: for (k = 0; k < dof; ++k) {
6074: a[k] += values[off+k];
6075: }
6076: } else {
6077: for (k = 0; k < dof; ++k) {
6078: a[k] += values[off+dof-k-1];
6079: }
6080: }
6081: } else {
6082: PetscSectionGetConstraintIndices(section, cp, &cdofs);
6083: if (o >= 0) {
6084: for (k = 0; k < dof; ++k) {
6085: if ((cind < cdof) && (k == cdofs[cind])) {++cind; continue;}
6086: a[k] += values[off+k];
6087: }
6088: } else {
6089: for (k = 0; k < dof; ++k) {
6090: if ((cind < cdof) && (k == cdofs[cind])) {++cind; continue;}
6091: a[k] += values[off+dof-k-1];
6092: }
6093: }
6094: }
6095: }
6096: }
6097: VecRestoreArray(v, &array);
6098: return(0);
6099: }
6101: /*@C
6102: DMPlexVecSetClosure - Set an array of the values on the closure of 'point'
6104: Not collective
6106: Input Parameters:
6107: + dm - The DM
6108: . section - The section describing the layout in v, or NULL to use the default section
6109: . v - The local vector
6110: . point - The point in the DM
6111: . values - The array of values
6112: - mode - The insert mode. One of INSERT_ALL_VALUES, ADD_ALL_VALUES, INSERT_VALUES, ADD_VALUES, INSERT_BC_VALUES, and ADD_BC_VALUES,
6113: where INSERT_ALL_VALUES and ADD_ALL_VALUES also overwrite boundary conditions.
6115: Fortran Notes:
6116: This routine is only available in Fortran 90, and you must include petsc.h90 in your code.
6118: Level: intermediate
6120: .seealso DMPlexVecGetClosure(), DMPlexMatSetClosure()
6121: @*/
6122: PetscErrorCode DMPlexVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt point, const PetscScalar values[], InsertMode mode)
6123: {
6124: PetscSection clSection;
6125: IS clPoints;
6126: PetscScalar *array;
6127: PetscInt *points = NULL;
6128: const PetscInt *clp, *clperm = NULL;
6129: PetscInt depth, numFields, numPoints, p, clsize;
6130: PetscErrorCode ierr;
6134: if (!section) {DMGetLocalSection(dm, §ion);}
6137: DMPlexGetDepth(dm, &depth);
6138: PetscSectionGetNumFields(section, &numFields);
6139: if (depth == 1 && numFields < 2 && mode == ADD_VALUES) {
6140: DMPlexVecSetClosure_Depth1_Static(dm, section, v, point, values, mode);
6141: return(0);
6142: }
6143: /* Get points */
6144: DMPlexGetCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
6145: for (clsize=0,p=0; p<numPoints; p++) {
6146: PetscInt dof;
6147: PetscSectionGetDof(section, points[2*p], &dof);
6148: clsize += dof;
6149: }
6150: PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject) dm, depth, clsize, &clperm);
6151: /* Get array */
6152: VecGetArray(v, &array);
6153: /* Get values */
6154: if (numFields > 0) {
6155: PetscInt offset = 0, f;
6156: for (f = 0; f < numFields; ++f) {
6157: const PetscInt **perms = NULL;
6158: const PetscScalar **flips = NULL;
6160: PetscSectionGetFieldPointSyms(section,f,numPoints,points,&perms,&flips);
6161: switch (mode) {
6162: case INSERT_VALUES:
6163: for (p = 0; p < numPoints; p++) {
6164: const PetscInt point = points[2*p];
6165: const PetscInt *perm = perms ? perms[p] : NULL;
6166: const PetscScalar *flip = flips ? flips[p] : NULL;
6167: updatePointFields_private(section, point, perm, flip, f, insert, PETSC_FALSE, clperm, values, &offset, array);
6168: } break;
6169: case INSERT_ALL_VALUES:
6170: for (p = 0; p < numPoints; p++) {
6171: const PetscInt point = points[2*p];
6172: const PetscInt *perm = perms ? perms[p] : NULL;
6173: const PetscScalar *flip = flips ? flips[p] : NULL;
6174: updatePointFields_private(section, point, perm, flip, f, insert, PETSC_TRUE, clperm, values, &offset, array);
6175: } break;
6176: case INSERT_BC_VALUES:
6177: for (p = 0; p < numPoints; p++) {
6178: const PetscInt point = points[2*p];
6179: const PetscInt *perm = perms ? perms[p] : NULL;
6180: const PetscScalar *flip = flips ? flips[p] : NULL;
6181: updatePointFieldsBC_private(section, point, perm, flip, f, -1, NULL, insert, clperm, values, &offset, array);
6182: } break;
6183: case ADD_VALUES:
6184: for (p = 0; p < numPoints; p++) {
6185: const PetscInt point = points[2*p];
6186: const PetscInt *perm = perms ? perms[p] : NULL;
6187: const PetscScalar *flip = flips ? flips[p] : NULL;
6188: updatePointFields_private(section, point, perm, flip, f, add, PETSC_FALSE, clperm, values, &offset, array);
6189: } break;
6190: case ADD_ALL_VALUES:
6191: for (p = 0; p < numPoints; p++) {
6192: const PetscInt point = points[2*p];
6193: const PetscInt *perm = perms ? perms[p] : NULL;
6194: const PetscScalar *flip = flips ? flips[p] : NULL;
6195: updatePointFields_private(section, point, perm, flip, f, add, PETSC_TRUE, clperm, values, &offset, array);
6196: } break;
6197: case ADD_BC_VALUES:
6198: for (p = 0; p < numPoints; p++) {
6199: const PetscInt point = points[2*p];
6200: const PetscInt *perm = perms ? perms[p] : NULL;
6201: const PetscScalar *flip = flips ? flips[p] : NULL;
6202: updatePointFieldsBC_private(section, point, perm, flip, f, -1, NULL, add, clperm, values, &offset, array);
6203: } break;
6204: default:
6205: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insert mode %d", mode);
6206: }
6207: PetscSectionRestoreFieldPointSyms(section,f,numPoints,points,&perms,&flips);
6208: }
6209: } else {
6210: PetscInt dof, off;
6211: const PetscInt **perms = NULL;
6212: const PetscScalar **flips = NULL;
6214: PetscSectionGetPointSyms(section,numPoints,points,&perms,&flips);
6215: switch (mode) {
6216: case INSERT_VALUES:
6217: for (p = 0, off = 0; p < numPoints; p++, off += dof) {
6218: const PetscInt point = points[2*p];
6219: const PetscInt *perm = perms ? perms[p] : NULL;
6220: const PetscScalar *flip = flips ? flips[p] : NULL;
6221: PetscSectionGetDof(section, point, &dof);
6222: updatePoint_private(section, point, dof, insert, PETSC_FALSE, perm, flip, clperm, values, off, array);
6223: } break;
6224: case INSERT_ALL_VALUES:
6225: for (p = 0, off = 0; p < numPoints; p++, off += dof) {
6226: const PetscInt point = points[2*p];
6227: const PetscInt *perm = perms ? perms[p] : NULL;
6228: const PetscScalar *flip = flips ? flips[p] : NULL;
6229: PetscSectionGetDof(section, point, &dof);
6230: updatePoint_private(section, point, dof, insert, PETSC_TRUE, perm, flip, clperm, values, off, array);
6231: } break;
6232: case INSERT_BC_VALUES:
6233: for (p = 0, off = 0; p < numPoints; p++, off += dof) {
6234: const PetscInt point = points[2*p];
6235: const PetscInt *perm = perms ? perms[p] : NULL;
6236: const PetscScalar *flip = flips ? flips[p] : NULL;
6237: PetscSectionGetDof(section, point, &dof);
6238: updatePointBC_private(section, point, dof, insert, perm, flip, clperm, values, off, array);
6239: } break;
6240: case ADD_VALUES:
6241: for (p = 0, off = 0; p < numPoints; p++, off += dof) {
6242: const PetscInt point = points[2*p];
6243: const PetscInt *perm = perms ? perms[p] : NULL;
6244: const PetscScalar *flip = flips ? flips[p] : NULL;
6245: PetscSectionGetDof(section, point, &dof);
6246: updatePoint_private(section, point, dof, add, PETSC_FALSE, perm, flip, clperm, values, off, array);
6247: } break;
6248: case ADD_ALL_VALUES:
6249: for (p = 0, off = 0; p < numPoints; p++, off += dof) {
6250: const PetscInt point = points[2*p];
6251: const PetscInt *perm = perms ? perms[p] : NULL;
6252: const PetscScalar *flip = flips ? flips[p] : NULL;
6253: PetscSectionGetDof(section, point, &dof);
6254: updatePoint_private(section, point, dof, add, PETSC_TRUE, perm, flip, clperm, values, off, array);
6255: } break;
6256: case ADD_BC_VALUES:
6257: for (p = 0, off = 0; p < numPoints; p++, off += dof) {
6258: const PetscInt point = points[2*p];
6259: const PetscInt *perm = perms ? perms[p] : NULL;
6260: const PetscScalar *flip = flips ? flips[p] : NULL;
6261: PetscSectionGetDof(section, point, &dof);
6262: updatePointBC_private(section, point, dof, add, perm, flip, clperm, values, off, array);
6263: } break;
6264: default:
6265: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insert mode %d", mode);
6266: }
6267: PetscSectionRestorePointSyms(section,numPoints,points,&perms,&flips);
6268: }
6269: /* Cleanup points */
6270: DMPlexRestoreCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
6271: /* Cleanup array */
6272: VecRestoreArray(v, &array);
6273: return(0);
6274: }
6276: /* Check whether the given point is in the label. If not, update the offset to skip this point */
6277: PETSC_STATIC_INLINE PetscErrorCode CheckPoint_Private(DMLabel label, PetscInt labelId, PetscSection section, PetscInt point, PetscInt f, PetscInt *offset)
6278: {
6280: if (label) {
6281: PetscInt val, fdof;
6284: /* There is a problem with this:
6285: Suppose we have two label values, defining surfaces, interecting along a line in 3D. When we add cells to the label, the cells that
6286: touch both surfaces must pick a label value. Thus we miss setting values for the surface with that other value intersecting that cell.
6287: Thus I am only going to check val != -1, not val != labelId
6288: */
6289: DMLabelGetValue(label, point, &val);
6290: if (val < 0) {
6291: PetscSectionGetFieldDof(section, point, f, &fdof);
6292: *offset += fdof;
6293: PetscFunctionReturn(1);
6294: }
6295: }
6296: return(0);
6297: }
6299: /* Unlike DMPlexVecSetClosure(), this uses plex-native closure permutation, not a user-specified permutation such as DMPlexSetClosurePermutationTensor(). */
6300: PetscErrorCode DMPlexVecSetFieldClosure_Internal(DM dm, PetscSection section, Vec v, PetscBool fieldActive[], PetscInt point, PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt labelId, const PetscScalar values[], InsertMode mode)
6301: {
6302: PetscSection clSection;
6303: IS clPoints;
6304: PetscScalar *array;
6305: PetscInt *points = NULL;
6306: const PetscInt *clp;
6307: PetscInt numFields, numPoints, p;
6308: PetscInt offset = 0, f;
6309: PetscErrorCode ierr;
6313: if (!section) {DMGetLocalSection(dm, §ion);}
6316: PetscSectionGetNumFields(section, &numFields);
6317: /* Get points */
6318: DMPlexGetCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
6319: /* Get array */
6320: VecGetArray(v, &array);
6321: /* Get values */
6322: for (f = 0; f < numFields; ++f) {
6323: const PetscInt **perms = NULL;
6324: const PetscScalar **flips = NULL;
6326: if (!fieldActive[f]) {
6327: for (p = 0; p < numPoints*2; p += 2) {
6328: PetscInt fdof;
6329: PetscSectionGetFieldDof(section, points[p], f, &fdof);
6330: offset += fdof;
6331: }
6332: continue;
6333: }
6334: PetscSectionGetFieldPointSyms(section,f,numPoints,points,&perms,&flips);
6335: switch (mode) {
6336: case INSERT_VALUES:
6337: for (p = 0; p < numPoints; p++) {
6338: const PetscInt point = points[2*p];
6339: const PetscInt *perm = perms ? perms[p] : NULL;
6340: const PetscScalar *flip = flips ? flips[p] : NULL;
6341: CheckPoint_Private(label, labelId, section, point, f, &offset); if (ierr) continue;
6342: updatePointFields_private(section, point, perm, flip, f, insert, PETSC_FALSE, NULL, values, &offset, array);
6343: } break;
6344: case INSERT_ALL_VALUES:
6345: for (p = 0; p < numPoints; p++) {
6346: const PetscInt point = points[2*p];
6347: const PetscInt *perm = perms ? perms[p] : NULL;
6348: const PetscScalar *flip = flips ? flips[p] : NULL;
6349: CheckPoint_Private(label, labelId, section, point, f, &offset); if (ierr) continue;
6350: updatePointFields_private(section, point, perm, flip, f, insert, PETSC_TRUE, NULL, values, &offset, array);
6351: } break;
6352: case INSERT_BC_VALUES:
6353: for (p = 0; p < numPoints; p++) {
6354: const PetscInt point = points[2*p];
6355: const PetscInt *perm = perms ? perms[p] : NULL;
6356: const PetscScalar *flip = flips ? flips[p] : NULL;
6357: CheckPoint_Private(label, labelId, section, point, f, &offset); if (ierr) continue;
6358: updatePointFieldsBC_private(section, point, perm, flip, f, Ncc, comps, insert, NULL, values, &offset, array);
6359: } break;
6360: case ADD_VALUES:
6361: for (p = 0; p < numPoints; p++) {
6362: const PetscInt point = points[2*p];
6363: const PetscInt *perm = perms ? perms[p] : NULL;
6364: const PetscScalar *flip = flips ? flips[p] : NULL;
6365: CheckPoint_Private(label, labelId, section, point, f, &offset); if (ierr) continue;
6366: updatePointFields_private(section, point, perm, flip, f, add, PETSC_FALSE, NULL, values, &offset, array);
6367: } break;
6368: case ADD_ALL_VALUES:
6369: for (p = 0; p < numPoints; p++) {
6370: const PetscInt point = points[2*p];
6371: const PetscInt *perm = perms ? perms[p] : NULL;
6372: const PetscScalar *flip = flips ? flips[p] : NULL;
6373: CheckPoint_Private(label, labelId, section, point, f, &offset); if (ierr) continue;
6374: updatePointFields_private(section, point, perm, flip, f, add, PETSC_TRUE, NULL, values, &offset, array);
6375: } break;
6376: default:
6377: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insert mode %d", mode);
6378: }
6379: PetscSectionRestoreFieldPointSyms(section,f,numPoints,points,&perms,&flips);
6380: }
6381: /* Cleanup points */
6382: DMPlexRestoreCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
6383: /* Cleanup array */
6384: VecRestoreArray(v, &array);
6385: return(0);
6386: }
6388: static PetscErrorCode DMPlexPrintMatSetValues(PetscViewer viewer, Mat A, PetscInt point, PetscInt numRIndices, const PetscInt rindices[], PetscInt numCIndices, const PetscInt cindices[], const PetscScalar values[])
6389: {
6390: PetscMPIInt rank;
6391: PetscInt i, j;
6395: MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank);
6396: PetscViewerASCIIPrintf(viewer, "[%d]mat for point %D\n", rank, point);
6397: for (i = 0; i < numRIndices; i++) {PetscViewerASCIIPrintf(viewer, "[%d]mat row indices[%D] = %D\n", rank, i, rindices[i]);}
6398: for (i = 0; i < numCIndices; i++) {PetscViewerASCIIPrintf(viewer, "[%d]mat col indices[%D] = %D\n", rank, i, cindices[i]);}
6399: numCIndices = numCIndices ? numCIndices : numRIndices;
6400: if (!values) return(0);
6401: for (i = 0; i < numRIndices; i++) {
6402: PetscViewerASCIIPrintf(viewer, "[%d]", rank);
6403: for (j = 0; j < numCIndices; j++) {
6404: #if defined(PETSC_USE_COMPLEX)
6405: PetscViewerASCIIPrintf(viewer, " (%g,%g)", (double)PetscRealPart(values[i*numCIndices+j]), (double)PetscImaginaryPart(values[i*numCIndices+j]));
6406: #else
6407: PetscViewerASCIIPrintf(viewer, " %g", (double)values[i*numCIndices+j]);
6408: #endif
6409: }
6410: PetscViewerASCIIPrintf(viewer, "\n");
6411: }
6412: return(0);
6413: }
6415: /*
6416: DMPlexGetIndicesPoint_Internal - Add the indices for dofs on a point to an index array
6418: Input Parameters:
6419: + section - The section for this data layout
6420: . islocal - Is the section (and thus indices being requested) local or global?
6421: . point - The point contributing dofs with these indices
6422: . off - The global offset of this point
6423: . loff - The local offset of each field
6424: . setBC - The flag determining whether to include indices of boundary values
6425: . perm - A permutation of the dofs on this point, or NULL
6426: - indperm - A permutation of the entire indices array, or NULL
6428: Output Parameter:
6429: . indices - Indices for dofs on this point
6431: Level: developer
6433: Note: The indices could be local or global, depending on the value of 'off'.
6434: */
6435: PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection section, PetscBool islocal,PetscInt point, PetscInt off, PetscInt *loff, PetscBool setBC, const PetscInt perm[], const PetscInt indperm[], PetscInt indices[])
6436: {
6437: PetscInt dof; /* The number of unknowns on this point */
6438: PetscInt cdof; /* The number of constraints on this point */
6439: const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
6440: PetscInt cind = 0, k;
6441: PetscErrorCode ierr;
6444: if (!islocal && setBC) SETERRQ(PetscObjectComm((PetscObject)section),PETSC_ERR_ARG_INCOMP,"setBC incompatible with global indices; use a local section or disable setBC");
6445: PetscSectionGetDof(section, point, &dof);
6446: PetscSectionGetConstraintDof(section, point, &cdof);
6447: if (!cdof || setBC) {
6448: for (k = 0; k < dof; ++k) {
6449: const PetscInt preind = perm ? *loff+perm[k] : *loff+k;
6450: const PetscInt ind = indperm ? indperm[preind] : preind;
6452: indices[ind] = off + k;
6453: }
6454: } else {
6455: PetscSectionGetConstraintIndices(section, point, &cdofs);
6456: for (k = 0; k < dof; ++k) {
6457: const PetscInt preind = perm ? *loff+perm[k] : *loff+k;
6458: const PetscInt ind = indperm ? indperm[preind] : preind;
6460: if ((cind < cdof) && (k == cdofs[cind])) {
6461: /* Insert check for returning constrained indices */
6462: indices[ind] = -(off+k+1);
6463: ++cind;
6464: } else {
6465: indices[ind] = off + k - (islocal ? 0 : cind);
6466: }
6467: }
6468: }
6469: *loff += dof;
6470: return(0);
6471: }
6473: /*
6474: DMPlexGetIndicesPointFields_Internal - gets section indices for a point in its canonical ordering.
6476: Input Parameters:
6477: + section - a section (global or local)
6478: - islocal - PETSC_TRUE if requesting local indices (i.e., section is local); PETSC_FALSE for global
6479: . point - point within section
6480: . off - The offset of this point in the (local or global) indexed space - should match islocal and (usually) the section
6481: . foffs - array of length numFields containing the offset in canonical point ordering (the location in indices) of each field
6482: . setBC - identify constrained (boundary condition) points via involution.
6483: . perms - perms[f][permsoff][:] is a permutation of dofs within each field
6484: . permsoff - offset
6485: - indperm - index permutation
6487: Output Parameter:
6488: . foffs - each entry is incremented by the number of (unconstrained if setBC=FALSE) dofs in that field
6489: . indices - array to hold indices (as defined by section) of each dof associated with point
6491: Notes:
6492: If section is local and setBC=true, there is no distinction between constrained and unconstrained dofs.
6493: If section is local and setBC=false, the indices for constrained points are the involution -(i+1) of their position
6494: in the local vector.
6496: If section is global and setBC=false, the indices for constrained points are negative (and their value is not
6497: significant). It is invalid to call with a global section and setBC=true.
6499: Developer Note:
6500: The section is only used for field layout, so islocal is technically a statement about the offset (off). At some point
6501: in the future, global sections may have fields set, in which case we could pass the global section and obtain the
6502: offset could be obtained from the section instead of passing it explicitly as we do now.
6504: Example:
6505: Suppose a point contains one field with three components, and for which the unconstrained indices are {10, 11, 12}.
6506: When the middle component is constrained, we get the array {10, -12, 12} for (islocal=TRUE, setBC=FALSE).
6507: Note that -12 is the involution of 11, so the user can involute negative indices to recover local indices.
6508: The global vector does not store constrained dofs, so when this function returns global indices, say {110, -112, 111}, the value of -112 is an arbitrary flag that should not be interpreted beyond its sign.
6510: Level: developer
6511: */
6512: PetscErrorCode DMPlexGetIndicesPointFields_Internal(PetscSection section, PetscBool islocal, PetscInt point, PetscInt off, PetscInt foffs[], PetscBool setBC, const PetscInt ***perms, PetscInt permsoff, const PetscInt indperm[], PetscInt indices[])
6513: {
6514: PetscInt numFields, foff, f;
6518: if (!islocal && setBC) SETERRQ(PetscObjectComm((PetscObject)section),PETSC_ERR_ARG_INCOMP,"setBC incompatible with global indices; use a local section or disable setBC");
6519: PetscSectionGetNumFields(section, &numFields);
6520: for (f = 0, foff = 0; f < numFields; ++f) {
6521: PetscInt fdof, cfdof;
6522: const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
6523: PetscInt cind = 0, b;
6524: const PetscInt *perm = (perms && perms[f]) ? perms[f][permsoff] : NULL;
6526: PetscSectionGetFieldDof(section, point, f, &fdof);
6527: PetscSectionGetFieldConstraintDof(section, point, f, &cfdof);
6528: if (!cfdof || setBC) {
6529: for (b = 0; b < fdof; ++b) {
6530: const PetscInt preind = perm ? foffs[f]+perm[b] : foffs[f]+b;
6531: const PetscInt ind = indperm ? indperm[preind] : preind;
6533: indices[ind] = off+foff+b;
6534: }
6535: } else {
6536: PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs);
6537: for (b = 0; b < fdof; ++b) {
6538: const PetscInt preind = perm ? foffs[f]+perm[b] : foffs[f]+b;
6539: const PetscInt ind = indperm ? indperm[preind] : preind;
6541: if ((cind < cfdof) && (b == fcdofs[cind])) {
6542: indices[ind] = -(off+foff+b+1);
6543: ++cind;
6544: } else {
6545: indices[ind] = off + foff + b - (islocal ? 0 : cind);
6546: }
6547: }
6548: }
6549: foff += (setBC || islocal ? fdof : (fdof - cfdof));
6550: foffs[f] += fdof;
6551: }
6552: return(0);
6553: }
6555: /*
6556: This version believes the globalSection offsets for each field, rather than just the point offset
6558: . foffs - The offset into 'indices' for each field, since it is segregated by field
6560: Notes:
6561: The semantics of this function relate to that of setBC=FALSE in DMPlexGetIndicesPointFields_Internal.
6562: Since this function uses global indices, setBC=TRUE would be invalid, so no such argument exists.
6563: */
6564: static PetscErrorCode DMPlexGetIndicesPointFieldsSplit_Internal(PetscSection section, PetscSection globalSection, PetscInt point, PetscInt foffs[], const PetscInt ***perms, PetscInt permsoff, const PetscInt indperm[], PetscInt indices[])
6565: {
6566: PetscInt numFields, foff, f;
6570: PetscSectionGetNumFields(section, &numFields);
6571: for (f = 0; f < numFields; ++f) {
6572: PetscInt fdof, cfdof;
6573: const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
6574: PetscInt cind = 0, b;
6575: const PetscInt *perm = (perms && perms[f]) ? perms[f][permsoff] : NULL;
6577: PetscSectionGetFieldDof(section, point, f, &fdof);
6578: PetscSectionGetFieldConstraintDof(section, point, f, &cfdof);
6579: PetscSectionGetFieldOffset(globalSection, point, f, &foff);
6580: if (!cfdof) {
6581: for (b = 0; b < fdof; ++b) {
6582: const PetscInt preind = perm ? foffs[f]+perm[b] : foffs[f]+b;
6583: const PetscInt ind = indperm ? indperm[preind] : preind;
6585: indices[ind] = foff+b;
6586: }
6587: } else {
6588: PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs);
6589: for (b = 0; b < fdof; ++b) {
6590: const PetscInt preind = perm ? foffs[f]+perm[b] : foffs[f]+b;
6591: const PetscInt ind = indperm ? indperm[preind] : preind;
6593: if ((cind < cfdof) && (b == fcdofs[cind])) {
6594: indices[ind] = -(foff+b+1);
6595: ++cind;
6596: } else {
6597: indices[ind] = foff+b-cind;
6598: }
6599: }
6600: }
6601: foffs[f] += fdof;
6602: }
6603: return(0);
6604: }
6606: PetscErrorCode DMPlexAnchorsModifyMat(DM dm, PetscSection section, PetscInt numPoints, PetscInt numIndices, const PetscInt points[], const PetscInt ***perms, const PetscScalar values[], PetscInt *outNumPoints, PetscInt *outNumIndices, PetscInt *outPoints[], PetscScalar *outValues[], PetscInt offsets[], PetscBool multiplyLeft)
6607: {
6608: Mat cMat;
6609: PetscSection aSec, cSec;
6610: IS aIS;
6611: PetscInt aStart = -1, aEnd = -1;
6612: const PetscInt *anchors;
6613: PetscInt numFields, f, p, q, newP = 0;
6614: PetscInt newNumPoints = 0, newNumIndices = 0;
6615: PetscInt *newPoints, *indices, *newIndices;
6616: PetscInt maxAnchor, maxDof;
6617: PetscInt newOffsets[32];
6618: PetscInt *pointMatOffsets[32];
6619: PetscInt *newPointOffsets[32];
6620: PetscScalar *pointMat[32];
6621: PetscScalar *newValues=NULL,*tmpValues;
6622: PetscBool anyConstrained = PETSC_FALSE;
6623: PetscErrorCode ierr;
6628: PetscSectionGetNumFields(section, &numFields);
6630: DMPlexGetAnchors(dm,&aSec,&aIS);
6631: /* if there are point-to-point constraints */
6632: if (aSec) {
6633: PetscArrayzero(newOffsets, 32);
6634: ISGetIndices(aIS,&anchors);
6635: PetscSectionGetChart(aSec,&aStart,&aEnd);
6636: /* figure out how many points are going to be in the new element matrix
6637: * (we allow double counting, because it's all just going to be summed
6638: * into the global matrix anyway) */
6639: for (p = 0; p < 2*numPoints; p+=2) {
6640: PetscInt b = points[p];
6641: PetscInt bDof = 0, bSecDof;
6643: PetscSectionGetDof(section,b,&bSecDof);
6644: if (!bSecDof) {
6645: continue;
6646: }
6647: if (b >= aStart && b < aEnd) {
6648: PetscSectionGetDof(aSec,b,&bDof);
6649: }
6650: if (bDof) {
6651: /* this point is constrained */
6652: /* it is going to be replaced by its anchors */
6653: PetscInt bOff, q;
6655: anyConstrained = PETSC_TRUE;
6656: newNumPoints += bDof;
6657: PetscSectionGetOffset(aSec,b,&bOff);
6658: for (q = 0; q < bDof; q++) {
6659: PetscInt a = anchors[bOff + q];
6660: PetscInt aDof;
6662: PetscSectionGetDof(section,a,&aDof);
6663: newNumIndices += aDof;
6664: for (f = 0; f < numFields; ++f) {
6665: PetscInt fDof;
6667: PetscSectionGetFieldDof(section, a, f, &fDof);
6668: newOffsets[f+1] += fDof;
6669: }
6670: }
6671: }
6672: else {
6673: /* this point is not constrained */
6674: newNumPoints++;
6675: newNumIndices += bSecDof;
6676: for (f = 0; f < numFields; ++f) {
6677: PetscInt fDof;
6679: PetscSectionGetFieldDof(section, b, f, &fDof);
6680: newOffsets[f+1] += fDof;
6681: }
6682: }
6683: }
6684: }
6685: if (!anyConstrained) {
6686: if (outNumPoints) *outNumPoints = 0;
6687: if (outNumIndices) *outNumIndices = 0;
6688: if (outPoints) *outPoints = NULL;
6689: if (outValues) *outValues = NULL;
6690: if (aSec) {ISRestoreIndices(aIS,&anchors);}
6691: return(0);
6692: }
6694: if (outNumPoints) *outNumPoints = newNumPoints;
6695: if (outNumIndices) *outNumIndices = newNumIndices;
6697: for (f = 0; f < numFields; ++f) newOffsets[f+1] += newOffsets[f];
6699: if (!outPoints && !outValues) {
6700: if (offsets) {
6701: for (f = 0; f <= numFields; f++) {
6702: offsets[f] = newOffsets[f];
6703: }
6704: }
6705: if (aSec) {ISRestoreIndices(aIS,&anchors);}
6706: return(0);
6707: }
6709: if (numFields && newOffsets[numFields] != newNumIndices) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid size for closure %D should be %D", newOffsets[numFields], newNumIndices);
6711: DMGetDefaultConstraints(dm, &cSec, &cMat);
6713: /* workspaces */
6714: if (numFields) {
6715: for (f = 0; f < numFields; f++) {
6716: DMGetWorkArray(dm,numPoints+1,MPIU_INT,&pointMatOffsets[f]);
6717: DMGetWorkArray(dm,numPoints+1,MPIU_INT,&newPointOffsets[f]);
6718: }
6719: }
6720: else {
6721: DMGetWorkArray(dm,numPoints+1,MPIU_INT,&pointMatOffsets[0]);
6722: DMGetWorkArray(dm,numPoints,MPIU_INT,&newPointOffsets[0]);
6723: }
6725: /* get workspaces for the point-to-point matrices */
6726: if (numFields) {
6727: PetscInt totalOffset, totalMatOffset;
6729: for (p = 0; p < numPoints; p++) {
6730: PetscInt b = points[2*p];
6731: PetscInt bDof = 0, bSecDof;
6733: PetscSectionGetDof(section,b,&bSecDof);
6734: if (!bSecDof) {
6735: for (f = 0; f < numFields; f++) {
6736: newPointOffsets[f][p + 1] = 0;
6737: pointMatOffsets[f][p + 1] = 0;
6738: }
6739: continue;
6740: }
6741: if (b >= aStart && b < aEnd) {
6742: PetscSectionGetDof(aSec, b, &bDof);
6743: }
6744: if (bDof) {
6745: for (f = 0; f < numFields; f++) {
6746: PetscInt fDof, q, bOff, allFDof = 0;
6748: PetscSectionGetFieldDof(section, b, f, &fDof);
6749: PetscSectionGetOffset(aSec, b, &bOff);
6750: for (q = 0; q < bDof; q++) {
6751: PetscInt a = anchors[bOff + q];
6752: PetscInt aFDof;
6754: PetscSectionGetFieldDof(section, a, f, &aFDof);
6755: allFDof += aFDof;
6756: }
6757: newPointOffsets[f][p+1] = allFDof;
6758: pointMatOffsets[f][p+1] = fDof * allFDof;
6759: }
6760: }
6761: else {
6762: for (f = 0; f < numFields; f++) {
6763: PetscInt fDof;
6765: PetscSectionGetFieldDof(section, b, f, &fDof);
6766: newPointOffsets[f][p+1] = fDof;
6767: pointMatOffsets[f][p+1] = 0;
6768: }
6769: }
6770: }
6771: for (f = 0, totalOffset = 0, totalMatOffset = 0; f < numFields; f++) {
6772: newPointOffsets[f][0] = totalOffset;
6773: pointMatOffsets[f][0] = totalMatOffset;
6774: for (p = 0; p < numPoints; p++) {
6775: newPointOffsets[f][p+1] += newPointOffsets[f][p];
6776: pointMatOffsets[f][p+1] += pointMatOffsets[f][p];
6777: }
6778: totalOffset = newPointOffsets[f][numPoints];
6779: totalMatOffset = pointMatOffsets[f][numPoints];
6780: DMGetWorkArray(dm,pointMatOffsets[f][numPoints],MPIU_SCALAR,&pointMat[f]);
6781: }
6782: }
6783: else {
6784: for (p = 0; p < numPoints; p++) {
6785: PetscInt b = points[2*p];
6786: PetscInt bDof = 0, bSecDof;
6788: PetscSectionGetDof(section,b,&bSecDof);
6789: if (!bSecDof) {
6790: newPointOffsets[0][p + 1] = 0;
6791: pointMatOffsets[0][p + 1] = 0;
6792: continue;
6793: }
6794: if (b >= aStart && b < aEnd) {
6795: PetscSectionGetDof(aSec, b, &bDof);
6796: }
6797: if (bDof) {
6798: PetscInt bOff, q, allDof = 0;
6800: PetscSectionGetOffset(aSec, b, &bOff);
6801: for (q = 0; q < bDof; q++) {
6802: PetscInt a = anchors[bOff + q], aDof;
6804: PetscSectionGetDof(section, a, &aDof);
6805: allDof += aDof;
6806: }
6807: newPointOffsets[0][p+1] = allDof;
6808: pointMatOffsets[0][p+1] = bSecDof * allDof;
6809: }
6810: else {
6811: newPointOffsets[0][p+1] = bSecDof;
6812: pointMatOffsets[0][p+1] = 0;
6813: }
6814: }
6815: newPointOffsets[0][0] = 0;
6816: pointMatOffsets[0][0] = 0;
6817: for (p = 0; p < numPoints; p++) {
6818: newPointOffsets[0][p+1] += newPointOffsets[0][p];
6819: pointMatOffsets[0][p+1] += pointMatOffsets[0][p];
6820: }
6821: DMGetWorkArray(dm,pointMatOffsets[0][numPoints],MPIU_SCALAR,&pointMat[0]);
6822: }
6824: /* output arrays */
6825: DMGetWorkArray(dm,2*newNumPoints,MPIU_INT,&newPoints);
6827: /* get the point-to-point matrices; construct newPoints */
6828: PetscSectionGetMaxDof(aSec, &maxAnchor);
6829: PetscSectionGetMaxDof(section, &maxDof);
6830: DMGetWorkArray(dm,maxDof,MPIU_INT,&indices);
6831: DMGetWorkArray(dm,maxAnchor*maxDof,MPIU_INT,&newIndices);
6832: if (numFields) {
6833: for (p = 0, newP = 0; p < numPoints; p++) {
6834: PetscInt b = points[2*p];
6835: PetscInt o = points[2*p+1];
6836: PetscInt bDof = 0, bSecDof;
6838: PetscSectionGetDof(section, b, &bSecDof);
6839: if (!bSecDof) {
6840: continue;
6841: }
6842: if (b >= aStart && b < aEnd) {
6843: PetscSectionGetDof(aSec, b, &bDof);
6844: }
6845: if (bDof) {
6846: PetscInt fStart[32], fEnd[32], fAnchorStart[32], fAnchorEnd[32], bOff, q;
6848: fStart[0] = 0;
6849: fEnd[0] = 0;
6850: for (f = 0; f < numFields; f++) {
6851: PetscInt fDof;
6853: PetscSectionGetFieldDof(cSec, b, f, &fDof);
6854: fStart[f+1] = fStart[f] + fDof;
6855: fEnd[f+1] = fStart[f+1];
6856: }
6857: PetscSectionGetOffset(cSec, b, &bOff);
6858: DMPlexGetIndicesPointFields_Internal(cSec, PETSC_TRUE, b, bOff, fEnd, PETSC_TRUE, perms, p, NULL, indices);
6860: fAnchorStart[0] = 0;
6861: fAnchorEnd[0] = 0;
6862: for (f = 0; f < numFields; f++) {
6863: PetscInt fDof = newPointOffsets[f][p + 1] - newPointOffsets[f][p];
6865: fAnchorStart[f+1] = fAnchorStart[f] + fDof;
6866: fAnchorEnd[f+1] = fAnchorStart[f + 1];
6867: }
6868: PetscSectionGetOffset(aSec, b, &bOff);
6869: for (q = 0; q < bDof; q++) {
6870: PetscInt a = anchors[bOff + q], aOff;
6872: /* we take the orientation of ap into account in the order that we constructed the indices above: the newly added points have no orientation */
6873: newPoints[2*(newP + q)] = a;
6874: newPoints[2*(newP + q) + 1] = 0;
6875: PetscSectionGetOffset(section, a, &aOff);
6876: DMPlexGetIndicesPointFields_Internal(section, PETSC_TRUE, a, aOff, fAnchorEnd, PETSC_TRUE, NULL, -1, NULL, newIndices);
6877: }
6878: newP += bDof;
6880: if (outValues) {
6881: /* get the point-to-point submatrix */
6882: for (f = 0; f < numFields; f++) {
6883: MatGetValues(cMat,fEnd[f]-fStart[f],indices + fStart[f],fAnchorEnd[f] - fAnchorStart[f],newIndices + fAnchorStart[f],pointMat[f] + pointMatOffsets[f][p]);
6884: }
6885: }
6886: }
6887: else {
6888: newPoints[2 * newP] = b;
6889: newPoints[2 * newP + 1] = o;
6890: newP++;
6891: }
6892: }
6893: } else {
6894: for (p = 0; p < numPoints; p++) {
6895: PetscInt b = points[2*p];
6896: PetscInt o = points[2*p+1];
6897: PetscInt bDof = 0, bSecDof;
6899: PetscSectionGetDof(section, b, &bSecDof);
6900: if (!bSecDof) {
6901: continue;
6902: }
6903: if (b >= aStart && b < aEnd) {
6904: PetscSectionGetDof(aSec, b, &bDof);
6905: }
6906: if (bDof) {
6907: PetscInt bEnd = 0, bAnchorEnd = 0, bOff;
6909: PetscSectionGetOffset(cSec, b, &bOff);
6910: DMPlexGetIndicesPoint_Internal(cSec, PETSC_TRUE, b, bOff, &bEnd, PETSC_TRUE, (perms && perms[0]) ? perms[0][p] : NULL, NULL, indices);
6912: PetscSectionGetOffset (aSec, b, &bOff);
6913: for (q = 0; q < bDof; q++) {
6914: PetscInt a = anchors[bOff + q], aOff;
6916: /* we take the orientation of ap into account in the order that we constructed the indices above: the newly added points have no orientation */
6918: newPoints[2*(newP + q)] = a;
6919: newPoints[2*(newP + q) + 1] = 0;
6920: PetscSectionGetOffset(section, a, &aOff);
6921: DMPlexGetIndicesPoint_Internal(section, PETSC_TRUE, a, aOff, &bAnchorEnd, PETSC_TRUE, NULL, NULL, newIndices);
6922: }
6923: newP += bDof;
6925: /* get the point-to-point submatrix */
6926: if (outValues) {
6927: MatGetValues(cMat,bEnd,indices,bAnchorEnd,newIndices,pointMat[0] + pointMatOffsets[0][p]);
6928: }
6929: }
6930: else {
6931: newPoints[2 * newP] = b;
6932: newPoints[2 * newP + 1] = o;
6933: newP++;
6934: }
6935: }
6936: }
6938: if (outValues) {
6939: DMGetWorkArray(dm,newNumIndices*numIndices,MPIU_SCALAR,&tmpValues);
6940: PetscArrayzero(tmpValues,newNumIndices*numIndices);
6941: /* multiply constraints on the right */
6942: if (numFields) {
6943: for (f = 0; f < numFields; f++) {
6944: PetscInt oldOff = offsets[f];
6946: for (p = 0; p < numPoints; p++) {
6947: PetscInt cStart = newPointOffsets[f][p];
6948: PetscInt b = points[2 * p];
6949: PetscInt c, r, k;
6950: PetscInt dof;
6952: PetscSectionGetFieldDof(section,b,f,&dof);
6953: if (!dof) {
6954: continue;
6955: }
6956: if (pointMatOffsets[f][p] < pointMatOffsets[f][p + 1]) {
6957: PetscInt nCols = newPointOffsets[f][p+1]-cStart;
6958: const PetscScalar *mat = pointMat[f] + pointMatOffsets[f][p];
6960: for (r = 0; r < numIndices; r++) {
6961: for (c = 0; c < nCols; c++) {
6962: for (k = 0; k < dof; k++) {
6963: tmpValues[r * newNumIndices + cStart + c] += values[r * numIndices + oldOff + k] * mat[k * nCols + c];
6964: }
6965: }
6966: }
6967: }
6968: else {
6969: /* copy this column as is */
6970: for (r = 0; r < numIndices; r++) {
6971: for (c = 0; c < dof; c++) {
6972: tmpValues[r * newNumIndices + cStart + c] = values[r * numIndices + oldOff + c];
6973: }
6974: }
6975: }
6976: oldOff += dof;
6977: }
6978: }
6979: }
6980: else {
6981: PetscInt oldOff = 0;
6982: for (p = 0; p < numPoints; p++) {
6983: PetscInt cStart = newPointOffsets[0][p];
6984: PetscInt b = points[2 * p];
6985: PetscInt c, r, k;
6986: PetscInt dof;
6988: PetscSectionGetDof(section,b,&dof);
6989: if (!dof) {
6990: continue;
6991: }
6992: if (pointMatOffsets[0][p] < pointMatOffsets[0][p + 1]) {
6993: PetscInt nCols = newPointOffsets[0][p+1]-cStart;
6994: const PetscScalar *mat = pointMat[0] + pointMatOffsets[0][p];
6996: for (r = 0; r < numIndices; r++) {
6997: for (c = 0; c < nCols; c++) {
6998: for (k = 0; k < dof; k++) {
6999: tmpValues[r * newNumIndices + cStart + c] += mat[k * nCols + c] * values[r * numIndices + oldOff + k];
7000: }
7001: }
7002: }
7003: }
7004: else {
7005: /* copy this column as is */
7006: for (r = 0; r < numIndices; r++) {
7007: for (c = 0; c < dof; c++) {
7008: tmpValues[r * newNumIndices + cStart + c] = values[r * numIndices + oldOff + c];
7009: }
7010: }
7011: }
7012: oldOff += dof;
7013: }
7014: }
7016: if (multiplyLeft) {
7017: DMGetWorkArray(dm,newNumIndices*newNumIndices,MPIU_SCALAR,&newValues);
7018: PetscArrayzero(newValues,newNumIndices*newNumIndices);
7019: /* multiply constraints transpose on the left */
7020: if (numFields) {
7021: for (f = 0; f < numFields; f++) {
7022: PetscInt oldOff = offsets[f];
7024: for (p = 0; p < numPoints; p++) {
7025: PetscInt rStart = newPointOffsets[f][p];
7026: PetscInt b = points[2 * p];
7027: PetscInt c, r, k;
7028: PetscInt dof;
7030: PetscSectionGetFieldDof(section,b,f,&dof);
7031: if (pointMatOffsets[f][p] < pointMatOffsets[f][p + 1]) {
7032: PetscInt nRows = newPointOffsets[f][p+1]-rStart;
7033: const PetscScalar *PETSC_RESTRICT mat = pointMat[f] + pointMatOffsets[f][p];
7035: for (r = 0; r < nRows; r++) {
7036: for (c = 0; c < newNumIndices; c++) {
7037: for (k = 0; k < dof; k++) {
7038: newValues[(rStart + r) * newNumIndices + c] += mat[k * nRows + r] * tmpValues[(oldOff + k) * newNumIndices + c];
7039: }
7040: }
7041: }
7042: }
7043: else {
7044: /* copy this row as is */
7045: for (r = 0; r < dof; r++) {
7046: for (c = 0; c < newNumIndices; c++) {
7047: newValues[(rStart + r) * newNumIndices + c] = tmpValues[(oldOff + r) * newNumIndices + c];
7048: }
7049: }
7050: }
7051: oldOff += dof;
7052: }
7053: }
7054: }
7055: else {
7056: PetscInt oldOff = 0;
7058: for (p = 0; p < numPoints; p++) {
7059: PetscInt rStart = newPointOffsets[0][p];
7060: PetscInt b = points[2 * p];
7061: PetscInt c, r, k;
7062: PetscInt dof;
7064: PetscSectionGetDof(section,b,&dof);
7065: if (pointMatOffsets[0][p] < pointMatOffsets[0][p + 1]) {
7066: PetscInt nRows = newPointOffsets[0][p+1]-rStart;
7067: const PetscScalar *PETSC_RESTRICT mat = pointMat[0] + pointMatOffsets[0][p];
7069: for (r = 0; r < nRows; r++) {
7070: for (c = 0; c < newNumIndices; c++) {
7071: for (k = 0; k < dof; k++) {
7072: newValues[(rStart + r) * newNumIndices + c] += mat[k * nRows + r] * tmpValues[(oldOff + k) * newNumIndices + c];
7073: }
7074: }
7075: }
7076: }
7077: else {
7078: /* copy this row as is */
7079: for (r = 0; r < dof; r++) {
7080: for (c = 0; c < newNumIndices; c++) {
7081: newValues[(rStart + r) * newNumIndices + c] = tmpValues[(oldOff + r) * newNumIndices + c];
7082: }
7083: }
7084: }
7085: oldOff += dof;
7086: }
7087: }
7089: DMRestoreWorkArray(dm,newNumIndices*numIndices,MPIU_SCALAR,&tmpValues);
7090: }
7091: else {
7092: newValues = tmpValues;
7093: }
7094: }
7096: /* clean up */
7097: DMRestoreWorkArray(dm,maxDof,MPIU_INT,&indices);
7098: DMRestoreWorkArray(dm,maxAnchor*maxDof,MPIU_INT,&newIndices);
7100: if (numFields) {
7101: for (f = 0; f < numFields; f++) {
7102: DMRestoreWorkArray(dm,pointMatOffsets[f][numPoints],MPIU_SCALAR,&pointMat[f]);
7103: DMRestoreWorkArray(dm,numPoints+1,MPIU_INT,&pointMatOffsets[f]);
7104: DMRestoreWorkArray(dm,numPoints+1,MPIU_INT,&newPointOffsets[f]);
7105: }
7106: }
7107: else {
7108: DMRestoreWorkArray(dm,pointMatOffsets[0][numPoints],MPIU_SCALAR,&pointMat[0]);
7109: DMRestoreWorkArray(dm,numPoints+1,MPIU_INT,&pointMatOffsets[0]);
7110: DMRestoreWorkArray(dm,numPoints+1,MPIU_INT,&newPointOffsets[0]);
7111: }
7112: ISRestoreIndices(aIS,&anchors);
7114: /* output */
7115: if (outPoints) {
7116: *outPoints = newPoints;
7117: }
7118: else {
7119: DMRestoreWorkArray(dm,2*newNumPoints,MPIU_INT,&newPoints);
7120: }
7121: if (outValues) {
7122: *outValues = newValues;
7123: }
7124: for (f = 0; f <= numFields; f++) {
7125: offsets[f] = newOffsets[f];
7126: }
7127: return(0);
7128: }
7130: /*@C
7131: DMPlexGetClosureIndices - Gets the global dof indices associated with the closure of the given point within the provided sections.
7133: Not collective
7135: Input Parameters:
7136: + dm - The DM
7137: . section - The PetscSection describing the points (a local section)
7138: . idxSection - The PetscSection from which to obtain indices (may be local or global)
7139: . point - The point defining the closure
7140: - useClPerm - Use the closure point permutation if available
7142: Output Parameters:
7143: + numIndices - The number of dof indices in the closure of point with the input sections
7144: . indices - The dof indices
7145: . outOffsets - Array to write the field offsets into, or NULL
7146: - values - The input values, which may be modified if sign flips are induced by the point symmetries, or NULL
7148: Notes:
7149: Must call DMPlexRestoreClosureIndices() to free allocated memory
7151: If idxSection is global, any constrained dofs (see DMAddBoundary(), for example) will get negative indices. The value
7152: of those indices is not significant. If idxSection is local, the constrained dofs will yield the involution -(idx+1)
7153: of their index in a local vector. A caller who does not wish to distinguish those points may recover the nonnegative
7154: indices via involution, -(-(idx+1)+1)==idx. Local indices are provided when idxSection == section, otherwise global
7155: indices (with the above semantics) are implied.
7157: Level: advanced
7159: .seealso DMPlexRestoreClosureIndices(), DMPlexVecGetClosure(), DMPlexMatSetClosure(), DMGetLocalSection(), DMGetGlobalSection()
7160: @*/
7161: PetscErrorCode DMPlexGetClosureIndices(DM dm, PetscSection section, PetscSection idxSection, PetscInt point, PetscBool useClPerm,
7162: PetscInt *numIndices, PetscInt *indices[], PetscInt outOffsets[], PetscScalar *values[])
7163: {
7164: /* Closure ordering */
7165: PetscSection clSection;
7166: IS clPoints;
7167: const PetscInt *clp;
7168: PetscInt *points;
7169: const PetscInt *clperm = NULL;
7170: /* Dof permutation and sign flips */
7171: const PetscInt **perms[32] = {NULL};
7172: const PetscScalar **flips[32] = {NULL};
7173: PetscScalar *valCopy = NULL;
7174: /* Hanging node constraints */
7175: PetscInt *pointsC = NULL;
7176: PetscScalar *valuesC = NULL;
7177: PetscInt NclC, NiC;
7179: PetscInt *idx;
7180: PetscInt Nf, Ncl, Ni = 0, offsets[32], p, f;
7181: PetscBool isLocal = (section == idxSection) ? PETSC_TRUE : PETSC_FALSE;
7182: PetscErrorCode ierr;
7192: PetscSectionGetNumFields(section, &Nf);
7193: if (Nf > 31) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Number of fields %D limited to 31", Nf);
7194: PetscArrayzero(offsets, 32);
7195: /* 1) Get points in closure */
7196: DMPlexGetCompressedClosure(dm, section, point, &Ncl, &points, &clSection, &clPoints, &clp);
7197: if (useClPerm) {
7198: PetscInt depth, clsize;
7199: DMPlexGetPointDepth(dm, point, &depth);
7200: for (clsize=0,p=0; p<Ncl; p++) {
7201: PetscInt dof;
7202: PetscSectionGetDof(section, points[2*p], &dof);
7203: clsize += dof;
7204: }
7205: PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject) dm, depth, clsize, &clperm);
7206: }
7207: /* 2) Get number of indices on these points and field offsets from section */
7208: for (p = 0; p < Ncl*2; p += 2) {
7209: PetscInt dof, fdof;
7211: PetscSectionGetDof(section, points[p], &dof);
7212: for (f = 0; f < Nf; ++f) {
7213: PetscSectionGetFieldDof(section, points[p], f, &fdof);
7214: offsets[f+1] += fdof;
7215: }
7216: Ni += dof;
7217: }
7218: for (f = 1; f < Nf; ++f) offsets[f+1] += offsets[f];
7219: if (Nf && offsets[Nf] != Ni) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Invalid size for closure %D should be %D", offsets[Nf], Ni);
7220: /* 3) Get symmetries and sign flips. Apply sign flips to values if passed in (only works for square values matrix) */
7221: for (f = 0; f < PetscMax(1, Nf); ++f) {
7222: if (Nf) {PetscSectionGetFieldPointSyms(section, f, Ncl, points, &perms[f], &flips[f]);}
7223: else {PetscSectionGetPointSyms(section, Ncl, points, &perms[f], &flips[f]);}
7224: /* may need to apply sign changes to the element matrix */
7225: if (values && flips[f]) {
7226: PetscInt foffset = offsets[f];
7228: for (p = 0; p < Ncl; ++p) {
7229: PetscInt pnt = points[2*p], fdof;
7230: const PetscScalar *flip = flips[f] ? flips[f][p] : NULL;
7232: if (!Nf) {PetscSectionGetDof(section, pnt, &fdof);}
7233: else {PetscSectionGetFieldDof(section, pnt, f, &fdof);}
7234: if (flip) {
7235: PetscInt i, j, k;
7237: if (!valCopy) {
7238: DMGetWorkArray(dm, Ni*Ni, MPIU_SCALAR, &valCopy);
7239: for (j = 0; j < Ni * Ni; ++j) valCopy[j] = (*values)[j];
7240: *values = valCopy;
7241: }
7242: for (i = 0; i < fdof; ++i) {
7243: PetscScalar fval = flip[i];
7245: for (k = 0; k < Ni; ++k) {
7246: valCopy[Ni * (foffset + i) + k] *= fval;
7247: valCopy[Ni * k + (foffset + i)] *= fval;
7248: }
7249: }
7250: }
7251: foffset += fdof;
7252: }
7253: }
7254: }
7255: /* 4) Apply hanging node constraints. Get new symmetries and replace all storage with constrained storage */
7256: DMPlexAnchorsModifyMat(dm, section, Ncl, Ni, points, perms, values ? *values : NULL, &NclC, &NiC, &pointsC, values ? &valuesC : NULL, offsets, PETSC_TRUE);
7257: if (NclC) {
7258: if (valCopy) {DMRestoreWorkArray(dm, Ni*Ni, MPIU_SCALAR, &valCopy);}
7259: for (f = 0; f < PetscMax(1, Nf); ++f) {
7260: if (Nf) {PetscSectionRestoreFieldPointSyms(section, f, Ncl, points, &perms[f], &flips[f]);}
7261: else {PetscSectionRestorePointSyms(section, Ncl, points, &perms[f], &flips[f]);}
7262: }
7263: for (f = 0; f < PetscMax(1, Nf); ++f) {
7264: if (Nf) {PetscSectionGetFieldPointSyms(section, f, NclC, pointsC, &perms[f], &flips[f]);}
7265: else {PetscSectionGetPointSyms(section, NclC, pointsC, &perms[f], &flips[f]);}
7266: }
7267: DMPlexRestoreCompressedClosure(dm, section, point, &Ncl, &points, &clSection, &clPoints, &clp);
7268: Ncl = NclC;
7269: Ni = NiC;
7270: points = pointsC;
7271: if (values) *values = valuesC;
7272: }
7273: /* 5) Calculate indices */
7274: DMGetWorkArray(dm, Ni, MPIU_INT, &idx);
7275: if (Nf) {
7276: PetscInt idxOff;
7277: PetscBool useFieldOffsets;
7279: if (outOffsets) {for (f = 0; f <= Nf; f++) outOffsets[f] = offsets[f];}
7280: PetscSectionGetUseFieldOffsets(idxSection, &useFieldOffsets);
7281: if (useFieldOffsets) {
7282: for (p = 0; p < Ncl; ++p) {
7283: const PetscInt pnt = points[p*2];
7285: DMPlexGetIndicesPointFieldsSplit_Internal(section, idxSection, pnt, offsets, perms, p, clperm, idx);
7286: }
7287: } else {
7288: for (p = 0; p < Ncl; ++p) {
7289: const PetscInt pnt = points[p*2];
7291: PetscSectionGetOffset(idxSection, pnt, &idxOff);
7292: /* Note that we pass a local section even though we're using global offsets. This is because global sections do
7293: * not (at the time of this writing) have fields set. They probably should, in which case we would pass the
7294: * global section. */
7295: DMPlexGetIndicesPointFields_Internal(section, isLocal, pnt, idxOff < 0 ? -(idxOff+1) : idxOff, offsets, PETSC_FALSE, perms, p, clperm, idx);
7296: }
7297: }
7298: } else {
7299: PetscInt off = 0, idxOff;
7301: for (p = 0; p < Ncl; ++p) {
7302: const PetscInt pnt = points[p*2];
7303: const PetscInt *perm = perms[0] ? perms[0][p] : NULL;
7305: PetscSectionGetOffset(idxSection, pnt, &idxOff);
7306: /* Note that we pass a local section even though we're using global offsets. This is because global sections do
7307: * not (at the time of this writing) have fields set. They probably should, in which case we would pass the global section. */
7308: DMPlexGetIndicesPoint_Internal(section, isLocal, pnt, idxOff < 0 ? -(idxOff+1) : idxOff, &off, PETSC_FALSE, perm, clperm, idx);
7309: }
7310: }
7311: /* 6) Cleanup */
7312: for (f = 0; f < PetscMax(1, Nf); ++f) {
7313: if (Nf) {PetscSectionRestoreFieldPointSyms(section, f, Ncl, points, &perms[f], &flips[f]);}
7314: else {PetscSectionRestorePointSyms(section, Ncl, points, &perms[f], &flips[f]);}
7315: }
7316: if (NclC) {
7317: DMRestoreWorkArray(dm, NclC*2, MPIU_INT, &pointsC);
7318: } else {
7319: DMPlexRestoreCompressedClosure(dm, section, point, &Ncl, &points, &clSection, &clPoints, &clp);
7320: }
7322: if (numIndices) *numIndices = Ni;
7323: if (indices) *indices = idx;
7324: return(0);
7325: }
7327: /*@C
7328: DMPlexRestoreClosureIndices - Restores the global dof indices associated with the closure of the given point within the provided sections.
7330: Not collective
7332: Input Parameters:
7333: + dm - The DM
7334: . section - The PetscSection describing the points (a local section)
7335: . idxSection - The PetscSection from which to obtain indices (may be local or global)
7336: . point - The point defining the closure
7337: - useClPerm - Use the closure point permutation if available
7339: Output Parameters:
7340: + numIndices - The number of dof indices in the closure of point with the input sections
7341: . indices - The dof indices
7342: . outOffsets - Array to write the field offsets into, or NULL
7343: - values - The input values, which may be modified if sign flips are induced by the point symmetries, or NULL
7345: Notes:
7346: If values were modified, the user is responsible for calling DMRestoreWorkArray(dm, 0, MPIU_SCALAR, &values).
7348: If idxSection is global, any constrained dofs (see DMAddBoundary(), for example) will get negative indices. The value
7349: of those indices is not significant. If idxSection is local, the constrained dofs will yield the involution -(idx+1)
7350: of their index in a local vector. A caller who does not wish to distinguish those points may recover the nonnegative
7351: indices via involution, -(-(idx+1)+1)==idx. Local indices are provided when idxSection == section, otherwise global
7352: indices (with the above semantics) are implied.
7354: Level: advanced
7356: .seealso DMPlexGetClosureIndices(), DMPlexVecGetClosure(), DMPlexMatSetClosure(), DMGetLocalSection(), DMGetGlobalSection()
7357: @*/
7358: PetscErrorCode DMPlexRestoreClosureIndices(DM dm, PetscSection section, PetscSection idxSection, PetscInt point, PetscBool useClPerm,
7359: PetscInt *numIndices, PetscInt *indices[], PetscInt outOffsets[], PetscScalar *values[])
7360: {
7366: DMRestoreWorkArray(dm, 0, MPIU_INT, indices);
7367: return(0);
7368: }
7370: /*@C
7371: DMPlexMatSetClosure - Set an array of the values on the closure of 'point'
7373: Not collective
7375: Input Parameters:
7376: + dm - The DM
7377: . section - The section describing the layout in v, or NULL to use the default section
7378: . globalSection - The section describing the layout in v, or NULL to use the default global section
7379: . A - The matrix
7380: . point - The point in the DM
7381: . values - The array of values
7382: - mode - The insert mode, where INSERT_ALL_VALUES and ADD_ALL_VALUES also overwrite boundary conditions
7384: Fortran Notes:
7385: This routine is only available in Fortran 90, and you must include petsc.h90 in your code.
7387: Level: intermediate
7389: .seealso DMPlexMatSetClosureGeneral(), DMPlexVecGetClosure(), DMPlexVecSetClosure()
7390: @*/
7391: PetscErrorCode DMPlexMatSetClosure(DM dm, PetscSection section, PetscSection globalSection, Mat A, PetscInt point, const PetscScalar values[], InsertMode mode)
7392: {
7393: DM_Plex *mesh = (DM_Plex*) dm->data;
7394: PetscInt *indices;
7395: PetscInt numIndices;
7396: const PetscScalar *valuesOrig = values;
7397: PetscErrorCode ierr;
7401: if (!section) {DMGetLocalSection(dm, §ion);}
7403: if (!globalSection) {DMGetGlobalSection(dm, &globalSection);}
7407: DMPlexGetClosureIndices(dm, section, globalSection, point, PETSC_TRUE, &numIndices, &indices, NULL, (PetscScalar **) &values);
7409: if (mesh->printSetValues) {DMPlexPrintMatSetValues(PETSC_VIEWER_STDOUT_SELF, A, point, numIndices, indices, 0, NULL, values);}
7410: MatSetValues(A, numIndices, indices, numIndices, indices, values, mode);
7411: if (ierr) {
7412: PetscMPIInt rank;
7413: PetscErrorCode ierr2;
7415: ierr2 = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank);CHKERRMPI(ierr2);
7416: ierr2 = (*PetscErrorPrintf)("[%d]ERROR in DMPlexMatSetClosure\n", rank);CHKERRQ(ierr2);
7417: ierr2 = DMPlexPrintMatSetValues(PETSC_VIEWER_STDERR_SELF, A, point, numIndices, indices, 0, NULL, values);CHKERRQ(ierr2);
7418: ierr2 = DMPlexRestoreClosureIndices(dm, section, globalSection, point, PETSC_TRUE, &numIndices, &indices, NULL, (PetscScalar **) &values);CHKERRQ(ierr2);
7419: if (values != valuesOrig) {ierr2 = DMRestoreWorkArray(dm, 0, MPIU_SCALAR, &values);CHKERRQ(ierr2);}
7420:
7421: }
7422: if (mesh->printFEM > 1) {
7423: PetscInt i;
7424: PetscPrintf(PETSC_COMM_SELF, " Indices:");
7425: for (i = 0; i < numIndices; ++i) {PetscPrintf(PETSC_COMM_SELF, " %D", indices[i]);}
7426: PetscPrintf(PETSC_COMM_SELF, "\n");
7427: }
7429: DMPlexRestoreClosureIndices(dm, section, globalSection, point, PETSC_TRUE, &numIndices, &indices, NULL, (PetscScalar **) &values);
7430: if (values != valuesOrig) {DMRestoreWorkArray(dm, 0, MPIU_SCALAR, &values);}
7431: return(0);
7432: }
7434: /*@C
7435: DMPlexMatSetClosure - Set an array of the values on the closure of 'point' using a different row and column section
7437: Not collective
7439: Input Parameters:
7440: + dmRow - The DM for the row fields
7441: . sectionRow - The section describing the layout, or NULL to use the default section in dmRow
7442: . globalSectionRow - The section describing the layout, or NULL to use the default global section in dmRow
7443: . dmCol - The DM for the column fields
7444: . sectionCol - The section describing the layout, or NULL to use the default section in dmCol
7445: . globalSectionCol - The section describing the layout, or NULL to use the default global section in dmCol
7446: . A - The matrix
7447: . point - The point in the DMs
7448: . values - The array of values
7449: - mode - The insert mode, where INSERT_ALL_VALUES and ADD_ALL_VALUES also overwrite boundary conditions
7451: Level: intermediate
7453: .seealso DMPlexMatSetClosure(), DMPlexVecGetClosure(), DMPlexVecSetClosure()
7454: @*/
7455: PetscErrorCode DMPlexMatSetClosureGeneral(DM dmRow, PetscSection sectionRow, PetscSection globalSectionRow, DM dmCol, PetscSection sectionCol, PetscSection globalSectionCol, Mat A, PetscInt point, const PetscScalar values[], InsertMode mode)
7456: {
7457: DM_Plex *mesh = (DM_Plex*) dmRow->data;
7458: PetscInt *indicesRow, *indicesCol;
7459: PetscInt numIndicesRow, numIndicesCol;
7460: const PetscScalar *valuesOrig = values;
7461: PetscErrorCode ierr;
7465: if (!sectionRow) {DMGetLocalSection(dmRow, §ionRow);}
7467: if (!globalSectionRow) {DMGetGlobalSection(dmRow, &globalSectionRow);}
7470: if (!sectionCol) {DMGetLocalSection(dmCol, §ionCol);}
7472: if (!globalSectionCol) {DMGetGlobalSection(dmCol, &globalSectionCol);}
7476: DMPlexGetClosureIndices(dmRow, sectionRow, globalSectionRow, point, PETSC_TRUE, &numIndicesRow, &indicesRow, NULL, (PetscScalar **) &values);
7477: DMPlexGetClosureIndices(dmCol, sectionCol, globalSectionCol, point, PETSC_TRUE, &numIndicesCol, &indicesCol, NULL, (PetscScalar **) &values);
7479: if (mesh->printSetValues) {DMPlexPrintMatSetValues(PETSC_VIEWER_STDOUT_SELF, A, point, numIndicesRow, indicesRow, numIndicesCol, indicesCol, values);}
7480: MatSetValues(A, numIndicesRow, indicesRow, numIndicesCol, indicesCol, values, mode);
7481: if (ierr) {
7482: PetscMPIInt rank;
7483: PetscErrorCode ierr2;
7485: ierr2 = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank);CHKERRMPI(ierr2);
7486: ierr2 = (*PetscErrorPrintf)("[%d]ERROR in DMPlexMatSetClosure\n", rank);CHKERRQ(ierr2);
7487: ierr2 = DMPlexPrintMatSetValues(PETSC_VIEWER_STDERR_SELF, A, point, numIndicesRow, indicesRow, numIndicesCol, indicesCol, values);CHKERRQ(ierr2);
7488: ierr2 = DMPlexRestoreClosureIndices(dmRow, sectionRow, globalSectionRow, point, PETSC_TRUE, &numIndicesRow, &indicesRow, NULL, (PetscScalar **) &values);CHKERRQ(ierr2);
7489: ierr2 = DMPlexRestoreClosureIndices(dmCol, sectionCol, globalSectionCol, point, PETSC_TRUE, &numIndicesCol, &indicesRow, NULL, (PetscScalar **) &values);CHKERRQ(ierr2);
7490: if (values != valuesOrig) {ierr2 = DMRestoreWorkArray(dmRow, 0, MPIU_SCALAR, &values);CHKERRQ(ierr2);}
7491:
7492: }
7494: DMPlexRestoreClosureIndices(dmRow, sectionRow, globalSectionRow, point, PETSC_TRUE, &numIndicesRow, &indicesRow, NULL, (PetscScalar **) &values);
7495: DMPlexRestoreClosureIndices(dmCol, sectionCol, globalSectionCol, point, PETSC_TRUE, &numIndicesCol, &indicesCol, NULL, (PetscScalar **) &values);
7496: if (values != valuesOrig) {DMRestoreWorkArray(dmRow, 0, MPIU_SCALAR, &values);}
7497: return(0);
7498: }
7500: PetscErrorCode DMPlexMatSetClosureRefined(DM dmf, PetscSection fsection, PetscSection globalFSection, DM dmc, PetscSection csection, PetscSection globalCSection, Mat A, PetscInt point, const PetscScalar values[], InsertMode mode)
7501: {
7502: DM_Plex *mesh = (DM_Plex*) dmf->data;
7503: PetscInt *fpoints = NULL, *ftotpoints = NULL;
7504: PetscInt *cpoints = NULL;
7505: PetscInt *findices, *cindices;
7506: const PetscInt *fclperm = NULL, *cclperm = NULL; /* Closure permutations cannot work here */
7507: PetscInt foffsets[32], coffsets[32];
7508: DMPolytopeType ct;
7509: PetscInt numFields, numSubcells, maxFPoints, numFPoints, numCPoints, numFIndices, numCIndices, dof, off, globalOff, pStart, pEnd, p, q, r, s, f;
7510: PetscErrorCode ierr;
7515: if (!fsection) {DMGetLocalSection(dmf, &fsection);}
7517: if (!csection) {DMGetLocalSection(dmc, &csection);}
7519: if (!globalFSection) {DMGetGlobalSection(dmf, &globalFSection);}
7521: if (!globalCSection) {DMGetGlobalSection(dmc, &globalCSection);}
7524: PetscSectionGetNumFields(fsection, &numFields);
7525: if (numFields > 31) SETERRQ1(PetscObjectComm((PetscObject)dmf), PETSC_ERR_ARG_OUTOFRANGE, "Number of fields %D limited to 31", numFields);
7526: PetscArrayzero(foffsets, 32);
7527: PetscArrayzero(coffsets, 32);
7528: /* Column indices */
7529: DMPlexGetTransitiveClosure(dmc, point, PETSC_TRUE, &numCPoints, &cpoints);
7530: maxFPoints = numCPoints;
7531: /* Compress out points not in the section */
7532: /* TODO: Squeeze out points with 0 dof as well */
7533: PetscSectionGetChart(csection, &pStart, &pEnd);
7534: for (p = 0, q = 0; p < numCPoints*2; p += 2) {
7535: if ((cpoints[p] >= pStart) && (cpoints[p] < pEnd)) {
7536: cpoints[q*2] = cpoints[p];
7537: cpoints[q*2+1] = cpoints[p+1];
7538: ++q;
7539: }
7540: }
7541: numCPoints = q;
7542: for (p = 0, numCIndices = 0; p < numCPoints*2; p += 2) {
7543: PetscInt fdof;
7545: PetscSectionGetDof(csection, cpoints[p], &dof);
7546: if (!dof) continue;
7547: for (f = 0; f < numFields; ++f) {
7548: PetscSectionGetFieldDof(csection, cpoints[p], f, &fdof);
7549: coffsets[f+1] += fdof;
7550: }
7551: numCIndices += dof;
7552: }
7553: for (f = 1; f < numFields; ++f) coffsets[f+1] += coffsets[f];
7554: /* Row indices */
7555: DMPlexGetCellType(dmc, point, &ct);
7556: {
7557: DMPlexTransform tr;
7558: DMPolytopeType *rct;
7559: PetscInt *rsize, *rcone, *rornt, Nt;
7561: DMPlexTransformCreate(PETSC_COMM_SELF, &tr);
7562: DMPlexTransformSetType(tr, DMPLEXREFINEREGULAR);
7563: DMPlexTransformCellTransform(tr, ct, point, NULL, &Nt, &rct, &rsize, &rcone, &rornt);
7564: numSubcells = rsize[Nt-1];
7565: DMPlexTransformDestroy(&tr);
7566: }
7567: DMGetWorkArray(dmf, maxFPoints*2*numSubcells, MPIU_INT, &ftotpoints);
7568: for (r = 0, q = 0; r < numSubcells; ++r) {
7569: /* TODO Map from coarse to fine cells */
7570: DMPlexGetTransitiveClosure(dmf, point*numSubcells + r, PETSC_TRUE, &numFPoints, &fpoints);
7571: /* Compress out points not in the section */
7572: PetscSectionGetChart(fsection, &pStart, &pEnd);
7573: for (p = 0; p < numFPoints*2; p += 2) {
7574: if ((fpoints[p] >= pStart) && (fpoints[p] < pEnd)) {
7575: PetscSectionGetDof(fsection, fpoints[p], &dof);
7576: if (!dof) continue;
7577: for (s = 0; s < q; ++s) if (fpoints[p] == ftotpoints[s*2]) break;
7578: if (s < q) continue;
7579: ftotpoints[q*2] = fpoints[p];
7580: ftotpoints[q*2+1] = fpoints[p+1];
7581: ++q;
7582: }
7583: }
7584: DMPlexRestoreTransitiveClosure(dmf, point, PETSC_TRUE, &numFPoints, &fpoints);
7585: }
7586: numFPoints = q;
7587: for (p = 0, numFIndices = 0; p < numFPoints*2; p += 2) {
7588: PetscInt fdof;
7590: PetscSectionGetDof(fsection, ftotpoints[p], &dof);
7591: if (!dof) continue;
7592: for (f = 0; f < numFields; ++f) {
7593: PetscSectionGetFieldDof(fsection, ftotpoints[p], f, &fdof);
7594: foffsets[f+1] += fdof;
7595: }
7596: numFIndices += dof;
7597: }
7598: for (f = 1; f < numFields; ++f) foffsets[f+1] += foffsets[f];
7600: if (numFields && foffsets[numFields] != numFIndices) SETERRQ2(PetscObjectComm((PetscObject)dmf), PETSC_ERR_PLIB, "Invalid size for closure %D should be %D", foffsets[numFields], numFIndices);
7601: if (numFields && coffsets[numFields] != numCIndices) SETERRQ2(PetscObjectComm((PetscObject)dmc), PETSC_ERR_PLIB, "Invalid size for closure %D should be %D", coffsets[numFields], numCIndices);
7602: DMGetWorkArray(dmf, numFIndices, MPIU_INT, &findices);
7603: DMGetWorkArray(dmc, numCIndices, MPIU_INT, &cindices);
7604: if (numFields) {
7605: const PetscInt **permsF[32] = {NULL};
7606: const PetscInt **permsC[32] = {NULL};
7608: for (f = 0; f < numFields; f++) {
7609: PetscSectionGetFieldPointSyms(fsection,f,numFPoints,ftotpoints,&permsF[f],NULL);
7610: PetscSectionGetFieldPointSyms(csection,f,numCPoints,cpoints,&permsC[f],NULL);
7611: }
7612: for (p = 0; p < numFPoints; p++) {
7613: PetscSectionGetOffset(globalFSection, ftotpoints[2*p], &globalOff);
7614: DMPlexGetIndicesPointFields_Internal(fsection, PETSC_FALSE, ftotpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, foffsets, PETSC_FALSE, permsF, p, fclperm, findices);
7615: }
7616: for (p = 0; p < numCPoints; p++) {
7617: PetscSectionGetOffset(globalCSection, cpoints[2*p], &globalOff);
7618: DMPlexGetIndicesPointFields_Internal(csection, PETSC_FALSE, cpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, coffsets, PETSC_FALSE, permsC, p, cclperm, cindices);
7619: }
7620: for (f = 0; f < numFields; f++) {
7621: PetscSectionRestoreFieldPointSyms(fsection,f,numFPoints,ftotpoints,&permsF[f],NULL);
7622: PetscSectionRestoreFieldPointSyms(csection,f,numCPoints,cpoints,&permsC[f],NULL);
7623: }
7624: } else {
7625: const PetscInt **permsF = NULL;
7626: const PetscInt **permsC = NULL;
7628: PetscSectionGetPointSyms(fsection,numFPoints,ftotpoints,&permsF,NULL);
7629: PetscSectionGetPointSyms(csection,numCPoints,cpoints,&permsC,NULL);
7630: for (p = 0, off = 0; p < numFPoints; p++) {
7631: const PetscInt *perm = permsF ? permsF[p] : NULL;
7633: PetscSectionGetOffset(globalFSection, ftotpoints[2*p], &globalOff);
7634: DMPlexGetIndicesPoint_Internal(fsection, PETSC_FALSE, ftotpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, &off, PETSC_FALSE, perm, fclperm, findices);
7635: }
7636: for (p = 0, off = 0; p < numCPoints; p++) {
7637: const PetscInt *perm = permsC ? permsC[p] : NULL;
7639: PetscSectionGetOffset(globalCSection, cpoints[2*p], &globalOff);
7640: DMPlexGetIndicesPoint_Internal(csection, PETSC_FALSE, cpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, &off, PETSC_FALSE, perm, cclperm, cindices);
7641: }
7642: PetscSectionRestorePointSyms(fsection,numFPoints,ftotpoints,&permsF,NULL);
7643: PetscSectionRestorePointSyms(csection,numCPoints,cpoints,&permsC,NULL);
7644: }
7645: if (mesh->printSetValues) {DMPlexPrintMatSetValues(PETSC_VIEWER_STDOUT_SELF, A, point, numFIndices, findices, numCIndices, cindices, values);}
7646: /* TODO: flips */
7647: MatSetValues(A, numFIndices, findices, numCIndices, cindices, values, mode);
7648: if (ierr) {
7649: PetscMPIInt rank;
7650: PetscErrorCode ierr2;
7652: ierr2 = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank);CHKERRMPI(ierr2);
7653: ierr2 = (*PetscErrorPrintf)("[%d]ERROR in DMPlexMatSetClosure\n", rank);CHKERRQ(ierr2);
7654: ierr2 = DMPlexPrintMatSetValues(PETSC_VIEWER_STDERR_SELF, A, point, numFIndices, findices, numCIndices, cindices, values);CHKERRQ(ierr2);
7655: ierr2 = DMRestoreWorkArray(dmf, numFIndices, MPIU_INT, &findices);CHKERRQ(ierr2);
7656: ierr2 = DMRestoreWorkArray(dmc, numCIndices, MPIU_INT, &cindices);CHKERRQ(ierr2);
7657:
7658: }
7659: DMRestoreWorkArray(dmf, numCPoints*2*4, MPIU_INT, &ftotpoints);
7660: DMPlexRestoreTransitiveClosure(dmc, point, PETSC_TRUE, &numCPoints, &cpoints);
7661: DMRestoreWorkArray(dmf, numFIndices, MPIU_INT, &findices);
7662: DMRestoreWorkArray(dmc, numCIndices, MPIU_INT, &cindices);
7663: return(0);
7664: }
7666: PetscErrorCode DMPlexMatGetClosureIndicesRefined(DM dmf, PetscSection fsection, PetscSection globalFSection, DM dmc, PetscSection csection, PetscSection globalCSection, PetscInt point, PetscInt cindices[], PetscInt findices[])
7667: {
7668: PetscInt *fpoints = NULL, *ftotpoints = NULL;
7669: PetscInt *cpoints = NULL;
7670: PetscInt foffsets[32], coffsets[32];
7671: const PetscInt *fclperm = NULL, *cclperm = NULL; /* Closure permutations cannot work here */
7672: DMPolytopeType ct;
7673: PetscInt numFields, numSubcells, maxFPoints, numFPoints, numCPoints, numFIndices, numCIndices, dof, off, globalOff, pStart, pEnd, p, q, r, s, f;
7679: if (!fsection) {DMGetLocalSection(dmf, &fsection);}
7681: if (!csection) {DMGetLocalSection(dmc, &csection);}
7683: if (!globalFSection) {DMGetGlobalSection(dmf, &globalFSection);}
7685: if (!globalCSection) {DMGetGlobalSection(dmc, &globalCSection);}
7687: PetscSectionGetNumFields(fsection, &numFields);
7688: if (numFields > 31) SETERRQ1(PetscObjectComm((PetscObject)dmf), PETSC_ERR_ARG_OUTOFRANGE, "Number of fields %D limited to 31", numFields);
7689: PetscArrayzero(foffsets, 32);
7690: PetscArrayzero(coffsets, 32);
7691: /* Column indices */
7692: DMPlexGetTransitiveClosure(dmc, point, PETSC_TRUE, &numCPoints, &cpoints);
7693: maxFPoints = numCPoints;
7694: /* Compress out points not in the section */
7695: /* TODO: Squeeze out points with 0 dof as well */
7696: PetscSectionGetChart(csection, &pStart, &pEnd);
7697: for (p = 0, q = 0; p < numCPoints*2; p += 2) {
7698: if ((cpoints[p] >= pStart) && (cpoints[p] < pEnd)) {
7699: cpoints[q*2] = cpoints[p];
7700: cpoints[q*2+1] = cpoints[p+1];
7701: ++q;
7702: }
7703: }
7704: numCPoints = q;
7705: for (p = 0, numCIndices = 0; p < numCPoints*2; p += 2) {
7706: PetscInt fdof;
7708: PetscSectionGetDof(csection, cpoints[p], &dof);
7709: if (!dof) continue;
7710: for (f = 0; f < numFields; ++f) {
7711: PetscSectionGetFieldDof(csection, cpoints[p], f, &fdof);
7712: coffsets[f+1] += fdof;
7713: }
7714: numCIndices += dof;
7715: }
7716: for (f = 1; f < numFields; ++f) coffsets[f+1] += coffsets[f];
7717: /* Row indices */
7718: DMPlexGetCellType(dmc, point, &ct);
7719: {
7720: DMPlexTransform tr;
7721: DMPolytopeType *rct;
7722: PetscInt *rsize, *rcone, *rornt, Nt;
7724: DMPlexTransformCreate(PETSC_COMM_SELF, &tr);
7725: DMPlexTransformSetType(tr, DMPLEXREFINEREGULAR);
7726: DMPlexTransformCellTransform(tr, ct, point, NULL, &Nt, &rct, &rsize, &rcone, &rornt);
7727: numSubcells = rsize[Nt-1];
7728: DMPlexTransformDestroy(&tr);
7729: }
7730: DMGetWorkArray(dmf, maxFPoints*2*numSubcells, MPIU_INT, &ftotpoints);
7731: for (r = 0, q = 0; r < numSubcells; ++r) {
7732: /* TODO Map from coarse to fine cells */
7733: DMPlexGetTransitiveClosure(dmf, point*numSubcells + r, PETSC_TRUE, &numFPoints, &fpoints);
7734: /* Compress out points not in the section */
7735: PetscSectionGetChart(fsection, &pStart, &pEnd);
7736: for (p = 0; p < numFPoints*2; p += 2) {
7737: if ((fpoints[p] >= pStart) && (fpoints[p] < pEnd)) {
7738: PetscSectionGetDof(fsection, fpoints[p], &dof);
7739: if (!dof) continue;
7740: for (s = 0; s < q; ++s) if (fpoints[p] == ftotpoints[s*2]) break;
7741: if (s < q) continue;
7742: ftotpoints[q*2] = fpoints[p];
7743: ftotpoints[q*2+1] = fpoints[p+1];
7744: ++q;
7745: }
7746: }
7747: DMPlexRestoreTransitiveClosure(dmf, point, PETSC_TRUE, &numFPoints, &fpoints);
7748: }
7749: numFPoints = q;
7750: for (p = 0, numFIndices = 0; p < numFPoints*2; p += 2) {
7751: PetscInt fdof;
7753: PetscSectionGetDof(fsection, ftotpoints[p], &dof);
7754: if (!dof) continue;
7755: for (f = 0; f < numFields; ++f) {
7756: PetscSectionGetFieldDof(fsection, ftotpoints[p], f, &fdof);
7757: foffsets[f+1] += fdof;
7758: }
7759: numFIndices += dof;
7760: }
7761: for (f = 1; f < numFields; ++f) foffsets[f+1] += foffsets[f];
7763: if (numFields && foffsets[numFields] != numFIndices) SETERRQ2(PetscObjectComm((PetscObject)dmf), PETSC_ERR_PLIB, "Invalid size for closure %D should be %D", foffsets[numFields], numFIndices);
7764: if (numFields && coffsets[numFields] != numCIndices) SETERRQ2(PetscObjectComm((PetscObject)dmc), PETSC_ERR_PLIB, "Invalid size for closure %D should be %D", coffsets[numFields], numCIndices);
7765: if (numFields) {
7766: const PetscInt **permsF[32] = {NULL};
7767: const PetscInt **permsC[32] = {NULL};
7769: for (f = 0; f < numFields; f++) {
7770: PetscSectionGetFieldPointSyms(fsection,f,numFPoints,ftotpoints,&permsF[f],NULL);
7771: PetscSectionGetFieldPointSyms(csection,f,numCPoints,cpoints,&permsC[f],NULL);
7772: }
7773: for (p = 0; p < numFPoints; p++) {
7774: PetscSectionGetOffset(globalFSection, ftotpoints[2*p], &globalOff);
7775: DMPlexGetIndicesPointFields_Internal(fsection, PETSC_FALSE, ftotpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, foffsets, PETSC_FALSE, permsF, p, fclperm, findices);
7776: }
7777: for (p = 0; p < numCPoints; p++) {
7778: PetscSectionGetOffset(globalCSection, cpoints[2*p], &globalOff);
7779: DMPlexGetIndicesPointFields_Internal(csection, PETSC_FALSE, cpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, coffsets, PETSC_FALSE, permsC, p, cclperm, cindices);
7780: }
7781: for (f = 0; f < numFields; f++) {
7782: PetscSectionRestoreFieldPointSyms(fsection,f,numFPoints,ftotpoints,&permsF[f],NULL);
7783: PetscSectionRestoreFieldPointSyms(csection,f,numCPoints,cpoints,&permsC[f],NULL);
7784: }
7785: } else {
7786: const PetscInt **permsF = NULL;
7787: const PetscInt **permsC = NULL;
7789: PetscSectionGetPointSyms(fsection,numFPoints,ftotpoints,&permsF,NULL);
7790: PetscSectionGetPointSyms(csection,numCPoints,cpoints,&permsC,NULL);
7791: for (p = 0, off = 0; p < numFPoints; p++) {
7792: const PetscInt *perm = permsF ? permsF[p] : NULL;
7794: PetscSectionGetOffset(globalFSection, ftotpoints[2*p], &globalOff);
7795: DMPlexGetIndicesPoint_Internal(fsection, PETSC_FALSE, ftotpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, &off, PETSC_FALSE, perm, fclperm, findices);
7796: }
7797: for (p = 0, off = 0; p < numCPoints; p++) {
7798: const PetscInt *perm = permsC ? permsC[p] : NULL;
7800: PetscSectionGetOffset(globalCSection, cpoints[2*p], &globalOff);
7801: DMPlexGetIndicesPoint_Internal(csection, PETSC_FALSE, cpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, &off, PETSC_FALSE, perm, cclperm, cindices);
7802: }
7803: PetscSectionRestorePointSyms(fsection,numFPoints,ftotpoints,&permsF,NULL);
7804: PetscSectionRestorePointSyms(csection,numCPoints,cpoints,&permsC,NULL);
7805: }
7806: DMRestoreWorkArray(dmf, numCPoints*2*4, MPIU_INT, &ftotpoints);
7807: DMPlexRestoreTransitiveClosure(dmc, point, PETSC_TRUE, &numCPoints, &cpoints);
7808: return(0);
7809: }
7811: /*@C
7812: DMPlexGetVTKCellHeight - Returns the height in the DAG used to determine which points are cells (normally 0)
7814: Input Parameter:
7815: . dm - The DMPlex object
7817: Output Parameter:
7818: . cellHeight - The height of a cell
7820: Level: developer
7822: .seealso DMPlexSetVTKCellHeight()
7823: @*/
7824: PetscErrorCode DMPlexGetVTKCellHeight(DM dm, PetscInt *cellHeight)
7825: {
7826: DM_Plex *mesh = (DM_Plex*) dm->data;
7831: *cellHeight = mesh->vtkCellHeight;
7832: return(0);
7833: }
7835: /*@C
7836: DMPlexSetVTKCellHeight - Sets the height in the DAG used to determine which points are cells (normally 0)
7838: Input Parameters:
7839: + dm - The DMPlex object
7840: - cellHeight - The height of a cell
7842: Level: developer
7844: .seealso DMPlexGetVTKCellHeight()
7845: @*/
7846: PetscErrorCode DMPlexSetVTKCellHeight(DM dm, PetscInt cellHeight)
7847: {
7848: DM_Plex *mesh = (DM_Plex*) dm->data;
7852: mesh->vtkCellHeight = cellHeight;
7853: return(0);
7854: }
7856: /*@
7857: DMPlexGetGhostCellStratum - Get the range of cells which are used to enforce FV boundary conditions
7859: Input Parameter:
7860: . dm - The DMPlex object
7862: Output Parameters:
7863: + gcStart - The first ghost cell, or NULL
7864: - gcEnd - The upper bound on ghost cells, or NULL
7866: Level: advanced
7868: .seealso DMPlexConstructGhostCells(), DMPlexGetGhostCellStratum()
7869: @*/
7870: PetscErrorCode DMPlexGetGhostCellStratum(DM dm, PetscInt *gcStart, PetscInt *gcEnd)
7871: {
7872: DMLabel ctLabel;
7877: DMPlexGetCellTypeLabel(dm, &ctLabel);
7878: DMLabelGetStratumBounds(ctLabel, DM_POLYTOPE_FV_GHOST, gcStart, gcEnd);
7879: return(0);
7880: }
7882: /* We can easily have a form that takes an IS instead */
7883: PetscErrorCode DMPlexCreateNumbering_Plex(DM dm, PetscInt pStart, PetscInt pEnd, PetscInt shift, PetscInt *globalSize, PetscSF sf, IS *numbering)
7884: {
7885: PetscSection section, globalSection;
7886: PetscInt *numbers, p;
7890: PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion);
7891: PetscSectionSetChart(section, pStart, pEnd);
7892: for (p = pStart; p < pEnd; ++p) {
7893: PetscSectionSetDof(section, p, 1);
7894: }
7895: PetscSectionSetUp(section);
7896: PetscSectionCreateGlobalSection(section, sf, PETSC_FALSE, PETSC_FALSE, &globalSection);
7897: PetscMalloc1(pEnd - pStart, &numbers);
7898: for (p = pStart; p < pEnd; ++p) {
7899: PetscSectionGetOffset(globalSection, p, &numbers[p-pStart]);
7900: if (numbers[p-pStart] < 0) numbers[p-pStart] -= shift;
7901: else numbers[p-pStart] += shift;
7902: }
7903: ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd - pStart, numbers, PETSC_OWN_POINTER, numbering);
7904: if (globalSize) {
7905: PetscLayout layout;
7906: PetscSectionGetPointLayout(PetscObjectComm((PetscObject) dm), globalSection, &layout);
7907: PetscLayoutGetSize(layout, globalSize);
7908: PetscLayoutDestroy(&layout);
7909: }
7910: PetscSectionDestroy(§ion);
7911: PetscSectionDestroy(&globalSection);
7912: return(0);
7913: }
7915: PetscErrorCode DMPlexCreateCellNumbering_Internal(DM dm, PetscBool includeHybrid, IS *globalCellNumbers)
7916: {
7917: PetscInt cellHeight, cStart, cEnd;
7921: DMPlexGetVTKCellHeight(dm, &cellHeight);
7922: if (includeHybrid) {DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);}
7923: else {DMPlexGetSimplexOrBoxCells(dm, cellHeight, &cStart, &cEnd);}
7924: DMPlexCreateNumbering_Plex(dm, cStart, cEnd, 0, NULL, dm->sf, globalCellNumbers);
7925: return(0);
7926: }
7928: /*@
7929: DMPlexGetCellNumbering - Get a global cell numbering for all cells on this process
7931: Input Parameter:
7932: . dm - The DMPlex object
7934: Output Parameter:
7935: . globalCellNumbers - Global cell numbers for all cells on this process
7937: Level: developer
7939: .seealso DMPlexGetVertexNumbering()
7940: @*/
7941: PetscErrorCode DMPlexGetCellNumbering(DM dm, IS *globalCellNumbers)
7942: {
7943: DM_Plex *mesh = (DM_Plex*) dm->data;
7948: if (!mesh->globalCellNumbers) {DMPlexCreateCellNumbering_Internal(dm, PETSC_FALSE, &mesh->globalCellNumbers);}
7949: *globalCellNumbers = mesh->globalCellNumbers;
7950: return(0);
7951: }
7953: PetscErrorCode DMPlexCreateVertexNumbering_Internal(DM dm, PetscBool includeHybrid, IS *globalVertexNumbers)
7954: {
7955: PetscInt vStart, vEnd;
7960: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
7961: DMPlexCreateNumbering_Plex(dm, vStart, vEnd, 0, NULL, dm->sf, globalVertexNumbers);
7962: return(0);
7963: }
7965: /*@
7966: DMPlexGetVertexNumbering - Get a global vertex numbering for all vertices on this process
7968: Input Parameter:
7969: . dm - The DMPlex object
7971: Output Parameter:
7972: . globalVertexNumbers - Global vertex numbers for all vertices on this process
7974: Level: developer
7976: .seealso DMPlexGetCellNumbering()
7977: @*/
7978: PetscErrorCode DMPlexGetVertexNumbering(DM dm, IS *globalVertexNumbers)
7979: {
7980: DM_Plex *mesh = (DM_Plex*) dm->data;
7985: if (!mesh->globalVertexNumbers) {DMPlexCreateVertexNumbering_Internal(dm, PETSC_FALSE, &mesh->globalVertexNumbers);}
7986: *globalVertexNumbers = mesh->globalVertexNumbers;
7987: return(0);
7988: }
7990: /*@
7991: DMPlexCreatePointNumbering - Create a global numbering for all points on this process
7993: Input Parameter:
7994: . dm - The DMPlex object
7996: Output Parameter:
7997: . globalPointNumbers - Global numbers for all points on this process
7999: Level: developer
8001: .seealso DMPlexGetCellNumbering()
8002: @*/
8003: PetscErrorCode DMPlexCreatePointNumbering(DM dm, IS *globalPointNumbers)
8004: {
8005: IS nums[4];
8006: PetscInt depths[4], gdepths[4], starts[4];
8007: PetscInt depth, d, shift = 0;
8012: DMPlexGetDepth(dm, &depth);
8013: /* For unstratified meshes use dim instead of depth */
8014: if (depth < 0) {DMGetDimension(dm, &depth);}
8015: for (d = 0; d <= depth; ++d) {
8016: PetscInt end;
8018: depths[d] = depth-d;
8019: DMPlexGetDepthStratum(dm, depths[d], &starts[d], &end);
8020: if (!(starts[d]-end)) { starts[d] = depths[d] = -1; }
8021: }
8022: PetscSortIntWithArray(depth+1, starts, depths);
8023: MPIU_Allreduce(depths, gdepths, depth+1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject) dm));
8024: for (d = 0; d <= depth; ++d) {
8025: if (starts[d] >= 0 && depths[d] != gdepths[d]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected depth %D, found %D",depths[d],gdepths[d]);
8026: }
8027: for (d = 0; d <= depth; ++d) {
8028: PetscInt pStart, pEnd, gsize;
8030: DMPlexGetDepthStratum(dm, gdepths[d], &pStart, &pEnd);
8031: DMPlexCreateNumbering_Plex(dm, pStart, pEnd, shift, &gsize, dm->sf, &nums[d]);
8032: shift += gsize;
8033: }
8034: ISConcatenate(PetscObjectComm((PetscObject) dm), depth+1, nums, globalPointNumbers);
8035: for (d = 0; d <= depth; ++d) {ISDestroy(&nums[d]);}
8036: return(0);
8037: }
8039: /*@
8040: DMPlexCreateRankField - Create a cell field whose value is the rank of the owner
8042: Input Parameter:
8043: . dm - The DMPlex object
8045: Output Parameter:
8046: . ranks - The rank field
8048: Options Database Keys:
8049: . -dm_partition_view - Adds the rank field into the DM output from -dm_view using the same viewer
8051: Level: intermediate
8053: .seealso: DMView()
8054: @*/
8055: PetscErrorCode DMPlexCreateRankField(DM dm, Vec *ranks)
8056: {
8057: DM rdm;
8058: PetscFE fe;
8059: PetscScalar *r;
8060: PetscMPIInt rank;
8061: DMPolytopeType ct;
8062: PetscInt dim, cStart, cEnd, c;
8063: PetscBool simplex;
8069: MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
8070: DMClone(dm, &rdm);
8071: DMGetDimension(rdm, &dim);
8072: DMPlexGetHeightStratum(rdm, 0, &cStart, &cEnd);
8073: DMPlexGetCellType(dm, cStart, &ct);
8074: simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
8075: PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "PETSc___rank_", -1, &fe);
8076: PetscObjectSetName((PetscObject) fe, "rank");
8077: DMSetField(rdm, 0, NULL, (PetscObject) fe);
8078: PetscFEDestroy(&fe);
8079: DMCreateDS(rdm);
8080: DMCreateGlobalVector(rdm, ranks);
8081: PetscObjectSetName((PetscObject) *ranks, "partition");
8082: VecGetArray(*ranks, &r);
8083: for (c = cStart; c < cEnd; ++c) {
8084: PetscScalar *lr;
8086: DMPlexPointGlobalRef(rdm, c, r, &lr);
8087: if (lr) *lr = rank;
8088: }
8089: VecRestoreArray(*ranks, &r);
8090: DMDestroy(&rdm);
8091: return(0);
8092: }
8094: /*@
8095: DMPlexCreateLabelField - Create a cell field whose value is the label value for that cell
8097: Input Parameters:
8098: + dm - The DMPlex
8099: - label - The DMLabel
8101: Output Parameter:
8102: . val - The label value field
8104: Options Database Keys:
8105: . -dm_label_view - Adds the label value field into the DM output from -dm_view using the same viewer
8107: Level: intermediate
8109: .seealso: DMView()
8110: @*/
8111: PetscErrorCode DMPlexCreateLabelField(DM dm, DMLabel label, Vec *val)
8112: {
8113: DM rdm;
8114: PetscFE fe;
8115: PetscScalar *v;
8116: PetscInt dim, cStart, cEnd, c;
8123: DMClone(dm, &rdm);
8124: DMGetDimension(rdm, &dim);
8125: PetscFECreateDefault(PetscObjectComm((PetscObject) rdm), dim, 1, PETSC_TRUE, "PETSc___label_value_", -1, &fe);
8126: PetscObjectSetName((PetscObject) fe, "label_value");
8127: DMSetField(rdm, 0, NULL, (PetscObject) fe);
8128: PetscFEDestroy(&fe);
8129: DMCreateDS(rdm);
8130: DMPlexGetHeightStratum(rdm, 0, &cStart, &cEnd);
8131: DMCreateGlobalVector(rdm, val);
8132: PetscObjectSetName((PetscObject) *val, "label_value");
8133: VecGetArray(*val, &v);
8134: for (c = cStart; c < cEnd; ++c) {
8135: PetscScalar *lv;
8136: PetscInt cval;
8138: DMPlexPointGlobalRef(rdm, c, v, &lv);
8139: DMLabelGetValue(label, c, &cval);
8140: *lv = cval;
8141: }
8142: VecRestoreArray(*val, &v);
8143: DMDestroy(&rdm);
8144: return(0);
8145: }
8147: /*@
8148: DMPlexCheckSymmetry - Check that the adjacency information in the mesh is symmetric.
8150: Input Parameter:
8151: . dm - The DMPlex object
8153: Notes:
8154: This is a useful diagnostic when creating meshes programmatically.
8156: For the complete list of DMPlexCheck* functions, see DMSetFromOptions().
8158: Level: developer
8160: .seealso: DMCreate(), DMSetFromOptions()
8161: @*/
8162: PetscErrorCode DMPlexCheckSymmetry(DM dm)
8163: {
8164: PetscSection coneSection, supportSection;
8165: const PetscInt *cone, *support;
8166: PetscInt coneSize, c, supportSize, s;
8167: PetscInt pStart, pEnd, p, pp, csize, ssize;
8168: PetscBool storagecheck = PETSC_TRUE;
8169: PetscErrorCode ierr;
8173: DMViewFromOptions(dm, NULL, "-sym_dm_view");
8174: DMPlexGetConeSection(dm, &coneSection);
8175: DMPlexGetSupportSection(dm, &supportSection);
8176: /* Check that point p is found in the support of its cone points, and vice versa */
8177: DMPlexGetChart(dm, &pStart, &pEnd);
8178: for (p = pStart; p < pEnd; ++p) {
8179: DMPlexGetConeSize(dm, p, &coneSize);
8180: DMPlexGetCone(dm, p, &cone);
8181: for (c = 0; c < coneSize; ++c) {
8182: PetscBool dup = PETSC_FALSE;
8183: PetscInt d;
8184: for (d = c-1; d >= 0; --d) {
8185: if (cone[c] == cone[d]) {dup = PETSC_TRUE; break;}
8186: }
8187: DMPlexGetSupportSize(dm, cone[c], &supportSize);
8188: DMPlexGetSupport(dm, cone[c], &support);
8189: for (s = 0; s < supportSize; ++s) {
8190: if (support[s] == p) break;
8191: }
8192: if ((s >= supportSize) || (dup && (support[s+1] != p))) {
8193: PetscPrintf(PETSC_COMM_SELF, "p: %D cone: ", p);
8194: for (s = 0; s < coneSize; ++s) {
8195: PetscPrintf(PETSC_COMM_SELF, "%D, ", cone[s]);
8196: }
8197: PetscPrintf(PETSC_COMM_SELF, "\n");
8198: PetscPrintf(PETSC_COMM_SELF, "p: %D support: ", cone[c]);
8199: for (s = 0; s < supportSize; ++s) {
8200: PetscPrintf(PETSC_COMM_SELF, "%D, ", support[s]);
8201: }
8202: PetscPrintf(PETSC_COMM_SELF, "\n");
8203: if (dup) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D not repeatedly found in support of repeated cone point %D", p, cone[c]);
8204: else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D not found in support of cone point %D", p, cone[c]);
8205: }
8206: }
8207: DMPlexGetTreeParent(dm, p, &pp, NULL);
8208: if (p != pp) { storagecheck = PETSC_FALSE; continue; }
8209: DMPlexGetSupportSize(dm, p, &supportSize);
8210: DMPlexGetSupport(dm, p, &support);
8211: for (s = 0; s < supportSize; ++s) {
8212: DMPlexGetConeSize(dm, support[s], &coneSize);
8213: DMPlexGetCone(dm, support[s], &cone);
8214: for (c = 0; c < coneSize; ++c) {
8215: DMPlexGetTreeParent(dm, cone[c], &pp, NULL);
8216: if (cone[c] != pp) { c = 0; break; }
8217: if (cone[c] == p) break;
8218: }
8219: if (c >= coneSize) {
8220: PetscPrintf(PETSC_COMM_SELF, "p: %D support: ", p);
8221: for (c = 0; c < supportSize; ++c) {
8222: PetscPrintf(PETSC_COMM_SELF, "%D, ", support[c]);
8223: }
8224: PetscPrintf(PETSC_COMM_SELF, "\n");
8225: PetscPrintf(PETSC_COMM_SELF, "p: %D cone: ", support[s]);
8226: for (c = 0; c < coneSize; ++c) {
8227: PetscPrintf(PETSC_COMM_SELF, "%D, ", cone[c]);
8228: }
8229: PetscPrintf(PETSC_COMM_SELF, "\n");
8230: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D not found in cone of support point %D", p, support[s]);
8231: }
8232: }
8233: }
8234: if (storagecheck) {
8235: PetscSectionGetStorageSize(coneSection, &csize);
8236: PetscSectionGetStorageSize(supportSection, &ssize);
8237: if (csize != ssize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Total cone size %D != Total support size %D", csize, ssize);
8238: }
8239: return(0);
8240: }
8242: /*
8243: For submeshes with cohesive cells (see DMPlexConstructCohesiveCells()), we allow a special case where some of the boundary of a face (edges and vertices) are not duplicated. We call these special boundary points "unsplit", since the same edge or vertex appears in both copies of the face. These unsplit points throw off our counting, so we have to explicitly account for them here.
8244: */
8245: static PetscErrorCode DMPlexCellUnsplitVertices_Private(DM dm, PetscInt c, DMPolytopeType ct, PetscInt *unsplit)
8246: {
8247: DMPolytopeType cct;
8248: PetscInt ptpoints[4];
8249: const PetscInt *cone, *ccone, *ptcone;
8250: PetscInt coneSize, cp, cconeSize, ccp, npt = 0, pt;
8251: PetscErrorCode ierr;
8254: *unsplit = 0;
8255: switch (ct) {
8256: case DM_POLYTOPE_POINT_PRISM_TENSOR:
8257: ptpoints[npt++] = c;
8258: break;
8259: case DM_POLYTOPE_SEG_PRISM_TENSOR:
8260: DMPlexGetCone(dm, c, &cone);
8261: DMPlexGetConeSize(dm, c, &coneSize);
8262: for (cp = 0; cp < coneSize; ++cp) {
8263: DMPlexGetCellType(dm, cone[cp], &cct);
8264: if (cct == DM_POLYTOPE_POINT_PRISM_TENSOR) ptpoints[npt++] = cone[cp];
8265: }
8266: break;
8267: case DM_POLYTOPE_TRI_PRISM_TENSOR:
8268: case DM_POLYTOPE_QUAD_PRISM_TENSOR:
8269: DMPlexGetCone(dm, c, &cone);
8270: DMPlexGetConeSize(dm, c, &coneSize);
8271: for (cp = 0; cp < coneSize; ++cp) {
8272: DMPlexGetCone(dm, cone[cp], &ccone);
8273: DMPlexGetConeSize(dm, cone[cp], &cconeSize);
8274: for (ccp = 0; ccp < cconeSize; ++ccp) {
8275: DMPlexGetCellType(dm, ccone[ccp], &cct);
8276: if (cct == DM_POLYTOPE_POINT_PRISM_TENSOR) {
8277: PetscInt p;
8278: for (p = 0; p < npt; ++p) if (ptpoints[p] == ccone[ccp]) break;
8279: if (p == npt) ptpoints[npt++] = ccone[ccp];
8280: }
8281: }
8282: }
8283: break;
8284: default: break;
8285: }
8286: for (pt = 0; pt < npt; ++pt) {
8287: DMPlexGetCone(dm, ptpoints[pt], &ptcone);
8288: if (ptcone[0] == ptcone[1]) ++(*unsplit);
8289: }
8290: return(0);
8291: }
8293: /*@
8294: DMPlexCheckSkeleton - Check that each cell has the correct number of vertices
8296: Input Parameters:
8297: + dm - The DMPlex object
8298: - cellHeight - Normally 0
8300: Notes:
8301: This is a useful diagnostic when creating meshes programmatically.
8302: Currently applicable only to homogeneous simplex or tensor meshes.
8304: For the complete list of DMPlexCheck* functions, see DMSetFromOptions().
8306: Level: developer
8308: .seealso: DMCreate(), DMSetFromOptions()
8309: @*/
8310: PetscErrorCode DMPlexCheckSkeleton(DM dm, PetscInt cellHeight)
8311: {
8312: DMPlexInterpolatedFlag interp;
8313: DMPolytopeType ct;
8314: PetscInt vStart, vEnd, cStart, cEnd, c;
8315: PetscErrorCode ierr;
8319: DMPlexIsInterpolated(dm, &interp);
8320: DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);
8321: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
8322: for (c = cStart; c < cEnd; ++c) {
8323: PetscInt *closure = NULL;
8324: PetscInt coneSize, closureSize, cl, Nv = 0;
8326: DMPlexGetCellType(dm, c, &ct);
8327: if ((PetscInt) ct < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %D has no cell type", c);
8328: if (ct == DM_POLYTOPE_UNKNOWN) continue;
8329: if (interp == DMPLEX_INTERPOLATED_FULL) {
8330: DMPlexGetConeSize(dm, c, &coneSize);
8331: if (coneSize != DMPolytopeTypeGetConeSize(ct)) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %D of type %s has cone size %D != %D", c, DMPolytopeTypes[ct], coneSize, DMPolytopeTypeGetConeSize(ct));
8332: }
8333: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
8334: for (cl = 0; cl < closureSize*2; cl += 2) {
8335: const PetscInt p = closure[cl];
8336: if ((p >= vStart) && (p < vEnd)) ++Nv;
8337: }
8338: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
8339: /* Special Case: Tensor faces with identified vertices */
8340: if (Nv < DMPolytopeTypeGetNumVertices(ct)) {
8341: PetscInt unsplit;
8343: DMPlexCellUnsplitVertices_Private(dm, c, ct, &unsplit);
8344: if (Nv + unsplit == DMPolytopeTypeGetNumVertices(ct)) continue;
8345: }
8346: if (Nv != DMPolytopeTypeGetNumVertices(ct)) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %D of type %s has %D vertices != %D", c, DMPolytopeTypes[ct], Nv, DMPolytopeTypeGetNumVertices(ct));
8347: }
8348: return(0);
8349: }
8351: /*@
8352: DMPlexCheckFaces - Check that the faces of each cell give a vertex order this is consistent with what we expect from the cell type
8354: Not Collective
8356: Input Parameters:
8357: + dm - The DMPlex object
8358: - cellHeight - Normally 0
8360: Notes:
8361: This is a useful diagnostic when creating meshes programmatically.
8362: This routine is only relevant for meshes that are fully interpolated across all ranks.
8363: It will error out if a partially interpolated mesh is given on some rank.
8364: It will do nothing for locally uninterpolated mesh (as there is nothing to check).
8366: For the complete list of DMPlexCheck* functions, see DMSetFromOptions().
8368: Level: developer
8370: .seealso: DMCreate(), DMPlexGetVTKCellHeight(), DMSetFromOptions()
8371: @*/
8372: PetscErrorCode DMPlexCheckFaces(DM dm, PetscInt cellHeight)
8373: {
8374: PetscInt dim, depth, vStart, vEnd, cStart, cEnd, c, h;
8376: DMPlexInterpolatedFlag interpEnum;
8380: DMPlexIsInterpolated(dm, &interpEnum);
8381: if (interpEnum == DMPLEX_INTERPOLATED_NONE) return(0);
8382: if (interpEnum == DMPLEX_INTERPOLATED_PARTIAL) {
8383: PetscMPIInt rank;
8384: MPI_Comm comm;
8386: PetscObjectGetComm((PetscObject) dm, &comm);
8387: MPI_Comm_rank(comm, &rank);
8388: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Mesh is only partially interpolated on rank %d, this is currently not supported", rank);
8389: }
8391: DMGetDimension(dm, &dim);
8392: DMPlexGetDepth(dm, &depth);
8393: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
8394: for (h = cellHeight; h < PetscMin(depth, dim); ++h) {
8395: DMPlexGetHeightStratum(dm, h, &cStart, &cEnd);
8396: for (c = cStart; c < cEnd; ++c) {
8397: const PetscInt *cone, *ornt, *faceSizes, *faces;
8398: const DMPolytopeType *faceTypes;
8399: DMPolytopeType ct;
8400: PetscInt numFaces, coneSize, f;
8401: PetscInt *closure = NULL, closureSize, cl, numCorners = 0, fOff = 0, unsplit;
8403: DMPlexGetCellType(dm, c, &ct);
8404: DMPlexCellUnsplitVertices_Private(dm, c, ct, &unsplit);
8405: if (unsplit) continue;
8406: DMPlexGetConeSize(dm, c, &coneSize);
8407: DMPlexGetCone(dm, c, &cone);
8408: DMPlexGetConeOrientation(dm, c, &ornt);
8409: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
8410: for (cl = 0; cl < closureSize*2; cl += 2) {
8411: const PetscInt p = closure[cl];
8412: if ((p >= vStart) && (p < vEnd)) closure[numCorners++] = p;
8413: }
8414: DMPlexGetRawFaces_Internal(dm, ct, closure, &numFaces, &faceTypes, &faceSizes, &faces);
8415: if (coneSize != numFaces) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %D of type %s has %D faces but should have %D", c, DMPolytopeTypes[ct], coneSize, numFaces);
8416: for (f = 0; f < numFaces; ++f) {
8417: DMPolytopeType fct;
8418: PetscInt *fclosure = NULL, fclosureSize, cl, fnumCorners = 0, v;
8420: DMPlexGetCellType(dm, cone[f], &fct);
8421: DMPlexGetTransitiveClosure_Internal(dm, cone[f], ornt[f], PETSC_TRUE, &fclosureSize, &fclosure);
8422: for (cl = 0; cl < fclosureSize*2; cl += 2) {
8423: const PetscInt p = fclosure[cl];
8424: if ((p >= vStart) && (p < vEnd)) fclosure[fnumCorners++] = p;
8425: }
8426: if (fnumCorners != faceSizes[f]) SETERRQ7(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %D of type %s (cone idx %D) of cell %D of type %s has %D vertices but should have %D", cone[f], DMPolytopeTypes[fct], f, c, DMPolytopeTypes[ct], fnumCorners, faceSizes[f]);
8427: for (v = 0; v < fnumCorners; ++v) {
8428: if (fclosure[v] != faces[fOff+v]) {
8429: PetscInt v1;
8431: PetscPrintf(PETSC_COMM_SELF, "face closure:");
8432: for (v1 = 0; v1 < fnumCorners; ++v1) {PetscPrintf(PETSC_COMM_SELF, " %D", fclosure[v1]);}
8433: PetscPrintf(PETSC_COMM_SELF, "\ncell face:");
8434: for (v1 = 0; v1 < fnumCorners; ++v1) {PetscPrintf(PETSC_COMM_SELF, " %D", faces[fOff+v1]);}
8435: PetscPrintf(PETSC_COMM_SELF, "\n");
8436: SETERRQ9(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %D of type %s (cone idx %d, ornt %D) of cell %D of type %s vertex %D, %D != %D", cone[f], DMPolytopeTypes[fct], f, ornt[f], c, DMPolytopeTypes[ct], v, fclosure[v], faces[fOff+v]);
8437: }
8438: }
8439: DMPlexRestoreTransitiveClosure(dm, cone[f], PETSC_TRUE, &fclosureSize, &fclosure);
8440: fOff += faceSizes[f];
8441: }
8442: DMPlexRestoreRawFaces_Internal(dm, ct, closure, &numFaces, &faceTypes, &faceSizes, &faces);
8443: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
8444: }
8445: }
8446: return(0);
8447: }
8449: /*@
8450: DMPlexCheckGeometry - Check the geometry of mesh cells
8452: Input Parameter:
8453: . dm - The DMPlex object
8455: Notes:
8456: This is a useful diagnostic when creating meshes programmatically.
8458: For the complete list of DMPlexCheck* functions, see DMSetFromOptions().
8460: Level: developer
8462: .seealso: DMCreate(), DMSetFromOptions()
8463: @*/
8464: PetscErrorCode DMPlexCheckGeometry(DM dm)
8465: {
8466: Vec coordinates;
8467: PetscReal detJ, J[9], refVol = 1.0;
8468: PetscReal vol;
8469: PetscBool periodic;
8470: PetscInt dim, depth, dE, d, cStart, cEnd, c;
8474: DMGetDimension(dm, &dim);
8475: DMGetCoordinateDim(dm, &dE);
8476: if (dim != dE) return(0);
8477: DMPlexGetDepth(dm, &depth);
8478: DMGetPeriodicity(dm, &periodic, NULL, NULL, NULL);
8479: for (d = 0; d < dim; ++d) refVol *= 2.0;
8480: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
8481: /* Make sure local coordinates are created, because that step is collective */
8482: DMGetCoordinatesLocal(dm, &coordinates);
8483: for (c = cStart; c < cEnd; ++c) {
8484: DMPolytopeType ct;
8485: PetscInt unsplit;
8486: PetscBool ignoreZeroVol = PETSC_FALSE;
8488: DMPlexGetCellType(dm, c, &ct);
8489: switch (ct) {
8490: case DM_POLYTOPE_SEG_PRISM_TENSOR:
8491: case DM_POLYTOPE_TRI_PRISM_TENSOR:
8492: case DM_POLYTOPE_QUAD_PRISM_TENSOR:
8493: ignoreZeroVol = PETSC_TRUE; break;
8494: default: break;
8495: }
8496: switch (ct) {
8497: case DM_POLYTOPE_TRI_PRISM:
8498: case DM_POLYTOPE_TRI_PRISM_TENSOR:
8499: case DM_POLYTOPE_QUAD_PRISM_TENSOR:
8500: case DM_POLYTOPE_PYRAMID:
8501: continue;
8502: default: break;
8503: }
8504: DMPlexCellUnsplitVertices_Private(dm, c, ct, &unsplit);
8505: if (unsplit) continue;
8506: DMPlexComputeCellGeometryFEM(dm, c, NULL, NULL, J, NULL, &detJ);
8507: if (detJ < -PETSC_SMALL || (detJ <= 0.0 && !ignoreZeroVol)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh cell %D of type %s is inverted, |J| = %g", c, DMPolytopeTypes[ct], (double) detJ);
8508: PetscInfo2(dm, "Cell %D FEM Volume %g\n", c, (double) detJ*refVol);
8509: if (depth > 1 && !periodic) {
8510: DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);
8511: if (vol < -PETSC_SMALL || (vol <= 0.0 && !ignoreZeroVol)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh cell %D of type %s is inverted, vol = %g", c, DMPolytopeTypes[ct], (double) vol);
8512: PetscInfo2(dm, "Cell %D FVM Volume %g\n", c, (double) vol);
8513: }
8514: }
8515: return(0);
8516: }
8518: /*@
8519: DMPlexCheckPointSF - Check that several necessary conditions are met for the point SF of this plex.
8521: Input Parameters:
8522: . dm - The DMPlex object
8524: Notes:
8525: This is mainly intended for debugging/testing purposes.
8526: It currently checks only meshes with no partition overlapping.
8528: For the complete list of DMPlexCheck* functions, see DMSetFromOptions().
8530: Level: developer
8532: .seealso: DMGetPointSF(), DMSetFromOptions()
8533: @*/
8534: PetscErrorCode DMPlexCheckPointSF(DM dm)
8535: {
8536: PetscSF pointSF;
8537: PetscInt cellHeight, cStart, cEnd, l, nleaves, nroots, overlap;
8538: const PetscInt *locals, *rootdegree;
8539: PetscBool distributed;
8540: PetscErrorCode ierr;
8544: DMGetPointSF(dm, &pointSF);
8545: DMPlexIsDistributed(dm, &distributed);
8546: if (!distributed) return(0);
8547: DMPlexGetOverlap(dm, &overlap);
8548: if (overlap) {
8549: PetscPrintf(PetscObjectComm((PetscObject)dm), "Warning: DMPlexCheckPointSF() is currently not implemented for meshes with partition overlapping");
8550: return(0);
8551: }
8552: if (!pointSF) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This DMPlex is distributed but does not have PointSF attached");
8553: PetscSFGetGraph(pointSF, &nroots, &nleaves, &locals, NULL);
8554: if (nroots < 0) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This DMPlex is distributed but its PointSF has no graph set");
8555: PetscSFComputeDegreeBegin(pointSF, &rootdegree);
8556: PetscSFComputeDegreeEnd(pointSF, &rootdegree);
8558: /* 1) check there are no faces in 2D, cells in 3D, in interface */
8559: DMPlexGetVTKCellHeight(dm, &cellHeight);
8560: DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);
8561: for (l = 0; l < nleaves; ++l) {
8562: const PetscInt point = locals[l];
8564: if (point >= cStart && point < cEnd) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point SF contains %D which is a cell", point);
8565: }
8567: /* 2) if some point is in interface, then all its cone points must be also in interface (either as leaves or roots) */
8568: for (l = 0; l < nleaves; ++l) {
8569: const PetscInt point = locals[l];
8570: const PetscInt *cone;
8571: PetscInt coneSize, c, idx;
8573: DMPlexGetConeSize(dm, point, &coneSize);
8574: DMPlexGetCone(dm, point, &cone);
8575: for (c = 0; c < coneSize; ++c) {
8576: if (!rootdegree[cone[c]]) {
8577: PetscFindInt(cone[c], nleaves, locals, &idx);
8578: if (idx < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point SF contains %D but not %D from its cone", point, cone[c]);
8579: }
8580: }
8581: }
8582: return(0);
8583: }
8585: PetscErrorCode DMPlexCheckAll_Internal(DM dm, PetscInt cellHeight)
8586: {
8590: DMPlexCheckSymmetry(dm);
8591: DMPlexCheckSkeleton(dm, cellHeight);
8592: DMPlexCheckFaces(dm, cellHeight);
8593: DMPlexCheckGeometry(dm);
8594: DMPlexCheckPointSF(dm);
8595: DMPlexCheckInterfaceCones(dm);
8596: return(0);
8597: }
8599: typedef struct cell_stats
8600: {
8601: PetscReal min, max, sum, squaresum;
8602: PetscInt count;
8603: } cell_stats_t;
8605: static void MPIAPI cell_stats_reduce(void *a, void *b, int * len, MPI_Datatype *datatype)
8606: {
8607: PetscInt i, N = *len;
8609: for (i = 0; i < N; i++) {
8610: cell_stats_t *A = (cell_stats_t *) a;
8611: cell_stats_t *B = (cell_stats_t *) b;
8613: B->min = PetscMin(A->min,B->min);
8614: B->max = PetscMax(A->max,B->max);
8615: B->sum += A->sum;
8616: B->squaresum += A->squaresum;
8617: B->count += A->count;
8618: }
8619: }
8621: /*@
8622: DMPlexCheckCellShape - Checks the Jacobian of the mapping from reference to real cells and computes some minimal statistics.
8624: Collective on dm
8626: Input Parameters:
8627: + dm - The DMPlex object
8628: . output - If true, statistics will be displayed on stdout
8629: - condLimit - Display all cells above this condition number, or PETSC_DETERMINE for no cell output
8631: Notes:
8632: This is mainly intended for debugging/testing purposes.
8634: For the complete list of DMPlexCheck* functions, see DMSetFromOptions().
8636: Level: developer
8638: .seealso: DMSetFromOptions(), DMPlexComputeOrthogonalQuality()
8639: @*/
8640: PetscErrorCode DMPlexCheckCellShape(DM dm, PetscBool output, PetscReal condLimit)
8641: {
8642: DM dmCoarse;
8643: cell_stats_t stats, globalStats;
8644: MPI_Comm comm = PetscObjectComm((PetscObject)dm);
8645: PetscReal *J, *invJ, min = 0, max = 0, mean = 0, stdev = 0;
8646: PetscReal limit = condLimit > 0 ? condLimit : PETSC_MAX_REAL;
8647: PetscInt cdim, cStart, cEnd, c, eStart, eEnd, count = 0;
8648: PetscMPIInt rank,size;
8653: stats.min = PETSC_MAX_REAL;
8654: stats.max = PETSC_MIN_REAL;
8655: stats.sum = stats.squaresum = 0.;
8656: stats.count = 0;
8658: MPI_Comm_size(comm, &size);
8659: MPI_Comm_rank(comm, &rank);
8660: DMGetCoordinateDim(dm,&cdim);
8661: PetscMalloc2(PetscSqr(cdim), &J, PetscSqr(cdim), &invJ);
8662: DMPlexGetSimplexOrBoxCells(dm,0,&cStart,&cEnd);
8663: DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);
8664: for (c = cStart; c < cEnd; c++) {
8665: PetscInt i;
8666: PetscReal frobJ = 0., frobInvJ = 0., cond2, cond, detJ;
8668: DMPlexComputeCellGeometryAffineFEM(dm,c,NULL,J,invJ,&detJ);
8669: if (detJ < 0.0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh cell %D is inverted", c);
8670: for (i = 0; i < PetscSqr(cdim); ++i) {
8671: frobJ += J[i] * J[i];
8672: frobInvJ += invJ[i] * invJ[i];
8673: }
8674: cond2 = frobJ * frobInvJ;
8675: cond = PetscSqrtReal(cond2);
8677: stats.min = PetscMin(stats.min,cond);
8678: stats.max = PetscMax(stats.max,cond);
8679: stats.sum += cond;
8680: stats.squaresum += cond2;
8681: stats.count++;
8682: if (output && cond > limit) {
8683: PetscSection coordSection;
8684: Vec coordsLocal;
8685: PetscScalar *coords = NULL;
8686: PetscInt Nv, d, clSize, cl, *closure = NULL;
8688: DMGetCoordinatesLocal(dm, &coordsLocal);
8689: DMGetCoordinateSection(dm, &coordSection);
8690: DMPlexVecGetClosure(dm, coordSection, coordsLocal, c, &Nv, &coords);
8691: PetscSynchronizedPrintf(comm, "[%d] Cell %D cond %g\n", rank, c, (double) cond);
8692: for (i = 0; i < Nv/cdim; ++i) {
8693: PetscSynchronizedPrintf(comm, " Vertex %D: (", i);
8694: for (d = 0; d < cdim; ++d) {
8695: if (d > 0) {PetscSynchronizedPrintf(comm, ", ");}
8696: PetscSynchronizedPrintf(comm, "%g", (double) PetscRealPart(coords[i*cdim+d]));
8697: }
8698: PetscSynchronizedPrintf(comm, ")\n");
8699: }
8700: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);
8701: for (cl = 0; cl < clSize*2; cl += 2) {
8702: const PetscInt edge = closure[cl];
8704: if ((edge >= eStart) && (edge < eEnd)) {
8705: PetscReal len;
8707: DMPlexComputeCellGeometryFVM(dm, edge, &len, NULL, NULL);
8708: PetscSynchronizedPrintf(comm, " Edge %D: length %g\n", edge, (double) len);
8709: }
8710: }
8711: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);
8712: DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, c, &Nv, &coords);
8713: }
8714: }
8715: if (output) {PetscSynchronizedFlush(comm, NULL);}
8717: if (size > 1) {
8718: PetscMPIInt blockLengths[2] = {4,1};
8719: MPI_Aint blockOffsets[2] = {offsetof(cell_stats_t,min),offsetof(cell_stats_t,count)};
8720: MPI_Datatype blockTypes[2] = {MPIU_REAL,MPIU_INT}, statType;
8721: MPI_Op statReduce;
8723: MPI_Type_create_struct(2,blockLengths,blockOffsets,blockTypes,&statType);
8724: MPI_Type_commit(&statType);
8725: MPI_Op_create(cell_stats_reduce, PETSC_TRUE, &statReduce);
8726: MPI_Reduce(&stats,&globalStats,1,statType,statReduce,0,comm);
8727: MPI_Op_free(&statReduce);
8728: MPI_Type_free(&statType);
8729: } else {
8730: PetscArraycpy(&globalStats,&stats,1);
8731: }
8732: if (rank == 0) {
8733: count = globalStats.count;
8734: min = globalStats.min;
8735: max = globalStats.max;
8736: mean = globalStats.sum / globalStats.count;
8737: stdev = globalStats.count > 1 ? PetscSqrtReal(PetscMax((globalStats.squaresum - globalStats.count * mean * mean) / (globalStats.count - 1),0)) : 0.0;
8738: }
8740: if (output) {
8741: PetscPrintf(comm,"Mesh with %D cells, shape condition numbers: min = %g, max = %g, mean = %g, stddev = %g\n", count, (double) min, (double) max, (double) mean, (double) stdev);
8742: }
8743: PetscFree2(J,invJ);
8745: DMGetCoarseDM(dm,&dmCoarse);
8746: if (dmCoarse) {
8747: PetscBool isplex;
8749: PetscObjectTypeCompare((PetscObject)dmCoarse,DMPLEX,&isplex);
8750: if (isplex) {
8751: DMPlexCheckCellShape(dmCoarse,output,condLimit);
8752: }
8753: }
8754: return(0);
8755: }
8757: /*@
8758: DMPlexComputeOrthogonalQuality - Compute cell-wise orthogonal quality mesh statistic. Optionally tags all cells with
8759: orthogonal quality below given tolerance.
8761: Collective on dm
8763: Input Parameters:
8764: + dm - The DMPlex object
8765: . fv - Optional PetscFV object for pre-computed cell/face centroid information
8766: - atol - [0, 1] Absolute tolerance for tagging cells.
8768: Output Parameters:
8769: + OrthQual - Vec containing orthogonal quality per cell
8770: - OrthQualLabel - DMLabel tagging cells below atol with DM_ADAPT_REFINE
8772: Options Database Keys:
8773: + -dm_plex_orthogonal_quality_label_view - view OrthQualLabel if label is requested. Currently only PETSCVIEWERASCII is
8774: supported.
8775: - -dm_plex_orthogonal_quality_vec_view - view OrthQual vector.
8777: Notes:
8778: Orthogonal quality is given by the following formula:
8780: \min \left[ \frac{A_i \cdot f_i}{\|A_i\| \|f_i\|} , \frac{A_i \cdot c_i}{\|A_i\| \|c_i\|} \right]
8782: Where A_i is the i'th face-normal vector, f_i is the vector from the cell centroid to the i'th face centroid, and c_i
8783: is the vector from the current cells centroid to the centroid of its i'th neighbor (which shares a face with the
8784: current cell). This computes the vector similarity between each cell face and its corresponding neighbor centroid by
8785: calculating the cosine of the angle between these vectors.
8787: Orthogonal quality ranges from 1 (best) to 0 (worst).
8789: This routine is mainly useful for FVM, however is not restricted to only FVM. The PetscFV object is optionally used to check for
8790: pre-computed FVM cell data, but if it is not passed in then this data will be computed.
8792: Cells are tagged if they have an orthogonal quality less than or equal to the absolute tolerance.
8794: Level: intermediate
8796: .seealso: DMPlexCheckCellShape(), DMCreateLabel()
8797: @*/
8798: PetscErrorCode DMPlexComputeOrthogonalQuality(DM dm, PetscFV fv, PetscReal atol, Vec *OrthQual, DMLabel *OrthQualLabel)
8799: {
8800: PetscInt nc, cellHeight, cStart, cEnd, cell, cellIter = 0;
8801: PetscInt *idx;
8802: PetscScalar *oqVals;
8803: const PetscScalar *cellGeomArr, *faceGeomArr;
8804: PetscReal *ci, *fi, *Ai;
8805: MPI_Comm comm;
8806: Vec cellgeom, facegeom;
8807: DM dmFace, dmCell;
8808: IS glob;
8809: ISLocalToGlobalMapping ltog;
8810: PetscViewer vwr;
8811: PetscErrorCode ierr;
8817: if (PetscUnlikelyDebug(atol < 0.0 || atol > 1.0)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g not in [0,1]",(double)atol);
8818: PetscObjectGetComm((PetscObject) dm, &comm);
8819: DMGetDimension(dm, &nc);
8820: if (nc < 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DM must have dimension >= 2 (current %D)", nc);
8821: {
8822: DMPlexInterpolatedFlag interpFlag;
8824: DMPlexIsInterpolated(dm, &interpFlag);
8825: if (interpFlag != DMPLEX_INTERPOLATED_FULL) {
8826: PetscMPIInt rank;
8828: MPI_Comm_rank(comm, &rank);
8829: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DM must be fully interpolated, DM on rank %d is not fully interpolated", rank);
8830: }
8831: }
8832: if (OrthQualLabel) {
8834: DMCreateLabel(dm, "Orthogonal_Quality");
8835: DMGetLabel(dm, "Orthogonal_Quality", OrthQualLabel);
8836: } else {*OrthQualLabel = NULL;}
8837: DMPlexGetVTKCellHeight(dm, &cellHeight);
8838: DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);
8839: DMPlexCreateCellNumbering_Internal(dm, PETSC_TRUE, &glob);
8840: ISLocalToGlobalMappingCreateIS(glob, <og);
8841: ISLocalToGlobalMappingSetType(ltog, ISLOCALTOGLOBALMAPPINGHASH);
8842: VecCreate(comm, OrthQual);
8843: VecSetType(*OrthQual, VECSTANDARD);
8844: VecSetSizes(*OrthQual, cEnd-cStart, PETSC_DETERMINE);
8845: VecSetLocalToGlobalMapping(*OrthQual, ltog);
8846: VecSetUp(*OrthQual);
8847: ISDestroy(&glob);
8848: ISLocalToGlobalMappingDestroy(<og);
8849: DMPlexGetDataFVM(dm, fv, &cellgeom, &facegeom, NULL);
8850: VecGetArrayRead(cellgeom, &cellGeomArr);
8851: VecGetArrayRead(facegeom, &faceGeomArr);
8852: VecGetDM(cellgeom, &dmCell);
8853: VecGetDM(facegeom, &dmFace);
8854: PetscMalloc5(cEnd-cStart, &idx, cEnd-cStart, &oqVals, nc, &ci, nc, &fi, nc, &Ai);
8855: for (cell = cStart; cell < cEnd; cellIter++,cell++) {
8856: PetscInt cellneigh, cellneighiter = 0, adjSize = PETSC_DETERMINE;
8857: PetscInt cellarr[2], *adj = NULL;
8858: PetscScalar *cArr, *fArr;
8859: PetscReal minvalc = 1.0, minvalf = 1.0;
8860: PetscFVCellGeom *cg;
8862: idx[cellIter] = cell-cStart;
8863: cellarr[0] = cell;
8864: /* Make indexing into cellGeom easier */
8865: DMPlexPointLocalRead(dmCell, cell, cellGeomArr, &cg);
8866: DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &adjSize, &adj);
8867: /* Technically 1 too big, but easier than fiddling with empty adjacency array */
8868: PetscCalloc2(adjSize, &cArr, adjSize, &fArr);
8869: for (cellneigh = 0; cellneigh < adjSize; cellneighiter++,cellneigh++) {
8870: PetscInt i;
8871: const PetscInt neigh = adj[cellneigh];
8872: PetscReal normci = 0, normfi = 0, normai = 0;
8873: PetscFVCellGeom *cgneigh;
8874: PetscFVFaceGeom *fg;
8876: /* Don't count ourselves in the neighbor list */
8877: if (neigh == cell) continue;
8878: DMPlexPointLocalRead(dmCell, neigh, cellGeomArr, &cgneigh);
8879: cellarr[1] = neigh;
8880: {
8881: PetscInt numcovpts;
8882: const PetscInt *covpts;
8884: DMPlexGetMeet(dm, 2, cellarr, &numcovpts, &covpts);
8885: DMPlexPointLocalRead(dmFace, covpts[0], faceGeomArr, &fg);
8886: DMPlexRestoreMeet(dm, 2, cellarr, &numcovpts, &covpts);
8887: }
8889: /* Compute c_i, f_i and their norms */
8890: for (i = 0; i < nc; i++) {
8891: ci[i] = cgneigh->centroid[i] - cg->centroid[i];
8892: fi[i] = fg->centroid[i] - cg->centroid[i];
8893: Ai[i] = fg->normal[i];
8894: normci += PetscPowReal(ci[i], 2);
8895: normfi += PetscPowReal(fi[i], 2);
8896: normai += PetscPowReal(Ai[i], 2);
8897: }
8898: normci = PetscSqrtReal(normci);
8899: normfi = PetscSqrtReal(normfi);
8900: normai = PetscSqrtReal(normai);
8902: /* Normalize and compute for each face-cell-normal pair */
8903: for (i = 0; i < nc; i++) {
8904: ci[i] = ci[i]/normci;
8905: fi[i] = fi[i]/normfi;
8906: Ai[i] = Ai[i]/normai;
8907: /* PetscAbs because I don't know if normals are guaranteed to point out */
8908: cArr[cellneighiter] += PetscAbs(Ai[i]*ci[i]);
8909: fArr[cellneighiter] += PetscAbs(Ai[i]*fi[i]);
8910: }
8911: if (PetscRealPart(cArr[cellneighiter]) < minvalc) {
8912: minvalc = PetscRealPart(cArr[cellneighiter]);
8913: }
8914: if (PetscRealPart(fArr[cellneighiter]) < minvalf) {
8915: minvalf = PetscRealPart(fArr[cellneighiter]);
8916: }
8917: }
8918: PetscFree(adj);
8919: PetscFree2(cArr, fArr);
8920: /* Defer to cell if they're equal */
8921: oqVals[cellIter] = PetscMin(minvalf, minvalc);
8922: if (OrthQualLabel) {
8923: if (PetscRealPart(oqVals[cellIter]) <= atol) {DMLabelSetValue(*OrthQualLabel, cell, DM_ADAPT_REFINE);}
8924: }
8925: }
8926: VecSetValuesLocal(*OrthQual, cEnd-cStart, idx, oqVals, INSERT_VALUES);
8927: VecAssemblyBegin(*OrthQual);
8928: VecAssemblyEnd(*OrthQual);
8929: VecRestoreArrayRead(cellgeom, &cellGeomArr);
8930: VecRestoreArrayRead(facegeom, &faceGeomArr);
8931: PetscOptionsGetViewer(comm, NULL, NULL, "-dm_plex_orthogonal_quality_label_view", &vwr, NULL, NULL);
8932: if (OrthQualLabel) {
8933: if (vwr) {DMLabelView(*OrthQualLabel, vwr);}
8934: }
8935: PetscFree5(idx, oqVals, ci, fi, Ai);
8936: PetscViewerDestroy(&vwr);
8937: VecViewFromOptions(*OrthQual, NULL, "-dm_plex_orthogonal_quality_vec_view");
8938: return(0);
8939: }
8941: /* this is here insead of DMGetOutputDM because output DM still has constraints in the local indices that affect
8942: * interpolator construction */
8943: static PetscErrorCode DMGetFullDM(DM dm, DM *odm)
8944: {
8945: PetscSection section, newSection, gsection;
8946: PetscSF sf;
8947: PetscBool hasConstraints, ghasConstraints;
8953: DMGetLocalSection(dm, §ion);
8954: PetscSectionHasConstraints(section, &hasConstraints);
8955: MPI_Allreduce(&hasConstraints, &ghasConstraints, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject) dm));
8956: if (!ghasConstraints) {
8957: PetscObjectReference((PetscObject)dm);
8958: *odm = dm;
8959: return(0);
8960: }
8961: DMClone(dm, odm);
8962: DMCopyFields(dm, *odm);
8963: DMGetLocalSection(*odm, &newSection);
8964: DMGetPointSF(*odm, &sf);
8965: PetscSectionCreateGlobalSection(newSection, sf, PETSC_TRUE, PETSC_FALSE, &gsection);
8966: DMSetGlobalSection(*odm, gsection);
8967: PetscSectionDestroy(&gsection);
8968: return(0);
8969: }
8971: static PetscErrorCode DMCreateAffineInterpolationCorrection_Plex(DM dmc, DM dmf, Vec *shift)
8972: {
8973: DM dmco, dmfo;
8974: Mat interpo;
8975: Vec rscale;
8976: Vec cglobalo, clocal;
8977: Vec fglobal, fglobalo, flocal;
8978: PetscBool regular;
8982: DMGetFullDM(dmc, &dmco);
8983: DMGetFullDM(dmf, &dmfo);
8984: DMSetCoarseDM(dmfo, dmco);
8985: DMPlexGetRegularRefinement(dmf, ®ular);
8986: DMPlexSetRegularRefinement(dmfo, regular);
8987: DMCreateInterpolation(dmco, dmfo, &interpo, &rscale);
8988: DMCreateGlobalVector(dmco, &cglobalo);
8989: DMCreateLocalVector(dmc, &clocal);
8990: VecSet(cglobalo, 0.);
8991: VecSet(clocal, 0.);
8992: DMCreateGlobalVector(dmf, &fglobal);
8993: DMCreateGlobalVector(dmfo, &fglobalo);
8994: DMCreateLocalVector(dmf, &flocal);
8995: VecSet(fglobal, 0.);
8996: VecSet(fglobalo, 0.);
8997: VecSet(flocal, 0.);
8998: DMPlexInsertBoundaryValues(dmc, PETSC_TRUE, clocal, 0., NULL, NULL, NULL);
8999: DMLocalToGlobalBegin(dmco, clocal, INSERT_VALUES, cglobalo);
9000: DMLocalToGlobalEnd(dmco, clocal, INSERT_VALUES, cglobalo);
9001: MatMult(interpo, cglobalo, fglobalo);
9002: DMGlobalToLocalBegin(dmfo, fglobalo, INSERT_VALUES, flocal);
9003: DMGlobalToLocalEnd(dmfo, fglobalo, INSERT_VALUES, flocal);
9004: DMLocalToGlobalBegin(dmf, flocal, INSERT_VALUES, fglobal);
9005: DMLocalToGlobalEnd(dmf, flocal, INSERT_VALUES, fglobal);
9006: *shift = fglobal;
9007: VecDestroy(&flocal);
9008: VecDestroy(&fglobalo);
9009: VecDestroy(&clocal);
9010: VecDestroy(&cglobalo);
9011: VecDestroy(&rscale);
9012: MatDestroy(&interpo);
9013: DMDestroy(&dmfo);
9014: DMDestroy(&dmco);
9015: return(0);
9016: }
9018: PETSC_INTERN PetscErrorCode DMInterpolateSolution_Plex(DM coarse, DM fine, Mat interp, Vec coarseSol, Vec fineSol)
9019: {
9020: PetscObject shifto;
9021: Vec shift;
9026: if (!interp) {
9027: Vec rscale;
9029: DMCreateInterpolation(coarse, fine, &interp, &rscale);
9030: VecDestroy(&rscale);
9031: } else {
9032: PetscObjectReference((PetscObject)interp);
9033: }
9034: PetscObjectQuery((PetscObject)interp, "_DMInterpolateSolution_Plex_Vec", &shifto);
9035: if (!shifto) {
9036: DMCreateAffineInterpolationCorrection_Plex(coarse, fine, &shift);
9037: PetscObjectCompose((PetscObject)interp, "_DMInterpolateSolution_Plex_Vec", (PetscObject) shift);
9038: shifto = (PetscObject) shift;
9039: VecDestroy(&shift);
9040: }
9041: shift = (Vec) shifto;
9042: MatInterpolate(interp, coarseSol, fineSol);
9043: VecAXPY(fineSol, 1.0, shift);
9044: MatDestroy(&interp);
9045: return(0);
9046: }
9048: /* Pointwise interpolation
9049: Just code FEM for now
9050: u^f = I u^c
9051: sum_k u^f_k phi^f_k = I sum_j u^c_j phi^c_j
9052: u^f_i = sum_j psi^f_i I phi^c_j u^c_j
9053: I_{ij} = psi^f_i phi^c_j
9054: */
9055: PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling)
9056: {
9057: PetscSection gsc, gsf;
9058: PetscInt m, n;
9059: void *ctx;
9060: DM cdm;
9061: PetscBool regular, ismatis, isRefined = dmCoarse->data == dmFine->data ? PETSC_FALSE : PETSC_TRUE;
9065: DMGetGlobalSection(dmFine, &gsf);
9066: PetscSectionGetConstrainedStorageSize(gsf, &m);
9067: DMGetGlobalSection(dmCoarse, &gsc);
9068: PetscSectionGetConstrainedStorageSize(gsc, &n);
9070: PetscStrcmp(dmCoarse->mattype, MATIS, &ismatis);
9071: MatCreate(PetscObjectComm((PetscObject) dmCoarse), interpolation);
9072: MatSetSizes(*interpolation, m, n, PETSC_DETERMINE, PETSC_DETERMINE);
9073: MatSetType(*interpolation, ismatis ? MATAIJ : dmCoarse->mattype);
9074: DMGetApplicationContext(dmFine, &ctx);
9076: DMGetCoarseDM(dmFine, &cdm);
9077: DMPlexGetRegularRefinement(dmFine, ®ular);
9078: if (!isRefined || (regular && cdm == dmCoarse)) {DMPlexComputeInterpolatorNested(dmCoarse, dmFine, isRefined, *interpolation, ctx);}
9079: else {DMPlexComputeInterpolatorGeneral(dmCoarse, dmFine, *interpolation, ctx);}
9080: MatViewFromOptions(*interpolation, NULL, "-interp_mat_view");
9081: if (scaling) {
9082: /* Use naive scaling */
9083: DMCreateInterpolationScale(dmCoarse, dmFine, *interpolation, scaling);
9084: }
9085: return(0);
9086: }
9088: PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat)
9089: {
9091: VecScatter ctx;
9094: DMPlexComputeInjectorFEM(dmCoarse, dmFine, &ctx, NULL);
9095: MatCreateScatter(PetscObjectComm((PetscObject)ctx), ctx, mat);
9096: VecScatterDestroy(&ctx);
9097: return(0);
9098: }
9100: static void g0_identity_private(PetscInt dim, PetscInt Nf, PetscInt NfAux,
9101: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
9102: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
9103: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
9104: {
9105: g0[0] = 1.0;
9106: }
9108: PetscErrorCode DMCreateMassMatrix_Plex(DM dmCoarse, DM dmFine, Mat *mass)
9109: {
9110: PetscSection gsc, gsf;
9111: PetscInt m, n;
9112: void *ctx;
9113: DM cdm;
9114: PetscBool regular;
9118: if (dmFine == dmCoarse) {
9119: DM dmc;
9120: PetscDS ds;
9121: Vec u;
9122: IS cellIS;
9123: PetscFormKey key;
9124: PetscInt depth;
9126: DMClone(dmFine, &dmc);
9127: DMCopyDisc(dmFine, dmc);
9128: DMGetDS(dmc, &ds);
9129: PetscDSSetJacobian(ds, 0, 0, g0_identity_private, NULL, NULL, NULL);
9130: DMCreateMatrix(dmc, mass);
9131: DMGetGlobalVector(dmc, &u);
9132: DMPlexGetDepth(dmc, &depth);
9133: DMGetStratumIS(dmc, "depth", depth, &cellIS);
9134: MatZeroEntries(*mass);
9135: key.label = NULL;
9136: key.value = 0;
9137: key.field = 0;
9138: key.part = 0;
9139: DMPlexComputeJacobian_Internal(dmc, key, cellIS, 0.0, 0.0, u, NULL, *mass, *mass, NULL);
9140: ISDestroy(&cellIS);
9141: DMRestoreGlobalVector(dmc, &u);
9142: DMDestroy(&dmc);
9143: } else {
9144: DMGetGlobalSection(dmFine, &gsf);
9145: PetscSectionGetConstrainedStorageSize(gsf, &m);
9146: DMGetGlobalSection(dmCoarse, &gsc);
9147: PetscSectionGetConstrainedStorageSize(gsc, &n);
9149: MatCreate(PetscObjectComm((PetscObject) dmCoarse), mass);
9150: MatSetSizes(*mass, m, n, PETSC_DETERMINE, PETSC_DETERMINE);
9151: MatSetType(*mass, dmCoarse->mattype);
9152: DMGetApplicationContext(dmFine, &ctx);
9154: DMGetCoarseDM(dmFine, &cdm);
9155: DMPlexGetRegularRefinement(dmFine, ®ular);
9156: if (regular && cdm == dmCoarse) {DMPlexComputeMassMatrixNested(dmCoarse, dmFine, *mass, ctx);}
9157: else {DMPlexComputeMassMatrixGeneral(dmCoarse, dmFine, *mass, ctx);}
9158: }
9159: MatViewFromOptions(*mass, NULL, "-mass_mat_view");
9160: return(0);
9161: }
9163: /*@
9164: DMPlexGetRegularRefinement - Get the flag indicating that this mesh was obtained by regular refinement from its coarse mesh
9166: Input Parameter:
9167: . dm - The DMPlex object
9169: Output Parameter:
9170: . regular - The flag
9172: Level: intermediate
9174: .seealso: DMPlexSetRegularRefinement()
9175: @*/
9176: PetscErrorCode DMPlexGetRegularRefinement(DM dm, PetscBool *regular)
9177: {
9181: *regular = ((DM_Plex *) dm->data)->regularRefinement;
9182: return(0);
9183: }
9185: /*@
9186: DMPlexSetRegularRefinement - Set the flag indicating that this mesh was obtained by regular refinement from its coarse mesh
9188: Input Parameters:
9189: + dm - The DMPlex object
9190: - regular - The flag
9192: Level: intermediate
9194: .seealso: DMPlexGetRegularRefinement()
9195: @*/
9196: PetscErrorCode DMPlexSetRegularRefinement(DM dm, PetscBool regular)
9197: {
9200: ((DM_Plex *) dm->data)->regularRefinement = regular;
9201: return(0);
9202: }
9204: /* anchors */
9205: /*@
9206: DMPlexGetAnchors - Get the layout of the anchor (point-to-point) constraints. Typically, the user will not have to
9207: call DMPlexGetAnchors() directly: if there are anchors, then DMPlexGetAnchors() is called during DMGetConstraints().
9209: not collective
9211: Input Parameter:
9212: . dm - The DMPlex object
9214: Output Parameters:
9215: + anchorSection - If not NULL, set to the section describing which points anchor the constrained points.
9216: - anchorIS - If not NULL, set to the list of anchors indexed by anchorSection
9218: Level: intermediate
9220: .seealso: DMPlexSetAnchors(), DMGetConstraints(), DMSetConstraints()
9221: @*/
9222: PetscErrorCode DMPlexGetAnchors(DM dm, PetscSection *anchorSection, IS *anchorIS)
9223: {
9224: DM_Plex *plex = (DM_Plex *)dm->data;
9229: if (!plex->anchorSection && !plex->anchorIS && plex->createanchors) {(*plex->createanchors)(dm);}
9230: if (anchorSection) *anchorSection = plex->anchorSection;
9231: if (anchorIS) *anchorIS = plex->anchorIS;
9232: return(0);
9233: }
9235: /*@
9236: DMPlexSetAnchors - Set the layout of the local anchor (point-to-point) constraints. Unlike boundary conditions,
9237: when a point's degrees of freedom in a section are constrained to an outside value, the anchor constraints set a
9238: point's degrees of freedom to be a linear combination of other points' degrees of freedom.
9240: After specifying the layout of constraints with DMPlexSetAnchors(), one specifies the constraints by calling
9241: DMGetConstraints() and filling in the entries in the constraint matrix.
9243: collective on dm
9245: Input Parameters:
9246: + dm - The DMPlex object
9247: . anchorSection - The section that describes the mapping from constrained points to the anchor points listed in anchorIS. Must have a local communicator (PETSC_COMM_SELF or derivative).
9248: - anchorIS - The list of all anchor points. Must have a local communicator (PETSC_COMM_SELF or derivative).
9250: The reference counts of anchorSection and anchorIS are incremented.
9252: Level: intermediate
9254: .seealso: DMPlexGetAnchors(), DMGetConstraints(), DMSetConstraints()
9255: @*/
9256: PetscErrorCode DMPlexSetAnchors(DM dm, PetscSection anchorSection, IS anchorIS)
9257: {
9258: DM_Plex *plex = (DM_Plex *)dm->data;
9259: PetscMPIInt result;
9264: if (anchorSection) {
9266: MPI_Comm_compare(PETSC_COMM_SELF,PetscObjectComm((PetscObject)anchorSection),&result);
9267: if (result != MPI_CONGRUENT && result != MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMECOMM,"anchor section must have local communicator");
9268: }
9269: if (anchorIS) {
9271: MPI_Comm_compare(PETSC_COMM_SELF,PetscObjectComm((PetscObject)anchorIS),&result);
9272: if (result != MPI_CONGRUENT && result != MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMECOMM,"anchor IS must have local communicator");
9273: }
9275: PetscObjectReference((PetscObject)anchorSection);
9276: PetscSectionDestroy(&plex->anchorSection);
9277: plex->anchorSection = anchorSection;
9279: PetscObjectReference((PetscObject)anchorIS);
9280: ISDestroy(&plex->anchorIS);
9281: plex->anchorIS = anchorIS;
9283: if (PetscUnlikelyDebug(anchorIS && anchorSection)) {
9284: PetscInt size, a, pStart, pEnd;
9285: const PetscInt *anchors;
9287: PetscSectionGetChart(anchorSection,&pStart,&pEnd);
9288: ISGetLocalSize(anchorIS,&size);
9289: ISGetIndices(anchorIS,&anchors);
9290: for (a = 0; a < size; a++) {
9291: PetscInt p;
9293: p = anchors[a];
9294: if (p >= pStart && p < pEnd) {
9295: PetscInt dof;
9297: PetscSectionGetDof(anchorSection,p,&dof);
9298: if (dof) {
9299: PetscErrorCode ierr2;
9301: ierr2 = ISRestoreIndices(anchorIS,&anchors);CHKERRQ(ierr2);
9302: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Point %D cannot be constrained and an anchor",p);
9303: }
9304: }
9305: }
9306: ISRestoreIndices(anchorIS,&anchors);
9307: }
9308: /* reset the generic constraints */
9309: DMSetDefaultConstraints(dm,NULL,NULL);
9310: return(0);
9311: }
9313: static PetscErrorCode DMPlexCreateConstraintSection_Anchors(DM dm, PetscSection section, PetscSection *cSec)
9314: {
9315: PetscSection anchorSection;
9316: PetscInt pStart, pEnd, sStart, sEnd, p, dof, numFields, f;
9321: DMPlexGetAnchors(dm,&anchorSection,NULL);
9322: PetscSectionCreate(PETSC_COMM_SELF,cSec);
9323: PetscSectionGetNumFields(section,&numFields);
9324: if (numFields) {
9325: PetscInt f;
9326: PetscSectionSetNumFields(*cSec,numFields);
9328: for (f = 0; f < numFields; f++) {
9329: PetscInt numComp;
9331: PetscSectionGetFieldComponents(section,f,&numComp);
9332: PetscSectionSetFieldComponents(*cSec,f,numComp);
9333: }
9334: }
9335: PetscSectionGetChart(anchorSection,&pStart,&pEnd);
9336: PetscSectionGetChart(section,&sStart,&sEnd);
9337: pStart = PetscMax(pStart,sStart);
9338: pEnd = PetscMin(pEnd,sEnd);
9339: pEnd = PetscMax(pStart,pEnd);
9340: PetscSectionSetChart(*cSec,pStart,pEnd);
9341: for (p = pStart; p < pEnd; p++) {
9342: PetscSectionGetDof(anchorSection,p,&dof);
9343: if (dof) {
9344: PetscSectionGetDof(section,p,&dof);
9345: PetscSectionSetDof(*cSec,p,dof);
9346: for (f = 0; f < numFields; f++) {
9347: PetscSectionGetFieldDof(section,p,f,&dof);
9348: PetscSectionSetFieldDof(*cSec,p,f,dof);
9349: }
9350: }
9351: }
9352: PetscSectionSetUp(*cSec);
9353: PetscObjectSetName((PetscObject) *cSec, "Constraint Section");
9354: return(0);
9355: }
9357: static PetscErrorCode DMPlexCreateConstraintMatrix_Anchors(DM dm, PetscSection section, PetscSection cSec, Mat *cMat)
9358: {
9359: PetscSection aSec;
9360: PetscInt pStart, pEnd, p, sStart, sEnd, dof, aDof, aOff, off, nnz, annz, m, n, q, a, offset, *i, *j;
9361: const PetscInt *anchors;
9362: PetscInt numFields, f;
9363: IS aIS;
9365: MatType mtype;
9366: PetscBool iscuda,iskokkos;
9370: PetscSectionGetStorageSize(cSec, &m);
9371: PetscSectionGetStorageSize(section, &n);
9372: MatCreate(PETSC_COMM_SELF,cMat);
9373: MatSetSizes(*cMat,m,n,m,n);
9374: PetscStrcmp(dm->mattype,MATSEQAIJCUSPARSE,&iscuda);
9375: if (!iscuda) { PetscStrcmp(dm->mattype,MATMPIAIJCUSPARSE,&iscuda); }
9376: PetscStrcmp(dm->mattype,MATSEQAIJKOKKOS,&iskokkos);
9377: if (!iskokkos) { PetscStrcmp(dm->mattype,MATMPIAIJKOKKOS,&iskokkos); }
9378: if (iscuda) mtype = MATSEQAIJCUSPARSE;
9379: else if (iskokkos) mtype = MATSEQAIJKOKKOS;
9380: else mtype = MATSEQAIJ;
9381: MatSetType(*cMat,mtype);
9382: DMPlexGetAnchors(dm,&aSec,&aIS);
9383: ISGetIndices(aIS,&anchors);
9384: /* cSec will be a subset of aSec and section */
9385: PetscSectionGetChart(cSec,&pStart,&pEnd);
9386: PetscSectionGetChart(section,&sStart,&sEnd);
9387: PetscMalloc1(m+1,&i);
9388: i[0] = 0;
9389: PetscSectionGetNumFields(section,&numFields);
9390: for (p = pStart; p < pEnd; p++) {
9391: PetscInt rDof, rOff, r;
9393: PetscSectionGetDof(aSec,p,&rDof);
9394: if (!rDof) continue;
9395: PetscSectionGetOffset(aSec,p,&rOff);
9396: if (numFields) {
9397: for (f = 0; f < numFields; f++) {
9398: annz = 0;
9399: for (r = 0; r < rDof; r++) {
9400: a = anchors[rOff + r];
9401: if (a < sStart || a >= sEnd) continue;
9402: PetscSectionGetFieldDof(section,a,f,&aDof);
9403: annz += aDof;
9404: }
9405: PetscSectionGetFieldDof(cSec,p,f,&dof);
9406: PetscSectionGetFieldOffset(cSec,p,f,&off);
9407: for (q = 0; q < dof; q++) {
9408: i[off + q + 1] = i[off + q] + annz;
9409: }
9410: }
9411: }
9412: else {
9413: annz = 0;
9414: PetscSectionGetDof(cSec,p,&dof);
9415: for (q = 0; q < dof; q++) {
9416: a = anchors[rOff + q];
9417: if (a < sStart || a >= sEnd) continue;
9418: PetscSectionGetDof(section,a,&aDof);
9419: annz += aDof;
9420: }
9421: PetscSectionGetDof(cSec,p,&dof);
9422: PetscSectionGetOffset(cSec,p,&off);
9423: for (q = 0; q < dof; q++) {
9424: i[off + q + 1] = i[off + q] + annz;
9425: }
9426: }
9427: }
9428: nnz = i[m];
9429: PetscMalloc1(nnz,&j);
9430: offset = 0;
9431: for (p = pStart; p < pEnd; p++) {
9432: if (numFields) {
9433: for (f = 0; f < numFields; f++) {
9434: PetscSectionGetFieldDof(cSec,p,f,&dof);
9435: for (q = 0; q < dof; q++) {
9436: PetscInt rDof, rOff, r;
9437: PetscSectionGetDof(aSec,p,&rDof);
9438: PetscSectionGetOffset(aSec,p,&rOff);
9439: for (r = 0; r < rDof; r++) {
9440: PetscInt s;
9442: a = anchors[rOff + r];
9443: if (a < sStart || a >= sEnd) continue;
9444: PetscSectionGetFieldDof(section,a,f,&aDof);
9445: PetscSectionGetFieldOffset(section,a,f,&aOff);
9446: for (s = 0; s < aDof; s++) {
9447: j[offset++] = aOff + s;
9448: }
9449: }
9450: }
9451: }
9452: }
9453: else {
9454: PetscSectionGetDof(cSec,p,&dof);
9455: for (q = 0; q < dof; q++) {
9456: PetscInt rDof, rOff, r;
9457: PetscSectionGetDof(aSec,p,&rDof);
9458: PetscSectionGetOffset(aSec,p,&rOff);
9459: for (r = 0; r < rDof; r++) {
9460: PetscInt s;
9462: a = anchors[rOff + r];
9463: if (a < sStart || a >= sEnd) continue;
9464: PetscSectionGetDof(section,a,&aDof);
9465: PetscSectionGetOffset(section,a,&aOff);
9466: for (s = 0; s < aDof; s++) {
9467: j[offset++] = aOff + s;
9468: }
9469: }
9470: }
9471: }
9472: }
9473: MatSeqAIJSetPreallocationCSR(*cMat,i,j,NULL);
9474: PetscFree(i);
9475: PetscFree(j);
9476: ISRestoreIndices(aIS,&anchors);
9477: return(0);
9478: }
9480: PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm)
9481: {
9482: DM_Plex *plex = (DM_Plex *)dm->data;
9483: PetscSection anchorSection, section, cSec;
9484: Mat cMat;
9489: DMPlexGetAnchors(dm,&anchorSection,NULL);
9490: if (anchorSection) {
9491: PetscInt Nf;
9493: DMGetLocalSection(dm,§ion);
9494: DMPlexCreateConstraintSection_Anchors(dm,section,&cSec);
9495: DMPlexCreateConstraintMatrix_Anchors(dm,section,cSec,&cMat);
9496: DMGetNumFields(dm,&Nf);
9497: if (Nf && plex->computeanchormatrix) {(*plex->computeanchormatrix)(dm,section,cSec,cMat);}
9498: DMSetDefaultConstraints(dm,cSec,cMat);
9499: PetscSectionDestroy(&cSec);
9500: MatDestroy(&cMat);
9501: }
9502: return(0);
9503: }
9505: PetscErrorCode DMCreateSubDomainDM_Plex(DM dm, DMLabel label, PetscInt value, IS *is, DM *subdm)
9506: {
9507: IS subis;
9508: PetscSection section, subsection;
9512: DMGetLocalSection(dm, §ion);
9513: if (!section) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting subdomain");
9514: if (!subdm) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Must set output subDM for splitting subdomain");
9515: /* Create subdomain */
9516: DMPlexFilter(dm, label, value, subdm);
9517: /* Create submodel */
9518: DMPlexGetSubpointIS(*subdm, &subis);
9519: PetscSectionCreateSubmeshSection(section, subis, &subsection);
9520: DMSetLocalSection(*subdm, subsection);
9521: PetscSectionDestroy(&subsection);
9522: DMCopyDisc(dm, *subdm);
9523: /* Create map from submodel to global model */
9524: if (is) {
9525: PetscSection sectionGlobal, subsectionGlobal;
9526: IS spIS;
9527: const PetscInt *spmap;
9528: PetscInt *subIndices;
9529: PetscInt subSize = 0, subOff = 0, pStart, pEnd, p;
9530: PetscInt Nf, f, bs = -1, bsLocal[2], bsMinMax[2];
9532: DMPlexGetSubpointIS(*subdm, &spIS);
9533: ISGetIndices(spIS, &spmap);
9534: PetscSectionGetNumFields(section, &Nf);
9535: DMGetGlobalSection(dm, §ionGlobal);
9536: DMGetGlobalSection(*subdm, &subsectionGlobal);
9537: PetscSectionGetChart(subsection, &pStart, &pEnd);
9538: for (p = pStart; p < pEnd; ++p) {
9539: PetscInt gdof, pSubSize = 0;
9541: PetscSectionGetDof(sectionGlobal, p, &gdof);
9542: if (gdof > 0) {
9543: for (f = 0; f < Nf; ++f) {
9544: PetscInt fdof, fcdof;
9546: PetscSectionGetFieldDof(subsection, p, f, &fdof);
9547: PetscSectionGetFieldConstraintDof(subsection, p, f, &fcdof);
9548: pSubSize += fdof-fcdof;
9549: }
9550: subSize += pSubSize;
9551: if (pSubSize) {
9552: if (bs < 0) {
9553: bs = pSubSize;
9554: } else if (bs != pSubSize) {
9555: /* Layout does not admit a pointwise block size */
9556: bs = 1;
9557: }
9558: }
9559: }
9560: }
9561: /* Must have same blocksize on all procs (some might have no points) */
9562: bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; bsLocal[1] = bs;
9563: PetscGlobalMinMaxInt(PetscObjectComm((PetscObject) dm), bsLocal, bsMinMax);
9564: if (bsMinMax[0] != bsMinMax[1]) {bs = 1;}
9565: else {bs = bsMinMax[0];}
9566: PetscMalloc1(subSize, &subIndices);
9567: for (p = pStart; p < pEnd; ++p) {
9568: PetscInt gdof, goff;
9570: PetscSectionGetDof(subsectionGlobal, p, &gdof);
9571: if (gdof > 0) {
9572: const PetscInt point = spmap[p];
9574: PetscSectionGetOffset(sectionGlobal, point, &goff);
9575: for (f = 0; f < Nf; ++f) {
9576: PetscInt fdof, fcdof, fc, f2, poff = 0;
9578: /* Can get rid of this loop by storing field information in the global section */
9579: for (f2 = 0; f2 < f; ++f2) {
9580: PetscSectionGetFieldDof(section, p, f2, &fdof);
9581: PetscSectionGetFieldConstraintDof(section, p, f2, &fcdof);
9582: poff += fdof-fcdof;
9583: }
9584: PetscSectionGetFieldDof(section, p, f, &fdof);
9585: PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);
9586: for (fc = 0; fc < fdof-fcdof; ++fc, ++subOff) {
9587: subIndices[subOff] = goff+poff+fc;
9588: }
9589: }
9590: }
9591: }
9592: ISRestoreIndices(spIS, &spmap);
9593: ISCreateGeneral(PetscObjectComm((PetscObject)dm), subSize, subIndices, PETSC_OWN_POINTER, is);
9594: if (bs > 1) {
9595: /* We need to check that the block size does not come from non-contiguous fields */
9596: PetscInt i, j, set = 1;
9597: for (i = 0; i < subSize; i += bs) {
9598: for (j = 0; j < bs; ++j) {
9599: if (subIndices[i+j] != subIndices[i]+j) {set = 0; break;}
9600: }
9601: }
9602: if (set) {ISSetBlockSize(*is, bs);}
9603: }
9604: /* Attach nullspace */
9605: for (f = 0; f < Nf; ++f) {
9606: (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[f];
9607: if ((*subdm)->nullspaceConstructors[f]) break;
9608: }
9609: if (f < Nf) {
9610: MatNullSpace nullSpace;
9611: (*(*subdm)->nullspaceConstructors[f])(*subdm, f, f, &nullSpace);
9613: PetscObjectCompose((PetscObject) *is, "nullspace", (PetscObject) nullSpace);
9614: MatNullSpaceDestroy(&nullSpace);
9615: }
9616: }
9617: return(0);
9618: }
9620: /*@
9621: DMPlexMonitorThroughput - Report the cell throughput of FE integration
9623: Input Parameter:
9624: - dm - The DM
9626: Level: developer
9628: Options Database Keys:
9629: . -dm_plex_monitor_throughput - Activate the monitor
9631: .seealso: DMSetFromOptions(), DMPlexCreate()
9632: @*/
9633: PetscErrorCode DMPlexMonitorThroughput(DM dm, void *dummy)
9634: {
9635: #if defined(PETSC_USE_LOG)
9636: PetscStageLog stageLog;
9637: PetscLogEvent event;
9638: PetscLogStage stage;
9639: PetscEventPerfInfo eventInfo;
9640: PetscReal cellRate, flopRate;
9641: PetscInt cStart, cEnd, Nf, N;
9642: const char *name;
9643: PetscErrorCode ierr;
9644: #endif
9648: #if defined(PETSC_USE_LOG)
9649: PetscObjectGetName((PetscObject) dm, &name);
9650: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
9651: DMGetNumFields(dm, &Nf);
9652: PetscLogGetStageLog(&stageLog);
9653: PetscStageLogGetCurrent(stageLog, &stage);
9654: PetscLogEventGetId("DMPlexResidualFE", &event);
9655: PetscLogEventGetPerfInfo(stage, event, &eventInfo);
9656: N = (cEnd - cStart)*Nf*eventInfo.count;
9657: flopRate = eventInfo.flops/eventInfo.time;
9658: cellRate = N/eventInfo.time;
9659: PetscPrintf(PetscObjectComm((PetscObject) dm), "DM (%s) FE Residual Integration: %D integrals %D reps\n Cell rate: %.2g/s flop rate: %.2g MF/s\n", name ? name : "unknown", N, eventInfo.count, (double) cellRate, (double) (flopRate/1.e6));
9660: #else
9661: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Plex Throughput Monitor is not supported if logging is turned off. Reconfigure using --with-log.");
9662: #endif
9663: return(0);
9664: }