Actual source code: network.c

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

  3: /*
  4:   Creates the component header and value objects for a network point
  5: */
  6: static PetscErrorCode SetUpNetworkHeaderComponentValue(DM dm,DMNetworkComponentHeader header,DMNetworkComponentValue cvalue)
  7: {

 11:   /* Allocate arrays for component information */
 12:   PetscCalloc5(header->maxcomps,&header->size,header->maxcomps,&header->key,header->maxcomps,&header->offset,header->maxcomps,&header->nvar,header->maxcomps,&header->offsetvarrel);
 13:   PetscCalloc1(header->maxcomps,&cvalue->data);

 15:   /* 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.
 16:    If the data header struct changes then this header size calculation needs to be updated. */
 17:   header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5*header->maxcomps*sizeof(PetscInt);
 18:   header->hsize /= sizeof(DMNetworkComponentGenericDataType);
 19:   return(0);
 20: }

 22: /*@
 23:   DMNetworkGetPlex - Gets the Plex DM associated with this network DM

 25:   Not collective

 27:   Input Parameters:
 28: . dm - the dm object

 30:   Output Parameters:
 31: . plexdm - the plex dm object

 33:   Level: Advanced

 35: .seealso: DMNetworkCreate()
 36: @*/
 37: PetscErrorCode DMNetworkGetPlex(DM dm,DM *plexdm)
 38: {
 39:   DM_Network *network = (DM_Network*)dm->data;

 42:   *plexdm = network->plex;
 43:   return(0);
 44: }

 46: /*@
 47:   DMNetworkGetNumSubNetworks - Gets the the number of subnetworks

 49:   Not collective

 51:   Input Parameter:
 52: . dm - the dm object

 54:   Output Parameters:
 55: + nsubnet - local number of subnetworks
 56: - Nsubnet - global number of subnetworks

 58:   Level: beginner

 60: .seealso: DMNetworkCreate(), DMNetworkSetNumSubNetworks()
 61: @*/
 62: PetscErrorCode DMNetworkGetNumSubNetworks(DM dm,PetscInt *nsubnet,PetscInt *Nsubnet)
 63: {
 64:   DM_Network *network = (DM_Network*)dm->data;

 67:   if (nsubnet) *nsubnet = network->nsubnet;
 68:   if (Nsubnet) *Nsubnet = network->Nsubnet;
 69:   return(0);
 70: }

 72: /*@
 73:   DMNetworkSetNumSubNetworks - Sets the number of subnetworks

 75:   Collective on dm

 77:   Input Parameters:
 78: + dm - the dm object
 79: . nsubnet - local number of subnetworks
 80: - Nsubnet - global number of subnetworks

 82:    Level: beginner

 84: .seealso: DMNetworkCreate(), DMNetworkGetNumSubNetworks()
 85: @*/
 86: PetscErrorCode DMNetworkSetNumSubNetworks(DM dm,PetscInt nsubnet,PetscInt Nsubnet)
 87: {
 89:   DM_Network     *network = (DM_Network*)dm->data;

 92:   if (network->Nsubnet != 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network");


 98:   if (Nsubnet == PETSC_DECIDE) {
 99:     if (nsubnet < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of local subnetworks %D cannot be less than 0",nsubnet);
100:     MPIU_Allreduce(&nsubnet,&Nsubnet,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));
101:   }
102:   if (Nsubnet < 1) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Number of global subnetworks %D cannot be less than 1",Nsubnet);

104:   network->Nsubnet  = Nsubnet;
105:   network->nsubnet  = 0;       /* initia value; will be determind by DMNetworkAddSubnetwork() */
106:   PetscCalloc1(Nsubnet,&network->subnet);

108:   /* num of shared vertices */
109:   network->nsvtx = 0;
110:   network->Nsvtx = 0;
111:   return(0);
112: }

114: /*@
115:   DMNetworkAddSubnetwork - Add a subnetwork

117:   Collective on dm

119:   Input Parameters:
120: + dm - the dm object
121: . name - name of the subnetwork
122: . ne - number of local edges of this subnetwork
123: - edgelist - list of edges for this subnetwork

125:   Output Parameters:
126: . netnum - global index of the subnetwork

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

132:   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, 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.

134:   Level: beginner

136:   Example usage:
137:   Consider the following networks:
138:   1) A sigle subnetwork:
139: .vb
140:  network 0:
141:  rank[0]:
142:    v0 --> v2; v1 --> v2
143:  rank[1]:
144:    v3 --> v5; v4 --> v5
145: .ve

147:  The resulting input
148:  network 0:
149:  rank[0]:
150:    ne = 2
151:    edgelist = [0 2 | 1 2]
152:  rank[1]:
153:    ne = 2
154:    edgelist = [3 5 | 4 5]

156:   2) Two subnetworks:
157: .vb
158:  subnetwork 0:
159:  rank[0]:
160:    v0 --> v2; v2 --> v1; v1 --> v3;
161:  subnetwork 1:
162:  rank[1]:
163:    v0 --> v3; v3 --> v2; v2 --> v1;
164: .ve

166:  The resulting input
167:  subnetwork 0:
168:  rank[0]:
169:    ne = 3
170:    edgelist = [0 2 | 2 1 | 1 3]
171:  rank[1]:
172:    ne = 0
173:    edgelist = NULL

175:  subnetwork 1:
176:  rank[0]:
177:    ne = 0
178:    edgelist = NULL
179:  rank[1]:
180:    edgelist = [0 3 | 3 2 | 2 1]

182: .seealso: DMNetworkCreate(), DMNetworkSetNumSubnetworks()
183: @*/
184: PetscErrorCode DMNetworkAddSubnetwork(DM dm,const char* name,PetscInt ne,PetscInt edgelist[],PetscInt *netnum)
185: {
187:   DM_Network     *network = (DM_Network*)dm->data;
188:   PetscInt       i,Nedge,j,Nvtx,nvtx;
189:   PetscBT        table;

192:   /* Get global total Nvtx = max(edgelist[])+1 for this subnet */
193:   nvtx = -1; i = 0;
194:   for (j=0; j<ne; j++) {
195:     nvtx = PetscMax(nvtx, edgelist[i]); i++;
196:     nvtx = PetscMax(nvtx, edgelist[i]); i++;
197:   }
198:   nvtx++;
199:   MPIU_Allreduce(&nvtx,&Nvtx,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)dm));

201:   /* Get local nvtx for this subnet */
202:   PetscBTCreate(Nvtx,&table);
203:   PetscBTMemzero(Nvtx,table);
204:   i = 0;
205:   for (j=0; j<ne; j++) {
206:     PetscBTSet(table,edgelist[i]);
207:     i++;
208:     PetscBTSet(table,edgelist[i]);
209:     i++;
210:   }
211:   nvtx = 0;
212:   for (j=0; j<Nvtx; j++) {
213:     if (PetscBTLookup(table,j)) nvtx++;
214:   }
215:   PetscBTDestroy(&table);

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

220:   i = network->nsubnet;
221:   if (name) {
222:     PetscStrcpy(network->subnet[i].name,name);
223:   }
224:   network->subnet[i].nvtx     = nvtx; /* include ghost vertices */
225:   network->subnet[i].nedge    = ne;
226:   network->subnet[i].edgelist = edgelist;
227:   network->subnet[i].Nvtx     = Nvtx;
228:   network->subnet[i].Nedge    = Nedge;

230:   /* ----------------------------------------------------------
231:    p=v or e;
232:    subnet[0].pStart   = 0
233:    subnet[i+1].pStart = subnet[i].pEnd = subnet[i].pStart + (nE[i] or NV[i])
234:    ----------------------------------------------------------------------- */
235:   /* GLOBAL subnet[].vStart and vEnd, used by DMNetworkLayoutSetUp() */
236:   network->subnet[i].vStart = network->NVertices;
237:   network->subnet[i].vEnd   = network->subnet[i].vStart + network->subnet[i].Nvtx; /* global vEnd of subnet[i] */

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

242:   /* LOCAL subnet[].eStart and eEnd, used by DMNetworkLayoutSetUp() */
243:   network->subnet[i].eStart = network->nEdges;
244:   network->subnet[i].eEnd   = network->subnet[i].eStart + ne;
245:   network->nEdges += ne;
246:   network->NEdges += network->subnet[i].Nedge;

248:   PetscStrcpy(network->subnet[i].name,name);
249:   if (netnum) *netnum = network->nsubnet;
250:   network->nsubnet++;
251:   return(0);
252: }

254: /*
255:   SetUp a single svtx struct. See SVtx defined in dmnetworkimpl.h
256:   Set gidx and type if input v=(net,idx) is a from_vertex;
257:   Get gid, type and index in the svtx array if input v=(net,idx) is a to_vertex.

259:   Input:  Nsvtx, svtx, net, idx, gidx
260:   Output: gidx, svtype, svtx_idx
261:  */
262: static PetscErrorCode SVtxSetUp(PetscInt Nsvtx,SVtx *svtx,PetscInt net,PetscInt idx,PetscInt *gidx,SVtxType *svtype,PetscInt *svtx_idx)
263: {
264:   PetscInt i,j,*svto;
265:   SVtxType vtype;

268:   if (!Nsvtx) return(0);

270:   vtype = SVNONE;
271:   for (i=0; i<Nsvtx; i++) {
272:     if (net == svtx[i].sv[0] && idx == svtx[i].sv[1]) {
273:       /* (1) input vertex net.idx is a shared from_vertex, set its global index and output its svtype */
274:       svtx[i].gidx = *gidx; /* set gidx */
275:       vtype        = SVFROM;
276:     } else { /* loop over svtx[i].n */
277:       for (j=1; j<svtx[i].n; j++) {
278:         svto = svtx[i].sv + 2*j;
279:         if (net == svto[0] && idx == svto[1]) {
280:           /* input vertex net.idx is a shared to_vertex, output its global index and its svtype */
281:           *gidx = svtx[i].gidx; /* output gidx for to_vertex */
282:           vtype = SVTO;
283:         }
284:       }
285:     }
286:     if (vtype != SVNONE) break;
287:   }
288:   if (svtype)   *svtype   = vtype;
289:   if (svtx_idx) *svtx_idx = i;
290:   return(0);
291: }

293: /*
294:  Add a new shared vertex sv=(net,idx) to table svtas[ita]
295: */
296: static PetscErrorCode TableAddSVtx(PetscTable *svtas,PetscInt ita,PetscInt* tdata,PetscInt *sv_wk,PetscInt *ii,PetscInt *sedgelist,PetscInt k,DM_Network *network,PetscInt **ta2sv)
297: {
298:   PetscInt       net,idx,gidx,i=*ii;

302:   net = sv_wk[2*i]   = sedgelist[k];
303:   idx = sv_wk[2*i+1] = sedgelist[k+1];
304:   gidx = network->subnet[net].vStart + idx;
305:   PetscTableAdd(svtas[ita],gidx+1,tdata[ita]+1,INSERT_VALUES);
306:   *(ta2sv[ita] + tdata[ita]) = i; /* maps tdata to index of sv_wk; sv_wk keeps (net,idx) info */
307:   tdata[ita]++; (*ii)++;
308:   return(0);
309: }

311: /*
312:   Create an array of shared vertices. See SVtx and SVtxType in dmnetworkimpl.h

314:   Input:  dm, Nsedgelist, sedgelist
315:   Output: Nsvtx,svtx

317:   Note: Output svtx is organized as
318:         sv(net[0],idx[0]) --> sv(net[1],idx[1])
319:                           --> sv(net[1],idx[1])
320:                           ...
321:                           --> sv(net[n-1],idx[n-1])
322:         and net[0] < net[1] < ... < net[n-1]
323:         where sv[0] has SVFROM type, sv[i], i>0, has SVTO type.
324:  */
325: static PetscErrorCode SVtxCreate(DM dm,PetscInt Nsedgelist,PetscInt *sedgelist,PetscInt *Nsvtx,SVtx **svtx)
326: {
328:   SVtx           *sedges = NULL;
329:   PetscInt       *sv,k,j,nsv,*tdata,**ta2sv;
330:   PetscTable     *svtas;
331:   PetscInt       gidx,net,idx,i,nta,ita,idx_from,idx_to,n,*sv_wk;
332:   DM_Network     *network = (DM_Network*)dm->data;
333:   PetscTablePosition ppos;

336:   /* (1) Crete ctables svtas */
337:   PetscCalloc4(Nsedgelist,&svtas,Nsedgelist,&tdata,4*Nsedgelist,&sv_wk,2*Nsedgelist,&ta2sv);

339:   j   = 0;   /* sedgelist counter */
340:   k   = 0;   /* sedgelist vertex counter j = 4*k */
341:   i   = 0;   /* sv_wk (vertices added to the ctables) counter */
342:   nta = 0;   /* num of sv tables created */

344:   /* for j=0 */
345:   PetscTableCreate(2*Nsedgelist,network->NVertices+1,&svtas[nta]);
346:   PetscMalloc1(2*Nsedgelist,&ta2sv[nta]);

348:   TableAddSVtx(svtas,nta,tdata,sv_wk,&i,sedgelist,k,network,ta2sv);
349:   TableAddSVtx(svtas,nta,tdata,sv_wk,&i,sedgelist,k+2,network,ta2sv);
350:   nta++; k += 4;

352:   for (j = 1; j < Nsedgelist; j++) {
353:     for (ita = 0; ita < nta; ita++) {
354:       /* vfrom */
355:       net = sedgelist[k]; idx = sedgelist[k+1];
356:       gidx = network->subnet[net].vStart + idx; /* global index of the vertex net.idx before merging shared vertices */
357:       PetscTableFind(svtas[ita],gidx+1,&idx_from);

359:       /* vto */
360:       net = sedgelist[k+2]; idx = sedgelist[k+3];
361:       gidx = network->subnet[net].vStart + idx;
362:       PetscTableFind(svtas[ita],gidx+1,&idx_to);

364:       if (idx_from || idx_to) { /* vfrom or vto is on table svtas[ita] */
365:         idx_from--; idx_to--;
366:         if (idx_from < 0) { /* vto is on svtas[ita] */
367:           TableAddSVtx(svtas,ita,tdata,sv_wk,&i,sedgelist,k,network,ta2sv);
368:           break;
369:         } else if (idx_to < 0) {
370:           TableAddSVtx(svtas,ita,tdata,sv_wk,&i,sedgelist,k+2,network,ta2sv);
371:           break;
372:         }
373:       }
374:     }

376:     if (ita == nta) {
377:       PetscTableCreate(2*Nsedgelist,network->NVertices+1,&svtas[nta]);
378:       PetscMalloc1(2*Nsedgelist, &ta2sv[nta]);

380:       TableAddSVtx(svtas,nta,tdata,sv_wk,&i,sedgelist,k,network,ta2sv);
381:       TableAddSVtx(svtas,nta,tdata,sv_wk,&i,sedgelist,k+2,network,ta2sv);
382:       nta++;
383:     }
384:     k += 4;
385:   }

387:   /* (2) Construct sedges from ctable
388:      sedges: edges connect vertex sv[0]=(net[0],idx[0]) to vertices sv[k], k=1,...,n-1;
389:      net[k], k=0, ...,n-1, are in ascending order */
390:   PetscMalloc1(nta,&sedges);
391:   for (nsv = 0; nsv < nta; nsv++) {
392:     /* for a single svtx, put shared vertices in ascending order of gidx */
393:     PetscTableGetCount(svtas[nsv],&n);
394:     PetscCalloc1(2*n,&sv);
395:     sedges[nsv].sv   = sv;
396:     sedges[nsv].n    = n;
397:     sedges[nsv].gidx = -1; /* initialization */

399:     PetscTableGetHeadPosition(svtas[nsv],&ppos);
400:     for (k=0; k<n; k++) { /* gidx is sorted in ascending order */
401:       PetscTableGetNext(svtas[nsv],&ppos,&gidx,&i);
402:       gidx--; i--;

404:       j = ta2sv[nsv][i]; /* maps i to index of sv_wk */
405:       sv[2*k]   = sv_wk[2*j];
406:       sv[2*k+1] = sv_wk[2*j + 1];
407:     }
408:   }

410:   for (j=0; j<nta; j++) {
411:     PetscTableDestroy(&svtas[j]);
412:     PetscFree(ta2sv[j]);
413:   }
414:   PetscFree4(svtas,tdata,sv_wk,ta2sv);

416:   *Nsvtx = nta;
417:   *svtx  = sedges;
418:   return(0);
419: }

