Actual source code: network.c

  1: #include <petsc/private/dmnetworkimpl.h>

  3: PetscLogEvent DMNetwork_LayoutSetUp;
  4: PetscLogEvent DMNetwork_SetUpNetwork;
  5: PetscLogEvent DMNetwork_Distribute;

  7: /*
  8:   Creates the component header and value objects for a network point
  9: */
 10: static PetscErrorCode SetUpNetworkHeaderComponentValue(DM dm, DMNetworkComponentHeader header, DMNetworkComponentValue cvalue)
 11: {
 12:   /* Allocate arrays for component information */
 13:   PetscCalloc5(header->maxcomps, &header->size, header->maxcomps, &header->key, header->maxcomps, &header->offset, header->maxcomps, &header->nvar, header->maxcomps, &header->offsetvarrel);
 14:   PetscCalloc1(header->maxcomps, &cvalue->data);

 16:   /* The size of the header is the size of struct _p_DMNetworkComponentHeader. Since the struct contains PetscInt pointers we cannot use sizeof(struct). So, we need to explicitly calculate the size.
 17:    If the data header struct changes then this header size calculation needs to be updated. */
 18:   header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5 * header->maxcomps * sizeof(PetscInt);
 19:   header->hsize /= sizeof(DMNetworkComponentGenericDataType);
 20:   return 0;
 21: }

 23: PetscErrorCode DMNetworkInitializeHeaderComponentData(DM dm)
 24: {
 25:   DM_Network *network = (DM_Network *)dm->data;
 26:   PetscInt    np, p, defaultnumcomp = DMNETWORK_MAX_COMP_AT_POINT_DEFAULT;

 28:   np = network->cloneshared->pEnd - network->cloneshared->pStart;
 29:   if (network->header)
 30:     for (p = 0; p < np; p++) {
 31:       network->header[p].maxcomps = defaultnumcomp;
 32:       SetUpNetworkHeaderComponentValue(dm, &network->header[p], &network->cvalue[p]);
 33:     }

 35:   return 0;
 36: }
 37: /*@
 38:   DMNetworkGetPlex - Gets the Plex DM associated with this network DM

 40:   Not collective

 42:   Input Parameters:
 43: . dm - the dm object

 45:   Output Parameters:
 46: . plexdm - the plex dm object

 48:   Level: Advanced

 50: .seealso: `DMNetworkCreate()`
 51: @*/
 52: PetscErrorCode DMNetworkGetPlex(DM dm, DM *plexdm)
 53: {
 54:   DM_Network *network = (DM_Network *)dm->data;

 56:   *plexdm = network->plex;
 57:   return 0;
 58: }

 60: /*@
 61:   DMNetworkGetNumSubNetworks - Gets the the number of subnetworks

 63:   Not collective

 65:   Input Parameter:
 66: . dm - the dm object

 68:   Output Parameters:
 69: + nsubnet - local number of subnetworks
 70: - Nsubnet - global number of subnetworks

 72:   Level: beginner

 74: .seealso: `DMNetworkCreate()`, `DMNetworkSetNumSubNetworks()`
 75: @*/
 76: PetscErrorCode DMNetworkGetNumSubNetworks(DM dm, PetscInt *nsubnet, PetscInt *Nsubnet)
 77: {
 78:   DM_Network *network = (DM_Network *)dm->data;

 80:   if (nsubnet) *nsubnet = network->cloneshared->nsubnet;
 81:   if (Nsubnet) *Nsubnet = network->cloneshared->Nsubnet;
 82:   return 0;
 83: }

 85: /*@
 86:   DMNetworkSetNumSubNetworks - Sets the number of subnetworks

 88:   Collective on dm

 90:   Input Parameters:
 91: + dm - the dm object
 92: . nsubnet - local number of subnetworks
 93: - Nsubnet - global number of subnetworks

 95:    Level: beginner

 97: .seealso: `DMNetworkCreate()`, `DMNetworkGetNumSubNetworks()`
 98: @*/
 99: PetscErrorCode DMNetworkSetNumSubNetworks(DM dm, PetscInt nsubnet, PetscInt Nsubnet)
100: {
101:   DM_Network *network = (DM_Network *)dm->data;



109:   if (Nsubnet == PETSC_DECIDE) {
111:     MPIU_Allreduce(&nsubnet, &Nsubnet, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm));
112:   }

115:   network->cloneshared->Nsubnet = Nsubnet;
116:   network->cloneshared->nsubnet = 0; /* initial value; will be determind by DMNetworkAddSubnetwork() */
117:   PetscCalloc1(Nsubnet, &network->cloneshared->subnet);

119:   /* num of shared vertices */
120:   network->cloneshared->nsvtx = 0;
121:   network->cloneshared->Nsvtx = 0;
122:   return 0;
123: }

125: /*@
126:   DMNetworkAddSubnetwork - Add a subnetwork

128:   Collective on dm

130:   Input Parameters:
131: + dm - the dm object
132: . name - name of the subnetwork
133: . ne - number of local edges of this subnetwork
134: - edgelist - list of edges for this subnetwork, this is a one dimensional array with pairs of entries being the two vertices (in global numbering of the vertices) of each edge,
135: $            [first vertex of first edge, second vertex of first edge, first vertex of second edge, second vertex of second edge, etc]

137:   Output Parameters:
138: . netnum - global index of the subnetwork

140:   Notes:
141:   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
142:   not be destroyed before the call to DMNetworkLayoutSetUp()

144:   A network can comprise of a single subnetwork OR multiple subnetworks. For a single subnetwork, the subnetwork can be read either in serial or parallel. For a multiple subnetworks,
145:   each subnetwork topology needs to be set on a unique rank and the communicator size needs to be at least equal to the number of subnetworks.

147:   Level: beginner

149:   Example usage:
150:   Consider the following networks:
151:   1) A single subnetwork:
152: .vb
153:  network 0:
154:  rank[0]:
155:    v0 --> v2; v1 --> v2
156:  rank[1]:
157:    v3 --> v5; v4 --> v5
158: .ve

160:  The resulting input
161:  network 0:
162:  rank[0]:
163:    ne = 2
164:    edgelist = [0 2 | 1 2]
165:  rank[1]:
166:    ne = 2
167:    edgelist = [3 5 | 4 5]

169:   2) Two subnetworks:
170: .vb
171:  subnetwork 0:
172:  rank[0]:
173:    v0 --> v2; v2 --> v1; v1 --> v3;
174:  subnetwork 1:
175:  rank[1]:
176:    v0 --> v3; v3 --> v2; v2 --> v1;
177: .ve

179:  The resulting input
180:  subnetwork 0:
181:  rank[0]:
182:    ne = 3
183:    edgelist = [0 2 | 2 1 | 1 3]
184:  rank[1]:
185:    ne = 0
186:    edgelist = NULL

188:  subnetwork 1:
189:  rank[0]:
190:    ne = 0
191:    edgelist = NULL
192:  rank[1]:
193:    edgelist = [0 3 | 3 2 | 2 1]

195: .seealso: `DMNetworkCreate()`, `DMNetworkSetNumSubnetworks()`
196: @*/
197: PetscErrorCode DMNetworkAddSubnetwork(DM dm, const char *name, PetscInt ne, PetscInt edgelist[], PetscInt *netnum)
198: {
199:   DM_Network *network = (DM_Network *)dm->data;
200:   PetscInt    i, Nedge, j, Nvtx, nvtx, nvtx_min = -1, nvtx_max = 0;
201:   PetscBT     table;


205:   i = 0;
206:   if (ne) nvtx_min = nvtx_max = edgelist[0];
207:   for (j = 0; j < ne; j++) {
208:     nvtx_min = PetscMin(nvtx_min, edgelist[i]);
209:     nvtx_max = PetscMax(nvtx_max, edgelist[i]);
210:     i++;
211:     nvtx_min = PetscMin(nvtx_min, edgelist[i]);
212:     nvtx_max = PetscMax(nvtx_max, edgelist[i]);
213:     i++;
214:   }
215:   Nvtx = nvtx_max - nvtx_min + 1; /* approximated total local nvtx for this subnet */

217:   /* Get exact local nvtx for this subnet: counting local values between nvtx_min and nvtx_max */
218:   PetscBTCreate(Nvtx, &table);
219:   PetscBTMemzero(Nvtx, table);
220:   i = 0;
221:   for (j = 0; j < ne; j++) {
222:     PetscBTSet(table, edgelist[i++] - nvtx_min);
223:     PetscBTSet(table, edgelist[i++] - nvtx_min);
224:   }
225:   nvtx = 0;
226:   for (j = 0; j < Nvtx; j++) {
227:     if (PetscBTLookup(table, j)) nvtx++;
228:   }

230:   /* Get global total Nvtx = max(edgelist[])+1 for this subnet */
231:   MPIU_Allreduce(&nvtx_max, &Nvtx, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm));
232:   Nvtx++;
233:   PetscBTDestroy(&table);

235:   /* Get global total Nedge for this subnet */
236:   MPIU_Allreduce(&ne, &Nedge, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm));

238:   i = network->cloneshared->nsubnet;
239:   if (name) PetscStrcpy(network->cloneshared->subnet[i].name, name);
240:   network->cloneshared->subnet[i].nvtx     = nvtx; /* include ghost vertices */
241:   network->cloneshared->subnet[i].nedge    = ne;
242:   network->cloneshared->subnet[i].edgelist = edgelist;
243:   network->cloneshared->subnet[i].Nvtx     = Nvtx;
244:   network->cloneshared->subnet[i].Nedge    = Nedge;

246:   /* ----------------------------------------------------------
247:    p=v or e;
248:    subnet[0].pStart   = 0
249:    subnet[i+1].pStart = subnet[i].pEnd = subnet[i].pStart + (nE[i] or NV[i])
250:    ----------------------------------------------------------------------- */
251:   /* GLOBAL subnet[].vStart and vEnd, used by DMNetworkLayoutSetUp() */
252:   network->cloneshared->subnet[i].vStart = network->cloneshared->NVertices;
253:   network->cloneshared->subnet[i].vEnd   = network->cloneshared->subnet[i].vStart + network->cloneshared->subnet[i].Nvtx; /* global vEnd of subnet[i] */

255:   network->cloneshared->nVertices += nvtx; /* include ghost vertices */
256:   network->cloneshared->NVertices += network->cloneshared->subnet[i].Nvtx;

258:   /* LOCAL subnet[].eStart and eEnd, used by DMNetworkLayoutSetUp() */
259:   network->cloneshared->subnet[i].eStart = network->cloneshared->nEdges;
260:   network->cloneshared->subnet[i].eEnd   = network->cloneshared->subnet[i].eStart + ne;
261:   network->cloneshared->nEdges += ne;
262:   network->cloneshared->NEdges += network->cloneshared->subnet[i].Nedge;

264:   PetscStrcpy(network->cloneshared->subnet[i].name, name);
265:   if (netnum) *netnum = network->cloneshared->nsubnet;
266:   network->cloneshared->nsubnet++;
267:   return 0;
268: }

270: /*@C
271:   DMNetworkSharedVertexGetInfo - Get info of a shared vertex struct, see petsc/private/dmnetworkimpl.h

273:   Not collective

275:   Input Parameters:
276: + dm - the DM object
277: - v - vertex point

279:   Output Parameters:
280: + gidx - global number of this shared vertex in the internal dmplex
281: . n - number of subnetworks that share this vertex
282: - sv - array of size n: sv[2*i,2*i+1]=(net[i], idx[i]), i=0,...,n-1

284:   Level: intermediate

286: .seealso: `DMNetworkGetSharedVertices()`
287: @*/
288: PetscErrorCode DMNetworkSharedVertexGetInfo(DM dm, PetscInt v, PetscInt *gidx, PetscInt *n, const PetscInt **sv)
289: {
290:   DM_Network *network = (DM_Network *)dm->data;
291:   SVtx       *svtx    = network->cloneshared->svtx;
292:   PetscInt    i, gidx_tmp;

294:   DMNetworkGetGlobalVertexIndex(dm, v, &gidx_tmp);
295:   PetscTableFind(network->cloneshared->svtable, gidx_tmp + 1, &i);

298:   i--;
299:   if (gidx) *gidx = gidx_tmp;
300:   if (n) *n = svtx[i].n;
301:   if (sv) *sv = svtx[i].sv;
302:   return 0;
303: }

305: /*
306:   VtxGetInfo - Get info of an input vertex=(net,idx)

308:   Input Parameters:
309: + Nsvtx - global num of shared vertices
310: . svtx - array of shared vertices (global)
311: - (net,idx) - subnet number and local index for a vertex

313:   Output Parameters:
314: + gidx - global index of (net,idx)
315: . svtype - see petsc/private/dmnetworkimpl.h
316: - svtx_idx - ordering in the svtx array
317: */
318: static inline PetscErrorCode VtxGetInfo(PetscInt Nsvtx, SVtx *svtx, PetscInt net, PetscInt idx, PetscInt *gidx, SVtxType *svtype, PetscInt *svtx_idx)
319: {
320:   PetscInt i, j, *svto, g_idx;
321:   SVtxType vtype;

323:   if (!Nsvtx) return 0;

325:   g_idx = -1;
326:   vtype = SVNONE;

328:   for (i = 0; i < Nsvtx; i++) {
329:     if (net == svtx[i].sv[0] && idx == svtx[i].sv[1]) {
330:       g_idx = svtx[i].gidx;
331:       vtype = SVFROM;
332:     } else { /* loop over svtx[i].n */
333:       for (j = 1; j < svtx[i].n; j++) {
334:         svto = svtx[i].sv + 2 * j;
335:         if (net == svto[0] && idx == svto[1]) {
336:           /* input vertex net.idx is a shared to_vertex, output its global index and its svtype */
337:           g_idx = svtx[i].gidx; /* output gidx for to_vertex */
338:           vtype = SVTO;
339:         }
340:       }
341:     }
342:     if (vtype != SVNONE) break;
343:   }
344:   if (gidx) *gidx = g_idx;
345:   if (svtype) *svtype = vtype;
346:   if (svtx_idx) *svtx_idx = i;
347:   return 0;
348: }

350: /*
351:   TableAddSVtx - Add a new shared vertice from sedgelist[k] to a ctable svta

353:   Input:  network, sedgelist, k, svta
354:   Output: svta, tdata, ta2sv
355: */
356: static inline PetscErrorCode TableAddSVtx(DM_Network *network, PetscInt *sedgelist, PetscInt k, PetscTable svta, PetscInt *tdata, PetscInt *ta2sv)
357: {
358:   PetscInt net, idx, gidx;

360:   net  = sedgelist[k];
361:   idx  = sedgelist[k + 1];
362:   gidx = network->cloneshared->subnet[net].vStart + idx;
363:   PetscTableAdd(svta, gidx + 1, *tdata + 1, INSERT_VALUES);

365:   ta2sv[*tdata] = k; /* maps tdata to index of sedgelist */
366:   (*tdata)++;
367:   return 0;
368: }

