xref: /petsc/src/vec/is/sf/utils/sfutils.c (revision b0c7db2251b14c2bf1fcf12ecb441159a8187d90)
1*b0c7db22SLisandro Dalcin #include <petsc/private/sfimpl.h>       /*I  "petscsf.h"   I*/
2*b0c7db22SLisandro Dalcin #include <petsc/private/sectionimpl.h>
3*b0c7db22SLisandro Dalcin 
4*b0c7db22SLisandro Dalcin /*@C
5*b0c7db22SLisandro Dalcin    PetscSFSetGraphLayout - Set a parallel star forest via global indices and a PetscLayout
6*b0c7db22SLisandro Dalcin 
7*b0c7db22SLisandro Dalcin    Collective
8*b0c7db22SLisandro Dalcin 
9*b0c7db22SLisandro Dalcin    Input Arguments:
10*b0c7db22SLisandro Dalcin +  sf - star forest
11*b0c7db22SLisandro Dalcin .  layout - PetscLayout defining the global space
12*b0c7db22SLisandro Dalcin .  nleaves - number of leaf vertices on the current process, each of these references a root on any process
13*b0c7db22SLisandro Dalcin .  ilocal - locations of leaves in leafdata buffers, pass NULL for contiguous storage
14*b0c7db22SLisandro Dalcin .  localmode - copy mode for ilocal
15*b0c7db22SLisandro Dalcin -  iremote - remote locations of root vertices for each leaf on the current process
16*b0c7db22SLisandro Dalcin 
17*b0c7db22SLisandro Dalcin    Level: intermediate
18*b0c7db22SLisandro Dalcin 
19*b0c7db22SLisandro Dalcin    Developers Note: Local indices which are the identity permutation in the range [0,nleaves) are discarded as they
20*b0c7db22SLisandro Dalcin    encode contiguous storage. In such case, if localmode is PETSC_OWN_POINTER, the memory is deallocated as it is not
21*b0c7db22SLisandro Dalcin    needed
22*b0c7db22SLisandro Dalcin 
23*b0c7db22SLisandro Dalcin .seealso: PetscSFCreate(), PetscSFView(), PetscSFSetGraph(), PetscSFGetGraph()
24*b0c7db22SLisandro Dalcin @*/
25*b0c7db22SLisandro Dalcin PetscErrorCode PetscSFSetGraphLayout(PetscSF sf,PetscLayout layout,PetscInt nleaves,const PetscInt *ilocal,PetscCopyMode localmode,const PetscInt *iremote)
26*b0c7db22SLisandro Dalcin {
27*b0c7db22SLisandro Dalcin   PetscErrorCode ierr;
28*b0c7db22SLisandro Dalcin   PetscInt       i,nroots;
29*b0c7db22SLisandro Dalcin   PetscSFNode    *remote;
30*b0c7db22SLisandro Dalcin 
31*b0c7db22SLisandro Dalcin   PetscFunctionBegin;
32*b0c7db22SLisandro Dalcin   ierr = PetscLayoutGetLocalSize(layout,&nroots);CHKERRQ(ierr);
33*b0c7db22SLisandro Dalcin   ierr = PetscMalloc1(nleaves,&remote);CHKERRQ(ierr);
34*b0c7db22SLisandro Dalcin   for (i=0; i<nleaves; i++) {
35*b0c7db22SLisandro Dalcin     PetscMPIInt owner = -1;
36*b0c7db22SLisandro Dalcin     ierr = PetscLayoutFindOwner(layout,iremote[i],&owner);CHKERRQ(ierr);
37*b0c7db22SLisandro Dalcin     remote[i].rank  = owner;
38*b0c7db22SLisandro Dalcin     remote[i].index = iremote[i] - layout->range[owner];
39*b0c7db22SLisandro Dalcin   }
40*b0c7db22SLisandro Dalcin   ierr = PetscSFSetGraph(sf,nroots,nleaves,ilocal,localmode,remote,PETSC_OWN_POINTER);CHKERRQ(ierr);
41*b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
42*b0c7db22SLisandro Dalcin }
43*b0c7db22SLisandro Dalcin 
44*b0c7db22SLisandro Dalcin /*@
45*b0c7db22SLisandro Dalcin   PetscSFSetGraphSection - Sets the PetscSF graph encoding the parallel dof overlap based upon the PetscSections
46*b0c7db22SLisandro Dalcin   describing the data layout.
47*b0c7db22SLisandro Dalcin 
48*b0c7db22SLisandro Dalcin   Input Parameters:
49*b0c7db22SLisandro Dalcin + sf - The SF
50*b0c7db22SLisandro Dalcin . localSection - PetscSection describing the local data layout
51*b0c7db22SLisandro Dalcin - globalSection - PetscSection describing the global data layout
52*b0c7db22SLisandro Dalcin 
53*b0c7db22SLisandro Dalcin   Level: developer
54*b0c7db22SLisandro Dalcin 
55*b0c7db22SLisandro Dalcin .seealso: PetscSFSetGraph(), PetscSFSetGraphLayout()
56*b0c7db22SLisandro Dalcin @*/
57*b0c7db22SLisandro Dalcin PetscErrorCode PetscSFSetGraphSection(PetscSF sf, PetscSection localSection, PetscSection globalSection)
58*b0c7db22SLisandro Dalcin {
59*b0c7db22SLisandro Dalcin   MPI_Comm       comm;
60*b0c7db22SLisandro Dalcin   PetscLayout    layout;
61*b0c7db22SLisandro Dalcin   const PetscInt *ranges;
62*b0c7db22SLisandro Dalcin   PetscInt       *local;
63*b0c7db22SLisandro Dalcin   PetscSFNode    *remote;
64*b0c7db22SLisandro Dalcin   PetscInt       pStart, pEnd, p, nroots, nleaves = 0, l;
65*b0c7db22SLisandro Dalcin   PetscMPIInt    size, rank;
66*b0c7db22SLisandro Dalcin   PetscErrorCode ierr;
67*b0c7db22SLisandro Dalcin 
68*b0c7db22SLisandro Dalcin   PetscFunctionBegin;
69*b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1);
70*b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(localSection, PETSC_SECTION_CLASSID, 2);
71*b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(globalSection, PETSC_SECTION_CLASSID, 3);
72*b0c7db22SLisandro Dalcin 
73*b0c7db22SLisandro Dalcin   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
74*b0c7db22SLisandro Dalcin   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
75*b0c7db22SLisandro Dalcin   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
76*b0c7db22SLisandro Dalcin   ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr);
77*b0c7db22SLisandro Dalcin   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr);
78*b0c7db22SLisandro Dalcin   ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
79*b0c7db22SLisandro Dalcin   ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
80*b0c7db22SLisandro Dalcin   ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr);
81*b0c7db22SLisandro Dalcin   ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
82*b0c7db22SLisandro Dalcin   ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
83*b0c7db22SLisandro Dalcin   for (p = pStart; p < pEnd; ++p) {
84*b0c7db22SLisandro Dalcin     PetscInt gdof, gcdof;
85*b0c7db22SLisandro Dalcin 
86*b0c7db22SLisandro Dalcin     ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr);
87*b0c7db22SLisandro Dalcin     ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr);
88*b0c7db22SLisandro Dalcin     if (gcdof > (gdof < 0 ? -(gdof+1) : gdof)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D has %D constraints > %D dof", p, gcdof, (gdof < 0 ? -(gdof+1) : gdof));
89*b0c7db22SLisandro Dalcin     nleaves += gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
90*b0c7db22SLisandro Dalcin   }
91*b0c7db22SLisandro Dalcin   ierr = PetscMalloc1(nleaves, &local);CHKERRQ(ierr);
92*b0c7db22SLisandro Dalcin   ierr = PetscMalloc1(nleaves, &remote);CHKERRQ(ierr);
93*b0c7db22SLisandro Dalcin   for (p = pStart, l = 0; p < pEnd; ++p) {
94*b0c7db22SLisandro Dalcin     const PetscInt *cind;
95*b0c7db22SLisandro Dalcin     PetscInt       dof, cdof, off, gdof, gcdof, goff, gsize, d, c;
96*b0c7db22SLisandro Dalcin 
97*b0c7db22SLisandro Dalcin     ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr);
98*b0c7db22SLisandro Dalcin     ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr);
99*b0c7db22SLisandro Dalcin     ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr);
100*b0c7db22SLisandro Dalcin     ierr = PetscSectionGetConstraintIndices(localSection, p, &cind);CHKERRQ(ierr);
101*b0c7db22SLisandro Dalcin     ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr);
102*b0c7db22SLisandro Dalcin     ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr);
103*b0c7db22SLisandro Dalcin     ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr);
104*b0c7db22SLisandro Dalcin     if (!gdof) continue; /* Censored point */
105*b0c7db22SLisandro Dalcin     gsize = gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof;
106*b0c7db22SLisandro Dalcin     if (gsize != dof-cdof) {
107*b0c7db22SLisandro Dalcin       if (gsize != dof) SETERRQ4(comm, PETSC_ERR_ARG_WRONG, "Global dof %D for point %D is neither the constrained size %D, nor the unconstrained %D", gsize, p, dof-cdof, dof);
108*b0c7db22SLisandro Dalcin       cdof = 0; /* Ignore constraints */
109*b0c7db22SLisandro Dalcin     }
110*b0c7db22SLisandro Dalcin     for (d = 0, c = 0; d < dof; ++d) {
111*b0c7db22SLisandro Dalcin       if ((c < cdof) && (cind[c] == d)) {++c; continue;}
112*b0c7db22SLisandro Dalcin       local[l+d-c] = off+d;
113*b0c7db22SLisandro Dalcin     }
114*b0c7db22SLisandro Dalcin     if (gdof < 0) {
115*b0c7db22SLisandro Dalcin       for (d = 0; d < gsize; ++d, ++l) {
116*b0c7db22SLisandro Dalcin         PetscInt offset = -(goff+1) + d, r;
117*b0c7db22SLisandro Dalcin 
118*b0c7db22SLisandro Dalcin         ierr = PetscFindInt(offset,size+1,ranges,&r);CHKERRQ(ierr);
119*b0c7db22SLisandro Dalcin         if (r < 0) r = -(r+2);
120*b0c7db22SLisandro Dalcin         if ((r < 0) || (r >= size)) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D mapped to invalid process %D (%D, %D)", p, r, gdof, goff);
121*b0c7db22SLisandro Dalcin         remote[l].rank  = r;
122*b0c7db22SLisandro Dalcin         remote[l].index = offset - ranges[r];
123*b0c7db22SLisandro Dalcin       }
124*b0c7db22SLisandro Dalcin     } else {
125*b0c7db22SLisandro Dalcin       for (d = 0; d < gsize; ++d, ++l) {
126*b0c7db22SLisandro Dalcin         remote[l].rank  = rank;
127*b0c7db22SLisandro Dalcin         remote[l].index = goff+d - ranges[rank];
128*b0c7db22SLisandro Dalcin       }
129*b0c7db22SLisandro Dalcin     }
130*b0c7db22SLisandro Dalcin   }
131*b0c7db22SLisandro Dalcin   if (l != nleaves) SETERRQ2(comm, PETSC_ERR_PLIB, "Iteration error, l %D != nleaves %D", l, nleaves);
132*b0c7db22SLisandro Dalcin   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
133*b0c7db22SLisandro Dalcin   ierr = PetscSFSetGraph(sf, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr);
134*b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
135*b0c7db22SLisandro Dalcin }
136*b0c7db22SLisandro Dalcin 
137*b0c7db22SLisandro Dalcin static PetscErrorCode PetscSectionCheckConstraints_Static(PetscSection s)
138*b0c7db22SLisandro Dalcin {
139*b0c7db22SLisandro Dalcin   PetscErrorCode ierr;
140*b0c7db22SLisandro Dalcin 
141*b0c7db22SLisandro Dalcin   PetscFunctionBegin;
142*b0c7db22SLisandro Dalcin   if (!s->bc) {
143*b0c7db22SLisandro Dalcin     ierr = PetscSectionCreate(PETSC_COMM_SELF, &s->bc);CHKERRQ(ierr);
144*b0c7db22SLisandro Dalcin     ierr = PetscSectionSetChart(s->bc, s->pStart, s->pEnd);CHKERRQ(ierr);
145*b0c7db22SLisandro Dalcin   }
146*b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
147*b0c7db22SLisandro Dalcin }
148*b0c7db22SLisandro Dalcin 
149*b0c7db22SLisandro Dalcin /*@C
150*b0c7db22SLisandro Dalcin   PetscSFDistributeSection - Create a new PetscSection reorganized, moving from the root to the leaves of the SF
151*b0c7db22SLisandro Dalcin 
152*b0c7db22SLisandro Dalcin   Collective on sf
153*b0c7db22SLisandro Dalcin 
154*b0c7db22SLisandro Dalcin   Input Parameters:
155*b0c7db22SLisandro Dalcin + sf - The SF
156*b0c7db22SLisandro Dalcin - rootSection - Section defined on root space
157*b0c7db22SLisandro Dalcin 
158*b0c7db22SLisandro Dalcin   Output Parameters:
159*b0c7db22SLisandro Dalcin + remoteOffsets - root offsets in leaf storage, or NULL
160*b0c7db22SLisandro Dalcin - leafSection - Section defined on the leaf space
161*b0c7db22SLisandro Dalcin 
162*b0c7db22SLisandro Dalcin   Level: advanced
163*b0c7db22SLisandro Dalcin 
164*b0c7db22SLisandro Dalcin .seealso: PetscSFCreate()
165*b0c7db22SLisandro Dalcin @*/
166*b0c7db22SLisandro Dalcin PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection)
167*b0c7db22SLisandro Dalcin {
168*b0c7db22SLisandro Dalcin   PetscSF        embedSF;
169*b0c7db22SLisandro Dalcin   const PetscInt *indices;
170*b0c7db22SLisandro Dalcin   IS             selected;
171*b0c7db22SLisandro Dalcin   PetscInt       numFields, nroots, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, f, c;
172*b0c7db22SLisandro Dalcin   PetscBool      *sub, hasc;
173*b0c7db22SLisandro Dalcin   PetscErrorCode ierr;
174*b0c7db22SLisandro Dalcin 
175*b0c7db22SLisandro Dalcin   PetscFunctionBegin;
176*b0c7db22SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_DistSect,sf,0,0,0);CHKERRQ(ierr);
177*b0c7db22SLisandro Dalcin   ierr = PetscSectionGetNumFields(rootSection, &numFields);CHKERRQ(ierr);
178*b0c7db22SLisandro Dalcin   if (numFields) {ierr = PetscSectionSetNumFields(leafSection, numFields);CHKERRQ(ierr);}
179*b0c7db22SLisandro Dalcin   ierr = PetscMalloc1(numFields+2, &sub);CHKERRQ(ierr);
180*b0c7db22SLisandro Dalcin   sub[1] = rootSection->bc ? PETSC_TRUE : PETSC_FALSE;
181*b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) {
182*b0c7db22SLisandro Dalcin     PetscSectionSym sym;
183*b0c7db22SLisandro Dalcin     const char      *name   = NULL;
184*b0c7db22SLisandro Dalcin     PetscInt        numComp = 0;
185*b0c7db22SLisandro Dalcin 
186*b0c7db22SLisandro Dalcin     sub[2 + f] = rootSection->field[f]->bc ? PETSC_TRUE : PETSC_FALSE;
187*b0c7db22SLisandro Dalcin     ierr = PetscSectionGetFieldComponents(rootSection, f, &numComp);CHKERRQ(ierr);
188*b0c7db22SLisandro Dalcin     ierr = PetscSectionGetFieldName(rootSection, f, &name);CHKERRQ(ierr);
189*b0c7db22SLisandro Dalcin     ierr = PetscSectionGetFieldSym(rootSection, f, &sym);CHKERRQ(ierr);
190*b0c7db22SLisandro Dalcin     ierr = PetscSectionSetFieldComponents(leafSection, f, numComp);CHKERRQ(ierr);
191*b0c7db22SLisandro Dalcin     ierr = PetscSectionSetFieldName(leafSection, f, name);CHKERRQ(ierr);
192*b0c7db22SLisandro Dalcin     ierr = PetscSectionSetFieldSym(leafSection, f, sym);CHKERRQ(ierr);
193*b0c7db22SLisandro Dalcin     for (c = 0; c < rootSection->numFieldComponents[f]; ++c) {
194*b0c7db22SLisandro Dalcin       ierr = PetscSectionGetComponentName(rootSection, f, c, &name);CHKERRQ(ierr);
195*b0c7db22SLisandro Dalcin       ierr = PetscSectionSetComponentName(leafSection, f, c, name);CHKERRQ(ierr);
196*b0c7db22SLisandro Dalcin     }
197*b0c7db22SLisandro Dalcin   }
198*b0c7db22SLisandro Dalcin   ierr = PetscSectionGetChart(rootSection, &rpStart, &rpEnd);CHKERRQ(ierr);
199*b0c7db22SLisandro Dalcin   ierr = PetscSFGetGraph(sf,&nroots,NULL,NULL,NULL);CHKERRQ(ierr);
200*b0c7db22SLisandro Dalcin   rpEnd = PetscMin(rpEnd,nroots);
201*b0c7db22SLisandro Dalcin   rpEnd = PetscMax(rpStart,rpEnd);
202*b0c7db22SLisandro Dalcin   /* see if we can avoid creating the embedded SF, since it can cost more than an allreduce */
203*b0c7db22SLisandro Dalcin   sub[0] = (PetscBool)(nroots != rpEnd - rpStart);
204*b0c7db22SLisandro Dalcin   ierr = MPIU_Allreduce(MPI_IN_PLACE, sub, 2+numFields, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)sf));CHKERRQ(ierr);
205*b0c7db22SLisandro Dalcin   if (sub[0]) {
206*b0c7db22SLisandro Dalcin     ierr = ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected);CHKERRQ(ierr);
207*b0c7db22SLisandro Dalcin     ierr = ISGetIndices(selected, &indices);CHKERRQ(ierr);
208*b0c7db22SLisandro Dalcin     ierr = PetscSFCreateEmbeddedSF(sf, rpEnd - rpStart, indices, &embedSF);CHKERRQ(ierr);
209*b0c7db22SLisandro Dalcin     ierr = ISRestoreIndices(selected, &indices);CHKERRQ(ierr);
210*b0c7db22SLisandro Dalcin     ierr = ISDestroy(&selected);CHKERRQ(ierr);
211*b0c7db22SLisandro Dalcin   } else {
212*b0c7db22SLisandro Dalcin     ierr = PetscObjectReference((PetscObject)sf);CHKERRQ(ierr);
213*b0c7db22SLisandro Dalcin     embedSF = sf;
214*b0c7db22SLisandro Dalcin   }
215*b0c7db22SLisandro Dalcin   ierr = PetscSFGetLeafRange(embedSF, &lpStart, &lpEnd);CHKERRQ(ierr);
216*b0c7db22SLisandro Dalcin   lpEnd++;
217*b0c7db22SLisandro Dalcin 
218*b0c7db22SLisandro Dalcin   ierr = PetscSectionSetChart(leafSection, lpStart, lpEnd);CHKERRQ(ierr);
219*b0c7db22SLisandro Dalcin 
220*b0c7db22SLisandro Dalcin   /* Constrained dof section */
221*b0c7db22SLisandro Dalcin   hasc = sub[1];
222*b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) hasc = (PetscBool)(hasc || sub[2+f]);
223*b0c7db22SLisandro Dalcin 
224*b0c7db22SLisandro Dalcin   /* Could fuse these at the cost of copies and extra allocation */
225*b0c7db22SLisandro Dalcin   ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart]);CHKERRQ(ierr);
226*b0c7db22SLisandro Dalcin   ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart]);CHKERRQ(ierr);
227*b0c7db22SLisandro Dalcin   if (sub[1]) {
228*b0c7db22SLisandro Dalcin     ierr = PetscSectionCheckConstraints_Static(rootSection);CHKERRQ(ierr);
229*b0c7db22SLisandro Dalcin     ierr = PetscSectionCheckConstraints_Static(leafSection);CHKERRQ(ierr);
230*b0c7db22SLisandro Dalcin     ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart]);CHKERRQ(ierr);
231*b0c7db22SLisandro Dalcin     ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasDof[-rpStart], &leafSection->bc->atlasDof[-lpStart]);CHKERRQ(ierr);
232*b0c7db22SLisandro Dalcin   }
233*b0c7db22SLisandro Dalcin   for (f = 0; f < numFields; ++f) {
234*b0c7db22SLisandro Dalcin     ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart]);CHKERRQ(ierr);
235*b0c7db22SLisandro Dalcin     ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart]);CHKERRQ(ierr);
236*b0c7db22SLisandro Dalcin     if (sub[2+f]) {
237*b0c7db22SLisandro Dalcin       ierr = PetscSectionCheckConstraints_Static(rootSection->field[f]);CHKERRQ(ierr);
238*b0c7db22SLisandro Dalcin       ierr = PetscSectionCheckConstraints_Static(leafSection->field[f]);CHKERRQ(ierr);
239*b0c7db22SLisandro Dalcin       ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart]);CHKERRQ(ierr);
240*b0c7db22SLisandro Dalcin       ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasDof[-rpStart], &leafSection->field[f]->bc->atlasDof[-lpStart]);CHKERRQ(ierr);
241*b0c7db22SLisandro Dalcin     }
242*b0c7db22SLisandro Dalcin   }
243*b0c7db22SLisandro Dalcin   if (remoteOffsets) {
244*b0c7db22SLisandro Dalcin     ierr = PetscMalloc1(lpEnd - lpStart, remoteOffsets);CHKERRQ(ierr);
245*b0c7db22SLisandro Dalcin     ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart]);CHKERRQ(ierr);
246*b0c7db22SLisandro Dalcin     ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart]);CHKERRQ(ierr);
247*b0c7db22SLisandro Dalcin   }
248*b0c7db22SLisandro Dalcin   ierr = PetscSectionSetUp(leafSection);CHKERRQ(ierr);
249*b0c7db22SLisandro Dalcin   if (hasc) { /* need to communicate bcIndices */
250*b0c7db22SLisandro Dalcin     PetscSF  bcSF;
251*b0c7db22SLisandro Dalcin     PetscInt *rOffBc;
252*b0c7db22SLisandro Dalcin 
253*b0c7db22SLisandro Dalcin     ierr = PetscMalloc1(lpEnd - lpStart, &rOffBc);CHKERRQ(ierr);
254*b0c7db22SLisandro Dalcin     if (sub[1]) {
255*b0c7db22SLisandro Dalcin       ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart]);CHKERRQ(ierr);
256*b0c7db22SLisandro Dalcin       ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->bc->atlasOff[-rpStart], &rOffBc[-lpStart]);CHKERRQ(ierr);
257*b0c7db22SLisandro Dalcin       ierr = PetscSFCreateSectionSF(embedSF, rootSection->bc, rOffBc, leafSection->bc, &bcSF);CHKERRQ(ierr);
258*b0c7db22SLisandro Dalcin       ierr = PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices);CHKERRQ(ierr);
259*b0c7db22SLisandro Dalcin       ierr = PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->bcIndices, leafSection->bcIndices);CHKERRQ(ierr);
260*b0c7db22SLisandro Dalcin       ierr = PetscSFDestroy(&bcSF);CHKERRQ(ierr);
261*b0c7db22SLisandro Dalcin     }
262*b0c7db22SLisandro Dalcin     for (f = 0; f < numFields; ++f) {
263*b0c7db22SLisandro Dalcin       if (sub[2+f]) {
264*b0c7db22SLisandro Dalcin         ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart]);CHKERRQ(ierr);
265*b0c7db22SLisandro Dalcin         ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->bc->atlasOff[-rpStart], &rOffBc[-lpStart]);CHKERRQ(ierr);
266*b0c7db22SLisandro Dalcin         ierr = PetscSFCreateSectionSF(embedSF, rootSection->field[f]->bc, rOffBc, leafSection->field[f]->bc, &bcSF);CHKERRQ(ierr);
267*b0c7db22SLisandro Dalcin         ierr = PetscSFBcastBegin(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices);CHKERRQ(ierr);
268*b0c7db22SLisandro Dalcin         ierr = PetscSFBcastEnd(bcSF, MPIU_INT, rootSection->field[f]->bcIndices, leafSection->field[f]->bcIndices);CHKERRQ(ierr);
269*b0c7db22SLisandro Dalcin         ierr = PetscSFDestroy(&bcSF);CHKERRQ(ierr);
270*b0c7db22SLisandro Dalcin       }
271*b0c7db22SLisandro Dalcin     }
272*b0c7db22SLisandro Dalcin     ierr = PetscFree(rOffBc);CHKERRQ(ierr);
273*b0c7db22SLisandro Dalcin   }
274*b0c7db22SLisandro Dalcin   ierr = PetscSFDestroy(&embedSF);CHKERRQ(ierr);
275*b0c7db22SLisandro Dalcin   ierr = PetscFree(sub);CHKERRQ(ierr);
276*b0c7db22SLisandro Dalcin   ierr = PetscLogEventEnd(PETSCSF_DistSect,sf,0,0,0);CHKERRQ(ierr);
277*b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
278*b0c7db22SLisandro Dalcin }
279*b0c7db22SLisandro Dalcin 
280*b0c7db22SLisandro Dalcin /*@C
281*b0c7db22SLisandro Dalcin   PetscSFCreateRemoteOffsets - Create offsets for point data on remote processes
282*b0c7db22SLisandro Dalcin 
283*b0c7db22SLisandro Dalcin   Collective on sf
284*b0c7db22SLisandro Dalcin 
285*b0c7db22SLisandro Dalcin   Input Parameters:
286*b0c7db22SLisandro Dalcin + sf          - The SF
287*b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is layout for SF roots)
288*b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data  (this is layout for SF leaves)
289*b0c7db22SLisandro Dalcin 
290*b0c7db22SLisandro Dalcin   Output Parameter:
291*b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
292*b0c7db22SLisandro Dalcin 
293*b0c7db22SLisandro Dalcin   Level: developer
294*b0c7db22SLisandro Dalcin 
295*b0c7db22SLisandro Dalcin .seealso: PetscSFCreate()
296*b0c7db22SLisandro Dalcin @*/
297*b0c7db22SLisandro Dalcin PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets)
298*b0c7db22SLisandro Dalcin {
299*b0c7db22SLisandro Dalcin   PetscSF         embedSF;
300*b0c7db22SLisandro Dalcin   const PetscInt *indices;
301*b0c7db22SLisandro Dalcin   IS              selected;
302*b0c7db22SLisandro Dalcin   PetscInt        numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0;
303*b0c7db22SLisandro Dalcin   PetscErrorCode  ierr;
304*b0c7db22SLisandro Dalcin 
305*b0c7db22SLisandro Dalcin   PetscFunctionBegin;
306*b0c7db22SLisandro Dalcin   *remoteOffsets = NULL;
307*b0c7db22SLisandro Dalcin   ierr = PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL);CHKERRQ(ierr);
308*b0c7db22SLisandro Dalcin   if (numRoots < 0) PetscFunctionReturn(0);
309*b0c7db22SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_RemoteOff,sf,0,0,0);CHKERRQ(ierr);
310*b0c7db22SLisandro Dalcin   ierr = PetscSectionGetChart(rootSection, &rpStart, &rpEnd);CHKERRQ(ierr);
311*b0c7db22SLisandro Dalcin   ierr = PetscSectionGetChart(leafSection, &lpStart, &lpEnd);CHKERRQ(ierr);
312*b0c7db22SLisandro Dalcin   ierr = ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected);CHKERRQ(ierr);
313*b0c7db22SLisandro Dalcin   ierr = ISGetIndices(selected, &indices);CHKERRQ(ierr);
314*b0c7db22SLisandro Dalcin   ierr = PetscSFCreateEmbeddedSF(sf, rpEnd - rpStart, indices, &embedSF);CHKERRQ(ierr);
315*b0c7db22SLisandro Dalcin   ierr = ISRestoreIndices(selected, &indices);CHKERRQ(ierr);
316*b0c7db22SLisandro Dalcin   ierr = ISDestroy(&selected);CHKERRQ(ierr);
317*b0c7db22SLisandro Dalcin   ierr = PetscCalloc1(lpEnd - lpStart, remoteOffsets);CHKERRQ(ierr);
318*b0c7db22SLisandro Dalcin   ierr = PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart]);CHKERRQ(ierr);
319*b0c7db22SLisandro Dalcin   ierr = PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart]);CHKERRQ(ierr);
320*b0c7db22SLisandro Dalcin   ierr = PetscSFDestroy(&embedSF);CHKERRQ(ierr);
321*b0c7db22SLisandro Dalcin   ierr = PetscLogEventEnd(PETSCSF_RemoteOff,sf,0,0,0);CHKERRQ(ierr);
322*b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
323*b0c7db22SLisandro Dalcin }
324*b0c7db22SLisandro Dalcin 
325*b0c7db22SLisandro Dalcin /*@C
326*b0c7db22SLisandro Dalcin   PetscSFCreateSectionSF - Create an expanded SF of dofs, assuming the input SF relates points
327*b0c7db22SLisandro Dalcin 
328*b0c7db22SLisandro Dalcin   Collective on sf
329*b0c7db22SLisandro Dalcin 
330*b0c7db22SLisandro Dalcin   Input Parameters:
331*b0c7db22SLisandro Dalcin + sf - The SF
332*b0c7db22SLisandro Dalcin . rootSection - Data layout of remote points for outgoing data (this is usually the serial section)
333*b0c7db22SLisandro Dalcin . remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
334*b0c7db22SLisandro Dalcin - leafSection - Data layout of local points for incoming data  (this is the distributed section)
335*b0c7db22SLisandro Dalcin 
336*b0c7db22SLisandro Dalcin   Output Parameters:
337*b0c7db22SLisandro Dalcin - sectionSF - The new SF
338*b0c7db22SLisandro Dalcin 
339*b0c7db22SLisandro Dalcin   Note: Either rootSection or remoteOffsets can be specified
340*b0c7db22SLisandro Dalcin 
341*b0c7db22SLisandro Dalcin   Level: advanced
342*b0c7db22SLisandro Dalcin 
343*b0c7db22SLisandro Dalcin .seealso: PetscSFCreate()
344*b0c7db22SLisandro Dalcin @*/
345*b0c7db22SLisandro Dalcin PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF)
346*b0c7db22SLisandro Dalcin {
347*b0c7db22SLisandro Dalcin   MPI_Comm          comm;
348*b0c7db22SLisandro Dalcin   const PetscInt    *localPoints;
349*b0c7db22SLisandro Dalcin   const PetscSFNode *remotePoints;
350*b0c7db22SLisandro Dalcin   PetscInt          lpStart, lpEnd;
351*b0c7db22SLisandro Dalcin   PetscInt          numRoots, numSectionRoots, numPoints, numIndices = 0;
352*b0c7db22SLisandro Dalcin   PetscInt          *localIndices;
353*b0c7db22SLisandro Dalcin   PetscSFNode       *remoteIndices;
354*b0c7db22SLisandro Dalcin   PetscInt          i, ind;
355*b0c7db22SLisandro Dalcin   PetscErrorCode    ierr;
356*b0c7db22SLisandro Dalcin 
357*b0c7db22SLisandro Dalcin   PetscFunctionBegin;
358*b0c7db22SLisandro Dalcin   PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1);
359*b0c7db22SLisandro Dalcin   PetscValidPointer(rootSection,2);
360*b0c7db22SLisandro Dalcin   /* Cannot check PetscValidIntPointer(remoteOffsets,3) because it can be NULL if sf does not reference any points in leafSection */
361*b0c7db22SLisandro Dalcin   PetscValidPointer(leafSection,4);
362*b0c7db22SLisandro Dalcin   PetscValidPointer(sectionSF,5);
363*b0c7db22SLisandro Dalcin   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
364*b0c7db22SLisandro Dalcin   ierr = PetscSFCreate(comm, sectionSF);CHKERRQ(ierr);
365*b0c7db22SLisandro Dalcin   ierr = PetscSectionGetChart(leafSection, &lpStart, &lpEnd);CHKERRQ(ierr);
366*b0c7db22SLisandro Dalcin   ierr = PetscSectionGetStorageSize(rootSection, &numSectionRoots);CHKERRQ(ierr);
367*b0c7db22SLisandro Dalcin   ierr = PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints);CHKERRQ(ierr);
368*b0c7db22SLisandro Dalcin   if (numRoots < 0) PetscFunctionReturn(0);
369*b0c7db22SLisandro Dalcin   ierr = PetscLogEventBegin(PETSCSF_SectSF,sf,0,0,0);CHKERRQ(ierr);
370*b0c7db22SLisandro Dalcin   for (i = 0; i < numPoints; ++i) {
371*b0c7db22SLisandro Dalcin     PetscInt localPoint = localPoints ? localPoints[i] : i;
372*b0c7db22SLisandro Dalcin     PetscInt dof;
373*b0c7db22SLisandro Dalcin 
374*b0c7db22SLisandro Dalcin     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
375*b0c7db22SLisandro Dalcin       ierr = PetscSectionGetDof(leafSection, localPoint, &dof);CHKERRQ(ierr);
376*b0c7db22SLisandro Dalcin       numIndices += dof;
377*b0c7db22SLisandro Dalcin     }
378*b0c7db22SLisandro Dalcin   }
379*b0c7db22SLisandro Dalcin   ierr = PetscMalloc1(numIndices, &localIndices);CHKERRQ(ierr);
380*b0c7db22SLisandro Dalcin   ierr = PetscMalloc1(numIndices, &remoteIndices);CHKERRQ(ierr);
381*b0c7db22SLisandro Dalcin   /* Create new index graph */
382*b0c7db22SLisandro Dalcin   for (i = 0, ind = 0; i < numPoints; ++i) {
383*b0c7db22SLisandro Dalcin     PetscInt localPoint = localPoints ? localPoints[i] : i;
384*b0c7db22SLisandro Dalcin     PetscInt rank       = remotePoints[i].rank;
385*b0c7db22SLisandro Dalcin 
386*b0c7db22SLisandro Dalcin     if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
387*b0c7db22SLisandro Dalcin       PetscInt remoteOffset = remoteOffsets[localPoint-lpStart];
388*b0c7db22SLisandro Dalcin       PetscInt loff, dof, d;
389*b0c7db22SLisandro Dalcin 
390*b0c7db22SLisandro Dalcin       ierr = PetscSectionGetOffset(leafSection, localPoint, &loff);CHKERRQ(ierr);
391*b0c7db22SLisandro Dalcin       ierr = PetscSectionGetDof(leafSection, localPoint, &dof);CHKERRQ(ierr);
392*b0c7db22SLisandro Dalcin       for (d = 0; d < dof; ++d, ++ind) {
393*b0c7db22SLisandro Dalcin         localIndices[ind]        = loff+d;
394*b0c7db22SLisandro Dalcin         remoteIndices[ind].rank  = rank;
395*b0c7db22SLisandro Dalcin         remoteIndices[ind].index = remoteOffset+d;
396*b0c7db22SLisandro Dalcin       }
397*b0c7db22SLisandro Dalcin     }
398*b0c7db22SLisandro Dalcin   }
399*b0c7db22SLisandro Dalcin   if (numIndices != ind) SETERRQ2(comm, PETSC_ERR_PLIB, "Inconsistency in indices, %D should be %D", ind, numIndices);
400*b0c7db22SLisandro Dalcin   ierr = PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER);CHKERRQ(ierr);
401*b0c7db22SLisandro Dalcin   ierr = PetscSFSetUp(*sectionSF);CHKERRQ(ierr);
402*b0c7db22SLisandro Dalcin   ierr = PetscLogEventEnd(PETSCSF_SectSF,sf,0,0,0);CHKERRQ(ierr);
403*b0c7db22SLisandro Dalcin   PetscFunctionReturn(0);
404*b0c7db22SLisandro Dalcin }
405