421: static PetscErrorCode DMNetworkLayoutSetUp_Coupling(DM dm)
422: {
424:   DM_Network     *network = (DM_Network*)dm->data;
425:   PetscInt       i,j,ctr,np,*edges,*subnetvtx,*subnetedge,vStart;
426:   PetscInt       k,*vidxlTog,Nsv=0,Nsubnet=network->Nsubnet;
427:   PetscInt       *sedgelist=network->sedgelist;
428:   const PetscInt *cone;
429:   MPI_Comm       comm;
430:   PetscMPIInt    size,rank,*recvcounts=NULL,*displs=NULL;
431:   PetscInt       net,idx,gidx,nmerged,e,v,vfrom,vto,*vrange,*eowners,gidx_from,net_from,sv_idx;
432:   SVtxType       svtype = SVNONE;
433:   SVtx           *svtx=NULL;
434:   PetscSection   sectiong;

437:   /* This implementation requires user input each subnet by a single processor, thus subnet[net].nvtx=subnet[net].Nvtx */
438:   for (net=0; net<Nsubnet; net++) {
439:     if (network->subnet[net].nvtx && network->subnet[net].nvtx != network->subnet[net].Nvtx) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"subnetwork %D local num of vertices %D != %D global num",net,network->subnet[net].nvtx,network->subnet[net].Nvtx);
440:   }

442:   PetscObjectGetComm((PetscObject)dm,&comm);
443:   MPI_Comm_rank(comm,&rank);
444:   MPI_Comm_size(comm,&size);

446:   /* (1) Create svtx[] from sedgelist */
447:   /* -------------------------------- */
448:   /* Nsv: global number of SVtx; svtx: shared vertices, see SVtx in dmnetworkimpl.h */
449:   SVtxCreate(dm,network->Nsvtx,sedgelist,&Nsv,&svtx);

451:   /* (2) Setup svtx; Shared vto vertices are merged to their vfrom vertex with same global vetex index (gidx) */
452:   /* -------------------------------------------------------------------------------------------------------- */
453:   /* (2.1) compute vrage[rank]: global index of 1st local vertex in proc[rank] */
454:   PetscMalloc3(size+1,&vrange,size,&displs,size,&recvcounts);
455:   for (i=0; i<size; i++) {displs[i] = i; recvcounts[i] = 1;}

457:   vrange[0] = 0;
458:   MPI_Allgatherv(&network->nVertices,1,MPIU_INT,vrange+1,recvcounts,displs,MPIU_INT,comm);
459:   for (i=2; i<size+1; i++) {
460:     vrange[i] += vrange[i-1];
461:   }

463:   /* (2.2) Create vidxlTog: maps UN-MERGED local vertex index i to global index gidx (plex, excluding ghost vertices) */
464:   PetscMalloc1(network->nVertices,&vidxlTog);
465:   i = 0; gidx = 0;
466:   nmerged = 0; /* local num of merged vertices */
467:   network->nsvtx = 0;
468:   for (net=0; net<Nsubnet; net++) {
469:     for (idx=0; idx<network->subnet[net].Nvtx; idx++) {
470:       gidx_from = gidx;
471:       sv_idx    = -1;

473:       SVtxSetUp(Nsv,svtx,net,idx,&gidx_from,&svtype,&sv_idx);
474:       if (svtype == SVTO) {
475:         if (network->subnet[net].nvtx) {/* this proc owns sv_to */
476:           net_from = svtx[sv_idx].sv[0]; /* subnet num of its shared vertex */
477:           if (network->subnet[net_from].nvtx == 0) {
478:             /* this proc does not own v_from, thus a new local coupling vertex */
479:             network->nsvtx++;
480:           }
481:           vidxlTog[i++] = gidx_from;
482:           nmerged++; /* a coupling vertex -- merged */
483:         }
484:       } else {
485:         if (svtype == SVFROM) {
486:           if (network->subnet[net].nvtx) {
487:             /* this proc owns this v_from, a new local coupling vertex */
488:             network->nsvtx++;
489:           }
490:         }
491:         if (network->subnet[net].nvtx) vidxlTog[i++] = gidx;
492:         gidx++;
493:       }
494:     }
495:   }
496: #if defined(PETSC_USE_DEBUG)
497:   if (i != network->nVertices) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"%D != %D nVertices",i,network->nVertices);
498: #endif

500:   /* (2.3) Setup svtable for querry shared vertices */
501:   for (v=0; v<Nsv; v++) {
502:     gidx = svtx[v].gidx;
503:     PetscTableAdd(network->svtable,gidx+1,v+1,INSERT_VALUES);
504:   }

506:   /* (2.4) Shared vertices in the subnetworks are merged, update global NVertices: np = sum(local nmerged) */
507:   MPI_Allreduce(&nmerged,&np,1,MPIU_INT,MPI_SUM,comm);
508:   network->NVertices -= np;

510:   PetscCalloc1(2*network->nEdges,&edges);

512:   ctr = 0;
513:   for (net=0; net<Nsubnet; net++) {
514:     for (j = 0; j < network->subnet[net].nedge; j++) {
515:       /* vfrom: */
516:       i = network->subnet[net].edgelist[2*j] + (network->subnet[net].vStart - vrange[rank]);
517:       edges[2*ctr] = vidxlTog[i];

519:       /* vto */
520:       i = network->subnet[net].edgelist[2*j+1] + (network->subnet[net].vStart - vrange[rank]);
521:       edges[2*ctr+1] = vidxlTog[i];
522:       ctr++;
523:     }
524:   }
525:   PetscFree3(vrange,displs,recvcounts);
526:   PetscFree(vidxlTog);

528:   /* (3) Create network->plex */
529:   DMCreate(comm,&network->plex);
530:   DMSetType(network->plex,DMPLEX);
531:   DMSetDimension(network->plex,1);
532:   if (size == 1) {
533:     DMPlexBuildFromCellList(network->plex,network->nEdges,network->nVertices-nmerged,2,edges);
534:   } else {
535:     DMPlexBuildFromCellListParallel(network->plex,network->nEdges,network->nVertices-nmerged,PETSC_DECIDE,2,edges,NULL);
536:   }

538:   DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);
539:   DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);
540:   DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);
541:   vStart = network->vStart;

543:   PetscSectionCreate(comm,&network->DataSection);
544:   PetscSectionCreate(comm,&network->DofSection);
545:   PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);
546:   PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);

548:   np = network->pEnd - network->pStart;
549:   PetscCalloc2(np,&network->header,np,&network->cvalue);
550:   for (i=0; i<np; i++) {
551:     network->header[i].maxcomps = 1;
552:     SetUpNetworkHeaderComponentValue(dm,&network->header[i],&network->cvalue[i]);
553:   }

555:   /* (4) Create edge and vertex arrays for the subnetworks */
556:   PetscCalloc2(network->nEdges,&subnetedge,network->nVertices+network->nsvtx,&subnetvtx); /* Maps local edge/vertex to local subnetwork's edge/vertex */
557:   network->subnetedge = subnetedge;
558:   network->subnetvtx  = subnetvtx;
559:   for (j=0; j < Nsubnet; j++) {
560:     network->subnet[j].edges = subnetedge;
561:     subnetedge              += network->subnet[j].nedge;

563:     network->subnet[j].vertices = subnetvtx;
564:     subnetvtx                  += network->subnet[j].nvtx;
565:   }
566:   network->svertices = subnetvtx;

568:   /* Get edge ownership */
569:   PetscMalloc1(size+1,&eowners);
570:   np = network->eEnd - network->eStart; /* num of local edges */
571:   MPI_Allgather(&np,1,MPIU_INT,eowners+1,1,MPIU_INT,comm);
572:   eowners[0] = 0;
573:   for (i=2; i<=size; i++) eowners[i] += eowners[i-1];

575:   /* Setup edge and vertex arrays for subnetworks */
576:   e = 0;
577:   for (i=0; i < Nsubnet; i++) {
578:     ctr = 0;
579:     for (j = 0; j < network->subnet[i].nedge; j++) {
580:       /* edge e */
581:       network->header[e].index    = e + eowners[rank]; /* Global edge index */
582:       network->header[e].subnetid = i;                 /* Subnetwork id */
583:       network->subnet[i].edges[j] = e;
584:       network->header[e].ndata           = 0;
585:       network->header[e].offset[0]       = 0;
586:       network->header[e].offsetvarrel[0] = 0;
587:       PetscSectionAddDof(network->DataSection,e,network->header[e].hsize);

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

592:       /* vertex cone[0] */
593:       v = cone[0];
594:       network->header[v].index    = edges[2*e];   /* Global vertex index */
595:       network->header[v].subnetid = i;            /* Subnetwork id */
596:       vfrom = network->subnet[i].edgelist[2*ctr];
597:       network->subnet[i].vertices[vfrom] = v;     /* user's subnet[].dix = petsc's v */

599:       /* vertex cone[1] */
600:       v = cone[1];
601:       network->header[v].index    = edges[2*e+1]; /* Global vertex index */
602:       network->header[v].subnetid = i;
603:       vto = network->subnet[i].edgelist[2*ctr+1];
604:       network->subnet[i].vertices[vto]= v;

606:       e++; ctr++;
607:     }
608:   }
609:   PetscFree(eowners);
610:   PetscFree(edges);

612:   /* Set vertex array for the subnetworks */
613:   k = 0;
614:   for (v=vStart; v<network->vEnd; v++) { /* local vertices, including ghosts */
615:     network->header[v].ndata           = 0;
616:     network->header[v].offset[0]       = 0;
617:     network->header[v].offsetvarrel[0] = 0;
618:     PetscSectionAddDof(network->DataSection,v,network->header[v].hsize);

620:     /* shared vertex */
621:     PetscTableFind(network->svtable,network->header[v].index+1,&i);
622:     if (i) network->svertices[k++] = v;
623:   }
624: #if defined(PETSC_USE_DEBUG)
625:   if (k != network->nsvtx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"k %D != %D nsvtx",k,network->nsvtx);
626: #endif

628:   network->svtx  = svtx;
629:   network->Nsvtx = Nsv;
630:   PetscFree(sedgelist);

632:   /* Create a global section to be used by DMNetworkIsGhostVertex() which is a non-collective routine */
633:   /* see snes_tutorials_network-ex1_4 */
634:   DMGetGlobalSection(network->plex,&sectiong);
635:   return(0);
636: }

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

641:   Collective on dm

643:   Input Parameters:
644: . dm - the dmnetwork object

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

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

652:   Level: beginner

654: .seealso: DMNetworkSetNumSubNetworks(), DMNetworkAddSubnetwork()
655: @*/
656: PetscErrorCode DMNetworkLayoutSetUp(DM dm)
657: {
659:   DM_Network     *network = (DM_Network*)dm->data;
660:   PetscInt       i,j,ctr,Nsubnet=network->Nsubnet,*eowners,np,*edges,*subnetvtx,*subnetedge,e,v,vfrom,vto;
661:   const PetscInt *cone;
662:   MPI_Comm       comm;
663:   PetscMPIInt    size,rank;

666:   if (network->nsubnet != network->Nsubnet) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Must call DMNetworkAddSubnetwork() %D times",network->Nsubnet);

668:   /* Create svtable for querry shared vertices */
669:   PetscTableCreate(network->Nsvtx,network->NVertices+1,&network->svtable);

671:   if (network->Nsvtx) {
672:     DMNetworkLayoutSetUp_Coupling(dm);
673:     return(0);
674:   }

676:   PetscObjectGetComm((PetscObject)dm,&comm);
677:   MPI_Comm_rank(comm,&rank);
678:   MPI_Comm_size(comm,&size);

680:   /* Create LOCAL edgelist in global vertex ordering for the network by concatenating local input edgelists of the subnetworks */
681:   PetscCalloc1(2*network->nEdges,&edges);
682:   ctr = 0;
683:   for (i=0; i < Nsubnet; i++) {
684:     for (j = 0; j < network->subnet[i].nedge; j++) {
685:       edges[2*ctr]   = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
686:       edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
687:       ctr++;
688:     }
689:   }

691:   /* Create network->plex; One dimensional network, numCorners=2 */
692:   DMCreate(comm,&network->plex);
693:   DMSetType(network->plex,DMPLEX);
694:   DMSetDimension(network->plex,1);
695:   if (size == 1) {
696:     DMPlexBuildFromCellList(network->plex,network->nEdges,network->nVertices,2,edges);
697:   } else {
698:     DMPlexBuildFromCellListParallel(network->plex,network->nEdges,network->nVertices,PETSC_DECIDE,2,edges,NULL);
699:   }

701:   DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);
702:   DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);
703:   DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);