370: /*
371:   SharedVtxCreate - Create an array of global shared vertices. See SVtx and SVtxType in dmnetworkimpl.h

373:   Input:  dm, Nsedgelist, sedgelist

375:   Note: Output svtx is organized as
376:         sv(net[0],idx[0]) --> sv(net[1],idx[1])
377:                           --> sv(net[1],idx[1])
378:                           ...
379:                           --> sv(net[n-1],idx[n-1])
380:         and net[0] < net[1] < ... < net[n-1]
381:         where sv[0] has SVFROM type, sv[i], i>0, has SVTO type.
382:  */
383: static PetscErrorCode SharedVtxCreate(DM dm, PetscInt Nsedgelist, PetscInt *sedgelist)
384: {
385:   SVtx              *svtx = NULL;
386:   PetscInt          *sv, k, j, nsv, *tdata, **ta2sv;
387:   PetscTable        *svtas;
388:   PetscInt           gidx, net, idx, i, nta, ita, idx_from, idx_to, n, *net_tmp, *idx_tmp, *gidx_tmp;
389:   DM_Network        *network = (DM_Network *)dm->data;
390:   PetscTablePosition ppos;

392:   /* (1) Crete an array of ctables svtas to map (net,idx) -> gidx; a svtas[] for a shared/merged vertex */
393:   PetscCalloc3(Nsedgelist, &svtas, Nsedgelist, &tdata, 2 * Nsedgelist, &ta2sv);

395:   k   = 0; /* sedgelist vertex counter j = 4*k */
396:   nta = 0; /* num of svta tables created = num of groups of shared vertices */

398:   /* for j=0 */
399:   PetscTableCreate(2 * Nsedgelist, network->cloneshared->NVertices + 1, &svtas[nta]);
400:   PetscMalloc1(2 * Nsedgelist, &ta2sv[nta]);

402:   TableAddSVtx(network, sedgelist, k, svtas[nta], &tdata[nta], ta2sv[nta]);
403:   TableAddSVtx(network, sedgelist, k + 2, svtas[nta], &tdata[nta], ta2sv[nta]);
404:   nta++;
405:   k += 4;

407:   for (j = 1; j < Nsedgelist; j++) { /* j: sedgelist counter */
408:     for (ita = 0; ita < nta; ita++) {
409:       /* vfrom */
410:       net  = sedgelist[k];
411:       idx  = sedgelist[k + 1];
412:       gidx = network->cloneshared->subnet[net].vStart + idx; /* global index of the vertex net.idx before merging shared vertices */
413:       PetscTableFind(svtas[ita], gidx + 1, &idx_from);

415:       /* vto */
416:       net  = sedgelist[k + 2];
417:       idx  = sedgelist[k + 3];
418:       gidx = network->cloneshared->subnet[net].vStart + idx;
419:       PetscTableFind(svtas[ita], gidx + 1, &idx_to);

421:       if (idx_from || idx_to) { /* vfrom or vto is on table svtas[ita] */
422:         idx_from--;
423:         idx_to--;
424:         if (idx_from < 0) { /* vto is on svtas[ita] */
425:           TableAddSVtx(network, sedgelist, k, svtas[ita], &tdata[ita], ta2sv[ita]);
426:           break;
427:         } else if (idx_to < 0) {
428:           TableAddSVtx(network, sedgelist, k + 2, svtas[ita], &tdata[ita], ta2sv[ita]);
429:           break;
430:         }
431:       }
432:     }

434:     if (ita == nta) {
435:       PetscTableCreate(2 * Nsedgelist, network->cloneshared->NVertices + 1, &svtas[nta]);
436:       PetscMalloc1(2 * Nsedgelist, &ta2sv[nta]);

438:       TableAddSVtx(network, sedgelist, k, svtas[nta], &tdata[nta], ta2sv[nta]);
439:       TableAddSVtx(network, sedgelist, k + 2, svtas[nta], &tdata[nta], ta2sv[nta]);
440:       nta++;
441:     }
442:     k += 4;
443:   }

445:   /* (2) Create svtable for query shared vertices using gidx */
446:   PetscTableCreate(nta, network->cloneshared->NVertices + 1, &network->cloneshared->svtable);

448:   /* (3) Construct svtx from svtas
449:    svtx: array of SVtx: sv[0]=(net[0],idx[0]) to vertices sv[k], k=1,...,n-1. */
450:   PetscMalloc1(nta, &svtx);
451:   for (nsv = 0; nsv < nta; nsv++) {
452:     /* for a single svtx, put shared vertices in ascending order of gidx */
453:     PetscTableGetCount(svtas[nsv], &n);
454:     PetscCalloc1(2 * n, &sv);
455:     PetscMalloc3(n, &gidx_tmp, n, &net_tmp, n, &idx_tmp);
456:     svtx[nsv].sv   = sv;
457:     svtx[nsv].n    = n;
458:     svtx[nsv].gidx = network->cloneshared->NVertices; /* initialization */

460:     PetscTableGetHeadPosition(svtas[nsv], &ppos);
461:     for (k = 0; k < n; k++) { /* gidx is sorted in ascending order */
462:       PetscTableGetNext(svtas[nsv], &ppos, &gidx, &i);
463:       gidx--;
464:       i--;
465:       j           = ta2sv[nsv][i];    /* maps i to index of sedgelist */
466:       net_tmp[k]  = sedgelist[j];     /* subnet number */
467:       idx_tmp[k]  = sedgelist[j + 1]; /* index on the subnet */
468:       gidx_tmp[k] = gidx;             /* gidx in un-merged dmnetwork */
469:     }

471:     /* current implementation requires sv[]=[net,idx] in ascending order of its gidx in un-merged dmnetwork */
472:     PetscSortIntWithArrayPair(n, gidx_tmp, net_tmp, idx_tmp);
473:     svtx[nsv].gidx = gidx_tmp[0]; /* = min(gidx) */
474:     for (k = 0; k < n; k++) {
475:       sv[2 * k]     = net_tmp[k];
476:       sv[2 * k + 1] = idx_tmp[k];
477:     }
478:     PetscFree3(gidx_tmp, net_tmp, idx_tmp);

480:     /* Setup svtable for query shared vertices */
481:     PetscTableAdd(network->cloneshared->svtable, svtx[nsv].gidx + 1, nsv + 1, INSERT_VALUES);
482:   }

484:   for (j = 0; j < nta; j++) {
485:     PetscTableDestroy(&svtas[j]);
486:     PetscFree(ta2sv[j]);
487:   }
488:   PetscFree3(svtas, tdata, ta2sv);

490:   network->cloneshared->Nsvtx = nta;
491:   network->cloneshared->svtx  = svtx;
492:   return 0;
493: }

495: /*
496:   GetEdgelist_Coupling - Get an integrated edgelist for dmplex from user-provided subnet[].edgelist when subnets are coupled by shared vertices

498:   Input Parameters:
499: . dm - the dmnetwork object

501:    Output Parameters:
502: +  edges - the integrated edgelist for dmplex
503: -  nmerged_ptr - num of vertices being merged
504: */
505: static PetscErrorCode GetEdgelist_Coupling(DM dm, PetscInt *edges, PetscInt *nmerged_ptr)
506: {
507:   MPI_Comm    comm;
508:   PetscMPIInt size, rank, *recvcounts = NULL, *displs = NULL;
509:   DM_Network *network = (DM_Network *)dm->data;
510:   PetscInt    i, j, ctr, np;
511:   PetscInt   *vidxlTog, Nsv, Nsubnet = network->cloneshared->Nsubnet;
512:   PetscInt   *sedgelist = network->cloneshared->sedgelist;
513:   PetscInt    net, idx, gidx, nmerged, *vrange, gidx_from, net_from, sv_idx;
514:   SVtxType    svtype = SVNONE;
515:   SVtx       *svtx;

517:   PetscObjectGetComm((PetscObject)dm, &comm);
518:   MPI_Comm_rank(comm, &rank);
519:   MPI_Comm_size(comm, &size);

521:   /* (1) Create global svtx[] from sedgelist */
522:   /* --------------------------------------- */
523:   SharedVtxCreate(dm, network->cloneshared->Nsvtx, sedgelist);
524:   Nsv  = network->cloneshared->Nsvtx;
525:   svtx = network->cloneshared->svtx;

527:   /* (2) Merge shared vto vertices to their vfrom vertex with same global vertex index (gidx) */
528:   /* --------------------------------------------------------------------------------------- */
529:   /* (2.1) compute vrage[rank]: global index of 1st local vertex in proc[rank] */
530:   PetscMalloc4(size + 1, &vrange, size, &displs, size, &recvcounts, network->cloneshared->nVertices, &vidxlTog);
531:   for (i = 0; i < size; i++) {
532:     displs[i]     = i;
533:     recvcounts[i] = 1;
534:   }

536:   vrange[0] = 0;
537:   MPI_Allgatherv(&network->cloneshared->nVertices, 1, MPIU_INT, vrange + 1, recvcounts, displs, MPIU_INT, comm);
538:   for (i = 2; i < size + 1; i++) vrange[i] += vrange[i - 1];

540:   /* (2.2) Create vidxlTog: maps UN-MERGED local vertex index i to global index gidx (plex, excluding ghost vertices) */
541:   i                           = 0;
542:   gidx                        = 0;
543:   nmerged                     = 0; /* local num of merged vertices */
544:   network->cloneshared->nsvtx = 0; /* local num of SVtx structs, including ghosts */
545:   for (net = 0; net < Nsubnet; net++) {
546:     for (idx = 0; idx < network->cloneshared->subnet[net].Nvtx; idx++) { /* Note: global subnet[net].Nvtx */
547:       VtxGetInfo(Nsv, svtx, net, idx, &gidx_from, &svtype, &sv_idx);
548:       if (svtype == SVTO) {
549:         if (network->cloneshared->subnet[net].nvtx) { /* this proc owns sv_to */
550:           net_from = svtx[sv_idx].sv[0];              /* subnet number of its shared vertex */
551:           if (network->cloneshared->subnet[net_from].nvtx == 0) {
552:             /* this proc does not own v_from, thus a ghost local vertex */
553:             network->cloneshared->nsvtx++;
554:           }
555:           vidxlTog[i++] = gidx_from; /* gidx before merging! Bug??? */
556:           nmerged++;                 /* a shared vertex -- merged */
557:         }
558:       } else {
559:         if (svtype == SVFROM && network->cloneshared->subnet[net].nvtx) {
560:           /* this proc owns this v_from, a new local shared vertex */
561:           network->cloneshared->nsvtx++;
562:         }
563:         if (network->cloneshared->subnet[net].nvtx) vidxlTog[i++] = gidx;
564:         gidx++;
565:       }
566:     }
567:   }
568: #if defined(PETSC_USE_DEBUG)
570: #endif

572:   /* (2.3) Shared vertices in the subnetworks are merged, update global NVertices: np = sum(local nmerged) */
573:   MPI_Allreduce(&nmerged, &np, 1, MPIU_INT, MPI_SUM, comm);
574:   network->cloneshared->NVertices -= np;

576:   ctr = 0;
577:   for (net = 0; net < Nsubnet; net++) {
578:     for (j = 0; j < network->cloneshared->subnet[net].nedge; j++) {
579:       /* vfrom: */
580:       i              = network->cloneshared->subnet[net].edgelist[2 * j] + (network->cloneshared->subnet[net].vStart - vrange[rank]);
581:       edges[2 * ctr] = vidxlTog[i];

583:       /* vto */
584:       i                  = network->cloneshared->subnet[net].edgelist[2 * j + 1] + (network->cloneshared->subnet[net].vStart - vrange[rank]);
585:       edges[2 * ctr + 1] = vidxlTog[i];
586:       ctr++;
587:     }
588:   }
589:   PetscFree4(vrange, displs, recvcounts, vidxlTog);
590:   PetscFree(sedgelist); /* created in DMNetworkAddSharedVertices() */

592:   *nmerged_ptr = nmerged;
593:   return 0;
594: }

596: PetscErrorCode DMNetworkInitializeNonTopological(DM dm)
597: {
598:   DM_Network *network = (DM_Network *)dm->data;
599:   PetscInt    p, pStart = network->cloneshared->pStart, pEnd = network->cloneshared->pEnd;
600:   MPI_Comm    comm;

602:   PetscObjectGetComm((PetscObject)dm, &comm);

604:   PetscSectionCreate(comm, &network->DataSection);
605:   PetscSectionCreate(comm, &network->DofSection);
606:   PetscSectionSetChart(network->DataSection, pStart, pEnd);
607:   PetscSectionSetChart(network->DofSection, pStart, pEnd);

609:   DMNetworkInitializeHeaderComponentData(dm);

611:   for (p = 0; p < pEnd - pStart; p++) {
612:     network->header[p].ndata           = 0;
613:     network->header[p].offset[0]       = 0;
614:     network->header[p].offsetvarrel[0] = 0;
615:     PetscSectionAddDof(network->DataSection, p, network->header[p].hsize);
616:   }
617:   return 0;
618: }

620: /*@
621:   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network

623:   Not Collective

625:   Input Parameters:
626: . dm - the dmnetwork object

628:   Notes:
629:   This routine should be called after the network sizes and edgelists have been provided. It creates
630:   the bare layout of the network and sets up the network to begin insertion of components.

632:   All the components should be registered before calling this routine.

634:   Level: beginner

636: .seealso: `DMNetworkSetNumSubNetworks()`, `DMNetworkAddSubnetwork()`
637: @*/
638: PetscErrorCode DMNetworkLayoutSetUp(DM dm)
639: {
640:   DM_Network     *network = (DM_Network *)dm->data;
641:   PetscInt        i, j, ctr, Nsubnet = network->cloneshared->Nsubnet, np, *edges, *subnetvtx, *subnetedge, e, v, vfrom, vto, net, globaledgeoff;
642:   const PetscInt *cone;
643:   MPI_Comm        comm;
644:   PetscMPIInt     size;
645:   PetscSection    sectiong;
646:   PetscInt        nmerged = 0;

648:   PetscLogEventBegin(DMNetwork_LayoutSetUp, dm, 0, 0, 0);

651:   /* This implementation requires user input each subnet by a single processor when Nsubnet>1, thus subnet[net].nvtx=subnet[net].Nvtx when net>0 */
652:   for (net = 1; net < Nsubnet; net++) {
653:     if (network->cloneshared->subnet[net].nvtx)
655:                  network->cloneshared->subnet[net].nvtx, network->cloneshared->subnet[net].Nvtx);
656:   }

658:   PetscObjectGetComm((PetscObject)dm, &comm);
659:   MPI_Comm_size(comm, &size);

661:   /* Create LOCAL edgelist in global vertex ordering for the network by concatenating local input edgelists of the subnetworks */
662:   PetscCalloc1(2 * network->cloneshared->nEdges, &edges);

664:   if (network->cloneshared->Nsvtx) { /* subnetworks are coupled via shared vertices */
665:     GetEdgelist_Coupling(dm, edges, &nmerged);
666:   } else { /* subnetworks are not coupled */
667:     /* Create a 0-size svtable for query shared vertices */
668:     PetscTableCreate(0, network->cloneshared->NVertices + 1, &network->cloneshared->svtable);
669:     ctr = 0;
670:     for (i = 0; i < Nsubnet; i++) {
671:       for (j = 0; j < network->cloneshared->subnet[i].nedge; j++) {
672:         edges[2 * ctr]     = network->cloneshared->subnet[i].vStart + network->cloneshared->subnet[i].edgelist[2 * j];
673:         edges[2 * ctr + 1] = network->cloneshared->subnet[i].vStart + network->cloneshared->subnet[i].edgelist[2 * j + 1];
674:         ctr++;
675:       }
676:     }
677:   }

679:   /* Create network->plex; One dimensional network, numCorners=2 */
680:   DMCreate(comm, &network->plex);
681:   DMSetType(network->plex, DMPLEX);
682:   DMSetDimension(network->plex, 1);

684:   if (size == 1) DMPlexBuildFromCellList(network->plex, network->cloneshared->nEdges, PETSC_DECIDE, 2, edges);
685:   else DMPlexBuildFromCellListParallel(network->plex, network->cloneshared->nEdges, PETSC_DECIDE, PETSC_DECIDE, 2, edges, NULL, NULL);

687:   DMPlexGetChart(network->plex, &network->cloneshared->pStart, &network->cloneshared->pEnd);
688:   DMPlexGetHeightStratum(network->plex, 0, &network->cloneshared->eStart, &network->cloneshared->eEnd);
689:   DMPlexGetHeightStratum(network->plex, 1, &network->cloneshared->vStart, &network->cloneshared->vEnd);
690:   np = network->cloneshared->pEnd - network->cloneshared->pStart;
691:   PetscCalloc2(np, &network->header, np, &network->cvalue);

693:   /* Create edge and vertex arrays for the subnetworks
694:      This implementation assumes that DMNetwork reads
695:      (1) a single subnetwork in parallel; or
696:      (2) n subnetworks using n processors, one subnetwork/processor.
697:   */
698:   PetscCalloc2(network->cloneshared->nEdges, &subnetedge, network->cloneshared->nVertices + network->cloneshared->nsvtx, &subnetvtx); /* Maps local edge/vertex to local subnetwork's edge/vertex */
699:   network->cloneshared->subnetedge = subnetedge;
700:   network->cloneshared->subnetvtx  = subnetvtx;
701:   for (j = 0; j < Nsubnet; j++) {
702:     network->cloneshared->subnet[j].edges = subnetedge;
703:     subnetedge += network->cloneshared->subnet[j].nedge;

705:     network->cloneshared->subnet[j].vertices = subnetvtx;
706:     subnetvtx += network->cloneshared->subnet[j].nvtx;
707:   }
708:   network->cloneshared->svertices = subnetvtx;

710:   /* Get edge ownership */
711:   np = network->cloneshared->eEnd - network->cloneshared->eStart;
712:   MPI_Scan(&np, &globaledgeoff, 1, MPIU_INT, MPI_SUM, comm);
713:   globaledgeoff -= np;

715:   /* Setup local edge and vertex arrays for subnetworks */
716:   e = 0;
717:   for (i = 0; i < Nsubnet; i++) {
718:     ctr = 0;
719:     for (j = 0; j < network->cloneshared->subnet[i].nedge; j++) {
720:       /* edge e */
721:       network->header[e].index                 = e + globaledgeoff; /* Global edge index */
722:       network->header[e].subnetid              = i;
723:       network->cloneshared->subnet[i].edges[j] = e;

725:       /* connected vertices */
726:       DMPlexGetCone(network->plex, e, &cone);

728:       /* vertex cone[0] */
729:       v                           = cone[0];
730:       network->header[v].index    = edges[2 * e]; /* Global vertex index */
731:       network->header[v].subnetid = i;            /* Subnetwork id */
732:       if (Nsubnet == 1) {
733:         network->cloneshared->subnet[i].vertices[v - network->cloneshared->vStart] = v; /* user's subnet[].idx = petsc's v */
734:       } else {
735:         vfrom                                           = network->cloneshared->subnet[i].edgelist[2 * ctr]; /* =subnet[i].idx, Global index! */
736:         network->cloneshared->subnet[i].vertices[vfrom] = v;                                                 /* user's subnet[].dix = petsc's v */
737:       }

739:       /* vertex cone[1] */
740:       v                           = cone[1];
741:       network->header[v].index    = edges[2 * e + 1]; /* Global vertex index */
742:       network->header[v].subnetid = i;                /* Subnetwork id */
743:       if (Nsubnet == 1) {
744:         network->cloneshared->subnet[i].vertices[v - network->cloneshared->vStart] = v; /* user's subnet[].idx = petsc's v */
745:       } else {
746:         vto                                           = network->cloneshared->subnet[i].edgelist[2 * ctr + 1]; /* =subnet[i].idx, Global index! */
747:         network->cloneshared->subnet[i].vertices[vto] = v;                                                     /* user's subnet[].dix = petsc's v */
748:       }

750:       e++;
751:       ctr++;
752:     }
753:   }
754:   PetscFree(edges);

