xref: /petsc/src/dm/impls/plex/plexindices.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2d9917b9dSMatthew G. Knepley 
3d9917b9dSMatthew G. Knepley /*@
4d9917b9dSMatthew G. Knepley   DMPlexCreateClosureIndex - Calculate an index for the given PetscSection for the closure operation on the DM
5d9917b9dSMatthew G. Knepley 
6d9917b9dSMatthew G. Knepley   Not collective
7d9917b9dSMatthew G. Knepley 
8d9917b9dSMatthew G. Knepley   Input Parameters:
9d9917b9dSMatthew G. Knepley + dm - The DM
103e03ebd8SStefano Zampini - section - The section describing the layout in the local vector, or NULL to use the default section
11d9917b9dSMatthew G. Knepley 
12d9917b9dSMatthew G. Knepley   Note:
13d9917b9dSMatthew G. Knepley   This should greatly improve the performance of the closure operations, at the cost of additional memory.
14d9917b9dSMatthew G. Knepley 
15d9917b9dSMatthew G. Knepley   Level: intermediate
16d9917b9dSMatthew G. Knepley 
17d9917b9dSMatthew G. Knepley .seealso DMPlexVecGetClosure(), DMPlexVecRestoreClosure(), DMPlexVecSetClosure(), DMPlexMatSetClosure()
18d9917b9dSMatthew G. Knepley @*/
19d9917b9dSMatthew G. Knepley PetscErrorCode DMPlexCreateClosureIndex(DM dm, PetscSection section)
20d9917b9dSMatthew G. Knepley {
21d9917b9dSMatthew G. Knepley   PetscSection   closureSection;
22d9917b9dSMatthew G. Knepley   IS             closureIS;
23d9917b9dSMatthew G. Knepley   PetscInt      *clPoints;
24d9917b9dSMatthew G. Knepley   PetscInt       pStart, pEnd, sStart, sEnd, point, clSize;
25d9917b9dSMatthew G. Knepley 
26d9917b9dSMatthew G. Knepley   PetscFunctionBegin;
27d9917b9dSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
289566063dSJacob Faibussowitsch   if (!section) PetscCall(DMGetLocalSection(dm, &section));
291e0781a8SMatthew G. Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 2);
309566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &sStart, &sEnd));
319566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
329566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) section), &closureSection));
339566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(closureSection, pStart, pEnd));
34d9917b9dSMatthew G. Knepley   for (point = pStart; point < pEnd; ++point) {
35d9917b9dSMatthew G. Knepley     PetscInt *points = NULL, numPoints, p, dof, cldof = 0;
36d9917b9dSMatthew G. Knepley 
379566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &numPoints, &points));
38d9917b9dSMatthew G. Knepley     for (p = 0; p < numPoints*2; p += 2) {
39d9917b9dSMatthew G. Knepley       if ((points[p] >= sStart) && (points[p] < sEnd)) {
409566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, points[p], &dof));
41d9917b9dSMatthew G. Knepley         if (dof) cldof += 2;
42d9917b9dSMatthew G. Knepley       }
43d9917b9dSMatthew G. Knepley     }
449566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &numPoints, &points));
459566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(closureSection, point, cldof));
46d9917b9dSMatthew G. Knepley   }
479566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(closureSection));
489566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(closureSection, &clSize));
499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(clSize, &clPoints));
50d9917b9dSMatthew G. Knepley   for (point = pStart; point < pEnd; ++point) {
51d9917b9dSMatthew G. Knepley     PetscInt *points = NULL, numPoints, p, q, dof, cldof, cloff;
52d9917b9dSMatthew G. Knepley 
539566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(closureSection, point, &cldof));
549566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(closureSection, point, &cloff));
559566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &numPoints, &points));
56d9917b9dSMatthew G. Knepley     for (p = 0, q = 0; p < numPoints*2; p += 2) {
57d9917b9dSMatthew G. Knepley       if ((points[p] >= sStart) && (points[p] < sEnd)) {
589566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, points[p], &dof));
59d9917b9dSMatthew G. Knepley         if (dof) {
60d9917b9dSMatthew G. Knepley           clPoints[cloff+q*2]   = points[p];
61d9917b9dSMatthew G. Knepley           clPoints[cloff+q*2+1] = points[p+1];
62d9917b9dSMatthew G. Knepley           ++q;
63d9917b9dSMatthew G. Knepley         }
64d9917b9dSMatthew G. Knepley       }
65d9917b9dSMatthew G. Knepley     }
669566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &numPoints, &points));
67*63a3b9bcSJacob Faibussowitsch     PetscCheck(q*2 == cldof,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Invalid size for closure %" PetscInt_FMT " should be %" PetscInt_FMT, q*2, cldof);
68d9917b9dSMatthew G. Knepley   }
699566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, clSize, clPoints, PETSC_OWN_POINTER, &closureIS));
709566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetClosureIndex(section, (PetscObject) dm, closureSection, closureIS));
719566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&closureSection));
729566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&closureIS));
73d9917b9dSMatthew G. Knepley   PetscFunctionReturn(0);
74d9917b9dSMatthew G. Knepley }
75