705:   PetscSectionCreate(comm,&network->DataSection);
706:   PetscSectionCreate(comm,&network->DofSection);
707:   PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);
708:   PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);

710:   np = network->pEnd - network->pStart;
711:   PetscCalloc2(np,&network->header,np,&network->cvalue);
712:   for (i=0; i < np; i++) {
713:     network->header[i].maxcomps = 1;
714:     SetUpNetworkHeaderComponentValue(dm,&network->header[i],&network->cvalue[i]);
715:   }

717:   /* Create edge and vertex arrays for the subnetworks
718:      This implementation assumes that DMNetwork reads
719:      (1) a single subnetwork in parallel; or
720:      (2) n subnetworks using n processors, one subnetwork/processor.
721:    */
722:   PetscCalloc2(network->nEdges,&subnetedge,network->nVertices,&subnetvtx); /* Maps local edge/vertex to local subnetwork's edge/vertex */

724:   network->subnetedge = subnetedge;
725:   network->subnetvtx  = subnetvtx;
726:   for (j=0; j < network->Nsubnet; j++) {
727:     network->subnet[j].edges = subnetedge;
728:     subnetedge              += network->subnet[j].nedge;

730:     network->subnet[j].vertices = subnetvtx;
731:     subnetvtx                  += network->subnet[j].nvtx;
732:   }

734:   /* Get edge ownership */
735:   PetscMalloc1(size+1,&eowners);
736:   np = network->eEnd - network->eStart;
737:   MPI_Allgather(&np,1,MPIU_INT,eowners+1,1,MPIU_INT,comm);
738:   eowners[0] = 0;
739:   for (i=2; i<=size; i++) eowners[i] += eowners[i-1];

741:   /* Setup edge and vertex arrays for subnetworks */
742:   e = 0;
743:   for (i=0; i < Nsubnet; i++) {
744:     ctr = 0;
745:     for (j = 0; j < network->subnet[i].nedge; j++) {
746:       /* edge e */
747:       network->header[e].index    = e + eowners[rank];   /* Global edge index */
748:       network->header[e].subnetid = i;
749:       network->subnet[i].edges[j] = e;

751:       network->header[e].ndata           = 0;
752:       network->header[e].offset[0]       = 0;
753:       network->header[e].offsetvarrel[0] = 0;
754:       PetscSectionAddDof(network->DataSection,e,network->header[e].hsize);

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

759:       /* vertex cone[0] */
760:       v = cone[0];
761:       network->header[v].index     = edges[2*e];  /* Global vertex index */
762:       network->header[v].subnetid  = i;           /* Subnetwork id */
763:       if (Nsubnet == 1) {
764:         network->subnet[i].vertices[v - network->vStart] = v; /* user's subnet[].idx = petsc's v */
765:       } else {
766:         vfrom = network->subnet[i].edgelist[2*ctr];     /* =subnet[i].idx, Global index! */
767:         network->subnet[i].vertices[vfrom] = v; /* user's subnet[].dix = petsc's v */
768:       }

770:       /* vertex cone[1] */
771:       v = cone[1];
772:       network->header[v].index    = edges[2*e+1];   /* Global vertex index */
773:       network->header[v].subnetid = i;              /* Subnetwork id */
774:       if (Nsubnet == 1) {
775:         network->subnet[i].vertices[v - network->vStart] = v; /* user's subnet[].idx = petsc's v */
776:       } else {
777:         vto = network->subnet[i].edgelist[2*ctr+1];     /* =subnet[i].idx, Global index! */
778:         network->subnet[i].vertices[vto] = v; /* user's subnet[].dix = petsc's v */
779:       }

781:       e++; ctr++;
782:     }
783:   }
784:   PetscFree(eowners);
785:   PetscFree(edges); /* local edge list with global idx used by DMPlexBuildFromCellList() */

787:   for (v = network->vStart; v < network->vEnd; v++) {
788:     network->header[v].ndata           = 0;
789:     network->header[v].offset[0]       = 0;
790:     network->header[v].offsetvarrel[0] = 0;
791:     PetscSectionAddDof(network->DataSection,v,network->header[e].hsize);
792:   }
793:   return(0);
794: }

796: /*@C
797:   DMNetworkGetSubnetwork - Returns the information about a requested subnetwork

799:   Not collective

801:   Input Parameters:
802: + dm - the DM object
803: - netnum - the global index of the subnetwork

805:   Output Parameters:
806: + nv - number of vertices (local)
807: . ne - number of edges (local)
808: . vtx - local vertices of the subnetwork
809: - edge - local edges of the subnetwork

811:   Notes:
812:   Cannot call this routine before DMNetworkLayoutSetup()

814:   Level: intermediate

816: .seealso: DMNetworkCreate(), DMNetworkAddSubnetwork(), DMNetworkLayoutSetUp()
817: @*/
818: PetscErrorCode DMNetworkGetSubnetwork(DM dm,PetscInt netnum,PetscInt *nv,PetscInt *ne,const PetscInt **vtx,const PetscInt **edge)
819: {
820:   DM_Network *network = (DM_Network*)dm->data;

823:   if (netnum >= network->Nsubnet) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Subnet index %D exceeds the num of subnets %D",netnum,network->Nsubnet);
824:   if (nv) *nv     = network->subnet[netnum].nvtx;
825:   if (ne) *ne     = network->subnet[netnum].nedge;
826:   if (vtx) *vtx   = network->subnet[netnum].vertices;
827:   if (edge) *edge = network->subnet[netnum].edges;
828:   return(0);
829: }

831: /*@
832:   DMNetworkAddSharedVertices - Add shared vertices that connect two given subnetworks

834:   Collective on dm

836:   Input Parameters:
837: + dm - the dm object
838: . anetnum - first subnetwork global numbering returned by DMNetworkAddSubnetwork()
839: . bnetnum - second subnetwork global numbering returned by DMNetworkAddSubnetwork()
840: . nsvtx - number of vertices that are shared by the two subnetworks
841: . asvtx - vertex index in the first subnetwork
842: - bsvtx - vertex index in the second subnetwork

844:   Level: beginner

846: .seealso: DMNetworkCreate(), DMNetworkAddSubnetwork(), DMNetworkGetSharedVertices()
847: @*/
848: PetscErrorCode DMNetworkAddSharedVertices(DM dm,PetscInt anetnum,PetscInt bnetnum,PetscInt nsvtx,PetscInt asvtx[],PetscInt bsvtx[])
849: {
851:   DM_Network     *network = (DM_Network*)dm->data;
852:   PetscInt       i,nsubnet = network->Nsubnet,*sedgelist,Nsvtx=network->Nsvtx;

855:   if (anetnum == bnetnum) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Subnetworks must have different netnum");
856:   if (anetnum < 0 || bnetnum < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"netnum cannot be negative");
857:   if (!Nsvtx) {
858:     /* allocate network->sedgelist to hold at most 2*nsubnet pairs of shared vertices */
859:     PetscMalloc1(2*4*nsubnet,&network->sedgelist);
860:   }

862:   sedgelist = network->sedgelist;
863:   for (i=0; i<nsvtx; i++) {
864:     sedgelist[4*Nsvtx]   = anetnum; sedgelist[4*Nsvtx+1] = asvtx[i];
865:     sedgelist[4*Nsvtx+2] = bnetnum; sedgelist[4*Nsvtx+3] = bsvtx[i];
866:     Nsvtx++;
867:   }
868:   if (Nsvtx > 2*nsubnet) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"allocate more space for coupling edgelist");
869:   network->Nsvtx = Nsvtx;
870:   return(0);
871: }

873: /*@C
874:   DMNetworkGetSharedVertices - Returns the info for the shared vertices

876:   Not collective

878:   Input Parameter:
879: . dm - the DM object

881:   Output Parameters:
882: + nsv - number of local shared vertices
883: - svtx - local shared vertices

885:   Notes:
886:   Cannot call this routine before DMNetworkLayoutSetup()

888:   Level: intermediate

890: .seealso: DMNetworkGetSubnetwork(), DMNetworkLayoutSetUp(), DMNetworkAddSharedVertices()
891: @*/
892: PetscErrorCode DMNetworkGetSharedVertices(DM dm,PetscInt *nsv,const PetscInt **svtx)
893: {
894:   DM_Network *net = (DM_Network*)dm->data;

897:   if (net->Nsvtx) {
898:     *nsv  = net->nsvtx;
899:     *svtx = net->svertices;
900:   } else {
901:     *nsv  = 0;
902:     *svtx = NULL;
903:   }
904:   return(0);
905: }

907: /*@C
908:   DMNetworkRegisterComponent - Registers the network component

910:   Logically collective on dm

912:   Input Parameters:
913: + dm - the network object
914: . name - the component name
915: - size - the storage size in bytes for this component data

917:    Output Parameters:
918: .  key - an integer key that defines the component

920:    Notes
921:    This routine should be called by all processors before calling DMNetworkLayoutSetup().

923:    Level: beginner

925: .seealso: DMNetworkCreate(), DMNetworkLayoutSetUp()
926: @*/
927: PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,size_t size,PetscInt *key)
928: {
929:   PetscErrorCode        ierr;
930:   DM_Network            *network = (DM_Network*) dm->data;
931:   DMNetworkComponent    *component=NULL,*newcomponent=NULL;
932:   PetscBool             flg=PETSC_FALSE;
933:   PetscInt              i;

936:   if (!network->component) {
937:     PetscCalloc1(network->max_comps_registered,&network->component);
938:   }

940:   for (i=0; i < network->ncomponent; i++) {
941:     PetscStrcmp(network->component[i].name,name,&flg);
942:     if (flg) {
943:       *key = i;
944:       return(0);
945:     }
946:   }

948:   if (network->ncomponent == network->max_comps_registered) {
949:     /* Reached max allowed so resize component */
950:     network->max_comps_registered += 2;
951:     PetscCalloc1(network->max_comps_registered,&newcomponent);
952:     /* Copy over the previous component info */
953:     for (i=0; i < network->ncomponent; i++) {
954:       PetscStrcpy(newcomponent[i].name,network->component[i].name);
955:       newcomponent[i].size = network->component[i].size;
956:     }
957:     /* Free old one */
958:     PetscFree(network->component);
959:     /* Update pointer */
960:     network->component = newcomponent;
961:   }

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

965:   PetscStrcpy(component->name,name);
966:   component->size = size/sizeof(DMNetworkComponentGenericDataType);
967:   *key = network->ncomponent;
968:   network->ncomponent++;
969:   return(0);
970: }

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

975:   Not Collective

977:   Input Parameter:
978: . dm - the DMNetwork object

980:   Output Parameters:
981: + vStart - the first vertex point
982: - vEnd - one beyond the last vertex point

984:   Level: beginner

986: .seealso: DMNetworkGetEdgeRange()
987: @*/
988: PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
989: {
990:   DM_Network *network = (DM_Network*)dm->data;

993:   if (vStart) *vStart = network->vStart;
994:   if (vEnd) *vEnd = network->vEnd;
995:   return(0);
996: }

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

1001:   Not Collective

1003:   Input Parameter:
1004: . dm - the DMNetwork object

1006:   Output Parameters:
1007: + eStart - The first edge point
1008: - eEnd - One beyond the last edge point

1010:   Level: beginner

1012: .seealso: DMNetworkGetVertexRange()
1013: @*/
1014: PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
1015: {
1016:   DM_Network *network = (DM_Network*)dm->data;

1019:   if (eStart) *eStart = network->eStart;
1020:   if (eEnd) *eEnd = network->eEnd;
1021:   return(0);
1022: }

1024: static PetscErrorCode DMNetworkGetIndex(DM dm,PetscInt p,PetscInt *index)
1025: {
1026:   PetscErrorCode    ierr;
1027:   DM_Network        *network = (DM_Network*)dm->data;
1028:   PetscInt          offsetp;
1029:   DMNetworkComponentHeader header;

1032:   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first");
1033:   PetscSectionGetOffset(network->DataSection,p,&offsetp);
1034:   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
1035:   *index = header->index;
1036:   return(0);
1037: }

1039: /*@
1040:   DMNetworkGetGlobalEdgeIndex - Get the global numbering for the edge on the network

1042:   Not Collective

1044:   Input Parameters:
1045: + dm - DMNetwork object
1046: - p - edge point

1048:   Output Parameters:
1049: . index - the global numbering for the edge

1051:   Level: intermediate

1053: .seealso: DMNetworkGetGlobalVertexIndex()
1054: @*/
1055: PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
1056: {

1060:   DMNetworkGetIndex(dm,p,index);
1061:   return(0);
1062: }

1064: /*@
1065:   DMNetworkGetGlobalVertexIndex - Get the global numbering for the vertex on the network

1067:   Not Collective

1069:   Input Parameters:
1070: + dm - DMNetwork object
1071: - p  - vertex point

1073:   Output Parameters:
1074: . index - the global numbering for the vertex

1076:   Level: intermediate

1078: .seealso: DMNetworkGetGlobalEdgeIndex()
1079: @*/
1080: PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
1081: {

1085:   DMNetworkGetIndex(dm,p,index);
1086:   return(0);
1087: }

1089: /*@
1090:   DMNetworkGetNumComponents - Get the number of components at a vertex/edge

1092:   Not Collective

1094:   Input Parameters:
1095: + dm - the DMNetwork object
1096: - p - vertex/edge point

1098:   Output Parameters:
1099: . numcomponents - Number of components at the vertex/edge

1101:   Level: beginner

1103: .seealso: DMNetworkRegisterComponent(), DMNetworkAddComponent()
1104: @*/
1105: PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
1106: {
1108:   PetscInt       offset;
1109:   DM_Network     *network = (DM_Network*)dm->data;

1112:   PetscSectionGetOffset(network->DataSection,p,&offset);
1113:   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
1114:   return(0);
1115: }

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

1120:   Not Collective

1122:   Input Parameters:
1123: + dm - the DMNetwork object
1124: . p - the edge/vertex point
1125: - compnum - component number; use ALL_COMPONENTS if no specific component is requested