756:   /* Set local vertex array for the subnetworks */
757:   j = 0;
758:   for (v = network->cloneshared->vStart; v < network->cloneshared->vEnd; v++) {
759:     /* local shared vertex */
760:     PetscTableFind(network->cloneshared->svtable, network->header[v].index + 1, &i);
761:     if (i) network->cloneshared->svertices[j++] = v;
762:   }

764:   /* Create a global section to be used by DMNetworkIsGhostVertex() which is a non-collective routine */
765:   /* see snes_tutorials_network-ex1_4 */
766:   DMGetGlobalSection(network->plex, &sectiong);
767:   /* Initialize non-topological data structures  */
768:   DMNetworkInitializeNonTopological(dm);
769:   PetscLogEventEnd(DMNetwork_LayoutSetUp, dm, 0, 0, 0);
770:   return 0;
771: }

773: /*@C
774:   DMNetworkGetSubnetwork - Returns the information about a requested subnetwork

776:   Not collective

778:   Input Parameters:
779: + dm - the DM object
780: - netnum - the global index of the subnetwork

782:   Output Parameters:
783: + nv - number of vertices (local)
784: . ne - number of edges (local)
785: . vtx - local vertices of the subnetwork
786: - edge - local edges of the subnetwork

788:   Notes:
789:     Cannot call this routine before DMNetworkLayoutSetup()

791:     The local vertices returned on each rank are determined by DMNetwork. The user does not have any control over what vertices are local.

793:   Level: intermediate

795: .seealso: `DMNetworkCreate()`, `DMNetworkAddSubnetwork()`, `DMNetworkLayoutSetUp()`
796: @*/
797: PetscErrorCode DMNetworkGetSubnetwork(DM dm, PetscInt netnum, PetscInt *nv, PetscInt *ne, const PetscInt **vtx, const PetscInt **edge)
798: {
799:   DM_Network *network = (DM_Network *)dm->data;

802:   if (nv) *nv = network->cloneshared->subnet[netnum].nvtx;
803:   if (ne) *ne = network->cloneshared->subnet[netnum].nedge;
804:   if (vtx) *vtx = network->cloneshared->subnet[netnum].vertices;
805:   if (edge) *edge = network->cloneshared->subnet[netnum].edges;
806:   return 0;
807: }

809: /*@
810:   DMNetworkAddSharedVertices - Add shared vertices that connect two given subnetworks

812:   Collective on dm

814:   Input Parameters:
815: + dm - the dm object
816: . anetnum - first subnetwork global numbering returned by DMNetworkAddSubnetwork()
817: . bnetnum - second subnetwork global numbering returned by DMNetworkAddSubnetwork()
818: . nsvtx - number of vertices that are shared by the two subnetworks
819: . asvtx - vertex index in the first subnetwork
820: - bsvtx - vertex index in the second subnetwork

822:   Level: beginner

824: .seealso: `DMNetworkCreate()`, `DMNetworkAddSubnetwork()`, `DMNetworkGetSharedVertices()`
825: @*/
826: PetscErrorCode DMNetworkAddSharedVertices(DM dm, PetscInt anetnum, PetscInt bnetnum, PetscInt nsvtx, PetscInt asvtx[], PetscInt bsvtx[])
827: {
828:   DM_Network *network = (DM_Network *)dm->data;
829:   PetscInt    i, nsubnet = network->cloneshared->Nsubnet, *sedgelist, Nsvtx = network->cloneshared->Nsvtx;

833:   if (!Nsvtx) {
834:     /* allocate network->sedgelist to hold at most 2*nsubnet pairs of shared vertices */
835:     PetscMalloc1(2 * 4 * nsubnet, &network->cloneshared->sedgelist);
836:   }

838:   sedgelist = network->cloneshared->sedgelist;
839:   for (i = 0; i < nsvtx; i++) {
840:     sedgelist[4 * Nsvtx]     = anetnum;
841:     sedgelist[4 * Nsvtx + 1] = asvtx[i];
842:     sedgelist[4 * Nsvtx + 2] = bnetnum;
843:     sedgelist[4 * Nsvtx + 3] = bsvtx[i];
844:     Nsvtx++;
845:   }
847:   network->cloneshared->Nsvtx = Nsvtx;
848:   return 0;
849: }

851: /*@C
852:   DMNetworkGetSharedVertices - Returns the info for the shared vertices

854:   Not collective

856:   Input Parameter:
857: . dm - the DM object

859:   Output Parameters:
860: + nsv - number of local shared vertices
861: - svtx - local shared vertices

863:   Notes:
864:   Cannot call this routine before DMNetworkLayoutSetup()

866:   Level: intermediate

868: .seealso: `DMNetworkGetSubnetwork()`, `DMNetworkLayoutSetUp()`, `DMNetworkAddSharedVertices()`
869: @*/
870: PetscErrorCode DMNetworkGetSharedVertices(DM dm, PetscInt *nsv, const PetscInt **svtx)
871: {
872:   DM_Network *net = (DM_Network *)dm->data;

875:   if (nsv) *nsv = net->cloneshared->nsvtx;
876:   if (svtx) *svtx = net->cloneshared->svertices;
877:   return 0;
878: }

880: /*@C
881:   DMNetworkRegisterComponent - Registers the network component

883:   Logically collective on dm

885:   Input Parameters:
886: + dm - the network object
887: . name - the component name
888: - size - the storage size in bytes for this component data

890:    Output Parameters:
891: .  key - an integer key that defines the component

893:    Note:
894:    This routine should be called by all processors before calling `DMNetworkLayoutSetup()`.

896:    Level: beginner

898: .seealso: `DMNetworkCreate()`, `DMNetworkLayoutSetUp()`
899: @*/
900: PetscErrorCode DMNetworkRegisterComponent(DM dm, const char *name, size_t size, PetscInt *key)
901: {
902:   DM_Network         *network   = (DM_Network *)dm->data;
903:   DMNetworkComponent *component = NULL, *newcomponent = NULL;
904:   PetscBool           flg = PETSC_FALSE;
905:   PetscInt            i;

907:   if (!network->component) PetscCalloc1(network->max_comps_registered, &network->component);

909:   for (i = 0; i < network->ncomponent; i++) {
910:     PetscStrcmp(network->component[i].name, name, &flg);
911:     if (flg) {
912:       *key = i;
913:       return 0;
914:     }
915:   }

917:   if (network->ncomponent == network->max_comps_registered) {
918:     /* Reached max allowed so resize component */
919:     network->max_comps_registered += 2;
920:     PetscCalloc1(network->max_comps_registered, &newcomponent);
921:     /* Copy over the previous component info */
922:     for (i = 0; i < network->ncomponent; i++) {
923:       PetscStrcpy(newcomponent[i].name, network->component[i].name);
924:       newcomponent[i].size = network->component[i].size;
925:     }
926:     /* Free old one */
927:     PetscFree(network->component);
928:     /* Update pointer */
929:     network->component = newcomponent;
930:   }

932:   component = &network->component[network->ncomponent];

934:   PetscStrcpy(component->name, name);
935:   component->size = size / sizeof(DMNetworkComponentGenericDataType);
936:   *key            = network->ncomponent;
937:   network->ncomponent++;
938:   return 0;
939: }

941: /*@
942:   DMNetworkGetVertexRange - Get the bounds [start, end) for the local vertices

944:   Not Collective

946:   Input Parameter:
947: . dm - the DMNetwork object

949:   Output Parameters:
950: + vStart - the first vertex point
951: - vEnd - one beyond the last vertex point

953:   Level: beginner

955: .seealso: `DMNetworkGetEdgeRange()`
956: @*/
957: PetscErrorCode DMNetworkGetVertexRange(DM dm, PetscInt *vStart, PetscInt *vEnd)
958: {
959:   DM_Network *network = (DM_Network *)dm->data;

961:   if (vStart) *vStart = network->cloneshared->vStart;
962:   if (vEnd) *vEnd = network->cloneshared->vEnd;
963:   return 0;
964: }

966: /*@
967:   DMNetworkGetEdgeRange - Get the bounds [start, end) for the local edges

969:   Not Collective

971:   Input Parameter:
972: . dm - the DMNetwork object

974:   Output Parameters:
975: + eStart - The first edge point
976: - eEnd - One beyond the last edge point

978:   Level: beginner

980: .seealso: `DMNetworkGetVertexRange()`
981: @*/
982: PetscErrorCode DMNetworkGetEdgeRange(DM dm, PetscInt *eStart, PetscInt *eEnd)
983: {
984:   DM_Network *network = (DM_Network *)dm->data;

987:   if (eStart) *eStart = network->cloneshared->eStart;
988:   if (eEnd) *eEnd = network->cloneshared->eEnd;
989:   return 0;
990: }

992: PetscErrorCode DMNetworkGetIndex(DM dm, PetscInt p, PetscInt *index)
993: {
994:   DM_Network *network = (DM_Network *)dm->data;

996:   if (network->header) {
997:     *index = network->header[p].index;
998:   } else {
999:     PetscInt                 offsetp;
1000:     DMNetworkComponentHeader header;

1002:     PetscSectionGetOffset(network->DataSection, p, &offsetp);
1003:     header = (DMNetworkComponentHeader)(network->componentdataarray + offsetp);
1004:     *index = header->index;
1005:   }
1006:   return 0;
1007: }

1009: PetscErrorCode DMNetworkGetSubnetID(DM dm, PetscInt p, PetscInt *subnetid)
1010: {
1011:   DM_Network *network = (DM_Network *)dm->data;

1013:   if (network->header) {
1014:     *subnetid = network->header[p].subnetid;
1015:   } else {
1016:     PetscInt                 offsetp;
1017:     DMNetworkComponentHeader header;

1019:     PetscSectionGetOffset(network->DataSection, p, &offsetp);
1020:     header    = (DMNetworkComponentHeader)(network->componentdataarray + offsetp);
1021:     *subnetid = header->subnetid;
1022:   }
1023:   return 0;
1024: }

1026: /*@
1027:   DMNetworkGetGlobalEdgeIndex - Get the global numbering for the edge on the network

1029:   Not Collective

1031:   Input Parameters:
1032: + dm - DMNetwork object
1033: - p - edge point

1035:   Output Parameters:
1036: . index - the global numbering for the edge

1038:   Level: intermediate

1040: .seealso: `DMNetworkGetGlobalVertexIndex()`
1041: @*/
1042: PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm, PetscInt p, PetscInt *index)
1043: {
1044:   DMNetworkGetIndex(dm, p, index);
1045:   return 0;
1046: }

1048: /*@
1049:   DMNetworkGetGlobalVertexIndex - Get the global numbering for the vertex on the network

1051:   Not Collective

1053:   Input Parameters:
1054: + dm - DMNetwork object
1055: - p  - vertex point

1057:   Output Parameters:
1058: . index - the global numbering for the vertex

1060:   Level: intermediate

1062: .seealso: `DMNetworkGetGlobalEdgeIndex()`, `DMNetworkGetLocalVertexIndex()`
1063: @*/
1064: PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm, PetscInt p, PetscInt *index)
1065: {
1066:   DMNetworkGetIndex(dm, p, index);
1067:   return 0;
1068: }

1070: /*@
1071:   DMNetworkGetNumComponents - Get the number of components at a vertex/edge

1073:   Not Collective

1075:   Input Parameters:
1076: + dm - the DMNetwork object
1077: - p - vertex/edge point

1079:   Output Parameters:
1080: . numcomponents - Number of components at the vertex/edge

1082:   Level: beginner

1084: .seealso: `DMNetworkRegisterComponent()`, `DMNetworkAddComponent()`
1085: @*/
1086: PetscErrorCode DMNetworkGetNumComponents(DM dm, PetscInt p, PetscInt *numcomponents)
1087: {
1088:   PetscInt    offset;
1089:   DM_Network *network = (DM_Network *)dm->data;

1091:   PetscSectionGetOffset(network->DataSection, p, &offset);
1092:   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray + offset))->ndata;
1093:   return 0;
1094: }

1096: /*@
1097:   DMNetworkGetLocalVecOffset - Get the offset for accessing the variables associated with a component at the given vertex/edge from the local vector

1099:   Not Collective

1101:   Input Parameters:
1102: + dm - the DMNetwork object
1103: . p - the edge or vertex point
1104: - compnum - component number; use ALL_COMPONENTS if no specific component is requested

1106:   Output Parameters:
1107: . offset - the local offset

1109:   Notes:
1110:     These offsets can be passed to MatSetValuesLocal() for matrices obtained with DMCreateMatrix().

1112:     For vectors obtained with DMCreateLocalVector() the offsets can be used with VecSetValues().

1114:     For vectors obtained with DMCreateLocalVector() and the array obtained with VecGetArray(vec,&array) you can access or set
1115:     the vector values with array[offset].

1117:     For vectors obtained with DMCreateGlobalVector() the offsets can be used with VecSetValuesLocal().

1119:   Level: intermediate

1121: .seealso: `DMGetLocalVector()`, `DMNetworkGetComponent()`, `DMNetworkGetGlobalVecOffset()`, `DMCreateGlobalVector()`, `VecGetArray()`, `VecSetValuesLocal()`, `MatSetValuesLocal()`
1122: @*/
1123: PetscErrorCode DMNetworkGetLocalVecOffset(DM dm, PetscInt p, PetscInt compnum, PetscInt *offset)
1124: {
1125:   DM_Network              *network = (DM_Network *)dm->data;
1126:   PetscInt                 offsetp, offsetd;
1127:   DMNetworkComponentHeader header;

1129:   PetscSectionGetOffset(network->plex->localSection, p, &offsetp);
1130:   if (compnum == ALL_COMPONENTS) {
1131:     *offset = offsetp;
1132:     return 0;
1133:   }

1135:   PetscSectionGetOffset(network->DataSection, p, &offsetd);
1136:   header  = (DMNetworkComponentHeader)(network->componentdataarray + offsetd);
1137:   *offset = offsetp + header->offsetvarrel[compnum];
1138:   return 0;
1139: }

1141: /*@
1142:   DMNetworkGetGlobalVecOffset - Get the global offset for accessing the variables associated with a component for the given vertex/edge from the global vector

1144:   Not Collective

1146:   Input Parameters:
1147: + dm - the DMNetwork object
1148: . p - the edge or vertex point
1149: - compnum - component number; use ALL_COMPONENTS if no specific component is requested

1151:   Output Parameters:
1152: . offsetg - the global offset

1154:   Notes:
1155:     These offsets can be passed to MatSetValues() for matrices obtained with DMCreateMatrix().

1157:     For vectors obtained with DMCreateGlobalVector() the offsets can be used with VecSetValues().

1159:     For vectors obtained with DMCreateGlobalVector() and the array obtained with VecGetArray(vec,&array) you can access or set
1160:     the vector values with array[offset - rstart] where restart is obtained with VecGetOwnershipRange(v,&rstart,NULL);

1162:   Level: intermediate

1164: .seealso: `DMNetworkGetLocalVecOffset()`, `DMGetGlobalVector()`, `DMNetworkGetComponent()`, `DMCreateGlobalVector()`, `VecGetArray()`, `VecSetValues()`, `MatSetValues()`
1165: @*/
1166: PetscErrorCode DMNetworkGetGlobalVecOffset(DM dm, PetscInt p, PetscInt compnum, PetscInt *offsetg)
1167: {
1168:   DM_Network              *network = (DM_Network *)dm->data;
1169:   PetscInt                 offsetp, offsetd;
1170:   DMNetworkComponentHeader header;

1172:   PetscSectionGetOffset(network->plex->globalSection, p, &offsetp);
1173:   if (offsetp < 0) offsetp = -(offsetp + 1); /* Convert to actual global offset for ghost vertex */

1175:   if (compnum == ALL_COMPONENTS) {
1176:     *offsetg = offsetp;
1177:     return 0;
1178:   }
1179:   PetscSectionGetOffset(network->DataSection, p, &offsetd);
1180:   header   = (DMNetworkComponentHeader)(network->componentdataarray + offsetd);
1181:   *offsetg = offsetp + header->offsetvarrel[compnum];
1182:   return 0;
1183: }

