xref: /petsc/src/vec/is/sf/utils/sfutils.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1b0c7db22SLisandro Dalcin #include <petsc/private/sfimpl.h>       /*I  "petscsf.h"   I*/
2b0c7db22SLisandro Dalcin #include <petsc/private/sectionimpl.h>
3b0c7db22SLisandro Dalcin 
4b0c7db22SLisandro Dalcin /*@C
5b0c7db22SLisandro Dalcin    PetscSFSetGraphLayout - Set a parallel star forest via global indices and a PetscLayout
6b0c7db22SLisandro Dalcin 
7b0c7db22SLisandro Dalcin    Collective
8b0c7db22SLisandro Dalcin 
94165533cSJose E. Roman    Input Parameters:
10b0c7db22SLisandro Dalcin +  sf - star forest
11ddea5d60SJunchao Zhang .  layout - PetscLayout defining the global space for roots
12b0c7db22SLisandro Dalcin .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
13b0c7db22SLisandro Dalcin .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage
14b0c7db22SLisandro Dalcin .  localmode - copy mode for ilocal
1586c56f52SVaclav Hapla -  iremote - root vertices in global numbering corresponding to leaves in ilocal
16b0c7db22SLisandro Dalcin 
17b0c7db22SLisandro Dalcin    Level: intermediate
18b0c7db22SLisandro Dalcin 
1986c56f52SVaclav Hapla    Notes:
2086c56f52SVaclav Hapla    Global indices must lie in [0, N) where N is the global size of layout.
21db2b9530SVaclav Hapla    Leaf indices in ilocal get sorted; this means the user-provided array gets sorted if localmode is PETSC_OWN_POINTER.
2286c56f52SVaclav Hapla 
23db2b9530SVaclav Hapla    Developer Notes:
2486c56f52SVaclav Hapla    Local indices which are the identity permutation in the range [0,nleaves) are discarded as they
25b0c7db22SLisandro Dalcin    encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not
26b0c7db22SLisandro Dalcin    needed
27b0c7db22SLisandro Dalcin 
28b0c7db22SLisandro Dalcin .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
29b0c7db22SLisandro Dalcin @*/
30db2b9530SVaclav Hapla PetscErrorCode PetscSFSetGraphLayout(PetscSF sf,PetscLayout layout,PetscInt nleaves,PetscInt *ilocal,PetscCopyMode localmode,const PetscInt *iremote)
31b0c7db22SLisandro Dalcin {
3238a25198SStefano Zampini   const PetscInt *range;
3338a25198SStefano Zampini   PetscInt       i, nroots, ls = -1, ln = -1;
3438a25198SStefano Zampini   PetscMPIInt    lr = -1;
35b0c7db22SLisandro Dalcin   PetscSFNode    *remote;
36b0c7db22SLisandro Dalcin 
37b0c7db22SLisandro Dalcin   PetscFunctionBegin;
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetLocalSize(layout,&nroots));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetRanges(layout,&range));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nleaves,&remote));
4138a25198SStefano Zampini   if (nleaves) { ls = iremote[0] + 1; }
42b0c7db22SLisandro Dalcin   for (i=0; i<nleaves; i++) {
4338a25198SStefano Zampini     const PetscInt idx = iremote[i] - ls;
4438a25198SStefano Zampini     if (idx < 0 || idx >= ln) { /* short-circuit the search */
455f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLayoutFindOwnerIndex(layout,iremote[i],&lr,&remote[i].index));
4638a25198SStefano Zampini       remote[i].rank = lr;
4738a25198SStefano Zampini       ls = range[lr];
4838a25198SStefano Zampini       ln = range[lr+1] - ls;
4938a25198SStefano Zampini     } else {
5038a25198SStefano Zampini       remote[i].rank  = lr;
5138a25198SStefano Zampini       remote[i].index = idx;
5238a25198SStefano Zampini     }
53b0c7db22SLisandro Dalcin   }
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetGraph(sf,nroots,nleaves,ilocal,localmode,remote,PETSC_OWN_POINTER));
55b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
56b0c7db22SLisandro Dalcin }
57b0c7db22SLisandro Dalcin 
58b0c7db22SLisandro Dalcin /*@
59b0c7db22SLisandro Dalcin   PetscSFSetGraphSection - Sets the PetscSF graph encoding the parallel dof overlap based upon the PetscSections
60b0c7db22SLisandro Dalcin   describing the data layout.
61b0c7db22SLisandro Dalcin 
62b0c7db22SLisandro Dalcin   Input Parameters:
63b0c7db22SLisandro Dalcin + sf - The SF
64b0c7db22SLisandro Dalcin . localSection - PetscSection describing the local data layout
65b0c7db22SLisandro Dalcin - globalSection - PetscSection describing the global data layout
66b0c7db22SLisandro Dalcin 
67b0c7db22SLisandro Dalcin   Level: developer
68b0c7db22SLisandro Dalcin 
69b0c7db22SLisandro Dalcin .seealso: PetscSFSetGraph(), PetscSFSetGraphLayout()
70b0c7db22SLisandro Dalcin @*/
71b0c7db22SLisandro Dalcin PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection)
72b0c7db22SLisandro Dalcin {
73b0c7db22SLisandro Dalcin   MPI_Comm       comm;
74b0c7db22SLisandro Dalcin   PetscLayout    layout;
75b0c7db22SLisandro Dalcin   const PetscInt *ranges;
76b0c7db22SLisandro Dalcin   PetscInt       *local;
77b0c7db22SLisandro Dalcin   PetscSFNode    *remote;
78b0c7db22SLisandro Dalcin   PetscInt       pStart, pEnd, p, nroots, nleaves = 0, l;
79b0c7db22SLisandro Dalcin   PetscMPIInt    size, rank;
80b0c7db22SLisandro Dalcin 
81b0c7db22SLisandro Dalcin   PetscFunctionBegin;
82b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
83b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2);
84b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3);
85b0c7db22SLisandro Dalcin 
865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)sf,&comm));
875f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm, &size));
885f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetChart(globalSection, &pStart, &pEnd));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetConstrainedStorageSize(globalSection, &nroots));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutCreate(comm, &layout));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetBlockSize(layout, 1));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetLocalSize(layout, nroots));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetUp(layout));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetRanges(layout, &ranges));
96b0c7db22SLisandro Dalcin   for (p = pStart; p < pEnd; ++p) {
97b0c7db22SLisandro Dalcin     PetscInt gdof, gcdof;
98b0c7db22SLisandro Dalcin 
995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(globalSection, p, &gdof));
1005f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
1012c71b3e2SJacob Faibussowitsch     PetscCheckFalse(gcdof > (gdof < 0 ? -(gdof+1) : gdof),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " has %" PetscInt_FMT " constraints > %" PetscInt_FMT " dof", p, gcdof, (gdof < 0 ? -(gdof+1) : gdof));
102b0c7db22SLisandro Dalcin     nleaves += gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
103b0c7db22SLisandro Dalcin   }
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nleaves, &local));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nleaves, &remote));
106b0c7db22SLisandro Dalcin   for (p = pStart, l = 0; p < pEnd; ++p) {
107b0c7db22SLisandro Dalcin     const PetscInt *cind;
108b0c7db22SLisandro Dalcin     PetscInt       dof, cdof, off, gdof, gcdof, goff, gsize, d, c;
109b0c7db22SLisandro Dalcin 
1105f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(localSection, p, &dof));
1115f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(localSection, p, &off));
1125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetConstraintDof(localSection, p, &cdof));
1135f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetConstraintIndices(localSection, p, &cind));
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(globalSection, p, &gdof));
1155f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetConstraintDof(globalSection, p, &gcdof));
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(globalSection, p, &goff));
117b0c7db22SLisandro Dalcin     if (!gdof) continue; /* Censored point */
118b0c7db22SLisandro Dalcin     gsize = gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
119b0c7db22SLisandro Dalcin     if (gsize != dof-cdof) {
1202c71b3e2SJacob Faibussowitsch       PetscCheckFalse(gsize != dof,comm, PETSC_ERR_ARG_WRONG, "Global dof %" PetscInt_FMT " for point %" PetscInt_FMT " is neither the constrained size %" PetscInt_FMT ", nor the unconstrained %" PetscInt_FMT, gsize, p, dof-cdof, dof);
121b0c7db22SLisandro Dalcin       cdof = 0; /* Ignore constraints */
122b0c7db22SLisandro Dalcin     }
123b0c7db22SLisandro Dalcin     for (d = 0, c = 0; d < dof; ++d) {
124b0c7db22SLisandro Dalcin       if ((c < cdof) && (cind[c] == d)) {++c; continue;}
125b0c7db22SLisandro Dalcin       local[l+d-c] = off+d;
126b0c7db22SLisandro Dalcin     }
1272c71b3e2SJacob Faibussowitsch     PetscCheckFalse(d-c != gsize,comm, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT ": Global dof %" PetscInt_FMT " != %" PetscInt_FMT " size - number of constraints", p, gsize, d-c);
128b0c7db22SLisandro Dalcin     if (gdof < 0) {
129b0c7db22SLisandro Dalcin       for (d = 0; d < gsize; ++d, ++l) {
130b0c7db22SLisandro Dalcin         PetscInt offset = -(goff+1) + d, r;
131b0c7db22SLisandro Dalcin 
1325f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFindInt(offset,size+1,ranges,&r));
133b0c7db22SLisandro Dalcin         if (r < 0) r = -(r+2);
1342c71b3e2SJacob Faibussowitsch         PetscCheckFalse((r < 0) || (r >= size),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " mapped to invalid process %" PetscInt_FMT " (%" PetscInt_FMT ", %" PetscInt_FMT ")", p, r, gdof, goff);
135b0c7db22SLisandro Dalcin         remote[l].rank  = r;
136b0c7db22SLisandro Dalcin         remote[l].index = offset - ranges[r];
137b0c7db22SLisandro Dalcin       }
138b0c7db22SLisandro Dalcin     } else {
139b0c7db22SLisandro Dalcin       for (d = 0; d < gsize; ++d, ++l) {
140b0c7db22SLisandro Dalcin         remote[l].rank  = rank;
141b0c7db22SLisandro Dalcin         remote[l].index = goff+d - ranges[rank];
142b0c7db22SLisandro Dalcin       }
143b0c7db22SLisandro Dalcin     }
144b0c7db22SLisandro Dalcin   }
1452c71b3e2SJacob Faibussowitsch   PetscCheckFalse(l != nleaves,comm, PETSC_ERR_PLIB, "Iteration error, l %" PetscInt_FMT " != nleaves %" PetscInt_FMT, l, nleaves);
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutDestroy(&layout));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER));
148b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
149b0c7db22SLisandro Dalcin }
150b0c7db22SLisandro Dalcin 
151b0c7db22SLisandro Dalcin static PetscErrorCode PetscSectionCheckConstraints_Static(PetscSection s)
152b0c7db22SLisandro Dalcin {
153b0c7db22SLisandro Dalcin   PetscFunctionBegin;
154b0c7db22SLisandro Dalcin   if (!s->bc) {
1555f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCreate(PETSC_COMM_SELF, &s->bc));
1565f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetChart(s->bc, s->pStart, s->pEnd));
157b0c7db22SLisandro Dalcin   }
158b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
159b0c7db22SLisandro Dalcin }
160b0c7db22SLisandro Dalcin 
161b0c7db22SLisandro Dalcin /*@C
162b0c7db22SLisandro Dalcin   PetscSFDistributeSection - Create a new PetscSection reorganized, moving from the root to the leaves of the SF
163b0c7db22SLisandro Dalcin 
164b0c7db22SLisandro Dalcin   Collective on sf
165b0c7db22SLisandro Dalcin 
166b0c7db22SLisandro Dalcin   Input Parameters:
167b0c7db22SLisandro Dalcin + sf - The SF
168b0c7db22SLisandro Dalcin - rootSection - Section defined on root space
169b0c7db22SLisandro Dalcin 
170b0c7db22SLisandro Dalcin   Output Parameters:
171b0c7db22SLisandro Dalcin + remoteOffsets - root offsets in leaf storage, or NULL
172b0c7db22SLisandro Dalcin - leafSection - Section defined on the leaf space
173b0c7db22SLisandro Dalcin 
174b0c7db22SLisandro Dalcin   Level: advanced
175b0c7db22SLisandro Dalcin 
176b0c7db22SLisandro Dalcin .seealso: PetscSFCreate()
177b0c7db22SLisandro Dalcin @*/
178b0c7db22SLisandro Dalcin PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection)
179b0c7db22SLisandro Dalcin {
180b0c7db22SLisandro Dalcin   PetscSF        embedSF;
181b0c7db22SLisandro Dalcin   const PetscInt *indices;
182b0c7db22SLisandro Dalcin   IS             selected;
183b0c7db22SLisandro Dalcin   PetscInt       numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c;
184b0c7db22SLisandro Dalcin   PetscBool      *sub, hasc;
185b0c7db22SLisandro Dalcin 
186b0c7db22SLisandro Dalcin   PetscFunctionBegin;
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(PETSCSF_DistSect,sf,0,0,0));
1885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetNumFields(rootSection, &numFields));
189029e8818Sksagiyam   if (numFields) {
190029e8818Sksagiyam     IS perm;
191029e8818Sksagiyam 
192029e8818Sksagiyam     /* PetscSectionSetNumFields() calls PetscSectionReset(), which destroys
193029e8818Sksagiyam        leafSection->perm. To keep this permutation set by the user, we grab
194029e8818Sksagiyam        the reference before calling PetscSectionSetNumFields() and set it
195029e8818Sksagiyam        back after. */
1965f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetPermutation(leafSection, &perm));
1975f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectReference((PetscObject)perm));
1985f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetNumFields(leafSection, numFields));
1995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetPermutation(leafSection, perm));
2005f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&perm));
201029e8818Sksagiyam   }
2025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(numFields+2, &sub));
203b0c7db22SLisandro Dalcin   sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE;
204b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) {
2052ee82556SMatthew G. Knepley     PetscSectionSym sym, dsym = NULL;
206b0c7db22SLisandro Dalcin     const char      *name   = NULL;
207b0c7db22SLisandro Dalcin     PetscInt        numComp = 0;
208b0c7db22SLisandro Dalcin 
209b0c7db22SLisandro Dalcin     sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE;
2105f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetFieldComponents(rootSection, f, &numComp));
2115f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetFieldName(rootSection, f, &name));
2125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetFieldSym(rootSection, f, &sym));
2135f80ce2aSJacob Faibussowitsch     if (sym) CHKERRQ(PetscSectionSymDistribute(sym, sf, &dsym));
2145f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetFieldComponents(leafSection, f, numComp));
2155f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetFieldName(leafSection, f, name));
2165f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetFieldSym(leafSection, f, dsym));
2175f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSymDestroy(&dsym));
218b0c7db22SLisandro Dalcin     for (c = 0; c < rootSection->numFieldComponents[f]; ++c) {
2195f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetComponentName(rootSection, f, c, &name));
2205f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionSetComponentName(leafSection, f, c, name));
221b0c7db22SLisandro Dalcin     }
222b0c7db22SLisandro Dalcin   }
2235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetChart(rootSection, &rpStart, &rpEnd));
2245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetGraph(sf,&nroots,NULL,NULL,NULL));
225b0c7db22SLisandro Dalcin   rpEnd = PetscMin(rpEnd,nroots);
226b0c7db22SLisandro Dalcin   rpEnd = PetscMax(rpStart,rpEnd);
227b0c7db22SLisandro Dalcin   /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */
228b0c7db22SLisandro Dalcin   sub[0] = (PetscBool)(nroots != rpEnd - rpStart);
2295f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPIU_Allreduce(MPI_IN_PLACE, sub, 2+numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf)));
230b0c7db22SLisandro Dalcin   if (sub[0]) {
2315f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected));
2325f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(selected, &indices));
2335f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF));
2345f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(selected, &indices));
2355f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&selected));
236b0c7db22SLisandro Dalcin   } else {
2375f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectReference((PetscObject)sf));
238b0c7db22SLisandro Dalcin     embedSF = sf;
239b0c7db22SLisandro Dalcin   }
2405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd));
241b0c7db22SLisandro Dalcin   lpEnd++;
242b0c7db22SLisandro Dalcin 
2435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetChart(leafSection, lpStart, lpEnd));
244b0c7db22SLisandro Dalcin 
245b0c7db22SLisandro Dalcin   /* Constrained dof section */
246b0c7db22SLisandro Dalcin   hasc = sub[1];
247b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2+f]);
248b0c7db22SLisandro Dalcin 
249b0c7db22SLisandro Dalcin   /* Could fuse these at the cost of copies and extra allocation */
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart],MPI_REPLACE));
2515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart],MPI_REPLACE));
252b0c7db22SLisandro Dalcin   if (sub[1]) {
2535f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCheckConstraints_Static(rootSection));
2545f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCheckConstraints_Static(leafSection));
2555f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart],MPI_REPLACE));
2565f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart],MPI_REPLACE));
257b0c7db22SLisandro Dalcin   }
258b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) {
2595f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart],MPI_REPLACE));
2605f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart],MPI_REPLACE));
261b0c7db22SLisandro Dalcin     if (sub[2+f]) {
2625f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionCheckConstraints_Static(rootSection->field[f]));
2635f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionCheckConstraints_Static(leafSection->field[f]));
2645f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart],MPI_REPLACE));
2655f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart],MPI_REPLACE));
266b0c7db22SLisandro Dalcin     }
267b0c7db22SLisandro Dalcin   }
268b0c7db22SLisandro Dalcin   if (remoteOffsets) {
2695f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(lpEnd - lpStart, remoteOffsets));
2705f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE));
2715f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE));
272b0c7db22SLisandro Dalcin   }
2735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetUp(leafSection));
274b0c7db22SLisandro Dalcin   if (hasc) { /* need to communicate bcIndices */
275b0c7db22SLisandro Dalcin     PetscSF  bcSF;
276b0c7db22SLisandro Dalcin     PetscInt *rOffBc;
277b0c7db22SLisandro Dalcin 
2785f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(lpEnd - lpStart, &rOffBc));
279b0c7db22SLisandro Dalcin     if (sub[1]) {
2805f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE));
2815f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE));
2825f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF));
2835f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices,MPI_REPLACE));
2845f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices,MPI_REPLACE));
2855f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFDestroy(&bcSF));
286b0c7db22SLisandro Dalcin     }
287b0c7db22SLisandro Dalcin     for (f = 0; f < numFields; ++f) {
288b0c7db22SLisandro Dalcin       if (sub[2+f]) {
2895f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE));
2905f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart],MPI_REPLACE));
2915f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF));
2925f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices,MPI_REPLACE));
2935f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices,MPI_REPLACE));
2945f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSFDestroy(&bcSF));
295b0c7db22SLisandro Dalcin       }
296b0c7db22SLisandro Dalcin     }
2975f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(rOffBc));
298b0c7db22SLisandro Dalcin   }
2995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&embedSF));
3005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(sub));
3015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(PETSCSF_DistSect,sf,0,0,0));
302b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
303b0c7db22SLisandro Dalcin }
304b0c7db22SLisandro Dalcin 
305b0c7db22SLisandro Dalcin /*@C
306b0c7db22SLisandro Dalcin   PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes
307b0c7db22SLisandro Dalcin 
308b0c7db22SLisandro Dalcin   Collective on sf
309b0c7db22SLisandro Dalcin 
310b0c7db22SLisandro Dalcin   Input Parameters:
311b0c7db22SLisandro Dalcin + sf          - The SF
312b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots)
313b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data  (this is layout for SF leaves)
314b0c7db22SLisandro Dalcin 
315b0c7db22SLisandro Dalcin   Output Parameter:
316b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
317b0c7db22SLisandro Dalcin 
318b0c7db22SLisandro Dalcin   Level: developer
319b0c7db22SLisandro Dalcin 
320b0c7db22SLisandro Dalcin .seealso: PetscSFCreate()
321b0c7db22SLisandro Dalcin @*/
322b0c7db22SLisandro Dalcin PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets)
323b0c7db22SLisandro Dalcin {
324b0c7db22SLisandro Dalcin   PetscSF         embedSF;
325b0c7db22SLisandro Dalcin   const PetscInt *indices;
326b0c7db22SLisandro Dalcin   IS              selected;
327b0c7db22SLisandro Dalcin   PetscInt        numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0;
328b0c7db22SLisandro Dalcin 
329b0c7db22SLisandro Dalcin   PetscFunctionBegin;
330b0c7db22SLisandro Dalcin   *remoteOffsets = NULL;
3315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL));
332b0c7db22SLisandro Dalcin   if (numRoots < 0) PetscFunctionReturn(0);
3335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(PETSCSF_RemoteOff,sf,0,0,0));
3345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetChart(rootSection, &rpStart, &rpEnd));
3355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetChart(leafSection, &lpStart, &lpEnd));
3365f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected));
3375f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(selected, &indices));
3385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreateEmbeddedRootSF(sf, rpEnd - rpStart, indices, &embedSF));
3395f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(selected, &indices));
3405f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&selected));
3415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(lpEnd - lpStart, remoteOffsets));
3425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE));
3435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart],MPI_REPLACE));
3445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&embedSF));
3455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(PETSCSF_RemoteOff,sf,0,0,0));
346b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
347b0c7db22SLisandro Dalcin }
348b0c7db22SLisandro Dalcin 
349b0c7db22SLisandro Dalcin /*@C
350b0c7db22SLisandro Dalcin   PetscSFCreateSectionSF - Create an expanded SF of dofs, assuming the input SF relates points
351b0c7db22SLisandro Dalcin 
352b0c7db22SLisandro Dalcin   Collective on sf
353b0c7db22SLisandro Dalcin 
354b0c7db22SLisandro Dalcin   Input Parameters:
355b0c7db22SLisandro Dalcin + sf - The SF
356b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is usually the serial section)
357b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
358b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data  (this is the distributed section)
359b0c7db22SLisandro Dalcin 
360b0c7db22SLisandro Dalcin   Output Parameters:
361b0c7db22SLisandro Dalcin - sectionSF - The new SF
362b0c7db22SLisandro Dalcin 
363b0c7db22SLisandro Dalcin   Note: Either rootSection or remoteOffsets can be specified
364b0c7db22SLisandro Dalcin 
365b0c7db22SLisandro Dalcin   Level: advanced
366b0c7db22SLisandro Dalcin 
367b0c7db22SLisandro Dalcin .seealso: PetscSFCreate()
368b0c7db22SLisandro Dalcin @*/
369b0c7db22SLisandro Dalcin PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF)
370b0c7db22SLisandro Dalcin {
371b0c7db22SLisandro Dalcin   MPI_Comm          comm;
372b0c7db22SLisandro Dalcin   const PetscInt    *localPoints;
373b0c7db22SLisandro Dalcin   const PetscSFNode *remotePoints;
374b0c7db22SLisandro Dalcin   PetscInt          lpStart, lpEnd;
375b0c7db22SLisandro Dalcin   PetscInt          numRoots, numSectionRoots, numPoints, numIndices = 0;
376b0c7db22SLisandro Dalcin   PetscInt          *localIndices;
377b0c7db22SLisandro Dalcin   PetscSFNode       *remoteIndices;
378b0c7db22SLisandro Dalcin   PetscInt          i, ind;
379b0c7db22SLisandro Dalcin 
380b0c7db22SLisandro Dalcin   PetscFunctionBegin;
381b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
382b0c7db22SLisandro Dalcin   PetscValidPointer(rootSection,2);
383b0c7db22SLisandro Dalcin   /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */
384b0c7db22SLisandro Dalcin   PetscValidPointer(leafSection,4);
385b0c7db22SLisandro Dalcin   PetscValidPointer(sectionSF,5);
3865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)sf,&comm));
3875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreate(comm, sectionSF));
3885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetChart(leafSection, &lpStart, &lpEnd));
3895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetStorageSize(rootSection, &numSectionRoots));
3905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints));
391b0c7db22SLisandro Dalcin   if (numRoots < 0) PetscFunctionReturn(0);
3925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(PETSCSF_SectSF,sf,0,0,0));
393b0c7db22SLisandro Dalcin   for (i = 0; i < numPoints; ++i) {
394b0c7db22SLisandro Dalcin     PetscInt localPoint = localPoints ? localPoints[i] : i;
395b0c7db22SLisandro Dalcin     PetscInt dof;
396b0c7db22SLisandro Dalcin 
397b0c7db22SLisandro Dalcin     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
3985f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(leafSection, localPoint, &dof));
399b0c7db22SLisandro Dalcin       numIndices += dof;
400b0c7db22SLisandro Dalcin     }
401b0c7db22SLisandro Dalcin   }
4025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(numIndices, &localIndices));
4035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(numIndices, &remoteIndices));
404b0c7db22SLisandro Dalcin   /* Create new index graph */
405b0c7db22SLisandro Dalcin   for (i = 0, ind = 0; i < numPoints; ++i) {
406b0c7db22SLisandro Dalcin     PetscInt localPoint = localPoints ? localPoints[i] : i;
407b0c7db22SLisandro Dalcin     PetscInt rank       = remotePoints[i].rank;
408b0c7db22SLisandro Dalcin 
409b0c7db22SLisandro Dalcin     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
410b0c7db22SLisandro Dalcin       PetscInt remoteOffset = remoteOffsets[localPoint-lpStart];
411b0c7db22SLisandro Dalcin       PetscInt loff, dof, d;
412b0c7db22SLisandro Dalcin 
4135f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetOffset(leafSection, localPoint, &loff));
4145f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(leafSection, localPoint, &dof));
415b0c7db22SLisandro Dalcin       for (d = 0; d < dof; ++d, ++ind) {
416b0c7db22SLisandro Dalcin         localIndices[ind]        = loff+d;
417b0c7db22SLisandro Dalcin         remoteIndices[ind].rank  = rank;
418b0c7db22SLisandro Dalcin         remoteIndices[ind].index = remoteOffset+d;
419b0c7db22SLisandro Dalcin       }
420b0c7db22SLisandro Dalcin     }
421b0c7db22SLisandro Dalcin   }
4222c71b3e2SJacob Faibussowitsch   PetscCheckFalse(numIndices != ind,comm, PETSC_ERR_PLIB, "Inconsistency in indices, %" PetscInt_FMT " should be %" PetscInt_FMT, ind, numIndices);
4235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER));
4245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetUp(*sectionSF));
4255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(PETSCSF_SectSF,sf,0,0,0));
426b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
427b0c7db22SLisandro Dalcin }
4288fa9e22eSVaclav Hapla 
4298fa9e22eSVaclav Hapla /*@C
4308fa9e22eSVaclav Hapla    PetscSFCreateFromLayouts - Creates a parallel star forest mapping two PetscLayout objects
4318fa9e22eSVaclav Hapla 
4328fa9e22eSVaclav Hapla    Collective
4338fa9e22eSVaclav Hapla 
4344165533cSJose E. Roman    Input Parameters:
4358fa9e22eSVaclav Hapla +  rmap - PetscLayout defining the global root space
4368fa9e22eSVaclav Hapla -  lmap - PetscLayout defining the global leaf space
4378fa9e22eSVaclav Hapla 
4384165533cSJose E. Roman    Output Parameter:
4398fa9e22eSVaclav Hapla .  sf - The parallel star forest
4408fa9e22eSVaclav Hapla 
4418fa9e22eSVaclav Hapla    Level: intermediate
4428fa9e22eSVaclav Hapla 
4438fa9e22eSVaclav Hapla .seealso: PetscSFCreate(), PetscLayoutCreate(), PetscSFSetGraphLayout()
4448fa9e22eSVaclav Hapla @*/
4458fa9e22eSVaclav Hapla PetscErrorCode PetscSFCreateFromLayouts(PetscLayout rmap, PetscLayout lmap, PetscSF* sf)
4468fa9e22eSVaclav Hapla {
4478fa9e22eSVaclav Hapla   PetscInt       i,nroots,nleaves = 0;
4488fa9e22eSVaclav Hapla   PetscInt       rN, lst, len;
4498fa9e22eSVaclav Hapla   PetscMPIInt    owner = -1;
4508fa9e22eSVaclav Hapla   PetscSFNode    *remote;
4518fa9e22eSVaclav Hapla   MPI_Comm       rcomm = rmap->comm;
4528fa9e22eSVaclav Hapla   MPI_Comm       lcomm = lmap->comm;
4538fa9e22eSVaclav Hapla   PetscMPIInt    flg;
4548fa9e22eSVaclav Hapla 
4558fa9e22eSVaclav Hapla   PetscFunctionBegin;
4568fa9e22eSVaclav Hapla   PetscValidPointer(sf,3);
457*28b400f6SJacob Faibussowitsch   PetscCheck(rmap->setupcalled,rcomm,PETSC_ERR_ARG_WRONGSTATE,"Root layout not setup");
458*28b400f6SJacob Faibussowitsch   PetscCheck(lmap->setupcalled,lcomm,PETSC_ERR_ARG_WRONGSTATE,"Leaf layout not setup");
4595f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_compare(rcomm,lcomm,&flg));
4602c71b3e2SJacob Faibussowitsch   PetscCheckFalse(flg != MPI_CONGRUENT && flg != MPI_IDENT,rcomm,PETSC_ERR_SUP,"cannot map two layouts with non-matching communicators");
4615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreate(rcomm,sf));
4625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetLocalSize(rmap,&nroots));
4635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetSize(rmap,&rN));
4645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetRange(lmap,&lst,&len));
4655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(len-lst,&remote));
4668fa9e22eSVaclav Hapla   for (i = lst; i < len && i < rN; i++) {
4678fa9e22eSVaclav Hapla     if (owner < -1 || i >= rmap->range[owner+1]) {
4685f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLayoutFindOwner(rmap,i,&owner));
4698fa9e22eSVaclav Hapla     }
4708fa9e22eSVaclav Hapla     remote[nleaves].rank  = owner;
4718fa9e22eSVaclav Hapla     remote[nleaves].index = i - rmap->range[owner];
4728fa9e22eSVaclav Hapla     nleaves++;
4738fa9e22eSVaclav Hapla   }
4745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetGraph(*sf,nroots,nleaves,NULL,PETSC_OWN_POINTER,remote,PETSC_COPY_VALUES));
4755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(remote));
4768fa9e22eSVaclav Hapla   PetscFunctionReturn(0);
4778fa9e22eSVaclav Hapla }
4788fa9e22eSVaclav Hapla 
4798fa9e22eSVaclav Hapla /* TODO: handle nooffprocentries like MatZeroRowsMapLocal_Private, since this code is the same */
4808fa9e22eSVaclav Hapla PetscErrorCode PetscLayoutMapLocal(PetscLayout map,PetscInt N,const PetscInt idxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs)
4818fa9e22eSVaclav Hapla {
4828fa9e22eSVaclav Hapla   PetscInt      *owners = map->range;
4838fa9e22eSVaclav Hapla   PetscInt       n      = map->n;
4848fa9e22eSVaclav Hapla   PetscSF        sf;
4858fa9e22eSVaclav Hapla   PetscInt      *lidxs,*work = NULL;
4868fa9e22eSVaclav Hapla   PetscSFNode   *ridxs;
4878fa9e22eSVaclav Hapla   PetscMPIInt    rank, p = 0;
4888fa9e22eSVaclav Hapla   PetscInt       r, len = 0;
4898fa9e22eSVaclav Hapla 
4908fa9e22eSVaclav Hapla   PetscFunctionBegin;
4918fa9e22eSVaclav Hapla   if (on) *on = 0;              /* squelch -Wmaybe-uninitialized */
4928fa9e22eSVaclav Hapla   /* Create SF where leaves are input idxs and roots are owned idxs */
4935f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(map->comm,&rank));
4945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&lidxs));
4958fa9e22eSVaclav Hapla   for (r = 0; r < n; ++r) lidxs[r] = -1;
4965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(N,&ridxs));
4978fa9e22eSVaclav Hapla   for (r = 0; r < N; ++r) {
4988fa9e22eSVaclav Hapla     const PetscInt idx = idxs[r];
4992c71b3e2SJacob Faibussowitsch     PetscCheckFalse(idx < 0 || map->N <= idx,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %" PetscInt_FMT " out of range [0,%" PetscInt_FMT ")",idx,map->N);
5008fa9e22eSVaclav Hapla     if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */
5015f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscLayoutFindOwner(map,idx,&p));
5028fa9e22eSVaclav Hapla     }
5038fa9e22eSVaclav Hapla     ridxs[r].rank = p;
5048fa9e22eSVaclav Hapla     ridxs[r].index = idxs[r] - owners[p];
5058fa9e22eSVaclav Hapla   }
5065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreate(map->comm,&sf));
5075f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER));
5085f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR));
5095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR));
5108fa9e22eSVaclav Hapla   if (ogidxs) { /* communicate global idxs */
5118fa9e22eSVaclav Hapla     PetscInt cum = 0,start,*work2;
5128fa9e22eSVaclav Hapla 
5135f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(n,&work));
5145f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(N,&work2));
5158fa9e22eSVaclav Hapla     for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++;
5165f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm));
5178fa9e22eSVaclav Hapla     start -= cum;
5188fa9e22eSVaclav Hapla     cum = 0;
5198fa9e22eSVaclav Hapla     for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++;
5205f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPI_REPLACE));
5215f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPI_REPLACE));
5225f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(work2));
5238fa9e22eSVaclav Hapla   }
5245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sf));
5258fa9e22eSVaclav Hapla   /* Compress and put in indices */
5268fa9e22eSVaclav Hapla   for (r = 0; r < n; ++r)
5278fa9e22eSVaclav Hapla     if (lidxs[r] >= 0) {
5288fa9e22eSVaclav Hapla       if (work) work[len] = work[r];
5298fa9e22eSVaclav Hapla       lidxs[len++] = r;
5308fa9e22eSVaclav Hapla     }
5318fa9e22eSVaclav Hapla   if (on) *on = len;
5328fa9e22eSVaclav Hapla   if (oidxs) *oidxs = lidxs;
5338fa9e22eSVaclav Hapla   if (ogidxs) *ogidxs = work;
5348fa9e22eSVaclav Hapla   PetscFunctionReturn(0);
5358fa9e22eSVaclav Hapla }
536deffd5ebSksagiyam 
537deffd5ebSksagiyam /*@
538deffd5ebSksagiyam   PetscSFCreateByMatchingIndices - Create SF by matching root and leaf indices
539deffd5ebSksagiyam 
540deffd5ebSksagiyam   Collective
541deffd5ebSksagiyam 
542deffd5ebSksagiyam   Input Parameters:
543deffd5ebSksagiyam + layout - PetscLayout defining the global index space and the rank that brokers each index
544deffd5ebSksagiyam . numRootIndices - size of rootIndices
545deffd5ebSksagiyam . rootIndices - PetscInt array of global indices of which this process requests ownership
546deffd5ebSksagiyam . rootLocalIndices - root local index permutation (NULL if no permutation)
547deffd5ebSksagiyam . rootLocalOffset - offset to be added to root local indices
548deffd5ebSksagiyam . numLeafIndices - size of leafIndices
549deffd5ebSksagiyam . leafIndices - PetscInt array of global indices with which this process requires data associated
550deffd5ebSksagiyam . leafLocalIndices - leaf local index permutation (NULL if no permutation)
551deffd5ebSksagiyam - leafLocalOffset - offset to be added to leaf local indices
552deffd5ebSksagiyam 
553d8d19677SJose E. Roman   Output Parameters:
554deffd5ebSksagiyam + sfA - star forest representing the communication pattern from the layout space to the leaf space (NULL if not needed)
555deffd5ebSksagiyam - sf - star forest representing the communication pattern from the root space to the leaf space
556deffd5ebSksagiyam 
557deffd5ebSksagiyam   Example 1:
558deffd5ebSksagiyam $
559deffd5ebSksagiyam $  rank             : 0            1            2
560deffd5ebSksagiyam $  rootIndices      : [1 0 2]      [3]          [3]
561deffd5ebSksagiyam $  rootLocalOffset  : 100          200          300
562deffd5ebSksagiyam $  layout           : [0 1]        [2]          [3]
563deffd5ebSksagiyam $  leafIndices      : [0]          [2]          [0 3]
564deffd5ebSksagiyam $  leafLocalOffset  : 400          500          600
565deffd5ebSksagiyam $
5667ef5781fSVaclav Hapla would build the following SF
567deffd5ebSksagiyam $
568deffd5ebSksagiyam $  [0] 400 <- (0,101)
569deffd5ebSksagiyam $  [1] 500 <- (0,102)
570deffd5ebSksagiyam $  [2] 600 <- (0,101)
571deffd5ebSksagiyam $  [2] 601 <- (2,300)
572deffd5ebSksagiyam $
573deffd5ebSksagiyam   Example 2:
574deffd5ebSksagiyam $
575deffd5ebSksagiyam $  rank             : 0               1               2
576deffd5ebSksagiyam $  rootIndices      : [1 0 2]         [3]             [3]
577deffd5ebSksagiyam $  rootLocalOffset  : 100             200             300
578deffd5ebSksagiyam $  layout           : [0 1]           [2]             [3]
579deffd5ebSksagiyam $  leafIndices      : rootIndices     rootIndices     rootIndices
580deffd5ebSksagiyam $  leafLocalOffset  : rootLocalOffset rootLocalOffset rootLocalOffset
581deffd5ebSksagiyam $
5827ef5781fSVaclav Hapla would build the following SF
583deffd5ebSksagiyam $
584deffd5ebSksagiyam $  [1] 200 <- (2,300)
585deffd5ebSksagiyam $
586deffd5ebSksagiyam   Example 3:
587deffd5ebSksagiyam $
588deffd5ebSksagiyam $  No process requests ownership of global index 1, but no process needs it.
589deffd5ebSksagiyam $
590deffd5ebSksagiyam $  rank             : 0            1            2
591deffd5ebSksagiyam $  numRootIndices   : 2            1            1
592deffd5ebSksagiyam $  rootIndices      : [0 2]        [3]          [3]
593deffd5ebSksagiyam $  rootLocalOffset  : 100          200          300
594deffd5ebSksagiyam $  layout           : [0 1]        [2]          [3]
595deffd5ebSksagiyam $  numLeafIndices   : 1            1            2
596deffd5ebSksagiyam $  leafIndices      : [0]          [2]          [0 3]
597deffd5ebSksagiyam $  leafLocalOffset  : 400          500          600
598deffd5ebSksagiyam $
5997ef5781fSVaclav Hapla would build the following SF
600deffd5ebSksagiyam $
601deffd5ebSksagiyam $  [0] 400 <- (0,100)
602deffd5ebSksagiyam $  [1] 500 <- (0,101)
603deffd5ebSksagiyam $  [2] 600 <- (0,100)
604deffd5ebSksagiyam $  [2] 601 <- (2,300)
605deffd5ebSksagiyam $
606deffd5ebSksagiyam 
607deffd5ebSksagiyam   Notes:
608deffd5ebSksagiyam   The layout parameter represents any partitioning of [0, N), where N is the total number of global indices, and its
609deffd5ebSksagiyam   local size can be set to PETSC_DECIDE.
610deffd5ebSksagiyam   If a global index x lies in the partition owned by process i, each process whose rootIndices contains x requests
611deffd5ebSksagiyam   ownership of x and sends its own rank and the local index of x to process i.
612deffd5ebSksagiyam   If multiple processes request ownership of x, the one with the highest rank is to own x.
613deffd5ebSksagiyam   Process i then broadcasts the ownership information, so that each process whose leafIndices contains x knows the
614deffd5ebSksagiyam   ownership information of x.
615deffd5ebSksagiyam   The output sf is constructed by associating each leaf point to a root point in this way.
616deffd5ebSksagiyam 
617deffd5ebSksagiyam   Suppose there is point data ordered according to the global indices and partitioned according to the given layout.
618deffd5ebSksagiyam   The optional output PetscSF object sfA can be used to push such data to leaf points.
619deffd5ebSksagiyam 
620deffd5ebSksagiyam   All indices in rootIndices and leafIndices must lie in the layout range. The union (over all processes) of rootIndices
621deffd5ebSksagiyam   must cover that of leafIndices, but need not cover the entire layout.
622deffd5ebSksagiyam 
623deffd5ebSksagiyam   If (leafIndices, leafLocalIndices, leafLocalOffset) == (rootIndices, rootLocalIndices, rootLocalOffset), the output
624deffd5ebSksagiyam   star forest is almost identity, so will only include non-trivial part of the map.
625deffd5ebSksagiyam 
626deffd5ebSksagiyam   Developer Notes:
627deffd5ebSksagiyam   Current approach of a process of the highest rank gaining the ownership may cause load imbalance; consider using
628deffd5ebSksagiyam   hash(rank, root_local_index) as the bid for the ownership determination.
629deffd5ebSksagiyam 
630deffd5ebSksagiyam   Level: advanced
631deffd5ebSksagiyam 
632deffd5ebSksagiyam .seealso: PetscSFCreate()
633deffd5ebSksagiyam @*/
634deffd5ebSksagiyam PetscErrorCode PetscSFCreateByMatchingIndices(PetscLayout layout, PetscInt numRootIndices, const PetscInt *rootIndices, const PetscInt *rootLocalIndices, PetscInt rootLocalOffset, PetscInt numLeafIndices, const PetscInt *leafIndices, const PetscInt *leafLocalIndices, PetscInt leafLocalOffset, PetscSF *sfA, PetscSF *sf)
635deffd5ebSksagiyam {
636deffd5ebSksagiyam   MPI_Comm        comm = layout->comm;
637deffd5ebSksagiyam   PetscMPIInt     size, rank;
638deffd5ebSksagiyam   PetscSF         sf1;
639deffd5ebSksagiyam   PetscSFNode    *owners, *buffer, *iremote;
640deffd5ebSksagiyam   PetscInt       *ilocal, nleaves, N, n, i;
641deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG)
642deffd5ebSksagiyam   PetscInt        N1;
643deffd5ebSksagiyam #endif
644deffd5ebSksagiyam   PetscBool       flag;
645deffd5ebSksagiyam 
646deffd5ebSksagiyam   PetscFunctionBegin;
647deffd5ebSksagiyam   if (rootIndices)      PetscValidIntPointer(rootIndices,3);
648deffd5ebSksagiyam   if (rootLocalIndices) PetscValidIntPointer(rootLocalIndices,4);
649deffd5ebSksagiyam   if (leafIndices)      PetscValidIntPointer(leafIndices,7);
650deffd5ebSksagiyam   if (leafLocalIndices) PetscValidIntPointer(leafLocalIndices,8);
651deffd5ebSksagiyam   if (sfA)              PetscValidPointer(sfA,10);
652deffd5ebSksagiyam   PetscValidPointer(sf,11);
6532c71b3e2SJacob Faibussowitsch   PetscCheckFalse(numRootIndices < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numRootIndices (%" PetscInt_FMT ") must be non-negative", numRootIndices);
6542c71b3e2SJacob Faibussowitsch   PetscCheckFalse(numLeafIndices < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "numLeafIndices (%" PetscInt_FMT ") must be non-negative", numLeafIndices);
6555f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm, &size));
6565f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
6575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetUp(layout));
6585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetSize(layout, &N));
6595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetLocalSize(layout, &n));
660deffd5ebSksagiyam   flag = (PetscBool)(leafIndices == rootIndices);
6615f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPIU_Allreduce(MPI_IN_PLACE, &flag, 1, MPIU_BOOL, MPI_LAND, comm));
6622c71b3e2SJacob Faibussowitsch   PetscCheckFalse(flag && numLeafIndices != numRootIndices,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "leafIndices == rootIndices, but numLeafIndices (%" PetscInt_FMT ") != numRootIndices(%" PetscInt_FMT ")", numLeafIndices, numRootIndices);
663deffd5ebSksagiyam #if defined(PETSC_USE_DEBUG)
664deffd5ebSksagiyam   N1 = PETSC_MIN_INT;
665deffd5ebSksagiyam   for (i = 0; i < numRootIndices; i++) if (rootIndices[i] > N1) N1 = rootIndices[i];
6665f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm));
6672c71b3e2SJacob Faibussowitsch   PetscCheckFalse(N1 >= N,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. root index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N);
668deffd5ebSksagiyam   if (!flag) {
669deffd5ebSksagiyam     N1 = PETSC_MIN_INT;
670deffd5ebSksagiyam     for (i = 0; i < numLeafIndices; i++) if (leafIndices[i] > N1) N1 = leafIndices[i];
6715f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Allreduce(MPI_IN_PLACE, &N1, 1, MPIU_INT, MPI_MAX, comm));
6722c71b3e2SJacob Faibussowitsch     PetscCheckFalse(N1 >= N,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Max. leaf index (%" PetscInt_FMT ") out of layout range [0,%" PetscInt_FMT ")", N1, N);
673deffd5ebSksagiyam   }
674deffd5ebSksagiyam #endif
675deffd5ebSksagiyam   /* Reduce: owners -> buffer */
6765f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n, &buffer));
6775f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreate(comm, &sf1));
6785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetFromOptions(sf1));
6795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetGraphLayout(sf1, layout, numRootIndices, NULL, PETSC_OWN_POINTER, rootIndices));
6805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(numRootIndices, &owners));
681deffd5ebSksagiyam   for (i = 0; i < numRootIndices; ++i) {
682deffd5ebSksagiyam     owners[i].rank = rank;
683deffd5ebSksagiyam     owners[i].index = rootLocalOffset + (rootLocalIndices ? rootLocalIndices[i] : i);
684deffd5ebSksagiyam   }
685deffd5ebSksagiyam   for (i = 0; i < n; ++i) {
686deffd5ebSksagiyam     buffer[i].index = -1;
687deffd5ebSksagiyam     buffer[i].rank = -1;
688deffd5ebSksagiyam   }
6895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFReduceBegin(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC));
6905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFReduceEnd(sf1, MPIU_2INT, owners, buffer, MPI_MAXLOC));
691deffd5ebSksagiyam   /* Bcast: buffer -> owners */
692deffd5ebSksagiyam   if (!flag) {
693deffd5ebSksagiyam     /* leafIndices is different from rootIndices */
6945f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(owners));
6955f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFSetGraphLayout(sf1, layout, numLeafIndices, NULL, PETSC_OWN_POINTER, leafIndices));
6965f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(numLeafIndices, &owners));
697deffd5ebSksagiyam   }
6985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastBegin(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE));
6995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastEnd(sf1, MPIU_2INT, buffer, owners, MPI_REPLACE));
7002c71b3e2SJacob Faibussowitsch   for (i = 0; i < numLeafIndices; ++i) PetscCheckFalse(owners[i].rank < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global point %" PetscInt_FMT " was unclaimed", leafIndices[i]);
7015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(buffer));
702deffd5ebSksagiyam   if (sfA) {*sfA = sf1;}
7035f80ce2aSJacob Faibussowitsch   else     CHKERRQ(PetscSFDestroy(&sf1));
704deffd5ebSksagiyam   /* Create sf */
705deffd5ebSksagiyam   if (flag && rootLocalIndices == leafLocalIndices && leafLocalOffset == rootLocalOffset) {
706deffd5ebSksagiyam     /* leaf space == root space */
707deffd5ebSksagiyam     for (i = 0, nleaves = 0; i < numLeafIndices; ++i) if (owners[i].rank != rank) ++nleaves;
7085f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nleaves, &ilocal));
7095f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nleaves, &iremote));
710deffd5ebSksagiyam     for (i = 0, nleaves = 0; i < numLeafIndices; ++i) {
711deffd5ebSksagiyam       if (owners[i].rank != rank) {
712deffd5ebSksagiyam         ilocal[nleaves]        = leafLocalOffset + i;
713deffd5ebSksagiyam         iremote[nleaves].rank  = owners[i].rank;
714deffd5ebSksagiyam         iremote[nleaves].index = owners[i].index;
715deffd5ebSksagiyam         ++nleaves;
716deffd5ebSksagiyam       }
717deffd5ebSksagiyam     }
7185f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(owners));
719deffd5ebSksagiyam   } else {
720deffd5ebSksagiyam     nleaves = numLeafIndices;
7215f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(nleaves, &ilocal));
722deffd5ebSksagiyam     for (i = 0; i < nleaves; ++i) {ilocal[i] = leafLocalOffset + (leafLocalIndices ? leafLocalIndices[i] : i);}
723deffd5ebSksagiyam     iremote = owners;
724deffd5ebSksagiyam   }
7255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreate(comm, sf));
7265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetFromOptions(*sf));
7275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFSetGraph(*sf, rootLocalOffset + numRootIndices, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
728deffd5ebSksagiyam   PetscFunctionReturn(0);
729deffd5ebSksagiyam }
730