1127:   Output Parameters:
1128: . offset - the local offset

1130:   Level: intermediate

1132: .seealso: DMGetLocalVector(), DMNetworkGetComponent(), DMNetworkGetGlobalVecOffset()
1133: @*/
1134: PetscErrorCode DMNetworkGetLocalVecOffset(DM dm,PetscInt p,PetscInt compnum,PetscInt *offset)
1135: {
1137:   DM_Network     *network = (DM_Network*)dm->data;
1138:   PetscInt       offsetp,offsetd;
1139:   DMNetworkComponentHeader header;

1142:   PetscSectionGetOffset(network->plex->localSection,p,&offsetp);
1143:   if (compnum == ALL_COMPONENTS) {
1144:     *offset = offsetp;
1145:     return(0);
1146:   }

1148:   PetscSectionGetOffset(network->DataSection,p,&offsetd);
1149:   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetd);
1150:   *offset = offsetp + header->offsetvarrel[compnum];
1151:   return(0);
1152: }

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

1157:   Not Collective

1159:   Input Parameters:
1160: + dm - the DMNetwork object
1161: . p - the edge/vertex point
1162: - compnum - component number; use ALL_COMPONENTS if no specific component is requested

1164:   Output Parameters:
1165: . offsetg - the global offset

1167:   Level: intermediate

1169: .seealso: DMNetworkGetLocalVecOffset(), DMGetGlobalVector(), DMNetworkGetComponent()
1170: @*/
1171: PetscErrorCode DMNetworkGetGlobalVecOffset(DM dm,PetscInt p,PetscInt compnum,PetscInt *offsetg)
1172: {
1174:   DM_Network     *network = (DM_Network*)dm->data;
1175:   PetscInt       offsetp,offsetd;
1176:   DMNetworkComponentHeader header;

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

1182:   if (compnum == ALL_COMPONENTS) {
1183:     *offsetg = offsetp;
1184:     return(0);
1185:   }
1186:   PetscSectionGetOffset(network->DataSection,p,&offsetd);
1187:   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetd);
1188:   *offsetg = offsetp + header->offsetvarrel[compnum];
1189:   return(0);
1190: }

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

1195:   Not Collective

1197:   Input Parameters:
1198: + dm - the DMNetwork object
1199: - p - the edge point

1201:   Output Parameters:
1202: . offset - the offset

1204:   Level: intermediate

1206: .seealso: DMNetworkGetLocalVecOffset(), DMGetLocalVector()
1207: @*/
1208: PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
1209: {
1211:   DM_Network     *network = (DM_Network*)dm->data;

1214:   PetscSectionGetOffset(network->edge.DofSection,p,offset);
1215:   return(0);
1216: }

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

1221:   Not Collective

1223:   Input Parameters:
1224: + dm - the DMNetwork object
1225: - p - the vertex point

1227:   Output Parameters:
1228: . offset - the offset

1230:   Level: intermediate

1232: .seealso: DMNetworkGetEdgeOffset(), DMGetLocalVector()
1233: @*/
1234: PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
1235: {
1237:   DM_Network     *network = (DM_Network*)dm->data;

1240:   p -= network->vStart;
1241:   PetscSectionGetOffset(network->vertex.DofSection,p,offset);
1242:   return(0);
1243: }

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

1248:   Not Collective

1250:   Input Parameters:
1251: + dm - the DMNetwork
1252: . p - the vertex/edge point
1253: . componentkey - component key returned while registering the component; ignored if compvalue=NULL
1254: . compvalue - pointer to the data structure for the component, or NULL if not required.
1255: - nvar - number of variables for the component at the vertex/edge point

1257:   Level: beginner

1259: .seealso: DMNetworkGetComponent()
1260: @*/
1261: PetscErrorCode DMNetworkAddComponent(DM dm,PetscInt p,PetscInt componentkey,void* compvalue,PetscInt nvar)
1262: {
1263:   PetscErrorCode           ierr;
1264:   DM_Network               *network = (DM_Network*)dm->data;
1265:   DMNetworkComponent       *component = &network->component[componentkey];
1266:   DMNetworkComponentHeader header;
1267:   DMNetworkComponentValue  cvalue;
1268:   PetscBool                sharedv=PETSC_FALSE;
1269:   PetscInt                 compnum;
1270:   PetscInt                 *compsize,*compkey,*compoffset,*compnvar,*compoffsetvarrel;
1271:   void*                    *compdata;

1274:   PetscSectionAddDof(network->DofSection,p,nvar);
1275:   if (!compvalue) return(0);

1277:   DMNetworkIsSharedVertex(dm,p,&sharedv);
1278:   if (sharedv) {
1279:     PetscBool ghost;
1280:     DMNetworkIsGhostVertex(dm,p,&ghost);
1281:     if (ghost) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Adding a component at a leaf(ghost) shared vertex is not supported");
1282:   }

1284:   header = &network->header[p];
1285:   cvalue = &network->cvalue[p];
1286:   if (header->ndata == header->maxcomps) {
1287:     PetscInt additional_size;

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

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

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

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

1301:     /* Copy over component info */
1302:     PetscMemcpy(compsize,header->size,header->ndata*sizeof(PetscInt));
1303:     PetscMemcpy(compkey,header->key,header->ndata*sizeof(PetscInt));
1304:     PetscMemcpy(compoffset,header->offset,header->ndata*sizeof(PetscInt));
1305:     PetscMemcpy(compnvar,header->nvar,header->ndata*sizeof(PetscInt));
1306:     PetscMemcpy(compoffsetvarrel,header->offsetvarrel,header->ndata*sizeof(PetscInt));

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

1311:     /* Free old arrays */
1312:     PetscFree5(header->size,header->key,header->offset,header->nvar,header->offsetvarrel);
1313:     PetscFree(cvalue->data);

1315:     /* Update pointers */
1316:     header->size = compsize;
1317:     header->key  = compkey;
1318:     header->offset = compoffset;
1319:     header->nvar = compnvar;
1320:     header->offsetvarrel = compoffsetvarrel;

1322:     cvalue->data = compdata;

1324:     /* Update DataSection Dofs */
1325:     /* 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 */
1326:     additional_size = (5*(header->maxcomps - header->ndata)*sizeof(PetscInt))/sizeof(DMNetworkComponentGenericDataType);
1327:     PetscSectionAddDof(network->DataSection,p,additional_size);
1328:   }
1329:   header = &network->header[p];
1330:   cvalue = &network->cvalue[p];

1332:   compnum = header->ndata;

1334:   header->size[compnum] = component->size;
1335:   PetscSectionAddDof(network->DataSection,p,component->size);
1336:   header->key[compnum] = componentkey;
1337:   if (compnum != 0) header->offset[compnum] = header->offset[compnum-1] + header->size[compnum-1];
1338:   cvalue->data[compnum] = (void*)compvalue;

1340:   /* variables */
1341:   header->nvar[compnum] += nvar;
1342:   if (compnum != 0) header->offsetvarrel[compnum] = header->offsetvarrel[compnum-1] + header->nvar[compnum-1];

1344:   header->ndata++;
1345:   return(0);
1346: }

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

1351:   Not Collective

1353:   Input Parameters:
1354: + dm - the DMNetwork object
1355: . p - vertex/edge point
1356: - compnum - component number; use ALL_COMPONENTS if sum up all the components

1358:   Output Parameters:
1359: + compkey - the key obtained when registering the component (use NULL if not required)
1360: . component - the component data (use NULL if not required)
1361: - nvar  - number of variables (use NULL if not required)

1363:   Level: beginner

1365: .seealso: DMNetworkAddComponent(), DMNetworkGetNumComponents()
1366: @*/
1367: PetscErrorCode DMNetworkGetComponent(DM dm,PetscInt p,PetscInt compnum,PetscInt *compkey,void **component,PetscInt *nvar)
1368: {
1370:   DM_Network     *network = (DM_Network*)dm->data;
1371:   PetscInt       offset = 0;
1372:   DMNetworkComponentHeader header;

1375:   if (compnum == ALL_COMPONENTS) {
1376:     PetscSectionGetDof(network->DofSection,p,nvar);
1377:     return(0);
1378:   }

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

1383:   if (compnum >= 0) {
1384:     if (compkey) *compkey = header->key[compnum];
1385:     if (component) {
1386:       offset += header->hsize+header->offset[compnum];
1387:       *component = network->componentdataarray+offset;
1388:     }
1389:   }

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

1393:   return(0);
1394: }

1396: /*
1397:  Sets up the array that holds the data for all components and its associated section.
1398:  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.
1399: */
1400: PetscErrorCode DMNetworkComponentSetUp(DM dm)
1401: {
1402:   PetscErrorCode           ierr;
1403:   DM_Network               *network = (DM_Network*)dm->data;
1404:   PetscInt                 arr_size,p,offset,offsetp,ncomp,i,*headerarr;
1405:   DMNetworkComponentHeader header;
1406:   DMNetworkComponentValue  cvalue;
1407:   DMNetworkComponentHeader headerinfo;
1408:   DMNetworkComponentGenericDataType *componentdataarray;

1411:   PetscSectionSetUp(network->DataSection);
1412:   PetscSectionGetStorageSize(network->DataSection,&arr_size);
1413:   /* arr_size+1 fixes pipeline test of opensolaris-misc for src/dm/tests/ex10.c -- Do not know why */
1414:   PetscCalloc1(arr_size+1,&network->componentdataarray);
1415:   componentdataarray = network->componentdataarray;
1416:   for (p = network->pStart; p < network->pEnd; p++) {
1417:     PetscSectionGetOffset(network->DataSection,p,&offsetp);
1418:     /* Copy header */
1419:     header = &network->header[p];
1420:     headerinfo = (DMNetworkComponentHeader)(componentdataarray+offsetp);
1421:     PetscMemcpy(headerinfo,header,sizeof(struct _p_DMNetworkComponentHeader));
1422:     headerarr = (PetscInt*)(headerinfo+1);
1423:     PetscMemcpy(headerarr,header->size,header->maxcomps*sizeof(PetscInt));
1424:     headerarr += header->maxcomps;
1425:     PetscMemcpy(headerarr,header->key,header->maxcomps*sizeof(PetscInt));
1426:     headerarr += header->maxcomps;
1427:     PetscMemcpy(headerarr,header->offset,header->maxcomps*sizeof(PetscInt));
1428:     headerarr += header->maxcomps;
1429:     PetscMemcpy(headerarr,header->nvar,header->maxcomps*sizeof(PetscInt));
1430:     headerarr += header->maxcomps;
1431:     PetscMemcpy(headerarr,header->offsetvarrel,header->maxcomps*sizeof(PetscInt));

1433:     /* Copy data */
1434:     cvalue = &network->cvalue[p];
1435:     ncomp  = header->ndata;

1437:     for (i = 0; i < ncomp; i++) {
1438:       offset = offsetp + header->hsize + header->offset[i];
1439:       PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));
1440:     }
1441:   }
1442:   return(0);
1443: }

1445: /* Sets up the section for dofs. This routine is called during DMSetUp() */
1446: static PetscErrorCode DMNetworkVariablesSetUp(DM dm)
1447: {
1449:   DM_Network     *network = (DM_Network*)dm->data;
1450:   MPI_Comm       comm;
1451:   PetscMPIInt    size;

1454:   PetscObjectGetComm((PetscObject)dm,&comm);
1455:   MPI_Comm_size(comm,&size);

1457:   if (size > 1) { /* Sync nvar at shared vertices for all processes */
1458:     PetscSF           sf = network->plex->sf;
1459:     PetscInt          *local_nvar, *remote_nvar,nroots,nleaves,p=-1,i,nsv;
1460:     const PetscInt    *svtx;
1461:     PetscBool         ghost;
1462:     /*
1463:      PetscMPIInt       rank;
1464:      const PetscInt    *ilocal;
1465:      const PetscSFNode *iremote;
1466:      MPI_Comm_rank(comm,&rank);
1467:      PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,&iremote);
1468:     */
1469:     PetscSFGetGraph(sf,&nroots,&nleaves,NULL,NULL);
1470:     PetscCalloc2(nroots,&local_nvar,nroots,&remote_nvar);

1472:     /* Leaves copy user's nvar to local_nvar */
1473:     DMNetworkGetSharedVertices(dm,&nsv,&svtx);
1474:     for (i=0; i<nsv; i++) {
1475:       p = svtx[i];
1476:       DMNetworkIsGhostVertex(dm,p,&ghost);
1477:       if (!ghost) continue;
1478:       PetscSectionGetDof(network->DofSection,p,&local_nvar[p]);
1479:       /* printf("[%d] Before SFReduce: leaf local_nvar[%d] = %d\n",rank,p,local_nvar[p]); */
1480:     }

1482:     /* Leaves add local_nvar to root remote_nvar */
1483:     PetscSFReduceBegin(sf, MPIU_INT, local_nvar, remote_nvar, MPI_SUM);
1484:     PetscSFReduceEnd(sf, MPIU_INT, local_nvar, remote_nvar, MPI_SUM);
1485:     /* printf("[%d] remote_nvar[%d] = %d\n",rank,p,remote_nvar[p]); */

1487:     /* Update roots' local_nvar */
1488:     for (i=0; i<nsv; i++) {
1489:       p = svtx[i];
1490:       DMNetworkIsGhostVertex(dm,p,&ghost);
1491:       if (ghost) continue;

1493:       PetscSectionAddDof(network->DofSection,p,remote_nvar[p]);
1494:       PetscSectionGetDof(network->DofSection,p,&local_nvar[p]);
1495:       /* printf("[%d]  After SFReduce: root local_nvar[%d] = %d\n",rank,p,local_nvar[p]); */
1496:     }

1498:     /* Roots Bcast nvar to leaves */
1499:     PetscSFBcastBegin(sf, MPIU_INT, local_nvar, remote_nvar,MPI_REPLACE);
1500:     PetscSFBcastEnd(sf, MPIU_INT, local_nvar, remote_nvar,MPI_REPLACE);

1502:     /* Leaves reset receved/remote nvar to dm */
1503:     for (i=0; i<nsv; i++) {
1504:       DMNetworkIsGhostVertex(dm,p,&ghost);
1505:       if (!ghost) continue;
1506:       p = svtx[i];
1507:       /* printf("[%d] leaf reset nvar %d at p= %d \n",rank,remote_nvar[p],p); */
1508:       /* DMNetworkSetNumVariables(dm,p,remote_nvar[p]) */
1509:       PetscSectionSetDof(network->DofSection,p,remote_nvar[p]);
1510:     }

1512:     PetscFree2(local_nvar,remote_nvar);
1513:   }

1515:   PetscSectionSetUp(network->DofSection);
1516:   return(0);
1517: }