1185: /*@
1186:   DMNetworkGetEdgeOffset - Get the offset for accessing the variables associated with the given edge from the local subvector

1188:   Not Collective

1190:   Input Parameters:
1191: + dm - the DMNetwork object
1192: - p - the edge point

1194:   Output Parameters:
1195: . offset - the offset

1197:   Level: intermediate

1199: .seealso: `DMNetworkGetLocalVecOffset()`, `DMGetLocalVector()`
1200: @*/
1201: PetscErrorCode DMNetworkGetEdgeOffset(DM dm, PetscInt p, PetscInt *offset)
1202: {
1203:   DM_Network *network = (DM_Network *)dm->data;

1205:   PetscSectionGetOffset(network->edge.DofSection, p, offset);
1206:   return 0;
1207: }

1209: /*@
1210:   DMNetworkGetVertexOffset - Get the offset for accessing the variables associated with the given vertex from the local subvector

1212:   Not Collective

1214:   Input Parameters:
1215: + dm - the DMNetwork object
1216: - p - the vertex point

1218:   Output Parameters:
1219: . offset - the offset

1221:   Level: intermediate

1223: .seealso: `DMNetworkGetEdgeOffset()`, `DMGetLocalVector()`
1224: @*/
1225: PetscErrorCode DMNetworkGetVertexOffset(DM dm, PetscInt p, PetscInt *offset)
1226: {
1227:   DM_Network *network = (DM_Network *)dm->data;

1229:   p -= network->cloneshared->vStart;
1230:   PetscSectionGetOffset(network->vertex.DofSection, p, offset);
1231:   return 0;
1232: }

1234: /*@
1235:   DMNetworkAddComponent - Adds a network component and number of variables at the given point (vertex/edge)

1237:   Collective on dm

1239:   Input Parameters:
1240: + dm - the DMNetwork
1241: . p - the vertex/edge point. These points are local indices provided by DMNetworkGetSubnetwork()
1242: . componentkey - component key returned while registering the component with DMNetworkRegisterComponent()
1243: . compvalue - pointer to the data structure for the component, or NULL if the component does not require data, this data is not copied so you cannot
1244:               free this space until after DMSetUp() is called.
1245: - nvar - number of variables for the component at the vertex/edge point, zero if the component does not introduce any degrees of freedom at the point

1247:   Notes:
1248:     The owning rank and any other ranks that have this point as a ghost location must call this routine to add a component and number of variables in the same order at the given point.

1250:     DMNetworkLayoutSetUp() must be called before this routine.

1252:   Developer Notes:
1253:      The requirement that all the ranks with access to a vertex (as owner or as ghost) add all the components comes from a limitation of the underlying implementation based on DMPLEX.
1254:   Level: beginner

1256: .seealso: `DMNetworkGetComponent()`, `DMNetworkGetSubnetwork()`, `DMNetworkIsGhostVertex()`, `DMNetworkLayoutSetUp()`
1257: @*/
1258: PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p, PetscInt componentkey, void *compvalue, PetscInt nvar)
1259: {
1260:   DM_Network              *network   = (DM_Network *)dm->data;
1261:   DMNetworkComponent      *component = &network->component[componentkey];
1262:   DMNetworkComponentHeader header;
1263:   DMNetworkComponentValue  cvalue;
1264:   PetscInt                 compnum;
1265:   PetscInt                *compsize, *compkey, *compoffset, *compnvar, *compoffsetvarrel;
1266:   void                   **compdata;

1270:   /* The owning rank and all ghost ranks add nvar */
1271:   PetscSectionAddDof(network->DofSection, p, nvar);

1273:   /* The owning rank and all ghost ranks add a component, including compvalue=NULL */
1274:   header = &network->header[p];
1275:   cvalue = &network->cvalue[p];
1276:   if (header->ndata == header->maxcomps) {
1277:     PetscInt additional_size;

1279:     /* Reached limit so resize header component arrays */
1280:     header->maxcomps += 2;

1282:     /* Allocate arrays for component information and value */
1283:     PetscCalloc5(header->maxcomps, &compsize, header->maxcomps, &compkey, header->maxcomps, &compoffset, header->maxcomps, &compnvar, header->maxcomps, &compoffsetvarrel);
1284:     PetscMalloc1(header->maxcomps, &compdata);

1286:     /* Recalculate header size */
1287:     header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5 * header->maxcomps * sizeof(PetscInt);

1289:     header->hsize /= sizeof(DMNetworkComponentGenericDataType);

1291:     /* Copy over component info */
1292:     PetscMemcpy(compsize, header->size, header->ndata * sizeof(PetscInt));
1293:     PetscMemcpy(compkey, header->key, header->ndata * sizeof(PetscInt));
1294:     PetscMemcpy(compoffset, header->offset, header->ndata * sizeof(PetscInt));
1295:     PetscMemcpy(compnvar, header->nvar, header->ndata * sizeof(PetscInt));
1296:     PetscMemcpy(compoffsetvarrel, header->offsetvarrel, header->ndata * sizeof(PetscInt));

1298:     /* Copy over component data pointers */
1299:     PetscMemcpy(compdata, cvalue->data, header->ndata * sizeof(void *));

1301:     /* Free old arrays */
1302:     PetscFree5(header->size, header->key, header->offset, header->nvar, header->offsetvarrel);
1303:     PetscFree(cvalue->data);

1305:     /* Update pointers */
1306:     header->size         = compsize;
1307:     header->key          = compkey;
1308:     header->offset       = compoffset;
1309:     header->nvar         = compnvar;
1310:     header->offsetvarrel = compoffsetvarrel;

1312:     cvalue->data = compdata;

1314:     /* Update DataSection Dofs */
1315:     /* The dofs for datasection point p equals sizeof the header (i.e. header->hsize) + sizes of the components added at point p. With the resizing of the header, we need to update the dofs for point p. Hence, we add the extra size added for the header */
1316:     additional_size = (5 * (header->maxcomps - header->ndata) * sizeof(PetscInt)) / sizeof(DMNetworkComponentGenericDataType);
1317:     PetscSectionAddDof(network->DataSection, p, additional_size);
1318:   }
1319:   header = &network->header[p];
1320:   cvalue = &network->cvalue[p];

1322:   compnum = header->ndata;

1324:   header->size[compnum] = component->size;
1325:   PetscSectionAddDof(network->DataSection, p, component->size);
1326:   header->key[compnum] = componentkey;
1327:   if (compnum != 0) header->offset[compnum] = header->offset[compnum - 1] + header->size[compnum - 1];
1328:   cvalue->data[compnum] = (void *)compvalue;

1330:   /* variables */
1331:   header->nvar[compnum] += nvar;
1332:   if (compnum != 0) header->offsetvarrel[compnum] = header->offsetvarrel[compnum - 1] + header->nvar[compnum - 1];

1334:   header->ndata++;
1335:   return 0;
1336: }

1338: /*@
1339:   DMNetworkGetComponent - Gets the component key, the component data, and the number of variables at a given network point

1341:   Not Collective

1343:   Input Parameters:
1344: + dm - the DMNetwork object
1345: . p - vertex/edge point
1346: - compnum - component number; use ALL_COMPONENTS if sum up all the components

1348:   Output Parameters:
1349: + compkey - the key obtained when registering the component (use NULL if not required)
1350: . component - the component data (use NULL if not required)
1351: - nvar  - number of variables (use NULL if not required)

1353:   Level: beginner

1355: .seealso: `DMNetworkAddComponent()`, `DMNetworkGetNumComponents()`
1356: @*/
1357: PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *compkey, void **component, PetscInt *nvar)
1358: {
1359:   DM_Network              *network = (DM_Network *)dm->data;
1360:   PetscInt                 offset  = 0;
1361:   DMNetworkComponentHeader header;

1363:   if (compnum == ALL_COMPONENTS) {
1364:     PetscSectionGetDof(network->DofSection, p, nvar);
1365:     return 0;
1366:   }

1368:   PetscSectionGetOffset(network->DataSection, p, &offset);
1369:   header = (DMNetworkComponentHeader)(network->componentdataarray + offset);

1371:   if (compnum >= 0) {
1372:     if (compkey) *compkey = header->key[compnum];
1373:     if (component) {
1374:       offset += header->hsize + header->offset[compnum];
1375:       *component = network->componentdataarray + offset;
1376:     }
1377:   }

1379:   if (nvar) *nvar = header->nvar[compnum];

1381:   return 0;
1382: }

1384: /*
1385:  Sets up the array that holds the data for all components and its associated section.
1386:  It copies the data for all components in a contiguous array called componentdataarray. The component data is stored pointwise with an additional header (metadata) stored for each point. The header has metadata information such as number of components at each point, number of variables for each component, offsets for the components data, etc.
1387: */
1388: PetscErrorCode DMNetworkComponentSetUp(DM dm)
1389: {
1390:   DM_Network                        *network = (DM_Network *)dm->data;
1391:   PetscInt                           arr_size, p, offset, offsetp, ncomp, i, *headerarr;
1392:   DMNetworkComponentHeader           header;
1393:   DMNetworkComponentValue            cvalue;
1394:   DMNetworkComponentHeader           headerinfo;
1395:   DMNetworkComponentGenericDataType *componentdataarray;

1397:   PetscSectionSetUp(network->DataSection);
1398:   PetscSectionGetStorageSize(network->DataSection, &arr_size);
1399:   /* arr_size+1 fixes pipeline test of opensolaris-misc for src/dm/tests/ex10.c -- Do not know why */
1400:   PetscCalloc1(arr_size + 1, &network->componentdataarray);
1401:   componentdataarray = network->componentdataarray;
1402:   for (p = network->cloneshared->pStart; p < network->cloneshared->pEnd; p++) {
1403:     PetscSectionGetOffset(network->DataSection, p, &offsetp);
1404:     /* Copy header */
1405:     header     = &network->header[p];
1406:     headerinfo = (DMNetworkComponentHeader)(componentdataarray + offsetp);
1407:     PetscMemcpy(headerinfo, header, sizeof(struct _p_DMNetworkComponentHeader));
1408:     headerarr = (PetscInt *)(headerinfo + 1);
1409:     PetscMemcpy(headerarr, header->size, header->maxcomps * sizeof(PetscInt));
1410:     headerinfo->size = headerarr;
1411:     headerarr += header->maxcomps;
1412:     PetscMemcpy(headerarr, header->key, header->maxcomps * sizeof(PetscInt));
1413:     headerinfo->key = headerarr;
1414:     headerarr += header->maxcomps;
1415:     PetscMemcpy(headerarr, header->offset, header->maxcomps * sizeof(PetscInt));
1416:     headerinfo->offset = headerarr;
1417:     headerarr += header->maxcomps;
1418:     PetscMemcpy(headerarr, header->nvar, header->maxcomps * sizeof(PetscInt));
1419:     headerinfo->nvar = headerarr;
1420:     headerarr += header->maxcomps;
1421:     PetscMemcpy(headerarr, header->offsetvarrel, header->maxcomps * sizeof(PetscInt));
1422:     headerinfo->offsetvarrel = headerarr;

1424:     /* Copy data */
1425:     cvalue = &network->cvalue[p];
1426:     ncomp  = header->ndata;

1428:     for (i = 0; i < ncomp; i++) {
1429:       offset = offsetp + header->hsize + header->offset[i];
1430:       PetscMemcpy(componentdataarray + offset, cvalue->data[i], header->size[i] * sizeof(DMNetworkComponentGenericDataType));
1431:     }
1432:   }

1434:   for (i = network->cloneshared->pStart; i < network->cloneshared->pEnd; i++) {
1435:     PetscFree5(network->header[i].size, network->header[i].key, network->header[i].offset, network->header[i].nvar, network->header[i].offsetvarrel);
1436:     PetscFree(network->cvalue[i].data);
1437:   }
1438:   PetscFree2(network->header, network->cvalue);
1439:   return 0;
1440: }

1442: /* Sets up the section for dofs. This routine is called during DMSetUp() */
1443: static PetscErrorCode DMNetworkVariablesSetUp(DM dm)
1444: {
1445:   DM_Network *network = (DM_Network *)dm->data;

1447:   PetscSectionSetUp(network->DofSection);
1448:   return 0;
1449: }

1451: /* Get a subsection from a range of points */
1452: static PetscErrorCode DMNetworkGetSubSection_private(PetscSection main, PetscInt pstart, PetscInt pend, PetscSection *subsection)
1453: {
1454:   PetscInt i, nvar;

1456:   PetscSectionCreate(PetscObjectComm((PetscObject)main), subsection);
1457:   PetscSectionSetChart(*subsection, 0, pend - pstart);
1458:   for (i = pstart; i < pend; i++) {
1459:     PetscSectionGetDof(main, i, &nvar);
1460:     PetscSectionSetDof(*subsection, i - pstart, nvar);
1461:   }

1463:   PetscSectionSetUp(*subsection);
1464:   return 0;
1465: }

1467: /* Create a submap of points with a GlobalToLocal structure */
1468: static PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
1469: {
1470:   PetscInt i, *subpoints;

1472:   /* Create index sets to map from "points" to "subpoints" */
1473:   PetscMalloc1(pend - pstart, &subpoints);
1474:   for (i = pstart; i < pend; i++) subpoints[i - pstart] = i;
1475:   ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, pend - pstart, subpoints, PETSC_COPY_VALUES, map);
1476:   PetscFree(subpoints);
1477:   return 0;
1478: }

1480: /*@
1481:   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute

1483:   Collective on dm

1485:   Input Parameters:
1486: . dm - the DMNetworkObject

1488:   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:

1490:   points = [0 1 2 3 4 5 6]

1492:   where edges = [0,1,2,3] and vertices = [4,5,6]. The new orderings will be specific to the subset (i.e vertices = [0,1,2] <- [4,5,6]).

1494:   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.

1496:   Level: intermediate

1498: @*/
1499: PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
1500: {
1501:   MPI_Comm    comm;
1502:   PetscMPIInt size;
1503:   DM_Network *network = (DM_Network *)dm->data;

1505:   PetscObjectGetComm((PetscObject)dm, &comm);
1506:   MPI_Comm_size(comm, &size);

1508:   /* Create maps for vertices and edges */
1509:   DMNetworkSetSubMap_private(network->cloneshared->vStart, network->cloneshared->vEnd, &network->vertex.mapping);
1510:   DMNetworkSetSubMap_private(network->cloneshared->eStart, network->cloneshared->eEnd, &network->edge.mapping);

1512:   /* Create local sub-sections */
1513:   DMNetworkGetSubSection_private(network->DofSection, network->cloneshared->vStart, network->cloneshared->vEnd, &network->vertex.DofSection);
1514:   DMNetworkGetSubSection_private(network->DofSection, network->cloneshared->eStart, network->cloneshared->eEnd, &network->edge.DofSection);

1516:   if (size > 1) {
1517:     PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);

1519:     PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);
1520:     PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);
1521:     PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);
1522:   } else {
1523:     /* create structures for vertex */
1524:     PetscSectionClone(network->vertex.DofSection, &network->vertex.GlobalDofSection);
1525:     /* create structures for edge */
1526:     PetscSectionClone(network->edge.DofSection, &network->edge.GlobalDofSection);
1527:   }

1529:   /* Add viewers */
1530:   PetscObjectSetName((PetscObject)network->edge.GlobalDofSection, "Global edge dof section");
1531:   PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection, "Global vertex dof section");
1532:   PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");
1533:   PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");
1534:   return 0;
1535: }

1537: /*
1538:    Setup a lookup btable for the input v's owning subnetworks
1539:    - add all owing subnetworks that connect to this v to the btable
1540:      vertex_subnetid = supportingedge_subnetid
1541: */
1542: static inline PetscErrorCode SetSubnetIdLookupBT(DM dm, PetscInt v, PetscInt Nsubnet, PetscBT btable)
1543: {
1544:   PetscInt                 e, nedges, offset;
1545:   const PetscInt          *edges;
1546:   DM_Network              *newDMnetwork = (DM_Network *)dm->data;
1547:   DMNetworkComponentHeader header;

1549:   PetscBTMemzero(Nsubnet, btable);
1550:   DMNetworkGetSupportingEdges(dm, v, &nedges, &edges);
1551:   for (e = 0; e < nedges; e++) {
1552:     PetscSectionGetOffset(newDMnetwork->DataSection, edges[e], &offset);
1553:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);
1554:     PetscBTSet(btable, header->subnetid);
1555:   }
1556:   return 0;
1557: }

