1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2af0996ceSBarry Smith #include <petsc/private/isimpl.h>
30c312b8eSJed Brown #include <petscsf.h>
4a9fb6443SMatthew G. Knepley #include <petscds.h>
5cb1e1211SMatthew G Knepley
676185916SToby Isaac /* get adjacencies due to point-to-point constraints that can't be found with DMPlexGetAdjacency() */
DMPlexComputeAnchorAdjacencies(DM dm,PetscBool useCone,PetscBool useClosure,PetscSection * anchorSectionAdj,PetscInt * anchorAdj[])7d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexComputeAnchorAdjacencies(DM dm, PetscBool useCone, PetscBool useClosure, PetscSection *anchorSectionAdj, PetscInt *anchorAdj[])
8d71ae5a4SJacob Faibussowitsch {
976185916SToby Isaac PetscInt pStart, pEnd;
10e101f074SMatthew G. Knepley PetscSection section, sectionGlobal, adjSec, aSec;
1176185916SToby Isaac IS aIS;
1276185916SToby Isaac
1376185916SToby Isaac PetscFunctionBegin;
149566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion));
159566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dm, §ionGlobal));
169566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)section), &adjSec));
179566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
189566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(adjSec, pStart, pEnd));
1976185916SToby Isaac
209566063dSJacob Faibussowitsch PetscCall(DMPlexGetAnchors(dm, &aSec, &aIS));
2176185916SToby Isaac if (aSec) {
2276185916SToby Isaac const PetscInt *anchors;
2376185916SToby Isaac PetscInt p, q, a, aSize, *offsets, aStart, aEnd, *inverse, iSize, *adj, adjSize;
2476185916SToby Isaac PetscInt *tmpAdjP = NULL, *tmpAdjQ = NULL;
2576185916SToby Isaac PetscSection inverseSec;
2676185916SToby Isaac
2776185916SToby Isaac /* invert the constraint-to-anchor map */
289566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)aSec), &inverseSec));
299566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(inverseSec, pStart, pEnd));
309566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(aIS, &aSize));
319566063dSJacob Faibussowitsch PetscCall(ISGetIndices(aIS, &anchors));
3276185916SToby Isaac
3376185916SToby Isaac for (p = 0; p < aSize; p++) {
3476185916SToby Isaac PetscInt a = anchors[p];
3576185916SToby Isaac
369566063dSJacob Faibussowitsch PetscCall(PetscSectionAddDof(inverseSec, a, 1));
3776185916SToby Isaac }
389566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(inverseSec));
399566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(inverseSec, &iSize));
409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(iSize, &inverse));
419566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(pEnd - pStart, &offsets));
429566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd));
4376185916SToby Isaac for (p = aStart; p < aEnd; p++) {
4476185916SToby Isaac PetscInt dof, off;
4576185916SToby Isaac
469566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(aSec, p, &dof));
479566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(aSec, p, &off));
4876185916SToby Isaac
4976185916SToby Isaac for (q = 0; q < dof; q++) {
5076185916SToby Isaac PetscInt iOff;
5176185916SToby Isaac
5276185916SToby Isaac a = anchors[off + q];
539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(inverseSec, a, &iOff));
5476185916SToby Isaac inverse[iOff + offsets[a - pStart]++] = p;
5576185916SToby Isaac }
5676185916SToby Isaac }
579566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(aIS, &anchors));
589566063dSJacob Faibussowitsch PetscCall(PetscFree(offsets));
5976185916SToby Isaac
6076185916SToby Isaac /* construct anchorAdj and adjSec
6176185916SToby Isaac *
6276185916SToby Isaac * loop over anchors:
6376185916SToby Isaac * construct anchor adjacency
6476185916SToby Isaac * loop over constrained:
6576185916SToby Isaac * construct constrained adjacency
6676185916SToby Isaac * if not in anchor adjacency, add to dofs
6776185916SToby Isaac * setup adjSec, allocate anchorAdj
6876185916SToby Isaac * loop over anchors:
6976185916SToby Isaac * construct anchor adjacency
7076185916SToby Isaac * loop over constrained:
7176185916SToby Isaac * construct constrained adjacency
7276185916SToby Isaac * if not in anchor adjacency
7376185916SToby Isaac * if not already in list, put in list
7476185916SToby Isaac * sort, unique, reduce dof count
7576185916SToby Isaac * optional: compactify
7676185916SToby Isaac */
7776185916SToby Isaac for (p = pStart; p < pEnd; p++) {
7876185916SToby Isaac PetscInt iDof, iOff, i, r, s, numAdjP = PETSC_DETERMINE;
7976185916SToby Isaac
809566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(inverseSec, p, &iDof));
8176185916SToby Isaac if (!iDof) continue;
829566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(inverseSec, p, &iOff));
839566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, PETSC_TRUE, &numAdjP, &tmpAdjP));
8476185916SToby Isaac for (i = 0; i < iDof; i++) {
8576185916SToby Isaac PetscInt iNew = 0, qAdj, qAdjDof, qAdjCDof, numAdjQ = PETSC_DETERMINE;
8676185916SToby Isaac
8776185916SToby Isaac q = inverse[iOff + i];
889566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, q, useCone, useClosure, PETSC_TRUE, &numAdjQ, &tmpAdjQ));
8976185916SToby Isaac for (r = 0; r < numAdjQ; r++) {
9076185916SToby Isaac qAdj = tmpAdjQ[r];
9176185916SToby Isaac if ((qAdj < pStart) || (qAdj >= pEnd)) continue;
9276185916SToby Isaac for (s = 0; s < numAdjP; s++) {
9376185916SToby Isaac if (qAdj == tmpAdjP[s]) break;
9476185916SToby Isaac }
9576185916SToby Isaac if (s < numAdjP) continue;
969566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, qAdj, &qAdjDof));
979566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, qAdj, &qAdjCDof));
9876185916SToby Isaac iNew += qAdjDof - qAdjCDof;
9976185916SToby Isaac }
1009566063dSJacob Faibussowitsch PetscCall(PetscSectionAddDof(adjSec, p, iNew));
10176185916SToby Isaac }
10276185916SToby Isaac }
10376185916SToby Isaac
1049566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(adjSec));
1059566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(adjSec, &adjSize));
1069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(adjSize, &adj));
10776185916SToby Isaac
10876185916SToby Isaac for (p = pStart; p < pEnd; p++) {
10976185916SToby Isaac PetscInt iDof, iOff, i, r, s, aOff, aOffOrig, aDof, numAdjP = PETSC_DETERMINE;
11076185916SToby Isaac
1119566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(inverseSec, p, &iDof));
11276185916SToby Isaac if (!iDof) continue;
1139566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(inverseSec, p, &iOff));
1149566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, PETSC_TRUE, &numAdjP, &tmpAdjP));
1159566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(adjSec, p, &aDof));
1169566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(adjSec, p, &aOff));
11776185916SToby Isaac aOffOrig = aOff;
11876185916SToby Isaac for (i = 0; i < iDof; i++) {
11976185916SToby Isaac PetscInt qAdj, qAdjDof, qAdjCDof, qAdjOff, nd, numAdjQ = PETSC_DETERMINE;
12076185916SToby Isaac
12176185916SToby Isaac q = inverse[iOff + i];
1229566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, q, useCone, useClosure, PETSC_TRUE, &numAdjQ, &tmpAdjQ));
12376185916SToby Isaac for (r = 0; r < numAdjQ; r++) {
12476185916SToby Isaac qAdj = tmpAdjQ[r];
12576185916SToby Isaac if ((qAdj < pStart) || (qAdj >= pEnd)) continue;
12676185916SToby Isaac for (s = 0; s < numAdjP; s++) {
12776185916SToby Isaac if (qAdj == tmpAdjP[s]) break;
12876185916SToby Isaac }
12976185916SToby Isaac if (s < numAdjP) continue;
1309566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, qAdj, &qAdjDof));
1319566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, qAdj, &qAdjCDof));
1329566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, qAdj, &qAdjOff));
133ad540459SPierre Jolivet for (nd = 0; nd < qAdjDof - qAdjCDof; ++nd) adj[aOff++] = (qAdjOff < 0 ? -(qAdjOff + 1) : qAdjOff) + nd;
13476185916SToby Isaac }
13576185916SToby Isaac }
1368e3a54c0SPierre Jolivet PetscCall(PetscSortRemoveDupsInt(&aDof, PetscSafePointerPlusOffset(adj, aOffOrig)));
1379566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(adjSec, p, aDof));
13876185916SToby Isaac }
13976185916SToby Isaac *anchorAdj = adj;
14076185916SToby Isaac
14176185916SToby Isaac /* clean up */
1429566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&inverseSec));
1439566063dSJacob Faibussowitsch PetscCall(PetscFree(inverse));
1449566063dSJacob Faibussowitsch PetscCall(PetscFree(tmpAdjP));
1459566063dSJacob Faibussowitsch PetscCall(PetscFree(tmpAdjQ));
1469371c9d4SSatish Balay } else {
14776185916SToby Isaac *anchorAdj = NULL;
1489566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(adjSec));
14976185916SToby Isaac }
15076185916SToby Isaac *anchorSectionAdj = adjSec;
1513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
15276185916SToby Isaac }
15376185916SToby Isaac
154f0b52e9eSJames Wright // Based off of `PetscIntViewNumColumns()`
PetscIntViewPairs(PetscInt N,PetscInt Ncol,const PetscInt idx1[],const PetscInt idx2[],PetscViewer viewer)155f0b52e9eSJames Wright static PetscErrorCode PetscIntViewPairs(PetscInt N, PetscInt Ncol, const PetscInt idx1[], const PetscInt idx2[], PetscViewer viewer)
156f0b52e9eSJames Wright {
157f0b52e9eSJames Wright PetscMPIInt rank, size;
158f0b52e9eSJames Wright PetscInt j, i, n = N / Ncol, p = N % Ncol;
159f0b52e9eSJames Wright PetscBool isascii;
160f0b52e9eSJames Wright MPI_Comm comm;
161f0b52e9eSJames Wright
162f0b52e9eSJames Wright PetscFunctionBegin;
163f0b52e9eSJames Wright if (!viewer) viewer = PETSC_VIEWER_STDOUT_SELF;
164f0b52e9eSJames Wright if (N) PetscAssertPointer(idx1, 3);
165f0b52e9eSJames Wright if (N) PetscAssertPointer(idx2, 4);
166f0b52e9eSJames Wright PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5);
167f0b52e9eSJames Wright PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
168f0b52e9eSJames Wright PetscCallMPI(MPI_Comm_size(comm, &size));
169f0b52e9eSJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank));
170f0b52e9eSJames Wright
171f0b52e9eSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
172f0b52e9eSJames Wright if (isascii) {
173f0b52e9eSJames Wright PetscCall(PetscViewerASCIIPushSynchronized(viewer));
174f0b52e9eSJames Wright for (i = 0; i < n; i++) {
175f0b52e9eSJames Wright if (size > 1) {
176f0b52e9eSJames Wright PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT ":", rank, Ncol * i));
177f0b52e9eSJames Wright } else {
178f0b52e9eSJames Wright PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%" PetscInt_FMT ":", Ncol * i));
179f0b52e9eSJames Wright }
180f0b52e9eSJames Wright for (j = 0; j < Ncol; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " (%" PetscInt_FMT ", %" PetscInt_FMT ")", idx1[i * Ncol + j], idx2[i * Ncol + j]));
181f0b52e9eSJames Wright PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
182f0b52e9eSJames Wright }
183f0b52e9eSJames Wright if (p) {
184f0b52e9eSJames Wright if (size > 1) {
185f0b52e9eSJames Wright PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT ":", rank, Ncol * n));
186f0b52e9eSJames Wright } else {
187f0b52e9eSJames Wright PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%" PetscInt_FMT ":", Ncol * n));
188f0b52e9eSJames Wright }
189f0b52e9eSJames Wright for (i = 0; i < p; i++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " (%" PetscInt_FMT ", %" PetscInt_FMT ")", idx1[Ncol * n + i], idx2[Ncol * n + i]));
190f0b52e9eSJames Wright PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
191f0b52e9eSJames Wright }
192f0b52e9eSJames Wright PetscCall(PetscViewerFlush(viewer));
193f0b52e9eSJames Wright PetscCall(PetscViewerASCIIPopSynchronized(viewer));
194f0b52e9eSJames Wright } else {
195f0b52e9eSJames Wright const char *tname;
196f0b52e9eSJames Wright PetscCall(PetscObjectGetName((PetscObject)viewer, &tname));
197f0b52e9eSJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle that PetscViewer of type %s", tname);
198f0b52e9eSJames Wright }
199f0b52e9eSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
200f0b52e9eSJames Wright }
201f0b52e9eSJames Wright
2028c5add6aSPierre Jolivet // Determine if any of the local adjacencies match a leaf and root of the pointSF.
203a38eeca9SJames Wright // When using isoperiodic boundary conditions, it is possible for periodic (leaf) and donor (root) pairs to be on the same rank.
2048c5add6aSPierre Jolivet // This check is done to ensure the adjacency in these cases is only counted for one of the mesh points rather than both.
AdjancencyContainsLeafRootPair(PetscSection myRankPairSection,const PetscInt leaves[],const PetscInt roots[],PetscInt numAdj,const PetscInt tmpAdj[],PetscInt * num_leaves_found,PetscInt leaves_found[])205*9262ad6eSJames Wright static inline PetscErrorCode AdjancencyContainsLeafRootPair(PetscSection myRankPairSection, const PetscInt leaves[], const PetscInt roots[], PetscInt numAdj, const PetscInt tmpAdj[], PetscInt *num_leaves_found, PetscInt leaves_found[])
206a38eeca9SJames Wright {
207*9262ad6eSJames Wright PetscInt root_index = -1, leaf_, num_roots, num_leaves;
208a38eeca9SJames Wright
209a38eeca9SJames Wright PetscFunctionBeginUser;
210*9262ad6eSJames Wright PetscCall(PetscSectionGetChart(myRankPairSection, NULL, &num_roots));
211*9262ad6eSJames Wright PetscCall(PetscSectionGetStorageSize(myRankPairSection, &num_leaves));
212a38eeca9SJames Wright *num_leaves_found = 0;
213a38eeca9SJames Wright for (PetscInt q = 0; q < numAdj; q++) {
214a38eeca9SJames Wright const PetscInt padj = tmpAdj[q];
215*9262ad6eSJames Wright PetscCall(PetscFindInt(padj, num_roots, roots, &root_index));
216a38eeca9SJames Wright if (root_index >= 0) {
217*9262ad6eSJames Wright PetscInt dof, offset;
218a38eeca9SJames Wright
219*9262ad6eSJames Wright PetscCall(PetscSectionGetDof(myRankPairSection, root_index, &dof));
220*9262ad6eSJames Wright PetscCall(PetscSectionGetOffset(myRankPairSection, root_index, &offset));
221*9262ad6eSJames Wright
222*9262ad6eSJames Wright for (PetscInt l = 0; l < dof; l++) {
223*9262ad6eSJames Wright leaf_ = leaves[offset + l];
224a38eeca9SJames Wright for (PetscInt q = 0; q < numAdj; q++) {
225a38eeca9SJames Wright if (tmpAdj[q] == leaf_) {
226a38eeca9SJames Wright leaves_found[*num_leaves_found] = leaf_;
227a38eeca9SJames Wright (*num_leaves_found)++;
228*9262ad6eSJames Wright break;
229*9262ad6eSJames Wright }
230*9262ad6eSJames Wright }
231a38eeca9SJames Wright }
232a38eeca9SJames Wright }
233a38eeca9SJames Wright }
234a38eeca9SJames Wright
235a38eeca9SJames Wright PetscCall(PetscIntSortSemiOrdered(*num_leaves_found, leaves_found));
236a38eeca9SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
237a38eeca9SJames Wright }
238a38eeca9SJames Wright
DMPlexCreateAdjacencySection_Static(DM dm,PetscInt bs,PetscSF sfDof,PetscBool useCone,PetscBool useClosure,PetscBool useAnchors,PetscSection * sA,PetscInt ** colIdx)239d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateAdjacencySection_Static(DM dm, PetscInt bs, PetscSF sfDof, PetscBool useCone, PetscBool useClosure, PetscBool useAnchors, PetscSection *sA, PetscInt **colIdx)
240d71ae5a4SJacob Faibussowitsch {
241cb1e1211SMatthew G Knepley MPI_Comm comm;
242a38eeca9SJames Wright PetscMPIInt myrank;
243a9fb6443SMatthew G. Knepley PetscBool doCommLocal, doComm, debug = PETSC_FALSE;
244a9fb6443SMatthew G. Knepley PetscSF sf, sfAdj;
245*9262ad6eSJames Wright PetscSection section, sectionGlobal, leafSectionAdj, rootSectionAdj, sectionAdj, anchorSectionAdj, myRankPairSection;
246a9fb6443SMatthew G. Knepley PetscInt nroots, nleaves, l, p, r;
247cb1e1211SMatthew G Knepley const PetscInt *leaves;
248cb1e1211SMatthew G Knepley const PetscSFNode *remotes;
249*9262ad6eSJames Wright PetscInt dim, pStart, pEnd, numDof, globalOffStart, globalOffEnd, numCols, adjSize;
250*9262ad6eSJames Wright PetscInt *tmpAdj = NULL, *adj, *rootAdj, *anchorAdj = NULL, *cols, *remoteOffsets, *myRankPairRoots = NULL, *myRankPairLeaves = NULL;
251cb1e1211SMatthew G Knepley
252cb1e1211SMatthew G Knepley PetscFunctionBegin;
2539566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_view_preallocation", &debug, NULL));
255a38eeca9SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &myrank));
2569566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
257a38eeca9SJames Wright PetscCall(DMGetIsoperiodicPointSF_Internal(dm, &sf));
2589566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion));
2599566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dm, §ionGlobal));
2609566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
261a38eeca9SJames Wright doCommLocal = nroots >= 0 ? PETSC_TRUE : PETSC_FALSE;
2625440e5dcSBarry Smith PetscCallMPI(MPIU_Allreduce(&doCommLocal, &doComm, 1, MPI_C_BOOL, MPI_LAND, comm));
263cb1e1211SMatthew G Knepley /* Create section for dof adjacency (dof ==> # adj dof) */
2649566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
2659566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &numDof));
2669566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &leafSectionAdj));
2679566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(leafSectionAdj, 0, numDof));
2689566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &rootSectionAdj));
2699566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(rootSectionAdj, 0, numDof));
2709566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, &leaves, &remotes));
271a38eeca9SJames Wright
272*9262ad6eSJames Wright // Store leaf-root pairs if the leaf and root are both located on current rank.
273*9262ad6eSJames Wright // The point in myRankPairSection is an index into myRankPairRoots.
274*9262ad6eSJames Wright // The section then defines the number of pairs for that root and the offset into myRankPairLeaves to access them.
275*9262ad6eSJames Wright {
276*9262ad6eSJames Wright PetscSegBuffer seg_roots, seg_leaves;
277*9262ad6eSJames Wright PetscCount buffer_size;
278*9262ad6eSJames Wright PetscInt *roots_with_dups, num_non_dups, num_pairs = 0;
279*9262ad6eSJames Wright
280*9262ad6eSJames Wright PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 32, &seg_roots));
281*9262ad6eSJames Wright PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 32, &seg_leaves));
282a38eeca9SJames Wright for (PetscInt l = 0; l < nleaves; l++) {
283a38eeca9SJames Wright if (remotes[l].rank == myrank) {
284*9262ad6eSJames Wright PetscInt *slot;
285*9262ad6eSJames Wright PetscCall(PetscSegBufferGetInts(seg_roots, 1, &slot));
286*9262ad6eSJames Wright *slot = remotes[l].index;
287*9262ad6eSJames Wright PetscCall(PetscSegBufferGetInts(seg_leaves, 1, &slot));
288*9262ad6eSJames Wright *slot = leaves[l];
289a38eeca9SJames Wright }
290a38eeca9SJames Wright }
291*9262ad6eSJames Wright PetscCall(PetscSegBufferGetSize(seg_roots, &buffer_size));
292*9262ad6eSJames Wright PetscCall(PetscIntCast(buffer_size, &num_pairs));
293*9262ad6eSJames Wright PetscCall(PetscSegBufferExtractAlloc(seg_roots, &roots_with_dups));
294*9262ad6eSJames Wright PetscCall(PetscSegBufferExtractAlloc(seg_leaves, &myRankPairLeaves));
295*9262ad6eSJames Wright PetscCall(PetscSegBufferDestroy(&seg_roots));
296*9262ad6eSJames Wright PetscCall(PetscSegBufferDestroy(&seg_leaves));
297*9262ad6eSJames Wright
298*9262ad6eSJames Wright PetscCall(PetscIntSortSemiOrderedWithArray(num_pairs, roots_with_dups, myRankPairLeaves));
299a38eeca9SJames Wright if (debug) {
300f0b52e9eSJames Wright PetscCall(PetscPrintf(comm, "Root/leaf pairs on the same rank:\n"));
301*9262ad6eSJames Wright PetscCall(PetscIntViewPairs(num_pairs, 5, roots_with_dups, myRankPairLeaves, NULL));
302a38eeca9SJames Wright }
303*9262ad6eSJames Wright PetscCall(PetscMalloc1(num_pairs, &myRankPairRoots));
304*9262ad6eSJames Wright PetscCall(PetscArraycpy(myRankPairRoots, roots_with_dups, num_pairs));
305*9262ad6eSJames Wright
306*9262ad6eSJames Wright num_non_dups = num_pairs;
307*9262ad6eSJames Wright PetscCall(PetscSortedRemoveDupsInt(&num_non_dups, myRankPairRoots));
308*9262ad6eSJames Wright PetscCall(PetscSectionCreate(comm, &myRankPairSection));
309*9262ad6eSJames Wright PetscCall(PetscSectionSetChart(myRankPairSection, 0, num_non_dups));
310*9262ad6eSJames Wright for (PetscInt p = 0; p < num_pairs; p++) {
311*9262ad6eSJames Wright PetscInt root = roots_with_dups[p], index;
312*9262ad6eSJames Wright PetscCall(PetscFindInt(root, num_non_dups, myRankPairRoots, &index));
313*9262ad6eSJames Wright PetscAssert(index >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Root not found after removal of duplicates");
314*9262ad6eSJames Wright PetscCall(PetscSectionAddDof(myRankPairSection, index, 1));
315*9262ad6eSJames Wright }
316*9262ad6eSJames Wright PetscCall(PetscSectionSetUp(myRankPairSection));
317*9262ad6eSJames Wright
318*9262ad6eSJames Wright if (debug) {
319*9262ad6eSJames Wright PetscCall(PetscPrintf(comm, "Root/leaf pair section on the same rank:\n"));
320*9262ad6eSJames Wright PetscCall(PetscIntView(num_non_dups, myRankPairRoots, NULL));
321*9262ad6eSJames Wright PetscCall(PetscSectionArrayView(myRankPairSection, myRankPairLeaves, PETSC_INT, NULL));
322*9262ad6eSJames Wright }
323*9262ad6eSJames Wright PetscCall(PetscFree(roots_with_dups));
324*9262ad6eSJames Wright }
325*9262ad6eSJames Wright
326cb1e1211SMatthew G Knepley /*
327964bf7afSMatthew G. Knepley section - maps points to (# dofs, local dofs)
328964bf7afSMatthew G. Knepley sectionGlobal - maps points to (# dofs, global dofs)
329964bf7afSMatthew G. Knepley leafSectionAdj - maps unowned local dofs to # adj dofs
330964bf7afSMatthew G. Knepley rootSectionAdj - maps owned local dofs to # adj dofs
331964bf7afSMatthew G. Knepley adj - adj global dofs indexed by leafSectionAdj
332964bf7afSMatthew G. Knepley rootAdj - adj global dofs indexed by rootSectionAdj
333964bf7afSMatthew G. Knepley sf - describes shared points across procs
334964bf7afSMatthew G. Knepley sfDof - describes shared dofs across procs
335964bf7afSMatthew G. Knepley sfAdj - describes shared adjacent dofs across procs
336cb1e1211SMatthew G Knepley ** The bootstrapping process involves six rounds with similar structure of visiting neighbors of each point.
33776185916SToby Isaac (0). If there are point-to-point constraints, add the adjacencies of constrained points to anchors in anchorAdj
33876185916SToby Isaac (This is done in DMPlexComputeAnchorAdjacencies())
339cb1e1211SMatthew G Knepley 1. Visit unowned points on interface, count adjacencies placing in leafSectionAdj
340cb1e1211SMatthew G Knepley Reduce those counts to rootSectionAdj (now redundantly counting some interface points)
341cb1e1211SMatthew G Knepley 2. Visit owned points on interface, count adjacencies placing in rootSectionAdj
342cb1e1211SMatthew G Knepley Create sfAdj connecting rootSectionAdj and leafSectionAdj
343cb1e1211SMatthew G Knepley 3. Visit unowned points on interface, write adjacencies to adj
344cb1e1211SMatthew G Knepley Gather adj to rootAdj (note that there is redundancy in rootAdj when multiple procs find the same adjacencies)
345cb1e1211SMatthew G Knepley 4. Visit owned points on interface, write adjacencies to rootAdj
346cb1e1211SMatthew G Knepley Remove redundancy in rootAdj
347cb1e1211SMatthew G Knepley ** The last two traversals use transitive closure
348cb1e1211SMatthew G Knepley 5. Visit all owned points in the subdomain, count dofs for each point (sectionAdj)
349cb1e1211SMatthew G Knepley Allocate memory addressed by sectionAdj (cols)
350cb1e1211SMatthew G Knepley 6. Visit all owned points in the subdomain, insert dof adjacencies into cols
351cb1e1211SMatthew G Knepley ** Knowing all the column adjacencies, check ownership and sum into dnz and onz
352cb1e1211SMatthew G Knepley */
3539566063dSJacob Faibussowitsch PetscCall(DMPlexComputeAnchorAdjacencies(dm, useCone, useClosure, &anchorSectionAdj, &anchorAdj));
354cb1e1211SMatthew G Knepley for (l = 0; l < nleaves; ++l) {
35576185916SToby Isaac PetscInt dof, off, d, q, anDof;
35670034214SMatthew G. Knepley PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
357cb1e1211SMatthew G Knepley
358fb47e2feSMatthew G. Knepley if ((p < pStart) || (p >= pEnd)) continue;
3599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof));
3609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off));
3619566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
362cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) {
363cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q];
364cb1e1211SMatthew G Knepley PetscInt ndof, ncdof;
365cb1e1211SMatthew G Knepley
366cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue;
3679566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, padj, &ndof));
3689566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
36948a46eb9SPierre Jolivet for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(leafSectionAdj, d, ndof - ncdof));
370cb1e1211SMatthew G Knepley }
3719566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
37276185916SToby Isaac if (anDof) {
37348a46eb9SPierre Jolivet for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(leafSectionAdj, d, anDof));
37476185916SToby Isaac }
375cb1e1211SMatthew G Knepley }
3769566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(leafSectionAdj));
377cb1e1211SMatthew G Knepley if (debug) {
3789566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n"));
3799566063dSJacob Faibussowitsch PetscCall(PetscSectionView(leafSectionAdj, NULL));
380cb1e1211SMatthew G Knepley }
381cb1e1211SMatthew G Knepley /* Get maximum remote adjacency sizes for owned dofs on interface (roots) */
38247a6104aSMatthew G. Knepley if (doComm) {
3839566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM));
3849566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM));
38569c11d05SVaclav Hapla PetscCall(PetscSectionInvalidateMaxDof_Internal(rootSectionAdj));
386cb1e1211SMatthew G Knepley }
387cb1e1211SMatthew G Knepley if (debug) {
388145b44c9SPierre Jolivet PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Roots:\n"));
3899566063dSJacob Faibussowitsch PetscCall(PetscSectionView(rootSectionAdj, NULL));
390cb1e1211SMatthew G Knepley }
391cb1e1211SMatthew G Knepley /* Add in local adjacency sizes for owned dofs on interface (roots) */
392cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) {
39376185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof;
394cb1e1211SMatthew G Knepley
3959566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof));
3969566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off));
397cb1e1211SMatthew G Knepley if (!dof) continue;
3989566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
399cb1e1211SMatthew G Knepley if (adof <= 0) continue;
4009566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
401cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) {
402cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q];
403cb1e1211SMatthew G Knepley PetscInt ndof, ncdof;
404cb1e1211SMatthew G Knepley
405cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue;
4069566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, padj, &ndof));
4079566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
40848a46eb9SPierre Jolivet for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(rootSectionAdj, d, ndof - ncdof));
409cb1e1211SMatthew G Knepley }
4109566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
41176185916SToby Isaac if (anDof) {
41248a46eb9SPierre Jolivet for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(rootSectionAdj, d, anDof));
41376185916SToby Isaac }
414cb1e1211SMatthew G Knepley }
4159566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(rootSectionAdj));
416cb1e1211SMatthew G Knepley if (debug) {
417145b44c9SPierre Jolivet PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Roots after local additions:\n"));
4189566063dSJacob Faibussowitsch PetscCall(PetscSectionView(rootSectionAdj, NULL));
419cb1e1211SMatthew G Knepley }
420cb1e1211SMatthew G Knepley /* Create adj SF based on dof SF */
4219566063dSJacob Faibussowitsch PetscCall(PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets));
4229566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj));
4239566063dSJacob Faibussowitsch PetscCall(PetscFree(remoteOffsets));
4246a1a2f7bSJames Wright if (debug) {
4259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Adjacency SF for Preallocation:\n"));
4269566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfAdj, NULL));
427cb1e1211SMatthew G Knepley }
428cb1e1211SMatthew G Knepley /* Create leaf adjacency */
4299566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(leafSectionAdj));
4309566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(leafSectionAdj, &adjSize));
4319566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(adjSize, &adj));
432cb1e1211SMatthew G Knepley for (l = 0; l < nleaves; ++l) {
43376185916SToby Isaac PetscInt dof, off, d, q, anDof, anOff;
43470034214SMatthew G. Knepley PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
435cb1e1211SMatthew G Knepley
436fb47e2feSMatthew G. Knepley if ((p < pStart) || (p >= pEnd)) continue;
4379566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof));
4389566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off));
4399566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
4409566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
4419566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
442cb1e1211SMatthew G Knepley for (d = off; d < off + dof; ++d) {
443cb1e1211SMatthew G Knepley PetscInt aoff, i = 0;
444cb1e1211SMatthew G Knepley
4459566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(leafSectionAdj, d, &aoff));
446cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) {
447cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q];
448cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, ngoff, nd;
449cb1e1211SMatthew G Knepley
450cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue;
4519566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, padj, &ndof));
4529566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
4539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
454cb1e1211SMatthew G Knepley for (nd = 0; nd < ndof - ncdof; ++nd) {
455cb1e1211SMatthew G Knepley adj[aoff + i] = (ngoff < 0 ? -(ngoff + 1) : ngoff) + nd;
456cb1e1211SMatthew G Knepley ++i;
457cb1e1211SMatthew G Knepley }
458cb1e1211SMatthew G Knepley }
45976185916SToby Isaac for (q = 0; q < anDof; q++) {
46076185916SToby Isaac adj[aoff + i] = anchorAdj[anOff + q];
46176185916SToby Isaac ++i;
46276185916SToby Isaac }
463cb1e1211SMatthew G Knepley }
464cb1e1211SMatthew G Knepley }
465cb1e1211SMatthew G Knepley /* Debugging */
466cb1e1211SMatthew G Knepley if (debug) {
4679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Leaf adjacency indices\n"));
4686a1a2f7bSJames Wright PetscCall(PetscSectionArrayView(leafSectionAdj, adj, PETSC_INT, NULL));
469cb1e1211SMatthew G Knepley }
470543482b8SMatthew G. Knepley /* Gather adjacent indices to root */
4719566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(rootSectionAdj, &adjSize));
4729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(adjSize, &rootAdj));
473cb1e1211SMatthew G Knepley for (r = 0; r < adjSize; ++r) rootAdj[r] = -1;
47447a6104aSMatthew G. Knepley if (doComm) {
475543482b8SMatthew G. Knepley const PetscInt *indegree;
476543482b8SMatthew G. Knepley PetscInt *remoteadj, radjsize = 0;
477543482b8SMatthew G. Knepley
4789566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sfAdj, &indegree));
4799566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sfAdj, &indegree));
480543482b8SMatthew G. Knepley for (p = 0; p < adjSize; ++p) radjsize += indegree[p];
4819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(radjsize, &remoteadj));
4829566063dSJacob Faibussowitsch PetscCall(PetscSFGatherBegin(sfAdj, MPIU_INT, adj, remoteadj));
4839566063dSJacob Faibussowitsch PetscCall(PetscSFGatherEnd(sfAdj, MPIU_INT, adj, remoteadj));
4846dba2905SMatthew G. Knepley for (p = 0, l = 0, r = 0; p < adjSize; ++p, l = PetscMax(p, l + indegree[p - 1])) {
485543482b8SMatthew G. Knepley PetscInt s;
4866dba2905SMatthew G. Knepley for (s = 0; s < indegree[p]; ++s, ++r) rootAdj[l + s] = remoteadj[r];
487543482b8SMatthew G. Knepley }
48863a3b9bcSJacob Faibussowitsch PetscCheck(r == radjsize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %" PetscInt_FMT " != %" PetscInt_FMT, r, radjsize);
48963a3b9bcSJacob Faibussowitsch PetscCheck(l == adjSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %" PetscInt_FMT " != %" PetscInt_FMT, l, adjSize);
4909566063dSJacob Faibussowitsch PetscCall(PetscFree(remoteadj));
491cb1e1211SMatthew G Knepley }
4929566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfAdj));
4939566063dSJacob Faibussowitsch PetscCall(PetscFree(adj));
494cb1e1211SMatthew G Knepley /* Debugging */
495cb1e1211SMatthew G Knepley if (debug) {
4969566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Root adjacency indices after gather\n"));
4976a1a2f7bSJames Wright PetscCall(PetscSectionArrayView(rootSectionAdj, rootAdj, PETSC_INT, NULL));
498cb1e1211SMatthew G Knepley }
499cb1e1211SMatthew G Knepley /* Add in local adjacency indices for owned dofs on interface (roots) */
500cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) {
50176185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof, anOff;
502cb1e1211SMatthew G Knepley
5039566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof));
5049566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off));
505cb1e1211SMatthew G Knepley if (!dof) continue;
5069566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
507cb1e1211SMatthew G Knepley if (adof <= 0) continue;
5089566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
5099566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
5109566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
511cb1e1211SMatthew G Knepley for (d = off; d < off + dof; ++d) {
512cb1e1211SMatthew G Knepley PetscInt adof, aoff, i;
513cb1e1211SMatthew G Knepley
5149566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, d, &adof));
5159566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(rootSectionAdj, d, &aoff));
516cb1e1211SMatthew G Knepley i = adof - 1;
51776185916SToby Isaac for (q = 0; q < anDof; q++) {
51876185916SToby Isaac rootAdj[aoff + i] = anchorAdj[anOff + q];
51976185916SToby Isaac --i;
52076185916SToby Isaac }
521cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) {
522cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q];
523cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, ngoff, nd;
524cb1e1211SMatthew G Knepley
525cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue;
5269566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, padj, &ndof));
5279566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
5289566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
529cb1e1211SMatthew G Knepley for (nd = 0; nd < ndof - ncdof; ++nd) {
530cb1e1211SMatthew G Knepley rootAdj[aoff + i] = ngoff < 0 ? -(ngoff + 1) + nd : ngoff + nd;
531cb1e1211SMatthew G Knepley --i;
532cb1e1211SMatthew G Knepley }
533cb1e1211SMatthew G Knepley }
534cb1e1211SMatthew G Knepley }
535cb1e1211SMatthew G Knepley }
536cb1e1211SMatthew G Knepley /* Debugging */
537cb1e1211SMatthew G Knepley if (debug) {
5386a1a2f7bSJames Wright PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Roots before compression:\n"));
5396a1a2f7bSJames Wright PetscCall(PetscSectionArrayView(rootSectionAdj, rootAdj, PETSC_INT, NULL));
540cb1e1211SMatthew G Knepley }
541cb1e1211SMatthew G Knepley /* Compress indices */
5429566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(rootSectionAdj));
543cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) {
544cb1e1211SMatthew G Knepley PetscInt dof, cdof, off, d;
545cb1e1211SMatthew G Knepley PetscInt adof, aoff;
546cb1e1211SMatthew G Knepley
5479566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof));
5489566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
5499566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off));
550cb1e1211SMatthew G Knepley if (!dof) continue;
5519566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
552cb1e1211SMatthew G Knepley if (adof <= 0) continue;
553cb1e1211SMatthew G Knepley for (d = off; d < off + dof - cdof; ++d) {
5549566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, d, &adof));
5559566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(rootSectionAdj, d, &aoff));
5569566063dSJacob Faibussowitsch PetscCall(PetscSortRemoveDupsInt(&adof, &rootAdj[aoff]));
5579566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(rootSectionAdj, d, adof));
558cb1e1211SMatthew G Knepley }
559cb1e1211SMatthew G Knepley }
560cb1e1211SMatthew G Knepley /* Debugging */
561cb1e1211SMatthew G Knepley if (debug) {
562145b44c9SPierre Jolivet PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Roots after compression:\n"));
5636a1a2f7bSJames Wright PetscCall(PetscSectionArrayView(rootSectionAdj, rootAdj, PETSC_INT, NULL));
564cb1e1211SMatthew G Knepley }
565cb1e1211SMatthew G Knepley /* Build adjacency section: Maps global indices to sets of adjacent global indices */
5669566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd));
5679566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, §ionAdj));
5689566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(sectionAdj, globalOffStart, globalOffEnd));
569a38eeca9SJames Wright
570*9262ad6eSJames Wright PetscInt *exclude_leaves, num_exclude_leaves = 0, max_adjacency_size = 0;
571a38eeca9SJames Wright PetscCall(DMPlexGetMaxAdjacencySize_Internal(dm, useAnchors, &max_adjacency_size));
572a38eeca9SJames Wright PetscCall(PetscMalloc1(max_adjacency_size, &exclude_leaves));
573a38eeca9SJames Wright
574cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) {
57576185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof;
576cb1e1211SMatthew G Knepley PetscBool found = PETSC_TRUE;
577cb1e1211SMatthew G Knepley
5789566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof));
5799566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
5809566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off));
5819566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
582cb1e1211SMatthew G Knepley for (d = 0; d < dof - cdof; ++d) {
583cb1e1211SMatthew G Knepley PetscInt ldof, rdof;
584cb1e1211SMatthew G Knepley
5859566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSectionAdj, off + d, &ldof));
5869566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, off + d, &rdof));
587cb1e1211SMatthew G Knepley if (ldof > 0) {
588cb1e1211SMatthew G Knepley /* We do not own this point */
589cb1e1211SMatthew G Knepley } else if (rdof > 0) {
5909566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(sectionAdj, goff + d, rdof));
591cb1e1211SMatthew G Knepley } else {
592cb1e1211SMatthew G Knepley found = PETSC_FALSE;
593cb1e1211SMatthew G Knepley }
594cb1e1211SMatthew G Knepley }
595cb1e1211SMatthew G Knepley if (found) continue;
5969566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof));
5979566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
5989566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
599*9262ad6eSJames Wright PetscCall(AdjancencyContainsLeafRootPair(myRankPairSection, myRankPairLeaves, myRankPairRoots, numAdj, tmpAdj, &num_exclude_leaves, exclude_leaves));
600cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) {
601cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q];
602a38eeca9SJames Wright PetscInt ndof, ncdof, noff, count;
603cb1e1211SMatthew G Knepley
604cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue;
6059566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, padj, &ndof));
6069566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
6079566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, padj, &noff));
608a38eeca9SJames Wright // If leaf-root pair are both on this rank, only count root
609a38eeca9SJames Wright PetscCall(PetscFindInt(padj, num_exclude_leaves, exclude_leaves, &count));
610a38eeca9SJames Wright if (count >= 0) continue;
61148a46eb9SPierre Jolivet for (d = goff; d < goff + dof - cdof; ++d) PetscCall(PetscSectionAddDof(sectionAdj, d, ndof - ncdof));
612cb1e1211SMatthew G Knepley }
6139566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
61476185916SToby Isaac if (anDof) {
61548a46eb9SPierre Jolivet for (d = goff; d < goff + dof - cdof; ++d) PetscCall(PetscSectionAddDof(sectionAdj, d, anDof));
61676185916SToby Isaac }
617cb1e1211SMatthew G Knepley }
6189566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(sectionAdj));
619cb1e1211SMatthew G Knepley if (debug) {
6209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation:\n"));
6219566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sectionAdj, NULL));
622cb1e1211SMatthew G Knepley }
623cb1e1211SMatthew G Knepley /* Get adjacent indices */
6249566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(sectionAdj, &numCols));
6259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCols, &cols));
626cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) {
62776185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof, anOff;
628cb1e1211SMatthew G Knepley PetscBool found = PETSC_TRUE;
629cb1e1211SMatthew G Knepley
6309566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof));
6319566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
6329566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off));
6339566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
634cb1e1211SMatthew G Knepley for (d = 0; d < dof - cdof; ++d) {
635cb1e1211SMatthew G Knepley PetscInt ldof, rdof;
636cb1e1211SMatthew G Knepley
6379566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSectionAdj, off + d, &ldof));
6389566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, off + d, &rdof));
639cb1e1211SMatthew G Knepley if (ldof > 0) {
640cb1e1211SMatthew G Knepley /* We do not own this point */
641cb1e1211SMatthew G Knepley } else if (rdof > 0) {
642cb1e1211SMatthew G Knepley PetscInt aoff, roff;
643cb1e1211SMatthew G Knepley
6449566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionAdj, goff + d, &aoff));
6459566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(rootSectionAdj, off + d, &roff));
6469566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&cols[aoff], &rootAdj[roff], rdof));
647cb1e1211SMatthew G Knepley } else {
648cb1e1211SMatthew G Knepley found = PETSC_FALSE;
649cb1e1211SMatthew G Knepley }
650cb1e1211SMatthew G Knepley }
651cb1e1211SMatthew G Knepley if (found) continue;
6529566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
653*9262ad6eSJames Wright PetscCall(AdjancencyContainsLeafRootPair(myRankPairSection, myRankPairLeaves, myRankPairRoots, numAdj, tmpAdj, &num_exclude_leaves, exclude_leaves));
6549566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
6559566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
656cb1e1211SMatthew G Knepley for (d = goff; d < goff + dof - cdof; ++d) {
657cb1e1211SMatthew G Knepley PetscInt adof, aoff, i = 0;
658cb1e1211SMatthew G Knepley
6599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionAdj, d, &adof));
6609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionAdj, d, &aoff));
661cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) {
662a38eeca9SJames Wright const PetscInt padj = tmpAdj[q], *ncind;
663a38eeca9SJames Wright PetscInt ndof, ncdof, ngoff, nd, count;
664cb1e1211SMatthew G Knepley
665cb1e1211SMatthew G Knepley /* Adjacent points may not be in the section chart */
666cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue;
6679566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, padj, &ndof));
6689566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
6699566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintIndices(section, padj, &ncind));
6709566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
671a38eeca9SJames Wright // If leaf-root pair are both on this rank, only count root
672a38eeca9SJames Wright PetscCall(PetscFindInt(padj, num_exclude_leaves, exclude_leaves, &count));
673a38eeca9SJames Wright if (count >= 0) continue;
674ad540459SPierre Jolivet for (nd = 0; nd < ndof - ncdof; ++nd, ++i) cols[aoff + i] = ngoff < 0 ? -(ngoff + 1) + nd : ngoff + nd;
675cb1e1211SMatthew G Knepley }
676ad540459SPierre Jolivet for (q = 0; q < anDof; q++, i++) cols[aoff + i] = anchorAdj[anOff + q];
67763a3b9bcSJacob Faibussowitsch PetscCheck(i == adof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of entries %" PetscInt_FMT " != %" PetscInt_FMT " for dof %" PetscInt_FMT " (point %" PetscInt_FMT ")", i, adof, d, p);
678cb1e1211SMatthew G Knepley }
679cb1e1211SMatthew G Knepley }
6809566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&anchorSectionAdj));
6819566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&leafSectionAdj));
6829566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&rootSectionAdj));
683*9262ad6eSJames Wright PetscCall(PetscSectionDestroy(&myRankPairSection));
684a38eeca9SJames Wright PetscCall(PetscFree(exclude_leaves));
685*9262ad6eSJames Wright PetscCall(PetscFree(myRankPairLeaves));
686*9262ad6eSJames Wright PetscCall(PetscFree(myRankPairRoots));
6879566063dSJacob Faibussowitsch PetscCall(PetscFree(anchorAdj));
6889566063dSJacob Faibussowitsch PetscCall(PetscFree(rootAdj));
6899566063dSJacob Faibussowitsch PetscCall(PetscFree(tmpAdj));
690cb1e1211SMatthew G Knepley /* Debugging */
691cb1e1211SMatthew G Knepley if (debug) {
6929566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Column indices\n"));
6936a1a2f7bSJames Wright PetscCall(PetscSectionArrayView(sectionAdj, cols, PETSC_INT, NULL));
694cb1e1211SMatthew G Knepley }
695a9fb6443SMatthew G. Knepley
696a9fb6443SMatthew G. Knepley *sA = sectionAdj;
697a9fb6443SMatthew G. Knepley *colIdx = cols;
6983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
699a9fb6443SMatthew G. Knepley }
700a9fb6443SMatthew G. Knepley
DMPlexUpdateAllocation_Static(DM dm,PetscLayout rLayout,PetscInt bs,PetscInt f,PetscSection sectionAdj,const PetscInt cols[],PetscInt dnz[],PetscInt onz[],PetscInt dnzu[],PetscInt onzu[])701d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexUpdateAllocation_Static(DM dm, PetscLayout rLayout, PetscInt bs, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[])
702d71ae5a4SJacob Faibussowitsch {
703e101f074SMatthew G. Knepley PetscSection section;
704a9fb6443SMatthew G. Knepley PetscInt rStart, rEnd, r, pStart, pEnd, p;
705a9fb6443SMatthew G. Knepley
706a9fb6443SMatthew G. Knepley PetscFunctionBegin;
707a9fb6443SMatthew G. Knepley /* This loop needs to change to a loop over points, then field dofs, which means we need to look both sections */
7089566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(rLayout, &rStart, &rEnd));
7091dca8a05SBarry Smith PetscCheck((rStart % bs) == 0 && (rEnd % bs) == 0, PetscObjectComm((PetscObject)rLayout), PETSC_ERR_ARG_WRONG, "Invalid layout [%" PetscInt_FMT ", %" PetscInt_FMT ") for matrix, must be divisible by block size %" PetscInt_FMT, rStart, rEnd, bs);
710a9fb6443SMatthew G. Knepley if (f >= 0 && bs == 1) {
7119566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion));
7129566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
713a9fb6443SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) {
714294b7009SMatthew G. Knepley PetscInt rS, rE;
715a9fb6443SMatthew G. Knepley
7169566063dSJacob Faibussowitsch PetscCall(DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE));
717a9fb6443SMatthew G. Knepley for (r = rS; r < rE; ++r) {
718a9fb6443SMatthew G. Knepley PetscInt numCols, cStart, c;
719a9fb6443SMatthew G. Knepley
7209566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
7219566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
722a9fb6443SMatthew G. Knepley for (c = cStart; c < cStart + numCols; ++c) {
723a9fb6443SMatthew G. Knepley if ((cols[c] >= rStart) && (cols[c] < rEnd)) {
724a9fb6443SMatthew G. Knepley ++dnz[r - rStart];
725a9fb6443SMatthew G. Knepley if (cols[c] >= r) ++dnzu[r - rStart];
726a9fb6443SMatthew G. Knepley } else {
727a9fb6443SMatthew G. Knepley ++onz[r - rStart];
728a9fb6443SMatthew G. Knepley if (cols[c] >= r) ++onzu[r - rStart];
729a9fb6443SMatthew G. Knepley }
730a9fb6443SMatthew G. Knepley }
731a9fb6443SMatthew G. Knepley }
732a9fb6443SMatthew G. Knepley }
733a9fb6443SMatthew G. Knepley } else {
734cb1e1211SMatthew G Knepley /* Only loop over blocks of rows */
735cb1e1211SMatthew G Knepley for (r = rStart / bs; r < rEnd / bs; ++r) {
736cb1e1211SMatthew G Knepley const PetscInt row = r * bs;
737cb1e1211SMatthew G Knepley PetscInt numCols, cStart, c;
738cb1e1211SMatthew G Knepley
7399566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionAdj, row, &numCols));
7409566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionAdj, row, &cStart));
741cb1e1211SMatthew G Knepley for (c = cStart; c < cStart + numCols; ++c) {
742e7bcfa36SMatthew G. Knepley if ((cols[c] >= rStart) && (cols[c] < rEnd)) {
743e7bcfa36SMatthew G. Knepley ++dnz[r - rStart / bs];
744e7bcfa36SMatthew G. Knepley if (cols[c] >= row) ++dnzu[r - rStart / bs];
745cb1e1211SMatthew G Knepley } else {
746e7bcfa36SMatthew G. Knepley ++onz[r - rStart / bs];
747e7bcfa36SMatthew G. Knepley if (cols[c] >= row) ++onzu[r - rStart / bs];
748cb1e1211SMatthew G Knepley }
749cb1e1211SMatthew G Knepley }
750cb1e1211SMatthew G Knepley }
751a9fb6443SMatthew G. Knepley for (r = 0; r < (rEnd - rStart) / bs; ++r) {
752cb1e1211SMatthew G Knepley dnz[r] /= bs;
753cb1e1211SMatthew G Knepley onz[r] /= bs;
754cb1e1211SMatthew G Knepley dnzu[r] /= bs;
755cb1e1211SMatthew G Knepley onzu[r] /= bs;
756cb1e1211SMatthew G Knepley }
757cb1e1211SMatthew G Knepley }
7583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
759a9fb6443SMatthew G. Knepley }
760a9fb6443SMatthew G. Knepley
DMPlexFillMatrix_Static(DM dm,PetscLayout rLayout,PetscInt bs,PetscInt f,PetscSection sectionAdj,const PetscInt cols[],Mat A)761d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexFillMatrix_Static(DM dm, PetscLayout rLayout, PetscInt bs, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], Mat A)
762d71ae5a4SJacob Faibussowitsch {
763e101f074SMatthew G. Knepley PetscSection section;
764a9fb6443SMatthew G. Knepley PetscScalar *values;
765a9fb6443SMatthew G. Knepley PetscInt rStart, rEnd, r, pStart, pEnd, p, len, maxRowLen = 0;
766a9fb6443SMatthew G. Knepley
767a9fb6443SMatthew G. Knepley PetscFunctionBegin;
7689566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(rLayout, &rStart, &rEnd));
769a9fb6443SMatthew G. Knepley for (r = rStart; r < rEnd; ++r) {
7709566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionAdj, r, &len));
771a9fb6443SMatthew G. Knepley maxRowLen = PetscMax(maxRowLen, len);
772a9fb6443SMatthew G. Knepley }
7739566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(maxRowLen, &values));
774a9fb6443SMatthew G. Knepley if (f >= 0 && bs == 1) {
7759566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion));
7769566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
777a9fb6443SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) {
778294b7009SMatthew G. Knepley PetscInt rS, rE;
779a9fb6443SMatthew G. Knepley
7809566063dSJacob Faibussowitsch PetscCall(DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE));
781a9fb6443SMatthew G. Knepley for (r = rS; r < rE; ++r) {
7827e01ee02SMatthew G. Knepley PetscInt numCols, cStart;
783a9fb6443SMatthew G. Knepley
7849566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
7859566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
7869566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES));
787a9fb6443SMatthew G. Knepley }
788a9fb6443SMatthew G. Knepley }
789a9fb6443SMatthew G. Knepley } else {
790a9fb6443SMatthew G. Knepley for (r = rStart; r < rEnd; ++r) {
791a9fb6443SMatthew G. Knepley PetscInt numCols, cStart;
792a9fb6443SMatthew G. Knepley
7939566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
7949566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
7959566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES));
796a9fb6443SMatthew G. Knepley }
797a9fb6443SMatthew G. Knepley }
7989566063dSJacob Faibussowitsch PetscCall(PetscFree(values));
7993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
800a9fb6443SMatthew G. Knepley }
801a9fb6443SMatthew G. Knepley
802cc4c1da9SBarry Smith /*@
803a1cb98faSBarry Smith DMPlexPreallocateOperator - Calculate the matrix nonzero pattern based upon the information in the `DM`,
804a1cb98faSBarry Smith the `PetscDS` it contains, and the default `PetscSection`.
80525e402d2SMatthew G. Knepley
80625e402d2SMatthew G. Knepley Collective
80725e402d2SMatthew G. Knepley
8084165533cSJose E. Roman Input Parameters:
809a1cb98faSBarry Smith + dm - The `DMPLEX`
81025e402d2SMatthew G. Knepley . bs - The matrix blocksize
81125e402d2SMatthew G. Knepley . dnz - An array to hold the number of nonzeros in the diagonal block
81225e402d2SMatthew G. Knepley . onz - An array to hold the number of nonzeros in the off-diagonal block
81325e402d2SMatthew G. Knepley . dnzu - An array to hold the number of nonzeros in the upper triangle of the diagonal block
81425e402d2SMatthew G. Knepley . onzu - An array to hold the number of nonzeros in the upper triangle of the off-diagonal block
815a1cb98faSBarry Smith - fillMatrix - If `PETSC_TRUE`, fill the matrix with zeros
81625e402d2SMatthew G. Knepley
8174165533cSJose E. Roman Output Parameter:
81825e402d2SMatthew G. Knepley . A - The preallocated matrix
81925e402d2SMatthew G. Knepley
82025e402d2SMatthew G. Knepley Level: advanced
82125e402d2SMatthew G. Knepley
8221cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMCreateMatrix()`
82325e402d2SMatthew G. Knepley @*/
DMPlexPreallocateOperator(DM dm,PetscInt bs,PetscInt dnz[],PetscInt onz[],PetscInt dnzu[],PetscInt onzu[],Mat A,PetscBool fillMatrix)824d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexPreallocateOperator(DM dm, PetscInt bs, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix)
825d71ae5a4SJacob Faibussowitsch {
826a9fb6443SMatthew G. Knepley MPI_Comm comm;
827a9fb6443SMatthew G. Knepley MatType mtype;
828a9fb6443SMatthew G. Knepley PetscSF sf, sfDof;
829e101f074SMatthew G. Knepley PetscSection section;
830a9fb6443SMatthew G. Knepley PetscInt *remoteOffsets;
831a9fb6443SMatthew G. Knepley PetscSection sectionAdj[4] = {NULL, NULL, NULL, NULL};
832a9fb6443SMatthew G. Knepley PetscInt *cols[4] = {NULL, NULL, NULL, NULL};
833a9fb6443SMatthew G. Knepley PetscBool useCone, useClosure;
8347e01ee02SMatthew G. Knepley PetscInt Nf, f, idx, locRows;
835a9fb6443SMatthew G. Knepley PetscLayout rLayout;
836a9fb6443SMatthew G. Knepley PetscBool isSymBlock, isSymSeqBlock, isSymMPIBlock, debug = PETSC_FALSE;
837a9fb6443SMatthew G. Knepley
838a9fb6443SMatthew G. Knepley PetscFunctionBegin;
839a9fb6443SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
840064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(A, MAT_CLASSID, 7);
8414f572ea9SToby Isaac if (dnz) PetscAssertPointer(dnz, 3);
8424f572ea9SToby Isaac if (onz) PetscAssertPointer(onz, 4);
8434f572ea9SToby Isaac if (dnzu) PetscAssertPointer(dnzu, 5);
8444f572ea9SToby Isaac if (onzu) PetscAssertPointer(onzu, 6);
845a38eeca9SJames Wright PetscCall(DMGetIsoperiodicPointSF_Internal(dm, &sf));
8469566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion));
8479566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_view_preallocation", &debug, NULL));
8489566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
8499566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_Preallocate, dm, 0, 0, 0));
850a9fb6443SMatthew G. Knepley /* Create dof SF based on point SF */
851a9fb6443SMatthew G. Knepley if (debug) {
852e101f074SMatthew G. Knepley PetscSection section, sectionGlobal;
853a9fb6443SMatthew G. Knepley PetscSF sf;
854a9fb6443SMatthew G. Knepley
855a38eeca9SJames Wright PetscCall(DMGetIsoperiodicPointSF_Internal(dm, &sf));
8569566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion));
8579566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dm, §ionGlobal));
8589566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Input Section for Preallocation:\n"));
8599566063dSJacob Faibussowitsch PetscCall(PetscSectionView(section, NULL));
8609566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Input Global Section for Preallocation:\n"));
8619566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sectionGlobal, NULL));
8629566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Input SF for Preallocation:\n"));
8639566063dSJacob Faibussowitsch PetscCall(PetscSFView(sf, NULL));
864a9fb6443SMatthew G. Knepley }
8659566063dSJacob Faibussowitsch PetscCall(PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets));
8669566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof));
8679566063dSJacob Faibussowitsch PetscCall(PetscFree(remoteOffsets));
8686a1a2f7bSJames Wright if (debug) {
8699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Dof SF for Preallocation:\n"));
8709566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfDof, NULL));
871a9fb6443SMatthew G. Knepley }
872a9fb6443SMatthew G. Knepley /* Create allocation vectors from adjacency graph */
8739566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &locRows, NULL));
8749566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &rLayout));
8759566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
8769566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
8779566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(rLayout));
878a9fb6443SMatthew G. Knepley /* There are 4 types of adjacency */
8799566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &Nf));
880a9fb6443SMatthew G. Knepley if (Nf < 1 || bs > 1) {
8819566063dSJacob Faibussowitsch PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure));
882a9fb6443SMatthew G. Knepley idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
8839566063dSJacob Faibussowitsch PetscCall(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, §ionAdj[idx], &cols[idx]));
8849566063dSJacob Faibussowitsch PetscCall(DMPlexUpdateAllocation_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu));
885a9fb6443SMatthew G. Knepley } else {
886a9fb6443SMatthew G. Knepley for (f = 0; f < Nf; ++f) {
8879566063dSJacob Faibussowitsch PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure));
888a9fb6443SMatthew G. Knepley idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
8899566063dSJacob Faibussowitsch if (!sectionAdj[idx]) PetscCall(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, §ionAdj[idx], &cols[idx]));
8909566063dSJacob Faibussowitsch PetscCall(DMPlexUpdateAllocation_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu));
891a9fb6443SMatthew G. Knepley }
892a9fb6443SMatthew G. Knepley }
8939566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfDof));
894cb1e1211SMatthew G Knepley /* Set matrix pattern */
8959566063dSJacob Faibussowitsch PetscCall(MatXAIJSetPreallocation(A, bs, dnz, onz, dnzu, onzu));
8969566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
89789545effSMatthew G. Knepley /* Check for symmetric storage */
8989566063dSJacob Faibussowitsch PetscCall(MatGetType(A, &mtype));
8999566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(mtype, MATSBAIJ, &isSymBlock));
9009566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock));
9019566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock));
9029566063dSJacob Faibussowitsch if (isSymBlock || isSymSeqBlock || isSymMPIBlock) PetscCall(MatSetOption(A, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
903cb1e1211SMatthew G Knepley /* Fill matrix with zeros */
904cb1e1211SMatthew G Knepley if (fillMatrix) {
905a9fb6443SMatthew G. Knepley if (Nf < 1 || bs > 1) {
9069566063dSJacob Faibussowitsch PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure));
907a9fb6443SMatthew G. Knepley idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
9089566063dSJacob Faibussowitsch PetscCall(DMPlexFillMatrix_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], A));
909a9fb6443SMatthew G. Knepley } else {
910a9fb6443SMatthew G. Knepley for (f = 0; f < Nf; ++f) {
9119566063dSJacob Faibussowitsch PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure));
912a9fb6443SMatthew G. Knepley idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
9139566063dSJacob Faibussowitsch PetscCall(DMPlexFillMatrix_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], A));
914cb1e1211SMatthew G Knepley }
915cb1e1211SMatthew G Knepley }
9169566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
9179566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
918cb1e1211SMatthew G Knepley }
9199566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&rLayout));
9209371c9d4SSatish Balay for (idx = 0; idx < 4; ++idx) {
9219371c9d4SSatish Balay PetscCall(PetscSectionDestroy(§ionAdj[idx]));
9229371c9d4SSatish Balay PetscCall(PetscFree(cols[idx]));
9239371c9d4SSatish Balay }
9249566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_Preallocate, dm, 0, 0, 0));
9253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
926cb1e1211SMatthew G Knepley }
927cb1e1211SMatthew G Knepley
928cb1e1211SMatthew G Knepley #if 0
929cb1e1211SMatthew G Knepley PetscErrorCode DMPlexPreallocateOperator_2(DM dm, PetscInt bs, PetscSection section, PetscSection sectionGlobal, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix)
930cb1e1211SMatthew G Knepley {
931cb1e1211SMatthew G Knepley PetscInt *tmpClosure,*tmpAdj,*visits;
932cb1e1211SMatthew G Knepley PetscInt c,cStart,cEnd,pStart,pEnd;
933cb1e1211SMatthew G Knepley
934cb1e1211SMatthew G Knepley PetscFunctionBegin;
9359566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
9369566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth));
9379566063dSJacob Faibussowitsch PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize));
938cb1e1211SMatthew G Knepley
939e7b80042SMatthew G. Knepley maxClosureSize = 2*PetscMax(PetscPowInt(mesh->maxConeSize,depth+1),PetscPowInt(mesh->maxSupportSize,depth+1));
940cb1e1211SMatthew G Knepley
9419566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
942cb1e1211SMatthew G Knepley npoints = pEnd - pStart;
943cb1e1211SMatthew G Knepley
9449566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(maxClosureSize,&tmpClosure,npoints,&lvisits,npoints,&visits));
9459566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(lvisits,pEnd-pStart));
9469566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(visits,pEnd-pStart));
9479566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
948cb1e1211SMatthew G Knepley for (c=cStart; c<cEnd; c++) {
949cb1e1211SMatthew G Knepley PetscInt *support = tmpClosure;
9509566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_FALSE, &supportSize, (PetscInt**)&support));
951cb1e1211SMatthew G Knepley for (p=0; p<supportSize; p++) lvisits[support[p]]++;
952cb1e1211SMatthew G Knepley }
9539566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf,MPIU_INT,lvisits,visits,MPI_SUM));
9549566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd (sf,MPIU_INT,lvisits,visits,MPI_SUM));
9559566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf,MPIU_INT,visits,lvisits,MPI_REPLACE));
9569566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd (sf,MPIU_INT,visits,lvisits));
957cb1e1211SMatthew G Knepley
9589566063dSJacob Faibussowitsch PetscCall(PetscSFGetRootRanks());
959cb1e1211SMatthew G Knepley
9609566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(maxClosureSize*maxClosureSize,&cellmat,npoints,&owner));
961cb1e1211SMatthew G Knepley for (c=cStart; c<cEnd; c++) {
9629566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(cellmat,maxClosureSize*maxClosureSize));
963cb1e1211SMatthew G Knepley /*
964cb1e1211SMatthew G Knepley Depth-first walk of transitive closure.
965cb1e1211SMatthew G Knepley At each leaf frame f of transitive closure that we see, add 1/visits[f] to each pair (p,q) not marked as done in cellmat.
966cb1e1211SMatthew G Knepley This contribution is added to dnz if owning ranks of p and q match, to onz otherwise.
967cb1e1211SMatthew G Knepley */
968cb1e1211SMatthew G Knepley }
969cb1e1211SMatthew G Knepley
9709566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf,MPIU_INT,ldnz,dnz,MPI_SUM));
9719566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd (sf,MPIU_INT,lonz,onz,MPI_SUM));
9723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
973cb1e1211SMatthew G Knepley }
974cb1e1211SMatthew G Knepley #endif
975