1519: /* Get a subsection from a range of points */
1520: static PetscErrorCode DMNetworkGetSubSection_private(PetscSection main,PetscInt pstart,PetscInt pend,PetscSection *subsection)
1521: {
1523:   PetscInt       i, nvar;

1526:   PetscSectionCreate(PetscObjectComm((PetscObject)main), subsection);
1527:   PetscSectionSetChart(*subsection, 0, pend - pstart);
1528:   for (i = pstart; i < pend; i++) {
1529:     PetscSectionGetDof(main,i,&nvar);
1530:     PetscSectionSetDof(*subsection, i - pstart, nvar);
1531:   }

1533:   PetscSectionSetUp(*subsection);
1534:   return(0);
1535: }

1537: /* Create a submap of points with a GlobalToLocal structure */
1538: static PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
1539: {
1541:   PetscInt       i, *subpoints;

1544:   /* Create index sets to map from "points" to "subpoints" */
1545:   PetscMalloc1(pend - pstart, &subpoints);
1546:   for (i = pstart; i < pend; i++) {
1547:     subpoints[i - pstart] = i;
1548:   }
1549:   ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);
1550:   PetscFree(subpoints);
1551:   return(0);
1552: }

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

1557:   Collective on dm

1559:   Input Parameters:
1560: . dm - the DMNetworkObject

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

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

1566:   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]).

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

1570:   Level: intermediate

1572: @*/
1573: PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
1574: {
1576:   MPI_Comm       comm;
1577:   PetscMPIInt    size;
1578:   DM_Network     *network = (DM_Network*)dm->data;

1581:   PetscObjectGetComm((PetscObject)dm,&comm);
1582:   MPI_Comm_size(comm, &size);

1584:   /* Create maps for vertices and edges */
1585:   DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);
1586:   DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);

1588:   /* Create local sub-sections */
1589:   DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);
1590:   DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);

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

1595:     PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);
1596:     PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);
1597:     PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);
1598:   } else {
1599:     /* create structures for vertex */
1600:     PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);
1601:     /* create structures for edge */
1602:     PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);
1603:   }

1605:   /* Add viewers */
1606:   PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");
1607:   PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");
1608:   PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");
1609:   PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");
1610:   return(0);
1611: }

1613: /*
1614:    Add all subnetid for the input vertex v in this process to the btable
1615:    vertex_subnetid = supportingedge_subnetid
1616: */
1617: PETSC_STATIC_INLINE PetscErrorCode SetSubnetIdLookupBT(DM dm,PetscInt v,PetscInt Nsubnet,PetscBT btable)
1618: {
1620:   PetscInt       e,nedges,offset;
1621:   const PetscInt *edges;
1622:   DM_Network     *newDMnetwork = (DM_Network*)dm->data;
1623:   DMNetworkComponentHeader header;

1626:   DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);
1627:   PetscBTMemzero(Nsubnet,btable);
1628:   for (e=0; e<nedges; e++) {
1629:     PetscSectionGetOffset(newDMnetwork->DataSection,edges[e],&offset);
1630:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1631:     PetscBTSet(btable,header->subnetid);
1632:   }
1633:   return(0);
1634: }

1636: /*@
1637:   DMNetworkDistribute - Distributes the network and moves associated component data

1639:   Collective

1641:   Input Parameters:
1642: + DM - the DMNetwork object
1643: - overlap - the overlap of partitions, 0 is the default

1645:   Options Database Key:
1646: + -dmnetwork_view - Calls DMView() at the conclusion of DMSetUp()
1647: - -dmnetwork_view_distributed - Calls DMView() at the conclusion of DMNetworkDistribute()

1649:   Notes:
1650:   Distributes the network with <overlap>-overlapping partitioning of the edges.

1652:   Level: intermediate

1654: .seealso: DMNetworkCreate()
1655: @*/
1656: PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
1657: {
1658:   MPI_Comm       comm;
1660:   PetscMPIInt    size;
1661:   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
1662:   DM_Network     *newDMnetwork;
1663:   PetscSF        pointsf=NULL;
1664:   DM             newDM;
1665:   PetscInt       j,e,v,offset,*subnetvtx,*subnetedge,Nsubnet,gidx,svtx_idx,nv;
1666:   PetscInt       to_net,from_net,*svto;
1667:   PetscBT        btable;
1668:   PetscPartitioner         part;
1669:   DMNetworkComponentHeader header;

1672:   PetscObjectGetComm((PetscObject)*dm,&comm);
1673:   MPI_Comm_size(comm, &size);
1674:   if (size == 1) return(0);

1676:   if (overlap) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"overlap %D != 0 is not supported yet",overlap);

1678:   /* 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. */
1679:   DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);
1680:   newDMnetwork = (DM_Network*)newDM->data;
1681:   newDMnetwork->max_comps_registered = oldDMnetwork->max_comps_registered;
1682:   PetscMalloc1(newDMnetwork->max_comps_registered,&newDMnetwork->component);

1684:   /* Enable runtime options for petscpartitioner */
1685:   DMPlexGetPartitioner(oldDMnetwork->plex,&part);
1686:   PetscPartitionerSetFromOptions(part);

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

1691:   /* Distribute dof section */
1692:   PetscSectionCreate(comm,&newDMnetwork->DofSection);
1693:   PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);

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

1699:   PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);
1700:   DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);
1701:   DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);
1702:   newDMnetwork->nEdges    = newDMnetwork->eEnd - newDMnetwork->eStart;
1703:   newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
1704:   newDMnetwork->NVertices = oldDMnetwork->NVertices;
1705:   newDMnetwork->NEdges    = oldDMnetwork->NEdges;
1706:   newDMnetwork->svtable   = oldDMnetwork->svtable; /* global table! */
1707:   oldDMnetwork->svtable   = NULL;

1709:   /* Set Dof section as the section for dm */
1710:   DMSetLocalSection(newDMnetwork->plex,newDMnetwork->DofSection);
1711:   DMGetGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);

1713:   /* Setup subnetwork info in the newDM */
1714:   newDMnetwork->Nsubnet = oldDMnetwork->Nsubnet;
1715:   newDMnetwork->Nsvtx   = oldDMnetwork->Nsvtx;
1716:   oldDMnetwork->Nsvtx   = 0;
1717:   newDMnetwork->svtx    = oldDMnetwork->svtx; /* global vertices! */
1718:   oldDMnetwork->svtx    = NULL;
1719:   PetscCalloc1(newDMnetwork->Nsubnet,&newDMnetwork->subnet);

1721:   /* Copy over the global number of vertices and edges in each subnetwork.
1722:      Note: these are calculated in DMNetworkLayoutSetUp()
1723:   */
1724:   Nsubnet = newDMnetwork->Nsubnet;
1725:   for (j = 0; j < Nsubnet; j++) {
1726:     newDMnetwork->subnet[j].Nvtx  = oldDMnetwork->subnet[j].Nvtx;
1727:     newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge;
1728:   }

1730:   /* Count local nedges for subnetworks */
1731:   for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++) {
1732:     PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);
1733:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1734:     /* Update pointers */
1735:     header->size          = (PetscInt*)(header + 1);
1736:     header->key           = header->size   + header->maxcomps;
1737:     header->offset        = header->key    + header->maxcomps;
1738:     header->nvar          = header->offset + header->maxcomps;
1739:     header->offsetvarrel  = header->nvar   + header->maxcomps;

1741:     newDMnetwork->subnet[header->subnetid].nedge++;
1742:   }

1744:   /* Count local nvtx for subnetworks */
1745:   PetscBTCreate(Nsubnet,&btable);
1746:   for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++) {
1747:     PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);
1748:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1749:     /* Update pointers */
1750:     header->size          = (PetscInt*)(header + 1);
1751:     header->key           = header->size   + header->maxcomps;
1752:     header->offset        = header->key    + header->maxcomps;
1753:     header->nvar          = header->offset + header->maxcomps;
1754:     header->offsetvarrel  = header->nvar   + header->maxcomps;

1756:     /* shared vertices: use gidx=header->index to check if v is a shared vertex */
1757:     gidx = header->index;
1758:     PetscTableFind(newDMnetwork->svtable,gidx+1,&svtx_idx);
1759:     svtx_idx--;

1761:     if (svtx_idx < 0) { /* not a shared vertex */
1762:       newDMnetwork->subnet[header->subnetid].nvtx++;
1763:     } else { /* a shared vertex belongs to more than one subnetworks, it is being counted by multiple subnets */
1764:       SetSubnetIdLookupBT(newDM,v,Nsubnet,btable);

1766:       from_net = newDMnetwork->svtx[svtx_idx].sv[0];
1767:       if (PetscBTLookup(btable,from_net)) newDMnetwork->subnet[from_net].nvtx++; /* sv is on from_net */

1769:       for (j=1; j<newDMnetwork->svtx[svtx_idx].n; j++) {
1770:         svto   = newDMnetwork->svtx[svtx_idx].sv + 2*j;
1771:         to_net = svto[0];
1772:         if (PetscBTLookup(btable,to_net)) newDMnetwork->subnet[to_net].nvtx++; /* sv is on to_net */
1773:       }
1774:     }
1775:   }

1777:   /* Get total local nvtx for subnetworks */
1778:   nv = 0;
1779:   for (j=0; j<Nsubnet; j++) nv += newDMnetwork->subnet[j].nvtx;
1780:   nv += newDMnetwork->Nsvtx;

1782:   /* Now create the vertices and edge arrays for the subnetworks */
1783:   PetscCalloc2(newDMnetwork->nEdges,&subnetedge,nv,&subnetvtx); /* Maps local vertex to local subnetwork's vertex */
1784:   newDMnetwork->subnetedge = subnetedge;
1785:   newDMnetwork->subnetvtx  = subnetvtx;
1786:   for (j=0; j < newDMnetwork->Nsubnet; j++) {
1787:     newDMnetwork->subnet[j].edges = subnetedge;
1788:     subnetedge                   += newDMnetwork->subnet[j].nedge;

1790:     newDMnetwork->subnet[j].vertices = subnetvtx;
1791:     subnetvtx                       += newDMnetwork->subnet[j].nvtx;

1793:     /* 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. */
1794:     newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0;
1795:   }
1796:   newDMnetwork->svertices = subnetvtx;

1798:   /* Set the edges and vertices in each subnetwork */
1799:   for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++) {
1800:     PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);
1801:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1802:     newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++] = e;
1803:   }

1805:   nv = 0;
1806:   for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++) {
1807:     PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);
1808:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);

1810:     /* coupling vertices: use gidx = header->index to check if v is a coupling vertex */
1811:     PetscTableFind(newDMnetwork->svtable,header->index+1,&svtx_idx);
1812:     svtx_idx--;
1813:     if (svtx_idx < 0) {
1814:       newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++] = v;
1815:     } else { /* a shared vertex */
1816:       newDMnetwork->svertices[nv++] = v;

1818:       /* add all subnetid for this shared vertex in this process to btable */
1819:       SetSubnetIdLookupBT(newDM,v,Nsubnet,btable);

1821:       from_net = newDMnetwork->svtx[svtx_idx].sv[0];
1822:       if (PetscBTLookup(btable,from_net))
1823:         newDMnetwork->subnet[from_net].vertices[newDMnetwork->subnet[from_net].nvtx++] = v;

1825:       for (j=1; j<newDMnetwork->svtx[svtx_idx].n; j++) {
1826:         svto   = newDMnetwork->svtx[svtx_idx].sv + 2*j;
1827:         to_net = svto[0];
1828:         if (PetscBTLookup(btable,to_net))
1829:           newDMnetwork->subnet[to_net].vertices[newDMnetwork->subnet[to_net].nvtx++] = v;
1830:       }
1831:     }
1832:   }
1833:   newDMnetwork->nsvtx = nv;   /* num of local shared vertices */

1835:   newDM->setupcalled = (*dm)->setupcalled;
1836:   newDMnetwork->distributecalled = PETSC_TRUE;

1838:   /* Free spaces */
1839:   PetscSFDestroy(&pointsf);
1840:   DMDestroy(dm);
1841:   PetscBTDestroy(&btable);

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

1846:   *dm  = newDM;
1847:   return(0);
1848: }

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

1853:  Collective

1855:   Input Parameters:
1856: + mainSF - the original SF structure
1857: - map - a ISLocalToGlobal mapping that contains the subset of points

1859:   Output Parameter:
1860: . subSF - a subset of the mainSF for the desired subset.

1862:   Level: intermediate
1863: @*/
1864: PetscErrorCode PetscSFGetSubSF(PetscSF mainsf,ISLocalToGlobalMapping map,PetscSF *subSF)
1865: {
1866:   PetscErrorCode        ierr;
1867:   PetscInt              nroots, nleaves, *ilocal_sub;
1868:   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
1869:   PetscInt              *local_points, *remote_points;
1870:   PetscSFNode           *iremote_sub;
1871:   const PetscInt        *ilocal;
1872:   const PetscSFNode     *iremote;

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

1877:   /* Look for leaves that pertain to the subset of points. Get the local ordering */
1878:   PetscMalloc1(nleaves,&ilocal_map);
1879:   ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);
1880:   for (i = 0; i < nleaves; i++) {
1881:     if (ilocal_map[i] != -1) nleaves_sub += 1;
1882:   }
1883:   /* Re-number ilocal with subset numbering. Need information from roots */
1884:   PetscMalloc2(nroots,&local_points,nroots,&remote_points);
1885:   for (i = 0; i < nroots; i++) local_points[i] = i;
1886:   ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);
1887:   PetscSFBcastBegin(mainsf, MPIU_INT, local_points, remote_points,MPI_REPLACE);
1888:   PetscSFBcastEnd(mainsf, MPIU_INT, local_points, remote_points,MPI_REPLACE);
1889:   /* Fill up graph using local (that is, local to the subset) numbering. */
1890:   PetscMalloc1(nleaves_sub,&ilocal_sub);
1891:   PetscMalloc1(nleaves_sub,&iremote_sub);
1892:   nleaves_sub = 0;
1893:   for (i = 0; i < nleaves; i++) {
1894:     if (ilocal_map[i] != -1) {
1895:       ilocal_sub[nleaves_sub] = ilocal_map[i];
1896:       iremote_sub[nleaves_sub].rank = iremote[i].rank;
1897:       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
1898:       nleaves_sub += 1;
1899:     }
1900:   }
1901:   PetscFree2(local_points,remote_points);
1902:   ISLocalToGlobalMappingGetSize(map,&nroots_sub);