1559: /*@
1560:   DMNetworkDistribute - Distributes the network and moves associated component data

1562:   Collective

1564:   Input Parameters:
1565: + DM - the DMNetwork object
1566: - overlap - the overlap of partitions, 0 is the default

1568:   Options Database Key:
1569: + -dmnetwork_view - Calls DMView() at the conclusion of DMSetUp()
1570: - -dmnetwork_view_distributed - Calls DMView() at the conclusion of DMNetworkDistribute()

1572:   Notes:
1573:   Distributes the network with <overlap>-overlapping partitioning of the edges.

1575:   Level: intermediate

1577: .seealso: `DMNetworkCreate()`
1578: @*/
1579: PetscErrorCode DMNetworkDistribute(DM *dm, PetscInt overlap)
1580: {
1581:   MPI_Comm                 comm;
1582:   PetscMPIInt              size;
1583:   DM_Network              *oldDMnetwork = (DM_Network *)((*dm)->data);
1584:   DM_Network              *newDMnetwork;
1585:   PetscSF                  pointsf = NULL;
1586:   DM                       newDM;
1587:   PetscInt                 j, e, v, offset, *subnetvtx, *subnetedge, Nsubnet, gidx, svtx_idx, nv;
1588:   PetscInt                 net, *sv;
1589:   PetscBT                  btable;
1590:   PetscPartitioner         part;
1591:   DMNetworkComponentHeader header;

1595:   PetscObjectGetComm((PetscObject)*dm, &comm);
1596:   MPI_Comm_size(comm, &size);
1597:   if (size == 1) return 0;


1601:   /* This routine moves the component data to the appropriate processors. It makes use of the DataSection and the componentdataarray to move the component data to appropriate processors and returns a new DataSection and new componentdataarray. */
1602:   PetscLogEventBegin(DMNetwork_Distribute, dm, 0, 0, 0);
1603:   DMNetworkCreate(PetscObjectComm((PetscObject)*dm), &newDM);
1604:   newDMnetwork                       = (DM_Network *)newDM->data;
1605:   newDMnetwork->max_comps_registered = oldDMnetwork->max_comps_registered;
1606:   PetscMalloc1(newDMnetwork->max_comps_registered, &newDMnetwork->component);

1608:   /* Enable runtime options for petscpartitioner */
1609:   DMPlexGetPartitioner(oldDMnetwork->plex, &part);
1610:   PetscPartitionerSetFromOptions(part);

1612:   /* Distribute plex dm */
1613:   DMPlexDistribute(oldDMnetwork->plex, overlap, &pointsf, &newDMnetwork->plex);

1615:   /* Distribute dof section */
1616:   PetscSectionCreate(comm, &newDMnetwork->DofSection);
1617:   PetscSFDistributeSection(pointsf, oldDMnetwork->DofSection, NULL, newDMnetwork->DofSection);

1619:   /* Distribute data and associated section */
1620:   PetscSectionCreate(comm, &newDMnetwork->DataSection);
1621:   DMPlexDistributeData(newDMnetwork->plex, pointsf, oldDMnetwork->DataSection, MPIU_INT, (void *)oldDMnetwork->componentdataarray, newDMnetwork->DataSection, (void **)&newDMnetwork->componentdataarray);

1623:   PetscSectionGetChart(newDMnetwork->DataSection, &newDMnetwork->cloneshared->pStart, &newDMnetwork->cloneshared->pEnd);
1624:   DMPlexGetHeightStratum(newDMnetwork->plex, 0, &newDMnetwork->cloneshared->eStart, &newDMnetwork->cloneshared->eEnd);
1625:   DMPlexGetHeightStratum(newDMnetwork->plex, 1, &newDMnetwork->cloneshared->vStart, &newDMnetwork->cloneshared->vEnd);
1626:   newDMnetwork->cloneshared->nEdges    = newDMnetwork->cloneshared->eEnd - newDMnetwork->cloneshared->eStart;
1627:   newDMnetwork->cloneshared->nVertices = newDMnetwork->cloneshared->vEnd - newDMnetwork->cloneshared->vStart;
1628:   newDMnetwork->cloneshared->NVertices = oldDMnetwork->cloneshared->NVertices;
1629:   newDMnetwork->cloneshared->NEdges    = oldDMnetwork->cloneshared->NEdges;
1630:   newDMnetwork->cloneshared->svtable   = oldDMnetwork->cloneshared->svtable; /* global table! */
1631:   oldDMnetwork->cloneshared->svtable   = NULL;

1633:   /* Set Dof section as the section for dm */
1634:   DMSetLocalSection(newDMnetwork->plex, newDMnetwork->DofSection);
1635:   DMGetGlobalSection(newDMnetwork->plex, &newDMnetwork->GlobalDofSection);

1637:   /* Setup subnetwork info in the newDM */
1638:   newDMnetwork->cloneshared->Nsubnet = oldDMnetwork->cloneshared->Nsubnet;
1639:   newDMnetwork->cloneshared->Nsvtx   = oldDMnetwork->cloneshared->Nsvtx;
1640:   oldDMnetwork->cloneshared->Nsvtx   = 0;
1641:   newDMnetwork->cloneshared->svtx    = oldDMnetwork->cloneshared->svtx; /* global vertices! */
1642:   oldDMnetwork->cloneshared->svtx    = NULL;
1643:   PetscCalloc1(newDMnetwork->cloneshared->Nsubnet, &newDMnetwork->cloneshared->subnet);

1645:   /* Copy over the global number of vertices and edges in each subnetwork.
1646:      Note: these are calculated in DMNetworkLayoutSetUp()
1647:   */
1648:   Nsubnet = newDMnetwork->cloneshared->Nsubnet;
1649:   for (j = 0; j < Nsubnet; j++) {
1650:     newDMnetwork->cloneshared->subnet[j].Nvtx  = oldDMnetwork->cloneshared->subnet[j].Nvtx;
1651:     newDMnetwork->cloneshared->subnet[j].Nedge = oldDMnetwork->cloneshared->subnet[j].Nedge;
1652:   }

1654:   /* Count local nedges for subnetworks */
1655:   for (e = newDMnetwork->cloneshared->eStart; e < newDMnetwork->cloneshared->eEnd; e++) {
1656:     PetscSectionGetOffset(newDMnetwork->DataSection, e, &offset);
1657:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);

1659:     /* Update pointers */
1660:     header->size         = (PetscInt *)(header + 1);
1661:     header->key          = header->size + header->maxcomps;
1662:     header->offset       = header->key + header->maxcomps;
1663:     header->nvar         = header->offset + header->maxcomps;
1664:     header->offsetvarrel = header->nvar + header->maxcomps;

1666:     newDMnetwork->cloneshared->subnet[header->subnetid].nedge++;
1667:   }

1669:   /* Setup a btable to keep track subnetworks owned by this process at a shared vertex */
1670:   if (newDMnetwork->cloneshared->Nsvtx) PetscBTCreate(Nsubnet, &btable);

1672:   /* Count local nvtx for subnetworks */
1673:   for (v = newDMnetwork->cloneshared->vStart; v < newDMnetwork->cloneshared->vEnd; v++) {
1674:     PetscSectionGetOffset(newDMnetwork->DataSection, v, &offset);
1675:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);

1677:     /* Update pointers */
1678:     header->size         = (PetscInt *)(header + 1);
1679:     header->key          = header->size + header->maxcomps;
1680:     header->offset       = header->key + header->maxcomps;
1681:     header->nvar         = header->offset + header->maxcomps;
1682:     header->offsetvarrel = header->nvar + header->maxcomps;

1684:     /* shared vertices: use gidx=header->index to check if v is a shared vertex */
1685:     gidx = header->index;
1686:     PetscTableFind(newDMnetwork->cloneshared->svtable, gidx + 1, &svtx_idx);
1687:     svtx_idx--;

1689:     if (svtx_idx < 0) { /* not a shared vertex */
1690:       newDMnetwork->cloneshared->subnet[header->subnetid].nvtx++;
1691:     } else { /* a shared vertex belongs to more than one subnetworks, it is being counted by multiple subnets */
1692:       /* Setup a lookup btable for this v's owning subnetworks */
1693:       SetSubnetIdLookupBT(newDM, v, Nsubnet, btable);

1695:       for (j = 0; j < newDMnetwork->cloneshared->svtx[svtx_idx].n; j++) {
1696:         sv  = newDMnetwork->cloneshared->svtx[svtx_idx].sv + 2 * j;
1697:         net = sv[0];
1698:         if (PetscBTLookup(btable, net)) newDMnetwork->cloneshared->subnet[net].nvtx++; /* sv is on net owned by this proces */
1699:       }
1700:     }
1701:   }

1703:   /* Get total local nvtx for subnetworks */
1704:   nv = 0;
1705:   for (j = 0; j < Nsubnet; j++) nv += newDMnetwork->cloneshared->subnet[j].nvtx;
1706:   nv += newDMnetwork->cloneshared->Nsvtx;

1708:   /* Now create the vertices and edge arrays for the subnetworks */
1709:   PetscCalloc2(newDMnetwork->cloneshared->nEdges, &subnetedge, nv, &subnetvtx); /* Maps local vertex to local subnetwork's vertex */
1710:   newDMnetwork->cloneshared->subnetedge = subnetedge;
1711:   newDMnetwork->cloneshared->subnetvtx  = subnetvtx;
1712:   for (j = 0; j < newDMnetwork->cloneshared->Nsubnet; j++) {
1713:     newDMnetwork->cloneshared->subnet[j].edges = subnetedge;
1714:     subnetedge += newDMnetwork->cloneshared->subnet[j].nedge;

1716:     newDMnetwork->cloneshared->subnet[j].vertices = subnetvtx;
1717:     subnetvtx += newDMnetwork->cloneshared->subnet[j].nvtx;

1719:     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. These get updated when the vertices and edges are added. */
1720:     newDMnetwork->cloneshared->subnet[j].nvtx = newDMnetwork->cloneshared->subnet[j].nedge = 0;
1721:   }
1722:   newDMnetwork->cloneshared->svertices = subnetvtx;

1724:   /* Set the edges and vertices in each subnetwork */
1725:   for (e = newDMnetwork->cloneshared->eStart; e < newDMnetwork->cloneshared->eEnd; e++) {
1726:     PetscSectionGetOffset(newDMnetwork->DataSection, e, &offset);
1727:     header                                                                                                                 = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);
1728:     newDMnetwork->cloneshared->subnet[header->subnetid].edges[newDMnetwork->cloneshared->subnet[header->subnetid].nedge++] = e;
1729:   }

1731:   nv = 0;
1732:   for (v = newDMnetwork->cloneshared->vStart; v < newDMnetwork->cloneshared->vEnd; v++) {
1733:     PetscSectionGetOffset(newDMnetwork->DataSection, v, &offset);
1734:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);

1736:     /* coupling vertices: use gidx = header->index to check if v is a coupling vertex */
1737:     PetscTableFind(newDMnetwork->cloneshared->svtable, header->index + 1, &svtx_idx);
1738:     svtx_idx--;
1739:     if (svtx_idx < 0) {
1740:       newDMnetwork->cloneshared->subnet[header->subnetid].vertices[newDMnetwork->cloneshared->subnet[header->subnetid].nvtx++] = v;
1741:     } else { /* a shared vertex */
1742:       newDMnetwork->cloneshared->svertices[nv++] = v;

1744:       /* Setup a lookup btable for this v's owning subnetworks */
1745:       SetSubnetIdLookupBT(newDM, v, Nsubnet, btable);

1747:       for (j = 0; j < newDMnetwork->cloneshared->svtx[svtx_idx].n; j++) {
1748:         sv  = newDMnetwork->cloneshared->svtx[svtx_idx].sv + 2 * j;
1749:         net = sv[0];
1750:         if (PetscBTLookup(btable, net)) newDMnetwork->cloneshared->subnet[net].vertices[newDMnetwork->cloneshared->subnet[net].nvtx++] = v;
1751:       }
1752:     }
1753:   }
1754:   newDMnetwork->cloneshared->nsvtx = nv; /* num of local shared vertices */

1756:   newDM->setupcalled                          = (*dm)->setupcalled;
1757:   newDMnetwork->cloneshared->distributecalled = PETSC_TRUE;

1759:   /* Free spaces */
1760:   PetscSFDestroy(&pointsf);
1761:   DMDestroy(dm);
1762:   if (newDMnetwork->cloneshared->Nsvtx) PetscBTDestroy(&btable);

1764:   /* View distributed dmnetwork */
1765:   DMViewFromOptions(newDM, NULL, "-dmnetwork_view_distributed");

1767:   *dm = newDM;
1768:   PetscLogEventEnd(DMNetwork_Distribute, dm, 0, 0, 0);
1769:   return 0;
1770: }

1772: /*@C
1773:   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering

1775:  Collective

1777:   Input Parameters:
1778: + mainSF - the original SF structure
1779: - map - a ISLocalToGlobal mapping that contains the subset of points

1781:   Output Parameter:
1782: . subSF - a subset of the mainSF for the desired subset.

1784:   Level: intermediate
1785: @*/
1786: PetscErrorCode PetscSFGetSubSF(PetscSF mainsf, ISLocalToGlobalMapping map, PetscSF *subSF)
1787: {
1788:   PetscInt           nroots, nleaves, *ilocal_sub;
1789:   PetscInt           i, *ilocal_map, nroots_sub, nleaves_sub = 0;
1790:   PetscInt          *local_points, *remote_points;
1791:   PetscSFNode       *iremote_sub;
1792:   const PetscInt    *ilocal;
1793:   const PetscSFNode *iremote;

1795:   PetscSFGetGraph(mainsf, &nroots, &nleaves, &ilocal, &iremote);

1797:   /* Look for leaves that pertain to the subset of points. Get the local ordering */
1798:   PetscMalloc1(nleaves, &ilocal_map);
1799:   ISGlobalToLocalMappingApply(map, IS_GTOLM_MASK, nleaves, ilocal, NULL, ilocal_map);
1800:   for (i = 0; i < nleaves; i++) {
1801:     if (ilocal_map[i] != -1) nleaves_sub += 1;
1802:   }
1803:   /* Re-number ilocal with subset numbering. Need information from roots */
1804:   PetscMalloc2(nroots, &local_points, nroots, &remote_points);
1805:   for (i = 0; i < nroots; i++) local_points[i] = i;
1806:   ISGlobalToLocalMappingApply(map, IS_GTOLM_MASK, nroots, local_points, NULL, local_points);
1807:   PetscSFBcastBegin(mainsf, MPIU_INT, local_points, remote_points, MPI_REPLACE);
1808:   PetscSFBcastEnd(mainsf, MPIU_INT, local_points, remote_points, MPI_REPLACE);
1809:   /* Fill up graph using local (that is, local to the subset) numbering. */
1810:   PetscMalloc1(nleaves_sub, &ilocal_sub);
1811:   PetscMalloc1(nleaves_sub, &iremote_sub);
1812:   nleaves_sub = 0;
1813:   for (i = 0; i < nleaves; i++) {
1814:     if (ilocal_map[i] != -1) {
1815:       ilocal_sub[nleaves_sub]        = ilocal_map[i];
1816:       iremote_sub[nleaves_sub].rank  = iremote[i].rank;
1817:       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
1818:       nleaves_sub += 1;
1819:     }
1820:   }
1821:   PetscFree2(local_points, remote_points);
1822:   ISLocalToGlobalMappingGetSize(map, &nroots_sub);

1824:   /* Create new subSF */
1825:   PetscSFCreate(PETSC_COMM_WORLD, subSF);
1826:   PetscSFSetFromOptions(*subSF);
1827:   PetscSFSetGraph(*subSF, nroots_sub, nleaves_sub, ilocal_sub, PETSC_OWN_POINTER, iremote_sub, PETSC_COPY_VALUES);
1828:   PetscFree(ilocal_map);
1829:   PetscFree(iremote_sub);
1830:   return 0;
1831: }

1833: /*@C
1834:   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point

1836:   Not Collective

1838:   Input Parameters:
1839: + dm - the DMNetwork object
1840: - p  - the vertex point

1842:   Output Parameters:
1843: + nedges - number of edges connected to this vertex point
1844: - edges  - list of edge points

1846:   Level: beginner

1848:   Fortran Notes:
1849:   Since it returns an array, this routine is only available in Fortran 90, and you must
1850:   include petsc.h90 in your code.

1852: .seealso: `DMNetworkCreate()`, `DMNetworkGetConnectedVertices()`
1853: @*/
1854: PetscErrorCode DMNetworkGetSupportingEdges(DM dm, PetscInt vertex, PetscInt *nedges, const PetscInt *edges[])
1855: {
1856:   DM_Network *network = (DM_Network *)dm->data;

1858:   DMPlexGetSupportSize(network->plex, vertex, nedges);
1859:   if (edges) DMPlexGetSupport(network->plex, vertex, edges);
1860:   return 0;
1861: }

1863: /*@C
1864:   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point

1866:   Not Collective

1868:   Input Parameters:
1869: + dm - the DMNetwork object
1870: - p - the edge point

1872:   Output Parameters:
1873: . vertices - vertices connected to this edge

1875:   Level: beginner

1877:   Fortran Notes:
1878:   Since it returns an array, this routine is only available in Fortran 90, and you must
1879:   include petsc.h90 in your code.

1881: .seealso: `DMNetworkCreate()`, `DMNetworkGetSupportingEdges()`
1882: @*/
1883: PetscErrorCode DMNetworkGetConnectedVertices(DM dm, PetscInt edge, const PetscInt *vertices[])
1884: {
1885:   DM_Network *network = (DM_Network *)dm->data;

1887:   DMPlexGetCone(network->plex, edge, vertices);
1888:   return 0;
1889: }

