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