1904:   /* Create new subSF */
1905:   PetscSFCreate(PETSC_COMM_WORLD,subSF);
1906:   PetscSFSetFromOptions(*subSF);
1907:   PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);
1908:   PetscFree(ilocal_map);
1909:   PetscFree(iremote_sub);
1910:   return(0);
1911: }

1913: /*@C
1914:   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point

1916:   Not Collective

1918:   Input Parameters:
1919: + dm - the DMNetwork object
1920: - p  - the vertex point

1922:   Output Parameters:
1923: + nedges - number of edges connected to this vertex point
1924: - edges  - list of edge points

1926:   Level: beginner

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

1932: .seealso: DMNetworkCreate(), DMNetworkGetConnectedVertices()
1933: @*/
1934: PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
1935: {
1937:   DM_Network     *network = (DM_Network*)dm->data;

1940:   DMPlexGetSupportSize(network->plex,vertex,nedges);
1941:   if (edges) {DMPlexGetSupport(network->plex,vertex,edges);}
1942:   return(0);
1943: }

1945: /*@C
1946:   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point

1948:   Not Collective

1950:   Input Parameters:
1951: + dm - the DMNetwork object
1952: - p - the edge point

1954:   Output Parameters:
1955: . vertices - vertices connected to this edge

1957:   Level: beginner

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

1963: .seealso: DMNetworkCreate(), DMNetworkGetSupportingEdges()
1964: @*/
1965: PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
1966: {
1968:   DM_Network     *network = (DM_Network*)dm->data;

1971:   DMPlexGetCone(network->plex,edge,vertices);
1972:   return(0);
1973: }

1975: /*@
1976:   DMNetworkIsSharedVertex - Returns TRUE if the vertex is shared by subnetworks

1978:   Not Collective

1980:   Input Parameters:
1981: + dm - the DMNetwork object
1982: - p - the vertex point

1984:   Output Parameter:
1985: . flag - TRUE if the vertex is shared by subnetworks

1987:   Level: beginner

1989: .seealso: DMNetworkAddSharedVertices(), DMNetworkIsGhostVertex()
1990: @*/
1991: PetscErrorCode DMNetworkIsSharedVertex(DM dm,PetscInt p,PetscBool *flag)
1992: {
1994:   PetscInt       i;

1997:   *flag = PETSC_FALSE;

1999:   if (dm->setupcalled) { /* DMNetworkGetGlobalVertexIndex() requires DMSetUp() be called */
2000:     DM_Network     *network = (DM_Network*)dm->data;
2001:     PetscInt       gidx;
2002:     DMNetworkGetGlobalVertexIndex(dm,p,&gidx);
2003:     PetscTableFind(network->svtable,gidx+1,&i);
2004:     if (i) *flag = PETSC_TRUE;
2005:   } else { /* would be removed? */
2006:     PetscInt       nv;
2007:     const PetscInt *vtx;
2008:     DMNetworkGetSharedVertices(dm,&nv,&vtx);
2009:     for (i=0; i<nv; i++) {
2010:       if (p == vtx[i]) {
2011:         *flag = PETSC_TRUE;
2012:         break;
2013:       }
2014:     }
2015:   }
2016:   return(0);
2017: }

2019: /*@
2020:   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex

2022:   Not Collective

2024:   Input Parameters:
2025: + dm - the DMNetwork object
2026: - p - the vertex point

2028:   Output Parameter:
2029: . isghost - TRUE if the vertex is a ghost point

2031:   Level: beginner

2033: .seealso: DMNetworkGetConnectedVertices(), DMNetworkGetVertexRange(), DMNetworkIsSharedVertex()
2034: @*/
2035: PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
2036: {
2038:   DM_Network     *network = (DM_Network*)dm->data;
2039:   PetscInt       offsetg;
2040:   PetscSection   sectiong;

2043:   *isghost = PETSC_FALSE;
2044:   DMGetGlobalSection(network->plex,&sectiong);
2045:   PetscSectionGetOffset(sectiong,p,&offsetg);
2046:   if (offsetg < 0) *isghost = PETSC_TRUE;
2047:   return(0);
2048: }

2050: PetscErrorCode DMSetUp_Network(DM dm)
2051: {
2053:   DM_Network     *network=(DM_Network*)dm->data;

2056:   DMNetworkComponentSetUp(dm);
2057:   DMNetworkVariablesSetUp(dm);

2059:   DMSetLocalSection(network->plex,network->DofSection);
2060:   DMGetGlobalSection(network->plex,&network->GlobalDofSection);

2062:   dm->setupcalled = PETSC_TRUE;

2064:   /* View dmnetwork */
2065:   DMViewFromOptions(dm,NULL,"-dmnetwork_view");
2066:   return(0);
2067: }

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

2073:   Collective

2075:   Input Parameters:
2076: + dm - the DMNetwork object
2077: . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
2078: - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices

2080:  Level: intermediate

2082: @*/
2083: PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
2084: {
2085:   DM_Network     *network=(DM_Network*)dm->data;
2087:   PetscInt       nVertices = network->nVertices;

2090:   network->userEdgeJacobian   = eflg;
2091:   network->userVertexJacobian = vflg;

2093:   if (eflg && !network->Je) {
2094:     PetscCalloc1(3*network->nEdges,&network->Je);
2095:   }

2097:   if (vflg && !network->Jv && nVertices) {
2098:     PetscInt       i,*vptr,nedges,vStart=network->vStart;
2099:     PetscInt       nedges_total;
2100:     const PetscInt *edges;

2102:     /* count nvertex_total */
2103:     nedges_total = 0;
2104:     PetscMalloc1(nVertices+1,&vptr);

2106:     vptr[0] = 0;
2107:     for (i=0; i<nVertices; i++) {
2108:       DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);
2109:       nedges_total += nedges;
2110:       vptr[i+1] = vptr[i] + 2*nedges + 1;
2111:     }

2113:     PetscCalloc1(2*nedges_total+nVertices,&network->Jv);
2114:     network->Jvptr = vptr;
2115:   }
2116:   return(0);
2117: }

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

2122:   Not Collective

2124:   Input Parameters:
2125: + dm - the DMNetwork object
2126: . p - the edge point
2127: - J - array (size = 3) of Jacobian submatrices for this edge point:
2128:         J[0]: this edge
2129:         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()

2131:   Level: advanced

2133: .seealso: DMNetworkVertexSetMatrix()
2134: @*/
2135: PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
2136: {
2137:   DM_Network *network=(DM_Network*)dm->data;

2140:   if (!network->Je) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");

2142:   if (J) {
2143:     network->Je[3*p]   = J[0];
2144:     network->Je[3*p+1] = J[1];
2145:     network->Je[3*p+2] = J[2];
2146:   }
2147:   return(0);
2148: }

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

2153:   Not Collective

2155:   Input Parameters:
2156: + dm - The DMNetwork object
2157: . p - the vertex point
2158: - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
2159:         J[0]:       this vertex
2160:         J[1+2*i]:   i-th supporting edge
2161:         J[1+2*i+1]: i-th connected vertex

2163:   Level: advanced

2165: .seealso: DMNetworkEdgeSetMatrix()
2166: @*/
2167: PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
2168: {
2170:   DM_Network     *network=(DM_Network*)dm->data;
2171:   PetscInt       i,*vptr,nedges,vStart=network->vStart;
2172:   const PetscInt *edges;

2175:   if (!network->Jv) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");

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

2181:     /* Set Jacobian for each supporting edge and connected vertex */
2182:     DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);
2183:     for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
2184:   }
2185:   return(0);
2186: }

2188: PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
2189: {
2191:   PetscInt       j;
2192:   PetscScalar    val=(PetscScalar)ncols;

2195:   if (!ghost) {
2196:     for (j=0; j<nrows; j++) {
2197:       VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);
2198:     }
2199:   } else {
2200:     for (j=0; j<nrows; j++) {
2201:       VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);
2202:     }
2203:   }
2204:   return(0);
2205: }

2207: PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
2208: {
2210:   PetscInt       j,ncols_u;
2211:   PetscScalar    val;

2214:   if (!ghost) {
2215:     for (j=0; j<nrows; j++) {
2216:       MatGetRow(Ju,j,&ncols_u,NULL,NULL);
2217:       val = (PetscScalar)ncols_u;
2218:       VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);
2219:       MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);
2220:     }
2221:   } else {
2222:     for (j=0; j<nrows; j++) {
2223:       MatGetRow(Ju,j,&ncols_u,NULL,NULL);
2224:       val = (PetscScalar)ncols_u;
2225:       VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);
2226:       MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);
2227:     }
2228:   }
2229:   return(0);
2230: }

2232: PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
2233: {

2237:   if (Ju) {
2238:     MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);
2239:   } else {
2240:     MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);
2241:   }
2242:   return(0);
2243: }

2245: PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
2246: {
2248:   PetscInt       j,*cols;
2249:   PetscScalar    *zeros;

2252:   PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);
2253:   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
2254:   MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);
2255:   PetscFree2(cols,zeros);
2256:   return(0);
2257: }

2259: PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
2260: {
2262:   PetscInt       j,M,N,row,col,ncols_u;
2263:   const PetscInt *cols;
2264:   PetscScalar    zero=0.0;

2267:   MatGetSize(Ju,&M,&N);
2268:   if (nrows != M || ncols != N) SETERRQ4(PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%D by %D must equal %D by %D",nrows,ncols,M,N);

2270:   for (row=0; row<nrows; row++) {
2271:     MatGetRow(Ju,row,&ncols_u,&cols,NULL);
2272:     for (j=0; j<ncols_u; j++) {
2273:       col = cols[j] + cstart;
2274:       MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);
2275:     }
2276:     MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);
2277:   }
2278:   return(0);
2279: }

2281: PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
2282: {

2286:   if (Ju) {
2287:     MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);
2288:   } else {
2289:     MatSetDenseblock_private(nrows,rows,ncols,cstart,J);
2290:   }
2291:   return(0);
2292: }

2294: /* 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.
2295: */
2296: PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
2297: {
2299:   PetscInt       i,size,dof;
2300:   PetscInt       *glob2loc;

2303:   PetscSectionGetStorageSize(localsec,&size);
2304:   PetscMalloc1(size,&glob2loc);

2306:   for (i = 0; i < size; i++) {
2307:     PetscSectionGetOffset(globalsec,i,&dof);
2308:     dof = (dof >= 0) ? dof : -(dof + 1);
2309:     glob2loc[i] = dof;
2310:   }

2312:   ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);
2313: #if 0
2314:   PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);
2315: #endif
2316:   return(0);
2317: }

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

2321: PetscErrorCode DMCreateMatrix_Network_Nest(DM dm,Mat *J)
2322: {
2324:   DM_Network     *network = (DM_Network*)dm->data;
2325:   PetscInt       eDof,vDof;
2326:   Mat            j11,j12,j21,j22,bA[2][2];
2327:   MPI_Comm       comm;
2328:   ISLocalToGlobalMapping eISMap,vISMap;

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

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

2336:   MatCreate(comm, &j11);
2337:   MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);
2338:   MatSetType(j11, MATMPIAIJ);

2340:   MatCreate(comm, &j12);
2341:   MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);
2342:   MatSetType(j12, MATMPIAIJ);

2344:   MatCreate(comm, &j21);
2345:   MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);
2346:   MatSetType(j21, MATMPIAIJ);

2348:   MatCreate(comm, &j22);
2349:   MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);
2350:   MatSetType(j22, MATMPIAIJ);

2352:   bA[0][0] = j11;
2353:   bA[0][1] = j12;
2354:   bA[1][0] = j21;
2355:   bA[1][1] = j22;

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

2360:   MatSetLocalToGlobalMapping(j11,eISMap,eISMap);
2361:   MatSetLocalToGlobalMapping(j12,eISMap,vISMap);
2362:   MatSetLocalToGlobalMapping(j21,vISMap,eISMap);
2363:   MatSetLocalToGlobalMapping(j22,vISMap,vISMap);

2365:   MatSetUp(j11);
2366:   MatSetUp(j12);
2367:   MatSetUp(j21);
2368:   MatSetUp(j22);

2370:   MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);
2371:   MatSetUp(*J);
2372:   MatNestSetVecType(*J,VECNEST);
2373:   MatDestroy(&j11);
2374:   MatDestroy(&j12);
2375:   MatDestroy(&j21);
2376:   MatDestroy(&j22);

2378:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
2379:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
2380:   MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);

2382:   /* Free structures */
2383:   ISLocalToGlobalMappingDestroy(&eISMap);
2384:   ISLocalToGlobalMappingDestroy(&vISMap);
2385:   return(0);
2386: }

2388: PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
2389: {
2391:   DM_Network     *network = (DM_Network*)dm->data;
2392:   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
2393:   PetscInt       cstart,ncols,j,e,v;
2394:   PetscBool      ghost,ghost_vc,ghost2,isNest;
2395:   Mat            Juser;
2396:   PetscSection   sectionGlobal;
2397:   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
2398:   const PetscInt *edges,*cone;
2399:   MPI_Comm       comm;
2400:   MatType        mtype;
2401:   Vec            vd_nz,vo_nz;
2402:   PetscInt       *dnnz,*onnz;
2403:   PetscScalar    *vdnz,*vonz;

2406:   mtype = dm->mattype;
2407:   PetscStrcmp(mtype,MATNEST,&isNest);
2408:   if (isNest) {
2409:     DMCreateMatrix_Network_Nest(dm,J);
2410:     MatSetDM(*J,dm);
2411:     return(0);
2412:   }

2414:   if (!network->userEdgeJacobian && !network->userVertexJacobian) {
2415:     /* user does not provide Jacobian blocks */
2416:     DMCreateMatrix_Plex(network->plex,J);
2417:     MatSetDM(*J,dm);
2418:     return(0);
2419:   }

2421:   MatCreate(PetscObjectComm((PetscObject)dm),J);
2422:   DMGetGlobalSection(network->plex,&sectionGlobal);
2423:   PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);
2424:   MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);