1891: /*@
1892:   DMNetworkIsSharedVertex - Returns TRUE if the vertex is shared by subnetworks

1894:   Not Collective

1896:   Input Parameters:
1897: + dm - the DMNetwork object
1898: - p - the vertex point

1900:   Output Parameter:
1901: . flag - TRUE if the vertex is shared by subnetworks

1903:   Level: beginner

1905: .seealso: `DMNetworkAddSharedVertices()`, `DMNetworkIsGhostVertex()`
1906: @*/
1907: PetscErrorCode DMNetworkIsSharedVertex(DM dm, PetscInt p, PetscBool *flag)
1908: {
1909:   PetscInt i;

1911:   *flag = PETSC_FALSE;

1913:   if (dm->setupcalled) { /* DMNetworkGetGlobalVertexIndex() requires DMSetUp() be called */
1914:     DM_Network *network = (DM_Network *)dm->data;
1915:     PetscInt    gidx;
1916:     DMNetworkGetGlobalVertexIndex(dm, p, &gidx);
1917:     PetscTableFind(network->cloneshared->svtable, gidx + 1, &i);
1918:     if (i) *flag = PETSC_TRUE;
1919:   } else { /* would be removed? */
1920:     PetscInt        nv;
1921:     const PetscInt *vtx;
1922:     DMNetworkGetSharedVertices(dm, &nv, &vtx);
1923:     for (i = 0; i < nv; i++) {
1924:       if (p == vtx[i]) {
1925:         *flag = PETSC_TRUE;
1926:         break;
1927:       }
1928:     }
1929:   }
1930:   return 0;
1931: }

1933: /*@
1934:   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex

1936:   Not Collective

1938:   Input Parameters:
1939: + dm - the DMNetwork object
1940: - p - the vertex point

1942:   Output Parameter:
1943: . isghost - TRUE if the vertex is a ghost point

1945:   Level: beginner

1947: .seealso: `DMNetworkGetConnectedVertices()`, `DMNetworkGetVertexRange()`, `DMNetworkIsSharedVertex()`
1948: @*/
1949: PetscErrorCode DMNetworkIsGhostVertex(DM dm, PetscInt p, PetscBool *isghost)
1950: {
1951:   DM_Network  *network = (DM_Network *)dm->data;
1952:   PetscInt     offsetg;
1953:   PetscSection sectiong;

1955:   *isghost = PETSC_FALSE;
1956:   DMGetGlobalSection(network->plex, &sectiong);
1957:   PetscSectionGetOffset(sectiong, p, &offsetg);
1958:   if (offsetg < 0) *isghost = PETSC_TRUE;
1959:   return 0;
1960: }

1962: PetscErrorCode DMSetUp_Network(DM dm)
1963: {
1964:   PetscLogEventBegin(DMNetwork_SetUpNetwork, dm, 0, 0, 0);
1965:   DMNetworkFinalizeComponents(dm);
1966:   /* View dmnetwork */
1967:   DMViewFromOptions(dm, NULL, "-dmnetwork_view");
1968:   PetscLogEventEnd(DMNetwork_SetUpNetwork, dm, 0, 0, 0);
1969:   return 0;
1970: }

1972: /*@
1973:   DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
1974:       -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?

1976:   Collective

1978:   Input Parameters:
1979: + dm - the DMNetwork object
1980: . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
1981: - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices

1983:  Level: intermediate

1985: @*/
1986: PetscErrorCode DMNetworkHasJacobian(DM dm, PetscBool eflg, PetscBool vflg)
1987: {
1988:   DM_Network *network   = (DM_Network *)dm->data;
1989:   PetscInt    nVertices = network->cloneshared->nVertices;

1991:   network->userEdgeJacobian   = eflg;
1992:   network->userVertexJacobian = vflg;

1994:   if (eflg && !network->Je) PetscCalloc1(3 * network->cloneshared->nEdges, &network->Je);

1996:   if (vflg && !network->Jv && nVertices) {
1997:     PetscInt        i, *vptr, nedges, vStart = network->cloneshared->vStart;
1998:     PetscInt        nedges_total;
1999:     const PetscInt *edges;

2001:     /* count nvertex_total */
2002:     nedges_total = 0;
2003:     PetscMalloc1(nVertices + 1, &vptr);

2005:     vptr[0] = 0;
2006:     for (i = 0; i < nVertices; i++) {
2007:       DMNetworkGetSupportingEdges(dm, i + vStart, &nedges, &edges);
2008:       nedges_total += nedges;
2009:       vptr[i + 1] = vptr[i] + 2 * nedges + 1;
2010:     }

2012:     PetscCalloc1(2 * nedges_total + nVertices, &network->Jv);
2013:     network->Jvptr = vptr;
2014:   }
2015:   return 0;
2016: }

2018: /*@
2019:   DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network

2021:   Not Collective

2023:   Input Parameters:
2024: + dm - the DMNetwork object
2025: . p - the edge point
2026: - J - array (size = 3) of Jacobian submatrices for this edge point:
2027:         J[0]: this edge
2028:         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()

2030:   Level: advanced

2032: .seealso: `DMNetworkVertexSetMatrix()`
2033: @*/
2034: PetscErrorCode DMNetworkEdgeSetMatrix(DM dm, PetscInt p, Mat J[])
2035: {
2036:   DM_Network *network = (DM_Network *)dm->data;


2040:   if (J) {
2041:     network->Je[3 * p]     = J[0];
2042:     network->Je[3 * p + 1] = J[1];
2043:     network->Je[3 * p + 2] = J[2];
2044:   }
2045:   return 0;
2046: }

2048: /*@
2049:   DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network

2051:   Not Collective

2053:   Input Parameters:
2054: + dm - The DMNetwork object
2055: . p - the vertex point
2056: - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
2057:         J[0]:       this vertex
2058:         J[1+2*i]:   i-th supporting edge
2059:         J[1+2*i+1]: i-th connected vertex

2061:   Level: advanced

2063: .seealso: `DMNetworkEdgeSetMatrix()`
2064: @*/
2065: PetscErrorCode DMNetworkVertexSetMatrix(DM dm, PetscInt p, Mat J[])
2066: {
2067:   DM_Network     *network = (DM_Network *)dm->data;
2068:   PetscInt        i, *vptr, nedges, vStart = network->cloneshared->vStart;
2069:   const PetscInt *edges;


2073:   if (J) {
2074:     vptr                          = network->Jvptr;
2075:     network->Jv[vptr[p - vStart]] = J[0]; /* Set Jacobian for this vertex */

2077:     /* Set Jacobian for each supporting edge and connected vertex */
2078:     DMNetworkGetSupportingEdges(dm, p, &nedges, &edges);
2079:     for (i = 1; i <= 2 * nedges; i++) network->Jv[vptr[p - vStart] + i] = J[i];
2080:   }
2081:   return 0;
2082: }

2084: static inline PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscBool ghost, Vec vdnz, Vec vonz)
2085: {
2086:   PetscInt    j;
2087:   PetscScalar val = (PetscScalar)ncols;

2089:   if (!ghost) {
2090:     for (j = 0; j < nrows; j++) VecSetValues(vdnz, 1, &rows[j], &val, ADD_VALUES);
2091:   } else {
2092:     for (j = 0; j < nrows; j++) VecSetValues(vonz, 1, &rows[j], &val, ADD_VALUES);
2093:   }
2094:   return 0;
2095: }

2097: static inline PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju, PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscBool ghost, Vec vdnz, Vec vonz)
2098: {
2099:   PetscInt    j, ncols_u;
2100:   PetscScalar val;

2102:   if (!ghost) {
2103:     for (j = 0; j < nrows; j++) {
2104:       MatGetRow(Ju, j, &ncols_u, NULL, NULL);
2105:       val = (PetscScalar)ncols_u;
2106:       VecSetValues(vdnz, 1, &rows[j], &val, ADD_VALUES);
2107:       MatRestoreRow(Ju, j, &ncols_u, NULL, NULL);
2108:     }
2109:   } else {
2110:     for (j = 0; j < nrows; j++) {
2111:       MatGetRow(Ju, j, &ncols_u, NULL, NULL);
2112:       val = (PetscScalar)ncols_u;
2113:       VecSetValues(vonz, 1, &rows[j], &val, ADD_VALUES);
2114:       MatRestoreRow(Ju, j, &ncols_u, NULL, NULL);
2115:     }
2116:   }
2117:   return 0;
2118: }

2120: static inline PetscErrorCode MatSetPreallocationblock_private(Mat Ju, PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscBool ghost, Vec vdnz, Vec vonz)
2121: {
2122:   if (Ju) {
2123:     MatSetPreallocationUserblock_private(Ju, nrows, rows, ncols, ghost, vdnz, vonz);
2124:   } else {
2125:     MatSetPreallocationDenseblock_private(nrows, rows, ncols, ghost, vdnz, vonz);
2126:   }
2127:   return 0;
2128: }

2130: static inline PetscErrorCode MatSetDenseblock_private(PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscInt cstart, Mat *J)
2131: {
2132:   PetscInt     j, *cols;
2133:   PetscScalar *zeros;

2135:   PetscCalloc2(ncols, &cols, nrows * ncols, &zeros);
2136:   for (j = 0; j < ncols; j++) cols[j] = j + cstart;
2137:   MatSetValues(*J, nrows, rows, ncols, cols, zeros, INSERT_VALUES);
2138:   PetscFree2(cols, zeros);
2139:   return 0;
2140: }

2142: static inline PetscErrorCode MatSetUserblock_private(Mat Ju, PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscInt cstart, Mat *J)
2143: {
2144:   PetscInt        j, M, N, row, col, ncols_u;
2145:   const PetscInt *cols;
2146:   PetscScalar     zero = 0.0;

2148:   MatGetSize(Ju, &M, &N);

2151:   for (row = 0; row < nrows; row++) {
2152:     MatGetRow(Ju, row, &ncols_u, &cols, NULL);
2153:     for (j = 0; j < ncols_u; j++) {
2154:       col = cols[j] + cstart;
2155:       MatSetValues(*J, 1, &rows[row], 1, &col, &zero, INSERT_VALUES);
2156:     }
2157:     MatRestoreRow(Ju, row, &ncols_u, &cols, NULL);
2158:   }
2159:   return 0;
2160: }

2162: static inline PetscErrorCode MatSetblock_private(Mat Ju, PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscInt cstart, Mat *J)
2163: {
2164:   if (Ju) {
2165:     MatSetUserblock_private(Ju, nrows, rows, ncols, cstart, J);
2166:   } else {
2167:     MatSetDenseblock_private(nrows, rows, ncols, cstart, J);
2168:   }
2169:   return 0;
2170: }

2172: /* Creates a GlobalToLocal mapping with a Local and Global section. This is akin to the routine DMGetLocalToGlobalMapping but without the need of providing a dm.
2173: */
2174: PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
2175: {
2176:   PetscInt  i, size, dof;
2177:   PetscInt *glob2loc;

2179:   PetscSectionGetStorageSize(localsec, &size);
2180:   PetscMalloc1(size, &glob2loc);

2182:   for (i = 0; i < size; i++) {
2183:     PetscSectionGetOffset(globalsec, i, &dof);
2184:     dof         = (dof >= 0) ? dof : -(dof + 1);
2185:     glob2loc[i] = dof;
2186:   }

2188:   ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, size, glob2loc, PETSC_OWN_POINTER, ltog);
2189: #if 0
2190:   PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);
2191: #endif
2192:   return 0;
2193: }

2195: #include <petsc/private/matimpl.h>

2197: PetscErrorCode DMCreateMatrix_Network_Nest(DM dm, Mat *J)
2198: {
2199:   DM_Network            *network = (DM_Network *)dm->data;
2200:   PetscInt               eDof, vDof;
2201:   Mat                    j11, j12, j21, j22, bA[2][2];
2202:   MPI_Comm               comm;
2203:   ISLocalToGlobalMapping eISMap, vISMap;

2205:   PetscObjectGetComm((PetscObject)dm, &comm);

2207:   PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection, &eDof);
2208:   PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection, &vDof);

2210:   MatCreate(comm, &j11);
2211:   MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);
2212:   MatSetType(j11, MATMPIAIJ);

2214:   MatCreate(comm, &j12);
2215:   MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);
2216:   MatSetType(j12, MATMPIAIJ);

2218:   MatCreate(comm, &j21);
2219:   MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);
2220:   MatSetType(j21, MATMPIAIJ);

2222:   MatCreate(comm, &j22);
2223:   MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);
2224:   MatSetType(j22, MATMPIAIJ);

2226:   bA[0][0] = j11;
2227:   bA[0][1] = j12;
2228:   bA[1][0] = j21;
2229:   bA[1][1] = j22;

2231:   CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection, network->edge.DofSection, &eISMap);
2232:   CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection, network->vertex.DofSection, &vISMap);

2234:   MatSetLocalToGlobalMapping(j11, eISMap, eISMap);
2235:   MatSetLocalToGlobalMapping(j12, eISMap, vISMap);
2236:   MatSetLocalToGlobalMapping(j21, vISMap, eISMap);
2237:   MatSetLocalToGlobalMapping(j22, vISMap, vISMap);

2239:   MatSetUp(j11);
2240:   MatSetUp(j12);
2241:   MatSetUp(j21);
2242:   MatSetUp(j22);

2244:   MatCreateNest(comm, 2, NULL, 2, NULL, &bA[0][0], J);
2245:   MatSetUp(*J);
2246:   MatNestSetVecType(*J, VECNEST);
2247:   MatDestroy(&j11);
2248:   MatDestroy(&j12);
2249:   MatDestroy(&j21);
2250:   MatDestroy(&j22);

2252:   MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);
2253:   MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);
2254:   MatSetOption(*J, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);

2256:   /* Free structures */
2257:   ISLocalToGlobalMappingDestroy(&eISMap);
2258:   ISLocalToGlobalMappingDestroy(&vISMap);
2259:   return 0;
2260: }

2262: PetscErrorCode DMCreateMatrix_Network(DM dm, Mat *J)
2263: {
2264:   DM_Network     *network = (DM_Network *)dm->data;
2265:   PetscInt        eStart, eEnd, vStart, vEnd, rstart, nrows, *rows, localSize;
2266:   PetscInt        cstart, ncols, j, e, v;
2267:   PetscBool       ghost, ghost_vc, ghost2, isNest;
2268:   Mat             Juser;
2269:   PetscSection    sectionGlobal;
2270:   PetscInt        nedges, *vptr = NULL, vc, *rows_v; /* suppress maybe-uninitialized warning */
2271:   const PetscInt *edges, *cone;
2272:   MPI_Comm        comm;
2273:   MatType         mtype;
2274:   Vec             vd_nz, vo_nz;
2275:   PetscInt       *dnnz, *onnz;
2276:   PetscScalar    *vdnz, *vonz;

2278:   mtype = dm->mattype;
2279:   PetscStrcmp(mtype, MATNEST, &isNest);
2280:   if (isNest) {
2281:     DMCreateMatrix_Network_Nest(dm, J);
2282:     MatSetDM(*J, dm);
2283:     return 0;
2284:   }

2286:   if (!network->userEdgeJacobian && !network->userVertexJacobian) {
2287:     /* user does not provide Jacobian blocks */
2288:     DMCreateMatrix_Plex(network->plex, J);
2289:     MatSetDM(*J, dm);
2290:     return 0;
2291:   }

2293:   MatCreate(PetscObjectComm((PetscObject)dm), J);
2294:   DMGetGlobalSection(network->plex, &sectionGlobal);
2295:   PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize);
2296:   MatSetSizes(*J, localSize, localSize, PETSC_DETERMINE, PETSC_DETERMINE);

2298:   MatSetType(*J, MATAIJ);
2299:   MatSetFromOptions(*J);

2301:   /* (1) Set matrix preallocation */
2302:   /*------------------------------*/
2303:   PetscObjectGetComm((PetscObject)dm, &comm);
2304:   VecCreate(comm, &vd_nz);
2305:   VecSetSizes(vd_nz, localSize, PETSC_DECIDE);
2306:   VecSetFromOptions(vd_nz);
2307:   VecSet(vd_nz, 0.0);
2308:   VecDuplicate(vd_nz, &vo_nz);

2310:   /* Set preallocation for edges */
2311:   /*-----------------------------*/
2312:   DMNetworkGetEdgeRange(dm, &eStart, &eEnd);

2314:   PetscMalloc1(localSize, &rows);
2315:   for (e = eStart; e < eEnd; e++) {
2316:     /* Get row indices */
2317:     DMNetworkGetGlobalVecOffset(dm, e, ALL_COMPONENTS, &rstart);
2318:     PetscSectionGetDof(network->DofSection, e, &nrows);
2319:     if (nrows) {
2320:       for (j = 0; j < nrows; j++) rows[j] = j + rstart;

2322:       /* Set preallocation for connected vertices */
2323:       DMNetworkGetConnectedVertices(dm, e, &cone);
2324:       for (v = 0; v < 2; v++) {
2325:         PetscSectionGetDof(network->DofSection, cone[v], &ncols);

2327:         if (network->Je) {
2328:           Juser = network->Je[3 * e + 1 + v]; /* Jacobian(e,v) */
2329:         } else Juser = NULL;
2330:         DMNetworkIsGhostVertex(dm, cone[v], &ghost);
2331:         MatSetPreallocationblock_private(Juser, nrows, rows, ncols, ghost, vd_nz, vo_nz);
2332:       }

2334:       /* Set preallocation for edge self */
2335:       cstart = rstart;
2336:       if (network->Je) {
2337:         Juser = network->Je[3 * e]; /* Jacobian(e,e) */
2338:       } else Juser = NULL;
2339:       MatSetPreallocationblock_private(Juser, nrows, rows, nrows, PETSC_FALSE, vd_nz, vo_nz);
2340:     }
2341:   }

