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