2426:   MatSetType(*J,MATAIJ);
2427:   MatSetFromOptions(*J);

2429:   /* (1) Set matrix preallocation */
2430:   /*------------------------------*/
2431:   PetscObjectGetComm((PetscObject)dm,&comm);
2432:   VecCreate(comm,&vd_nz);
2433:   VecSetSizes(vd_nz,localSize,PETSC_DECIDE);
2434:   VecSetFromOptions(vd_nz);
2435:   VecSet(vd_nz,0.0);
2436:   VecDuplicate(vd_nz,&vo_nz);

2438:   /* Set preallocation for edges */
2439:   /*-----------------------------*/
2440:   DMNetworkGetEdgeRange(dm,&eStart,&eEnd);

2442:   PetscMalloc1(localSize,&rows);
2443:   for (e=eStart; e<eEnd; e++) {
2444:     /* Get row indices */
2445:     DMNetworkGetGlobalVecOffset(dm,e,ALL_COMPONENTS,&rstart);
2446:     PetscSectionGetDof(network->DofSection,e,&nrows);
2447:     if (nrows) {
2448:       for (j=0; j<nrows; j++) rows[j] = j + rstart;

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

2455:         if (network->Je) {
2456:           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
2457:         } else Juser = NULL;
2458:         DMNetworkIsGhostVertex(dm,cone[v],&ghost);
2459:         MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);
2460:       }

2462:       /* Set preallocation for edge self */
2463:       cstart = rstart;
2464:       if (network->Je) {
2465:         Juser = network->Je[3*e]; /* Jacobian(e,e) */
2466:       } else Juser = NULL;
2467:       MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);
2468:     }
2469:   }

2471:   /* Set preallocation for vertices */
2472:   /*--------------------------------*/
2473:   DMNetworkGetVertexRange(dm,&vStart,&vEnd);
2474:   if (vEnd - vStart) vptr = network->Jvptr;

2476:   for (v=vStart; v<vEnd; v++) {
2477:     /* Get row indices */
2478:     DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&rstart);
2479:     PetscSectionGetDof(network->DofSection,v,&nrows);
2480:     if (!nrows) continue;

2482:     DMNetworkIsGhostVertex(dm,v,&ghost);
2483:     if (ghost) {
2484:       PetscMalloc1(nrows,&rows_v);
2485:     } else {
2486:       rows_v = rows;
2487:     }

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

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

2494:     for (e=0; e<nedges; e++) {
2495:       /* Supporting edges */
2496:       DMNetworkGetGlobalVecOffset(dm,edges[e],ALL_COMPONENTS,&cstart);
2497:       PetscSectionGetDof(network->DofSection,edges[e],&ncols);

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

2504:       /* Connected vertices */
2505:       DMNetworkGetConnectedVertices(dm,edges[e],&cone);
2506:       vc = (v == cone[0]) ? cone[1]:cone[0];
2507:       DMNetworkIsGhostVertex(dm,vc,&ghost_vc);

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

2511:       if (network->Jv) {
2512:         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
2513:       } else Juser = NULL;
2514:       if (ghost_vc||ghost) {
2515:         ghost2 = PETSC_TRUE;
2516:       } else {
2517:         ghost2 = PETSC_FALSE;
2518:       }
2519:       MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);
2520:     }

2522:     /* Set preallocation for vertex self */
2523:     DMNetworkIsGhostVertex(dm,v,&ghost);
2524:     if (!ghost) {
2525:       DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&cstart);
2526:       if (network->Jv) {
2527:         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
2528:       } else Juser = NULL;
2529:       MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);
2530:     }
2531:     if (ghost) {
2532:       PetscFree(rows_v);
2533:     }
2534:   }

2536:   VecAssemblyBegin(vd_nz);
2537:   VecAssemblyBegin(vo_nz);

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

2541:   VecAssemblyEnd(vd_nz);
2542:   VecAssemblyEnd(vo_nz);

2544:   VecGetArray(vd_nz,&vdnz);
2545:   VecGetArray(vo_nz,&vonz);
2546:   for (j=0; j<localSize; j++) {
2547:     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
2548:     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
2549:   }
2550:   VecRestoreArray(vd_nz,&vdnz);
2551:   VecRestoreArray(vo_nz,&vonz);
2552:   VecDestroy(&vd_nz);
2553:   VecDestroy(&vo_nz);

2555:   MatSeqAIJSetPreallocation(*J,0,dnnz);
2556:   MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);
2557:   MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);

2559:   PetscFree2(dnnz,onnz);

2561:   /* (2) Set matrix entries for edges */
2562:   /*----------------------------------*/
2563:   for (e=eStart; e<eEnd; e++) {
2564:     /* Get row indices */
2565:     DMNetworkGetGlobalVecOffset(dm,e,ALL_COMPONENTS,&rstart);
2566:     PetscSectionGetDof(network->DofSection,e,&nrows);
2567:     if (nrows) {
2568:       for (j=0; j<nrows; j++) rows[j] = j + rstart;

2570:       /* Set matrix entries for connected vertices */
2571:       DMNetworkGetConnectedVertices(dm,e,&cone);
2572:       for (v=0; v<2; v++) {
2573:         DMNetworkGetGlobalVecOffset(dm,cone[v],ALL_COMPONENTS,&cstart);
2574:         PetscSectionGetDof(network->DofSection,cone[v],&ncols);

2576:         if (network->Je) {
2577:           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
2578:         } else Juser = NULL;
2579:         MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);
2580:       }

2582:       /* Set matrix entries for edge self */
2583:       cstart = rstart;
2584:       if (network->Je) {
2585:         Juser = network->Je[3*e]; /* Jacobian(e,e) */
2586:       } else Juser = NULL;
2587:       MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);
2588:     }
2589:   }

2591:   /* Set matrix entries for vertices */
2592:   /*---------------------------------*/
2593:   for (v=vStart; v<vEnd; v++) {
2594:     /* Get row indices */
2595:     DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&rstart);
2596:     PetscSectionGetDof(network->DofSection,v,&nrows);
2597:     if (!nrows) continue;

2599:     DMNetworkIsGhostVertex(dm,v,&ghost);
2600:     if (ghost) {
2601:       PetscMalloc1(nrows,&rows_v);
2602:     } else {
2603:       rows_v = rows;
2604:     }
2605:     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;

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

2610:     for (e=0; e<nedges; e++) {
2611:       /* Supporting edges */
2612:       DMNetworkGetGlobalVecOffset(dm,edges[e],ALL_COMPONENTS,&cstart);
2613:       PetscSectionGetDof(network->DofSection,edges[e],&ncols);

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

2620:       /* Connected vertices */
2621:       DMNetworkGetConnectedVertices(dm,edges[e],&cone);
2622:       vc = (v == cone[0]) ? cone[1]:cone[0];

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

2627:       if (network->Jv) {
2628:         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
2629:       } else Juser = NULL;
2630:       MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);
2631:     }

2633:     /* Set matrix entries for vertex self */
2634:     if (!ghost) {
2635:       DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&cstart);
2636:       if (network->Jv) {
2637:         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
2638:       } else Juser = NULL;
2639:       MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);
2640:     }
2641:     if (ghost) {
2642:       PetscFree(rows_v);
2643:     }
2644:   }
2645:   PetscFree(rows);

2647:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
2648:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);

2650:   MatSetDM(*J,dm);
2651:   return(0);
2652: }

2654: PetscErrorCode DMDestroy_Network(DM dm)
2655: {
2657:   DM_Network     *network = (DM_Network*)dm->data;
2658:   PetscInt       j,np;

2661:   if (--network->refct > 0) return(0);
2662:   if (network->Je) {
2663:     PetscFree(network->Je);
2664:   }
2665:   if (network->Jv) {
2666:     PetscFree(network->Jvptr);
2667:     PetscFree(network->Jv);
2668:   }

2670:   ISLocalToGlobalMappingDestroy(&network->vertex.mapping);
2671:   PetscSectionDestroy(&network->vertex.DofSection);
2672:   PetscSectionDestroy(&network->vertex.GlobalDofSection);
2673:   if (network->vltog) {
2674:     PetscFree(network->vltog);
2675:   }
2676:   if (network->vertex.sf) {
2677:     PetscSFDestroy(&network->vertex.sf);
2678:   }
2679:   /* edge */
2680:   ISLocalToGlobalMappingDestroy(&network->edge.mapping);
2681:   PetscSectionDestroy(&network->edge.DofSection);
2682:   PetscSectionDestroy(&network->edge.GlobalDofSection);
2683:   if (network->edge.sf) {
2684:     PetscSFDestroy(&network->edge.sf);
2685:   }
2686:   DMDestroy(&network->plex);
2687:   PetscSectionDestroy(&network->DataSection);
2688:   PetscSectionDestroy(&network->DofSection);

2690:   for (j=0; j<network->Nsvtx; j++) {
2691:     PetscFree(network->svtx[j].sv);
2692:   }
2693:   if (network->svtx) {PetscFree(network->svtx);}
2694:   PetscFree2(network->subnetedge,network->subnetvtx);

2696:   PetscTableDestroy(&network->svtable);
2697:   PetscFree(network->subnet);
2698:   PetscFree(network->component);
2699:   PetscFree(network->componentdataarray);

2701:   if (network->header) {
2702:     np = network->pEnd - network->pStart;
2703:     for (j=0; j < np; j++) {
2704:       PetscFree5(network->header[j].size,network->header[j].key,network->header[j].offset,network->header[j].nvar,network->header[j].offsetvarrel);
2705:       PetscFree(network->cvalue[j].data);
2706:     }
2707:     PetscFree2(network->header,network->cvalue);
2708:   }
2709:   PetscFree(network);
2710:   return(0);
2711: }

2713: PetscErrorCode DMView_Network(DM dm,PetscViewer viewer)
2714: {
2716:   PetscBool      iascii;
2717:   PetscMPIInt    rank;

2720:   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first");
2721:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
2724:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
2725:   if (iascii) {
2726:     const PetscInt *cone,*vtx,*edges;
2727:     PetscInt       vfrom,vto,i,j,nv,ne,ncv,p,nsubnet;
2728:     DM_Network     *network = (DM_Network*)dm->data;

2730:     nsubnet = network->Nsubnet; /* num of subnetworks */
2731:     if (rank == 0) {
2732:       PetscPrintf(PETSC_COMM_SELF,"  NSubnets: %D; NEdges: %D; NVertices: %D; NSharedVertices: %D.\n",nsubnet,network->NEdges,network->NVertices,network->Nsvtx);
2733:     }

2735:     DMNetworkGetSharedVertices(dm,&ncv,&vtx);
2736:     PetscViewerASCIIPushSynchronized(viewer);
2737:     PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] nEdges: %D; nVertices: %D; nSharedVertices: %D\n",rank,network->nEdges,network->nVertices,ncv);

2739:     for (i=0; i<nsubnet; i++) {
2740:       DMNetworkGetSubnetwork(dm,i,&nv,&ne,&vtx,&edges);
2741:       if (ne) {
2742:         PetscViewerASCIISynchronizedPrintf(viewer, "     Subnet %D: nEdges %D, nVertices(include shared vertices) %D\n",i,ne,nv);
2743:         for (j=0; j<ne; j++) {
2744:           p = edges[j];
2745:           DMNetworkGetConnectedVertices(dm,p,&cone);
2746:           DMNetworkGetGlobalVertexIndex(dm,cone[0],&vfrom);
2747:           DMNetworkGetGlobalVertexIndex(dm,cone[1],&vto);
2748:           DMNetworkGetGlobalEdgeIndex(dm,edges[j],&p);
2749:           PetscViewerASCIISynchronizedPrintf(viewer, "       edge %D: %D ----> %D\n",p,vfrom,vto);
2750:         }
2751:       }
2752:     }

2754:     /* Shared vertices */
2755:     DMNetworkGetSharedVertices(dm,&ncv,&vtx);
2756:     if (ncv) {
2757:       SVtx       *svtx = network->svtx;
2758:       PetscInt    gidx,svtx_idx,nvto,vfrom_net,vfrom_idx,*svto;
2759:       PetscBool   ghost;
2760:       PetscViewerASCIISynchronizedPrintf(viewer, "     SharedVertices:\n");
2761:       for (i=0; i<ncv; i++) {
2762:         DMNetworkIsGhostVertex(dm,vtx[i],&ghost);
2763:         if (ghost) continue;

2765:         DMNetworkGetGlobalVertexIndex(dm,vtx[i],&gidx);
2766:         PetscTableFind(network->svtable,gidx+1,&svtx_idx);
2767:         svtx_idx--;
2768:         nvto = svtx[svtx_idx].n;

2770:         vfrom_net = svtx[svtx_idx].sv[0];
2771:         vfrom_idx = svtx[svtx_idx].sv[1];
2772:         PetscViewerASCIISynchronizedPrintf(viewer, "       svtx %D: global index %D, subnet[%D].%D ---->\n",i,gidx,vfrom_net,vfrom_idx);
2773:         for (j=1; j<nvto; j++) {
2774:           svto = svtx[svtx_idx].sv + 2*j;
2775:           PetscViewerASCIISynchronizedPrintf(viewer, "                                           ----> subnet[%D].%D\n",svto[0],svto[1]);
2776:         }
2777:       }
2778:     }
2779:     PetscViewerFlush(viewer);
2780:     PetscViewerASCIIPopSynchronized(viewer);
2781:   } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMNetwork writing", ((PetscObject)viewer)->type_name);
2782:   return(0);
2783: }

2785: PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
2786: {
2788:   DM_Network     *network = (DM_Network*)dm->data;

2791:   DMGlobalToLocalBegin(network->plex,g,mode,l);
2792:   return(0);
2793: }

2795: PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
2796: {
2798:   DM_Network     *network = (DM_Network*)dm->data;

2801:   DMGlobalToLocalEnd(network->plex,g,mode,l);
2802:   return(0);
2803: }

2805: PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
2806: {
2808:   DM_Network     *network = (DM_Network*)dm->data;

2811:   DMLocalToGlobalBegin(network->plex,l,mode,g);
2812:   return(0);
2813: }