2343:   /* Set preallocation for vertices */
2344:   /*--------------------------------*/
2345:   DMNetworkGetVertexRange(dm, &vStart, &vEnd);
2346:   if (vEnd - vStart) vptr = network->Jvptr;

2348:   for (v = vStart; v < vEnd; v++) {
2349:     /* Get row indices */
2350:     DMNetworkGetGlobalVecOffset(dm, v, ALL_COMPONENTS, &rstart);
2351:     PetscSectionGetDof(network->DofSection, v, &nrows);
2352:     if (!nrows) continue;

2354:     DMNetworkIsGhostVertex(dm, v, &ghost);
2355:     if (ghost) {
2356:       PetscMalloc1(nrows, &rows_v);
2357:     } else {
2358:       rows_v = rows;
2359:     }

2361:     for (j = 0; j < nrows; j++) rows_v[j] = j + rstart;

2363:     /* Get supporting edges and connected vertices */
2364:     DMNetworkGetSupportingEdges(dm, v, &nedges, &edges);

2366:     for (e = 0; e < nedges; e++) {
2367:       /* Supporting edges */
2368:       DMNetworkGetGlobalVecOffset(dm, edges[e], ALL_COMPONENTS, &cstart);
2369:       PetscSectionGetDof(network->DofSection, edges[e], &ncols);

2371:       if (network->Jv) {
2372:         Juser = network->Jv[vptr[v - vStart] + 2 * e + 1]; /* Jacobian(v,e) */
2373:       } else Juser = NULL;
2374:       MatSetPreallocationblock_private(Juser, nrows, rows_v, ncols, ghost, vd_nz, vo_nz);

2376:       /* Connected vertices */
2377:       DMNetworkGetConnectedVertices(dm, edges[e], &cone);
2378:       vc = (v == cone[0]) ? cone[1] : cone[0];
2379:       DMNetworkIsGhostVertex(dm, vc, &ghost_vc);

2381:       PetscSectionGetDof(network->DofSection, vc, &ncols);

2383:       if (network->Jv) {
2384:         Juser = network->Jv[vptr[v - vStart] + 2 * e + 2]; /* Jacobian(v,vc) */
2385:       } else Juser = NULL;
2386:       if (ghost_vc || ghost) {
2387:         ghost2 = PETSC_TRUE;
2388:       } else {
2389:         ghost2 = PETSC_FALSE;
2390:       }
2391:       MatSetPreallocationblock_private(Juser, nrows, rows_v, ncols, ghost2, vd_nz, vo_nz);
2392:     }

2394:     /* Set preallocation for vertex self */
2395:     DMNetworkIsGhostVertex(dm, v, &ghost);
2396:     if (!ghost) {
2397:       DMNetworkGetGlobalVecOffset(dm, v, ALL_COMPONENTS, &cstart);
2398:       if (network->Jv) {
2399:         Juser = network->Jv[vptr[v - vStart]]; /* Jacobian(v,v) */
2400:       } else Juser = NULL;
2401:       MatSetPreallocationblock_private(Juser, nrows, rows_v, nrows, PETSC_FALSE, vd_nz, vo_nz);
2402:     }
2403:     if (ghost) PetscFree(rows_v);
2404:   }

2406:   VecAssemblyBegin(vd_nz);
2407:   VecAssemblyBegin(vo_nz);

2409:   PetscMalloc2(localSize, &dnnz, localSize, &onnz);

2411:   VecAssemblyEnd(vd_nz);
2412:   VecAssemblyEnd(vo_nz);

2414:   VecGetArray(vd_nz, &vdnz);
2415:   VecGetArray(vo_nz, &vonz);
2416:   for (j = 0; j < localSize; j++) {
2417:     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
2418:     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
2419:   }
2420:   VecRestoreArray(vd_nz, &vdnz);
2421:   VecRestoreArray(vo_nz, &vonz);
2422:   VecDestroy(&vd_nz);
2423:   VecDestroy(&vo_nz);

2425:   MatSeqAIJSetPreallocation(*J, 0, dnnz);
2426:   MatMPIAIJSetPreallocation(*J, 0, dnnz, 0, onnz);
2427:   MatSetOption(*J, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE);

2429:   PetscFree2(dnnz, onnz);

2431:   /* (2) Set matrix entries for edges */
2432:   /*----------------------------------*/
2433:   for (e = eStart; e < eEnd; e++) {
2434:     /* Get row indices */
2435:     DMNetworkGetGlobalVecOffset(dm, e, ALL_COMPONENTS, &rstart);
2436:     PetscSectionGetDof(network->DofSection, e, &nrows);
2437:     if (nrows) {
2438:       for (j = 0; j < nrows; j++) rows[j] = j + rstart;

2440:       /* Set matrix entries for connected vertices */
2441:       DMNetworkGetConnectedVertices(dm, e, &cone);
2442:       for (v = 0; v < 2; v++) {
2443:         DMNetworkGetGlobalVecOffset(dm, cone[v], ALL_COMPONENTS, &cstart);
2444:         PetscSectionGetDof(network->DofSection, cone[v], &ncols);

2446:         if (network->Je) {
2447:           Juser = network->Je[3 * e + 1 + v]; /* Jacobian(e,v) */
2448:         } else Juser = NULL;
2449:         MatSetblock_private(Juser, nrows, rows, ncols, cstart, J);
2450:       }

2452:       /* Set matrix entries for edge self */
2453:       cstart = rstart;
2454:       if (network->Je) {
2455:         Juser = network->Je[3 * e]; /* Jacobian(e,e) */
2456:       } else Juser = NULL;
2457:       MatSetblock_private(Juser, nrows, rows, nrows, cstart, J);
2458:     }
2459:   }

2461:   /* Set matrix entries for vertices */
2462:   /*---------------------------------*/
2463:   for (v = vStart; v < vEnd; v++) {
2464:     /* Get row indices */
2465:     DMNetworkGetGlobalVecOffset(dm, v, ALL_COMPONENTS, &rstart);
2466:     PetscSectionGetDof(network->DofSection, v, &nrows);
2467:     if (!nrows) continue;

2469:     DMNetworkIsGhostVertex(dm, v, &ghost);
2470:     if (ghost) {
2471:       PetscMalloc1(nrows, &rows_v);
2472:     } else {
2473:       rows_v = rows;
2474:     }
2475:     for (j = 0; j < nrows; j++) rows_v[j] = j + rstart;

2477:     /* Get supporting edges and connected vertices */
2478:     DMNetworkGetSupportingEdges(dm, v, &nedges, &edges);

2480:     for (e = 0; e < nedges; e++) {
2481:       /* Supporting edges */
2482:       DMNetworkGetGlobalVecOffset(dm, edges[e], ALL_COMPONENTS, &cstart);
2483:       PetscSectionGetDof(network->DofSection, edges[e], &ncols);

2485:       if (network->Jv) {
2486:         Juser = network->Jv[vptr[v - vStart] + 2 * e + 1]; /* Jacobian(v,e) */
2487:       } else Juser = NULL;
2488:       MatSetblock_private(Juser, nrows, rows_v, ncols, cstart, J);

2490:       /* Connected vertices */
2491:       DMNetworkGetConnectedVertices(dm, edges[e], &cone);
2492:       vc = (v == cone[0]) ? cone[1] : cone[0];

2494:       DMNetworkGetGlobalVecOffset(dm, vc, ALL_COMPONENTS, &cstart);
2495:       PetscSectionGetDof(network->DofSection, vc, &ncols);

2497:       if (network->Jv) {
2498:         Juser = network->Jv[vptr[v - vStart] + 2 * e + 2]; /* Jacobian(v,vc) */
2499:       } else Juser = NULL;
2500:       MatSetblock_private(Juser, nrows, rows_v, ncols, cstart, J);
2501:     }

2503:     /* Set matrix entries for vertex self */
2504:     if (!ghost) {
2505:       DMNetworkGetGlobalVecOffset(dm, v, ALL_COMPONENTS, &cstart);
2506:       if (network->Jv) {
2507:         Juser = network->Jv[vptr[v - vStart]]; /* Jacobian(v,v) */
2508:       } else Juser = NULL;
2509:       MatSetblock_private(Juser, nrows, rows_v, nrows, cstart, J);
2510:     }
2511:     if (ghost) PetscFree(rows_v);
2512:   }
2513:   PetscFree(rows);

2515:   MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);
2516:   MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);

2518:   MatSetDM(*J, dm);
2519:   return 0;
2520: }

2522: static PetscErrorCode DMNetworkDestroyComponentData(DM dm)
2523: {
2524:   DM_Network *network = (DM_Network *)dm->data;
2525:   PetscInt    j, np;

2527:   if (network->header) {
2528:     np = network->cloneshared->pEnd - network->cloneshared->pStart;
2529:     for (j = 0; j < np; j++) {
2530:       PetscFree5(network->header[j].size, network->header[j].key, network->header[j].offset, network->header[j].nvar, network->header[j].offsetvarrel);
2531:       PetscFree(network->cvalue[j].data);
2532:     }
2533:     PetscFree2(network->header, network->cvalue);
2534:   }
2535:   return 0;
2536: }

2538: PetscErrorCode DMDestroy_Network(DM dm)
2539: {
2540:   DM_Network *network = (DM_Network *)dm->data;
2541:   PetscInt    j;

2543:   /*
2544:     Developers Note: Due to the mixed topological definition of DMNetwork and data defined ON the
2545:     network like DofSection, DataSection, *componentdataarray, and friends, when cloning, we share
2546:     only the true topological data, and make our own data ON the network. Thus refct only refers
2547:     to the number of references to topological data, and data ON the network is always destroyed.
2548:     It is understood this is atypical for a DM, but is very intentional.
2549:   */

2551:   /* Always destroy data ON the network */
2552:   PetscFree(network->Je);
2553:   if (network->Jv) {
2554:     PetscFree(network->Jvptr);
2555:     PetscFree(network->Jv);
2556:   }
2557:   PetscSectionDestroy(&network->DataSection);
2558:   PetscSectionDestroy(&network->DofSection);
2559:   PetscFree(network->component);
2560:   PetscFree(network->componentdataarray);
2561:   DMNetworkDestroyComponentData(dm);

2563:   DMDestroy(&network->plex); /* this is cloned in DMClone_Network, so safe to destroy */

2565:   /*
2566:   Developers Note: The DMNetworkVertexInfo and DMNetworkEdgeInfo data structures are completely
2567:   destroyed as they are again a mix of topological data:
2568:     ISLocalToGlobalMapping            mapping;
2569:     PetscSF                           sf;
2570:   and data ON the network
2571:     PetscSection                      DofSection;
2572:     PetscSection                      GlobalDofSection;
2573:   And the only way to access them currently is through DMNetworkAssembleGraphStructures which assembles
2574:   everything. So we must destroy everything and require DMNetworkAssembleGraphStructures is called again
2575:   for any clone.
2576:   */
2577:   ISLocalToGlobalMappingDestroy(&network->vertex.mapping);
2578:   PetscSectionDestroy(&network->vertex.DofSection);
2579:   PetscSectionDestroy(&network->vertex.GlobalDofSection);
2580:   PetscSFDestroy(&network->vertex.sf);
2581:   /* edge */
2582:   ISLocalToGlobalMappingDestroy(&network->edge.mapping);
2583:   PetscSectionDestroy(&network->edge.DofSection);
2584:   PetscSectionDestroy(&network->edge.GlobalDofSection);
2585:   PetscSFDestroy(&network->edge.sf);
2586:   /* Destroy the potentially cloneshared data */
2587:   if (--network->cloneshared->refct <= 0) {
2588:     /* Developer Note: I'm not sure if vltog can be reused or not, as I'm not sure what it's purpose is. I
2589:      naively think it can be reused. */
2590:     PetscFree(network->cloneshared->vltog);
2591:     for (j = 0; j < network->cloneshared->Nsvtx; j++) PetscFree(network->cloneshared->svtx[j].sv);
2592:     PetscFree(network->cloneshared->svtx);
2593:     PetscFree2(network->cloneshared->subnetedge, network->cloneshared->subnetvtx);
2594:     PetscTableDestroy(&network->cloneshared->svtable);
2595:     PetscFree(network->cloneshared->subnet);
2596:     PetscFree(network->cloneshared);
2597:   }
2598:   PetscFree(network); /* Always freed as this structure is copied in a clone, not cloneshared */
2599:   return 0;
2600: }

2602: PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
2603: {
2604:   PetscBool   iascii;
2605:   PetscMPIInt rank;

2609:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
2610:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);
2611:   if (iascii) {
2612:     const PetscInt *cone, *vtx, *edges;
2613:     PetscInt        vfrom, vto, i, j, nv, ne, nsv, p, nsubnet;
2614:     DM_Network     *network = (DM_Network *)dm->data;

2616:     nsubnet = network->cloneshared->Nsubnet; /* num of subnetworks */
2617:     if (rank == 0) {
2618:       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  NSubnets: %" PetscInt_FMT "; NEdges: %" PetscInt_FMT "; NVertices: %" PetscInt_FMT "; NSharedVertices: %" PetscInt_FMT ".\n", nsubnet, network->cloneshared->NEdges, network->cloneshared->NVertices,
2619:                             network->cloneshared->Nsvtx));
2620:     }

2622:     DMNetworkGetSharedVertices(dm, &nsv, NULL);
2623:     PetscViewerASCIIPushSynchronized(viewer);
2624:     PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] nEdges: %" PetscInt_FMT "; nVertices: %" PetscInt_FMT "; nSharedVertices: %" PetscInt_FMT "\n", rank, network->cloneshared->nEdges, network->cloneshared->nVertices, nsv);

2626:     for (i = 0; i < nsubnet; i++) {
2627:       DMNetworkGetSubnetwork(dm, i, &nv, &ne, &vtx, &edges);
2628:       if (ne) {
2629:         PetscViewerASCIISynchronizedPrintf(viewer, "     Subnet %" PetscInt_FMT ": nEdges %" PetscInt_FMT ", nVertices(include shared vertices) %" PetscInt_FMT "\n", i, ne, nv);
2630:         for (j = 0; j < ne; j++) {
2631:           p = edges[j];
2632:           DMNetworkGetConnectedVertices(dm, p, &cone);
2633:           DMNetworkGetGlobalVertexIndex(dm, cone[0], &vfrom);
2634:           DMNetworkGetGlobalVertexIndex(dm, cone[1], &vto);
2635:           DMNetworkGetGlobalEdgeIndex(dm, edges[j], &p);
2636:           PetscViewerASCIISynchronizedPrintf(viewer, "       edge %" PetscInt_FMT ": %" PetscInt_FMT " ----> %" PetscInt_FMT "\n", p, vfrom, vto);
2637:         }
2638:       }
2639:     }

2641:     /* Shared vertices */
2642:     DMNetworkGetSharedVertices(dm, NULL, &vtx);
2643:     if (nsv) {
2644:       PetscInt        gidx;
2645:       PetscBool       ghost;
2646:       const PetscInt *sv = NULL;

2648:       PetscViewerASCIISynchronizedPrintf(viewer, "     SharedVertices:\n");
2649:       for (i = 0; i < nsv; i++) {
2650:         DMNetworkIsGhostVertex(dm, vtx[i], &ghost);
2651:         if (ghost) continue;

2653:         DMNetworkSharedVertexGetInfo(dm, vtx[i], &gidx, &nv, &sv);
2654:         PetscViewerASCIISynchronizedPrintf(viewer, "       svtx %" PetscInt_FMT ": global index %" PetscInt_FMT ", subnet[%" PetscInt_FMT "].%" PetscInt_FMT " ---->\n", i, gidx, sv[0], sv[1]);
2655:         for (j = 1; j < nv; j++) PetscViewerASCIISynchronizedPrintf(viewer, "                                           ----> subnet[%" PetscInt_FMT "].%" PetscInt_FMT "\n", sv[2 * j], sv[2 * j + 1]);
2656:       }
2657:     }
2658:     PetscViewerFlush(viewer);
2659:     PetscViewerASCIIPopSynchronized(viewer);
2661:   return 0;
2662: }

2664: PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
2665: {
2666:   DM_Network *network = (DM_Network *)dm->data;

2668:   DMGlobalToLocalBegin(network->plex, g, mode, l);
2669:   return 0;
2670: }

2672: PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
2673: {
2674:   DM_Network *network = (DM_Network *)dm->data;

2676:   DMGlobalToLocalEnd(network->plex, g, mode, l);
2677:   return 0;
2678: }

2680: PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
2681: {
2682:   DM_Network *network = (DM_Network *)dm->data;

2684:   DMLocalToGlobalBegin(network->plex, l, mode, g);
2685:   return 0;
2686: }

