Actual source code: plex.c
petsc-3.8.4 2018-03-24
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>
9: /* Logging support */
10: PetscLogEvent DMPLEX_Interpolate, PETSCPARTITIONER_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_Preallocate, DMPLEX_ResidualFEM, DMPLEX_JacobianFEM, DMPLEX_InterpolatorFEM, DMPLEX_InjectorFEM, DMPLEX_IntegralFEM, DMPLEX_CreateGmsh;
12: PETSC_EXTERN PetscErrorCode VecView_MPI(Vec, PetscViewer);
14: /*@
15: DMPlexRefineSimplexToTensor - Uniformly refines simplicial cells into tensor product cells.
16: 3 quadrilaterals per triangle in 2D and 4 hexahedra per tetrahedron in 3D.
18: Collective
20: Input Parameters:
21: . dm - The DMPlex object
23: Output Parameters:
24: . dmRefined - The refined DMPlex object
26: Note: Returns NULL if the mesh is already a tensor product mesh.
28: Level: intermediate
30: .seealso: DMPlexCreate(), DMPlexSetRefinementUniform()
31: @*/
32: PetscErrorCode DMPlexRefineSimplexToTensor(DM dm, DM *dmRefined)
33: {
34: PetscInt dim, cMax, fMax, cStart, cEnd, coneSize;
35: CellRefiner cellRefiner;
36: PetscBool lop, allnoop, localized;
37: PetscErrorCode ierr;
41: DMGetDimension(dm, &dim);
42: DMPlexGetHybridBounds(dm,&cMax,&fMax,NULL,NULL);
43: if (cMax >= 0 || fMax >= 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot handle hybrid meshes yet");
44: DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);
45: if (!(cEnd - cStart)) cellRefiner = REFINER_NOOP;
46: else {
47: DMPlexGetConeSize(dm,cStart,&coneSize);
48: switch (dim) {
49: case 1:
50: cellRefiner = REFINER_NOOP;
51: break;
52: case 2:
53: switch (coneSize) {
54: case 3:
55: cellRefiner = REFINER_SIMPLEX_TO_HEX_2D;
56: break;
57: case 4:
58: cellRefiner = REFINER_NOOP;
59: break;
60: default: SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot handle coneSize %D with dimension %D",coneSize,dim);
61: }
62: break;
63: case 3:
64: switch (coneSize) {
65: case 4:
66: cellRefiner = REFINER_SIMPLEX_TO_HEX_3D;
67: break;
68: case 6:
69: cellRefiner = REFINER_NOOP;
70: break;
71: default: SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot handle coneSize %D with dimension %D",coneSize,dim);
72: }
73: break;
74: default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot handle dimension %D",dim);
75: }
76: }
77: /* return if we don't need to refine */
78: lop = (cellRefiner == REFINER_NOOP) ? PETSC_TRUE : PETSC_FALSE;
79: MPIU_Allreduce(&lop,&allnoop,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)dm));
80: if (allnoop) {
81: *dmRefined = NULL;
82: return(0);
83: }
84: DMPlexRefineUniform_Internal(dm, cellRefiner, dmRefined);
85: DMCopyBoundary(dm, *dmRefined);
86: DMGetCoordinatesLocalized(dm, &localized);
87: if (localized) {
88: DMLocalizeCoordinates(*dmRefined);
89: }
90: return(0);
91: }
93: PetscErrorCode DMPlexGetFieldType_Internal(DM dm, PetscSection section, PetscInt field, PetscInt *sStart, PetscInt *sEnd, PetscViewerVTKFieldType *ft)
94: {
95: PetscInt dim, pStart, pEnd, vStart, vEnd, cStart, cEnd, cEndInterior, vdof = 0, cdof = 0;
99: *ft = PETSC_VTK_POINT_FIELD;
100: DMGetDimension(dm, &dim);
101: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
102: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
103: DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
104: cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
105: PetscSectionGetChart(section, &pStart, &pEnd);
106: if (field >= 0) {
107: if ((vStart >= pStart) && (vStart < pEnd)) {PetscSectionGetFieldDof(section, vStart, field, &vdof);}
108: if ((cStart >= pStart) && (cStart < pEnd)) {PetscSectionGetFieldDof(section, cStart, field, &cdof);}
109: } else {
110: if ((vStart >= pStart) && (vStart < pEnd)) {PetscSectionGetDof(section, vStart, &vdof);}
111: if ((cStart >= pStart) && (cStart < pEnd)) {PetscSectionGetDof(section, cStart, &cdof);}
112: }
113: if (vdof) {
114: *sStart = vStart;
115: *sEnd = vEnd;
116: if (vdof == dim) *ft = PETSC_VTK_POINT_VECTOR_FIELD;
117: else *ft = PETSC_VTK_POINT_FIELD;
118: } else if (cdof) {
119: *sStart = cStart;
120: *sEnd = cEnd;
121: if (cdof == dim) *ft = PETSC_VTK_CELL_VECTOR_FIELD;
122: else *ft = PETSC_VTK_CELL_FIELD;
123: } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Could not classify input Vec for VTK");
124: return(0);
125: }
127: static PetscErrorCode VecView_Plex_Local_Draw(Vec v, PetscViewer viewer)
128: {
129: DM dm;
130: PetscSection s;
131: PetscDraw draw, popup;
132: DM cdm;
133: PetscSection coordSection;
134: Vec coordinates;
135: const PetscScalar *coords, *array;
136: PetscReal bound[4] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
137: PetscReal vbound[2], time;
138: PetscBool isnull, flg;
139: PetscInt dim, Nf, f, Nc, comp, vStart, vEnd, cStart, cEnd, c, N, level, step, w = 0;
140: const char *name;
141: char title[PETSC_MAX_PATH_LEN];
142: PetscErrorCode ierr;
145: PetscViewerDrawGetDraw(viewer, 0, &draw);
146: PetscDrawIsNull(draw, &isnull);
147: if (isnull) return(0);
149: VecGetDM(v, &dm);
150: DMGetCoordinateDim(dm, &dim);
151: if (dim != 2) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Cannot draw meshes of dimension %D. Use PETSCVIEWERGLVIS", dim);
152: DMGetDefaultSection(dm, &s);
153: PetscSectionGetNumFields(s, &Nf);
154: DMGetCoarsenLevel(dm, &level);
155: DMGetCoordinateDM(dm, &cdm);
156: DMGetDefaultSection(cdm, &coordSection);
157: DMGetCoordinatesLocal(dm, &coordinates);
158: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
159: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
161: PetscObjectGetName((PetscObject) v, &name);
162: DMGetOutputSequenceNumber(dm, &step, &time);
164: VecGetLocalSize(coordinates, &N);
165: VecGetArrayRead(coordinates, &coords);
166: for (c = 0; c < N; c += dim) {
167: bound[0] = PetscMin(bound[0], PetscRealPart(coords[c])); bound[2] = PetscMax(bound[2], PetscRealPart(coords[c]));
168: bound[1] = PetscMin(bound[1], PetscRealPart(coords[c+1])); bound[3] = PetscMax(bound[3], PetscRealPart(coords[c+1]));
169: }
170: VecRestoreArrayRead(coordinates, &coords);
171: PetscDrawClear(draw);
173: /* Could implement something like DMDASelectFields() */
174: for (f = 0; f < Nf; ++f) {
175: DM fdm = dm;
176: Vec fv = v;
177: IS fis;
178: char prefix[PETSC_MAX_PATH_LEN];
179: const char *fname;
181: PetscSectionGetFieldComponents(s, f, &Nc);
182: PetscSectionGetFieldName(s, f, &fname);
184: if (v->hdr.prefix) {PetscStrcpy(prefix, v->hdr.prefix);}
185: else {prefix[0] = '\0';}
186: if (Nf > 1) {
187: DMCreateSubDM(dm, 1, &f, &fis, &fdm);
188: VecGetSubVector(v, fis, &fv);
189: PetscStrcat(prefix, fname);
190: PetscStrcat(prefix, "_");
191: }
192: for (comp = 0; comp < Nc; ++comp, ++w) {
193: PetscInt nmax = 2;
195: PetscViewerDrawGetDraw(viewer, w, &draw);
196: if (Nc > 1) {PetscSNPrintf(title, sizeof(title), "%s:%s_%D Step: %D Time: %.4g", name, fname, comp, step, time);}
197: else {PetscSNPrintf(title, sizeof(title), "%s:%s Step: %D Time: %.4g", name, fname, step, time);}
198: PetscDrawSetTitle(draw, title);
200: /* TODO Get max and min only for this component */
201: PetscOptionsGetRealArray(NULL, prefix, "-vec_view_bounds", vbound, &nmax, &flg);
202: if (!flg) {
203: VecMin(fv, NULL, &vbound[0]);
204: VecMax(fv, NULL, &vbound[1]);
205: if (vbound[1] <= vbound[0]) vbound[1] = vbound[0] + 1.0;
206: }
207: PetscDrawGetPopup(draw, &popup);
208: PetscDrawScalePopup(popup, vbound[0], vbound[1]);
209: PetscDrawSetCoordinates(draw, bound[0], bound[1], bound[2], bound[3]);
211: VecGetArrayRead(fv, &array);
212: for (c = cStart; c < cEnd; ++c) {
213: PetscScalar *coords = NULL, *a = NULL;
214: PetscInt numCoords, color[4] = {-1,-1,-1,-1};
216: DMPlexPointLocalRead(fdm, c, array, &a);
217: if (a) {
218: color[0] = PetscDrawRealToColor(PetscRealPart(a[comp]), vbound[0], vbound[1]);
219: color[1] = color[2] = color[3] = color[0];
220: } else {
221: PetscScalar *vals = NULL;
222: PetscInt numVals, va;
224: DMPlexVecGetClosure(fdm, NULL, fv, c, &numVals, &vals);
225: 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);
226: switch (numVals/Nc) {
227: case 3: /* P1 Triangle */
228: case 4: /* P1 Quadrangle */
229: for (va = 0; va < numVals/Nc; ++va) color[va] = PetscDrawRealToColor(PetscRealPart(vals[va*Nc+comp]), vbound[0], vbound[1]);
230: break;
231: case 6: /* P2 Triangle */
232: case 8: /* P2 Quadrangle */
233: for (va = 0; va < numVals/(Nc*2); ++va) color[va] = PetscDrawRealToColor(PetscRealPart(vals[va*Nc+comp + numVals/(Nc*2)]), vbound[0], vbound[1]);
234: break;
235: default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of values for cell closure %D cannot be handled", numVals/Nc);
236: }
237: DMPlexVecRestoreClosure(fdm, NULL, fv, c, &numVals, &vals);
238: }
239: DMPlexVecGetClosure(dm, coordSection, coordinates, c, &numCoords, &coords);
240: switch (numCoords) {
241: case 6:
242: 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]);
243: break;
244: case 8:
245: 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]);
246: 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]);
247: break;
248: default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot draw cells with %D coordinates", numCoords);
249: }
250: DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, &numCoords, &coords);
251: }
252: VecRestoreArrayRead(fv, &array);
253: PetscDrawFlush(draw);
254: PetscDrawPause(draw);
255: PetscDrawSave(draw);
256: }
257: if (Nf > 1) {
258: VecRestoreSubVector(v, fis, &fv);
259: ISDestroy(&fis);
260: DMDestroy(&fdm);
261: }
262: }
263: return(0);
264: }
266: PetscErrorCode VecView_Plex_Local(Vec v, PetscViewer viewer)
267: {
268: DM dm;
269: PetscBool isvtk, ishdf5, isdraw, isseq, isglvis;
273: VecGetDM(v, &dm);
274: if (!dm) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
275: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);
276: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
277: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW, &isdraw);
278: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERGLVIS, &isglvis);
279: PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);
280: if (isvtk || ishdf5 || isglvis) {
281: PetscInt numFields;
282: PetscBool fem = PETSC_FALSE;
284: DMGetNumFields(dm, &numFields);
285: if (numFields) {
286: PetscObject fe;
288: DMGetField(dm, 0, &fe);
289: if (fe->classid == PETSCFE_CLASSID) fem = PETSC_TRUE;
290: }
291: if (fem) {DMPlexInsertBoundaryValues(dm, PETSC_TRUE, v, 0.0, NULL, NULL, NULL);}
292: }
293: if (isvtk) {
294: PetscSection section;
295: PetscViewerVTKFieldType ft;
296: PetscInt pStart, pEnd;
298: DMGetDefaultSection(dm, §ion);
299: DMPlexGetFieldType_Internal(dm, section, PETSC_DETERMINE, &pStart, &pEnd, &ft);
300: PetscObjectReference((PetscObject) v); /* viewer drops reference */
301: PetscViewerVTKAddField(viewer, (PetscObject) dm, DMPlexVTKWriteAll, ft, (PetscObject) v);
302: } else if (ishdf5) {
303: #if defined(PETSC_HAVE_HDF5)
304: VecView_Plex_Local_HDF5_Internal(v, viewer);
305: #else
306: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
307: #endif
308: } else if (isglvis) {
309: VecView_GLVis(v, viewer);
310: } else if (isdraw) {
311: VecView_Plex_Local_Draw(v, viewer);
312: } else {
313: if (isseq) {VecView_Seq(v, viewer);}
314: else {VecView_MPI(v, viewer);}
315: }
316: return(0);
317: }
319: PetscErrorCode VecView_Plex(Vec v, PetscViewer viewer)
320: {
321: DM dm;
322: PetscReal time = 0.0;
323: PetscBool isvtk, ishdf5, isdraw, isseq, 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, PETSCVIEWERDRAW, &isdraw);
333: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERGLVIS, &isglvis);
334: PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);
335: if (isvtk || isglvis) {
336: PetscInt num;
337: Vec locv;
338: const char *name;
340: DMGetLocalVector(dm, &locv);
341: PetscObjectGetName((PetscObject) v, &name);
342: PetscObjectSetName((PetscObject) locv, name);
343: DMGlobalToLocalBegin(dm, v, INSERT_VALUES, locv);
344: DMGlobalToLocalEnd(dm, v, INSERT_VALUES, locv);
345: DMGetOutputSequenceNumber(dm, &num, &time);
346: DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locv, time, NULL, NULL, NULL);
347: PetscViewerGLVisSetSnapId(viewer, num);
348: VecView_Plex_Local(locv, viewer);
349: DMRestoreLocalVector(dm, &locv);
350: } else if (ishdf5) {
351: #if defined(PETSC_HAVE_HDF5)
352: VecView_Plex_HDF5_Internal(v, viewer);
353: #else
354: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
355: #endif
356: } else if (isdraw) {
357: Vec locv;
358: const char *name;
360: DMGetLocalVector(dm, &locv);
361: PetscObjectGetName((PetscObject) v, &name);
362: PetscObjectSetName((PetscObject) locv, name);
363: DMGlobalToLocalBegin(dm, v, INSERT_VALUES, locv);
364: DMGlobalToLocalEnd(dm, v, INSERT_VALUES, locv);
365: DMGetOutputSequenceNumber(dm, NULL, &time);
366: DMPlexInsertBoundaryValues(dm, PETSC_TRUE, locv, time, NULL, NULL, NULL);
367: VecView_Plex_Local(locv, viewer);
368: DMRestoreLocalVector(dm, &locv);
369: } else {
370: if (isseq) {VecView_Seq(v, viewer);}
371: else {VecView_MPI(v, viewer);}
372: }
373: return(0);
374: }
376: PetscErrorCode VecView_Plex_Native(Vec originalv, PetscViewer viewer)
377: {
378: DM dm;
379: MPI_Comm comm;
380: PetscViewerFormat format;
381: Vec v;
382: PetscBool isvtk, ishdf5;
383: PetscErrorCode ierr;
386: VecGetDM(originalv, &dm);
387: PetscObjectGetComm((PetscObject) originalv, &comm);
388: if (!dm) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
389: PetscViewerGetFormat(viewer, &format);
390: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
391: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);
392: if (format == PETSC_VIEWER_NATIVE) {
393: const char *vecname;
394: PetscInt n, nroots;
396: if (dm->sfNatural) {
397: VecGetLocalSize(originalv, &n);
398: PetscSFGetGraph(dm->sfNatural, &nroots, NULL, NULL, NULL);
399: if (n == nroots) {
400: DMGetGlobalVector(dm, &v);
401: DMPlexGlobalToNaturalBegin(dm, originalv, v);
402: DMPlexGlobalToNaturalEnd(dm, originalv, v);
403: PetscObjectGetName((PetscObject) originalv, &vecname);
404: PetscObjectSetName((PetscObject) v, vecname);
405: } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "DM global to natural SF only handles global vectors");
406: } else SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "DM global to natural SF was not created");
407: } else {
408: /* we are viewing a natural DMPlex vec. */
409: v = originalv;
410: }
411: if (ishdf5) {
412: #if defined(PETSC_HAVE_HDF5)
413: VecView_Plex_HDF5_Native_Internal(v, viewer);
414: #else
415: SETERRQ(comm, PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
416: #endif
417: } else if (isvtk) {
418: SETERRQ(comm, PETSC_ERR_SUP, "VTK format does not support viewing in natural order. Please switch to HDF5.");
419: } else {
420: PetscBool isseq;
422: PetscObjectTypeCompare((PetscObject) v, VECSEQ, &isseq);
423: if (isseq) {VecView_Seq(v, viewer);}
424: else {VecView_MPI(v, viewer);}
425: }
426: if (format == PETSC_VIEWER_NATIVE) {DMRestoreGlobalVector(dm, &v);}
427: return(0);
428: }
430: PetscErrorCode VecLoad_Plex_Local(Vec v, PetscViewer viewer)
431: {
432: DM dm;
433: PetscBool ishdf5;
437: VecGetDM(v, &dm);
438: if (!dm) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
439: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
440: if (ishdf5) {
441: DM dmBC;
442: Vec gv;
443: const char *name;
445: DMGetOutputDM(dm, &dmBC);
446: DMGetGlobalVector(dmBC, &gv);
447: PetscObjectGetName((PetscObject) v, &name);
448: PetscObjectSetName((PetscObject) gv, name);
449: VecLoad_Default(gv, viewer);
450: DMGlobalToLocalBegin(dmBC, gv, INSERT_VALUES, v);
451: DMGlobalToLocalEnd(dmBC, gv, INSERT_VALUES, v);
452: DMRestoreGlobalVector(dmBC, &gv);
453: } else {
454: VecLoad_Default(v, viewer);
455: }
456: return(0);
457: }
459: PetscErrorCode VecLoad_Plex(Vec v, PetscViewer viewer)
460: {
461: DM dm;
462: PetscBool ishdf5;
466: VecGetDM(v, &dm);
467: if (!dm) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
468: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
469: if (ishdf5) {
470: #if defined(PETSC_HAVE_HDF5)
471: VecLoad_Plex_HDF5_Internal(v, viewer);
472: #else
473: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
474: #endif
475: } else {
476: VecLoad_Default(v, viewer);
477: }
478: return(0);
479: }
481: PetscErrorCode VecLoad_Plex_Native(Vec originalv, PetscViewer viewer)
482: {
483: DM dm;
484: PetscViewerFormat format;
485: PetscBool ishdf5;
486: PetscErrorCode ierr;
489: VecGetDM(originalv, &dm);
490: if (!dm) SETERRQ(PetscObjectComm((PetscObject) originalv), PETSC_ERR_ARG_WRONG, "Vector not generated from a DM");
491: PetscViewerGetFormat(viewer, &format);
492: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
493: if (format == PETSC_VIEWER_NATIVE) {
494: if (dm->sfNatural) {
495: if (ishdf5) {
496: #if defined(PETSC_HAVE_HDF5)
497: Vec v;
498: const char *vecname;
500: DMGetGlobalVector(dm, &v);
501: PetscObjectGetName((PetscObject) originalv, &vecname);
502: PetscObjectSetName((PetscObject) v, vecname);
503: VecLoad_Plex_HDF5_Native_Internal(v, viewer);
504: DMPlexNaturalToGlobalBegin(dm, v, originalv);
505: DMPlexNaturalToGlobalEnd(dm, v, originalv);
506: DMRestoreGlobalVector(dm, &v);
507: #else
508: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
509: #endif
510: } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Reading in natural order is not supported for anything but HDF5.");
511: }
512: }
513: return(0);
514: }
516: PETSC_UNUSED static PetscErrorCode DMPlexView_Ascii_Geometry(DM dm, PetscViewer viewer)
517: {
518: PetscSection coordSection;
519: Vec coordinates;
520: DMLabel depthLabel;
521: const char *name[4];
522: const PetscScalar *a;
523: PetscInt dim, pStart, pEnd, cStart, cEnd, c;
524: PetscErrorCode ierr;
527: DMGetDimension(dm, &dim);
528: DMGetCoordinatesLocal(dm, &coordinates);
529: DMGetCoordinateSection(dm, &coordSection);
530: DMPlexGetDepthLabel(dm, &depthLabel);
531: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
532: PetscSectionGetChart(coordSection, &pStart, &pEnd);
533: VecGetArrayRead(coordinates, &a);
534: name[0] = "vertex";
535: name[1] = "edge";
536: name[dim-1] = "face";
537: name[dim] = "cell";
538: for (c = cStart; c < cEnd; ++c) {
539: PetscInt *closure = NULL;
540: PetscInt closureSize, cl;
542: PetscViewerASCIIPrintf(viewer, "Geometry for cell %D:\n", c);
543: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
544: PetscViewerASCIIPushTab(viewer);
545: for (cl = 0; cl < closureSize*2; cl += 2) {
546: PetscInt point = closure[cl], depth, dof, off, d, p;
548: if ((point < pStart) || (point >= pEnd)) continue;
549: PetscSectionGetDof(coordSection, point, &dof);
550: if (!dof) continue;
551: DMLabelGetValue(depthLabel, point, &depth);
552: PetscSectionGetOffset(coordSection, point, &off);
553: PetscViewerASCIIPrintf(viewer, "%s %D coords:", name[depth], point);
554: for (p = 0; p < dof/dim; ++p) {
555: PetscViewerASCIIPrintf(viewer, " (");
556: for (d = 0; d < dim; ++d) {
557: if (d > 0) {PetscViewerASCIIPrintf(viewer, ", ");}
558: PetscViewerASCIIPrintf(viewer, "%g", PetscRealPart(a[off+p*dim+d]));
559: }
560: PetscViewerASCIIPrintf(viewer, ")");
561: }
562: PetscViewerASCIIPrintf(viewer, "\n");
563: }
564: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
565: PetscViewerASCIIPopTab(viewer);
566: }
567: VecRestoreArrayRead(coordinates, &a);
568: return(0);
569: }
571: static PetscErrorCode DMPlexView_Ascii(DM dm, PetscViewer viewer)
572: {
573: DM_Plex *mesh = (DM_Plex*) dm->data;
574: DM cdm;
575: DMLabel markers;
576: PetscSection coordSection;
577: Vec coordinates;
578: PetscViewerFormat format;
579: PetscErrorCode ierr;
582: DMGetCoordinateDM(dm, &cdm);
583: DMGetDefaultSection(cdm, &coordSection);
584: DMGetCoordinatesLocal(dm, &coordinates);
585: PetscViewerGetFormat(viewer, &format);
586: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
587: const char *name;
588: PetscInt maxConeSize, maxSupportSize;
589: PetscInt pStart, pEnd, p;
590: PetscMPIInt rank, size;
592: MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
593: MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
594: PetscObjectGetName((PetscObject) dm, &name);
595: DMPlexGetChart(dm, &pStart, &pEnd);
596: DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);
597: PetscViewerASCIIPrintf(viewer, "Mesh '%s':\n", name);
598: PetscViewerASCIIPrintf(viewer, "orientation is missing\n", name);
599: PetscViewerASCIIPrintf(viewer, "cap --> base:\n", name);
600: PetscViewerASCIIPushSynchronized(viewer);
601: PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Max sizes cone: %D support: %D\n", rank,maxConeSize, maxSupportSize);
602: for (p = pStart; p < pEnd; ++p) {
603: PetscInt dof, off, s;
605: PetscSectionGetDof(mesh->supportSection, p, &dof);
606: PetscSectionGetOffset(mesh->supportSection, p, &off);
607: for (s = off; s < off+dof; ++s) {
608: PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D ----> %D\n", rank, p, mesh->supports[s]);
609: }
610: }
611: PetscViewerFlush(viewer);
612: PetscViewerASCIIPrintf(viewer, "base <-- cap:\n", name);
613: for (p = pStart; p < pEnd; ++p) {
614: PetscInt dof, off, c;
616: PetscSectionGetDof(mesh->coneSection, p, &dof);
617: PetscSectionGetOffset(mesh->coneSection, p, &off);
618: for (c = off; c < off+dof; ++c) {
619: PetscViewerASCIISynchronizedPrintf(viewer, "[%d]: %D <---- %D (%D)\n", rank, p, mesh->cones[c], mesh->coneOrientations[c]);
620: }
621: }
622: PetscViewerFlush(viewer);
623: PetscViewerASCIIPopSynchronized(viewer);
624: PetscSectionGetChart(coordSection, &pStart, NULL);
625: if (pStart >= 0) {PetscSectionVecView(coordSection, coordinates, viewer);}
626: DMGetLabel(dm, "marker", &markers);
627: DMLabelView(markers,viewer);
628: if (size > 1) {
629: PetscSF sf;
631: DMGetPointSF(dm, &sf);
632: PetscSFView(sf, viewer);
633: }
634: PetscViewerFlush(viewer);
635: } else if (format == PETSC_VIEWER_ASCII_LATEX) {
636: const char *name, *color;
637: const char *defcolors[3] = {"gray", "orange", "green"};
638: const char *deflcolors[4] = {"blue", "cyan", "red", "magenta"};
639: PetscReal scale = 2.0;
640: PetscBool useNumbers = PETSC_TRUE, useLabels, useColors;
641: double tcoords[3];
642: PetscScalar *coords;
643: PetscInt numLabels, l, numColors, numLColors, dim, depth, cStart, cEnd, c, vStart, vEnd, v, eStart = 0, eEnd = 0, e, p;
644: PetscMPIInt rank, size;
645: char **names, **colors, **lcolors;
647: DMGetDimension(dm, &dim);
648: DMPlexGetDepth(dm, &depth);
649: DMGetNumLabels(dm, &numLabels);
650: numLabels = PetscMax(numLabels, 10);
651: numColors = 10;
652: numLColors = 10;
653: PetscCalloc3(numLabels, &names, numColors, &colors, numLColors, &lcolors);
654: PetscOptionsGetReal(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_scale", &scale, NULL);
655: PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_numbers", &useNumbers, NULL);
656: PetscOptionsGetStringArray(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_labels", names, &numLabels, &useLabels);
657: if (!useLabels) numLabels = 0;
658: PetscOptionsGetStringArray(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_colors", colors, &numColors, &useColors);
659: if (!useColors) {
660: numColors = 3;
661: for (c = 0; c < numColors; ++c) {PetscStrallocpy(defcolors[c], &colors[c]);}
662: }
663: PetscOptionsGetStringArray(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_view_lcolors", lcolors, &numLColors, &useColors);
664: if (!useColors) {
665: numLColors = 4;
666: for (c = 0; c < numLColors; ++c) {PetscStrallocpy(deflcolors[c], &lcolors[c]);}
667: }
668: MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
669: MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
670: PetscObjectGetName((PetscObject) dm, &name);
671: PetscViewerASCIIPrintf(viewer, "\
672: \\documentclass[tikz]{standalone}\n\n\
673: \\usepackage{pgflibraryshapes}\n\
674: \\usetikzlibrary{backgrounds}\n\
675: \\usetikzlibrary{arrows}\n\
676: \\begin{document}\n");
677: if (size > 1) {
678: PetscViewerASCIIPrintf(viewer, "%s for process ", name);
679: for (p = 0; p < size; ++p) {
680: if (p > 0 && p == size-1) {
681: PetscViewerASCIIPrintf(viewer, ", and ", colors[p%numColors], p);
682: } else if (p > 0) {
683: PetscViewerASCIIPrintf(viewer, ", ", colors[p%numColors], p);
684: }
685: PetscViewerASCIIPrintf(viewer, "{\\textcolor{%s}%D}", colors[p%numColors], p);
686: }
687: PetscViewerASCIIPrintf(viewer, ".\n\n\n");
688: }
689: PetscViewerASCIIPrintf(viewer, "\\begin{tikzpicture}[scale = %g,font=\\fontsize{8}{8}\\selectfont]\n", 1.0);
690: /* Plot vertices */
691: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
692: VecGetArray(coordinates, &coords);
693: PetscViewerASCIIPushSynchronized(viewer);
694: for (v = vStart; v < vEnd; ++v) {
695: PetscInt off, dof, d;
696: PetscBool isLabeled = PETSC_FALSE;
698: PetscSectionGetDof(coordSection, v, &dof);
699: PetscSectionGetOffset(coordSection, v, &off);
700: PetscViewerASCIISynchronizedPrintf(viewer, "\\path (");
701: if (PetscUnlikely(dof > 3)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"coordSection vertex %D has dof %D > 3",v,dof);
702: for (d = 0; d < dof; ++d) {
703: tcoords[d] = (double) (scale*PetscRealPart(coords[off+d]));
704: tcoords[d] = PetscAbsReal(tcoords[d]) < 1e-10 ? 0.0 : tcoords[d];
705: }
706: /* Rotate coordinates since PGF makes z point out of the page instead of up */
707: if (dim == 3) {PetscReal tmp = tcoords[1]; tcoords[1] = tcoords[2]; tcoords[2] = -tmp;}
708: for (d = 0; d < dof; ++d) {
709: if (d > 0) {PetscViewerASCIISynchronizedPrintf(viewer, ",");}
710: PetscViewerASCIISynchronizedPrintf(viewer, "%g", tcoords[d]);
711: }
712: color = colors[rank%numColors];
713: for (l = 0; l < numLabels; ++l) {
714: PetscInt val;
715: DMGetLabelValue(dm, names[l], v, &val);
716: if (val >= 0) {color = lcolors[l%numLColors]; isLabeled = PETSC_TRUE; break;}
717: }
718: if (useNumbers) {
719: PetscViewerASCIISynchronizedPrintf(viewer, ") node(%D_%d) [draw,shape=circle,color=%s] {%D};\n", v, rank, color, v);
720: } else {
721: PetscViewerASCIISynchronizedPrintf(viewer, ") node(%D_%d) [fill,inner sep=%dpt,shape=circle,color=%s] {};\n", v, rank, !isLabeled ? 1 : 2, color);
722: }
723: }
724: VecRestoreArray(coordinates, &coords);
725: PetscViewerFlush(viewer);
726: /* Plot edges */
727: if (depth > 1) {DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);}
728: if (dim < 3 && useNumbers) {
729: VecGetArray(coordinates, &coords);
730: PetscViewerASCIIPrintf(viewer, "\\path\n");
731: for (e = eStart; e < eEnd; ++e) {
732: const PetscInt *cone;
733: PetscInt coneSize, offA, offB, dof, d;
735: DMPlexGetConeSize(dm, e, &coneSize);
736: if (coneSize != 2) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Edge %D cone should have two vertices, not %D", e, coneSize);
737: DMPlexGetCone(dm, e, &cone);
738: PetscSectionGetDof(coordSection, cone[0], &dof);
739: PetscSectionGetOffset(coordSection, cone[0], &offA);
740: PetscSectionGetOffset(coordSection, cone[1], &offB);
741: PetscViewerASCIISynchronizedPrintf(viewer, "(");
742: for (d = 0; d < dof; ++d) {
743: tcoords[d] = (double) (scale*PetscRealPart(coords[offA+d]+coords[offB+d]));
744: tcoords[d] = PetscAbsReal(tcoords[d]) < 1e-10 ? 0.0 : tcoords[d];
745: }
746: /* Rotate coordinates since PGF makes z point out of the page instead of up */
747: if (dim == 3) {PetscReal tmp = tcoords[1]; tcoords[1] = tcoords[2]; tcoords[2] = -tmp;}
748: for (d = 0; d < dof; ++d) {
749: if (d > 0) {PetscViewerASCIISynchronizedPrintf(viewer, ",");}
750: PetscViewerASCIISynchronizedPrintf(viewer, "%g", (double)tcoords[d]);
751: }
752: color = colors[rank%numColors];
753: for (l = 0; l < numLabels; ++l) {
754: PetscInt val;
755: DMGetLabelValue(dm, names[l], v, &val);
756: if (val >= 0) {color = lcolors[l%numLColors]; break;}
757: }
758: PetscViewerASCIISynchronizedPrintf(viewer, ") node(%D_%d) [draw,shape=circle,color=%s] {%D} --\n", e, rank, color, e);
759: }
760: VecRestoreArray(coordinates, &coords);
761: PetscViewerFlush(viewer);
762: PetscViewerASCIIPrintf(viewer, "(0,0);\n");
763: }
764: /* Plot cells */
765: if (dim == 3 || !useNumbers) {
766: for (e = eStart; e < eEnd; ++e) {
767: const PetscInt *cone;
769: color = colors[rank%numColors];
770: for (l = 0; l < numLabels; ++l) {
771: PetscInt val;
772: DMGetLabelValue(dm, names[l], e, &val);
773: if (val >= 0) {color = lcolors[l%numLColors]; break;}
774: }
775: DMPlexGetCone(dm, e, &cone);
776: PetscViewerASCIISynchronizedPrintf(viewer, "\\draw[color=%s] (%D_%d) -- (%D_%d);\n", color, cone[0], rank, cone[1], rank);
777: }
778: } else {
779: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
780: for (c = cStart; c < cEnd; ++c) {
781: PetscInt *closure = NULL;
782: PetscInt closureSize, firstPoint = -1;
784: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
785: PetscViewerASCIISynchronizedPrintf(viewer, "\\draw[color=%s] ", colors[rank%numColors]);
786: for (p = 0; p < closureSize*2; p += 2) {
787: const PetscInt point = closure[p];
789: if ((point < vStart) || (point >= vEnd)) continue;
790: if (firstPoint >= 0) {PetscViewerASCIISynchronizedPrintf(viewer, " -- ");}
791: PetscViewerASCIISynchronizedPrintf(viewer, "(%D_%d)", point, rank);
792: if (firstPoint < 0) firstPoint = point;
793: }
794: /* Why doesn't this work? PetscViewerASCIISynchronizedPrintf(viewer, " -- cycle;\n"); */
795: PetscViewerASCIISynchronizedPrintf(viewer, " -- (%D_%d);\n", firstPoint, rank);
796: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
797: }
798: }
799: PetscViewerFlush(viewer);
800: PetscViewerASCIIPopSynchronized(viewer);
801: PetscViewerASCIIPrintf(viewer, "\\end{tikzpicture}\n");
802: PetscViewerASCIIPrintf(viewer, "\\end{document}\n", name);
803: for (l = 0; l < numLabels; ++l) {PetscFree(names[l]);}
804: for (c = 0; c < numColors; ++c) {PetscFree(colors[c]);}
805: for (c = 0; c < numLColors; ++c) {PetscFree(lcolors[c]);}
806: PetscFree3(names, colors, lcolors);
807: } else {
808: MPI_Comm comm;
809: PetscInt *sizes, *hybsizes;
810: PetscInt locDepth, depth, dim, d, pMax[4];
811: PetscInt pStart, pEnd, p;
812: PetscInt numLabels, l;
813: const char *name;
814: PetscMPIInt size;
816: PetscObjectGetComm((PetscObject)dm,&comm);
817: MPI_Comm_size(comm, &size);
818: DMGetDimension(dm, &dim);
819: PetscObjectGetName((PetscObject) dm, &name);
820: if (name) {PetscViewerASCIIPrintf(viewer, "%s in %D dimensions:\n", name, dim);}
821: else {PetscViewerASCIIPrintf(viewer, "Mesh in %D dimensions:\n", dim);}
822: DMPlexGetDepth(dm, &locDepth);
823: MPIU_Allreduce(&locDepth, &depth, 1, MPIU_INT, MPI_MAX, comm);
824: DMPlexGetHybridBounds(dm, &pMax[depth], depth > 0 ? &pMax[depth-1] : NULL, &pMax[1], &pMax[0]);
825: PetscMalloc2(size,&sizes,size,&hybsizes);
826: if (depth == 1) {
827: DMPlexGetDepthStratum(dm, 0, &pStart, &pEnd);
828: pEnd = pEnd - pStart;
829: MPI_Gather(&pEnd, 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm);
830: PetscViewerASCIIPrintf(viewer, " %d-cells:", 0);
831: for (p = 0; p < size; ++p) {PetscViewerASCIIPrintf(viewer, " %D", sizes[p]);}
832: PetscViewerASCIIPrintf(viewer, "\n");
833: DMPlexGetHeightStratum(dm, 0, &pStart, &pEnd);
834: pEnd = pEnd - pStart;
835: MPI_Gather(&pEnd, 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm);
836: PetscViewerASCIIPrintf(viewer, " %D-cells:", dim);
837: for (p = 0; p < size; ++p) {PetscViewerASCIIPrintf(viewer, " %D", sizes[p]);}
838: PetscViewerASCIIPrintf(viewer, "\n");
839: } else {
840: PetscMPIInt rank;
841: MPI_Comm_rank(comm, &rank);
842: for (d = 0; d <= dim; d++) {
843: DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
844: pEnd -= pStart;
845: pMax[d] -= pStart;
846: MPI_Gather(&pEnd, 1, MPIU_INT, sizes, 1, MPIU_INT, 0, comm);
847: MPI_Gather(&pMax[d], 1, MPIU_INT, hybsizes, 1, MPIU_INT, 0, comm);
848: PetscViewerASCIIPrintf(viewer, " %D-cells:", d);
849: for (p = 0; p < size; ++p) {
850: if (!rank) {
851: if (hybsizes[p] >= 0) {PetscViewerASCIIPrintf(viewer, " %D (%D)", sizes[p], sizes[p] - hybsizes[p]);}
852: else {PetscViewerASCIIPrintf(viewer, " %D", sizes[p]);}
853: }
854: }
855: PetscViewerASCIIPrintf(viewer, "\n");
856: }
857: }
858: PetscFree2(sizes,hybsizes);
859: DMGetNumLabels(dm, &numLabels);
860: if (numLabels) {PetscViewerASCIIPrintf(viewer, "Labels:\n");}
861: for (l = 0; l < numLabels; ++l) {
862: DMLabel label;
863: const char *name;
864: IS valueIS;
865: const PetscInt *values;
866: PetscInt numValues, v;
868: DMGetLabelName(dm, l, &name);
869: DMGetLabel(dm, name, &label);
870: DMLabelGetNumValues(label, &numValues);
871: PetscViewerASCIIPrintf(viewer, " %s: %D strata with value/size (", name, numValues);
872: DMLabelGetValueIS(label, &valueIS);
873: ISGetIndices(valueIS, &values);
874: PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);
875: for (v = 0; v < numValues; ++v) {
876: PetscInt size;
878: DMLabelGetStratumSize(label, values[v], &size);
879: if (v > 0) {PetscViewerASCIIPrintf(viewer, ", ");}
880: PetscViewerASCIIPrintf(viewer, "%D (%D)", values[v], size);
881: }
882: PetscViewerASCIIPrintf(viewer, ")\n");
883: PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);
884: ISRestoreIndices(valueIS, &values);
885: ISDestroy(&valueIS);
886: }
887: DMGetCoarseDM(dm, &cdm);
888: if (cdm) {
889: PetscViewerASCIIPushTab(viewer);
890: DMPlexView_Ascii(cdm, viewer);
891: PetscViewerASCIIPopTab(viewer);
892: }
893: }
894: return(0);
895: }
897: static PetscErrorCode DMPlexView_Draw(DM dm, PetscViewer viewer)
898: {
899: PetscDraw draw;
900: DM cdm;
901: PetscSection coordSection;
902: Vec coordinates;
903: const PetscScalar *coords;
904: PetscReal xyl[2],xyr[2],bound[4] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
905: PetscBool isnull;
906: PetscInt dim, vStart, vEnd, cStart, cEnd, c, N;
907: PetscErrorCode ierr;
910: DMGetCoordinateDim(dm, &dim);
911: if (dim != 2) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Cannot draw meshes of dimension %D", dim);
912: DMGetCoordinateDM(dm, &cdm);
913: DMGetDefaultSection(cdm, &coordSection);
914: DMGetCoordinatesLocal(dm, &coordinates);
915: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
916: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
918: PetscViewerDrawGetDraw(viewer, 0, &draw);
919: PetscDrawIsNull(draw, &isnull);
920: if (isnull) return(0);
921: PetscDrawSetTitle(draw, "Mesh");
923: VecGetLocalSize(coordinates, &N);
924: VecGetArrayRead(coordinates, &coords);
925: for (c = 0; c < N; c += dim) {
926: bound[0] = PetscMin(bound[0], PetscRealPart(coords[c])); bound[2] = PetscMax(bound[2], PetscRealPart(coords[c]));
927: bound[1] = PetscMin(bound[1], PetscRealPart(coords[c+1])); bound[3] = PetscMax(bound[3], PetscRealPart(coords[c+1]));
928: }
929: VecRestoreArrayRead(coordinates, &coords);
930: MPIU_Allreduce(&bound[0],xyl,2,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)dm));
931: MPIU_Allreduce(&bound[2],xyr,2,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)dm));
932: PetscDrawSetCoordinates(draw, xyl[0], xyl[1], xyr[0], xyr[1]);
933: PetscDrawClear(draw);
935: for (c = cStart; c < cEnd; ++c) {
936: PetscScalar *coords = NULL;
937: PetscInt numCoords,coneSize;
939: DMPlexGetConeSize(dm, c, &coneSize);
940: DMPlexVecGetClosure(dm, coordSection, coordinates, c, &numCoords, &coords);
941: switch (coneSize) {
942: case 3:
943: PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PETSC_DRAW_BLACK);
944: PetscDrawLine(draw, PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_BLACK);
945: PetscDrawLine(draw, PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[0]), PetscRealPart(coords[1]), PETSC_DRAW_BLACK);
946: break;
947: case 4:
948: PetscDrawLine(draw, PetscRealPart(coords[0]), PetscRealPart(coords[1]), PetscRealPart(coords[2]), PetscRealPart(coords[3]), PETSC_DRAW_BLACK);
949: PetscDrawLine(draw, PetscRealPart(coords[2]), PetscRealPart(coords[3]), PetscRealPart(coords[4]), PetscRealPart(coords[5]), PETSC_DRAW_BLACK);
950: PetscDrawLine(draw, PetscRealPart(coords[4]), PetscRealPart(coords[5]), PetscRealPart(coords[6]), PetscRealPart(coords[7]), PETSC_DRAW_BLACK);
951: PetscDrawLine(draw, PetscRealPart(coords[6]), PetscRealPart(coords[7]), PetscRealPart(coords[0]), PetscRealPart(coords[1]), PETSC_DRAW_BLACK);
952: break;
953: default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot draw cells with %D facets", coneSize);
954: }
955: DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, &numCoords, &coords);
956: }
957: PetscDrawFlush(draw);
958: PetscDrawPause(draw);
959: PetscDrawSave(draw);
960: return(0);
961: }
963: PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer)
964: {
965: PetscBool iascii, ishdf5, isvtk, isdraw, flg, isglvis;
966: PetscErrorCode ierr;
971: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
972: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERVTK, &isvtk);
973: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
974: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW, &isdraw);
975: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERGLVIS, &isglvis);
976: if (iascii) {
977: PetscViewerFormat format;
978: PetscViewerGetFormat(viewer, &format);
979: if (format == PETSC_VIEWER_ASCII_GLVIS) {
980: DMPlexView_GLVis(dm, viewer);
981: } else {
982: DMPlexView_Ascii(dm, viewer);
983: }
984: } else if (ishdf5) {
985: #if defined(PETSC_HAVE_HDF5)
986: PetscViewerPushFormat(viewer, PETSC_VIEWER_HDF5_VIZ);
987: DMPlexView_HDF5_Internal(dm, viewer);
988: PetscViewerPopFormat(viewer);
989: #else
990: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
991: #endif
992: } else if (isvtk) {
993: DMPlexVTKWriteAll((PetscObject) dm,viewer);
994: } else if (isdraw) {
995: DMPlexView_Draw(dm, viewer);
996: } else if (isglvis) {
997: DMPlexView_GLVis(dm, viewer);
998: }
999: /* Optionally view the partition */
1000: PetscOptionsHasName(((PetscObject) dm)->options, ((PetscObject) dm)->prefix, "-dm_partition_view", &flg);
1001: if (flg) {
1002: Vec ranks;
1003: DMPlexCreateRankField(dm, &ranks);
1004: VecView(ranks, viewer);
1005: VecDestroy(&ranks);
1006: }
1007: return(0);
1008: }
1010: PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer)
1011: {
1012: PetscBool isbinary, ishdf5;
1018: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERBINARY, &isbinary);
1019: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);
1020: if (isbinary) {SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Do not yet support binary viewers");}
1021: else if (ishdf5) {
1022: #if defined(PETSC_HAVE_HDF5)
1023: DMPlexLoad_HDF5_Internal(dm, viewer);
1024: #else
1025: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5");
1026: #endif
1027: }
1028: return(0);
1029: }
1031: PetscErrorCode DMDestroy_Plex(DM dm)
1032: {
1033: DM_Plex *mesh = (DM_Plex*) dm->data;
1037: PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",NULL);
1038: PetscObjectComposeFunction((PetscObject)dm,"DMPlexInsertBoundaryValues_C", NULL);
1039: if (--mesh->refct > 0) return(0);
1040: PetscSectionDestroy(&mesh->coneSection);
1041: PetscFree(mesh->cones);
1042: PetscFree(mesh->coneOrientations);
1043: PetscSectionDestroy(&mesh->supportSection);
1044: PetscSectionDestroy(&mesh->subdomainSection);
1045: PetscFree(mesh->supports);
1046: PetscFree(mesh->facesTmp);
1047: PetscFree(mesh->tetgenOpts);
1048: PetscFree(mesh->triangleOpts);
1049: PetscPartitionerDestroy(&mesh->partitioner);
1050: DMLabelDestroy(&mesh->subpointMap);
1051: ISDestroy(&mesh->globalVertexNumbers);
1052: ISDestroy(&mesh->globalCellNumbers);
1053: PetscSectionDestroy(&mesh->anchorSection);
1054: ISDestroy(&mesh->anchorIS);
1055: PetscSectionDestroy(&mesh->parentSection);
1056: PetscFree(mesh->parents);
1057: PetscFree(mesh->childIDs);
1058: PetscSectionDestroy(&mesh->childSection);
1059: PetscFree(mesh->children);
1060: DMDestroy(&mesh->referenceTree);
1061: PetscGridHashDestroy(&mesh->lbox);
1062: /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
1063: PetscFree(mesh);
1064: return(0);
1065: }
1067: PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J)
1068: {
1069: PetscSection sectionGlobal;
1070: PetscInt bs = -1, mbs;
1071: PetscInt localSize;
1072: PetscBool isShell, isBlock, isSeqBlock, isMPIBlock, isSymBlock, isSymSeqBlock, isSymMPIBlock, isMatIS;
1073: PetscErrorCode ierr;
1074: MatType mtype;
1075: ISLocalToGlobalMapping ltog;
1078: MatInitializePackage();
1079: mtype = dm->mattype;
1080: DMGetDefaultGlobalSection(dm, §ionGlobal);
1081: /* PetscSectionGetStorageSize(sectionGlobal, &localSize); */
1082: PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize);
1083: MatCreate(PetscObjectComm((PetscObject)dm), J);
1084: MatSetSizes(*J, localSize, localSize, PETSC_DETERMINE, PETSC_DETERMINE);
1085: MatSetType(*J, mtype);
1086: MatSetFromOptions(*J);
1087: MatGetBlockSize(*J, &mbs);
1088: if (mbs > 1) bs = mbs;
1089: PetscStrcmp(mtype, MATSHELL, &isShell);
1090: PetscStrcmp(mtype, MATBAIJ, &isBlock);
1091: PetscStrcmp(mtype, MATSEQBAIJ, &isSeqBlock);
1092: PetscStrcmp(mtype, MATMPIBAIJ, &isMPIBlock);
1093: PetscStrcmp(mtype, MATSBAIJ, &isSymBlock);
1094: PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock);
1095: PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock);
1096: PetscStrcmp(mtype, MATIS, &isMatIS);
1097: if (!isShell) {
1098: PetscSection subSection;
1099: PetscBool fillMatrix = (PetscBool)(!dm->prealloc_only && !isMatIS);
1100: PetscInt *dnz, *onz, *dnzu, *onzu, bsLocal, bsMax, bsMin, *ltogidx, lsize;
1101: PetscInt pStart, pEnd, p, dof, cdof;
1103: /* Set localtoglobalmapping on the matrix for MatSetValuesLocal() to work (it also creates the local matrices in case of MATIS) */
1104: if (isMatIS) { /* need a different l2g map than the one computed by DMGetLocalToGlobalMapping */
1105: PetscSection section;
1106: PetscInt size;
1108: DMGetDefaultSection(dm, §ion);
1109: PetscSectionGetStorageSize(section, &size);
1110: PetscMalloc1(size,<ogidx);
1111: DMPlexGetSubdomainSection(dm, &subSection);
1112: } else {
1113: DMGetLocalToGlobalMapping(dm,<og);
1114: }
1115: PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);
1116: for (p = pStart, lsize = 0; p < pEnd; ++p) {
1117: PetscInt bdof;
1119: PetscSectionGetDof(sectionGlobal, p, &dof);
1120: PetscSectionGetConstraintDof(sectionGlobal, p, &cdof);
1121: dof = dof < 0 ? -(dof+1) : dof;
1122: bdof = cdof && (dof-cdof) ? 1 : dof;
1123: if (dof) {
1124: if (bs < 0) {bs = bdof;}
1125: else if (bs != bdof) {bs = 1; if (!isMatIS) break;}
1126: }
1127: if (isMatIS) {
1128: PetscInt loff,c,off;
1129: PetscSectionGetOffset(subSection, p, &loff);
1130: PetscSectionGetOffset(sectionGlobal, p, &off);
1131: for (c = 0; c < dof-cdof; ++c, ++lsize) ltogidx[loff+c] = off > -1 ? off+c : -(off+1)+c;
1132: }
1133: }
1134: /* Must have same blocksize on all procs (some might have no points) */
1135: bsLocal = bs;
1136: MPIU_Allreduce(&bsLocal, &bsMax, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));
1137: bsLocal = bs < 0 ? bsMax : bs;
1138: MPIU_Allreduce(&bsLocal, &bsMin, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));
1139: if (bsMin != bsMax) {bs = 1;}
1140: else {bs = bsMax;}
1141: bs = bs < 0 ? 1 : bs;
1142: if (isMatIS) {
1143: PetscInt l;
1144: /* Must reduce indices by blocksize */
1145: if (bs > 1) for (l = 0; l < lsize; ++l) ltogidx[l] /= bs;
1146: ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), bs, lsize, ltogidx, PETSC_OWN_POINTER, <og);
1147: }
1148: MatSetLocalToGlobalMapping(*J,ltog,ltog);
1149: if (isMatIS) {
1150: ISLocalToGlobalMappingDestroy(<og);
1151: }
1152: PetscCalloc4(localSize/bs, &dnz, localSize/bs, &onz, localSize/bs, &dnzu, localSize/bs, &onzu);
1153: DMPlexPreallocateOperator(dm, bs, dnz, onz, dnzu, onzu, *J, fillMatrix);
1154: PetscFree4(dnz, onz, dnzu, onzu);
1155: }
1156: MatSetDM(*J, dm);
1157: return(0);
1158: }
1160: /*@
1161: DMPlexGetSubdomainSection - Returns the section associated with the subdomain
1163: Not collective
1165: Input Parameter:
1166: . mesh - The DMPlex
1168: Output Parameters:
1169: . subsection - The subdomain section
1171: Level: developer
1173: .seealso:
1174: @*/
1175: PetscErrorCode DMPlexGetSubdomainSection(DM dm, PetscSection *subsection)
1176: {
1177: DM_Plex *mesh = (DM_Plex*) dm->data;
1182: if (!mesh->subdomainSection) {
1183: PetscSection section;
1184: PetscSF sf;
1186: PetscSFCreate(PETSC_COMM_SELF,&sf);
1187: DMGetDefaultSection(dm,§ion);
1188: PetscSectionCreateGlobalSection(section,sf,PETSC_FALSE,PETSC_TRUE,&mesh->subdomainSection);
1189: PetscSFDestroy(&sf);
1190: }
1191: *subsection = mesh->subdomainSection;
1192: return(0);
1193: }
1195: /*@
1196: DMPlexGetChart - Return the interval for all mesh points [pStart, pEnd)
1198: Not collective
1200: Input Parameter:
1201: . mesh - The DMPlex
1203: Output Parameters:
1204: + pStart - The first mesh point
1205: - pEnd - The upper bound for mesh points
1207: Level: beginner
1209: .seealso: DMPlexCreate(), DMPlexSetChart()
1210: @*/
1211: PetscErrorCode DMPlexGetChart(DM dm, PetscInt *pStart, PetscInt *pEnd)
1212: {
1213: DM_Plex *mesh = (DM_Plex*) dm->data;
1218: PetscSectionGetChart(mesh->coneSection, pStart, pEnd);
1219: return(0);
1220: }
1222: /*@
1223: DMPlexSetChart - Set the interval for all mesh points [pStart, pEnd)
1225: Not collective
1227: Input Parameters:
1228: + mesh - The DMPlex
1229: . pStart - The first mesh point
1230: - pEnd - The upper bound for mesh points
1232: Output Parameters:
1234: Level: beginner
1236: .seealso: DMPlexCreate(), DMPlexGetChart()
1237: @*/
1238: PetscErrorCode DMPlexSetChart(DM dm, PetscInt pStart, PetscInt pEnd)
1239: {
1240: DM_Plex *mesh = (DM_Plex*) dm->data;
1245: PetscSectionSetChart(mesh->coneSection, pStart, pEnd);
1246: PetscSectionSetChart(mesh->supportSection, pStart, pEnd);
1247: return(0);
1248: }
1250: /*@
1251: DMPlexGetConeSize - Return the number of in-edges for this point in the DAG
1253: Not collective
1255: Input Parameters:
1256: + mesh - The DMPlex
1257: - p - The point, which must lie in the chart set with DMPlexSetChart()
1259: Output Parameter:
1260: . size - The cone size for point p
1262: Level: beginner
1264: .seealso: DMPlexCreate(), DMPlexSetConeSize(), DMPlexSetChart()
1265: @*/
1266: PetscErrorCode DMPlexGetConeSize(DM dm, PetscInt p, PetscInt *size)
1267: {
1268: DM_Plex *mesh = (DM_Plex*) dm->data;
1274: PetscSectionGetDof(mesh->coneSection, p, size);
1275: return(0);
1276: }
1278: /*@
1279: DMPlexSetConeSize - Set the number of in-edges for this point in the DAG
1281: Not collective
1283: Input Parameters:
1284: + mesh - The DMPlex
1285: . p - The point, which must lie in the chart set with DMPlexSetChart()
1286: - size - The cone size for point p
1288: Output Parameter:
1290: Note:
1291: This should be called after DMPlexSetChart().
1293: Level: beginner
1295: .seealso: DMPlexCreate(), DMPlexGetConeSize(), DMPlexSetChart()
1296: @*/
1297: PetscErrorCode DMPlexSetConeSize(DM dm, PetscInt p, PetscInt size)
1298: {
1299: DM_Plex *mesh = (DM_Plex*) dm->data;
1304: PetscSectionSetDof(mesh->coneSection, p, size);
1306: mesh->maxConeSize = PetscMax(mesh->maxConeSize, size);
1307: return(0);
1308: }
1310: /*@
1311: DMPlexAddConeSize - Add the given number of in-edges to this point in the DAG
1313: Not collective
1315: Input Parameters:
1316: + mesh - The DMPlex
1317: . p - The point, which must lie in the chart set with DMPlexSetChart()
1318: - size - The additional cone size for point p
1320: Output Parameter:
1322: Note:
1323: This should be called after DMPlexSetChart().
1325: Level: beginner
1327: .seealso: DMPlexCreate(), DMPlexSetConeSize(), DMPlexGetConeSize(), DMPlexSetChart()
1328: @*/
1329: PetscErrorCode DMPlexAddConeSize(DM dm, PetscInt p, PetscInt size)
1330: {
1331: DM_Plex *mesh = (DM_Plex*) dm->data;
1332: PetscInt csize;
1337: PetscSectionAddDof(mesh->coneSection, p, size);
1338: PetscSectionGetDof(mesh->coneSection, p, &csize);
1340: mesh->maxConeSize = PetscMax(mesh->maxConeSize, csize);
1341: return(0);
1342: }
1344: /*@C
1345: DMPlexGetCone - Return the points on the in-edges for this point in the DAG
1347: Not collective
1349: Input Parameters:
1350: + mesh - The DMPlex
1351: - p - The point, which must lie in the chart set with DMPlexSetChart()
1353: Output Parameter:
1354: . cone - An array of points which are on the in-edges for point p
1356: Level: beginner
1358: Fortran Notes:
1359: Since it returns an array, this routine is only available in Fortran 90, and you must
1360: include petsc.h90 in your code.
1362: You must also call DMPlexRestoreCone() after you finish using the returned array.
1364: .seealso: DMPlexCreate(), DMPlexSetCone(), DMPlexSetChart()
1365: @*/
1366: PetscErrorCode DMPlexGetCone(DM dm, PetscInt p, const PetscInt *cone[])
1367: {
1368: DM_Plex *mesh = (DM_Plex*) dm->data;
1369: PetscInt off;
1375: PetscSectionGetOffset(mesh->coneSection, p, &off);
1376: *cone = &mesh->cones[off];
1377: return(0);
1378: }
1380: /*@
1381: DMPlexSetCone - Set the points on the in-edges for this point in the DAG
1383: Not collective
1385: Input Parameters:
1386: + mesh - The DMPlex
1387: . p - The point, which must lie in the chart set with DMPlexSetChart()
1388: - cone - An array of points which are on the in-edges for point p
1390: Output Parameter:
1392: Note:
1393: This should be called after all calls to DMPlexSetConeSize() and DMSetUp().
1395: Level: beginner
1397: .seealso: DMPlexCreate(), DMPlexGetCone(), DMPlexSetChart(), DMPlexSetConeSize(), DMSetUp()
1398: @*/
1399: PetscErrorCode DMPlexSetCone(DM dm, PetscInt p, const PetscInt cone[])
1400: {
1401: DM_Plex *mesh = (DM_Plex*) dm->data;
1402: PetscInt pStart, pEnd;
1403: PetscInt dof, off, c;
1408: PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd);
1409: PetscSectionGetDof(mesh->coneSection, p, &dof);
1411: PetscSectionGetOffset(mesh->coneSection, p, &off);
1412: 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);
1413: for (c = 0; c < dof; ++c) {
1414: 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);
1415: mesh->cones[off+c] = cone[c];
1416: }
1417: return(0);
1418: }
1420: /*@C
1421: DMPlexGetConeOrientation - Return the orientations on the in-edges for this point in the DAG
1423: Not collective
1425: Input Parameters:
1426: + mesh - The DMPlex
1427: - p - The point, which must lie in the chart set with DMPlexSetChart()
1429: Output Parameter:
1430: . coneOrientation - An array of orientations which are on the in-edges for point p. An orientation is an
1431: integer giving the prescription for cone traversal. If it is negative, the cone is
1432: traversed in the opposite direction. Its value 'o', or if negative '-(o+1)', gives
1433: the index of the cone point on which to start.
1435: Level: beginner
1437: Fortran Notes:
1438: Since it returns an array, this routine is only available in Fortran 90, and you must
1439: include petsc.h90 in your code.
1441: You must also call DMPlexRestoreConeOrientation() after you finish using the returned array.
1443: .seealso: DMPlexCreate(), DMPlexGetCone(), DMPlexSetCone(), DMPlexSetChart()
1444: @*/
1445: PetscErrorCode DMPlexGetConeOrientation(DM dm, PetscInt p, const PetscInt *coneOrientation[])
1446: {
1447: DM_Plex *mesh = (DM_Plex*) dm->data;
1448: PetscInt off;
1453: #if defined(PETSC_USE_DEBUG)
1454: {
1455: PetscInt dof;
1456: PetscSectionGetDof(mesh->coneSection, p, &dof);
1458: }
1459: #endif
1460: PetscSectionGetOffset(mesh->coneSection, p, &off);
1462: *coneOrientation = &mesh->coneOrientations[off];
1463: return(0);
1464: }
1466: /*@
1467: DMPlexSetConeOrientation - Set the orientations on the in-edges for this point in the DAG
1469: Not collective
1471: Input Parameters:
1472: + mesh - The DMPlex
1473: . p - The point, which must lie in the chart set with DMPlexSetChart()
1474: - coneOrientation - An array of orientations which are on the in-edges for point p. An orientation is an
1475: integer giving the prescription for cone traversal. If it is negative, the cone is
1476: traversed in the opposite direction. Its value 'o', or if negative '-(o+1)', gives
1477: the index of the cone point on which to start.
1479: Output Parameter:
1481: Note:
1482: This should be called after all calls to DMPlexSetConeSize() and DMSetUp().
1484: Level: beginner
1486: .seealso: DMPlexCreate(), DMPlexGetConeOrientation(), DMPlexSetCone(), DMPlexSetChart(), DMPlexSetConeSize(), DMSetUp()
1487: @*/
1488: PetscErrorCode DMPlexSetConeOrientation(DM dm, PetscInt p, const PetscInt coneOrientation[])
1489: {
1490: DM_Plex *mesh = (DM_Plex*) dm->data;
1491: PetscInt pStart, pEnd;
1492: PetscInt dof, off, c;
1497: PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd);
1498: PetscSectionGetDof(mesh->coneSection, p, &dof);
1500: PetscSectionGetOffset(mesh->coneSection, p, &off);
1501: 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);
1502: for (c = 0; c < dof; ++c) {
1503: PetscInt cdof, o = coneOrientation[c];
1505: PetscSectionGetDof(mesh->coneSection, mesh->cones[off+c], &cdof);
1506: 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);
1507: mesh->coneOrientations[off+c] = o;
1508: }
1509: return(0);
1510: }
1512: /*@
1513: DMPlexInsertCone - Insert a point into the in-edges for the point p in the DAG
1515: Not collective
1517: Input Parameters:
1518: + mesh - The DMPlex
1519: . p - The point, which must lie in the chart set with DMPlexSetChart()
1520: . conePos - The local index in the cone where the point should be put
1521: - conePoint - The mesh point to insert
1523: Level: beginner
1525: .seealso: DMPlexCreate(), DMPlexGetCone(), DMPlexSetChart(), DMPlexSetConeSize(), DMSetUp()
1526: @*/
1527: PetscErrorCode DMPlexInsertCone(DM dm, PetscInt p, PetscInt conePos, PetscInt conePoint)
1528: {
1529: DM_Plex *mesh = (DM_Plex*) dm->data;
1530: PetscInt pStart, pEnd;
1531: PetscInt dof, off;
1536: PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd);
1537: 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);
1538: 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);
1539: PetscSectionGetDof(mesh->coneSection, p, &dof);
1540: PetscSectionGetOffset(mesh->coneSection, p, &off);
1541: 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);
1542: mesh->cones[off+conePos] = conePoint;
1543: return(0);
1544: }
1546: /*@
1547: DMPlexInsertConeOrientation - Insert a point orientation for the in-edge for the point p in the DAG
1549: Not collective
1551: Input Parameters:
1552: + mesh - The DMPlex
1553: . p - The point, which must lie in the chart set with DMPlexSetChart()
1554: . conePos - The local index in the cone where the point should be put
1555: - coneOrientation - The point orientation to insert
1557: Level: beginner
1559: .seealso: DMPlexCreate(), DMPlexGetCone(), DMPlexSetChart(), DMPlexSetConeSize(), DMSetUp()
1560: @*/
1561: PetscErrorCode DMPlexInsertConeOrientation(DM dm, PetscInt p, PetscInt conePos, PetscInt coneOrientation)
1562: {
1563: DM_Plex *mesh = (DM_Plex*) dm->data;
1564: PetscInt pStart, pEnd;
1565: PetscInt dof, off;
1570: PetscSectionGetChart(mesh->coneSection, &pStart, &pEnd);
1571: 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);
1572: PetscSectionGetDof(mesh->coneSection, p, &dof);
1573: PetscSectionGetOffset(mesh->coneSection, p, &off);
1574: 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);
1575: mesh->coneOrientations[off+conePos] = coneOrientation;
1576: return(0);
1577: }
1579: /*@
1580: DMPlexGetSupportSize - Return the number of out-edges for this point in the DAG
1582: Not collective
1584: Input Parameters:
1585: + mesh - The DMPlex
1586: - p - The point, which must lie in the chart set with DMPlexSetChart()
1588: Output Parameter:
1589: . size - The support size for point p
1591: Level: beginner
1593: .seealso: DMPlexCreate(), DMPlexSetConeSize(), DMPlexSetChart(), DMPlexGetConeSize()
1594: @*/
1595: PetscErrorCode DMPlexGetSupportSize(DM dm, PetscInt p, PetscInt *size)
1596: {
1597: DM_Plex *mesh = (DM_Plex*) dm->data;
1603: PetscSectionGetDof(mesh->supportSection, p, size);
1604: return(0);
1605: }
1607: /*@
1608: DMPlexSetSupportSize - Set the number of out-edges for this point in the DAG
1610: Not collective
1612: Input Parameters:
1613: + mesh - The DMPlex
1614: . p - The point, which must lie in the chart set with DMPlexSetChart()
1615: - size - The support size for point p
1617: Output Parameter:
1619: Note:
1620: This should be called after DMPlexSetChart().
1622: Level: beginner
1624: .seealso: DMPlexCreate(), DMPlexGetSupportSize(), DMPlexSetChart()
1625: @*/
1626: PetscErrorCode DMPlexSetSupportSize(DM dm, PetscInt p, PetscInt size)
1627: {
1628: DM_Plex *mesh = (DM_Plex*) dm->data;
1633: PetscSectionSetDof(mesh->supportSection, p, size);
1635: mesh->maxSupportSize = PetscMax(mesh->maxSupportSize, size);
1636: return(0);
1637: }
1639: /*@C
1640: DMPlexGetSupport - Return the points on the out-edges for this point in the DAG
1642: Not collective
1644: Input Parameters:
1645: + mesh - The DMPlex
1646: - p - The point, which must lie in the chart set with DMPlexSetChart()
1648: Output Parameter:
1649: . support - An array of points which are on the out-edges for point p
1651: Level: beginner
1653: Fortran Notes:
1654: Since it returns an array, this routine is only available in Fortran 90, and you must
1655: include petsc.h90 in your code.
1657: You must also call DMPlexRestoreSupport() after you finish using the returned array.
1659: .seealso: DMPlexCreate(), DMPlexSetCone(), DMPlexSetChart(), DMPlexGetCone()
1660: @*/
1661: PetscErrorCode DMPlexGetSupport(DM dm, PetscInt p, const PetscInt *support[])
1662: {
1663: DM_Plex *mesh = (DM_Plex*) dm->data;
1664: PetscInt off;
1670: PetscSectionGetOffset(mesh->supportSection, p, &off);
1671: *support = &mesh->supports[off];
1672: return(0);
1673: }
1675: /*@
1676: DMPlexSetSupport - Set the points on the out-edges for this point in the DAG
1678: Not collective
1680: Input Parameters:
1681: + mesh - The DMPlex
1682: . p - The point, which must lie in the chart set with DMPlexSetChart()
1683: - support - An array of points which are on the in-edges for point p
1685: Output Parameter:
1687: Note:
1688: This should be called after all calls to DMPlexSetSupportSize() and DMSetUp().
1690: Level: beginner
1692: .seealso: DMPlexCreate(), DMPlexGetSupport(), DMPlexSetChart(), DMPlexSetSupportSize(), DMSetUp()
1693: @*/
1694: PetscErrorCode DMPlexSetSupport(DM dm, PetscInt p, const PetscInt support[])
1695: {
1696: DM_Plex *mesh = (DM_Plex*) dm->data;
1697: PetscInt pStart, pEnd;
1698: PetscInt dof, off, c;
1703: PetscSectionGetChart(mesh->supportSection, &pStart, &pEnd);
1704: PetscSectionGetDof(mesh->supportSection, p, &dof);
1706: PetscSectionGetOffset(mesh->supportSection, p, &off);
1707: 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);
1708: for (c = 0; c < dof; ++c) {
1709: 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);
1710: mesh->supports[off+c] = support[c];
1711: }
1712: return(0);
1713: }
1715: /*@
1716: DMPlexInsertSupport - Insert a point into the out-edges for the point p in the DAG
1718: Not collective
1720: Input Parameters:
1721: + mesh - The DMPlex
1722: . p - The point, which must lie in the chart set with DMPlexSetChart()
1723: . supportPos - The local index in the cone where the point should be put
1724: - supportPoint - The mesh point to insert
1726: Level: beginner
1728: .seealso: DMPlexCreate(), DMPlexGetCone(), DMPlexSetChart(), DMPlexSetConeSize(), DMSetUp()
1729: @*/
1730: PetscErrorCode DMPlexInsertSupport(DM dm, PetscInt p, PetscInt supportPos, PetscInt supportPoint)
1731: {
1732: DM_Plex *mesh = (DM_Plex*) dm->data;
1733: PetscInt pStart, pEnd;
1734: PetscInt dof, off;
1739: PetscSectionGetChart(mesh->supportSection, &pStart, &pEnd);
1740: PetscSectionGetDof(mesh->supportSection, p, &dof);
1741: PetscSectionGetOffset(mesh->supportSection, p, &off);
1742: 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);
1743: 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);
1744: 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);
1745: mesh->supports[off+supportPos] = supportPoint;
1746: return(0);
1747: }
1749: /*@C
1750: DMPlexGetTransitiveClosure - Return the points on the transitive closure of the in-edges or out-edges for this point in the DAG
1752: Not collective
1754: Input Parameters:
1755: + mesh - The DMPlex
1756: . p - The point, which must lie in the chart set with DMPlexSetChart()
1757: . useCone - PETSC_TRUE for in-edges, otherwise use out-edges
1758: - points - If points is NULL on input, internal storage will be returned, otherwise the provided array is used
1760: Output Parameters:
1761: + numPoints - The number of points in the closure, so points[] is of size 2*numPoints
1762: - points - The points and point orientations, interleaved as pairs [p0, o0, p1, o1, ...]
1764: Note:
1765: If using internal storage (points is NULL on input), each call overwrites the last output.
1767: Fortran Notes:
1768: Since it returns an array, this routine is only available in Fortran 90, and you must
1769: include petsc.h90 in your code.
1771: The numPoints argument is not present in the Fortran 90 binding since it is internal to the array.
1773: Level: beginner
1775: .seealso: DMPlexRestoreTransitiveClosure(), DMPlexCreate(), DMPlexSetCone(), DMPlexSetChart(), DMPlexGetCone()
1776: @*/
1777: PetscErrorCode DMPlexGetTransitiveClosure(DM dm, PetscInt p, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
1778: {
1779: DM_Plex *mesh = (DM_Plex*) dm->data;
1780: PetscInt *closure, *fifo;
1781: const PetscInt *tmp = NULL, *tmpO = NULL;
1782: PetscInt tmpSize, t;
1783: PetscInt depth = 0, maxSize;
1784: PetscInt closureSize = 2, fifoSize = 0, fifoStart = 0;
1785: PetscErrorCode ierr;
1789: DMPlexGetDepth(dm, &depth);
1790: /* This is only 1-level */
1791: if (useCone) {
1792: DMPlexGetConeSize(dm, p, &tmpSize);
1793: DMPlexGetCone(dm, p, &tmp);
1794: DMPlexGetConeOrientation(dm, p, &tmpO);
1795: } else {
1796: DMPlexGetSupportSize(dm, p, &tmpSize);
1797: DMPlexGetSupport(dm, p, &tmp);
1798: }
1799: if (depth == 1) {
1800: if (*points) {
1801: closure = *points;
1802: } else {
1803: maxSize = 2*(PetscMax(mesh->maxConeSize, mesh->maxSupportSize)+1);
1804: DMGetWorkArray(dm, maxSize, PETSC_INT, &closure);
1805: }
1806: closure[0] = p; closure[1] = 0;
1807: for (t = 0; t < tmpSize; ++t, closureSize += 2) {
1808: closure[closureSize] = tmp[t];
1809: closure[closureSize+1] = tmpO ? tmpO[t] : 0;
1810: }
1811: if (numPoints) *numPoints = closureSize/2;
1812: if (points) *points = closure;
1813: return(0);
1814: }
1815: {
1816: PetscInt c, coneSeries, s,supportSeries;
1818: c = mesh->maxConeSize;
1819: coneSeries = (c > 1) ? ((PetscPowInt(c,depth+1)-1)/(c-1)) : depth+1;
1820: s = mesh->maxSupportSize;
1821: supportSeries = (s > 1) ? ((PetscPowInt(s,depth+1)-1)/(s-1)) : depth+1;
1822: maxSize = 2*PetscMax(coneSeries,supportSeries);
1823: }
1824: DMGetWorkArray(dm, maxSize, PETSC_INT, &fifo);
1825: if (*points) {
1826: closure = *points;
1827: } else {
1828: DMGetWorkArray(dm, maxSize, PETSC_INT, &closure);
1829: }
1830: closure[0] = p; closure[1] = 0;
1831: for (t = 0; t < tmpSize; ++t, closureSize += 2, fifoSize += 2) {
1832: const PetscInt cp = tmp[t];
1833: const PetscInt co = tmpO ? tmpO[t] : 0;
1835: closure[closureSize] = cp;
1836: closure[closureSize+1] = co;
1837: fifo[fifoSize] = cp;
1838: fifo[fifoSize+1] = co;
1839: }
1840: /* Should kick out early when depth is reached, rather than checking all vertices for empty cones */
1841: while (fifoSize - fifoStart) {
1842: const PetscInt q = fifo[fifoStart];
1843: const PetscInt o = fifo[fifoStart+1];
1844: const PetscInt rev = o >= 0 ? 0 : 1;
1845: const PetscInt off = rev ? -(o+1) : o;
1847: if (useCone) {
1848: DMPlexGetConeSize(dm, q, &tmpSize);
1849: DMPlexGetCone(dm, q, &tmp);
1850: DMPlexGetConeOrientation(dm, q, &tmpO);
1851: } else {
1852: DMPlexGetSupportSize(dm, q, &tmpSize);
1853: DMPlexGetSupport(dm, q, &tmp);
1854: tmpO = NULL;
1855: }
1856: for (t = 0; t < tmpSize; ++t) {
1857: const PetscInt i = ((rev ? tmpSize-t : t) + off)%tmpSize;
1858: const PetscInt cp = tmp[i];
1859: /* Must propogate orientation: When we reverse orientation, we both reverse the direction of iteration and start at the other end of the chain. */
1860: /* HACK: It is worse to get the size here, than to change the interpretation of -(*+1)
1861: const PetscInt co = tmpO ? (rev ? -(tmpO[i]+1) : tmpO[i]) : 0; */
1862: PetscInt co = tmpO ? tmpO[i] : 0;
1863: PetscInt c;
1865: if (rev) {
1866: PetscInt childSize, coff;
1867: DMPlexGetConeSize(dm, cp, &childSize);
1868: coff = tmpO[i] < 0 ? -(tmpO[i]+1) : tmpO[i];
1869: co = childSize ? -(((coff+childSize-1)%childSize)+1) : 0;
1870: }
1871: /* Check for duplicate */
1872: for (c = 0; c < closureSize; c += 2) {
1873: if (closure[c] == cp) break;
1874: }
1875: if (c == closureSize) {
1876: closure[closureSize] = cp;
1877: closure[closureSize+1] = co;
1878: fifo[fifoSize] = cp;
1879: fifo[fifoSize+1] = co;
1880: closureSize += 2;
1881: fifoSize += 2;
1882: }
1883: }
1884: fifoStart += 2;
1885: }
1886: if (numPoints) *numPoints = closureSize/2;
1887: if (points) *points = closure;
1888: DMRestoreWorkArray(dm, maxSize, PETSC_INT, &fifo);
1889: return(0);
1890: }
1892: /*@C
1893: DMPlexGetTransitiveClosure_Internal - Return the points on the transitive closure of the in-edges or out-edges for this point in the DAG with a specified initial orientation
1895: Not collective
1897: Input Parameters:
1898: + mesh - The DMPlex
1899: . p - The point, which must lie in the chart set with DMPlexSetChart()
1900: . orientation - The orientation of the point
1901: . useCone - PETSC_TRUE for in-edges, otherwise use out-edges
1902: - points - If points is NULL on input, internal storage will be returned, otherwise the provided array is used
1904: Output Parameters:
1905: + numPoints - The number of points in the closure, so points[] is of size 2*numPoints
1906: - points - The points and point orientations, interleaved as pairs [p0, o0, p1, o1, ...]
1908: Note:
1909: If using internal storage (points is NULL on input), each call overwrites the last output.
1911: Fortran Notes:
1912: Since it returns an array, this routine is only available in Fortran 90, and you must
1913: include petsc.h90 in your code.
1915: The numPoints argument is not present in the Fortran 90 binding since it is internal to the array.
1917: Level: beginner
1919: .seealso: DMPlexRestoreTransitiveClosure(), DMPlexCreate(), DMPlexSetCone(), DMPlexSetChart(), DMPlexGetCone()
1920: @*/
1921: PetscErrorCode DMPlexGetTransitiveClosure_Internal(DM dm, PetscInt p, PetscInt ornt, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
1922: {
1923: DM_Plex *mesh = (DM_Plex*) dm->data;
1924: PetscInt *closure, *fifo;
1925: const PetscInt *tmp = NULL, *tmpO = NULL;
1926: PetscInt tmpSize, t;
1927: PetscInt depth = 0, maxSize;
1928: PetscInt closureSize = 2, fifoSize = 0, fifoStart = 0;
1929: PetscErrorCode ierr;
1933: DMPlexGetDepth(dm, &depth);
1934: /* This is only 1-level */
1935: if (useCone) {
1936: DMPlexGetConeSize(dm, p, &tmpSize);
1937: DMPlexGetCone(dm, p, &tmp);
1938: DMPlexGetConeOrientation(dm, p, &tmpO);
1939: } else {
1940: DMPlexGetSupportSize(dm, p, &tmpSize);
1941: DMPlexGetSupport(dm, p, &tmp);
1942: }
1943: if (depth == 1) {
1944: if (*points) {
1945: closure = *points;
1946: } else {
1947: maxSize = 2*(PetscMax(mesh->maxConeSize, mesh->maxSupportSize)+1);
1948: DMGetWorkArray(dm, maxSize, PETSC_INT, &closure);
1949: }
1950: closure[0] = p; closure[1] = ornt;
1951: for (t = 0; t < tmpSize; ++t, closureSize += 2) {
1952: const PetscInt i = ornt >= 0 ? (t+ornt)%tmpSize : (-(ornt+1) + tmpSize-t)%tmpSize;
1953: closure[closureSize] = tmp[i];
1954: closure[closureSize+1] = tmpO ? tmpO[i] : 0;
1955: }
1956: if (numPoints) *numPoints = closureSize/2;
1957: if (points) *points = closure;
1958: return(0);
1959: }
1960: {
1961: PetscInt c, coneSeries, s,supportSeries;
1963: c = mesh->maxConeSize;
1964: coneSeries = (c > 1) ? ((PetscPowInt(c,depth+1)-1)/(c-1)) : depth+1;
1965: s = mesh->maxSupportSize;
1966: supportSeries = (s > 1) ? ((PetscPowInt(s,depth+1)-1)/(s-1)) : depth+1;
1967: maxSize = 2*PetscMax(coneSeries,supportSeries);
1968: }
1969: DMGetWorkArray(dm, maxSize, PETSC_INT, &fifo);
1970: if (*points) {
1971: closure = *points;
1972: } else {
1973: DMGetWorkArray(dm, maxSize, PETSC_INT, &closure);
1974: }
1975: closure[0] = p; closure[1] = ornt;
1976: for (t = 0; t < tmpSize; ++t, closureSize += 2, fifoSize += 2) {
1977: const PetscInt i = ornt >= 0 ? (t+ornt)%tmpSize : (-(ornt+1) + tmpSize-t)%tmpSize;
1978: const PetscInt cp = tmp[i];
1979: PetscInt co = tmpO ? tmpO[i] : 0;
1981: if (ornt < 0) {
1982: PetscInt childSize, coff;
1983: DMPlexGetConeSize(dm, cp, &childSize);
1984: coff = co < 0 ? -(tmpO[i]+1) : tmpO[i];
1985: co = childSize ? -(((coff+childSize-1)%childSize)+1) : 0;
1986: }
1987: closure[closureSize] = cp;
1988: closure[closureSize+1] = co;
1989: fifo[fifoSize] = cp;
1990: fifo[fifoSize+1] = co;
1991: }
1992: /* Should kick out early when depth is reached, rather than checking all vertices for empty cones */
1993: while (fifoSize - fifoStart) {
1994: const PetscInt q = fifo[fifoStart];
1995: const PetscInt o = fifo[fifoStart+1];
1996: const PetscInt rev = o >= 0 ? 0 : 1;
1997: const PetscInt off = rev ? -(o+1) : o;
1999: if (useCone) {
2000: DMPlexGetConeSize(dm, q, &tmpSize);
2001: DMPlexGetCone(dm, q, &tmp);
2002: DMPlexGetConeOrientation(dm, q, &tmpO);
2003: } else {
2004: DMPlexGetSupportSize(dm, q, &tmpSize);
2005: DMPlexGetSupport(dm, q, &tmp);
2006: tmpO = NULL;
2007: }
2008: for (t = 0; t < tmpSize; ++t) {
2009: const PetscInt i = ((rev ? tmpSize-t : t) + off)%tmpSize;
2010: const PetscInt cp = tmp[i];
2011: /* Must propogate orientation: When we reverse orientation, we both reverse the direction of iteration and start at the other end of the chain. */
2012: /* HACK: It is worse to get the size here, than to change the interpretation of -(*+1)
2013: const PetscInt co = tmpO ? (rev ? -(tmpO[i]+1) : tmpO[i]) : 0; */
2014: PetscInt co = tmpO ? tmpO[i] : 0;
2015: PetscInt c;
2017: if (rev) {
2018: PetscInt childSize, coff;
2019: DMPlexGetConeSize(dm, cp, &childSize);
2020: coff = tmpO[i] < 0 ? -(tmpO[i]+1) : tmpO[i];
2021: co = childSize ? -(((coff+childSize-1)%childSize)+1) : 0;
2022: }
2023: /* Check for duplicate */
2024: for (c = 0; c < closureSize; c += 2) {
2025: if (closure[c] == cp) break;
2026: }
2027: if (c == closureSize) {
2028: closure[closureSize] = cp;
2029: closure[closureSize+1] = co;
2030: fifo[fifoSize] = cp;
2031: fifo[fifoSize+1] = co;
2032: closureSize += 2;
2033: fifoSize += 2;
2034: }
2035: }
2036: fifoStart += 2;
2037: }
2038: if (numPoints) *numPoints = closureSize/2;
2039: if (points) *points = closure;
2040: DMRestoreWorkArray(dm, maxSize, PETSC_INT, &fifo);
2041: return(0);
2042: }
2044: /*@C
2045: DMPlexRestoreTransitiveClosure - Restore the array of points on the transitive closure of the in-edges or out-edges for this point in the DAG
2047: Not collective
2049: Input Parameters:
2050: + mesh - The DMPlex
2051: . p - The point, which must lie in the chart set with DMPlexSetChart()
2052: . useCone - PETSC_TRUE for in-edges, otherwise use out-edges
2053: . numPoints - The number of points in the closure, so points[] is of size 2*numPoints, zeroed on exit
2054: - points - The points and point orientations, interleaved as pairs [p0, o0, p1, o1, ...], zeroed on exit
2056: Note:
2057: If not using internal storage (points is not NULL on input), this call is unnecessary
2059: Fortran Notes:
2060: Since it returns an array, this routine is only available in Fortran 90, and you must
2061: include petsc.h90 in your code.
2063: The numPoints argument is not present in the Fortran 90 binding since it is internal to the array.
2065: Level: beginner
2067: .seealso: DMPlexGetTransitiveClosure(), DMPlexCreate(), DMPlexSetCone(), DMPlexSetChart(), DMPlexGetCone()
2068: @*/
2069: PetscErrorCode DMPlexRestoreTransitiveClosure(DM dm, PetscInt p, PetscBool useCone, PetscInt *numPoints, PetscInt *points[])
2070: {
2077: DMRestoreWorkArray(dm, 0, PETSC_INT, points);
2078: if (numPoints) *numPoints = 0;
2079: return(0);
2080: }
2082: /*@
2083: DMPlexGetMaxSizes - Return the maximum number of in-edges (cone) and out-edges (support) for any point in the DAG
2085: Not collective
2087: Input Parameter:
2088: . mesh - The DMPlex
2090: Output Parameters:
2091: + maxConeSize - The maximum number of in-edges
2092: - maxSupportSize - The maximum number of out-edges
2094: Level: beginner
2096: .seealso: DMPlexCreate(), DMPlexSetConeSize(), DMPlexSetChart()
2097: @*/
2098: PetscErrorCode DMPlexGetMaxSizes(DM dm, PetscInt *maxConeSize, PetscInt *maxSupportSize)
2099: {
2100: DM_Plex *mesh = (DM_Plex*) dm->data;
2104: if (maxConeSize) *maxConeSize = mesh->maxConeSize;
2105: if (maxSupportSize) *maxSupportSize = mesh->maxSupportSize;
2106: return(0);
2107: }
2109: PetscErrorCode DMSetUp_Plex(DM dm)
2110: {
2111: DM_Plex *mesh = (DM_Plex*) dm->data;
2112: PetscInt size;
2117: PetscSectionSetUp(mesh->coneSection);
2118: PetscSectionGetStorageSize(mesh->coneSection, &size);
2119: PetscMalloc1(size, &mesh->cones);
2120: PetscCalloc1(size, &mesh->coneOrientations);
2121: if (mesh->maxSupportSize) {
2122: PetscSectionSetUp(mesh->supportSection);
2123: PetscSectionGetStorageSize(mesh->supportSection, &size);
2124: PetscMalloc1(size, &mesh->supports);
2125: }
2126: return(0);
2127: }
2129: PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
2130: {
2134: if (subdm) {DMClone(dm, subdm);}
2135: DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);
2136: return(0);
2137: }
2139: /*@
2140: DMPlexSymmetrize - Create support (out-edge) information from cone (in-edge) information
2142: Not collective
2144: Input Parameter:
2145: . mesh - The DMPlex
2147: Output Parameter:
2149: Note:
2150: This should be called after all calls to DMPlexSetCone()
2152: Level: beginner
2154: .seealso: DMPlexCreate(), DMPlexSetChart(), DMPlexSetConeSize(), DMPlexSetCone()
2155: @*/
2156: PetscErrorCode DMPlexSymmetrize(DM dm)
2157: {
2158: DM_Plex *mesh = (DM_Plex*) dm->data;
2159: PetscInt *offsets;
2160: PetscInt supportSize;
2161: PetscInt pStart, pEnd, p;
2166: if (mesh->supports) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Supports were already setup in this DMPlex");
2167: /* Calculate support sizes */
2168: DMPlexGetChart(dm, &pStart, &pEnd);
2169: for (p = pStart; p < pEnd; ++p) {
2170: PetscInt dof, off, c;
2172: PetscSectionGetDof(mesh->coneSection, p, &dof);
2173: PetscSectionGetOffset(mesh->coneSection, p, &off);
2174: for (c = off; c < off+dof; ++c) {
2175: PetscSectionAddDof(mesh->supportSection, mesh->cones[c], 1);
2176: }
2177: }
2178: for (p = pStart; p < pEnd; ++p) {
2179: PetscInt dof;
2181: PetscSectionGetDof(mesh->supportSection, p, &dof);
2183: mesh->maxSupportSize = PetscMax(mesh->maxSupportSize, dof);
2184: }
2185: PetscSectionSetUp(mesh->supportSection);
2186: /* Calculate supports */
2187: PetscSectionGetStorageSize(mesh->supportSection, &supportSize);
2188: PetscMalloc1(supportSize, &mesh->supports);
2189: PetscCalloc1(pEnd - pStart, &offsets);
2190: for (p = pStart; p < pEnd; ++p) {
2191: PetscInt dof, off, c;
2193: PetscSectionGetDof(mesh->coneSection, p, &dof);
2194: PetscSectionGetOffset(mesh->coneSection, p, &off);
2195: for (c = off; c < off+dof; ++c) {
2196: const PetscInt q = mesh->cones[c];
2197: PetscInt offS;
2199: PetscSectionGetOffset(mesh->supportSection, q, &offS);
2201: mesh->supports[offS+offsets[q]] = p;
2202: ++offsets[q];
2203: }
2204: }
2205: PetscFree(offsets);
2206: return(0);
2207: }
2209: /*@
2210: DMPlexStratify - The DAG for most topologies is a graded poset (http://en.wikipedia.org/wiki/Graded_poset), and
2211: can be illustrated by a Hasse Diagram (a http://en.wikipedia.org/wiki/Hasse_diagram). The strata group all points of the
2212: same grade, and this function calculates the strata. This grade can be seen as the height (or depth) of the point in
2213: the DAG.
2215: Collective on dm
2217: Input Parameter:
2218: . mesh - The DMPlex
2220: Output Parameter:
2222: Notes:
2223: Concretely, DMPlexStratify() creates a new label named "depth" containing the dimension of each element: 0 for vertices,
2224: 1 for edges, and so on. The depth label can be accessed through DMPlexGetDepthLabel() or DMPlexGetDepthStratum(), or
2225: manually via DMGetLabel(). The height is defined implicitly by height = maxDimension - depth, and can be accessed
2226: via DMPlexGetHeightStratum(). For example, cells have height 0 and faces have height 1.
2228: DMPlexStratify() should be called after all calls to DMPlexSymmetrize()
2230: Level: beginner
2232: .seealso: DMPlexCreate(), DMPlexSymmetrize()
2233: @*/
2234: PetscErrorCode DMPlexStratify(DM dm)
2235: {
2236: DM_Plex *mesh = (DM_Plex*) dm->data;
2237: DMLabel label;
2238: PetscInt pStart, pEnd, p;
2239: PetscInt numRoots = 0, numLeaves = 0;
2244: PetscLogEventBegin(DMPLEX_Stratify,dm,0,0,0);
2245: /* Calculate depth */
2246: DMPlexGetChart(dm, &pStart, &pEnd);
2247: DMCreateLabel(dm, "depth");
2248: DMPlexGetDepthLabel(dm, &label);
2249: /* Initialize roots and count leaves */
2250: for (p = pStart; p < pEnd; ++p) {
2251: PetscInt coneSize, supportSize;
2253: DMPlexGetConeSize(dm, p, &coneSize);
2254: DMPlexGetSupportSize(dm, p, &supportSize);
2255: if (!coneSize && supportSize) {
2256: ++numRoots;
2257: DMLabelSetValue(label, p, 0);
2258: } else if (!supportSize && coneSize) {
2259: ++numLeaves;
2260: } else if (!supportSize && !coneSize) {
2261: /* Isolated points */
2262: DMLabelSetValue(label, p, 0);
2263: }
2264: }
2265: if (numRoots + numLeaves == (pEnd - pStart)) {
2266: for (p = pStart; p < pEnd; ++p) {
2267: PetscInt coneSize, supportSize;
2269: DMPlexGetConeSize(dm, p, &coneSize);
2270: DMPlexGetSupportSize(dm, p, &supportSize);
2271: if (!supportSize && coneSize) {
2272: DMLabelSetValue(label, p, 1);
2273: }
2274: }
2275: } else {
2276: IS pointIS;
2277: PetscInt numPoints = 0, level = 0;
2279: DMLabelGetStratumIS(label, level, &pointIS);
2280: if (pointIS) {ISGetLocalSize(pointIS, &numPoints);}
2281: while (numPoints) {
2282: const PetscInt *points;
2283: const PetscInt newLevel = level+1;
2285: ISGetIndices(pointIS, &points);
2286: for (p = 0; p < numPoints; ++p) {
2287: const PetscInt point = points[p];
2288: const PetscInt *support;
2289: PetscInt supportSize, s;
2291: DMPlexGetSupportSize(dm, point, &supportSize);
2292: DMPlexGetSupport(dm, point, &support);
2293: for (s = 0; s < supportSize; ++s) {
2294: DMLabelSetValue(label, support[s], newLevel);
2295: }
2296: }
2297: ISRestoreIndices(pointIS, &points);
2298: ++level;
2299: ISDestroy(&pointIS);
2300: DMLabelGetStratumIS(label, level, &pointIS);
2301: if (pointIS) {ISGetLocalSize(pointIS, &numPoints);}
2302: else {numPoints = 0;}
2303: }
2304: ISDestroy(&pointIS);
2305: }
2306: { /* just in case there is an empty process */
2307: PetscInt numValues, maxValues = 0, v;
2309: DMLabelGetNumValues(label,&numValues);
2310: for (v = 0; v < numValues; v++) {
2311: IS pointIS;
2313: DMLabelGetStratumIS(label, v, &pointIS);
2314: if (pointIS) {
2315: PetscInt min, max, numPoints;
2316: PetscInt start;
2317: PetscBool contig;
2319: ISGetLocalSize(pointIS, &numPoints);
2320: ISGetMinMax(pointIS, &min, &max);
2321: ISContiguousLocal(pointIS,min,max+1,&start,&contig);
2322: if (start == 0 && contig) {
2323: ISDestroy(&pointIS);
2324: ISCreateStride(PETSC_COMM_SELF,numPoints,min,1,&pointIS);
2325: DMLabelSetStratumIS(label, v, pointIS);
2326: }
2327: }
2328: ISDestroy(&pointIS);
2329: }
2330: MPI_Allreduce(&numValues,&maxValues,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)dm));
2331: for (v = numValues; v < maxValues; v++) {
2332: DMLabelAddStratum(label,v);
2333: }
2334: }
2336: DMLabelGetState(label, &mesh->depthState);
2337: PetscLogEventEnd(DMPLEX_Stratify,dm,0,0,0);
2338: return(0);
2339: }
2341: /*@C
2342: DMPlexGetJoin - Get an array for the join of the set of points
2344: Not Collective
2346: Input Parameters:
2347: + dm - The DMPlex object
2348: . numPoints - The number of input points for the join
2349: - points - The input points
2351: Output Parameters:
2352: + numCoveredPoints - The number of points in the join
2353: - coveredPoints - The points in the join
2355: Level: intermediate
2357: Note: Currently, this is restricted to a single level join
2359: Fortran Notes:
2360: Since it returns an array, this routine is only available in Fortran 90, and you must
2361: include petsc.h90 in your code.
2363: The numCoveredPoints argument is not present in the Fortran 90 binding since it is internal to the array.
2365: .keywords: mesh
2366: .seealso: DMPlexRestoreJoin(), DMPlexGetMeet()
2367: @*/
2368: PetscErrorCode DMPlexGetJoin(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt **coveredPoints)
2369: {
2370: DM_Plex *mesh = (DM_Plex*) dm->data;
2371: PetscInt *join[2];
2372: PetscInt joinSize, i = 0;
2373: PetscInt dof, off, p, c, m;
2381: DMGetWorkArray(dm, mesh->maxSupportSize, PETSC_INT, &join[0]);
2382: DMGetWorkArray(dm, mesh->maxSupportSize, PETSC_INT, &join[1]);
2383: /* Copy in support of first point */
2384: PetscSectionGetDof(mesh->supportSection, points[0], &dof);
2385: PetscSectionGetOffset(mesh->supportSection, points[0], &off);
2386: for (joinSize = 0; joinSize < dof; ++joinSize) {
2387: join[i][joinSize] = mesh->supports[off+joinSize];
2388: }
2389: /* Check each successive support */
2390: for (p = 1; p < numPoints; ++p) {
2391: PetscInt newJoinSize = 0;
2393: PetscSectionGetDof(mesh->supportSection, points[p], &dof);
2394: PetscSectionGetOffset(mesh->supportSection, points[p], &off);
2395: for (c = 0; c < dof; ++c) {
2396: const PetscInt point = mesh->supports[off+c];
2398: for (m = 0; m < joinSize; ++m) {
2399: if (point == join[i][m]) {
2400: join[1-i][newJoinSize++] = point;
2401: break;
2402: }
2403: }
2404: }
2405: joinSize = newJoinSize;
2406: i = 1-i;
2407: }
2408: *numCoveredPoints = joinSize;
2409: *coveredPoints = join[i];
2410: DMRestoreWorkArray(dm, mesh->maxSupportSize, PETSC_INT, &join[1-i]);
2411: return(0);
2412: }
2414: /*@C
2415: DMPlexRestoreJoin - Restore an array for the join of the set of points
2417: Not Collective
2419: Input Parameters:
2420: + dm - The DMPlex object
2421: . numPoints - The number of input points for the join
2422: - points - The input points
2424: Output Parameters:
2425: + numCoveredPoints - The number of points in the join
2426: - coveredPoints - The points in the join
2428: Fortran Notes:
2429: Since it returns an array, this routine is only available in Fortran 90, and you must
2430: include petsc.h90 in your code.
2432: The numCoveredPoints argument is not present in the Fortran 90 binding since it is internal to the array.
2434: Level: intermediate
2436: .keywords: mesh
2437: .seealso: DMPlexGetJoin(), DMPlexGetFullJoin(), DMPlexGetMeet()
2438: @*/
2439: PetscErrorCode DMPlexRestoreJoin(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt **coveredPoints)
2440: {
2448: DMRestoreWorkArray(dm, 0, PETSC_INT, (void*) coveredPoints);
2449: if (numCoveredPoints) *numCoveredPoints = 0;
2450: return(0);
2451: }
2453: /*@C
2454: DMPlexGetFullJoin - Get an array for the join of the set of points
2456: Not Collective
2458: Input Parameters:
2459: + dm - The DMPlex object
2460: . numPoints - The number of input points for the join
2461: - points - The input points
2463: Output Parameters:
2464: + numCoveredPoints - The number of points in the join
2465: - coveredPoints - The points in the join
2467: Fortran Notes:
2468: Since it returns an array, this routine is only available in Fortran 90, and you must
2469: include petsc.h90 in your code.
2471: The numCoveredPoints argument is not present in the Fortran 90 binding since it is internal to the array.
2473: Level: intermediate
2475: .keywords: mesh
2476: .seealso: DMPlexGetJoin(), DMPlexRestoreJoin(), DMPlexGetMeet()
2477: @*/
2478: PetscErrorCode DMPlexGetFullJoin(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt **coveredPoints)
2479: {
2480: DM_Plex *mesh = (DM_Plex*) dm->data;
2481: PetscInt *offsets, **closures;
2482: PetscInt *join[2];
2483: PetscInt depth = 0, maxSize, joinSize = 0, i = 0;
2484: PetscInt p, d, c, m, ms;
2493: DMPlexGetDepth(dm, &depth);
2494: PetscCalloc1(numPoints, &closures);
2495: DMGetWorkArray(dm, numPoints*(depth+2), PETSC_INT, &offsets);
2496: ms = mesh->maxSupportSize;
2497: maxSize = (ms > 1) ? ((PetscPowInt(ms,depth+1)-1)/(ms-1)) : depth + 1;
2498: DMGetWorkArray(dm, maxSize, PETSC_INT, &join[0]);
2499: DMGetWorkArray(dm, maxSize, PETSC_INT, &join[1]);
2501: for (p = 0; p < numPoints; ++p) {
2502: PetscInt closureSize;
2504: DMPlexGetTransitiveClosure(dm, points[p], PETSC_FALSE, &closureSize, &closures[p]);
2506: offsets[p*(depth+2)+0] = 0;
2507: for (d = 0; d < depth+1; ++d) {
2508: PetscInt pStart, pEnd, i;
2510: DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
2511: for (i = offsets[p*(depth+2)+d]; i < closureSize; ++i) {
2512: if ((pStart > closures[p][i*2]) || (pEnd <= closures[p][i*2])) {
2513: offsets[p*(depth+2)+d+1] = i;
2514: break;
2515: }
2516: }
2517: if (i == closureSize) offsets[p*(depth+2)+d+1] = i;
2518: }
2519: 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);
2520: }
2521: for (d = 0; d < depth+1; ++d) {
2522: PetscInt dof;
2524: /* Copy in support of first point */
2525: dof = offsets[d+1] - offsets[d];
2526: for (joinSize = 0; joinSize < dof; ++joinSize) {
2527: join[i][joinSize] = closures[0][(offsets[d]+joinSize)*2];
2528: }
2529: /* Check each successive cone */
2530: for (p = 1; p < numPoints && joinSize; ++p) {
2531: PetscInt newJoinSize = 0;
2533: dof = offsets[p*(depth+2)+d+1] - offsets[p*(depth+2)+d];
2534: for (c = 0; c < dof; ++c) {
2535: const PetscInt point = closures[p][(offsets[p*(depth+2)+d]+c)*2];
2537: for (m = 0; m < joinSize; ++m) {
2538: if (point == join[i][m]) {
2539: join[1-i][newJoinSize++] = point;
2540: break;
2541: }
2542: }
2543: }
2544: joinSize = newJoinSize;
2545: i = 1-i;
2546: }
2547: if (joinSize) break;
2548: }
2549: *numCoveredPoints = joinSize;
2550: *coveredPoints = join[i];
2551: for (p = 0; p < numPoints; ++p) {
2552: DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_FALSE, NULL, &closures[p]);
2553: }
2554: PetscFree(closures);
2555: DMRestoreWorkArray(dm, numPoints*(depth+2), PETSC_INT, &offsets);
2556: DMRestoreWorkArray(dm, mesh->maxSupportSize, PETSC_INT, &join[1-i]);
2557: return(0);
2558: }
2560: /*@C
2561: DMPlexGetMeet - Get an array for the meet of the set of points
2563: Not Collective
2565: Input Parameters:
2566: + dm - The DMPlex object
2567: . numPoints - The number of input points for the meet
2568: - points - The input points
2570: Output Parameters:
2571: + numCoveredPoints - The number of points in the meet
2572: - coveredPoints - The points in the meet
2574: Level: intermediate
2576: Note: Currently, this is restricted to a single level meet
2578: Fortran Notes:
2579: Since it returns an array, this routine is only available in Fortran 90, and you must
2580: include petsc.h90 in your code.
2582: The numCoveredPoints argument is not present in the Fortran 90 binding since it is internal to the array.
2584: .keywords: mesh
2585: .seealso: DMPlexRestoreMeet(), DMPlexGetJoin()
2586: @*/
2587: PetscErrorCode DMPlexGetMeet(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveringPoints, const PetscInt **coveringPoints)
2588: {
2589: DM_Plex *mesh = (DM_Plex*) dm->data;
2590: PetscInt *meet[2];
2591: PetscInt meetSize, i = 0;
2592: PetscInt dof, off, p, c, m;
2600: DMGetWorkArray(dm, mesh->maxConeSize, PETSC_INT, &meet[0]);
2601: DMGetWorkArray(dm, mesh->maxConeSize, PETSC_INT, &meet[1]);
2602: /* Copy in cone of first point */
2603: PetscSectionGetDof(mesh->coneSection, points[0], &dof);
2604: PetscSectionGetOffset(mesh->coneSection, points[0], &off);
2605: for (meetSize = 0; meetSize < dof; ++meetSize) {
2606: meet[i][meetSize] = mesh->cones[off+meetSize];
2607: }
2608: /* Check each successive cone */
2609: for (p = 1; p < numPoints; ++p) {
2610: PetscInt newMeetSize = 0;
2612: PetscSectionGetDof(mesh->coneSection, points[p], &dof);
2613: PetscSectionGetOffset(mesh->coneSection, points[p], &off);
2614: for (c = 0; c < dof; ++c) {
2615: const PetscInt point = mesh->cones[off+c];
2617: for (m = 0; m < meetSize; ++m) {
2618: if (point == meet[i][m]) {
2619: meet[1-i][newMeetSize++] = point;
2620: break;
2621: }
2622: }
2623: }
2624: meetSize = newMeetSize;
2625: i = 1-i;
2626: }
2627: *numCoveringPoints = meetSize;
2628: *coveringPoints = meet[i];
2629: DMRestoreWorkArray(dm, mesh->maxConeSize, PETSC_INT, &meet[1-i]);
2630: return(0);
2631: }
2633: /*@C
2634: DMPlexRestoreMeet - Restore an array for the meet of the set of points
2636: Not Collective
2638: Input Parameters:
2639: + dm - The DMPlex object
2640: . numPoints - The number of input points for the meet
2641: - points - The input points
2643: Output Parameters:
2644: + numCoveredPoints - The number of points in the meet
2645: - coveredPoints - The points in the meet
2647: Level: intermediate
2649: Fortran Notes:
2650: Since it returns an array, this routine is only available in Fortran 90, and you must
2651: include petsc.h90 in your code.
2653: The numCoveredPoints argument is not present in the Fortran 90 binding since it is internal to the array.
2655: .keywords: mesh
2656: .seealso: DMPlexGetMeet(), DMPlexGetFullMeet(), DMPlexGetJoin()
2657: @*/
2658: PetscErrorCode DMPlexRestoreMeet(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt **coveredPoints)
2659: {
2667: DMRestoreWorkArray(dm, 0, PETSC_INT, (void*) coveredPoints);
2668: if (numCoveredPoints) *numCoveredPoints = 0;
2669: return(0);
2670: }
2672: /*@C
2673: DMPlexGetFullMeet - Get an array for the meet of the set of points
2675: Not Collective
2677: Input Parameters:
2678: + dm - The DMPlex object
2679: . numPoints - The number of input points for the meet
2680: - points - The input points
2682: Output Parameters:
2683: + numCoveredPoints - The number of points in the meet
2684: - coveredPoints - The points in the meet
2686: Level: intermediate
2688: Fortran Notes:
2689: Since it returns an array, this routine is only available in Fortran 90, and you must
2690: include petsc.h90 in your code.
2692: The numCoveredPoints argument is not present in the Fortran 90 binding since it is internal to the array.
2694: .keywords: mesh
2695: .seealso: DMPlexGetMeet(), DMPlexRestoreMeet(), DMPlexGetJoin()
2696: @*/
2697: PetscErrorCode DMPlexGetFullMeet(DM dm, PetscInt numPoints, const PetscInt points[], PetscInt *numCoveredPoints, const PetscInt **coveredPoints)
2698: {
2699: DM_Plex *mesh = (DM_Plex*) dm->data;
2700: PetscInt *offsets, **closures;
2701: PetscInt *meet[2];
2702: PetscInt height = 0, maxSize, meetSize = 0, i = 0;
2703: PetscInt p, h, c, m, mc;
2712: DMPlexGetDepth(dm, &height);
2713: PetscMalloc1(numPoints, &closures);
2714: DMGetWorkArray(dm, numPoints*(height+2), PETSC_INT, &offsets);
2715: mc = mesh->maxConeSize;
2716: maxSize = (mc > 1) ? ((PetscPowInt(mc,height+1)-1)/(mc-1)) : height + 1;
2717: DMGetWorkArray(dm, maxSize, PETSC_INT, &meet[0]);
2718: DMGetWorkArray(dm, maxSize, PETSC_INT, &meet[1]);
2720: for (p = 0; p < numPoints; ++p) {
2721: PetscInt closureSize;
2723: DMPlexGetTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closures[p]);
2725: offsets[p*(height+2)+0] = 0;
2726: for (h = 0; h < height+1; ++h) {
2727: PetscInt pStart, pEnd, i;
2729: DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);
2730: for (i = offsets[p*(height+2)+h]; i < closureSize; ++i) {
2731: if ((pStart > closures[p][i*2]) || (pEnd <= closures[p][i*2])) {
2732: offsets[p*(height+2)+h+1] = i;
2733: break;
2734: }
2735: }
2736: if (i == closureSize) offsets[p*(height+2)+h+1] = i;
2737: }
2738: 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);
2739: }
2740: for (h = 0; h < height+1; ++h) {
2741: PetscInt dof;
2743: /* Copy in cone of first point */
2744: dof = offsets[h+1] - offsets[h];
2745: for (meetSize = 0; meetSize < dof; ++meetSize) {
2746: meet[i][meetSize] = closures[0][(offsets[h]+meetSize)*2];
2747: }
2748: /* Check each successive cone */
2749: for (p = 1; p < numPoints && meetSize; ++p) {
2750: PetscInt newMeetSize = 0;
2752: dof = offsets[p*(height+2)+h+1] - offsets[p*(height+2)+h];
2753: for (c = 0; c < dof; ++c) {
2754: const PetscInt point = closures[p][(offsets[p*(height+2)+h]+c)*2];
2756: for (m = 0; m < meetSize; ++m) {
2757: if (point == meet[i][m]) {
2758: meet[1-i][newMeetSize++] = point;
2759: break;
2760: }
2761: }
2762: }
2763: meetSize = newMeetSize;
2764: i = 1-i;
2765: }
2766: if (meetSize) break;
2767: }
2768: *numCoveredPoints = meetSize;
2769: *coveredPoints = meet[i];
2770: for (p = 0; p < numPoints; ++p) {
2771: DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_TRUE, NULL, &closures[p]);
2772: }
2773: PetscFree(closures);
2774: DMRestoreWorkArray(dm, numPoints*(height+2), PETSC_INT, &offsets);
2775: DMRestoreWorkArray(dm, mesh->maxConeSize, PETSC_INT, &meet[1-i]);
2776: return(0);
2777: }
2779: /*@C
2780: DMPlexEqual - Determine if two DMs have the same topology
2782: Not Collective
2784: Input Parameters:
2785: + dmA - A DMPlex object
2786: - dmB - A DMPlex object
2788: Output Parameters:
2789: . equal - PETSC_TRUE if the topologies are identical
2791: Level: intermediate
2793: Notes:
2794: We are not solving graph isomorphism, so we do not permutation.
2796: .keywords: mesh
2797: .seealso: DMPlexGetCone()
2798: @*/
2799: PetscErrorCode DMPlexEqual(DM dmA, DM dmB, PetscBool *equal)
2800: {
2801: PetscInt depth, depthB, pStart, pEnd, pStartB, pEndB, p;
2809: *equal = PETSC_FALSE;
2810: DMPlexGetDepth(dmA, &depth);
2811: DMPlexGetDepth(dmB, &depthB);
2812: if (depth != depthB) return(0);
2813: DMPlexGetChart(dmA, &pStart, &pEnd);
2814: DMPlexGetChart(dmB, &pStartB, &pEndB);
2815: if ((pStart != pStartB) || (pEnd != pEndB)) return(0);
2816: for (p = pStart; p < pEnd; ++p) {
2817: const PetscInt *cone, *coneB, *ornt, *orntB, *support, *supportB;
2818: PetscInt coneSize, coneSizeB, c, supportSize, supportSizeB, s;
2820: DMPlexGetConeSize(dmA, p, &coneSize);
2821: DMPlexGetCone(dmA, p, &cone);
2822: DMPlexGetConeOrientation(dmA, p, &ornt);
2823: DMPlexGetConeSize(dmB, p, &coneSizeB);
2824: DMPlexGetCone(dmB, p, &coneB);
2825: DMPlexGetConeOrientation(dmB, p, &orntB);
2826: if (coneSize != coneSizeB) return(0);
2827: for (c = 0; c < coneSize; ++c) {
2828: if (cone[c] != coneB[c]) return(0);
2829: if (ornt[c] != orntB[c]) return(0);
2830: }
2831: DMPlexGetSupportSize(dmA, p, &supportSize);
2832: DMPlexGetSupport(dmA, p, &support);
2833: DMPlexGetSupportSize(dmB, p, &supportSizeB);
2834: DMPlexGetSupport(dmB, p, &supportB);
2835: if (supportSize != supportSizeB) return(0);
2836: for (s = 0; s < supportSize; ++s) {
2837: if (support[s] != supportB[s]) return(0);
2838: }
2839: }
2840: *equal = PETSC_TRUE;
2841: return(0);
2842: }
2844: /*@C
2845: DMPlexGetNumFaceVertices - Returns the number of vertices on a face
2847: Not Collective
2849: Input Parameters:
2850: + dm - The DMPlex
2851: . cellDim - The cell dimension
2852: - numCorners - The number of vertices on a cell
2854: Output Parameters:
2855: . numFaceVertices - The number of vertices on a face
2857: Level: developer
2859: Notes:
2860: Of course this can only work for a restricted set of symmetric shapes
2862: .seealso: DMPlexGetCone()
2863: @*/
2864: PetscErrorCode DMPlexGetNumFaceVertices(DM dm, PetscInt cellDim, PetscInt numCorners, PetscInt *numFaceVertices)
2865: {
2866: MPI_Comm comm;
2870: PetscObjectGetComm((PetscObject)dm,&comm);
2872: switch (cellDim) {
2873: case 0:
2874: *numFaceVertices = 0;
2875: break;
2876: case 1:
2877: *numFaceVertices = 1;
2878: break;
2879: case 2:
2880: switch (numCorners) {
2881: case 3: /* triangle */
2882: *numFaceVertices = 2; /* Edge has 2 vertices */
2883: break;
2884: case 4: /* quadrilateral */
2885: *numFaceVertices = 2; /* Edge has 2 vertices */
2886: break;
2887: case 6: /* quadratic triangle, tri and quad cohesive Lagrange cells */
2888: *numFaceVertices = 3; /* Edge has 3 vertices */
2889: break;
2890: case 9: /* quadratic quadrilateral, quadratic quad cohesive Lagrange cells */
2891: *numFaceVertices = 3; /* Edge has 3 vertices */
2892: break;
2893: default:
2894: SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of face corners %D for dimension %D", numCorners, cellDim);
2895: }
2896: break;
2897: case 3:
2898: switch (numCorners) {
2899: case 4: /* tetradehdron */
2900: *numFaceVertices = 3; /* Face has 3 vertices */
2901: break;
2902: case 6: /* tet cohesive cells */
2903: *numFaceVertices = 4; /* Face has 4 vertices */
2904: break;
2905: case 8: /* hexahedron */
2906: *numFaceVertices = 4; /* Face has 4 vertices */
2907: break;
2908: case 9: /* tet cohesive Lagrange cells */
2909: *numFaceVertices = 6; /* Face has 6 vertices */
2910: break;
2911: case 10: /* quadratic tetrahedron */
2912: *numFaceVertices = 6; /* Face has 6 vertices */
2913: break;
2914: case 12: /* hex cohesive Lagrange cells */
2915: *numFaceVertices = 6; /* Face has 6 vertices */
2916: break;
2917: case 18: /* quadratic tet cohesive Lagrange cells */
2918: *numFaceVertices = 6; /* Face has 6 vertices */
2919: break;
2920: case 27: /* quadratic hexahedron, quadratic hex cohesive Lagrange cells */
2921: *numFaceVertices = 9; /* Face has 9 vertices */
2922: break;
2923: default:
2924: SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of face corners %D for dimension %D", numCorners, cellDim);
2925: }
2926: break;
2927: default:
2928: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid cell dimension %D", cellDim);
2929: }
2930: return(0);
2931: }
2933: /*@
2934: DMPlexGetDepthLabel - Get the DMLabel recording the depth of each point
2936: Not Collective
2938: Input Parameter:
2939: . dm - The DMPlex object
2941: Output Parameter:
2942: . depthLabel - The DMLabel recording point depth
2944: Level: developer
2946: .keywords: mesh, points
2947: .seealso: DMPlexGetDepth(), DMPlexGetHeightStratum(), DMPlexGetDepthStratum()
2948: @*/
2949: PetscErrorCode DMPlexGetDepthLabel(DM dm, DMLabel *depthLabel)
2950: {
2956: if (!dm->depthLabel) {DMGetLabel(dm, "depth", &dm->depthLabel);}
2957: *depthLabel = dm->depthLabel;
2958: return(0);
2959: }
2961: /*@
2962: DMPlexGetDepth - Get the depth of the DAG representing this mesh
2964: Not Collective
2966: Input Parameter:
2967: . dm - The DMPlex object
2969: Output Parameter:
2970: . depth - The number of strata (breadth first levels) in the DAG
2972: Level: developer
2974: .keywords: mesh, points
2975: .seealso: DMPlexGetDepthLabel(), DMPlexGetHeightStratum(), DMPlexGetDepthStratum()
2976: @*/
2977: PetscErrorCode DMPlexGetDepth(DM dm, PetscInt *depth)
2978: {
2979: DMLabel label;
2980: PetscInt d = 0;
2986: DMPlexGetDepthLabel(dm, &label);
2987: if (label) {DMLabelGetNumValues(label, &d);}
2988: *depth = d-1;
2989: return(0);
2990: }
2992: /*@
2993: DMPlexGetDepthStratum - Get the bounds [start, end) for all points at a certain depth.
2995: Not Collective
2997: Input Parameters:
2998: + dm - The DMPlex object
2999: - stratumValue - The requested depth
3001: Output Parameters:
3002: + start - The first point at this depth
3003: - end - One beyond the last point at this depth
3005: Level: developer
3007: .keywords: mesh, points
3008: .seealso: DMPlexGetHeightStratum(), DMPlexGetDepth()
3009: @*/
3010: PetscErrorCode DMPlexGetDepthStratum(DM dm, PetscInt stratumValue, PetscInt *start, PetscInt *end)
3011: {
3012: DMLabel label;
3013: PetscInt pStart, pEnd;
3020: DMPlexGetChart(dm, &pStart, &pEnd);
3021: if (pStart == pEnd) return(0);
3022: if (stratumValue < 0) {
3023: if (start) *start = pStart;
3024: if (end) *end = pEnd;
3025: return(0);
3026: }
3027: DMPlexGetDepthLabel(dm, &label);
3028: if (!label) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "No label named depth was found");
3029: DMLabelGetStratumBounds(label, stratumValue, start, end);
3030: return(0);
3031: }
3033: /*@
3034: DMPlexGetHeightStratum - Get the bounds [start, end) for all points at a certain height.
3036: Not Collective
3038: Input Parameters:
3039: + dm - The DMPlex object
3040: - stratumValue - The requested height
3042: Output Parameters:
3043: + start - The first point at this height
3044: - end - One beyond the last point at this height
3046: Level: developer
3048: .keywords: mesh, points
3049: .seealso: DMPlexGetDepthStratum(), DMPlexGetDepth()
3050: @*/
3051: PetscErrorCode DMPlexGetHeightStratum(DM dm, PetscInt stratumValue, PetscInt *start, PetscInt *end)
3052: {
3053: DMLabel label;
3054: PetscInt depth, pStart, pEnd;
3061: DMPlexGetChart(dm, &pStart, &pEnd);
3062: if (pStart == pEnd) return(0);
3063: if (stratumValue < 0) {
3064: if (start) *start = pStart;
3065: if (end) *end = pEnd;
3066: return(0);
3067: }
3068: DMPlexGetDepthLabel(dm, &label);
3069: if (!label) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "No label named depth was found");
3070: DMLabelGetNumValues(label, &depth);
3071: DMLabelGetStratumBounds(label, depth-1-stratumValue, start, end);
3072: return(0);
3073: }
3075: /* Set the number of dof on each point and separate by fields */
3076: static PetscErrorCode DMPlexCreateSectionInitial(DM dm, PetscInt dim, PetscInt numFields,const PetscInt numComp[],const PetscInt numDof[], PetscSection *section)
3077: {
3078: PetscInt *pMax;
3079: PetscInt depth, pStart = 0, pEnd = 0;
3080: PetscInt Nf, p, d, dep, f;
3081: PetscBool *isFE;
3085: PetscMalloc1(numFields, &isFE);
3086: DMGetNumFields(dm, &Nf);
3087: for (f = 0; f < numFields; ++f) {
3088: PetscObject obj;
3089: PetscClassId id;
3091: isFE[f] = PETSC_FALSE;
3092: if (f >= Nf) continue;
3093: DMGetField(dm, f, &obj);
3094: PetscObjectGetClassId(obj, &id);
3095: if (id == PETSCFE_CLASSID) {isFE[f] = PETSC_TRUE;}
3096: else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
3097: }
3098: PetscSectionCreate(PetscObjectComm((PetscObject)dm), section);
3099: if (numFields > 0) {
3100: PetscSectionSetNumFields(*section, numFields);
3101: if (numComp) {
3102: for (f = 0; f < numFields; ++f) {
3103: PetscSectionSetFieldComponents(*section, f, numComp[f]);
3104: if (isFE[f]) {
3105: PetscFE fe;
3106: PetscDualSpace dspace;
3107: const PetscInt ***perms;
3108: const PetscScalar ***flips;
3109: const PetscInt *numDof;
3111: DMGetField(dm,f,(PetscObject *) &fe);
3112: PetscFEGetDualSpace(fe,&dspace);
3113: PetscDualSpaceGetSymmetries(dspace,&perms,&flips);
3114: PetscDualSpaceGetNumDof(dspace,&numDof);
3115: if (perms || flips) {
3116: DM K;
3117: DMLabel depthLabel;
3118: PetscInt depth, h;
3119: PetscSectionSym sym;
3121: PetscDualSpaceGetDM(dspace,&K);
3122: DMPlexGetDepthLabel(dm,&depthLabel);
3123: DMPlexGetDepth(dm,&depth);
3124: PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)*section),depthLabel,&sym);
3125: for (h = 0; h <= depth; h++) {
3126: PetscDualSpace hspace;
3127: PetscInt kStart, kEnd;
3128: PetscInt kConeSize;
3129: const PetscInt **perms0 = NULL;
3130: const PetscScalar **flips0 = NULL;
3132: PetscDualSpaceGetHeightSubspace(dspace,h,&hspace);
3133: DMPlexGetHeightStratum(K,h,&kStart,&kEnd);
3134: if (!hspace) continue;
3135: PetscDualSpaceGetSymmetries(hspace,&perms,&flips);
3136: if (perms) perms0 = perms[0];
3137: if (flips) flips0 = flips[0];
3138: if (!(perms0 || flips0)) continue;
3139: DMPlexGetConeSize(K,kStart,&kConeSize);
3140: PetscSectionSymLabelSetStratum(sym,depth - h,numDof[depth - h],-kConeSize,kConeSize,PETSC_USE_POINTER,perms0 ? &perms0[-kConeSize] : NULL,flips0 ? &flips0[-kConeSize] : NULL);
3141: }
3142: PetscSectionSetFieldSym(*section,f,sym);
3143: PetscSectionSymDestroy(&sym);
3144: }
3145: }
3146: }
3147: }
3148: }
3149: DMPlexGetChart(dm, &pStart, &pEnd);
3150: PetscSectionSetChart(*section, pStart, pEnd);
3151: DMPlexGetDepth(dm, &depth);
3152: PetscMalloc1(depth+1,&pMax);
3153: DMPlexGetHybridBounds(dm, depth >= 0 ? &pMax[depth] : NULL, depth>1 ? &pMax[depth-1] : NULL, depth>2 ? &pMax[1] : NULL, &pMax[0]);
3154: for (dep = 0; dep <= depth; ++dep) {
3155: d = dim == depth ? dep : (!dep ? 0 : dim);
3156: DMPlexGetDepthStratum(dm, dep, &pStart, &pEnd);
3157: pMax[dep] = pMax[dep] < 0 ? pEnd : pMax[dep];
3158: for (p = pStart; p < pEnd; ++p) {
3159: PetscInt tot = 0;
3161: for (f = 0; f < numFields; ++f) {
3162: if (isFE[f] && p >= pMax[dep]) continue;
3163: PetscSectionSetFieldDof(*section, p, f, numDof[f*(dim+1)+d]);
3164: tot += numDof[f*(dim+1)+d];
3165: }
3166: PetscSectionSetDof(*section, p, tot);
3167: }
3168: }
3169: PetscFree(pMax);
3170: PetscFree(isFE);
3171: return(0);
3172: }
3174: /* Set the number of dof on each point and separate by fields
3175: If bcComps is NULL or the IS is NULL, constrain every dof on the point
3176: */
3177: static PetscErrorCode DMPlexCreateSectionBCDof(DM dm, PetscInt numBC, const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
3178: {
3179: PetscInt numFields;
3180: PetscInt bc;
3181: PetscSection aSec;
3185: PetscSectionGetNumFields(section, &numFields);
3186: for (bc = 0; bc < numBC; ++bc) {
3187: PetscInt field = 0;
3188: const PetscInt *comp;
3189: const PetscInt *idx;
3190: PetscInt Nc = -1, n, i;
3192: if (numFields) field = bcField[bc];
3193: if (bcComps && bcComps[bc]) {ISGetLocalSize(bcComps[bc], &Nc);}
3194: if (bcComps && bcComps[bc]) {ISGetIndices(bcComps[bc], &comp);}
3195: ISGetLocalSize(bcPoints[bc], &n);
3196: ISGetIndices(bcPoints[bc], &idx);
3197: for (i = 0; i < n; ++i) {
3198: const PetscInt p = idx[i];
3199: PetscInt numConst;
3201: if (numFields) {
3202: PetscSectionGetFieldDof(section, p, field, &numConst);
3203: } else {
3204: PetscSectionGetDof(section, p, &numConst);
3205: }
3206: /* If Nc < 0, constrain every dof on the point */
3207: if (Nc > 0) numConst = PetscMin(numConst, Nc);
3208: if (numFields) {PetscSectionAddFieldConstraintDof(section, p, field, numConst);}
3209: PetscSectionAddConstraintDof(section, p, numConst);
3210: }
3211: ISRestoreIndices(bcPoints[bc], &idx);
3212: if (bcComps && bcComps[bc]) {ISRestoreIndices(bcComps[bc], &comp);}
3213: }
3214: DMPlexGetAnchors(dm, &aSec, NULL);
3215: if (aSec) {
3216: PetscInt aStart, aEnd, a;
3218: PetscSectionGetChart(aSec, &aStart, &aEnd);
3219: for (a = aStart; a < aEnd; a++) {
3220: PetscInt dof, f;
3222: PetscSectionGetDof(aSec, a, &dof);
3223: if (dof) {
3224: /* if there are point-to-point constraints, then all dofs are constrained */
3225: PetscSectionGetDof(section, a, &dof);
3226: PetscSectionSetConstraintDof(section, a, dof);
3227: for (f = 0; f < numFields; f++) {
3228: PetscSectionGetFieldDof(section, a, f, &dof);
3229: PetscSectionSetFieldConstraintDof(section, a, f, dof);
3230: }
3231: }
3232: }
3233: }
3234: return(0);
3235: }
3237: /* Set the constrained field indices on each point
3238: If bcComps is NULL or the IS is NULL, constrain every dof on the point
3239: */
3240: static PetscErrorCode DMPlexCreateSectionBCIndicesField(DM dm, PetscInt numBC,const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], PetscSection section)
3241: {
3242: PetscSection aSec;
3243: PetscInt *indices;
3244: PetscInt numFields, cdof, maxDof = 0, pStart, pEnd, p, bc, f, d;
3248: PetscSectionGetNumFields(section, &numFields);
3249: if (!numFields) return(0);
3250: /* Initialize all field indices to -1 */
3251: PetscSectionGetChart(section, &pStart, &pEnd);
3252: for (p = pStart; p < pEnd; ++p) {PetscSectionGetConstraintDof(section, p, &cdof); maxDof = PetscMax(maxDof, cdof);}
3253: PetscMalloc1(maxDof, &indices);
3254: for (d = 0; d < maxDof; ++d) indices[d] = -1;
3255: for (p = pStart; p < pEnd; ++p) for (f = 0; f < numFields; ++f) {PetscSectionSetFieldConstraintIndices(section, p, f, indices);}
3256: /* Handle BC constraints */
3257: for (bc = 0; bc < numBC; ++bc) {
3258: const PetscInt field = bcField[bc];
3259: const PetscInt *comp, *idx;
3260: PetscInt Nc = -1, n, i;
3262: if (bcComps && bcComps[bc]) {ISGetLocalSize(bcComps[bc], &Nc);}
3263: if (bcComps && bcComps[bc]) {ISGetIndices(bcComps[bc], &comp);}
3264: ISGetLocalSize(bcPoints[bc], &n);
3265: ISGetIndices(bcPoints[bc], &idx);
3266: for (i = 0; i < n; ++i) {
3267: const PetscInt p = idx[i];
3268: const PetscInt *find;
3269: PetscInt fdof, fcdof, c;
3271: PetscSectionGetFieldDof(section, p, field, &fdof);
3272: if (!fdof) continue;
3273: if (Nc < 0) {
3274: for (d = 0; d < fdof; ++d) indices[d] = d;
3275: fcdof = fdof;
3276: } else {
3277: PetscSectionGetFieldConstraintDof(section, p, field, &fcdof);
3278: PetscSectionGetFieldConstraintIndices(section, p, field, &find);
3279: for (d = 0; d < fcdof; ++d) {if (find[d] < 0) break; indices[d] = find[d];}
3280: for (c = 0; c < Nc; ++c) indices[d++] = comp[c];
3281: PetscSortRemoveDupsInt(&d, indices);
3282: for (c = d; c < fcdof; ++c) indices[c] = -1;
3283: fcdof = d;
3284: }
3285: PetscSectionSetFieldConstraintDof(section, p, field, fcdof);
3286: PetscSectionSetFieldConstraintIndices(section, p, field, indices);
3287: }
3288: if (bcComps && bcComps[bc]) {ISRestoreIndices(bcComps[bc], &comp);}
3289: ISRestoreIndices(bcPoints[bc], &idx);
3290: }
3291: /* Handle anchors */
3292: DMPlexGetAnchors(dm, &aSec, NULL);
3293: if (aSec) {
3294: PetscInt aStart, aEnd, a;
3296: for (d = 0; d < maxDof; ++d) indices[d] = d;
3297: PetscSectionGetChart(aSec, &aStart, &aEnd);
3298: for (a = aStart; a < aEnd; a++) {
3299: PetscInt dof, f;
3301: PetscSectionGetDof(aSec, a, &dof);
3302: if (dof) {
3303: /* if there are point-to-point constraints, then all dofs are constrained */
3304: for (f = 0; f < numFields; f++) {
3305: PetscSectionSetFieldConstraintIndices(section, a, f, indices);
3306: }
3307: }
3308: }
3309: }
3310: PetscFree(indices);
3311: return(0);
3312: }
3314: /* Set the constrained indices on each point */
3315: static PetscErrorCode DMPlexCreateSectionBCIndices(DM dm, PetscSection section)
3316: {
3317: PetscInt *indices;
3318: PetscInt numFields, maxDof, pStart, pEnd, p, f, d;
3322: PetscSectionGetNumFields(section, &numFields);
3323: PetscSectionGetMaxDof(section, &maxDof);
3324: PetscSectionGetChart(section, &pStart, &pEnd);
3325: PetscMalloc1(maxDof, &indices);
3326: for (d = 0; d < maxDof; ++d) indices[d] = -1;
3327: for (p = pStart; p < pEnd; ++p) {
3328: PetscInt cdof, d;
3330: PetscSectionGetConstraintDof(section, p, &cdof);
3331: if (cdof) {
3332: if (numFields) {
3333: PetscInt numConst = 0, foff = 0;
3335: for (f = 0; f < numFields; ++f) {
3336: const PetscInt *find;
3337: PetscInt fcdof, fdof;
3339: PetscSectionGetFieldDof(section, p, f, &fdof);
3340: PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);
3341: /* Change constraint numbering from field component to local dof number */
3342: PetscSectionGetFieldConstraintIndices(section, p, f, &find);
3343: for (d = 0; d < fcdof; ++d) indices[numConst+d] = find[d] + foff;
3344: numConst += fcdof;
3345: foff += fdof;
3346: }
3347: if (cdof != numConst) {PetscSectionSetConstraintDof(section, p, numConst);}
3348: } else {
3349: for (d = 0; d < cdof; ++d) indices[d] = d;
3350: }
3351: PetscSectionSetConstraintIndices(section, p, indices);
3352: }
3353: }
3354: PetscFree(indices);
3355: return(0);
3356: }
3358: /*@C
3359: DMPlexCreateSection - Create a PetscSection based upon the dof layout specification provided.
3361: Not Collective
3363: Input Parameters:
3364: + dm - The DMPlex object
3365: . dim - The spatial dimension of the problem
3366: . numFields - The number of fields in the problem
3367: . numComp - An array of size numFields that holds the number of components for each field
3368: . numDof - An array of size numFields*(dim+1) which holds the number of dof for each field on a mesh piece of dimension d
3369: . numBC - The number of boundary conditions
3370: . bcField - An array of size numBC giving the field number for each boundry condition
3371: . bcComps - [Optional] An array of size numBC giving an IS holding the field components to which each boundary condition applies
3372: . bcPoints - An array of size numBC giving an IS holding the Plex points to which each boundary condition applies
3373: - perm - Optional permutation of the chart, or NULL
3375: Output Parameter:
3376: . section - The PetscSection object
3378: Notes: numDof[f*(dim+1)+d] gives the number of dof for field f on points of dimension d. For instance, numDof[1] is the
3379: number of dof for field 0 on each edge.
3381: The chart permutation is the same one set using PetscSectionSetPermutation()
3383: Level: developer
3385: Fortran Notes:
3386: A Fortran 90 version is available as DMPlexCreateSectionF90()
3388: .keywords: mesh, elements
3389: .seealso: DMPlexCreate(), PetscSectionCreate(), PetscSectionSetPermutation()
3390: @*/
3391: PetscErrorCode DMPlexCreateSection(DM dm, PetscInt dim, PetscInt numFields,const PetscInt numComp[],const PetscInt numDof[], PetscInt numBC,const PetscInt bcField[], const IS bcComps[], const IS bcPoints[], IS perm, PetscSection *section)
3392: {
3393: PetscSection aSec;
3397: DMPlexCreateSectionInitial(dm, dim, numFields, numComp, numDof, section);
3398: DMPlexCreateSectionBCDof(dm, numBC, bcField, bcComps, bcPoints, *section);
3399: if (perm) {PetscSectionSetPermutation(*section, perm);}
3400: PetscSectionSetUp(*section);
3401: DMPlexGetAnchors(dm,&aSec,NULL);
3402: if (numBC || aSec) {
3403: DMPlexCreateSectionBCIndicesField(dm, numBC, bcField, bcComps, bcPoints, *section);
3404: DMPlexCreateSectionBCIndices(dm, *section);
3405: }
3406: PetscSectionViewFromOptions(*section,NULL,"-section_view");
3407: return(0);
3408: }
3410: PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm)
3411: {
3412: PetscSection section, s;
3413: Mat m;
3414: PetscInt maxHeight;
3418: DMClone(dm, cdm);
3419: DMPlexGetMaxProjectionHeight(dm, &maxHeight);
3420: DMPlexSetMaxProjectionHeight(*cdm, maxHeight);
3421: PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion);
3422: DMSetDefaultSection(*cdm, section);
3423: PetscSectionDestroy(§ion);
3424: PetscSectionCreate(PETSC_COMM_SELF, &s);
3425: MatCreate(PETSC_COMM_SELF, &m);
3426: DMSetDefaultConstraints(*cdm, s, m);
3427: PetscSectionDestroy(&s);
3428: MatDestroy(&m);
3429: return(0);
3430: }
3432: /*@C
3433: DMPlexGetConeSection - Return a section which describes the layout of cone data
3435: Not Collective
3437: Input Parameters:
3438: . dm - The DMPlex object
3440: Output Parameter:
3441: . section - The PetscSection object
3443: Level: developer
3445: .seealso: DMPlexGetSupportSection(), DMPlexGetCones(), DMPlexGetConeOrientations()
3446: @*/
3447: PetscErrorCode DMPlexGetConeSection(DM dm, PetscSection *section)
3448: {
3449: DM_Plex *mesh = (DM_Plex*) dm->data;
3453: if (section) *section = mesh->coneSection;
3454: return(0);
3455: }
3457: /*@C
3458: DMPlexGetSupportSection - Return a section which describes the layout of support data
3460: Not Collective
3462: Input Parameters:
3463: . dm - The DMPlex object
3465: Output Parameter:
3466: . section - The PetscSection object
3468: Level: developer
3470: .seealso: DMPlexGetConeSection()
3471: @*/
3472: PetscErrorCode DMPlexGetSupportSection(DM dm, PetscSection *section)
3473: {
3474: DM_Plex *mesh = (DM_Plex*) dm->data;
3478: if (section) *section = mesh->supportSection;
3479: return(0);
3480: }
3482: /*@C
3483: DMPlexGetCones - Return cone data
3485: Not Collective
3487: Input Parameters:
3488: . dm - The DMPlex object
3490: Output Parameter:
3491: . cones - The cone for each point
3493: Level: developer
3495: .seealso: DMPlexGetConeSection()
3496: @*/
3497: PetscErrorCode DMPlexGetCones(DM dm, PetscInt *cones[])
3498: {
3499: DM_Plex *mesh = (DM_Plex*) dm->data;
3503: if (cones) *cones = mesh->cones;
3504: return(0);
3505: }
3507: /*@C
3508: DMPlexGetConeOrientations - Return cone orientation data
3510: Not Collective
3512: Input Parameters:
3513: . dm - The DMPlex object
3515: Output Parameter:
3516: . coneOrientations - The cone orientation for each point
3518: Level: developer
3520: .seealso: DMPlexGetConeSection()
3521: @*/
3522: PetscErrorCode DMPlexGetConeOrientations(DM dm, PetscInt *coneOrientations[])
3523: {
3524: DM_Plex *mesh = (DM_Plex*) dm->data;
3528: if (coneOrientations) *coneOrientations = mesh->coneOrientations;
3529: return(0);
3530: }
3532: /******************************** FEM Support **********************************/
3534: PetscErrorCode DMPlexCreateSpectralClosurePermutation(DM dm, PetscSection section)
3535: {
3536: PetscInt *perm;
3537: PetscInt dim, eStart, k, Nf, f, Nc, c, i, j, size = 0, offset = 0, foffset = 0;
3541: if (!section) {DMGetDefaultSection(dm, §ion);}
3542: DMGetDimension(dm, &dim);
3543: PetscSectionGetNumFields(section, &Nf);
3544: if (dim <= 1) return(0);
3545: for (f = 0; f < Nf; ++f) {
3546: /* An order k SEM disc has k-1 dofs on an edge */
3547: DMPlexGetDepthStratum(dm, 1, &eStart, NULL);
3548: PetscSectionGetFieldDof(section, eStart, f, &k);
3549: PetscSectionGetFieldComponents(section, f, &Nc);
3550: k = k/Nc + 1;
3551: size += PetscPowInt(k+1, dim)*Nc;
3552: }
3553: PetscMalloc1(size, &perm);
3554: for (f = 0; f < Nf; ++f) {
3555: switch (dim) {
3556: case 2:
3557: /* The original quad closure is oriented clockwise, {f, e_b, e_r, e_t, e_l, v_lb, v_rb, v_tr, v_tl} */
3558: DMPlexGetDepthStratum(dm, 1, &eStart, NULL);
3559: PetscSectionGetFieldDof(section, eStart, f, &k);
3560: PetscSectionGetFieldComponents(section, f, &Nc);
3561: k = k/Nc + 1;
3562: /* The SEM order is
3564: v_lb, {e_b}, v_rb,
3565: e^{(k-1)-i}_l, {f^{i*(k-1)}}, e^i_r,
3566: v_lt, reverse {e_t}, v_rt
3567: */
3568: {
3569: const PetscInt of = 0;
3570: const PetscInt oeb = of + PetscSqr(k-1);
3571: const PetscInt oer = oeb + (k-1);
3572: const PetscInt oet = oer + (k-1);
3573: const PetscInt oel = oet + (k-1);
3574: const PetscInt ovlb = oel + (k-1);
3575: const PetscInt ovrb = ovlb + 1;
3576: const PetscInt ovrt = ovrb + 1;
3577: const PetscInt ovlt = ovrt + 1;
3578: PetscInt o;
3580: /* bottom */
3581: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovlb*Nc + c + foffset;
3582: for (o = oeb; o < oer; ++o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
3583: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovrb*Nc + c + foffset;
3584: /* middle */
3585: for (i = 0; i < k-1; ++i) {
3586: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oel+(k-2)-i)*Nc + c + foffset;
3587: 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;
3588: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oer+i)*Nc + c + foffset;
3589: }
3590: /* top */
3591: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovlt*Nc + c + foffset;
3592: for (o = oel-1; o >= oet; --o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
3593: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovrt*Nc + c + foffset;
3594: foffset = offset;
3595: }
3596: break;
3597: case 3:
3598: /* The original hex closure is
3600: {c,
3601: f_b, f_t, f_f, f_b, f_r, f_l,
3602: e_bl, e_bb, e_br, e_bf, e_tf, e_tr, e_tb, e_tl, e_rf, e_lf, e_lb, e_rb,
3603: v_blf, v_blb, v_brb, v_brf, v_tlf, v_trf, v_trb, v_tlb}
3604: */
3605: DMPlexGetDepthStratum(dm, 1, &eStart, NULL);
3606: PetscSectionGetFieldDof(section, eStart, f, &k);
3607: PetscSectionGetFieldComponents(section, f, &Nc);
3608: k = k/Nc + 1;
3609: /* The SEM order is
3610: Bottom Slice
3611: v_blf, {e^{(k-1)-n}_bf}, v_brf,
3612: e^{i}_bl, f^{n*(k-1)+(k-1)-i}_b, e^{(k-1)-i}_br,
3613: v_blb, {e_bb}, v_brb,
3615: Middle Slice (j)
3616: {e^{(k-1)-j}_lf}, {f^{j*(k-1)+n}_f}, e^j_rf,
3617: f^{i*(k-1)+j}_l, {c^{(j*(k-1) + i)*(k-1)+n}_t}, f^{j*(k-1)+i}_r,
3618: e^j_lb, {f^{j*(k-1)+(k-1)-n}_b}, e^{(k-1)-j}_rb,
3620: Top Slice
3621: v_tlf, {e_tf}, v_trf,
3622: e^{(k-1)-i}_tl, {f^{i*(k-1)}_t}, e^{i}_tr,
3623: v_tlb, {e^{(k-1)-n}_tb}, v_trb,
3624: */
3625: {
3626: const PetscInt oc = 0;
3627: const PetscInt ofb = oc + PetscSqr(k-1)*(k-1);
3628: const PetscInt oft = ofb + PetscSqr(k-1);
3629: const PetscInt off = oft + PetscSqr(k-1);
3630: const PetscInt ofk = off + PetscSqr(k-1);
3631: const PetscInt ofr = ofk + PetscSqr(k-1);
3632: const PetscInt ofl = ofr + PetscSqr(k-1);
3633: const PetscInt oebl = ofl + PetscSqr(k-1);
3634: const PetscInt oebb = oebl + (k-1);
3635: const PetscInt oebr = oebb + (k-1);
3636: const PetscInt oebf = oebr + (k-1);
3637: const PetscInt oetf = oebf + (k-1);
3638: const PetscInt oetr = oetf + (k-1);
3639: const PetscInt oetb = oetr + (k-1);
3640: const PetscInt oetl = oetb + (k-1);
3641: const PetscInt oerf = oetl + (k-1);
3642: const PetscInt oelf = oerf + (k-1);
3643: const PetscInt oelb = oelf + (k-1);
3644: const PetscInt oerb = oelb + (k-1);
3645: const PetscInt ovblf = oerb + (k-1);
3646: const PetscInt ovblb = ovblf + 1;
3647: const PetscInt ovbrb = ovblb + 1;
3648: const PetscInt ovbrf = ovbrb + 1;
3649: const PetscInt ovtlf = ovbrf + 1;
3650: const PetscInt ovtrf = ovtlf + 1;
3651: const PetscInt ovtrb = ovtrf + 1;
3652: const PetscInt ovtlb = ovtrb + 1;
3653: PetscInt o, n;
3655: /* Bottom Slice */
3656: /* bottom */
3657: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovblf*Nc + c + foffset;
3658: for (o = oetf-1; o >= oebf; --o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
3659: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovbrf*Nc + c + foffset;
3660: /* middle */
3661: for (i = 0; i < k-1; ++i) {
3662: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oebl+i)*Nc + c + foffset;
3663: 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;}
3664: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oebr+(k-2)-i)*Nc + c + foffset;
3665: }
3666: /* top */
3667: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovblb*Nc + c + foffset;
3668: for (o = oebb; o < oebr; ++o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
3669: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovbrb*Nc + c + foffset;
3671: /* Middle Slice */
3672: for (j = 0; j < k-1; ++j) {
3673: /* bottom */
3674: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oelf+(k-2)-j)*Nc + c + foffset;
3675: 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;
3676: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oerf+j)*Nc + c + foffset;
3677: /* middle */
3678: for (i = 0; i < k-1; ++i) {
3679: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (ofl+i*(k-1)+j)*Nc + c + foffset;
3680: 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;
3681: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (ofr+j*(k-1)+i)*Nc + c + foffset;
3682: }
3683: /* top */
3684: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oelb+j)*Nc + c + foffset;
3685: 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;
3686: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oerb+(k-2)-j)*Nc + c + foffset;
3687: }
3689: /* Top Slice */
3690: /* bottom */
3691: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtlf*Nc + c + foffset;
3692: for (o = oetf; o < oetr; ++o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
3693: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtrf*Nc + c + foffset;
3694: /* middle */
3695: for (i = 0; i < k-1; ++i) {
3696: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oetl+(k-2)-i)*Nc + c + foffset;
3697: for (n = 0; n < k-1; ++n) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oft+i*(k-1)+n)*Nc + c + foffset;
3698: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = (oetr+i)*Nc + c + foffset;
3699: }
3700: /* top */
3701: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtlb*Nc + c + foffset;
3702: for (o = oetl-1; o >= oetb; --o) for (c = 0; c < Nc; ++c, ++offset) perm[offset] = o*Nc + c + foffset;
3703: for (c = 0; c < Nc; ++c, ++offset) perm[offset] = ovtrb*Nc + c + foffset;
3705: foffset = offset;
3706: }
3707: break;
3708: default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No spectral ordering for dimension %D", dim);
3709: }
3710: }
3711: if (offset != size) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Number of permutation entries %D != %D", offset, size);
3712: /* Check permutation */
3713: {
3714: PetscInt *check;
3716: PetscMalloc1(size, &check);
3717: 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]);}
3718: for (i = 0; i < size; ++i) check[perm[i]] = i;
3719: for (i = 0; i < size; ++i) {if (check[i] < 0) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Missing permutation index %D", i);}
3720: PetscFree(check);
3721: }
3722: PetscSectionSetClosurePermutation_Internal(section, (PetscObject) dm, size, PETSC_OWN_POINTER, perm);
3723: return(0);
3724: }
3726: PetscErrorCode DMPlexGetPointDualSpaceFEM(DM dm, PetscInt point, PetscInt field, PetscDualSpace *dspace)
3727: {
3728: PetscDS prob;
3729: PetscInt depth, Nf, h;
3730: DMLabel label;
3734: prob = dm->prob;
3735: Nf = prob->Nf;
3736: label = dm->depthLabel;
3737: *dspace = NULL;
3738: if (field < Nf) {
3739: PetscObject disc = prob->disc[field];
3741: if (disc->classid == PETSCFE_CLASSID) {
3742: PetscDualSpace dsp;
3744: PetscFEGetDualSpace((PetscFE)disc,&dsp);
3745: DMLabelGetNumValues(label,&depth);
3746: DMLabelGetValue(label,point,&h);
3747: h = depth - 1 - h;
3748: if (h) {
3749: PetscDualSpaceGetHeightSubspace(dsp,h,dspace);
3750: } else {
3751: *dspace = dsp;
3752: }
3753: }
3754: }
3755: return(0);
3756: }
3759: PETSC_STATIC_INLINE PetscErrorCode DMPlexVecGetClosure_Depth1_Static(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt *csize, PetscScalar *values[])
3760: {
3761: PetscScalar *array, *vArray;
3762: const PetscInt *cone, *coneO;
3763: PetscInt pStart, pEnd, p, numPoints, size = 0, offset = 0;
3764: PetscErrorCode ierr;
3767: PetscSectionGetChart(section, &pStart, &pEnd);
3768: DMPlexGetConeSize(dm, point, &numPoints);
3769: DMPlexGetCone(dm, point, &cone);
3770: DMPlexGetConeOrientation(dm, point, &coneO);
3771: if (!values || !*values) {
3772: if ((point >= pStart) && (point < pEnd)) {
3773: PetscInt dof;
3775: PetscSectionGetDof(section, point, &dof);
3776: size += dof;
3777: }
3778: for (p = 0; p < numPoints; ++p) {
3779: const PetscInt cp = cone[p];
3780: PetscInt dof;
3782: if ((cp < pStart) || (cp >= pEnd)) continue;
3783: PetscSectionGetDof(section, cp, &dof);
3784: size += dof;
3785: }
3786: if (!values) {
3787: if (csize) *csize = size;
3788: return(0);
3789: }
3790: DMGetWorkArray(dm, size, PETSC_SCALAR, &array);
3791: } else {
3792: array = *values;
3793: }
3794: size = 0;
3795: VecGetArray(v, &vArray);
3796: if ((point >= pStart) && (point < pEnd)) {
3797: PetscInt dof, off, d;
3798: PetscScalar *varr;
3800: PetscSectionGetDof(section, point, &dof);
3801: PetscSectionGetOffset(section, point, &off);
3802: varr = &vArray[off];
3803: for (d = 0; d < dof; ++d, ++offset) {
3804: array[offset] = varr[d];
3805: }
3806: size += dof;
3807: }
3808: for (p = 0; p < numPoints; ++p) {
3809: const PetscInt cp = cone[p];
3810: PetscInt o = coneO[p];
3811: PetscInt dof, off, d;
3812: PetscScalar *varr;
3814: if ((cp < pStart) || (cp >= pEnd)) continue;
3815: PetscSectionGetDof(section, cp, &dof);
3816: PetscSectionGetOffset(section, cp, &off);
3817: varr = &vArray[off];
3818: if (o >= 0) {
3819: for (d = 0; d < dof; ++d, ++offset) {
3820: array[offset] = varr[d];
3821: }
3822: } else {
3823: for (d = dof-1; d >= 0; --d, ++offset) {
3824: array[offset] = varr[d];
3825: }
3826: }
3827: size += dof;
3828: }
3829: VecRestoreArray(v, &vArray);
3830: if (!*values) {
3831: if (csize) *csize = size;
3832: *values = array;
3833: } else {
3834: if (size > *csize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Size of input array %D < actual size %D", *csize, size);
3835: *csize = size;
3836: }
3837: return(0);
3838: }
3840: static PetscErrorCode DMPlexGetCompressedClosure(DM dm, PetscSection section, PetscInt point, PetscInt *numPoints, PetscInt **points, PetscSection *clSec, IS *clPoints, const PetscInt **clp)
3841: {
3842: const PetscInt *cla;
3843: PetscInt np, *pts = NULL;
3847: PetscSectionGetClosureIndex(section, (PetscObject) dm, clSec, clPoints);
3848: if (!*clPoints) {
3849: PetscInt pStart, pEnd, p, q;
3851: PetscSectionGetChart(section, &pStart, &pEnd);
3852: DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &np, &pts);
3853: /* Compress out points not in the section */
3854: for (p = 0, q = 0; p < np; p++) {
3855: PetscInt r = pts[2*p];
3856: if ((r >= pStart) && (r < pEnd)) {
3857: pts[q*2] = r;
3858: pts[q*2+1] = pts[2*p+1];
3859: ++q;
3860: }
3861: }
3862: np = q;
3863: cla = NULL;
3864: } else {
3865: PetscInt dof, off;
3867: PetscSectionGetDof(*clSec, point, &dof);
3868: PetscSectionGetOffset(*clSec, point, &off);
3869: ISGetIndices(*clPoints, &cla);
3870: np = dof/2;
3871: pts = (PetscInt *) &cla[off];
3872: }
3873: *numPoints = np;
3874: *points = pts;
3875: *clp = cla;
3877: return(0);
3878: }
3880: static PetscErrorCode DMPlexRestoreCompressedClosure(DM dm, PetscSection section, PetscInt point, PetscInt *numPoints, PetscInt **points, PetscSection *clSec, IS *clPoints, const PetscInt **clp)
3881: {
3885: if (!*clPoints) {
3886: DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, numPoints, points);
3887: } else {
3888: ISRestoreIndices(*clPoints, clp);
3889: }
3890: *numPoints = 0;
3891: *points = NULL;
3892: *clSec = NULL;
3893: *clPoints = NULL;
3894: *clp = NULL;
3895: return(0);
3896: }
3898: 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[])
3899: {
3900: PetscInt offset = 0, p;
3901: const PetscInt **perms = NULL;
3902: const PetscScalar **flips = NULL;
3903: PetscErrorCode ierr;
3906: *size = 0;
3907: PetscSectionGetPointSyms(section,numPoints,points,&perms,&flips);
3908: for (p = 0; p < numPoints; p++) {
3909: const PetscInt point = points[2*p];
3910: const PetscInt *perm = perms ? perms[p] : NULL;
3911: const PetscScalar *flip = flips ? flips[p] : NULL;
3912: PetscInt dof, off, d;
3913: const PetscScalar *varr;
3915: PetscSectionGetDof(section, point, &dof);
3916: PetscSectionGetOffset(section, point, &off);
3917: varr = &vArray[off];
3918: if (clperm) {
3919: if (perm) {
3920: for (d = 0; d < dof; d++) array[clperm[offset + perm[d]]] = varr[d];
3921: } else {
3922: for (d = 0; d < dof; d++) array[clperm[offset + d ]] = varr[d];
3923: }
3924: if (flip) {
3925: for (d = 0; d < dof; d++) array[clperm[offset + d ]] *= flip[d];
3926: }
3927: } else {
3928: if (perm) {
3929: for (d = 0; d < dof; d++) array[offset + perm[d]] = varr[d];
3930: } else {
3931: for (d = 0; d < dof; d++) array[offset + d ] = varr[d];
3932: }
3933: if (flip) {
3934: for (d = 0; d < dof; d++) array[offset + d ] *= flip[d];
3935: }
3936: }
3937: offset += dof;
3938: }
3939: PetscSectionRestorePointSyms(section,numPoints,points,&perms,&flips);
3940: *size = offset;
3941: return(0);
3942: }
3944: 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[])
3945: {
3946: PetscInt offset = 0, f;
3947: PetscErrorCode ierr;
3950: *size = 0;
3951: for (f = 0; f < numFields; ++f) {
3952: PetscInt p;
3953: const PetscInt **perms = NULL;
3954: const PetscScalar **flips = NULL;
3956: PetscSectionGetFieldPointSyms(section,f,numPoints,points,&perms,&flips);
3957: for (p = 0; p < numPoints; p++) {
3958: const PetscInt point = points[2*p];
3959: PetscInt fdof, foff, b;
3960: const PetscScalar *varr;
3961: const PetscInt *perm = perms ? perms[p] : NULL;
3962: const PetscScalar *flip = flips ? flips[p] : NULL;
3964: PetscSectionGetFieldDof(section, point, f, &fdof);
3965: PetscSectionGetFieldOffset(section, point, f, &foff);
3966: varr = &vArray[foff];
3967: if (clperm) {
3968: if (perm) {for (b = 0; b < fdof; b++) {array[clperm[offset + perm[b]]] = varr[b];}}
3969: else {for (b = 0; b < fdof; b++) {array[clperm[offset + b ]] = varr[b];}}
3970: if (flip) {for (b = 0; b < fdof; b++) {array[clperm[offset + b ]] *= flip[b];}}
3971: } else {
3972: if (perm) {for (b = 0; b < fdof; b++) {array[offset + perm[b]] = varr[b];}}
3973: else {for (b = 0; b < fdof; b++) {array[offset + b ] = varr[b];}}
3974: if (flip) {for (b = 0; b < fdof; b++) {array[offset + b ] *= flip[b];}}
3975: }
3976: offset += fdof;
3977: }
3978: PetscSectionRestoreFieldPointSyms(section,f,numPoints,points,&perms,&flips);
3979: }
3980: *size = offset;
3981: return(0);
3982: }
3984: /*@C
3985: DMPlexVecGetClosure - Get an array of the values on the closure of 'point'
3987: Not collective
3989: Input Parameters:
3990: + dm - The DM
3991: . section - The section describing the layout in v, or NULL to use the default section
3992: . v - The local vector
3993: - point - The point in the DM
3995: Output Parameters:
3996: + csize - The number of values in the closure, or NULL
3997: - values - The array of values, which is a borrowed array and should not be freed
3999: Fortran Notes:
4000: Since it returns an array, this routine is only available in Fortran 90, and you must
4001: include petsc.h90 in your code.
4003: The csize argument is not present in the Fortran 90 binding since it is internal to the array.
4005: Level: intermediate
4007: .seealso DMPlexVecRestoreClosure(), DMPlexVecSetClosure(), DMPlexMatSetClosure()
4008: @*/
4009: PetscErrorCode DMPlexVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt *csize, PetscScalar *values[])
4010: {
4011: PetscSection clSection;
4012: IS clPoints;
4013: PetscScalar *array;
4014: const PetscScalar *vArray;
4015: PetscInt *points = NULL;
4016: const PetscInt *clp, *perm;
4017: PetscInt depth, numFields, numPoints, size;
4018: PetscErrorCode ierr;
4022: if (!section) {DMGetDefaultSection(dm, §ion);}
4025: DMPlexGetDepth(dm, &depth);
4026: PetscSectionGetNumFields(section, &numFields);
4027: if (depth == 1 && numFields < 2) {
4028: DMPlexVecGetClosure_Depth1_Static(dm, section, v, point, csize, values);
4029: return(0);
4030: }
4031: /* Get points */
4032: DMPlexGetCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
4033: PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject) dm, NULL, &perm);
4034: /* Get array */
4035: if (!values || !*values) {
4036: PetscInt asize = 0, dof, p;
4038: for (p = 0; p < numPoints*2; p += 2) {
4039: PetscSectionGetDof(section, points[p], &dof);
4040: asize += dof;
4041: }
4042: if (!values) {
4043: DMPlexRestoreCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
4044: if (csize) *csize = asize;
4045: return(0);
4046: }
4047: DMGetWorkArray(dm, asize, PETSC_SCALAR, &array);
4048: } else {
4049: array = *values;
4050: }
4051: VecGetArrayRead(v, &vArray);
4052: /* Get values */
4053: if (numFields > 0) {DMPlexVecGetClosure_Fields_Static(dm, section, numPoints, points, numFields, perm, vArray, &size, array);}
4054: else {DMPlexVecGetClosure_Static(dm, section, numPoints, points, perm, vArray, &size, array);}
4055: /* Cleanup points */
4056: DMPlexRestoreCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
4057: /* Cleanup array */
4058: VecRestoreArrayRead(v, &vArray);
4059: if (!*values) {
4060: if (csize) *csize = size;
4061: *values = array;
4062: } else {
4063: if (size > *csize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Size of input array %D < actual size %D", *csize, size);
4064: *csize = size;
4065: }
4066: return(0);
4067: }
4069: /*@C
4070: DMPlexVecRestoreClosure - Restore the array of the values on the closure of 'point'
4072: Not collective
4074: Input Parameters:
4075: + dm - The DM
4076: . section - The section describing the layout in v, or NULL to use the default section
4077: . v - The local vector
4078: . point - The point in the DM
4079: . csize - The number of values in the closure, or NULL
4080: - values - The array of values, which is a borrowed array and should not be freed
4082: Fortran Notes:
4083: Since it returns an array, this routine is only available in Fortran 90, and you must
4084: include petsc.h90 in your code.
4086: The csize argument is not present in the Fortran 90 binding since it is internal to the array.
4088: Level: intermediate
4090: .seealso DMPlexVecGetClosure(), DMPlexVecSetClosure(), DMPlexMatSetClosure()
4091: @*/
4092: PetscErrorCode DMPlexVecRestoreClosure(DM dm, PetscSection section, Vec v, PetscInt point, PetscInt *csize, PetscScalar *values[])
4093: {
4094: PetscInt size = 0;
4098: /* Should work without recalculating size */
4099: DMRestoreWorkArray(dm, size, PETSC_SCALAR, (void*) values);
4100: return(0);
4101: }
4103: PETSC_STATIC_INLINE void add (PetscScalar *x, PetscScalar y) {*x += y;}
4104: PETSC_STATIC_INLINE void insert(PetscScalar *x, PetscScalar y) {*x = y;}
4106: 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[])
4107: {
4108: PetscInt cdof; /* The number of constraints on this point */
4109: const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
4110: PetscScalar *a;
4111: PetscInt off, cind = 0, k;
4112: PetscErrorCode ierr;
4115: PetscSectionGetConstraintDof(section, point, &cdof);
4116: PetscSectionGetOffset(section, point, &off);
4117: a = &array[off];
4118: if (!cdof || setBC) {
4119: if (clperm) {
4120: if (perm) {for (k = 0; k < dof; ++k) {fuse(&a[k], values[clperm[offset+perm[k]]] * (flip ? flip[perm[k]] : 1.));}}
4121: else {for (k = 0; k < dof; ++k) {fuse(&a[k], values[clperm[offset+ k ]] * (flip ? flip[ k ] : 1.));}}
4122: } else {
4123: if (perm) {for (k = 0; k < dof; ++k) {fuse(&a[k], values[offset+perm[k]] * (flip ? flip[perm[k]] : 1.));}}
4124: else {for (k = 0; k < dof; ++k) {fuse(&a[k], values[offset+ k ] * (flip ? flip[ k ] : 1.));}}
4125: }
4126: } else {
4127: PetscSectionGetConstraintIndices(section, point, &cdofs);
4128: if (clperm) {
4129: if (perm) {for (k = 0; k < dof; ++k) {
4130: if ((cind < cdof) && (k == cdofs[cind])) {++cind; continue;}
4131: fuse(&a[k], values[clperm[offset+perm[k]]] * (flip ? flip[perm[k]] : 1.));
4132: }
4133: } else {
4134: for (k = 0; k < dof; ++k) {
4135: if ((cind < cdof) && (k == cdofs[cind])) {++cind; continue;}
4136: fuse(&a[k], values[clperm[offset+ k ]] * (flip ? flip[ k ] : 1.));
4137: }
4138: }
4139: } else {
4140: if (perm) {
4141: for (k = 0; k < dof; ++k) {
4142: if ((cind < cdof) && (k == cdofs[cind])) {++cind; continue;}
4143: fuse(&a[k], values[offset+perm[k]] * (flip ? flip[perm[k]] : 1.));
4144: }
4145: } else {
4146: for (k = 0; k < dof; ++k) {
4147: if ((cind < cdof) && (k == cdofs[cind])) {++cind; continue;}
4148: fuse(&a[k], values[offset+ k ] * (flip ? flip[ k ] : 1.));
4149: }
4150: }
4151: }
4152: }
4153: return(0);
4154: }
4156: 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[])
4157: {
4158: PetscInt cdof; /* The number of constraints on this point */
4159: const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
4160: PetscScalar *a;
4161: PetscInt off, cind = 0, k;
4162: PetscErrorCode ierr;
4165: PetscSectionGetConstraintDof(section, point, &cdof);
4166: PetscSectionGetOffset(section, point, &off);
4167: a = &array[off];
4168: if (cdof) {
4169: PetscSectionGetConstraintIndices(section, point, &cdofs);
4170: if (clperm) {
4171: if (perm) {
4172: for (k = 0; k < dof; ++k) {
4173: if ((cind < cdof) && (k == cdofs[cind])) {
4174: fuse(&a[k], values[clperm[offset+perm[k]]] * (flip ? flip[perm[k]] : 1.));
4175: cind++;
4176: }
4177: }
4178: } else {
4179: for (k = 0; k < dof; ++k) {
4180: if ((cind < cdof) && (k == cdofs[cind])) {
4181: fuse(&a[k], values[clperm[offset+ k ]] * (flip ? flip[ k ] : 1.));
4182: cind++;
4183: }
4184: }
4185: }
4186: } else {
4187: if (perm) {
4188: for (k = 0; k < dof; ++k) {
4189: if ((cind < cdof) && (k == cdofs[cind])) {
4190: fuse(&a[k], values[offset+perm[k]] * (flip ? flip[perm[k]] : 1.));
4191: cind++;
4192: }
4193: }
4194: } else {
4195: for (k = 0; k < dof; ++k) {
4196: if ((cind < cdof) && (k == cdofs[cind])) {
4197: fuse(&a[k], values[offset+ k ] * (flip ? flip[ k ] : 1.));
4198: cind++;
4199: }
4200: }
4201: }
4202: }
4203: }
4204: return(0);
4205: }
4207: 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[])
4208: {
4209: PetscScalar *a;
4210: PetscInt fdof, foff, fcdof, foffset = *offset;
4211: const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
4212: PetscInt cind = 0, b;
4213: PetscErrorCode ierr;
4216: PetscSectionGetFieldDof(section, point, f, &fdof);
4217: PetscSectionGetFieldConstraintDof(section, point, f, &fcdof);
4218: PetscSectionGetFieldOffset(section, point, f, &foff);
4219: a = &array[foff];
4220: if (!fcdof || setBC) {
4221: if (clperm) {
4222: if (perm) {for (b = 0; b < fdof; b++) {fuse(&a[b], values[clperm[foffset+perm[b]]] * (flip ? flip[perm[b]] : 1.));}}
4223: else {for (b = 0; b < fdof; b++) {fuse(&a[b], values[clperm[foffset+ b ]] * (flip ? flip[ b ] : 1.));}}
4224: } else {
4225: if (perm) {for (b = 0; b < fdof; b++) {fuse(&a[b], values[foffset+perm[b]] * (flip ? flip[perm[b]] : 1.));}}
4226: else {for (b = 0; b < fdof; b++) {fuse(&a[b], values[foffset+ b ] * (flip ? flip[ b ] : 1.));}}
4227: }
4228: } else {
4229: PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs);
4230: if (clperm) {
4231: if (perm) {
4232: for (b = 0; b < fdof; b++) {
4233: if ((cind < fcdof) && (b == fcdofs[cind])) {++cind; continue;}
4234: fuse(&a[b], values[clperm[foffset+perm[b]]] * (flip ? flip[perm[b]] : 1.));
4235: }
4236: } else {
4237: for (b = 0; b < fdof; b++) {
4238: if ((cind < fcdof) && (b == fcdofs[cind])) {++cind; continue;}
4239: fuse(&a[b], values[clperm[foffset+ b ]] * (flip ? flip[ b ] : 1.));
4240: }
4241: }
4242: } else {
4243: if (perm) {
4244: for (b = 0; b < fdof; b++) {
4245: if ((cind < fcdof) && (b == fcdofs[cind])) {++cind; continue;}
4246: fuse(&a[b], values[foffset+perm[b]] * (flip ? flip[perm[b]] : 1.));
4247: }
4248: } else {
4249: for (b = 0; b < fdof; b++) {
4250: if ((cind < fcdof) && (b == fcdofs[cind])) {++cind; continue;}
4251: fuse(&a[b], values[foffset+ b ] * (flip ? flip[ b ] : 1.));
4252: }
4253: }
4254: }
4255: }
4256: *offset += fdof;
4257: return(0);
4258: }
4260: PETSC_STATIC_INLINE PetscErrorCode updatePointFieldsBC_private(PetscSection section, PetscInt point, const PetscInt perm[], const PetscScalar flip[], PetscInt f, void (*fuse)(PetscScalar*, PetscScalar), const PetscInt clperm[], const PetscScalar values[], PetscInt *offset, PetscScalar array[])
4261: {
4262: PetscScalar *a;
4263: PetscInt fdof, foff, fcdof, foffset = *offset;
4264: const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
4265: PetscInt cind = 0, b;
4266: PetscErrorCode ierr;
4269: PetscSectionGetFieldDof(section, point, f, &fdof);
4270: PetscSectionGetFieldConstraintDof(section, point, f, &fcdof);
4271: PetscSectionGetFieldOffset(section, point, f, &foff);
4272: a = &array[foff];
4273: if (fcdof) {
4274: PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs);
4275: if (clperm) {
4276: if (perm) {
4277: for (b = 0; b < fdof; b++) {
4278: if ((cind < fcdof) && (b == fcdofs[cind])) {
4279: fuse(&a[b], values[clperm[foffset+perm[b]]] * (flip ? flip[perm[b]] : 1.));
4280: ++cind;
4281: }
4282: }
4283: } else {
4284: for (b = 0; b < fdof; b++) {
4285: if ((cind < fcdof) && (b == fcdofs[cind])) {
4286: fuse(&a[b], values[clperm[foffset+ b ]] * (flip ? flip[ b ] : 1.));
4287: ++cind;
4288: }
4289: }
4290: }
4291: } else {
4292: if (perm) {
4293: for (b = 0; b < fdof; b++) {
4294: if ((cind < fcdof) && (b == fcdofs[cind])) {
4295: fuse(&a[b], values[foffset+perm[b]] * (flip ? flip[perm[b]] : 1.));
4296: ++cind;
4297: }
4298: }
4299: } else {
4300: for (b = 0; b < fdof; b++) {
4301: if ((cind < fcdof) && (b == fcdofs[cind])) {
4302: fuse(&a[b], values[foffset+ b ] * (flip ? flip[ b ] : 1.));
4303: ++cind;
4304: }
4305: }
4306: }
4307: }
4308: }
4309: *offset += fdof;
4310: return(0);
4311: }
4313: PETSC_STATIC_INLINE PetscErrorCode DMPlexVecSetClosure_Depth1_Static(DM dm, PetscSection section, Vec v, PetscInt point, const PetscScalar values[], InsertMode mode)
4314: {
4315: PetscScalar *array;
4316: const PetscInt *cone, *coneO;
4317: PetscInt pStart, pEnd, p, numPoints, off, dof;
4318: PetscErrorCode ierr;
4321: PetscSectionGetChart(section, &pStart, &pEnd);
4322: DMPlexGetConeSize(dm, point, &numPoints);
4323: DMPlexGetCone(dm, point, &cone);
4324: DMPlexGetConeOrientation(dm, point, &coneO);
4325: VecGetArray(v, &array);
4326: for (p = 0, off = 0; p <= numPoints; ++p, off += dof) {
4327: const PetscInt cp = !p ? point : cone[p-1];
4328: const PetscInt o = !p ? 0 : coneO[p-1];
4330: if ((cp < pStart) || (cp >= pEnd)) {dof = 0; continue;}
4331: PetscSectionGetDof(section, cp, &dof);
4332: /* ADD_VALUES */
4333: {
4334: const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
4335: PetscScalar *a;
4336: PetscInt cdof, coff, cind = 0, k;
4338: PetscSectionGetConstraintDof(section, cp, &cdof);
4339: PetscSectionGetOffset(section, cp, &coff);
4340: a = &array[coff];
4341: if (!cdof) {
4342: if (o >= 0) {
4343: for (k = 0; k < dof; ++k) {
4344: a[k] += values[off+k];
4345: }
4346: } else {
4347: for (k = 0; k < dof; ++k) {
4348: a[k] += values[off+dof-k-1];
4349: }
4350: }
4351: } else {
4352: PetscSectionGetConstraintIndices(section, cp, &cdofs);
4353: if (o >= 0) {
4354: for (k = 0; k < dof; ++k) {
4355: if ((cind < cdof) && (k == cdofs[cind])) {++cind; continue;}
4356: a[k] += values[off+k];
4357: }
4358: } else {
4359: for (k = 0; k < dof; ++k) {
4360: if ((cind < cdof) && (k == cdofs[cind])) {++cind; continue;}
4361: a[k] += values[off+dof-k-1];
4362: }
4363: }
4364: }
4365: }
4366: }
4367: VecRestoreArray(v, &array);
4368: return(0);
4369: }
4371: /*@C
4372: DMPlexVecSetClosure - Set an array of the values on the closure of 'point'
4374: Not collective
4376: Input Parameters:
4377: + dm - The DM
4378: . section - The section describing the layout in v, or NULL to use the default section
4379: . v - The local vector
4380: . point - The point in the DM
4381: . values - The array of values
4382: - mode - The insert mode, where INSERT_ALL_VALUES and ADD_ALL_VALUES also overwrite boundary conditions
4384: Fortran Notes:
4385: This routine is only available in Fortran 90, and you must include petsc.h90 in your code.
4387: Level: intermediate
4389: .seealso DMPlexVecGetClosure(), DMPlexMatSetClosure()
4390: @*/
4391: PetscErrorCode DMPlexVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt point, const PetscScalar values[], InsertMode mode)
4392: {
4393: PetscSection clSection;
4394: IS clPoints;
4395: PetscScalar *array;
4396: PetscInt *points = NULL;
4397: const PetscInt *clp, *clperm;
4398: PetscInt depth, numFields, numPoints, p;
4399: PetscErrorCode ierr;
4403: if (!section) {DMGetDefaultSection(dm, §ion);}
4406: DMPlexGetDepth(dm, &depth);
4407: PetscSectionGetNumFields(section, &numFields);
4408: if (depth == 1 && numFields < 2 && mode == ADD_VALUES) {
4409: DMPlexVecSetClosure_Depth1_Static(dm, section, v, point, values, mode);
4410: return(0);
4411: }
4412: /* Get points */
4413: PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject) dm, NULL, &clperm);
4414: DMPlexGetCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
4415: /* Get array */
4416: VecGetArray(v, &array);
4417: /* Get values */
4418: if (numFields > 0) {
4419: PetscInt offset = 0, f;
4420: for (f = 0; f < numFields; ++f) {
4421: const PetscInt **perms = NULL;
4422: const PetscScalar **flips = NULL;
4424: PetscSectionGetFieldPointSyms(section,f,numPoints,points,&perms,&flips);
4425: switch (mode) {
4426: case INSERT_VALUES:
4427: for (p = 0; p < numPoints; p++) {
4428: const PetscInt point = points[2*p];
4429: const PetscInt *perm = perms ? perms[p] : NULL;
4430: const PetscScalar *flip = flips ? flips[p] : NULL;
4431: updatePointFields_private(section, point, perm, flip, f, insert, PETSC_FALSE, clperm, values, &offset, array);
4432: } break;
4433: case INSERT_ALL_VALUES:
4434: for (p = 0; p < numPoints; p++) {
4435: const PetscInt point = points[2*p];
4436: const PetscInt *perm = perms ? perms[p] : NULL;
4437: const PetscScalar *flip = flips ? flips[p] : NULL;
4438: updatePointFields_private(section, point, perm, flip, f, insert, PETSC_TRUE, clperm, values, &offset, array);
4439: } break;
4440: case INSERT_BC_VALUES:
4441: for (p = 0; p < numPoints; p++) {
4442: const PetscInt point = points[2*p];
4443: const PetscInt *perm = perms ? perms[p] : NULL;
4444: const PetscScalar *flip = flips ? flips[p] : NULL;
4445: updatePointFieldsBC_private(section, point, perm, flip, f, insert, clperm, values, &offset, array);
4446: } break;
4447: case ADD_VALUES:
4448: for (p = 0; p < numPoints; p++) {
4449: const PetscInt point = points[2*p];
4450: const PetscInt *perm = perms ? perms[p] : NULL;
4451: const PetscScalar *flip = flips ? flips[p] : NULL;
4452: updatePointFields_private(section, point, perm, flip, f, add, PETSC_FALSE, clperm, values, &offset, array);
4453: } break;
4454: case ADD_ALL_VALUES:
4455: for (p = 0; p < numPoints; p++) {
4456: const PetscInt point = points[2*p];
4457: const PetscInt *perm = perms ? perms[p] : NULL;
4458: const PetscScalar *flip = flips ? flips[p] : NULL;
4459: updatePointFields_private(section, point, perm, flip, f, add, PETSC_TRUE, clperm, values, &offset, array);
4460: } break;
4461: case ADD_BC_VALUES:
4462: for (p = 0; p < numPoints; p++) {
4463: const PetscInt point = points[2*p];
4464: const PetscInt *perm = perms ? perms[p] : NULL;
4465: const PetscScalar *flip = flips ? flips[p] : NULL;
4466: updatePointFieldsBC_private(section, point, perm, flip, f, add, clperm, values, &offset, array);
4467: } break;
4468: default:
4469: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insert mode %d", mode);
4470: }
4471: PetscSectionRestoreFieldPointSyms(section,f,numPoints,points,&perms,&flips);
4472: }
4473: } else {
4474: PetscInt dof, off;
4475: const PetscInt **perms = NULL;
4476: const PetscScalar **flips = NULL;
4478: PetscSectionGetPointSyms(section,numPoints,points,&perms,&flips);
4479: switch (mode) {
4480: case INSERT_VALUES:
4481: for (p = 0, off = 0; p < numPoints; p++, off += dof) {
4482: const PetscInt point = points[2*p];
4483: const PetscInt *perm = perms ? perms[p] : NULL;
4484: const PetscScalar *flip = flips ? flips[p] : NULL;
4485: PetscSectionGetDof(section, point, &dof);
4486: updatePoint_private(section, point, dof, insert, PETSC_FALSE, perm, flip, clperm, values, off, array);
4487: } break;
4488: case INSERT_ALL_VALUES:
4489: for (p = 0, off = 0; p < numPoints; p++, off += dof) {
4490: const PetscInt point = points[2*p];
4491: const PetscInt *perm = perms ? perms[p] : NULL;
4492: const PetscScalar *flip = flips ? flips[p] : NULL;
4493: PetscSectionGetDof(section, point, &dof);
4494: updatePoint_private(section, point, dof, insert, PETSC_TRUE, perm, flip, clperm, values, off, array);
4495: } break;
4496: case INSERT_BC_VALUES:
4497: for (p = 0, off = 0; p < numPoints; p++, off += dof) {
4498: const PetscInt point = points[2*p];
4499: const PetscInt *perm = perms ? perms[p] : NULL;
4500: const PetscScalar *flip = flips ? flips[p] : NULL;
4501: PetscSectionGetDof(section, point, &dof);
4502: updatePointBC_private(section, point, dof, insert, perm, flip, clperm, values, off, array);
4503: } break;
4504: case ADD_VALUES:
4505: for (p = 0, off = 0; p < numPoints; p++, off += dof) {
4506: const PetscInt point = points[2*p];
4507: const PetscInt *perm = perms ? perms[p] : NULL;
4508: const PetscScalar *flip = flips ? flips[p] : NULL;
4509: PetscSectionGetDof(section, point, &dof);
4510: updatePoint_private(section, point, dof, add, PETSC_FALSE, perm, flip, clperm, values, off, array);
4511: } break;
4512: case ADD_ALL_VALUES:
4513: for (p = 0, off = 0; p < numPoints; p++, off += dof) {
4514: const PetscInt point = points[2*p];
4515: const PetscInt *perm = perms ? perms[p] : NULL;
4516: const PetscScalar *flip = flips ? flips[p] : NULL;
4517: PetscSectionGetDof(section, point, &dof);
4518: updatePoint_private(section, point, dof, add, PETSC_TRUE, perm, flip, clperm, values, off, array);
4519: } break;
4520: case ADD_BC_VALUES:
4521: for (p = 0, off = 0; p < numPoints; p++, off += dof) {
4522: const PetscInt point = points[2*p];
4523: const PetscInt *perm = perms ? perms[p] : NULL;
4524: const PetscScalar *flip = flips ? flips[p] : NULL;
4525: PetscSectionGetDof(section, point, &dof);
4526: updatePointBC_private(section, point, dof, add, perm, flip, clperm, values, off, array);
4527: } break;
4528: default:
4529: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insert mode %d", mode);
4530: }
4531: PetscSectionRestorePointSyms(section,numPoints,points,&perms,&flips);
4532: }
4533: /* Cleanup points */
4534: DMPlexRestoreCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
4535: /* Cleanup array */
4536: VecRestoreArray(v, &array);
4537: return(0);
4538: }
4540: PetscErrorCode DMPlexVecSetFieldClosure_Internal(DM dm, PetscSection section, Vec v, PetscBool fieldActive[], PetscInt point, const PetscScalar values[], InsertMode mode)
4541: {
4542: PetscSection clSection;
4543: IS clPoints;
4544: PetscScalar *array;
4545: PetscInt *points = NULL;
4546: const PetscInt *clp, *clperm;
4547: PetscInt numFields, numPoints, p;
4548: PetscInt offset = 0, f;
4549: PetscErrorCode ierr;
4553: if (!section) {DMGetDefaultSection(dm, §ion);}
4556: PetscSectionGetNumFields(section, &numFields);
4557: /* Get points */
4558: PetscSectionGetClosureInversePermutation_Internal(section, (PetscObject) dm, NULL, &clperm);
4559: DMPlexGetCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
4560: /* Get array */
4561: VecGetArray(v, &array);
4562: /* Get values */
4563: for (f = 0; f < numFields; ++f) {
4564: const PetscInt **perms = NULL;
4565: const PetscScalar **flips = NULL;
4567: if (!fieldActive[f]) {
4568: for (p = 0; p < numPoints*2; p += 2) {
4569: PetscInt fdof;
4570: PetscSectionGetFieldDof(section, points[p], f, &fdof);
4571: offset += fdof;
4572: }
4573: continue;
4574: }
4575: PetscSectionGetFieldPointSyms(section,f,numPoints,points,&perms,&flips);
4576: switch (mode) {
4577: case INSERT_VALUES:
4578: for (p = 0; p < numPoints; p++) {
4579: const PetscInt point = points[2*p];
4580: const PetscInt *perm = perms ? perms[p] : NULL;
4581: const PetscScalar *flip = flips ? flips[p] : NULL;
4582: updatePointFields_private(section, point, perm, flip, f, insert, PETSC_FALSE, clperm, values, &offset, array);
4583: } break;
4584: case INSERT_ALL_VALUES:
4585: for (p = 0; p < numPoints; p++) {
4586: const PetscInt point = points[2*p];
4587: const PetscInt *perm = perms ? perms[p] : NULL;
4588: const PetscScalar *flip = flips ? flips[p] : NULL;
4589: updatePointFields_private(section, point, perm, flip, f, insert, PETSC_TRUE, clperm, values, &offset, array);
4590: } break;
4591: case INSERT_BC_VALUES:
4592: for (p = 0; p < numPoints; p++) {
4593: const PetscInt point = points[2*p];
4594: const PetscInt *perm = perms ? perms[p] : NULL;
4595: const PetscScalar *flip = flips ? flips[p] : NULL;
4596: updatePointFieldsBC_private(section, point, perm, flip, f, insert, clperm, values, &offset, array);
4597: } break;
4598: case ADD_VALUES:
4599: for (p = 0; p < numPoints; p++) {
4600: const PetscInt point = points[2*p];
4601: const PetscInt *perm = perms ? perms[p] : NULL;
4602: const PetscScalar *flip = flips ? flips[p] : NULL;
4603: updatePointFields_private(section, point, perm, flip, f, add, PETSC_FALSE, clperm, values, &offset, array);
4604: } break;
4605: case ADD_ALL_VALUES:
4606: for (p = 0; p < numPoints; p++) {
4607: const PetscInt point = points[2*p];
4608: const PetscInt *perm = perms ? perms[p] : NULL;
4609: const PetscScalar *flip = flips ? flips[p] : NULL;
4610: updatePointFields_private(section, point, perm, flip, f, add, PETSC_TRUE, clperm, values, &offset, array);
4611: } break;
4612: default:
4613: SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insert mode %d", mode);
4614: }
4615: PetscSectionRestoreFieldPointSyms(section,f,numPoints,points,&perms,&flips);
4616: }
4617: /* Cleanup points */
4618: DMPlexRestoreCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
4619: /* Cleanup array */
4620: VecRestoreArray(v, &array);
4621: return(0);
4622: }
4624: static PetscErrorCode DMPlexPrintMatSetValues(PetscViewer viewer, Mat A, PetscInt point, PetscInt numRIndices, const PetscInt rindices[], PetscInt numCIndices, const PetscInt cindices[], const PetscScalar values[])
4625: {
4626: PetscMPIInt rank;
4627: PetscInt i, j;
4631: MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank);
4632: PetscViewerASCIIPrintf(viewer, "[%d]mat for point %D\n", rank, point);
4633: for (i = 0; i < numRIndices; i++) {PetscViewerASCIIPrintf(viewer, "[%d]mat row indices[%D] = %D\n", rank, i, rindices[i]);}
4634: for (i = 0; i < numCIndices; i++) {PetscViewerASCIIPrintf(viewer, "[%d]mat col indices[%D] = %D\n", rank, i, cindices[i]);}
4635: numCIndices = numCIndices ? numCIndices : numRIndices;
4636: for (i = 0; i < numRIndices; i++) {
4637: PetscViewerASCIIPrintf(viewer, "[%d]", rank);
4638: for (j = 0; j < numCIndices; j++) {
4639: #if defined(PETSC_USE_COMPLEX)
4640: PetscViewerASCIIPrintf(viewer, " (%g,%g)", (double)PetscRealPart(values[i*numCIndices+j]), (double)PetscImaginaryPart(values[i*numCIndices+j]));
4641: #else
4642: PetscViewerASCIIPrintf(viewer, " %g", (double)values[i*numCIndices+j]);
4643: #endif
4644: }
4645: PetscViewerASCIIPrintf(viewer, "\n");
4646: }
4647: return(0);
4648: }
4650: /* . off - The global offset of this point */
4651: PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection section, PetscInt point, PetscInt off, PetscInt *loff, PetscBool setBC, const PetscInt perm[], PetscInt indices[])
4652: {
4653: PetscInt dof; /* The number of unknowns on this point */
4654: PetscInt cdof; /* The number of constraints on this point */
4655: const PetscInt *cdofs; /* The indices of the constrained dofs on this point */
4656: PetscInt cind = 0, k;
4657: PetscErrorCode ierr;
4660: PetscSectionGetDof(section, point, &dof);
4661: PetscSectionGetConstraintDof(section, point, &cdof);
4662: if (!cdof || setBC) {
4663: if (perm) {
4664: for (k = 0; k < dof; k++) indices[*loff+perm[k]] = off + k;
4665: } else {
4666: for (k = 0; k < dof; k++) indices[*loff+k] = off + k;
4667: }
4668: } else {
4669: PetscSectionGetConstraintIndices(section, point, &cdofs);
4670: if (perm) {
4671: for (k = 0; k < dof; ++k) {
4672: if ((cind < cdof) && (k == cdofs[cind])) {
4673: /* Insert check for returning constrained indices */
4674: indices[*loff+perm[k]] = -(off+k+1);
4675: ++cind;
4676: } else {
4677: indices[*loff+perm[k]] = off+k-cind;
4678: }
4679: }
4680: } else {
4681: for (k = 0; k < dof; ++k) {
4682: if ((cind < cdof) && (k == cdofs[cind])) {
4683: /* Insert check for returning constrained indices */
4684: indices[*loff+k] = -(off+k+1);
4685: ++cind;
4686: } else {
4687: indices[*loff+k] = off+k-cind;
4688: }
4689: }
4690: }
4691: }
4692: *loff += dof;
4693: return(0);
4694: }
4696: /* . off - The global offset of this point */
4697: PetscErrorCode DMPlexGetIndicesPointFields_Internal(PetscSection section, PetscInt point, PetscInt off, PetscInt foffs[], PetscBool setBC, const PetscInt ***perms, PetscInt permsoff, PetscInt indices[])
4698: {
4699: PetscInt numFields, foff, f;
4703: PetscSectionGetNumFields(section, &numFields);
4704: for (f = 0, foff = 0; f < numFields; ++f) {
4705: PetscInt fdof, cfdof;
4706: const PetscInt *fcdofs; /* The indices of the constrained dofs for field f on this point */
4707: PetscInt cind = 0, b;
4708: const PetscInt *perm = (perms && perms[f]) ? perms[f][permsoff] : NULL;
4710: PetscSectionGetFieldDof(section, point, f, &fdof);
4711: PetscSectionGetFieldConstraintDof(section, point, f, &cfdof);
4712: if (!cfdof || setBC) {
4713: if (perm) {for (b = 0; b < fdof; b++) {indices[foffs[f]+perm[b]] = off+foff+b;}}
4714: else {for (b = 0; b < fdof; b++) {indices[foffs[f]+ b ] = off+foff+b;}}
4715: } else {
4716: PetscSectionGetFieldConstraintIndices(section, point, f, &fcdofs);
4717: if (perm) {
4718: for (b = 0; b < fdof; b++) {
4719: if ((cind < cfdof) && (b == fcdofs[cind])) {
4720: indices[foffs[f]+perm[b]] = -(off+foff+b+1);
4721: ++cind;
4722: } else {
4723: indices[foffs[f]+perm[b]] = off+foff+b-cind;
4724: }
4725: }
4726: } else {
4727: for (b = 0; b < fdof; b++) {
4728: if ((cind < cfdof) && (b == fcdofs[cind])) {
4729: indices[foffs[f]+b] = -(off+foff+b+1);
4730: ++cind;
4731: } else {
4732: indices[foffs[f]+b] = off+foff+b-cind;
4733: }
4734: }
4735: }
4736: }
4737: foff += (setBC ? fdof : (fdof - cfdof));
4738: foffs[f] += fdof;
4739: }
4740: return(0);
4741: }
4743: 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)
4744: {
4745: Mat cMat;
4746: PetscSection aSec, cSec;
4747: IS aIS;
4748: PetscInt aStart = -1, aEnd = -1;
4749: const PetscInt *anchors;
4750: PetscInt numFields, f, p, q, newP = 0;
4751: PetscInt newNumPoints = 0, newNumIndices = 0;
4752: PetscInt *newPoints, *indices, *newIndices;
4753: PetscInt maxAnchor, maxDof;
4754: PetscInt newOffsets[32];
4755: PetscInt *pointMatOffsets[32];
4756: PetscInt *newPointOffsets[32];
4757: PetscScalar *pointMat[32];
4758: PetscScalar *newValues=NULL,*tmpValues;
4759: PetscBool anyConstrained = PETSC_FALSE;
4760: PetscErrorCode ierr;
4765: PetscSectionGetNumFields(section, &numFields);
4767: DMPlexGetAnchors(dm,&aSec,&aIS);
4768: /* if there are point-to-point constraints */
4769: if (aSec) {
4770: PetscMemzero(newOffsets, 32 * sizeof(PetscInt));
4771: ISGetIndices(aIS,&anchors);
4772: PetscSectionGetChart(aSec,&aStart,&aEnd);
4773: /* figure out how many points are going to be in the new element matrix
4774: * (we allow double counting, because it's all just going to be summed
4775: * into the global matrix anyway) */
4776: for (p = 0; p < 2*numPoints; p+=2) {
4777: PetscInt b = points[p];
4778: PetscInt bDof = 0, bSecDof;
4780: PetscSectionGetDof(section,b,&bSecDof);
4781: if (!bSecDof) {
4782: continue;
4783: }
4784: if (b >= aStart && b < aEnd) {
4785: PetscSectionGetDof(aSec,b,&bDof);
4786: }
4787: if (bDof) {
4788: /* this point is constrained */
4789: /* it is going to be replaced by its anchors */
4790: PetscInt bOff, q;
4792: anyConstrained = PETSC_TRUE;
4793: newNumPoints += bDof;
4794: PetscSectionGetOffset(aSec,b,&bOff);
4795: for (q = 0; q < bDof; q++) {
4796: PetscInt a = anchors[bOff + q];
4797: PetscInt aDof;
4799: PetscSectionGetDof(section,a,&aDof);
4800: newNumIndices += aDof;
4801: for (f = 0; f < numFields; ++f) {
4802: PetscInt fDof;
4804: PetscSectionGetFieldDof(section, a, f, &fDof);
4805: newOffsets[f+1] += fDof;
4806: }
4807: }
4808: }
4809: else {
4810: /* this point is not constrained */
4811: newNumPoints++;
4812: newNumIndices += bSecDof;
4813: for (f = 0; f < numFields; ++f) {
4814: PetscInt fDof;
4816: PetscSectionGetFieldDof(section, b, f, &fDof);
4817: newOffsets[f+1] += fDof;
4818: }
4819: }
4820: }
4821: }
4822: if (!anyConstrained) {
4823: if (outNumPoints) *outNumPoints = 0;
4824: if (outNumIndices) *outNumIndices = 0;
4825: if (outPoints) *outPoints = NULL;
4826: if (outValues) *outValues = NULL;
4827: if (aSec) {ISRestoreIndices(aIS,&anchors);}
4828: return(0);
4829: }
4831: if (outNumPoints) *outNumPoints = newNumPoints;
4832: if (outNumIndices) *outNumIndices = newNumIndices;
4834: for (f = 0; f < numFields; ++f) newOffsets[f+1] += newOffsets[f];
4836: if (!outPoints && !outValues) {
4837: if (offsets) {
4838: for (f = 0; f <= numFields; f++) {
4839: offsets[f] = newOffsets[f];
4840: }
4841: }
4842: if (aSec) {ISRestoreIndices(aIS,&anchors);}
4843: return(0);
4844: }
4846: if (numFields && newOffsets[numFields] != newNumIndices) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid size for closure %D should be %D", newOffsets[numFields], newNumIndices);
4848: DMGetDefaultConstraints(dm, &cSec, &cMat);
4850: /* workspaces */
4851: if (numFields) {
4852: for (f = 0; f < numFields; f++) {
4853: DMGetWorkArray(dm,numPoints+1,PETSC_INT,&pointMatOffsets[f]);
4854: DMGetWorkArray(dm,numPoints+1,PETSC_INT,&newPointOffsets[f]);
4855: }
4856: }
4857: else {
4858: DMGetWorkArray(dm,numPoints+1,PETSC_INT,&pointMatOffsets[0]);
4859: DMGetWorkArray(dm,numPoints,PETSC_INT,&newPointOffsets[0]);
4860: }
4862: /* get workspaces for the point-to-point matrices */
4863: if (numFields) {
4864: PetscInt totalOffset, totalMatOffset;
4866: for (p = 0; p < numPoints; p++) {
4867: PetscInt b = points[2*p];
4868: PetscInt bDof = 0, bSecDof;
4870: PetscSectionGetDof(section,b,&bSecDof);
4871: if (!bSecDof) {
4872: for (f = 0; f < numFields; f++) {
4873: newPointOffsets[f][p + 1] = 0;
4874: pointMatOffsets[f][p + 1] = 0;
4875: }
4876: continue;
4877: }
4878: if (b >= aStart && b < aEnd) {
4879: PetscSectionGetDof(aSec, b, &bDof);
4880: }
4881: if (bDof) {
4882: for (f = 0; f < numFields; f++) {
4883: PetscInt fDof, q, bOff, allFDof = 0;
4885: PetscSectionGetFieldDof(section, b, f, &fDof);
4886: PetscSectionGetOffset(aSec, b, &bOff);
4887: for (q = 0; q < bDof; q++) {
4888: PetscInt a = anchors[bOff + q];
4889: PetscInt aFDof;
4891: PetscSectionGetFieldDof(section, a, f, &aFDof);
4892: allFDof += aFDof;
4893: }
4894: newPointOffsets[f][p+1] = allFDof;
4895: pointMatOffsets[f][p+1] = fDof * allFDof;
4896: }
4897: }
4898: else {
4899: for (f = 0; f < numFields; f++) {
4900: PetscInt fDof;
4902: PetscSectionGetFieldDof(section, b, f, &fDof);
4903: newPointOffsets[f][p+1] = fDof;
4904: pointMatOffsets[f][p+1] = 0;
4905: }
4906: }
4907: }
4908: for (f = 0, totalOffset = 0, totalMatOffset = 0; f < numFields; f++) {
4909: newPointOffsets[f][0] = totalOffset;
4910: pointMatOffsets[f][0] = totalMatOffset;
4911: for (p = 0; p < numPoints; p++) {
4912: newPointOffsets[f][p+1] += newPointOffsets[f][p];
4913: pointMatOffsets[f][p+1] += pointMatOffsets[f][p];
4914: }
4915: totalOffset = newPointOffsets[f][numPoints];
4916: totalMatOffset = pointMatOffsets[f][numPoints];
4917: DMGetWorkArray(dm,pointMatOffsets[f][numPoints],PETSC_SCALAR,&pointMat[f]);
4918: }
4919: }
4920: else {
4921: for (p = 0; p < numPoints; p++) {
4922: PetscInt b = points[2*p];
4923: PetscInt bDof = 0, bSecDof;
4925: PetscSectionGetDof(section,b,&bSecDof);
4926: if (!bSecDof) {
4927: newPointOffsets[0][p + 1] = 0;
4928: pointMatOffsets[0][p + 1] = 0;
4929: continue;
4930: }
4931: if (b >= aStart && b < aEnd) {
4932: PetscSectionGetDof(aSec, b, &bDof);
4933: }
4934: if (bDof) {
4935: PetscInt bOff, q, allDof = 0;
4937: PetscSectionGetOffset(aSec, b, &bOff);
4938: for (q = 0; q < bDof; q++) {
4939: PetscInt a = anchors[bOff + q], aDof;
4941: PetscSectionGetDof(section, a, &aDof);
4942: allDof += aDof;
4943: }
4944: newPointOffsets[0][p+1] = allDof;
4945: pointMatOffsets[0][p+1] = bSecDof * allDof;
4946: }
4947: else {
4948: newPointOffsets[0][p+1] = bSecDof;
4949: pointMatOffsets[0][p+1] = 0;
4950: }
4951: }
4952: newPointOffsets[0][0] = 0;
4953: pointMatOffsets[0][0] = 0;
4954: for (p = 0; p < numPoints; p++) {
4955: newPointOffsets[0][p+1] += newPointOffsets[0][p];
4956: pointMatOffsets[0][p+1] += pointMatOffsets[0][p];
4957: }
4958: DMGetWorkArray(dm,pointMatOffsets[0][numPoints],PETSC_SCALAR,&pointMat[0]);
4959: }
4961: /* output arrays */
4962: DMGetWorkArray(dm,2*newNumPoints,PETSC_INT,&newPoints);
4964: /* get the point-to-point matrices; construct newPoints */
4965: PetscSectionGetMaxDof(aSec, &maxAnchor);
4966: PetscSectionGetMaxDof(section, &maxDof);
4967: DMGetWorkArray(dm,maxDof,PETSC_INT,&indices);
4968: DMGetWorkArray(dm,maxAnchor*maxDof,PETSC_INT,&newIndices);
4969: if (numFields) {
4970: for (p = 0, newP = 0; p < numPoints; p++) {
4971: PetscInt b = points[2*p];
4972: PetscInt o = points[2*p+1];
4973: PetscInt bDof = 0, bSecDof;
4975: PetscSectionGetDof(section, b, &bSecDof);
4976: if (!bSecDof) {
4977: continue;
4978: }
4979: if (b >= aStart && b < aEnd) {
4980: PetscSectionGetDof(aSec, b, &bDof);
4981: }
4982: if (bDof) {
4983: PetscInt fStart[32], fEnd[32], fAnchorStart[32], fAnchorEnd[32], bOff, q;
4985: fStart[0] = 0;
4986: fEnd[0] = 0;
4987: for (f = 0; f < numFields; f++) {
4988: PetscInt fDof;
4990: PetscSectionGetFieldDof(cSec, b, f, &fDof);
4991: fStart[f+1] = fStart[f] + fDof;
4992: fEnd[f+1] = fStart[f+1];
4993: }
4994: PetscSectionGetOffset(cSec, b, &bOff);
4995: DMPlexGetIndicesPointFields_Internal(cSec, b, bOff, fEnd, PETSC_TRUE, perms, p, indices);
4997: fAnchorStart[0] = 0;
4998: fAnchorEnd[0] = 0;
4999: for (f = 0; f < numFields; f++) {
5000: PetscInt fDof = newPointOffsets[f][p + 1] - newPointOffsets[f][p];
5002: fAnchorStart[f+1] = fAnchorStart[f] + fDof;
5003: fAnchorEnd[f+1] = fAnchorStart[f + 1];
5004: }
5005: PetscSectionGetOffset(aSec, b, &bOff);
5006: for (q = 0; q < bDof; q++) {
5007: PetscInt a = anchors[bOff + q], aOff;
5009: /* we take the orientation of ap into account in the order that we constructed the indices above: the newly added points have no orientation */
5010: newPoints[2*(newP + q)] = a;
5011: newPoints[2*(newP + q) + 1] = 0;
5012: PetscSectionGetOffset(section, a, &aOff);
5013: DMPlexGetIndicesPointFields_Internal(section, a, aOff, fAnchorEnd, PETSC_TRUE, NULL, -1, newIndices);
5014: }
5015: newP += bDof;
5017: if (outValues) {
5018: /* get the point-to-point submatrix */
5019: for (f = 0; f < numFields; f++) {
5020: MatGetValues(cMat,fEnd[f]-fStart[f],indices + fStart[f],fAnchorEnd[f] - fAnchorStart[f],newIndices + fAnchorStart[f],pointMat[f] + pointMatOffsets[f][p]);
5021: }
5022: }
5023: }
5024: else {
5025: newPoints[2 * newP] = b;
5026: newPoints[2 * newP + 1] = o;
5027: newP++;
5028: }
5029: }
5030: } else {
5031: for (p = 0; p < numPoints; p++) {
5032: PetscInt b = points[2*p];
5033: PetscInt o = points[2*p+1];
5034: PetscInt bDof = 0, bSecDof;
5036: PetscSectionGetDof(section, b, &bSecDof);
5037: if (!bSecDof) {
5038: continue;
5039: }
5040: if (b >= aStart && b < aEnd) {
5041: PetscSectionGetDof(aSec, b, &bDof);
5042: }
5043: if (bDof) {
5044: PetscInt bEnd = 0, bAnchorEnd = 0, bOff;
5046: PetscSectionGetOffset(cSec, b, &bOff);
5047: DMPlexGetIndicesPoint_Internal(cSec, b, bOff, &bEnd, PETSC_TRUE, (perms && perms[0]) ? perms[0][p] : NULL, indices);
5049: PetscSectionGetOffset (aSec, b, &bOff);
5050: for (q = 0; q < bDof; q++) {
5051: PetscInt a = anchors[bOff + q], aOff;
5053: /* we take the orientation of ap into account in the order that we constructed the indices above: the newly added points have no orientation */
5055: newPoints[2*(newP + q)] = a;
5056: newPoints[2*(newP + q) + 1] = 0;
5057: PetscSectionGetOffset(section, a, &aOff);
5058: DMPlexGetIndicesPoint_Internal(section, a, aOff, &bAnchorEnd, PETSC_TRUE, NULL, newIndices);
5059: }
5060: newP += bDof;
5062: /* get the point-to-point submatrix */
5063: if (outValues) {
5064: MatGetValues(cMat,bEnd,indices,bAnchorEnd,newIndices,pointMat[0] + pointMatOffsets[0][p]);
5065: }
5066: }
5067: else {
5068: newPoints[2 * newP] = b;
5069: newPoints[2 * newP + 1] = o;
5070: newP++;
5071: }
5072: }
5073: }
5075: if (outValues) {
5076: DMGetWorkArray(dm,newNumIndices*numIndices,PETSC_SCALAR,&tmpValues);
5077: PetscMemzero(tmpValues,newNumIndices*numIndices*sizeof(*tmpValues));
5078: /* multiply constraints on the right */
5079: if (numFields) {
5080: for (f = 0; f < numFields; f++) {
5081: PetscInt oldOff = offsets[f];
5083: for (p = 0; p < numPoints; p++) {
5084: PetscInt cStart = newPointOffsets[f][p];
5085: PetscInt b = points[2 * p];
5086: PetscInt c, r, k;
5087: PetscInt dof;
5089: PetscSectionGetFieldDof(section,b,f,&dof);
5090: if (!dof) {
5091: continue;
5092: }
5093: if (pointMatOffsets[f][p] < pointMatOffsets[f][p + 1]) {
5094: PetscInt nCols = newPointOffsets[f][p+1]-cStart;
5095: const PetscScalar *mat = pointMat[f] + pointMatOffsets[f][p];
5097: for (r = 0; r < numIndices; r++) {
5098: for (c = 0; c < nCols; c++) {
5099: for (k = 0; k < dof; k++) {
5100: tmpValues[r * newNumIndices + cStart + c] += values[r * numIndices + oldOff + k] * mat[k * nCols + c];
5101: }
5102: }
5103: }
5104: }
5105: else {
5106: /* copy this column as is */
5107: for (r = 0; r < numIndices; r++) {
5108: for (c = 0; c < dof; c++) {
5109: tmpValues[r * newNumIndices + cStart + c] = values[r * numIndices + oldOff + c];
5110: }
5111: }
5112: }
5113: oldOff += dof;
5114: }
5115: }
5116: }
5117: else {
5118: PetscInt oldOff = 0;
5119: for (p = 0; p < numPoints; p++) {
5120: PetscInt cStart = newPointOffsets[0][p];
5121: PetscInt b = points[2 * p];
5122: PetscInt c, r, k;
5123: PetscInt dof;
5125: PetscSectionGetDof(section,b,&dof);
5126: if (!dof) {
5127: continue;
5128: }
5129: if (pointMatOffsets[0][p] < pointMatOffsets[0][p + 1]) {
5130: PetscInt nCols = newPointOffsets[0][p+1]-cStart;
5131: const PetscScalar *mat = pointMat[0] + pointMatOffsets[0][p];
5133: for (r = 0; r < numIndices; r++) {
5134: for (c = 0; c < nCols; c++) {
5135: for (k = 0; k < dof; k++) {
5136: tmpValues[r * newNumIndices + cStart + c] += mat[k * nCols + c] * values[r * numIndices + oldOff + k];
5137: }
5138: }
5139: }
5140: }
5141: else {
5142: /* copy this column as is */
5143: for (r = 0; r < numIndices; r++) {
5144: for (c = 0; c < dof; c++) {
5145: tmpValues[r * newNumIndices + cStart + c] = values[r * numIndices + oldOff + c];
5146: }
5147: }
5148: }
5149: oldOff += dof;
5150: }
5151: }
5153: if (multiplyLeft) {
5154: DMGetWorkArray(dm,newNumIndices*newNumIndices,PETSC_SCALAR,&newValues);
5155: PetscMemzero(newValues,newNumIndices*newNumIndices*sizeof(*newValues));
5156: /* multiply constraints transpose on the left */
5157: if (numFields) {
5158: for (f = 0; f < numFields; f++) {
5159: PetscInt oldOff = offsets[f];
5161: for (p = 0; p < numPoints; p++) {
5162: PetscInt rStart = newPointOffsets[f][p];
5163: PetscInt b = points[2 * p];
5164: PetscInt c, r, k;
5165: PetscInt dof;
5167: PetscSectionGetFieldDof(section,b,f,&dof);
5168: if (pointMatOffsets[f][p] < pointMatOffsets[f][p + 1]) {
5169: PetscInt nRows = newPointOffsets[f][p+1]-rStart;
5170: const PetscScalar *PETSC_RESTRICT mat = pointMat[f] + pointMatOffsets[f][p];
5172: for (r = 0; r < nRows; r++) {
5173: for (c = 0; c < newNumIndices; c++) {
5174: for (k = 0; k < dof; k++) {
5175: newValues[(rStart + r) * newNumIndices + c] += mat[k * nRows + r] * tmpValues[(oldOff + k) * newNumIndices + c];
5176: }
5177: }
5178: }
5179: }
5180: else {
5181: /* copy this row as is */
5182: for (r = 0; r < dof; r++) {
5183: for (c = 0; c < newNumIndices; c++) {
5184: newValues[(rStart + r) * newNumIndices + c] = tmpValues[(oldOff + r) * newNumIndices + c];
5185: }
5186: }
5187: }
5188: oldOff += dof;
5189: }
5190: }
5191: }
5192: else {
5193: PetscInt oldOff = 0;
5195: for (p = 0; p < numPoints; p++) {
5196: PetscInt rStart = newPointOffsets[0][p];
5197: PetscInt b = points[2 * p];
5198: PetscInt c, r, k;
5199: PetscInt dof;
5201: PetscSectionGetDof(section,b,&dof);
5202: if (pointMatOffsets[0][p] < pointMatOffsets[0][p + 1]) {
5203: PetscInt nRows = newPointOffsets[0][p+1]-rStart;
5204: const PetscScalar *PETSC_RESTRICT mat = pointMat[0] + pointMatOffsets[0][p];
5206: for (r = 0; r < nRows; r++) {
5207: for (c = 0; c < newNumIndices; c++) {
5208: for (k = 0; k < dof; k++) {
5209: newValues[(rStart + r) * newNumIndices + c] += mat[k * nRows + r] * tmpValues[(oldOff + k) * newNumIndices + c];
5210: }
5211: }
5212: }
5213: }
5214: else {
5215: /* copy this row as is */
5216: for (r = 0; r < dof; r++) {
5217: for (c = 0; c < newNumIndices; c++) {
5218: newValues[(rStart + r) * newNumIndices + c] = tmpValues[(oldOff + r) * newNumIndices + c];
5219: }
5220: }
5221: }
5222: oldOff += dof;
5223: }
5224: }
5226: DMRestoreWorkArray(dm,newNumIndices*numIndices,PETSC_SCALAR,&tmpValues);
5227: }
5228: else {
5229: newValues = tmpValues;
5230: }
5231: }
5233: /* clean up */
5234: DMRestoreWorkArray(dm,maxDof,PETSC_INT,&indices);
5235: DMRestoreWorkArray(dm,maxAnchor*maxDof,PETSC_INT,&newIndices);
5237: if (numFields) {
5238: for (f = 0; f < numFields; f++) {
5239: DMRestoreWorkArray(dm,pointMatOffsets[f][numPoints],PETSC_SCALAR,&pointMat[f]);
5240: DMRestoreWorkArray(dm,numPoints+1,PETSC_INT,&pointMatOffsets[f]);
5241: DMRestoreWorkArray(dm,numPoints+1,PETSC_INT,&newPointOffsets[f]);
5242: }
5243: }
5244: else {
5245: DMRestoreWorkArray(dm,pointMatOffsets[0][numPoints],PETSC_SCALAR,&pointMat[0]);
5246: DMRestoreWorkArray(dm,numPoints+1,PETSC_INT,&pointMatOffsets[0]);
5247: DMRestoreWorkArray(dm,numPoints+1,PETSC_INT,&newPointOffsets[0]);
5248: }
5249: ISRestoreIndices(aIS,&anchors);
5251: /* output */
5252: if (outPoints) {
5253: *outPoints = newPoints;
5254: }
5255: else {
5256: DMRestoreWorkArray(dm,2*newNumPoints,PETSC_INT,&newPoints);
5257: }
5258: if (outValues) {
5259: *outValues = newValues;
5260: }
5261: for (f = 0; f <= numFields; f++) {
5262: offsets[f] = newOffsets[f];
5263: }
5264: return(0);
5265: }
5267: /*@C
5268: DMPlexGetClosureIndices - Get the global indices in a vector v for all points in the closure of the given point
5270: Not collective
5272: Input Parameters:
5273: + dm - The DM
5274: . section - The section describing the layout in v, or NULL to use the default section
5275: . globalSection - The section describing the parallel layout in v, or NULL to use the default section
5276: - point - The mesh point
5278: Output parameters:
5279: + numIndices - The number of indices
5280: . indices - The indices
5281: - outOffsets - Field offset if not NULL
5283: Note: Must call DMPlexRestoreClosureIndices() to free allocated memory
5285: Level: advanced
5287: .seealso DMPlexRestoreClosureIndices(), DMPlexVecGetClosure(), DMPlexMatSetClosure()
5288: @*/
5289: PetscErrorCode DMPlexGetClosureIndices(DM dm, PetscSection section, PetscSection globalSection, PetscInt point, PetscInt *numIndices, PetscInt **indices, PetscInt *outOffsets)
5290: {
5291: PetscSection clSection;
5292: IS clPoints;
5293: const PetscInt *clp;
5294: const PetscInt **perms[32] = {NULL};
5295: PetscInt *points = NULL, *pointsNew;
5296: PetscInt numPoints, numPointsNew;
5297: PetscInt offsets[32];
5298: PetscInt Nf, Nind, NindNew, off, globalOff, f, p;
5299: PetscErrorCode ierr;
5307: PetscSectionGetNumFields(section, &Nf);
5308: if (Nf > 31) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Number of fields %D limited to 31", Nf);
5309: PetscMemzero(offsets, 32 * sizeof(PetscInt));
5310: /* Get points in closure */
5311: DMPlexGetCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
5312: /* Get number of indices and indices per field */
5313: for (p = 0, Nind = 0; p < numPoints*2; p += 2) {
5314: PetscInt dof, fdof;
5316: PetscSectionGetDof(section, points[p], &dof);
5317: for (f = 0; f < Nf; ++f) {
5318: PetscSectionGetFieldDof(section, points[p], f, &fdof);
5319: offsets[f+1] += fdof;
5320: }
5321: Nind += dof;
5322: }
5323: for (f = 1; f < Nf; ++f) offsets[f+1] += offsets[f];
5324: if (Nf && offsets[Nf] != Nind) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Invalid size for closure %D should be %D", offsets[Nf], Nind);
5325: if (!Nf) offsets[1] = Nind;
5326: /* Get dual space symmetries */
5327: for (f = 0; f < PetscMax(1,Nf); f++) {
5328: if (Nf) {PetscSectionGetFieldPointSyms(section,f,numPoints,points,&perms[f],NULL);}
5329: else {PetscSectionGetPointSyms(section,numPoints,points,&perms[f],NULL);}
5330: }
5331: /* Correct for hanging node constraints */
5332: {
5333: DMPlexAnchorsModifyMat(dm, section, numPoints, Nind, points, perms, NULL, &numPointsNew, &NindNew, &pointsNew, NULL, offsets, PETSC_TRUE);
5334: if (numPointsNew) {
5335: for (f = 0; f < PetscMax(1,Nf); f++) {
5336: if (Nf) {PetscSectionRestoreFieldPointSyms(section,f,numPoints,points,&perms[f],NULL);}
5337: else {PetscSectionRestorePointSyms(section,numPoints,points,&perms[f],NULL);}
5338: }
5339: for (f = 0; f < PetscMax(1,Nf); f++) {
5340: if (Nf) {PetscSectionGetFieldPointSyms(section,f,numPointsNew,pointsNew,&perms[f],NULL);}
5341: else {PetscSectionGetPointSyms(section,numPointsNew,pointsNew,&perms[f],NULL);}
5342: }
5343: DMPlexRestoreCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
5344: numPoints = numPointsNew;
5345: Nind = NindNew;
5346: points = pointsNew;
5347: }
5348: }
5349: /* Calculate indices */
5350: DMGetWorkArray(dm, Nind, PETSC_INT, indices);
5351: if (Nf) {
5352: if (outOffsets) {
5353: PetscInt f;
5355: for (f = 0; f <= Nf; f++) {
5356: outOffsets[f] = offsets[f];
5357: }
5358: }
5359: for (p = 0; p < numPoints; p++) {
5360: PetscSectionGetOffset(globalSection, points[2*p], &globalOff);
5361: DMPlexGetIndicesPointFields_Internal(section, points[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, offsets, PETSC_FALSE, perms, p, *indices);
5362: }
5363: } else {
5364: for (p = 0, off = 0; p < numPoints; p++) {
5365: const PetscInt *perm = perms[0] ? perms[0][p] : NULL;
5367: PetscSectionGetOffset(globalSection, points[2*p], &globalOff);
5368: DMPlexGetIndicesPoint_Internal(section, points[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, &off, PETSC_FALSE, perm, *indices);
5369: }
5370: }
5371: /* Cleanup points */
5372: for (f = 0; f < PetscMax(1,Nf); f++) {
5373: if (Nf) {PetscSectionRestoreFieldPointSyms(section,f,numPoints,points,&perms[f],NULL);}
5374: else {PetscSectionRestorePointSyms(section,numPoints,points,&perms[f],NULL);}
5375: }
5376: if (numPointsNew) {
5377: DMRestoreWorkArray(dm, 2*numPointsNew, PETSC_INT, &pointsNew);
5378: } else {
5379: DMPlexRestoreCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
5380: }
5381: if (numIndices) *numIndices = Nind;
5382: return(0);
5383: }
5385: /*@C
5386: DMPlexRestoreClosureIndices - Restore the indices in a vector v for all points in the closure of the given point
5388: Not collective
5390: Input Parameters:
5391: + dm - The DM
5392: . section - The section describing the layout in v, or NULL to use the default section
5393: . globalSection - The section describing the parallel layout in v, or NULL to use the default section
5394: . point - The mesh point
5395: . numIndices - The number of indices
5396: . indices - The indices
5397: - outOffsets - Field offset if not NULL
5399: Level: advanced
5401: .seealso DMPlexGetClosureIndices(), DMPlexVecGetClosure(), DMPlexMatSetClosure()
5402: @*/
5403: PetscErrorCode DMPlexRestoreClosureIndices(DM dm, PetscSection section, PetscSection globalSection, PetscInt point, PetscInt *numIndices, PetscInt **indices,PetscInt *outOffsets)
5404: {
5410: DMRestoreWorkArray(dm, 0, PETSC_INT, indices);
5411: return(0);
5412: }
5414: /*@C
5415: DMPlexMatSetClosure - Set an array of the values on the closure of 'point'
5417: Not collective
5419: Input Parameters:
5420: + dm - The DM
5421: . section - The section describing the layout in v, or NULL to use the default section
5422: . globalSection - The section describing the layout in v, or NULL to use the default global section
5423: . A - The matrix
5424: . point - The point in the DM
5425: . values - The array of values
5426: - mode - The insert mode, where INSERT_ALL_VALUES and ADD_ALL_VALUES also overwrite boundary conditions
5428: Fortran Notes:
5429: This routine is only available in Fortran 90, and you must include petsc.h90 in your code.
5431: Level: intermediate
5433: .seealso DMPlexVecGetClosure(), DMPlexVecSetClosure()
5434: @*/
5435: PetscErrorCode DMPlexMatSetClosure(DM dm, PetscSection section, PetscSection globalSection, Mat A, PetscInt point, const PetscScalar values[], InsertMode mode)
5436: {
5437: DM_Plex *mesh = (DM_Plex*) dm->data;
5438: PetscSection clSection;
5439: IS clPoints;
5440: PetscInt *points = NULL, *newPoints;
5441: const PetscInt *clp;
5442: PetscInt *indices;
5443: PetscInt offsets[32];
5444: const PetscInt **perms[32] = {NULL};
5445: const PetscScalar **flips[32] = {NULL};
5446: PetscInt numFields, numPoints, newNumPoints, numIndices, newNumIndices, dof, off, globalOff, p, f;
5447: PetscScalar *valCopy = NULL;
5448: PetscScalar *newValues;
5449: PetscErrorCode ierr;
5453: if (!section) {DMGetDefaultSection(dm, §ion);}
5455: if (!globalSection) {DMGetDefaultGlobalSection(dm, &globalSection);}
5458: PetscSectionGetNumFields(section, &numFields);
5459: if (numFields > 31) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Number of fields %D limited to 31", numFields);
5460: PetscMemzero(offsets, 32 * sizeof(PetscInt));
5461: DMPlexGetCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
5462: for (p = 0, numIndices = 0; p < numPoints*2; p += 2) {
5463: PetscInt fdof;
5465: PetscSectionGetDof(section, points[p], &dof);
5466: for (f = 0; f < numFields; ++f) {
5467: PetscSectionGetFieldDof(section, points[p], f, &fdof);
5468: offsets[f+1] += fdof;
5469: }
5470: numIndices += dof;
5471: }
5472: for (f = 1; f < numFields; ++f) offsets[f+1] += offsets[f];
5474: if (numFields && offsets[numFields] != numIndices) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid size for closure %D should be %D", offsets[numFields], numIndices);
5475: /* Get symmetries */
5476: for (f = 0; f < PetscMax(1,numFields); f++) {
5477: if (numFields) {PetscSectionGetFieldPointSyms(section,f,numPoints,points,&perms[f],&flips[f]);}
5478: else {PetscSectionGetPointSyms(section,numPoints,points,&perms[f],&flips[f]);}
5479: if (values && flips[f]) { /* may need to apply sign changes to the element matrix */
5480: PetscInt foffset = offsets[f];
5482: for (p = 0; p < numPoints; p++) {
5483: PetscInt point = points[2*p], fdof;
5484: const PetscScalar *flip = flips[f] ? flips[f][p] : NULL;
5486: if (!numFields) {
5487: PetscSectionGetDof(section,point,&fdof);
5488: } else {
5489: PetscSectionGetFieldDof(section,point,f,&fdof);
5490: }
5491: if (flip) {
5492: PetscInt i, j, k;
5494: if (!valCopy) {
5495: DMGetWorkArray(dm,numIndices*numIndices,PETSC_SCALAR,&valCopy);
5496: for (j = 0; j < numIndices * numIndices; j++) valCopy[j] = values[j];
5497: values = valCopy;
5498: }
5499: for (i = 0; i < fdof; i++) {
5500: PetscScalar fval = flip[i];
5502: for (k = 0; k < numIndices; k++) {
5503: valCopy[numIndices * (foffset + i) + k] *= fval;
5504: valCopy[numIndices * k + (foffset + i)] *= fval;
5505: }
5506: }
5507: }
5508: foffset += fdof;
5509: }
5510: }
5511: }
5512: DMPlexAnchorsModifyMat(dm,section,numPoints,numIndices,points,perms,values,&newNumPoints,&newNumIndices,&newPoints,&newValues,offsets,PETSC_TRUE);
5513: if (newNumPoints) {
5514: if (valCopy) {
5515: DMRestoreWorkArray(dm,numIndices*numIndices,PETSC_SCALAR,&valCopy);
5516: }
5517: for (f = 0; f < PetscMax(1,numFields); f++) {
5518: if (numFields) {PetscSectionRestoreFieldPointSyms(section,f,numPoints,points,&perms[f],&flips[f]);}
5519: else {PetscSectionRestorePointSyms(section,numPoints,points,&perms[f],&flips[f]);}
5520: }
5521: for (f = 0; f < PetscMax(1,numFields); f++) {
5522: if (numFields) {PetscSectionGetFieldPointSyms(section,f,newNumPoints,newPoints,&perms[f],&flips[f]);}
5523: else {PetscSectionGetPointSyms(section,newNumPoints,newPoints,&perms[f],&flips[f]);}
5524: }
5525: DMPlexRestoreCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
5526: numPoints = newNumPoints;
5527: numIndices = newNumIndices;
5528: points = newPoints;
5529: values = newValues;
5530: }
5531: DMGetWorkArray(dm, numIndices, PETSC_INT, &indices);
5532: if (numFields) {
5533: for (p = 0; p < numPoints; p++) {
5534: PetscSectionGetOffset(globalSection, points[2*p], &globalOff);
5535: DMPlexGetIndicesPointFields_Internal(section, points[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, offsets, PETSC_FALSE, perms, p, indices);
5536: }
5537: } else {
5538: for (p = 0, off = 0; p < numPoints; p++) {
5539: const PetscInt *perm = perms[0] ? perms[0][p] : NULL;
5540: PetscSectionGetOffset(globalSection, points[2*p], &globalOff);
5541: DMPlexGetIndicesPoint_Internal(section, points[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, &off, PETSC_FALSE, perm, indices);
5542: }
5543: }
5544: if (mesh->printSetValues) {DMPlexPrintMatSetValues(PETSC_VIEWER_STDOUT_SELF, A, point, numIndices, indices, 0, NULL, values);}
5545: MatSetValues(A, numIndices, indices, numIndices, indices, values, mode);
5546: if (mesh->printFEM > 1) {
5547: PetscInt i;
5548: PetscPrintf(PETSC_COMM_SELF, " Indices:");
5549: for (i = 0; i < numIndices; ++i) {PetscPrintf(PETSC_COMM_SELF, " %D", indices[i]);}
5550: PetscPrintf(PETSC_COMM_SELF, "\n");
5551: }
5552: if (ierr) {
5553: PetscMPIInt rank;
5554: PetscErrorCode ierr2;
5556: ierr2 = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank);CHKERRQ(ierr2);
5557: ierr2 = (*PetscErrorPrintf)("[%d]ERROR in DMPlexMatSetClosure\n", rank);CHKERRQ(ierr2);
5558: ierr2 = DMPlexPrintMatSetValues(PETSC_VIEWER_STDERR_SELF, A, point, numIndices, indices, 0, NULL, values);CHKERRQ(ierr2);
5559: ierr2 = DMRestoreWorkArray(dm, numIndices, PETSC_INT, &indices);CHKERRQ(ierr2);
5560:
5561: }
5562: for (f = 0; f < PetscMax(1,numFields); f++) {
5563: if (numFields) {PetscSectionRestoreFieldPointSyms(section,f,numPoints,points,&perms[f],&flips[f]);}
5564: else {PetscSectionRestorePointSyms(section,numPoints,points,&perms[f],&flips[f]);}
5565: }
5566: if (newNumPoints) {
5567: DMRestoreWorkArray(dm,newNumIndices*newNumIndices,PETSC_SCALAR,&newValues);
5568: DMRestoreWorkArray(dm,2*newNumPoints,PETSC_INT,&newPoints);
5569: }
5570: else {
5571: if (valCopy) {
5572: DMRestoreWorkArray(dm,numIndices*numIndices,PETSC_SCALAR,&valCopy);
5573: }
5574: DMPlexRestoreCompressedClosure(dm,section,point,&numPoints,&points,&clSection,&clPoints,&clp);
5575: }
5576: DMRestoreWorkArray(dm, numIndices, PETSC_INT, &indices);
5577: return(0);
5578: }
5580: PetscErrorCode DMPlexMatSetClosureRefined(DM dmf, PetscSection fsection, PetscSection globalFSection, DM dmc, PetscSection csection, PetscSection globalCSection, Mat A, PetscInt point, const PetscScalar values[], InsertMode mode)
5581: {
5582: DM_Plex *mesh = (DM_Plex*) dmf->data;
5583: PetscInt *fpoints = NULL, *ftotpoints = NULL;
5584: PetscInt *cpoints = NULL;
5585: PetscInt *findices, *cindices;
5586: PetscInt foffsets[32], coffsets[32];
5587: CellRefiner cellRefiner;
5588: PetscInt numFields, numSubcells, maxFPoints, numFPoints, numCPoints, numFIndices, numCIndices, dof, off, globalOff, pStart, pEnd, p, q, r, s, f;
5589: PetscErrorCode ierr;
5594: if (!fsection) {DMGetDefaultSection(dmf, &fsection);}
5596: if (!csection) {DMGetDefaultSection(dmc, &csection);}
5598: if (!globalFSection) {DMGetDefaultGlobalSection(dmf, &globalFSection);}
5600: if (!globalCSection) {DMGetDefaultGlobalSection(dmc, &globalCSection);}
5603: PetscSectionGetNumFields(fsection, &numFields);
5604: if (numFields > 31) SETERRQ1(PetscObjectComm((PetscObject)dmf), PETSC_ERR_ARG_OUTOFRANGE, "Number of fields %D limited to 31", numFields);
5605: PetscMemzero(foffsets, 32 * sizeof(PetscInt));
5606: PetscMemzero(coffsets, 32 * sizeof(PetscInt));
5607: /* Column indices */
5608: DMPlexGetTransitiveClosure(dmc, point, PETSC_TRUE, &numCPoints, &cpoints);
5609: maxFPoints = numCPoints;
5610: /* Compress out points not in the section */
5611: /* TODO: Squeeze out points with 0 dof as well */
5612: PetscSectionGetChart(csection, &pStart, &pEnd);
5613: for (p = 0, q = 0; p < numCPoints*2; p += 2) {
5614: if ((cpoints[p] >= pStart) && (cpoints[p] < pEnd)) {
5615: cpoints[q*2] = cpoints[p];
5616: cpoints[q*2+1] = cpoints[p+1];
5617: ++q;
5618: }
5619: }
5620: numCPoints = q;
5621: for (p = 0, numCIndices = 0; p < numCPoints*2; p += 2) {
5622: PetscInt fdof;
5624: PetscSectionGetDof(csection, cpoints[p], &dof);
5625: if (!dof) continue;
5626: for (f = 0; f < numFields; ++f) {
5627: PetscSectionGetFieldDof(csection, cpoints[p], f, &fdof);
5628: coffsets[f+1] += fdof;
5629: }
5630: numCIndices += dof;
5631: }
5632: for (f = 1; f < numFields; ++f) coffsets[f+1] += coffsets[f];
5633: /* Row indices */
5634: DMPlexGetCellRefiner_Internal(dmc, &cellRefiner);
5635: CellRefinerGetAffineTransforms_Internal(cellRefiner, &numSubcells, NULL, NULL, NULL);
5636: DMGetWorkArray(dmf, maxFPoints*2*numSubcells, PETSC_INT, &ftotpoints);
5637: for (r = 0, q = 0; r < numSubcells; ++r) {
5638: /* TODO Map from coarse to fine cells */
5639: DMPlexGetTransitiveClosure(dmf, point*numSubcells + r, PETSC_TRUE, &numFPoints, &fpoints);
5640: /* Compress out points not in the section */
5641: PetscSectionGetChart(fsection, &pStart, &pEnd);
5642: for (p = 0; p < numFPoints*2; p += 2) {
5643: if ((fpoints[p] >= pStart) && (fpoints[p] < pEnd)) {
5644: PetscSectionGetDof(fsection, fpoints[p], &dof);
5645: if (!dof) continue;
5646: for (s = 0; s < q; ++s) if (fpoints[p] == ftotpoints[s*2]) break;
5647: if (s < q) continue;
5648: ftotpoints[q*2] = fpoints[p];
5649: ftotpoints[q*2+1] = fpoints[p+1];
5650: ++q;
5651: }
5652: }
5653: DMPlexRestoreTransitiveClosure(dmf, point, PETSC_TRUE, &numFPoints, &fpoints);
5654: }
5655: numFPoints = q;
5656: for (p = 0, numFIndices = 0; p < numFPoints*2; p += 2) {
5657: PetscInt fdof;
5659: PetscSectionGetDof(fsection, ftotpoints[p], &dof);
5660: if (!dof) continue;
5661: for (f = 0; f < numFields; ++f) {
5662: PetscSectionGetFieldDof(fsection, ftotpoints[p], f, &fdof);
5663: foffsets[f+1] += fdof;
5664: }
5665: numFIndices += dof;
5666: }
5667: for (f = 1; f < numFields; ++f) foffsets[f+1] += foffsets[f];
5669: if (numFields && foffsets[numFields] != numFIndices) SETERRQ2(PetscObjectComm((PetscObject)dmf), PETSC_ERR_PLIB, "Invalid size for closure %D should be %D", foffsets[numFields], numFIndices);
5670: if (numFields && coffsets[numFields] != numCIndices) SETERRQ2(PetscObjectComm((PetscObject)dmc), PETSC_ERR_PLIB, "Invalid size for closure %D should be %D", coffsets[numFields], numCIndices);
5671: DMGetWorkArray(dmf, numFIndices, PETSC_INT, &findices);
5672: DMGetWorkArray(dmc, numCIndices, PETSC_INT, &cindices);
5673: if (numFields) {
5674: const PetscInt **permsF[32] = {NULL};
5675: const PetscInt **permsC[32] = {NULL};
5677: for (f = 0; f < numFields; f++) {
5678: PetscSectionGetFieldPointSyms(fsection,f,numFPoints,ftotpoints,&permsF[f],NULL);
5679: PetscSectionGetFieldPointSyms(csection,f,numCPoints,cpoints,&permsC[f],NULL);
5680: }
5681: for (p = 0; p < numFPoints; p++) {
5682: PetscSectionGetOffset(globalFSection, ftotpoints[2*p], &globalOff);
5683: DMPlexGetIndicesPointFields_Internal(fsection, ftotpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, foffsets, PETSC_FALSE, permsF, p, findices);
5684: }
5685: for (p = 0; p < numCPoints; p++) {
5686: PetscSectionGetOffset(globalCSection, cpoints[2*p], &globalOff);
5687: DMPlexGetIndicesPointFields_Internal(csection, cpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, coffsets, PETSC_FALSE, permsC, p, cindices);
5688: }
5689: for (f = 0; f < numFields; f++) {
5690: PetscSectionRestoreFieldPointSyms(fsection,f,numFPoints,ftotpoints,&permsF[f],NULL);
5691: PetscSectionRestoreFieldPointSyms(csection,f,numCPoints,cpoints,&permsC[f],NULL);
5692: }
5693: } else {
5694: const PetscInt **permsF = NULL;
5695: const PetscInt **permsC = NULL;
5697: PetscSectionGetPointSyms(fsection,numFPoints,ftotpoints,&permsF,NULL);
5698: PetscSectionGetPointSyms(csection,numCPoints,cpoints,&permsC,NULL);
5699: for (p = 0, off = 0; p < numFPoints; p++) {
5700: const PetscInt *perm = permsF ? permsF[p] : NULL;
5702: PetscSectionGetOffset(globalFSection, ftotpoints[2*p], &globalOff);
5703: DMPlexGetIndicesPoint_Internal(fsection, ftotpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, &off, PETSC_FALSE, perm, findices);
5704: }
5705: for (p = 0, off = 0; p < numCPoints; p++) {
5706: const PetscInt *perm = permsC ? permsC[p] : NULL;
5708: PetscSectionGetOffset(globalCSection, cpoints[2*p], &globalOff);
5709: DMPlexGetIndicesPoint_Internal(csection, cpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, &off, PETSC_FALSE, perm, cindices);
5710: }
5711: PetscSectionRestorePointSyms(fsection,numFPoints,ftotpoints,&permsF,NULL);
5712: PetscSectionRestorePointSyms(csection,numCPoints,cpoints,&permsC,NULL);
5713: }
5714: if (mesh->printSetValues) {DMPlexPrintMatSetValues(PETSC_VIEWER_STDOUT_SELF, A, point, numFIndices, findices, numCIndices, cindices, values);}
5715: /* TODO: flips */
5716: MatSetValues(A, numFIndices, findices, numCIndices, cindices, values, mode);
5717: if (ierr) {
5718: PetscMPIInt rank;
5719: PetscErrorCode ierr2;
5721: ierr2 = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank);CHKERRQ(ierr2);
5722: ierr2 = (*PetscErrorPrintf)("[%d]ERROR in DMPlexMatSetClosure\n", rank);CHKERRQ(ierr2);
5723: ierr2 = DMPlexPrintMatSetValues(PETSC_VIEWER_STDERR_SELF, A, point, numFIndices, findices, numCIndices, cindices, values);CHKERRQ(ierr2);
5724: ierr2 = DMRestoreWorkArray(dmf, numFIndices, PETSC_INT, &findices);CHKERRQ(ierr2);
5725: ierr2 = DMRestoreWorkArray(dmc, numCIndices, PETSC_INT, &cindices);CHKERRQ(ierr2);
5726:
5727: }
5728: DMRestoreWorkArray(dmf, numCPoints*2*4, PETSC_INT, &ftotpoints);
5729: DMPlexRestoreTransitiveClosure(dmc, point, PETSC_TRUE, &numCPoints, &cpoints);
5730: DMRestoreWorkArray(dmf, numFIndices, PETSC_INT, &findices);
5731: DMRestoreWorkArray(dmc, numCIndices, PETSC_INT, &cindices);
5732: return(0);
5733: }
5735: PetscErrorCode DMPlexMatGetClosureIndicesRefined(DM dmf, PetscSection fsection, PetscSection globalFSection, DM dmc, PetscSection csection, PetscSection globalCSection, PetscInt point, PetscInt cindices[], PetscInt findices[])
5736: {
5737: PetscInt *fpoints = NULL, *ftotpoints = NULL;
5738: PetscInt *cpoints = NULL;
5739: PetscInt foffsets[32], coffsets[32];
5740: CellRefiner cellRefiner;
5741: PetscInt numFields, numSubcells, maxFPoints, numFPoints, numCPoints, numFIndices, numCIndices, dof, off, globalOff, pStart, pEnd, p, q, r, s, f;
5747: if (!fsection) {DMGetDefaultSection(dmf, &fsection);}
5749: if (!csection) {DMGetDefaultSection(dmc, &csection);}
5751: if (!globalFSection) {DMGetDefaultGlobalSection(dmf, &globalFSection);}
5753: if (!globalCSection) {DMGetDefaultGlobalSection(dmc, &globalCSection);}
5755: PetscSectionGetNumFields(fsection, &numFields);
5756: if (numFields > 31) SETERRQ1(PetscObjectComm((PetscObject)dmf), PETSC_ERR_ARG_OUTOFRANGE, "Number of fields %D limited to 31", numFields);
5757: PetscMemzero(foffsets, 32 * sizeof(PetscInt));
5758: PetscMemzero(coffsets, 32 * sizeof(PetscInt));
5759: /* Column indices */
5760: DMPlexGetTransitiveClosure(dmc, point, PETSC_TRUE, &numCPoints, &cpoints);
5761: maxFPoints = numCPoints;
5762: /* Compress out points not in the section */
5763: /* TODO: Squeeze out points with 0 dof as well */
5764: PetscSectionGetChart(csection, &pStart, &pEnd);
5765: for (p = 0, q = 0; p < numCPoints*2; p += 2) {
5766: if ((cpoints[p] >= pStart) && (cpoints[p] < pEnd)) {
5767: cpoints[q*2] = cpoints[p];
5768: cpoints[q*2+1] = cpoints[p+1];
5769: ++q;
5770: }
5771: }
5772: numCPoints = q;
5773: for (p = 0, numCIndices = 0; p < numCPoints*2; p += 2) {
5774: PetscInt fdof;
5776: PetscSectionGetDof(csection, cpoints[p], &dof);
5777: if (!dof) continue;
5778: for (f = 0; f < numFields; ++f) {
5779: PetscSectionGetFieldDof(csection, cpoints[p], f, &fdof);
5780: coffsets[f+1] += fdof;
5781: }
5782: numCIndices += dof;
5783: }
5784: for (f = 1; f < numFields; ++f) coffsets[f+1] += coffsets[f];
5785: /* Row indices */
5786: DMPlexGetCellRefiner_Internal(dmc, &cellRefiner);
5787: CellRefinerGetAffineTransforms_Internal(cellRefiner, &numSubcells, NULL, NULL, NULL);
5788: DMGetWorkArray(dmf, maxFPoints*2*numSubcells, PETSC_INT, &ftotpoints);
5789: for (r = 0, q = 0; r < numSubcells; ++r) {
5790: /* TODO Map from coarse to fine cells */
5791: DMPlexGetTransitiveClosure(dmf, point*numSubcells + r, PETSC_TRUE, &numFPoints, &fpoints);
5792: /* Compress out points not in the section */
5793: PetscSectionGetChart(fsection, &pStart, &pEnd);
5794: for (p = 0; p < numFPoints*2; p += 2) {
5795: if ((fpoints[p] >= pStart) && (fpoints[p] < pEnd)) {
5796: PetscSectionGetDof(fsection, fpoints[p], &dof);
5797: if (!dof) continue;
5798: for (s = 0; s < q; ++s) if (fpoints[p] == ftotpoints[s*2]) break;
5799: if (s < q) continue;
5800: ftotpoints[q*2] = fpoints[p];
5801: ftotpoints[q*2+1] = fpoints[p+1];
5802: ++q;
5803: }
5804: }
5805: DMPlexRestoreTransitiveClosure(dmf, point, PETSC_TRUE, &numFPoints, &fpoints);
5806: }
5807: numFPoints = q;
5808: for (p = 0, numFIndices = 0; p < numFPoints*2; p += 2) {
5809: PetscInt fdof;
5811: PetscSectionGetDof(fsection, ftotpoints[p], &dof);
5812: if (!dof) continue;
5813: for (f = 0; f < numFields; ++f) {
5814: PetscSectionGetFieldDof(fsection, ftotpoints[p], f, &fdof);
5815: foffsets[f+1] += fdof;
5816: }
5817: numFIndices += dof;
5818: }
5819: for (f = 1; f < numFields; ++f) foffsets[f+1] += foffsets[f];
5821: if (numFields && foffsets[numFields] != numFIndices) SETERRQ2(PetscObjectComm((PetscObject)dmf), PETSC_ERR_PLIB, "Invalid size for closure %D should be %D", foffsets[numFields], numFIndices);
5822: if (numFields && coffsets[numFields] != numCIndices) SETERRQ2(PetscObjectComm((PetscObject)dmc), PETSC_ERR_PLIB, "Invalid size for closure %D should be %D", coffsets[numFields], numCIndices);
5823: if (numFields) {
5824: const PetscInt **permsF[32] = {NULL};
5825: const PetscInt **permsC[32] = {NULL};
5827: for (f = 0; f < numFields; f++) {
5828: PetscSectionGetFieldPointSyms(fsection,f,numFPoints,ftotpoints,&permsF[f],NULL);
5829: PetscSectionGetFieldPointSyms(csection,f,numCPoints,cpoints,&permsC[f],NULL);
5830: }
5831: for (p = 0; p < numFPoints; p++) {
5832: PetscSectionGetOffset(globalFSection, ftotpoints[2*p], &globalOff);
5833: DMPlexGetIndicesPointFields_Internal(fsection, ftotpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, foffsets, PETSC_FALSE, permsF, p, findices);
5834: }
5835: for (p = 0; p < numCPoints; p++) {
5836: PetscSectionGetOffset(globalCSection, cpoints[2*p], &globalOff);
5837: DMPlexGetIndicesPointFields_Internal(csection, cpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, coffsets, PETSC_FALSE, permsC, p, cindices);
5838: }
5839: for (f = 0; f < numFields; f++) {
5840: PetscSectionRestoreFieldPointSyms(fsection,f,numFPoints,ftotpoints,&permsF[f],NULL);
5841: PetscSectionRestoreFieldPointSyms(csection,f,numCPoints,cpoints,&permsC[f],NULL);
5842: }
5843: } else {
5844: const PetscInt **permsF = NULL;
5845: const PetscInt **permsC = NULL;
5847: PetscSectionGetPointSyms(fsection,numFPoints,ftotpoints,&permsF,NULL);
5848: PetscSectionGetPointSyms(csection,numCPoints,cpoints,&permsC,NULL);
5849: for (p = 0, off = 0; p < numFPoints; p++) {
5850: const PetscInt *perm = permsF ? permsF[p] : NULL;
5852: PetscSectionGetOffset(globalFSection, ftotpoints[2*p], &globalOff);
5853: DMPlexGetIndicesPoint_Internal(fsection, ftotpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, &off, PETSC_FALSE, perm, findices);
5854: }
5855: for (p = 0, off = 0; p < numCPoints; p++) {
5856: const PetscInt *perm = permsC ? permsC[p] : NULL;
5858: PetscSectionGetOffset(globalCSection, cpoints[2*p], &globalOff);
5859: DMPlexGetIndicesPoint_Internal(csection, cpoints[2*p], globalOff < 0 ? -(globalOff+1) : globalOff, &off, PETSC_FALSE, perm, cindices);
5860: }
5861: PetscSectionRestorePointSyms(fsection,numFPoints,ftotpoints,&permsF,NULL);
5862: PetscSectionRestorePointSyms(csection,numCPoints,cpoints,&permsC,NULL);
5863: }
5864: DMRestoreWorkArray(dmf, numCPoints*2*4, PETSC_INT, &ftotpoints);
5865: DMPlexRestoreTransitiveClosure(dmc, point, PETSC_TRUE, &numCPoints, &cpoints);
5866: return(0);
5867: }
5869: /*@
5870: DMPlexGetHybridBounds - Get the first mesh point of each dimension which is a hybrid
5872: Input Parameter:
5873: . dm - The DMPlex object
5875: Output Parameters:
5876: + cMax - The first hybrid cell
5877: . fMax - The first hybrid face
5878: . eMax - The first hybrid edge
5879: - vMax - The first hybrid vertex
5881: Level: developer
5883: .seealso DMPlexCreateHybridMesh(), DMPlexSetHybridBounds()
5884: @*/
5885: PetscErrorCode DMPlexGetHybridBounds(DM dm, PetscInt *cMax, PetscInt *fMax, PetscInt *eMax, PetscInt *vMax)
5886: {
5887: DM_Plex *mesh = (DM_Plex*) dm->data;
5888: PetscInt dim;
5893: DMGetDimension(dm, &dim);
5894: if (cMax) *cMax = mesh->hybridPointMax[dim];
5895: if (fMax) *fMax = mesh->hybridPointMax[dim-1];
5896: if (eMax) *eMax = mesh->hybridPointMax[1];
5897: if (vMax) *vMax = mesh->hybridPointMax[0];
5898: return(0);
5899: }
5901: /*@
5902: DMPlexSetHybridBounds - Set the first mesh point of each dimension which is a hybrid
5904: Input Parameters:
5905: . dm - The DMPlex object
5906: . cMax - The first hybrid cell
5907: . fMax - The first hybrid face
5908: . eMax - The first hybrid edge
5909: - vMax - The first hybrid vertex
5911: Level: developer
5913: .seealso DMPlexCreateHybridMesh(), DMPlexGetHybridBounds()
5914: @*/
5915: PetscErrorCode DMPlexSetHybridBounds(DM dm, PetscInt cMax, PetscInt fMax, PetscInt eMax, PetscInt vMax)
5916: {
5917: DM_Plex *mesh = (DM_Plex*) dm->data;
5918: PetscInt dim;
5923: DMGetDimension(dm, &dim);
5924: if (cMax >= 0) mesh->hybridPointMax[dim] = cMax;
5925: if (fMax >= 0) mesh->hybridPointMax[dim-1] = fMax;
5926: if (eMax >= 0) mesh->hybridPointMax[1] = eMax;
5927: if (vMax >= 0) mesh->hybridPointMax[0] = vMax;
5928: return(0);
5929: }
5931: /*@C
5932: DMPlexGetVTKCellHeight - Returns the height in the DAG used to determine which points are cells (normally 0)
5934: Input Parameter:
5935: . dm - The DMPlex object
5937: Output Parameter:
5938: . cellHeight - The height of a cell
5940: Level: developer
5942: .seealso DMPlexSetVTKCellHeight()
5943: @*/
5944: PetscErrorCode DMPlexGetVTKCellHeight(DM dm, PetscInt *cellHeight)
5945: {
5946: DM_Plex *mesh = (DM_Plex*) dm->data;
5951: *cellHeight = mesh->vtkCellHeight;
5952: return(0);
5953: }
5955: /*@C
5956: DMPlexSetVTKCellHeight - Sets the height in the DAG used to determine which points are cells (normally 0)
5958: Input Parameters:
5959: + dm - The DMPlex object
5960: - cellHeight - The height of a cell
5962: Level: developer
5964: .seealso DMPlexGetVTKCellHeight()
5965: @*/
5966: PetscErrorCode DMPlexSetVTKCellHeight(DM dm, PetscInt cellHeight)
5967: {
5968: DM_Plex *mesh = (DM_Plex*) dm->data;
5972: mesh->vtkCellHeight = cellHeight;
5973: return(0);
5974: }
5976: /* We can easily have a form that takes an IS instead */
5977: static PetscErrorCode DMPlexCreateNumbering_Private(DM dm, PetscInt pStart, PetscInt pEnd, PetscInt shift, PetscInt *globalSize, PetscSF sf, IS *numbering)
5978: {
5979: PetscSection section, globalSection;
5980: PetscInt *numbers, p;
5984: PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion);
5985: PetscSectionSetChart(section, pStart, pEnd);
5986: for (p = pStart; p < pEnd; ++p) {
5987: PetscSectionSetDof(section, p, 1);
5988: }
5989: PetscSectionSetUp(section);
5990: PetscSectionCreateGlobalSection(section, sf, PETSC_FALSE, PETSC_FALSE, &globalSection);
5991: PetscMalloc1(pEnd - pStart, &numbers);
5992: for (p = pStart; p < pEnd; ++p) {
5993: PetscSectionGetOffset(globalSection, p, &numbers[p-pStart]);
5994: if (numbers[p-pStart] < 0) numbers[p-pStart] -= shift;
5995: else numbers[p-pStart] += shift;
5996: }
5997: ISCreateGeneral(PetscObjectComm((PetscObject) dm), pEnd - pStart, numbers, PETSC_OWN_POINTER, numbering);
5998: if (globalSize) {
5999: PetscLayout layout;
6000: PetscSectionGetPointLayout(PetscObjectComm((PetscObject) dm), globalSection, &layout);
6001: PetscLayoutGetSize(layout, globalSize);
6002: PetscLayoutDestroy(&layout);
6003: }
6004: PetscSectionDestroy(§ion);
6005: PetscSectionDestroy(&globalSection);
6006: return(0);
6007: }
6009: PetscErrorCode DMPlexCreateCellNumbering_Internal(DM dm, PetscBool includeHybrid, IS *globalCellNumbers)
6010: {
6011: PetscInt cellHeight, cStart, cEnd, cMax;
6015: DMPlexGetVTKCellHeight(dm, &cellHeight);
6016: DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);
6017: DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);
6018: if (cMax >= 0 && !includeHybrid) cEnd = PetscMin(cEnd, cMax);
6019: DMPlexCreateNumbering_Private(dm, cStart, cEnd, 0, NULL, dm->sf, globalCellNumbers);
6020: return(0);
6021: }
6023: /*@C
6024: DMPlexGetCellNumbering - Get a global cell numbering for all cells on this process
6026: Input Parameter:
6027: . dm - The DMPlex object
6029: Output Parameter:
6030: . globalCellNumbers - Global cell numbers for all cells on this process
6032: Level: developer
6034: .seealso DMPlexGetVertexNumbering()
6035: @*/
6036: PetscErrorCode DMPlexGetCellNumbering(DM dm, IS *globalCellNumbers)
6037: {
6038: DM_Plex *mesh = (DM_Plex*) dm->data;
6043: if (!mesh->globalCellNumbers) {DMPlexCreateCellNumbering_Internal(dm, PETSC_FALSE, &mesh->globalCellNumbers);}
6044: *globalCellNumbers = mesh->globalCellNumbers;
6045: return(0);
6046: }
6048: PetscErrorCode DMPlexCreateVertexNumbering_Internal(DM dm, PetscBool includeHybrid, IS *globalVertexNumbers)
6049: {
6050: PetscInt vStart, vEnd, vMax;
6055: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
6056: DMPlexGetHybridBounds(dm, NULL, NULL, NULL, &vMax);
6057: if (vMax >= 0 && !includeHybrid) vEnd = PetscMin(vEnd, vMax);
6058: DMPlexCreateNumbering_Private(dm, vStart, vEnd, 0, NULL, dm->sf, globalVertexNumbers);
6059: return(0);
6060: }
6062: /*@C
6063: DMPlexGetVertexNumbering - Get a global certex numbering for all vertices on this process
6065: Input Parameter:
6066: . dm - The DMPlex object
6068: Output Parameter:
6069: . globalVertexNumbers - Global vertex numbers for all vertices on this process
6071: Level: developer
6073: .seealso DMPlexGetCellNumbering()
6074: @*/
6075: PetscErrorCode DMPlexGetVertexNumbering(DM dm, IS *globalVertexNumbers)
6076: {
6077: DM_Plex *mesh = (DM_Plex*) dm->data;
6082: if (!mesh->globalVertexNumbers) {DMPlexCreateVertexNumbering_Internal(dm, PETSC_FALSE, &mesh->globalVertexNumbers);}
6083: *globalVertexNumbers = mesh->globalVertexNumbers;
6084: return(0);
6085: }
6087: /*@C
6088: DMPlexCreatePointNumbering - Create a global numbering for all points on this process
6090: Input Parameter:
6091: . dm - The DMPlex object
6093: Output Parameter:
6094: . globalPointNumbers - Global numbers for all points on this process
6096: Level: developer
6098: .seealso DMPlexGetCellNumbering()
6099: @*/
6100: PetscErrorCode DMPlexCreatePointNumbering(DM dm, IS *globalPointNumbers)
6101: {
6102: IS nums[4];
6103: PetscInt depths[4];
6104: PetscInt depth, d, shift = 0;
6109: DMPlexGetDepth(dm, &depth);
6110: /* For unstratified meshes use dim instead of depth */
6111: if (depth < 0) {DMGetDimension(dm, &depth);}
6112: depths[0] = depth; depths[1] = 0;
6113: for (d = 2; d <= depth; ++d) depths[d] = depth-d+1;
6114: for (d = 0; d <= depth; ++d) {
6115: PetscInt pStart, pEnd, gsize;
6117: DMPlexGetDepthStratum(dm, depths[d], &pStart, &pEnd);
6118: DMPlexCreateNumbering_Private(dm, pStart, pEnd, shift, &gsize, dm->sf, &nums[d]);
6119: shift += gsize;
6120: }
6121: ISConcatenate(PetscObjectComm((PetscObject) dm), depth+1, nums, globalPointNumbers);
6122: for (d = 0; d <= depth; ++d) {ISDestroy(&nums[d]);}
6123: return(0);
6124: }
6127: /*@
6128: DMPlexCreateRankField - Create a cell field whose value is the rank of the owner
6130: Input Parameter:
6131: . dm - The DMPlex object
6133: Output Parameter:
6134: . ranks - The rank field
6136: Options Database Keys:
6137: . -dm_partition_view - Adds the rank field into the DM output from -dm_view using the same viewer
6139: Level: intermediate
6141: .seealso: DMView()
6142: @*/
6143: PetscErrorCode DMPlexCreateRankField(DM dm, Vec *ranks)
6144: {
6145: DM rdm;
6146: PetscDS prob;
6147: PetscFE fe;
6148: PetscScalar *r;
6149: PetscMPIInt rank;
6150: PetscInt dim, cStart, cEnd, c;
6154: MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
6155: DMClone(dm, &rdm);
6156: DMGetDimension(rdm, &dim);
6157: PetscFECreateDefault(rdm, dim, 1, PETSC_TRUE, NULL, -1, &fe);
6158: PetscObjectSetName((PetscObject) fe, "rank");
6159: DMGetDS(rdm, &prob);
6160: PetscDSSetDiscretization(prob, 0, (PetscObject) fe);
6161: PetscFEDestroy(&fe);
6162: DMPlexGetHeightStratum(rdm, 0, &cStart, &cEnd);
6163: DMCreateGlobalVector(rdm, ranks);
6164: PetscObjectSetName((PetscObject) *ranks, "partition");
6165: VecGetArray(*ranks, &r);
6166: for (c = cStart; c < cEnd; ++c) {
6167: PetscScalar *lr;
6169: DMPlexPointGlobalRef(rdm, c, r, &lr);
6170: *lr = rank;
6171: }
6172: VecRestoreArray(*ranks, &r);
6173: DMDestroy(&rdm);
6174: return(0);
6175: }
6177: /*@
6178: DMPlexCheckSymmetry - Check that the adjacency information in the mesh is symmetric.
6180: Input Parameter:
6181: . dm - The DMPlex object
6183: Note: This is a useful diagnostic when creating meshes programmatically.
6185: Level: developer
6187: .seealso: DMCreate(), DMCheckSkeleton(), DMCheckFaces()
6188: @*/
6189: PetscErrorCode DMPlexCheckSymmetry(DM dm)
6190: {
6191: PetscSection coneSection, supportSection;
6192: const PetscInt *cone, *support;
6193: PetscInt coneSize, c, supportSize, s;
6194: PetscInt pStart, pEnd, p, csize, ssize;
6195: PetscErrorCode ierr;
6199: DMPlexGetConeSection(dm, &coneSection);
6200: DMPlexGetSupportSection(dm, &supportSection);
6201: /* Check that point p is found in the support of its cone points, and vice versa */
6202: DMPlexGetChart(dm, &pStart, &pEnd);
6203: for (p = pStart; p < pEnd; ++p) {
6204: DMPlexGetConeSize(dm, p, &coneSize);
6205: DMPlexGetCone(dm, p, &cone);
6206: for (c = 0; c < coneSize; ++c) {
6207: PetscBool dup = PETSC_FALSE;
6208: PetscInt d;
6209: for (d = c-1; d >= 0; --d) {
6210: if (cone[c] == cone[d]) {dup = PETSC_TRUE; break;}
6211: }
6212: DMPlexGetSupportSize(dm, cone[c], &supportSize);
6213: DMPlexGetSupport(dm, cone[c], &support);
6214: for (s = 0; s < supportSize; ++s) {
6215: if (support[s] == p) break;
6216: }
6217: if ((s >= supportSize) || (dup && (support[s+1] != p))) {
6218: PetscPrintf(PETSC_COMM_SELF, "p: %D cone: ", p);
6219: for (s = 0; s < coneSize; ++s) {
6220: PetscPrintf(PETSC_COMM_SELF, "%D, ", cone[s]);
6221: }
6222: PetscPrintf(PETSC_COMM_SELF, "\n");
6223: PetscPrintf(PETSC_COMM_SELF, "p: %D support: ", cone[c]);
6224: for (s = 0; s < supportSize; ++s) {
6225: PetscPrintf(PETSC_COMM_SELF, "%D, ", support[s]);
6226: }
6227: PetscPrintf(PETSC_COMM_SELF, "\n");
6228: if (dup) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D not repeatedly found in support of repeated cone point %D", p, cone[c]);
6229: else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D not found in support of cone point %D", p, cone[c]);
6230: }
6231: }
6232: DMPlexGetSupportSize(dm, p, &supportSize);
6233: DMPlexGetSupport(dm, p, &support);
6234: for (s = 0; s < supportSize; ++s) {
6235: DMPlexGetConeSize(dm, support[s], &coneSize);
6236: DMPlexGetCone(dm, support[s], &cone);
6237: for (c = 0; c < coneSize; ++c) {
6238: if (cone[c] == p) break;
6239: }
6240: if (c >= coneSize) {
6241: PetscPrintf(PETSC_COMM_SELF, "p: %D support: ", p);
6242: for (c = 0; c < supportSize; ++c) {
6243: PetscPrintf(PETSC_COMM_SELF, "%D, ", support[c]);
6244: }
6245: PetscPrintf(PETSC_COMM_SELF, "\n");
6246: PetscPrintf(PETSC_COMM_SELF, "p: %D cone: ", support[s]);
6247: for (c = 0; c < coneSize; ++c) {
6248: PetscPrintf(PETSC_COMM_SELF, "%D, ", cone[c]);
6249: }
6250: PetscPrintf(PETSC_COMM_SELF, "\n");
6251: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D not found in cone of support point %D", p, support[s]);
6252: }
6253: }
6254: }
6255: PetscSectionGetStorageSize(coneSection, &csize);
6256: PetscSectionGetStorageSize(supportSection, &ssize);
6257: if (csize != ssize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Total cone size %D != Total support size %D", csize, ssize);
6258: return(0);
6259: }
6261: /*@
6262: DMPlexCheckSkeleton - Check that each cell has the correct number of vertices
6264: Input Parameters:
6265: + dm - The DMPlex object
6266: . isSimplex - Are the cells simplices or tensor products
6267: - cellHeight - Normally 0
6269: Note: This is a useful diagnostic when creating meshes programmatically.
6271: Level: developer
6273: .seealso: DMCreate(), DMCheckSymmetry(), DMCheckFaces()
6274: @*/
6275: PetscErrorCode DMPlexCheckSkeleton(DM dm, PetscBool isSimplex, PetscInt cellHeight)
6276: {
6277: PetscInt dim, numCorners, numHybridCorners, vStart, vEnd, cStart, cEnd, cMax, c;
6282: DMGetDimension(dm, &dim);
6283: switch (dim) {
6284: case 1: numCorners = isSimplex ? 2 : 2; numHybridCorners = isSimplex ? 2 : 2; break;
6285: case 2: numCorners = isSimplex ? 3 : 4; numHybridCorners = isSimplex ? 4 : 4; break;
6286: case 3: numCorners = isSimplex ? 4 : 8; numHybridCorners = isSimplex ? 6 : 8; break;
6287: default:
6288: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle meshes of dimension %D", dim);
6289: }
6290: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
6291: DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);
6292: DMPlexGetHybridBounds(dm, &cMax, NULL, NULL, NULL);
6293: cMax = cMax >= 0 ? cMax : cEnd;
6294: for (c = cStart; c < cMax; ++c) {
6295: PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
6297: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
6298: for (cl = 0; cl < closureSize*2; cl += 2) {
6299: const PetscInt p = closure[cl];
6300: if ((p >= vStart) && (p < vEnd)) ++coneSize;
6301: }
6302: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
6303: if (coneSize != numCorners) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %D has %D vertices != %D", c, coneSize, numCorners);
6304: }
6305: for (c = cMax; c < cEnd; ++c) {
6306: PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
6308: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
6309: for (cl = 0; cl < closureSize*2; cl += 2) {
6310: const PetscInt p = closure[cl];
6311: if ((p >= vStart) && (p < vEnd)) ++coneSize;
6312: }
6313: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
6314: if (coneSize > numHybridCorners) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Hybrid cell %D has %D vertices > %D", c, coneSize, numHybridCorners);
6315: }
6316: return(0);
6317: }
6319: /*@
6320: DMPlexCheckFaces - Check that the faces of each cell give a vertex order this is consistent with what we expect from the cell type
6322: Input Parameters:
6323: + dm - The DMPlex object
6324: . isSimplex - Are the cells simplices or tensor products
6325: - cellHeight - Normally 0
6327: Note: This is a useful diagnostic when creating meshes programmatically.
6329: Level: developer
6331: .seealso: DMCreate(), DMCheckSymmetry(), DMCheckSkeleton()
6332: @*/
6333: PetscErrorCode DMPlexCheckFaces(DM dm, PetscBool isSimplex, PetscInt cellHeight)
6334: {
6335: PetscInt pMax[4];
6336: PetscInt dim, vStart, vEnd, cStart, cEnd, c, h;
6341: DMGetDimension(dm, &dim);
6342: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
6343: DMPlexGetHybridBounds(dm, &pMax[dim], &pMax[dim-1], &pMax[1], &pMax[0]);
6344: for (h = cellHeight; h < dim; ++h) {
6345: DMPlexGetHeightStratum(dm, h, &cStart, &cEnd);
6346: for (c = cStart; c < cEnd; ++c) {
6347: const PetscInt *cone, *ornt, *faces;
6348: PetscInt numFaces, faceSize, coneSize,f;
6349: PetscInt *closure = NULL, closureSize, cl, numCorners = 0;
6351: if (pMax[dim-h] >= 0 && c >= pMax[dim-h]) continue;
6352: DMPlexGetConeSize(dm, c, &coneSize);
6353: DMPlexGetCone(dm, c, &cone);
6354: DMPlexGetConeOrientation(dm, c, &ornt);
6355: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
6356: for (cl = 0; cl < closureSize*2; cl += 2) {
6357: const PetscInt p = closure[cl];
6358: if ((p >= vStart) && (p < vEnd)) closure[numCorners++] = p;
6359: }
6360: DMPlexGetRawFaces_Internal(dm, dim-h, numCorners, closure, &numFaces, &faceSize, &faces);
6361: if (coneSize != numFaces) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %D has %D faces but should have %D", c, coneSize, numFaces);
6362: for (f = 0; f < numFaces; ++f) {
6363: PetscInt *fclosure = NULL, fclosureSize, cl, fnumCorners = 0, v;
6365: DMPlexGetTransitiveClosure_Internal(dm, cone[f], ornt[f], PETSC_TRUE, &fclosureSize, &fclosure);
6366: for (cl = 0; cl < fclosureSize*2; cl += 2) {
6367: const PetscInt p = fclosure[cl];
6368: if ((p >= vStart) && (p < vEnd)) fclosure[fnumCorners++] = p;
6369: }
6370: if (fnumCorners != faceSize) SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %D (%D) of cell %D has %D vertices but should have %D", cone[f], f, c, fnumCorners, faceSize);
6371: for (v = 0; v < fnumCorners; ++v) {
6372: if (fclosure[v] != faces[f*faceSize+v]) SETERRQ6(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %D (%d) of cell %D vertex %D, %D != %D", cone[f], f, c, v, fclosure[v], faces[f*faceSize+v]);
6373: }
6374: DMPlexRestoreTransitiveClosure(dm, cone[f], PETSC_TRUE, &fclosureSize, &fclosure);
6375: }
6376: DMPlexRestoreFaces_Internal(dm, dim, c, &numFaces, &faceSize, &faces);
6377: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
6378: }
6379: }
6380: return(0);
6381: }
6383: /* Pointwise interpolation
6384: Just code FEM for now
6385: u^f = I u^c
6386: sum_k u^f_k phi^f_k = I sum_j u^c_j phi^c_j
6387: u^f_i = sum_j psi^f_i I phi^c_j u^c_j
6388: I_{ij} = psi^f_i phi^c_j
6389: */
6390: PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling)
6391: {
6392: PetscSection gsc, gsf;
6393: PetscInt m, n;
6394: void *ctx;
6395: DM cdm;
6396: PetscBool regular;
6400: DMGetDefaultGlobalSection(dmFine, &gsf);
6401: PetscSectionGetConstrainedStorageSize(gsf, &m);
6402: DMGetDefaultGlobalSection(dmCoarse, &gsc);
6403: PetscSectionGetConstrainedStorageSize(gsc, &n);
6405: MatCreate(PetscObjectComm((PetscObject) dmCoarse), interpolation);
6406: MatSetSizes(*interpolation, m, n, PETSC_DETERMINE, PETSC_DETERMINE);
6407: MatSetType(*interpolation, dmCoarse->mattype);
6408: DMGetApplicationContext(dmFine, &ctx);
6410: DMGetCoarseDM(dmFine, &cdm);
6411: DMPlexGetRegularRefinement(dmFine, ®ular);
6412: if (regular && cdm == dmCoarse) {DMPlexComputeInterpolatorNested(dmCoarse, dmFine, *interpolation, ctx);}
6413: else {DMPlexComputeInterpolatorGeneral(dmCoarse, dmFine, *interpolation, ctx);}
6414: MatViewFromOptions(*interpolation, NULL, "-interp_mat_view");
6415: /* Use naive scaling */
6416: DMCreateInterpolationScale(dmCoarse, dmFine, *interpolation, scaling);
6417: return(0);
6418: }
6420: PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat)
6421: {
6423: VecScatter ctx;
6426: DMPlexComputeInjectorFEM(dmCoarse, dmFine, &ctx, NULL);
6427: MatCreateScatter(PetscObjectComm((PetscObject)ctx), ctx, mat);
6428: VecScatterDestroy(&ctx);
6429: return(0);
6430: }
6432: PetscErrorCode DMCreateDefaultSection_Plex(DM dm)
6433: {
6434: PetscSection section;
6435: IS *bcPoints, *bcComps;
6436: PetscBool *isFE;
6437: PetscInt *bcFields, *numComp, *numDof;
6438: PetscInt depth, dim, numBd, numBC = 0, numFields, bd, bc = 0, f;
6439: PetscInt cStart, cEnd, cEndInterior;
6443: DMGetNumFields(dm, &numFields);
6444: if (!numFields) return(0);
6445: /* FE and FV boundary conditions are handled slightly differently */
6446: PetscMalloc1(numFields, &isFE);
6447: for (f = 0; f < numFields; ++f) {
6448: PetscObject obj;
6449: PetscClassId id;
6451: DMGetField(dm, f, &obj);
6452: PetscObjectGetClassId(obj, &id);
6453: if (id == PETSCFE_CLASSID) {isFE[f] = PETSC_TRUE;}
6454: else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
6455: else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f);
6456: }
6457: /* Allocate boundary point storage for FEM boundaries */
6458: DMPlexGetDepth(dm, &depth);
6459: DMGetDimension(dm, &dim);
6460: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
6461: DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);
6462: PetscDSGetNumBoundary(dm->prob, &numBd);
6463: for (bd = 0; bd < numBd; ++bd) {
6464: PetscInt field;
6465: DMBoundaryConditionType type;
6466: const char *labelName;
6467: DMLabel label;
6469: PetscDSGetBoundary(dm->prob, bd, &type, NULL, &labelName, &field, NULL, NULL, NULL, NULL, NULL, NULL);
6470: DMGetLabel(dm,labelName,&label);
6471: if (label && isFE[field] && (type & DM_BC_ESSENTIAL)) ++numBC;
6472: }
6473: /* Add ghost cell boundaries for FVM */
6474: for (f = 0; f < numFields; ++f) if (!isFE[f] && cEndInterior >= 0) ++numBC;
6475: PetscCalloc3(numBC,&bcFields,numBC,&bcPoints,numBC,&bcComps);
6476: /* Constrain ghost cells for FV */
6477: for (f = 0; f < numFields; ++f) {
6478: PetscInt *newidx, c;
6480: if (isFE[f] || cEndInterior < 0) continue;
6481: PetscMalloc1(cEnd-cEndInterior,&newidx);
6482: for (c = cEndInterior; c < cEnd; ++c) newidx[c-cEndInterior] = c;
6483: bcFields[bc] = f;
6484: ISCreateGeneral(PetscObjectComm((PetscObject) dm), cEnd-cEndInterior, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]);
6485: }
6486: /* Handle FEM Dirichlet boundaries */
6487: for (bd = 0; bd < numBd; ++bd) {
6488: const char *bdLabel;
6489: DMLabel label;
6490: const PetscInt *comps;
6491: const PetscInt *values;
6492: PetscInt bd2, field, numComps, numValues;
6493: DMBoundaryConditionType type;
6494: PetscBool duplicate = PETSC_FALSE;
6496: PetscDSGetBoundary(dm->prob, bd, &type, NULL, &bdLabel, &field, &numComps, &comps, NULL, &numValues, &values, NULL);
6497: DMGetLabel(dm, bdLabel, &label);
6498: if (!isFE[field] || !label) continue;
6499: /* Only want to modify label once */
6500: for (bd2 = 0; bd2 < bd; ++bd2) {
6501: const char *bdname;
6502: PetscDSGetBoundary(dm->prob, bd2, NULL, NULL, &bdname, NULL, NULL, NULL, NULL, NULL, NULL, NULL);
6503: PetscStrcmp(bdname, bdLabel, &duplicate);
6504: if (duplicate) break;
6505: }
6506: if (!duplicate && (isFE[field])) {
6507: /* don't complete cells, which are just present to give orientation to the boundary */
6508: DMPlexLabelComplete(dm, label);
6509: }
6510: /* Filter out cells, if you actually want to constrain cells you need to do things by hand right now */
6511: if (type & DM_BC_ESSENTIAL) {
6512: PetscInt *newidx;
6513: PetscInt n, newn = 0, p, v;
6515: bcFields[bc] = field;
6516: if (numComps) {ISCreateGeneral(PetscObjectComm((PetscObject) dm), numComps, comps, PETSC_COPY_VALUES, &bcComps[bc]);}
6517: for (v = 0; v < numValues; ++v) {
6518: IS tmp;
6519: const PetscInt *idx;
6521: DMGetStratumIS(dm, bdLabel, values[v], &tmp);
6522: if (!tmp) continue;
6523: ISGetLocalSize(tmp, &n);
6524: ISGetIndices(tmp, &idx);
6525: if (isFE[field]) {
6526: for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) ++newn;
6527: } else {
6528: for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) ++newn;
6529: }
6530: ISRestoreIndices(tmp, &idx);
6531: ISDestroy(&tmp);
6532: }
6533: PetscMalloc1(newn,&newidx);
6534: newn = 0;
6535: for (v = 0; v < numValues; ++v) {
6536: IS tmp;
6537: const PetscInt *idx;
6539: DMGetStratumIS(dm, bdLabel, values[v], &tmp);
6540: if (!tmp) continue;
6541: ISGetLocalSize(tmp, &n);
6542: ISGetIndices(tmp, &idx);
6543: if (isFE[field]) {
6544: for (p = 0; p < n; ++p) if ((idx[p] < cStart) || (idx[p] >= cEnd)) newidx[newn++] = idx[p];
6545: } else {
6546: for (p = 0; p < n; ++p) if ((idx[p] >= cStart) || (idx[p] < cEnd)) newidx[newn++] = idx[p];
6547: }
6548: ISRestoreIndices(tmp, &idx);
6549: ISDestroy(&tmp);
6550: }
6551: ISCreateGeneral(PetscObjectComm((PetscObject) dm), newn, newidx, PETSC_OWN_POINTER, &bcPoints[bc++]);
6552: }
6553: }
6554: /* Handle discretization */
6555: PetscCalloc2(numFields,&numComp,numFields*(dim+1),&numDof);
6556: for (f = 0; f < numFields; ++f) {
6557: PetscObject obj;
6559: DMGetField(dm, f, &obj);
6560: if (isFE[f]) {
6561: PetscFE fe = (PetscFE) obj;
6562: const PetscInt *numFieldDof;
6563: PetscInt d;
6565: PetscFEGetNumComponents(fe, &numComp[f]);
6566: PetscFEGetNumDof(fe, &numFieldDof);
6567: for (d = 0; d < dim+1; ++d) numDof[f*(dim+1)+d] = numFieldDof[d];
6568: } else {
6569: PetscFV fv = (PetscFV) obj;
6571: PetscFVGetNumComponents(fv, &numComp[f]);
6572: numDof[f*(dim+1)+dim] = numComp[f];
6573: }
6574: }
6575: for (f = 0; f < numFields; ++f) {
6576: PetscInt d;
6577: for (d = 1; d < dim; ++d) {
6578: if ((numDof[f*(dim+1)+d] > 0) && (depth < dim)) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated when unknowns are specified on edges or faces.");
6579: }
6580: }
6581: DMPlexCreateSection(dm, dim, numFields, numComp, numDof, numBC, bcFields, bcComps, bcPoints, NULL, §ion);
6582: for (f = 0; f < numFields; ++f) {
6583: PetscFE fe;
6584: const char *name;
6586: DMGetField(dm, f, (PetscObject *) &fe);
6587: PetscObjectGetName((PetscObject) fe, &name);
6588: PetscSectionSetFieldName(section, f, name);
6589: }
6590: DMSetDefaultSection(dm, section);
6591: PetscSectionDestroy(§ion);
6592: for (bc = 0; bc < numBC; ++bc) {ISDestroy(&bcPoints[bc]);ISDestroy(&bcComps[bc]);}
6593: PetscFree3(bcFields,bcPoints,bcComps);
6594: PetscFree2(numComp,numDof);
6595: PetscFree(isFE);
6596: return(0);
6597: }
6599: /*@
6600: DMPlexGetRegularRefinement - Get the flag indicating that this mesh was obtained by regular refinement from its coarse mesh
6602: Input Parameter:
6603: . dm - The DMPlex object
6605: Output Parameter:
6606: . regular - The flag
6608: Level: intermediate
6610: .seealso: DMPlexSetRegularRefinement()
6611: @*/
6612: PetscErrorCode DMPlexGetRegularRefinement(DM dm, PetscBool *regular)
6613: {
6617: *regular = ((DM_Plex *) dm->data)->regularRefinement;
6618: return(0);
6619: }
6621: /*@
6622: DMPlexSetRegularRefinement - Set the flag indicating that this mesh was obtained by regular refinement from its coarse mesh
6624: Input Parameters:
6625: + dm - The DMPlex object
6626: - regular - The flag
6628: Level: intermediate
6630: .seealso: DMPlexGetRegularRefinement()
6631: @*/
6632: PetscErrorCode DMPlexSetRegularRefinement(DM dm, PetscBool regular)
6633: {
6636: ((DM_Plex *) dm->data)->regularRefinement = regular;
6637: return(0);
6638: }
6640: /* anchors */
6641: /*@
6642: DMPlexGetAnchors - Get the layout of the anchor (point-to-point) constraints. Typically, the user will not have to
6643: call DMPlexGetAnchors() directly: if there are anchors, then DMPlexGetAnchors() is called during DMGetConstraints().
6645: not collective
6647: Input Parameters:
6648: . dm - The DMPlex object
6650: Output Parameters:
6651: + anchorSection - If not NULL, set to the section describing which points anchor the constrained points.
6652: - anchorIS - If not NULL, set to the list of anchors indexed by anchorSection
6655: Level: intermediate
6657: .seealso: DMPlexSetAnchors(), DMGetConstraints(), DMSetConstraints()
6658: @*/
6659: PetscErrorCode DMPlexGetAnchors(DM dm, PetscSection *anchorSection, IS *anchorIS)
6660: {
6661: DM_Plex *plex = (DM_Plex *)dm->data;
6666: if (!plex->anchorSection && !plex->anchorIS && plex->createanchors) {(*plex->createanchors)(dm);}
6667: if (anchorSection) *anchorSection = plex->anchorSection;
6668: if (anchorIS) *anchorIS = plex->anchorIS;
6669: return(0);
6670: }
6672: /*@
6673: DMPlexSetAnchors - Set the layout of the local anchor (point-to-point) constraints. Unlike boundary conditions,
6674: when a point's degrees of freedom in a section are constrained to an outside value, the anchor constraints set a
6675: point's degrees of freedom to be a linear combination of other points' degrees of freedom.
6677: After specifying the layout of constraints with DMPlexSetAnchors(), one specifies the constraints by calling
6678: DMGetConstraints() and filling in the entries in the constraint matrix.
6680: collective on dm
6682: Input Parameters:
6683: + dm - The DMPlex object
6684: . 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).
6685: - anchorIS - The list of all anchor points. Must have a local communicator (PETSC_COMM_SELF or derivative).
6687: The reference counts of anchorSection and anchorIS are incremented.
6689: Level: intermediate
6691: .seealso: DMPlexGetAnchors(), DMGetConstraints(), DMSetConstraints()
6692: @*/
6693: PetscErrorCode DMPlexSetAnchors(DM dm, PetscSection anchorSection, IS anchorIS)
6694: {
6695: DM_Plex *plex = (DM_Plex *)dm->data;
6696: PetscMPIInt result;
6701: if (anchorSection) {
6703: MPI_Comm_compare(PETSC_COMM_SELF,PetscObjectComm((PetscObject)anchorSection),&result);
6704: if (result != MPI_CONGRUENT && result != MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMECOMM,"anchor section must have local communicator");
6705: }
6706: if (anchorIS) {
6708: MPI_Comm_compare(PETSC_COMM_SELF,PetscObjectComm((PetscObject)anchorIS),&result);
6709: if (result != MPI_CONGRUENT && result != MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMECOMM,"anchor IS must have local communicator");
6710: }
6712: PetscObjectReference((PetscObject)anchorSection);
6713: PetscSectionDestroy(&plex->anchorSection);
6714: plex->anchorSection = anchorSection;
6716: PetscObjectReference((PetscObject)anchorIS);
6717: ISDestroy(&plex->anchorIS);
6718: plex->anchorIS = anchorIS;
6720: #if defined(PETSC_USE_DEBUG)
6721: if (anchorIS && anchorSection) {
6722: PetscInt size, a, pStart, pEnd;
6723: const PetscInt *anchors;
6725: PetscSectionGetChart(anchorSection,&pStart,&pEnd);
6726: ISGetLocalSize(anchorIS,&size);
6727: ISGetIndices(anchorIS,&anchors);
6728: for (a = 0; a < size; a++) {
6729: PetscInt p;
6731: p = anchors[a];
6732: if (p >= pStart && p < pEnd) {
6733: PetscInt dof;
6735: PetscSectionGetDof(anchorSection,p,&dof);
6736: if (dof) {
6737: PetscErrorCode ierr2;
6739: ierr2 = ISRestoreIndices(anchorIS,&anchors);CHKERRQ(ierr2);
6740: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Point %D cannot be constrained and an anchor",p);
6741: }
6742: }
6743: }
6744: ISRestoreIndices(anchorIS,&anchors);
6745: }
6746: #endif
6747: /* reset the generic constraints */
6748: DMSetDefaultConstraints(dm,NULL,NULL);
6749: return(0);
6750: }
6752: static PetscErrorCode DMPlexCreateConstraintSection_Anchors(DM dm, PetscSection section, PetscSection *cSec)
6753: {
6754: PetscSection anchorSection;
6755: PetscInt pStart, pEnd, sStart, sEnd, p, dof, numFields, f;
6760: DMPlexGetAnchors(dm,&anchorSection,NULL);
6761: PetscSectionCreate(PETSC_COMM_SELF,cSec);
6762: PetscSectionGetNumFields(section,&numFields);
6763: if (numFields) {
6764: PetscInt f;
6765: PetscSectionSetNumFields(*cSec,numFields);
6767: for (f = 0; f < numFields; f++) {
6768: PetscInt numComp;
6770: PetscSectionGetFieldComponents(section,f,&numComp);
6771: PetscSectionSetFieldComponents(*cSec,f,numComp);
6772: }
6773: }
6774: PetscSectionGetChart(anchorSection,&pStart,&pEnd);
6775: PetscSectionGetChart(section,&sStart,&sEnd);
6776: pStart = PetscMax(pStart,sStart);
6777: pEnd = PetscMin(pEnd,sEnd);
6778: pEnd = PetscMax(pStart,pEnd);
6779: PetscSectionSetChart(*cSec,pStart,pEnd);
6780: for (p = pStart; p < pEnd; p++) {
6781: PetscSectionGetDof(anchorSection,p,&dof);
6782: if (dof) {
6783: PetscSectionGetDof(section,p,&dof);
6784: PetscSectionSetDof(*cSec,p,dof);
6785: for (f = 0; f < numFields; f++) {
6786: PetscSectionGetFieldDof(section,p,f,&dof);
6787: PetscSectionSetFieldDof(*cSec,p,f,dof);
6788: }
6789: }
6790: }
6791: PetscSectionSetUp(*cSec);
6792: return(0);
6793: }
6795: static PetscErrorCode DMPlexCreateConstraintMatrix_Anchors(DM dm, PetscSection section, PetscSection cSec, Mat *cMat)
6796: {
6797: PetscSection aSec;
6798: PetscInt pStart, pEnd, p, dof, aDof, aOff, off, nnz, annz, m, n, q, a, offset, *i, *j;
6799: const PetscInt *anchors;
6800: PetscInt numFields, f;
6801: IS aIS;
6806: PetscSectionGetStorageSize(cSec, &m);
6807: PetscSectionGetStorageSize(section, &n);
6808: MatCreate(PETSC_COMM_SELF,cMat);
6809: MatSetSizes(*cMat,m,n,m,n);
6810: MatSetType(*cMat,MATSEQAIJ);
6811: DMPlexGetAnchors(dm,&aSec,&aIS);
6812: ISGetIndices(aIS,&anchors);
6813: /* cSec will be a subset of aSec and section */
6814: PetscSectionGetChart(cSec,&pStart,&pEnd);
6815: PetscMalloc1(m+1,&i);
6816: i[0] = 0;
6817: PetscSectionGetNumFields(section,&numFields);
6818: for (p = pStart; p < pEnd; p++) {
6819: PetscInt rDof, rOff, r;
6821: PetscSectionGetDof(aSec,p,&rDof);
6822: if (!rDof) continue;
6823: PetscSectionGetOffset(aSec,p,&rOff);
6824: if (numFields) {
6825: for (f = 0; f < numFields; f++) {
6826: annz = 0;
6827: for (r = 0; r < rDof; r++) {
6828: a = anchors[rOff + r];
6829: PetscSectionGetFieldDof(section,a,f,&aDof);
6830: annz += aDof;
6831: }
6832: PetscSectionGetFieldDof(cSec,p,f,&dof);
6833: PetscSectionGetFieldOffset(cSec,p,f,&off);
6834: for (q = 0; q < dof; q++) {
6835: i[off + q + 1] = i[off + q] + annz;
6836: }
6837: }
6838: }
6839: else {
6840: annz = 0;
6841: for (q = 0; q < dof; q++) {
6842: a = anchors[off + q];
6843: PetscSectionGetDof(section,a,&aDof);
6844: annz += aDof;
6845: }
6846: PetscSectionGetDof(cSec,p,&dof);
6847: PetscSectionGetOffset(cSec,p,&off);
6848: for (q = 0; q < dof; q++) {
6849: i[off + q + 1] = i[off + q] + annz;
6850: }
6851: }
6852: }
6853: nnz = i[m];
6854: PetscMalloc1(nnz,&j);
6855: offset = 0;
6856: for (p = pStart; p < pEnd; p++) {
6857: if (numFields) {
6858: for (f = 0; f < numFields; f++) {
6859: PetscSectionGetFieldDof(cSec,p,f,&dof);
6860: for (q = 0; q < dof; q++) {
6861: PetscInt rDof, rOff, r;
6862: PetscSectionGetDof(aSec,p,&rDof);
6863: PetscSectionGetOffset(aSec,p,&rOff);
6864: for (r = 0; r < rDof; r++) {
6865: PetscInt s;
6867: a = anchors[rOff + r];
6868: PetscSectionGetFieldDof(section,a,f,&aDof);
6869: PetscSectionGetFieldOffset(section,a,f,&aOff);
6870: for (s = 0; s < aDof; s++) {
6871: j[offset++] = aOff + s;
6872: }
6873: }
6874: }
6875: }
6876: }
6877: else {
6878: PetscSectionGetDof(cSec,p,&dof);
6879: for (q = 0; q < dof; q++) {
6880: PetscInt rDof, rOff, r;
6881: PetscSectionGetDof(aSec,p,&rDof);
6882: PetscSectionGetOffset(aSec,p,&rOff);
6883: for (r = 0; r < rDof; r++) {
6884: PetscInt s;
6886: a = anchors[rOff + r];
6887: PetscSectionGetDof(section,a,&aDof);
6888: PetscSectionGetOffset(section,a,&aOff);
6889: for (s = 0; s < aDof; s++) {
6890: j[offset++] = aOff + s;
6891: }
6892: }
6893: }
6894: }
6895: }
6896: MatSeqAIJSetPreallocationCSR(*cMat,i,j,NULL);
6897: PetscFree(i);
6898: PetscFree(j);
6899: ISRestoreIndices(aIS,&anchors);
6900: return(0);
6901: }
6903: PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm)
6904: {
6905: DM_Plex *plex = (DM_Plex *)dm->data;
6906: PetscSection anchorSection, section, cSec;
6907: Mat cMat;
6912: DMPlexGetAnchors(dm,&anchorSection,NULL);
6913: if (anchorSection) {
6914: PetscDS ds;
6915: PetscInt nf;
6917: DMGetDefaultSection(dm,§ion);
6918: DMPlexCreateConstraintSection_Anchors(dm,section,&cSec);
6919: DMPlexCreateConstraintMatrix_Anchors(dm,section,cSec,&cMat);
6920: DMGetDS(dm,&ds);
6921: PetscDSGetNumFields(ds,&nf);
6922: if (nf && plex->computeanchormatrix) {(*plex->computeanchormatrix)(dm,section,cSec,cMat);}
6923: DMSetDefaultConstraints(dm,cSec,cMat);
6924: PetscSectionDestroy(&cSec);
6925: MatDestroy(&cMat);
6926: }
6927: return(0);
6928: }