2815: PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
2816: {
2818:   DM_Network     *network = (DM_Network*)dm->data;

2821:   DMLocalToGlobalEnd(network->plex,l,mode,g);
2822:   return(0);
2823: }

2825: /*@
2826:   DMNetworkGetVertexLocalToGlobalOrdering - Get vertex global index

2828:   Not collective

2830:   Input Parameters:
2831: + dm - the dm object
2832: - vloc - local vertex ordering, start from 0

2834:   Output Parameters:
2835: .  vg  - global vertex ordering, start from 0

2837:   Level: advanced

2839: .seealso: DMNetworkSetVertexLocalToGlobalOrdering()
2840: @*/
2841: PetscErrorCode DMNetworkGetVertexLocalToGlobalOrdering(DM dm,PetscInt vloc,PetscInt *vg)
2842: {
2843:   DM_Network  *network = (DM_Network*)dm->data;
2844:   PetscInt    *vltog = network->vltog;

2847:   if (!vltog) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Must call DMNetworkSetVertexLocalToGlobalOrdering() first");
2848:   *vg = vltog[vloc];
2849:   return(0);
2850: }

2852: /*@
2853:   DMNetworkSetVertexLocalToGlobalOrdering - Create and setup vertex local to global map

2855:   Collective

2857:   Input Parameters:
2858: . dm - the dm object

2860:   Level: advanced

2862: .seealso: DMNetworkGetGlobalVertexIndex()
2863: @*/
2864: PetscErrorCode DMNetworkSetVertexLocalToGlobalOrdering(DM dm)
2865: {
2866:   PetscErrorCode    ierr;
2867:   DM_Network        *network=(DM_Network*)dm->data;
2868:   MPI_Comm          comm;
2869:   PetscMPIInt       rank,size,*displs=NULL,*recvcounts=NULL,remoterank;
2870:   PetscBool         ghost;
2871:   PetscInt          *vltog,nroots,nleaves,i,*vrange,k,N,lidx;
2872:   const PetscSFNode *iremote;
2873:   PetscSF           vsf;
2874:   Vec               Vleaves,Vleaves_seq;
2875:   VecScatter        ctx;
2876:   PetscScalar       *varr,val;
2877:   const PetscScalar *varr_read;

2880:   PetscObjectGetComm((PetscObject)dm,&comm);
2881:   MPI_Comm_size(comm,&size);
2882:   MPI_Comm_rank(comm,&rank);

2884:   if (size == 1) {
2885:     nroots = network->vEnd - network->vStart;
2886:     PetscMalloc1(nroots, &vltog);
2887:     for (i=0; i<nroots; i++) vltog[i] = i;
2888:     network->vltog = vltog;
2889:     return(0);
2890:   }

2892:   if (!network->distributecalled) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE,"Must call DMNetworkDistribute() first");
2893:   if (network->vltog) {
2894:     PetscFree(network->vltog);
2895:   }

2897:   DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);
2898:   PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);
2899:   vsf = network->vertex.sf;

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

2904:   for (i=0; i<size; i++) { displs[i] = i; recvcounts[i] = 1;}

2906:   i         = nroots - nleaves; /* local number of vertices, excluding ghosts */
2907:   vrange[0] = 0;
2908:   MPI_Allgatherv(&i,1,MPIU_INT,vrange+1,recvcounts,displs,MPIU_INT,comm);
2909:   for (i=2; i<size+1; i++) {vrange[i] += vrange[i-1];}

2911:   PetscMalloc1(nroots, &vltog);
2912:   network->vltog = vltog;

2914:   /* Set vltog for non-ghost vertices */
2915:   k = 0;
2916:   for (i=0; i<nroots; i++) {
2917:     DMNetworkIsGhostVertex(dm,i+network->vStart,&ghost);
2918:     if (ghost) continue;
2919:     vltog[i] = vrange[rank] + k++;
2920:   }
2921:   PetscFree3(vrange,displs,recvcounts);

2923:   /* Set vltog for ghost vertices */
2924:   /* (a) create parallel Vleaves and sequential Vleaves_seq to convert local iremote[*].index to global index */
2925:   VecCreate(comm,&Vleaves);
2926:   VecSetSizes(Vleaves,2*nleaves,PETSC_DETERMINE);
2927:   VecSetFromOptions(Vleaves);
2928:   VecGetArray(Vleaves,&varr);
2929:   for (i=0; i<nleaves; i++) {
2930:     varr[2*i]   = (PetscScalar)(iremote[i].rank);  /* rank of remote process */
2931:     varr[2*i+1] = (PetscScalar)(iremote[i].index); /* local index in remote process */
2932:   }
2933:   VecRestoreArray(Vleaves,&varr);

2935:   /* (b) scatter local info to remote processes via VecScatter() */
2936:   VecScatterCreateToAll(Vleaves,&ctx,&Vleaves_seq);
2937:   VecScatterBegin(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD);
2938:   VecScatterEnd(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD);

2940:   /* (c) convert local indices to global indices in parallel vector Vleaves */
2941:   VecGetSize(Vleaves_seq,&N);
2942:   VecGetArrayRead(Vleaves_seq,&varr_read);
2943:   for (i=0; i<N; i+=2) {
2944:     remoterank = (PetscMPIInt)PetscRealPart(varr_read[i]);
2945:     if (remoterank == rank) {
2946:       k = i+1; /* row number */
2947:       lidx = (PetscInt)PetscRealPart(varr_read[i+1]);
2948:       val  = (PetscScalar)vltog[lidx]; /* global index for non-ghost vertex computed above */
2949:       VecSetValues(Vleaves,1,&k,&val,INSERT_VALUES);
2950:     }
2951:   }
2952:   VecRestoreArrayRead(Vleaves_seq,&varr_read);
2953:   VecAssemblyBegin(Vleaves);
2954:   VecAssemblyEnd(Vleaves);

2956:   /* (d) Set vltog for ghost vertices by copying local values of Vleaves */
2957:   VecGetArrayRead(Vleaves,&varr_read);
2958:   k = 0;
2959:   for (i=0; i<nroots; i++) {
2960:     DMNetworkIsGhostVertex(dm,i+network->vStart,&ghost);
2961:     if (!ghost) continue;
2962:     vltog[i] = (PetscInt)PetscRealPart(varr_read[2*k+1]); k++;
2963:   }
2964:   VecRestoreArrayRead(Vleaves,&varr_read);

2966:   VecDestroy(&Vleaves);
2967:   VecDestroy(&Vleaves_seq);
2968:   VecScatterDestroy(&ctx);
2969:   return(0);
2970: }

2972: PETSC_STATIC_INLINE PetscErrorCode DMISAddSize_private(DM_Network *network,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *nidx)
2973: {
2974:   PetscErrorCode           ierr;
2975:   PetscInt                 i,j,ncomps,nvar,key,offset=0;
2976:   DMNetworkComponentHeader header;

2979:   PetscSectionGetOffset(network->DataSection,p,&offset);
2980:   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
2981:   header = (DMNetworkComponentHeader)(network->componentdataarray+offset);

2983:   for (i=0; i<ncomps; i++) {
2984:     key  = header->key[i];
2985:     nvar = header->nvar[i];
2986:     for (j=0; j<numkeys; j++) {
2987:       if (key == keys[j]) {
2988:         if (!blocksize || blocksize[j] == -1) {
2989:           *nidx += nvar;
2990:         } else {
2991:           *nidx += nselectedvar[j]*nvar/blocksize[j];
2992:         }
2993:       }
2994:     }
2995:   }
2996:   return(0);
2997: }

2999: PETSC_STATIC_INLINE PetscErrorCode DMISComputeIdx_private(DM dm,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],PetscInt *ii,PetscInt *idx)
3000: {
3001:   PetscErrorCode           ierr;
3002:   PetscInt                 i,j,ncomps,nvar,key,offsetg,k,k1,offset=0;
3003:   DM_Network               *network = (DM_Network*)dm->data;
3004:   DMNetworkComponentHeader header;

3007:   PetscSectionGetOffset(network->DataSection,p,&offset);
3008:   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
3009:   header = (DMNetworkComponentHeader)(network->componentdataarray+offset);

3011:   for (i=0; i<ncomps; i++) {
3012:     key  = header->key[i];
3013:     nvar = header->nvar[i];
3014:     for (j=0; j<numkeys; j++) {
3015:       if (key != keys[j]) continue;

3017:       DMNetworkGetGlobalVecOffset(dm,p,i,&offsetg);
3018:       if (!blocksize || blocksize[j] == -1) {
3019:         for (k=0; k<nvar; k++) idx[(*ii)++] = offsetg + k;
3020:       } else {
3021:         for (k=0; k<nvar; k+=blocksize[j]) {
3022:           for (k1=0; k1<nselectedvar[j]; k1++) idx[(*ii)++] = offsetg + k + selectedvar[j][k1];
3023:         }
3024:       }
3025:     }
3026:   }
3027:   return(0);
3028: }

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

3033:   Collective

3035:   Input Parameters:
3036: + dm - DMNetwork object
3037: . numkeys - number of keys
3038: . keys - array of keys that define the components of the variables you wish to extract
3039: . blocksize - block size of the variables associated to the component
3040: . nselectedvar - number of variables in each block to select
3041: - selectedvar - the offset into the block of each variable in each block to select

3043:   Output Parameters:
3044: . is - the index set

3046:   Level: Advanced

3048:   Notes:
3049:     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.

3051: .seealso: DMNetworkCreate(), ISCreateGeneral(), DMNetworkCreateLocalIS()
3052: @*/
3053: PetscErrorCode DMNetworkCreateIS(DM dm,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],IS *is)
3054: {
3056:   MPI_Comm       comm;
3057:   DM_Network     *network = (DM_Network*)dm->data;
3058:   PetscInt       i,p,estart,eend,vstart,vend,nidx,*idx;
3059:   PetscBool      ghost;

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

3064:   /* Check input parameters */
3065:   for (i=0; i<numkeys; i++) {
3066:     if (!blocksize || blocksize[i] == -1) continue;
3067:     if (nselectedvar[i] > blocksize[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"number of selectedvariables %D cannot be larger than blocksize %D",nselectedvar[i],blocksize[i]);
3068:   }

3070:   DMNetworkGetEdgeRange(dm,&estart,&eend);
3071:   DMNetworkGetVertexRange(dm,&vstart,&vend);

3073:   /* Get local number of idx */
3074:   nidx = 0;
3075:   for (p=estart; p<eend; p++) {
3076:     DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx);
3077:   }
3078:   for (p=vstart; p<vend; p++) {
3079:     DMNetworkIsGhostVertex(dm,p,&ghost);
3080:     if (ghost) continue;
3081:     DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx);
3082:   }

3084:   /* Compute idx */
3085:   PetscMalloc1(nidx,&idx);
3086:   i = 0;
3087:   for (p=estart; p<eend; p++) {
3088:     DMISComputeIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx);
3089:   }
3090:   for (p=vstart; p<vend; p++) {
3091:     DMNetworkIsGhostVertex(dm,p,&ghost);
3092:     if (ghost) continue;
3093:     DMISComputeIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx);
3094:   }

3096:   /* Create is */
3097:   ISCreateGeneral(comm,nidx,idx,PETSC_COPY_VALUES,is);
3098:   PetscFree(idx);
3099:   return(0);
3100: }

3102: PETSC_STATIC_INLINE PetscErrorCode DMISComputeLocalIdx_private(DM dm,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],PetscInt *ii,PetscInt *idx)
3103: {
3104:   PetscErrorCode           ierr;
3105:   PetscInt                 i,j,ncomps,nvar,key,offsetl,k,k1,offset=0;
3106:   DM_Network               *network = (DM_Network*)dm->data;
3107:   DMNetworkComponentHeader header;

3110:   PetscSectionGetOffset(network->DataSection,p,&offset);
3111:   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
3112:   header = (DMNetworkComponentHeader)(network->componentdataarray+offset);

3114:   for (i=0; i<ncomps; i++) {
3115:     key  = header->key[i];
3116:     nvar = header->nvar[i];
3117:     for (j=0; j<numkeys; j++) {
3118:       if (key != keys[j]) continue;

3120:       DMNetworkGetLocalVecOffset(dm,p,i,&offsetl);
3121:       if (!blocksize || blocksize[j] == -1) {
3122:         for (k=0; k<nvar; k++) idx[(*ii)++] = offsetl + k;
3123:       } else {
3124:         for (k=0; k<nvar; k+=blocksize[j]) {
3125:           for (k1=0; k1<nselectedvar[j]; k1++) idx[(*ii)++] = offsetl + k + selectedvar[j][k1];
3126:         }
3127:       }
3128:     }
3129:   }
3130:   return(0);
3131: }

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

3136:   Not collective

3138:   Input Parameters:
3139: + dm - DMNetwork object
3140: . numkeys - number of keys
3141: . keys - array of keys that define the components of the variables you wish to extract
3142: . blocksize - block size of the variables associated to the component
3143: . nselectedvar - number of variables in each block to select
3144: - selectedvar - the offset into the block of each variable in each block to select

3146:   Output Parameters:
3147: . is - the index set

3149:   Level: Advanced

3151:   Notes:
3152:     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.

3154: .seealso: DMNetworkCreate(), DMNetworkCreateIS, ISCreateGeneral()
3155: @*/
3156: PetscErrorCode DMNetworkCreateLocalIS(DM dm,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],IS *is)
3157: {
3159:   DM_Network     *network = (DM_Network*)dm->data;
3160:   PetscInt       i,p,pstart,pend,nidx,*idx;

3163:   /* Check input parameters */
3164:   for (i=0; i<numkeys; i++) {
3165:     if (!blocksize || blocksize[i] == -1) continue;
3166:     if (nselectedvar[i] > blocksize[i]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"number of selectedvariables %D cannot be larger than blocksize %D",nselectedvar[i],blocksize[i]);
3167:   }

3169:   pstart = network->pStart;
3170:   pend   = network->pEnd;

3172:   /* Get local number of idx */
3173:   nidx = 0;
3174:   for (p=pstart; p<pend; p++) {
3175:     DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx);
3176:   }

3178:   /* Compute local idx */
3179:   PetscMalloc1(nidx,&idx);
3180:   i = 0;
3181:   for (p=pstart; p<pend; p++) {
3182:     DMISComputeLocalIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx);
3183:   }

3185:   /* Create is */
3186:   ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,is);
3187:   PetscFree(idx);
3188:   return(0);
3189: }