2688: PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
2689: {
2690:   DM_Network *network = (DM_Network *)dm->data;

2692:   DMLocalToGlobalEnd(network->plex, l, mode, g);
2693:   return 0;
2694: }

2696: /*@
2697:   DMNetworkGetVertexLocalToGlobalOrdering - Get vertex global index

2699:   Not collective

2701:   Input Parameters:
2702: + dm - the dm object
2703: - vloc - local vertex ordering, start from 0

2705:   Output Parameters:
2706: .  vg  - global vertex ordering, start from 0

2708:   Level: advanced

2710: .seealso: `DMNetworkSetVertexLocalToGlobalOrdering()`
2711: @*/
2712: PetscErrorCode DMNetworkGetVertexLocalToGlobalOrdering(DM dm, PetscInt vloc, PetscInt *vg)
2713: {
2714:   DM_Network *network = (DM_Network *)dm->data;
2715:   PetscInt   *vltog   = network->cloneshared->vltog;

2718:   *vg = vltog[vloc];
2719:   return 0;
2720: }

2722: /*@
2723:   DMNetworkSetVertexLocalToGlobalOrdering - Create and setup vertex local to global map

2725:   Collective

2727:   Input Parameters:
2728: . dm - the dm object

2730:   Level: advanced

2732: .seealso: `DMNetworkGetGlobalVertexIndex()`
2733: @*/
2734: PetscErrorCode DMNetworkSetVertexLocalToGlobalOrdering(DM dm)
2735: {
2736:   DM_Network        *network = (DM_Network *)dm->data;
2737:   MPI_Comm           comm;
2738:   PetscMPIInt        rank, size, *displs = NULL, *recvcounts = NULL, remoterank;
2739:   PetscBool          ghost;
2740:   PetscInt          *vltog, nroots, nleaves, i, *vrange, k, N, lidx;
2741:   const PetscSFNode *iremote;
2742:   PetscSF            vsf;
2743:   Vec                Vleaves, Vleaves_seq;
2744:   VecScatter         ctx;
2745:   PetscScalar       *varr, val;
2746:   const PetscScalar *varr_read;

2748:   PetscObjectGetComm((PetscObject)dm, &comm);
2749:   MPI_Comm_size(comm, &size);
2750:   MPI_Comm_rank(comm, &rank);

2752:   if (size == 1) {
2753:     nroots = network->cloneshared->vEnd - network->cloneshared->vStart;
2754:     PetscMalloc1(nroots, &vltog);
2755:     for (i = 0; i < nroots; i++) vltog[i] = i;
2756:     network->cloneshared->vltog = vltog;
2757:     return 0;
2758:   }

2761:   if (network->cloneshared->vltog) PetscFree(network->cloneshared->vltog);

2763:   DMNetworkSetSubMap_private(network->cloneshared->vStart, network->cloneshared->vEnd, &network->vertex.mapping);
2764:   PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);
2765:   vsf = network->vertex.sf;

2767:   PetscMalloc3(size + 1, &vrange, size, &displs, size, &recvcounts);
2768:   PetscSFGetGraph(vsf, &nroots, &nleaves, NULL, &iremote);

2770:   for (i = 0; i < size; i++) {
2771:     displs[i]     = i;
2772:     recvcounts[i] = 1;
2773:   }

2775:   i         = nroots - nleaves; /* local number of vertices, excluding ghosts */
2776:   vrange[0] = 0;
2777:   MPI_Allgatherv(&i, 1, MPIU_INT, vrange + 1, recvcounts, displs, MPIU_INT, comm);
2778:   for (i = 2; i < size + 1; i++) vrange[i] += vrange[i - 1];

2780:   PetscMalloc1(nroots, &vltog);
2781:   network->cloneshared->vltog = vltog;

2783:   /* Set vltog for non-ghost vertices */
2784:   k = 0;
2785:   for (i = 0; i < nroots; i++) {
2786:     DMNetworkIsGhostVertex(dm, i + network->cloneshared->vStart, &ghost);
2787:     if (ghost) continue;
2788:     vltog[i] = vrange[rank] + k++;
2789:   }
2790:   PetscFree3(vrange, displs, recvcounts);

2792:   /* Set vltog for ghost vertices */
2793:   /* (a) create parallel Vleaves and sequential Vleaves_seq to convert local iremote[*].index to global index */
2794:   VecCreate(comm, &Vleaves);
2795:   VecSetSizes(Vleaves, 2 * nleaves, PETSC_DETERMINE);
2796:   VecSetFromOptions(Vleaves);
2797:   VecGetArray(Vleaves, &varr);
2798:   for (i = 0; i < nleaves; i++) {
2799:     varr[2 * i]     = (PetscScalar)(iremote[i].rank);  /* rank of remote process */
2800:     varr[2 * i + 1] = (PetscScalar)(iremote[i].index); /* local index in remote process */
2801:   }
2802:   VecRestoreArray(Vleaves, &varr);

2804:   /* (b) scatter local info to remote processes via VecScatter() */
2805:   VecScatterCreateToAll(Vleaves, &ctx, &Vleaves_seq);
2806:   VecScatterBegin(ctx, Vleaves, Vleaves_seq, INSERT_VALUES, SCATTER_FORWARD);
2807:   VecScatterEnd(ctx, Vleaves, Vleaves_seq, INSERT_VALUES, SCATTER_FORWARD);

2809:   /* (c) convert local indices to global indices in parallel vector Vleaves */
2810:   VecGetSize(Vleaves_seq, &N);
2811:   VecGetArrayRead(Vleaves_seq, &varr_read);
2812:   for (i = 0; i < N; i += 2) {
2813:     remoterank = (PetscMPIInt)PetscRealPart(varr_read[i]);
2814:     if (remoterank == rank) {
2815:       k    = i + 1; /* row number */
2816:       lidx = (PetscInt)PetscRealPart(varr_read[i + 1]);
2817:       val  = (PetscScalar)vltog[lidx]; /* global index for non-ghost vertex computed above */
2818:       VecSetValues(Vleaves, 1, &k, &val, INSERT_VALUES);
2819:     }
2820:   }
2821:   VecRestoreArrayRead(Vleaves_seq, &varr_read);
2822:   VecAssemblyBegin(Vleaves);
2823:   VecAssemblyEnd(Vleaves);

2825:   /* (d) Set vltog for ghost vertices by copying local values of Vleaves */
2826:   VecGetArrayRead(Vleaves, &varr_read);
2827:   k = 0;
2828:   for (i = 0; i < nroots; i++) {
2829:     DMNetworkIsGhostVertex(dm, i + network->cloneshared->vStart, &ghost);
2830:     if (!ghost) continue;
2831:     vltog[i] = (PetscInt)PetscRealPart(varr_read[2 * k + 1]);
2832:     k++;
2833:   }
2834:   VecRestoreArrayRead(Vleaves, &varr_read);

2836:   VecDestroy(&Vleaves);
2837:   VecDestroy(&Vleaves_seq);
2838:   VecScatterDestroy(&ctx);
2839:   return 0;
2840: }

2842: static inline PetscErrorCode DMISAddSize_private(DM_Network *network, PetscInt p, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *nidx)
2843: {
2844:   PetscInt                 i, j, ncomps, nvar, key, offset = 0;
2845:   DMNetworkComponentHeader header;

2847:   PetscSectionGetOffset(network->DataSection, p, &offset);
2848:   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray + offset))->ndata;
2849:   header = (DMNetworkComponentHeader)(network->componentdataarray + offset);

2851:   for (i = 0; i < ncomps; i++) {
2852:     key  = header->key[i];
2853:     nvar = header->nvar[i];
2854:     for (j = 0; j < numkeys; j++) {
2855:       if (key == keys[j]) {
2856:         if (!blocksize || blocksize[j] == -1) {
2857:           *nidx += nvar;
2858:         } else {
2859:           *nidx += nselectedvar[j] * nvar / blocksize[j];
2860:         }
2861:       }
2862:     }
2863:   }
2864:   return 0;
2865: }

2867: static inline PetscErrorCode DMISComputeIdx_private(DM dm, PetscInt p, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *selectedvar[], PetscInt *ii, PetscInt *idx)
2868: {
2869:   PetscInt                 i, j, ncomps, nvar, key, offsetg, k, k1, offset = 0;
2870:   DM_Network              *network = (DM_Network *)dm->data;
2871:   DMNetworkComponentHeader header;

2873:   PetscSectionGetOffset(network->DataSection, p, &offset);
2874:   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray + offset))->ndata;
2875:   header = (DMNetworkComponentHeader)(network->componentdataarray + offset);

2877:   for (i = 0; i < ncomps; i++) {
2878:     key  = header->key[i];
2879:     nvar = header->nvar[i];
2880:     for (j = 0; j < numkeys; j++) {
2881:       if (key != keys[j]) continue;

2883:       DMNetworkGetGlobalVecOffset(dm, p, i, &offsetg);
2884:       if (!blocksize || blocksize[j] == -1) {
2885:         for (k = 0; k < nvar; k++) idx[(*ii)++] = offsetg + k;
2886:       } else {
2887:         for (k = 0; k < nvar; k += blocksize[j]) {
2888:           for (k1 = 0; k1 < nselectedvar[j]; k1++) idx[(*ii)++] = offsetg + k + selectedvar[j][k1];
2889:         }
2890:       }
2891:     }
2892:   }
2893:   return 0;
2894: }

2896: /*@
2897:   DMNetworkCreateIS - Create an index set object from the global vector of the network

2899:   Collective

2901:   Input Parameters:
2902: + dm - DMNetwork object
2903: . numkeys - number of keys
2904: . keys - array of keys that define the components of the variables you wish to extract
2905: . blocksize - block size of the variables associated to the component
2906: . nselectedvar - number of variables in each block to select
2907: - selectedvar - the offset into the block of each variable in each block to select

2909:   Output Parameters:
2910: . is - the index set

2912:   Level: Advanced

2914:   Notes:
2915:     Use blocksize[i] of -1 to indicate select all the variables of the i-th component, for which nselectvar[i] and selectedvar[i] are ignored. Use NULL, NULL, NULL to indicate for all selected components one wishes to obtain all the values of that component. For example, DMNetworkCreateIS(dm,1,&key,NULL,NULL,NULL,&is) will return an is that extracts all the variables for the 0-th component.

2917: .seealso: `DMNetworkCreate()`, `ISCreateGeneral()`, `DMNetworkCreateLocalIS()`
2918: @*/
2919: PetscErrorCode DMNetworkCreateIS(DM dm, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *selectedvar[], IS *is)
2920: {
2921:   MPI_Comm    comm;
2922:   DM_Network *network = (DM_Network *)dm->data;
2923:   PetscInt    i, p, estart, eend, vstart, vend, nidx, *idx;
2924:   PetscBool   ghost;

2926:   PetscObjectGetComm((PetscObject)dm, &comm);

2928:   /* Check input parameters */
2929:   for (i = 0; i < numkeys; i++) {
2930:     if (!blocksize || blocksize[i] == -1) continue;
2932:   }

2934:   DMNetworkGetEdgeRange(dm, &estart, &eend);
2935:   DMNetworkGetVertexRange(dm, &vstart, &vend);

2937:   /* Get local number of idx */
2938:   nidx = 0;
2939:   for (p = estart; p < eend; p++) DMISAddSize_private(network, p, numkeys, keys, blocksize, nselectedvar, &nidx);
2940:   for (p = vstart; p < vend; p++) {
2941:     DMNetworkIsGhostVertex(dm, p, &ghost);
2942:     if (ghost) continue;
2943:     DMISAddSize_private(network, p, numkeys, keys, blocksize, nselectedvar, &nidx);
2944:   }

2946:   /* Compute idx */
2947:   PetscMalloc1(nidx, &idx);
2948:   i = 0;
2949:   for (p = estart; p < eend; p++) DMISComputeIdx_private(dm, p, numkeys, keys, blocksize, nselectedvar, selectedvar, &i, idx);
2950:   for (p = vstart; p < vend; p++) {
2951:     DMNetworkIsGhostVertex(dm, p, &ghost);
2952:     if (ghost) continue;
2953:     DMISComputeIdx_private(dm, p, numkeys, keys, blocksize, nselectedvar, selectedvar, &i, idx);
2954:   }

2956:   /* Create is */
2957:   ISCreateGeneral(comm, nidx, idx, PETSC_COPY_VALUES, is);
2958:   PetscFree(idx);
2959:   return 0;
2960: }

2962: static inline PetscErrorCode DMISComputeLocalIdx_private(DM dm, PetscInt p, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *selectedvar[], PetscInt *ii, PetscInt *idx)
2963: {
2964:   PetscInt                 i, j, ncomps, nvar, key, offsetl, k, k1, offset = 0;
2965:   DM_Network              *network = (DM_Network *)dm->data;
2966:   DMNetworkComponentHeader header;

2968:   PetscSectionGetOffset(network->DataSection, p, &offset);
2969:   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray + offset))->ndata;
2970:   header = (DMNetworkComponentHeader)(network->componentdataarray + offset);

2972:   for (i = 0; i < ncomps; i++) {
2973:     key  = header->key[i];
2974:     nvar = header->nvar[i];
2975:     for (j = 0; j < numkeys; j++) {
2976:       if (key != keys[j]) continue;

2978:       DMNetworkGetLocalVecOffset(dm, p, i, &offsetl);
2979:       if (!blocksize || blocksize[j] == -1) {
2980:         for (k = 0; k < nvar; k++) idx[(*ii)++] = offsetl + k;
2981:       } else {
2982:         for (k = 0; k < nvar; k += blocksize[j]) {
2983:           for (k1 = 0; k1 < nselectedvar[j]; k1++) idx[(*ii)++] = offsetl + k + selectedvar[j][k1];
2984:         }
2985:       }
2986:     }
2987:   }
2988:   return 0;
2989: }

2991: /*@
2992:   DMNetworkCreateLocalIS - Create an index set object from the local vector of the network

2994:   Not collective

2996:   Input Parameters:
2997: + dm - DMNetwork object
2998: . numkeys - number of keys
2999: . keys - array of keys that define the components of the variables you wish to extract
3000: . blocksize - block size of the variables associated to the component
3001: . nselectedvar - number of variables in each block to select
3002: - selectedvar - the offset into the block of each variable in each block to select

3004:   Output Parameters:
3005: . is - the index set

3007:   Level: Advanced

3009:   Notes:
3010:     Use blocksize[i] of -1 to indicate select all the variables of the i-th component, for which nselectvar[i] and selectedvar[i] are ignored. Use NULL, NULL, NULL to indicate for all selected components one wishes to obtain all the values of that component. For example, DMNetworkCreateLocalIS(dm,1,&key,NULL,NULL,NULL,&is) will return an is that extracts all the variables for the 0-th component.

3012: .seealso: `DMNetworkCreate()`, `DMNetworkCreateIS`, `ISCreateGeneral()`
3013: @*/
3014: PetscErrorCode DMNetworkCreateLocalIS(DM dm, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *selectedvar[], IS *is)
3015: {
3016:   DM_Network *network = (DM_Network *)dm->data;
3017:   PetscInt    i, p, pstart, pend, nidx, *idx;

3019:   /* Check input parameters */
3020:   for (i = 0; i < numkeys; i++) {
3021:     if (!blocksize || blocksize[i] == -1) continue;
3023:   }

3025:   pstart = network->cloneshared->pStart;
3026:   pend   = network->cloneshared->pEnd;

3028:   /* Get local number of idx */
3029:   nidx = 0;
3030:   for (p = pstart; p < pend; p++) DMISAddSize_private(network, p, numkeys, keys, blocksize, nselectedvar, &nidx);

3032:   /* Compute local idx */
3033:   PetscMalloc1(nidx, &idx);
3034:   i = 0;
3035:   for (p = pstart; p < pend; p++) DMISComputeLocalIdx_private(dm, p, numkeys, keys, blocksize, nselectedvar, selectedvar, &i, idx);

3037:   /* Create is */
3038:   ISCreateGeneral(PETSC_COMM_SELF, nidx, idx, PETSC_COPY_VALUES, is);
3039:   PetscFree(idx);
3040:   return 0;
3041: }
3042: /*@
3043:   DMNetworkFinalizeComponents - Sets up internal data structures for the sections and components. It is called after registering new components and adding all components
3044:   to the cloned network. After calling this subroutine, no new components can be added to the network.

3046:   Collective

3048:   Input Parameters:
3049: . dm - the dmnetwork object

3051:   Level: beginner

3053: .seealso: `DMNetworkAddComponent()`, `DMNetworkRegisterComponent()`, `DMSetUp()`
3054: @*/
3055: PetscErrorCode DMNetworkFinalizeComponents(DM dm)
3056: {
3057:   DM_Network *network = (DM_Network *)dm->data;

3059:   if (network->componentsetup) return 0;
3060:   DMNetworkComponentSetUp(dm);
3061:   DMNetworkVariablesSetUp(dm);
3062:   DMSetLocalSection(network->plex, network->DofSection);
3063:   DMGetGlobalSection(network->plex, &network->GlobalDofSection);
3064:   network->componentsetup = PETSC_TRUE;
3065:   return 0;
3066: }