Actual source code: plexindices.c
petsc-3.7.3 2016-08-01
1: #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
5: /*@
6: DMPlexCreateClosureIndex - Calculate an index for the given PetscSection for the closure operation on the DM
8: Not collective
10: Input Parameters:
11: + dm - The DM
12: - section - The section describing the layout in v, or NULL to use the default section
14: Note:
15: This should greatly improve the performance of the closure operations, at the cost of additional memory.
17: Level: intermediate
19: .seealso DMPlexVecGetClosure(), DMPlexVecRestoreClosure(), DMPlexVecSetClosure(), DMPlexMatSetClosure()
20: @*/
21: PetscErrorCode DMPlexCreateClosureIndex(DM dm, PetscSection section)
22: {
23: PetscSection closureSection;
24: IS closureIS;
25: PetscInt *clPoints;
26: PetscInt pStart, pEnd, sStart, sEnd, point, clSize;
31: if (!section) {DMGetDefaultSection(dm, §ion);}
33: PetscSectionGetChart(section, &sStart, &sEnd);
34: DMPlexGetChart(dm, &pStart, &pEnd);
35: PetscSectionCreate(PetscObjectComm((PetscObject) section), &closureSection);
36: PetscSectionSetChart(closureSection, pStart, pEnd);
37: for (point = pStart; point < pEnd; ++point) {
38: PetscInt *points = NULL, numPoints, p, dof, cldof = 0;
40: DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &numPoints, &points);
41: for (p = 0; p < numPoints*2; p += 2) {
42: if ((points[p] >= sStart) && (points[p] < sEnd)) {
43: PetscSectionGetDof(section, points[p], &dof);
44: if (dof) cldof += 2;
45: }
46: }
47: DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &numPoints, &points);
48: PetscSectionSetDof(closureSection, point, cldof);
49: }
50: PetscSectionSetUp(closureSection);
51: PetscSectionGetStorageSize(closureSection, &clSize);
52: PetscMalloc1(clSize, &clPoints);
53: for (point = pStart; point < pEnd; ++point) {
54: PetscInt *points = NULL, numPoints, p, q, dof, cldof, cloff;
56: PetscSectionGetDof(closureSection, point, &cldof);
57: PetscSectionGetOffset(closureSection, point, &cloff);
58: DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &numPoints, &points);
59: for (p = 0, q = 0; p < numPoints*2; p += 2) {
60: if ((points[p] >= sStart) && (points[p] < sEnd)) {
61: PetscSectionGetDof(section, points[p], &dof);
62: if (dof) {
63: clPoints[cloff+q*2] = points[p];
64: clPoints[cloff+q*2+1] = points[p+1];
65: ++q;
66: }
67: }
68: }
69: DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &numPoints, &points);
70: if (q*2 != cldof) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Invalid size for closure %d should be %d", q*2, cldof);
71: }
72: ISCreateGeneral(PETSC_COMM_SELF, clSize, clPoints, PETSC_OWN_POINTER, &closureIS);
73: PetscSectionSetClosureIndex(section, (PetscObject) dm, closureSection, closureIS);
74: return(0);
75: }