Actual source code: network.c

petsc-master 2019-12-13
Report Typos and Errors
  1:  #include <petsc/private/dmnetworkimpl.h>

  3: /*@
  4:   DMNetworkGetPlex - Gets the Plex DM associated with this network DM

  6:   Not collective

  8:   Input Parameters:
  9: + netdm - the dm object
 10: - plexmdm - the plex dm object

 12:   Level: Advanced

 14: .seealso: DMNetworkCreate()
 15: @*/
 16: PetscErrorCode DMNetworkGetPlex(DM netdm, DM *plexdm)
 17: {
 18:   DM_Network *network = (DM_Network*) netdm->data;

 21:   *plexdm = network->plex;
 22:   return(0);
 23: }

 25: /*@
 26:   DMNetworkGetSizes - Gets the the number of subnetworks and coupling subnetworks

 28:   Collective on dm

 30:   Input Parameters:
 31: + dm - the dm object
 32: . Nsubnet - global number of subnetworks
 33: - NsubnetCouple - global number of coupling subnetworks

 35:   Level: Intermediate

 37: .seealso: DMNetworkCreate()
 38: @*/
 39: PetscErrorCode DMNetworkGetSizes(DM netdm, PetscInt *Nsubnet, PetscInt *Ncsubnet)
 40: {
 41:   DM_Network *network = (DM_Network*) netdm->data;

 44:   *Nsubnet = network->nsubnet;
 45:   *Ncsubnet = network->ncsubnet;
 46:   return(0);
 47: }

 49: /*@
 50:   DMNetworkSetSizes - Sets the number of subnetworks, local and global vertices and edges for each subnetwork.

 52:   Collective on dm

 54:   Input Parameters:
 55: + dm - the dm object
 56: . Nsubnet - global number of subnetworks
 57: . nV - number of local vertices for each subnetwork
 58: . nE - number of local edges for each subnetwork
 59: . NsubnetCouple - global number of coupling subnetworks
 60: - nec - number of local edges for each coupling subnetwork

 62:    You cannot change the sizes once they have been set. nV, nE are arrays of length Nsubnet, and nec is array of length NsubnetCouple.

 64:    Level: intermediate

 66: .seealso: DMNetworkCreate()
 67: @*/
 68: PetscErrorCode DMNetworkSetSizes(DM dm,PetscInt Nsubnet,PetscInt nV[], PetscInt nE[],PetscInt NsubnetCouple,PetscInt nec[])
 69: {
 71:   DM_Network     *network = (DM_Network*) dm->data;
 72:   PetscInt       a[2],b[2],i;

 76:   if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet);
 77:   if (NsubnetCouple < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of coupling subnetworks %D cannot be less than 0",NsubnetCouple);

 81:   if (network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network");

 83:   if (!nV || !nE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Local vertex size or edge size must be provided");

 85:   network->nsubnet  = Nsubnet + NsubnetCouple;
 86:   network->ncsubnet = NsubnetCouple;
 87:   PetscCalloc1(Nsubnet+NsubnetCouple,&network->subnet);

 89:   /* ----------------------------------------------------------
 90:    p=v or e; P=V or E
 91:    subnet[0].pStart   = 0
 92:    subnet[i+1].pStart = subnet[i].pEnd = subnet[i].pStart + (nE[i] or NV[i])
 93:    ----------------------------------------------------------------------- */
 94:   for (i=0; i < Nsubnet; i++) {
 95:     /* Get global number of vertices and edges for subnet[i] */
 96:     a[0] = nV[i]; a[1] = nE[i]; /* local number of vertices (excluding ghost) and edges */
 97:     MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));
 98:     network->subnet[i].Nvtx = b[0];
 99:     network->subnet[i].Nedge = b[1];

101:     network->subnet[i].nvtx   = nV[i]; /* local nvtx, without ghost */

103:     /* global subnet[].vStart and vEnd, used by DMNetworkLayoutSetUp() */
104:     network->subnet[i].vStart = network->NVertices;
105:     network->subnet[i].vEnd   = network->subnet[i].vStart + network->subnet[i].Nvtx;

107:     network->nVertices += nV[i];
108:     network->NVertices += network->subnet[i].Nvtx;

110:     network->subnet[i].nedge  = nE[i];
111:     network->subnet[i].eStart = network->nEdges;
112:     network->subnet[i].eEnd   = network->subnet[i].eStart + nE[i];
113:     network->nEdges += nE[i];
114:     network->NEdges += network->subnet[i].Nedge;
115:   }

117:   /* coupling subnetwork */
118:   for (; i < Nsubnet+NsubnetCouple; i++) {
119:     /* Get global number of coupling edges for subnet[i] */
120:     MPIU_Allreduce(nec+(i-Nsubnet),&network->subnet[i].Nedge,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));

122:     network->subnet[i].nvtx   = 0; /* We design coupling subnetwork such that it does not have its own vertices */
123:     network->subnet[i].vStart = network->nVertices;
124:     network->subnet[i].vEnd   = network->subnet[i].vStart;

126:     network->subnet[i].nedge  = nec[i-Nsubnet];
127:     network->subnet[i].eStart = network->nEdges;
128:     network->subnet[i].eEnd = network->subnet[i].eStart + nec[i-Nsubnet];
129:     network->nEdges += nec[i-Nsubnet];
130:     network->NEdges += network->subnet[i].Nedge;
131:   }
132:   return(0);
133: }

135: /*@
136:   DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network

138:   Logically collective on dm

140:   Input Parameters:
141: + dm - the dm object
142: . edgelist - list of edges for each subnetwork
143: - edgelistCouple - list of edges for each coupling subnetwork

145:   Notes:
146:   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
147:   not be destroyed before the call to DMNetworkLayoutSetUp

149:   Level: intermediate

151:   Example usage:
152:   Consider the following 2 separate networks and a coupling network:

154: .vb
155:  network 0: v0 -> v1 -> v2 -> v3
156:  network 1: v1 -> v2 -> v0
157:  coupling network: network 1: v2 -> network 0: v0
158: .ve

160:  The resulting input
161:    edgelist[0] = [0 1 | 1 2 | 2 3];
162:    edgelist[1] = [1 2 | 2 0]
163:    edgelistCouple[0] = [(network)1 (v)2 (network)0 (v)0].

165: .seealso: DMNetworkCreate, DMNetworkSetSizes
166: @*/
167: PetscErrorCode DMNetworkSetEdgeList(DM dm,PetscInt *edgelist[],PetscInt *edgelistCouple[])
168: {
169:   DM_Network *network = (DM_Network*) dm->data;
170:   PetscInt   i;

173:   for (i=0; i < (network->nsubnet-network->ncsubnet); i++) network->subnet[i].edgelist = edgelist[i];
174:   if (network->ncsubnet) {
175:     PetscInt j = 0;
176:     if (!edgelistCouple) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Must provide edgelist_couple");
177:     while (i < network->nsubnet) network->subnet[i++].edgelist = edgelistCouple[j++];
178:   }
179:   return(0);
180: }

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

185:   Collective on dm

187:   Input Parameters
188: . DM - the dmnetwork object

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

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

196:   Level: intermediate

198: .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList
199: @*/
200: PetscErrorCode DMNetworkLayoutSetUp(DM dm)
201: {
203:   DM_Network     *network = (DM_Network*)dm->data;
204:   PetscInt       numCorners=2,spacedim=2,dim = 1; /* One dimensional network */
205:   PetscReal      *vertexcoords=NULL;
206:   PetscInt       i,j,ctr,nsubnet,*eowners,np,*edges,*subnetvtx,vStart;
207:   PetscInt       k,netid,vid, *vidxlTog,*edgelist_couple=NULL;
208:   const PetscInt *cone;
209:   MPI_Comm       comm;
210:   PetscMPIInt    size,rank;

213:   PetscObjectGetComm((PetscObject)dm,&comm);
214:   MPI_Comm_rank(comm,&rank);
215:   MPI_Comm_size(comm,&size);

217:   /* Create the local edgelist for the network by concatenating local input edgelists of the subnetworks */
218:   PetscCalloc2(numCorners*network->nVertices,&vertexcoords,2*network->nEdges,&edges);
219:   nsubnet = network->nsubnet - network->ncsubnet;
220:   ctr = 0;
221:   for (i=0; i < nsubnet; i++) {
222:     for (j = 0; j < network->subnet[i].nedge; j++) {
223:       edges[2*ctr]   = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
224:       edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
225:       ctr++;
226:     }
227:   }

229:   /* Append local coupling edgelists of the subnetworks */
230:   i       = nsubnet; /* netid of coupling subnet */
231:   nsubnet = network->nsubnet;
232:   while (i < nsubnet) {
233:     edgelist_couple = network->subnet[i].edgelist;

235:     k = 0;
236:     for (j = 0; j < network->subnet[i].nedge; j++) {
237:       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
238:       edges[2*ctr] = network->subnet[netid].vStart + vid; k += 2;

240:       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
241:       edges[2*ctr+1] = network->subnet[netid].vStart + vid; k+=2;
242:       ctr++;
243:     }
244:     i++;
245:   }
246:   /*
247:   if (rank == 0) {
248:     PetscPrintf(PETSC_COMM_SELF,"[%d] edgelist:\n",rank);
249:     for(i=0; i < network->nEdges; i++) {
250:       PetscPrintf(PETSC_COMM_SELF,"[%D %D]",edges[2*i],edges[2*i+1]);
251:       printf("\n");
252:     }
253:   }
254:    */

256:   /* Create network->plex */
257: #if defined(PETSC_USE_64BIT_INDICES)
258:   {
259:     int *edges64;
260:     np = network->nEdges*numCorners;
261:     PetscMalloc1(np,&edges64);
262:     for (i=0; i<np; i++) edges64[i] = (int)edges[i];

264:     if (size == 1) {
265:       DMPlexCreateFromCellList(comm,dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,(const int*)edges64,spacedim,(const double*)vertexcoords,&network->plex);
266:     } else {
267:       DMPlexCreateFromCellListParallel(comm,dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,(const int*)edges64,spacedim,(const PetscReal*)vertexcoords,NULL,&network->plex);
268:     }
269:     PetscFree(edges64);
270:   }
271: #else
272:   if (size == 1) {
273:     DMPlexCreateFromCellList(comm,dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,(const int*)edges,spacedim,(const double*)vertexcoords,&network->plex);
274:   } else {
275:     DMPlexCreateFromCellListParallel(comm,dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,(const int*)edges,spacedim,(const PetscReal*)vertexcoords,NULL,&network->plex);
276:   }
277: #endif

279:   DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);
280:   DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);
281:   DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);
282:   vStart = network->vStart;

284:   PetscSectionCreate(comm,&network->DataSection);
285:   PetscSectionCreate(comm,&network->DofSection);
286:   PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);
287:   PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);

289:   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
290:   np = network->pEnd - network->pStart;
291:   PetscCalloc2(np,&network->header,np,&network->cvalue);

293:   /* Create vidxlTog: maps local vertex index to global index */
294:   np = network->vEnd - vStart;
295:   PetscMalloc2(np,&vidxlTog,size+1,&eowners);
296:   ctr = 0;
297:   for (i=network->eStart; i<network->eEnd; i++) {
298:     DMNetworkGetConnectedVertices(dm,i,&cone);
299:     vidxlTog[cone[0] - vStart] = edges[2*ctr];
300:     vidxlTog[cone[1] - vStart] = edges[2*ctr+1];
301:     ctr++;
302:   }
303:   PetscFree2(vertexcoords,edges);

305:   /* Create vertices and edges array for the subnetworks */
306:   for (j=0; j < network->nsubnet; j++) {
307:     PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);

309:     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
310:        These get updated when the vertices and edges are added. */
311:     network->subnet[j].nvtx  = 0;
312:     network->subnet[j].nedge = 0;
313:   }
314:   PetscCalloc1(np,&network->subnetvtx);


317:   /* Get edge ownership */
318:   np = network->eEnd - network->eStart;
319:   MPI_Allgather(&np,1,MPIU_INT,eowners+1,1,MPIU_INT,comm);
320:   eowners[0] = 0;
321:   for (i=2; i<=size; i++) eowners[i] += eowners[i-1];

323:   i = 0; j = 0;
324:   while (i < np) { /* local edges, including coupling edges */
325:     network->header[i].index = i + eowners[rank];   /* Global edge index */

327:     if (j < network->nsubnet && i < network->subnet[j].eEnd) {
328:       network->header[i].subnetid = j; /* Subnetwork id */
329:       network->subnet[j].edges[network->subnet[j].nedge++] = i;

331:       network->header[i].ndata = 0;
332:       PetscSectionAddDof(network->DataSection,i,network->dataheadersize);
333:       network->header[i].offset[0] = 0;
334:       i++;
335:     }
336:     if (i >= network->subnet[j].eEnd) j++;
337:   }

339:   /* Count network->subnet[*].nvtx */
340:   for (i=vStart; i<network->vEnd; i++) { /* local vertices, including ghosts */
341:     k = vidxlTog[i-vStart];
342:     for (j=0; j < network->nsubnet; j++) {
343:       if (network->subnet[j].vStart <= k && k < network->subnet[j].vEnd) {
344:         network->subnet[j].nvtx++;
345:         break;
346:       }
347:     }
348:   }

350:   /* Set network->subnet[*].vertices on array network->subnetvtx */
351:   subnetvtx = network->subnetvtx;
352:   for (j=0; j<network->nsubnet; j++) {
353:     network->subnet[j].vertices = subnetvtx;
354:     subnetvtx                  += network->subnet[j].nvtx;
355:     network->subnet[j].nvtx = 0;
356:   }

358:   /* Set vertex array for the subnetworks */
359:   for (i=vStart; i<network->vEnd; i++) { /* local vertices, including ghosts */
360:     network->header[i].index = vidxlTog[i-vStart]; /*  Global vertex index */

362:     k = vidxlTog[i-vStart];
363:     for (j=0; j < network->nsubnet; j++) {
364:       if (network->subnet[j].vStart <= k && k < network->subnet[j].vEnd) {
365:         network->header[i].subnetid = j;
366:         network->subnet[j].vertices[network->subnet[j].nvtx++] = i;
367:         break;
368:       }
369:     }

371:     network->header[i].ndata = 0;
372:     PetscSectionAddDof(network->DataSection,i,network->dataheadersize);
373:     network->header[i].offset[0] = 0;
374:   }

376:   PetscFree2(vidxlTog,eowners);
377:   return(0);
378: }

380: /*@C
381:   DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork

383:   Input Parameters
384: + dm - the DM object
385: - id   - the ID (integer) of the subnetwork

387:   Output Parameters
388: + nv    - number of vertices (local)
389: . ne    - number of edges (local)
390: . vtx   - local vertices for this subnetwork
391: - edge  - local edges for this subnetwork

393:   Notes:
394:   Cannot call this routine before DMNetworkLayoutSetup()

396:   Level: intermediate

398: .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
399: @*/
400: PetscErrorCode DMNetworkGetSubnetworkInfo(DM dm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge)
401: {
402:   DM_Network *network = (DM_Network*)dm->data;

405:   if (id >= network->nsubnet) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Subnet ID %D exceeds the num of subnets %D",id,network->nsubnet);
406:   *nv   = network->subnet[id].nvtx;
407:   *ne   = network->subnet[id].nedge;
408:   *vtx  = network->subnet[id].vertices;
409:   *edge = network->subnet[id].edges;
410:   return(0);
411: }

413: /*@C
414:   DMNetworkGetSubnetworkCoupleInfo - Returns the info for the coupling subnetwork

416:   Input Parameters
417: + dm - the DM object
418: - id   - the ID (integer) of the coupling subnetwork

420:   Output Parameters
421: + ne - number of edges (local)
422: - edge  - local edges for this coupling subnetwork

424:   Notes:
425:   Cannot call this routine before DMNetworkLayoutSetup()

427:   Level: intermediate

429: .seealso: DMNetworkGetSubnetworkInfo, DMNetworkLayoutSetUp, DMNetworkCreate
430: @*/
431: PetscErrorCode DMNetworkGetSubnetworkCoupleInfo(DM dm,PetscInt id,PetscInt *ne,const PetscInt **edge)
432: {
433:   DM_Network *net = (DM_Network*)dm->data;
434:   PetscInt   id1;

437:   if (net->ncsubnet) {
438:     if (id >= net->ncsubnet) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Subnet ID %D exceeds the num of coupling subnets %D",id,net->ncsubnet);

440:     id1   = id + net->nsubnet - net->ncsubnet;
441:     *ne   = net->subnet[id1].nedge;
442:     *edge = net->subnet[id1].edges;
443:   } else {
444:     *ne   = 0;
445:     *edge = NULL;
446:   }
447:   return(0);
448: }

450: /*@C
451:   DMNetworkRegisterComponent - Registers the network component

453:   Logically collective on dm

455:   Input Parameters
456: + dm   - the network object
457: . name - the component name
458: - size - the storage size in bytes for this component data

460:    Output Parameters
461: .   key - an integer key that defines the component

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

466:    Level: intermediate

468: .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
469: @*/
470: PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,size_t size,PetscInt *key)
471: {
472:   PetscErrorCode        ierr;
473:   DM_Network            *network = (DM_Network*) dm->data;
474:   DMNetworkComponent    *component=&network->component[network->ncomponent];
475:   PetscBool             flg=PETSC_FALSE;
476:   PetscInt              i;

479:   for (i=0; i < network->ncomponent; i++) {
480:     PetscStrcmp(component->name,name,&flg);
481:     if (flg) {
482:       *key = i;
483:       return(0);
484:     }
485:   }
486:   if(network->ncomponent == MAX_COMPONENTS) {
487:     SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components registered exceeds the max %D",MAX_COMPONENTS);
488:   }

490:   PetscStrcpy(component->name,name);
491:   component->size = size/sizeof(DMNetworkComponentGenericDataType);
492:   *key = network->ncomponent;
493:   network->ncomponent++;
494:   return(0);
495: }

497: /*@
498:   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.

500:   Not Collective

502:   Input Parameters:
503: . dm - The DMNetwork object

505:   Output Paramters:
506: + vStart - The first vertex point
507: - vEnd   - One beyond the last vertex point

509:   Level: intermediate

511: .seealso: DMNetworkGetEdgeRange
512: @*/
513: PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
514: {
515:   DM_Network     *network = (DM_Network*)dm->data;

518:   if (vStart) *vStart = network->vStart;
519:   if (vEnd) *vEnd = network->vEnd;
520:   return(0);
521: }

523: /*@
524:   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.

526:   Not Collective

528:   Input Parameters:
529: . dm - The DMNetwork object

531:   Output Paramters:
532: + eStart - The first edge point
533: - eEnd   - One beyond the last edge point

535:   Level: intermediate

537: .seealso: DMNetworkGetVertexRange
538: @*/
539: PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
540: {
541:   DM_Network     *network = (DM_Network*)dm->data;

544:   if (eStart) *eStart = network->eStart;
545:   if (eEnd) *eEnd = network->eEnd;
546:   return(0);
547: }

549: /*@
550:   DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge.

552:   Not Collective

554:   Input Parameters:
555: + dm - DMNetwork object
556: - p  - edge point

558:   Output Paramters:
559: . index - user global numbering for the edge

561:   Level: intermediate

563: .seealso: DMNetworkGetGlobalVertexIndex
564: @*/
565: PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
566: {
567:   PetscErrorCode    ierr;
568:   DM_Network        *network = (DM_Network*)dm->data;
569:   PetscInt          offsetp;
570:   DMNetworkComponentHeader header;

573:   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first");
574:   PetscSectionGetOffset(network->DataSection,p,&offsetp);
575:   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
576:   *index = header->index;
577:   return(0);
578: }

580: /*@
581:   DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex.

583:   Not Collective

585:   Input Parameters:
586: + dm - DMNetwork object
587: - p  - vertex point

589:   Output Paramters:
590: . index - user global numbering for the vertex

592:   Level: intermediate

594: .seealso: DMNetworkGetGlobalEdgeIndex
595: @*/
596: PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
597: {
598:   PetscErrorCode    ierr;
599:   DM_Network        *network = (DM_Network*)dm->data;
600:   PetscInt          offsetp;
601:   DMNetworkComponentHeader header;

604:   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first");
605:   PetscSectionGetOffset(network->DataSection,p,&offsetp);
606:   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
607:   *index = header->index;
608:   return(0);
609: }

611: /*
612:   DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the
613:                                     component value from the component data array

615:   Not Collective

617:   Input Parameters:
618: + dm      - The DMNetwork object
619: . p       - vertex/edge point
620: - compnum - component number

622:   Output Parameters:
623: + compkey - the key obtained when registering the component
624: - offset  - offset into the component data array associated with the vertex/edge point

626:   Notes:
627:   Typical usage:

629:   DMNetworkGetComponentDataArray(dm, &arr);
630:   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
631:   Loop over vertices or edges
632:     DMNetworkGetNumComponents(dm,v,&numcomps);
633:     Loop over numcomps
634:       DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset);
635:       compdata = (UserCompDataType)(arr+offset);

637:   Level: intermediate

639: .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
640: */
641: PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
642: {
643:   PetscErrorCode           ierr;
644:   PetscInt                 offsetp;
645:   DMNetworkComponentHeader header;
646:   DM_Network               *network = (DM_Network*)dm->data;

649:   PetscSectionGetOffset(network->DataSection,p,&offsetp);
650:   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
651:   if (compkey) *compkey = header->key[compnum];
652:   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
653:   return(0);
654: }

656: /*@
657:   DMNetworkGetComponent - Returns the network component and its key

659:   Not Collective

661:   Input Parameters
662: + dm - DMNetwork object
663: . p  - edge or vertex point
664: - compnum - component number

666:   Output Parameters:
667: + compkey - the key set for this computing during registration
668: - component - the component data

670:   Notes:
671:   Typical usage:

673:   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
674:   Loop over vertices or edges
675:     DMNetworkGetNumComponents(dm,v,&numcomps);
676:     Loop over numcomps
677:       DMNetworkGetComponent(dm,v,compnum,&key,&component);

679:   Level: intermediate

681: .seealso: DMNetworkGetNumComponents, DMNetworkGetVariableOffset
682: @*/
683: PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *key, void **component)
684: {
686:   DM_Network     *network = (DM_Network*)dm->data;
687:   PetscInt       offsetd = 0;

690:   DMNetworkGetComponentKeyOffset(dm,p,compnum,key,&offsetd);
691:   *component = network->componentdataarray+offsetd;
692:   return(0);
693: }

695: /*@
696:   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)

698:   Not Collective

700:   Input Parameters:
701: + dm           - The DMNetwork object
702: . p            - vertex/edge point
703: . componentkey - component key returned while registering the component
704: - compvalue    - pointer to the data structure for the component

706:   Level: intermediate

708: .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
709: @*/
710: PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
711: {
712:   DM_Network               *network = (DM_Network*)dm->data;
713:   DMNetworkComponent       *component = &network->component[componentkey];
714:   DMNetworkComponentHeader header = &network->header[p];
715:   DMNetworkComponentValue  cvalue = &network->cvalue[p];
716:   PetscErrorCode           ierr;

719:   if (header->ndata == MAX_DATA_AT_POINT) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components at a point exceeds the max %D",MAX_DATA_AT_POINT);

721:   header->size[header->ndata] = component->size;
722:   PetscSectionAddDof(network->DataSection,p,component->size);
723:   header->key[header->ndata] = componentkey;
724:   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];

726:   cvalue->data[header->ndata] = (void*)compvalue;
727:   header->ndata++;
728:   return(0);
729: }

731: /*@
732:   DMNetworkGetNumComponents - Get the number of components at a vertex/edge

734:   Not Collective

736:   Input Parameters:
737: + dm - The DMNetwork object
738: - p  - vertex/edge point

740:   Output Parameters:
741: . numcomponents - Number of components at the vertex/edge

743:   Level: intermediate

745: .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
746: @*/
747: PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
748: {
750:   PetscInt       offset;
751:   DM_Network     *network = (DM_Network*)dm->data;

754:   PetscSectionGetOffset(network->DataSection,p,&offset);
755:   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
756:   return(0);
757: }

759: /*@
760:   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.

762:   Not Collective

764:   Input Parameters:
765: + dm     - The DMNetwork object
766: - p      - the edge/vertex point

768:   Output Parameters:
769: . offset - the offset

771:   Level: intermediate

773: .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
774: @*/
775: PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
776: {
778:   DM_Network     *network = (DM_Network*)dm->data;

781:   PetscSectionGetOffset(network->plex->localSection,p,offset);
782:   return(0);
783: }

785: /*@
786:   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.

788:   Not Collective

790:   Input Parameters:
791: + dm      - The DMNetwork object
792: - p       - the edge/vertex point

794:   Output Parameters:
795: . offsetg - the offset

797:   Level: intermediate

799: .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
800: @*/
801: PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
802: {
804:   DM_Network     *network = (DM_Network*)dm->data;

807:   PetscSectionGetOffset(network->plex->globalSection,p,offsetg);
808:   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */
809:   return(0);
810: }

812: /*@
813:   DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector.

815:   Not Collective

817:   Input Parameters:
818: + dm     - The DMNetwork object
819: - p      - the edge point

821:   Output Parameters:
822: . offset - the offset

824:   Level: intermediate

826: .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
827: @*/
828: PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
829: {
831:   DM_Network     *network = (DM_Network*)dm->data;


835:   PetscSectionGetOffset(network->edge.DofSection,p,offset);
836:   return(0);
837: }

839: /*@
840:   DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector.

842:   Not Collective

844:   Input Parameters:
845: + dm     - The DMNetwork object
846: - p      - the vertex point

848:   Output Parameters:
849: . offset - the offset

851:   Level: intermediate

853: .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
854: @*/
855: PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
856: {
858:   DM_Network     *network = (DM_Network*)dm->data;


862:   p -= network->vStart;

864:   PetscSectionGetOffset(network->vertex.DofSection,p,offset);
865:   return(0);
866: }
867: /*@
868:   DMNetworkAddNumVariables - Add number of variables associated with a given point.

870:   Not Collective

872:   Input Parameters:
873: + dm   - The DMNetworkObject
874: . p    - the vertex/edge point
875: - nvar - number of additional variables

877:   Level: intermediate

879: .seealso: DMNetworkSetNumVariables
880: @*/
881: PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
882: {
884:   DM_Network     *network = (DM_Network*)dm->data;

887:   PetscSectionAddDof(network->DofSection,p,nvar);
888:   return(0);
889: }

891: /*@
892:   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.

894:   Not Collective

896:   Input Parameters:
897: + dm   - The DMNetworkObject
898: - p    - the vertex/edge point

900:   Output Parameters:
901: . nvar - number of variables

903:   Level: intermediate

905: .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
906: @*/
907: PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
908: {
910:   DM_Network     *network = (DM_Network*)dm->data;

913:   PetscSectionGetDof(network->DofSection,p,nvar);
914:   return(0);
915: }

917: /*@
918:   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.

920:   Not Collective

922:   Input Parameters:
923: + dm   - The DMNetworkObject
924: . p    - the vertex/edge point
925: - nvar - number of variables

927:   Level: intermediate

929: .seealso: DMNetworkAddNumVariables
930: @*/
931: PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
932: {
934:   DM_Network     *network = (DM_Network*)dm->data;

937:   PetscSectionSetDof(network->DofSection,p,nvar);
938:   return(0);
939: }

941: /* Sets up the array that holds the data for all components and its associated section. This
942:    function is called during DMSetUp() */
943: PetscErrorCode DMNetworkComponentSetUp(DM dm)
944: {
945:   PetscErrorCode           ierr;
946:   DM_Network               *network = (DM_Network*)dm->data;
947:   PetscInt                 arr_size,p,offset,offsetp,ncomp,i;
948:   DMNetworkComponentHeader header;
949:   DMNetworkComponentValue  cvalue;
950:   DMNetworkComponentGenericDataType *componentdataarray;

953:   PetscSectionSetUp(network->DataSection);
954:   PetscSectionGetStorageSize(network->DataSection,&arr_size);
955:   PetscMalloc1(arr_size,&network->componentdataarray);
956:   componentdataarray = network->componentdataarray;
957:   for (p = network->pStart; p < network->pEnd; p++) {
958:     PetscSectionGetOffset(network->DataSection,p,&offsetp);
959:     /* Copy header */
960:     header = &network->header[p];
961:     PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));
962:     /* Copy data */
963:     cvalue = &network->cvalue[p];
964:     ncomp = header->ndata;
965:     for (i = 0; i < ncomp; i++) {
966:       offset = offsetp + network->dataheadersize + header->offset[i];
967:       PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));
968:     }
969:   }
970:   return(0);
971: }

973: /* Sets up the section for dofs. This routine is called during DMSetUp() */
974: PetscErrorCode DMNetworkVariablesSetUp(DM dm)
975: {
977:   DM_Network     *network = (DM_Network*)dm->data;

980:   PetscSectionSetUp(network->DofSection);
981:   return(0);
982: }

984: /*@C
985:   DMNetworkGetComponentDataArray - Returns the component data array

987:   Not Collective

989:   Input Parameters:
990: . dm - The DMNetwork Object

992:   Output Parameters:
993: . componentdataarray - array that holds data for all components

995:   Level: intermediate

997: .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents
998: @*/
999: PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
1000: {
1001:   DM_Network     *network = (DM_Network*)dm->data;

1004:   *componentdataarray = network->componentdataarray;
1005:   return(0);
1006: }

1008: /* Get a subsection from a range of points */
1009: PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection)
1010: {
1012:   PetscInt       i, nvar;

1015:   PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);
1016:   PetscSectionSetChart(*subsection, 0, pend - pstart);
1017:   for (i = pstart; i < pend; i++) {
1018:     PetscSectionGetDof(master,i,&nvar);
1019:     PetscSectionSetDof(*subsection, i - pstart, nvar);
1020:   }

1022:   PetscSectionSetUp(*subsection);
1023:   return(0);
1024: }

1026: /* Create a submap of points with a GlobalToLocal structure */
1027: PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
1028: {
1030:   PetscInt       i, *subpoints;

1033:   /* Create index sets to map from "points" to "subpoints" */
1034:   PetscMalloc1(pend - pstart, &subpoints);
1035:   for (i = pstart; i < pend; i++) {
1036:     subpoints[i - pstart] = i;
1037:   }
1038:   ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);
1039:   PetscFree(subpoints);
1040:   return(0);
1041: }

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

1046:   Collective

1048:   Input Parameters:
1049: . dm   - The DMNetworkObject

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

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

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

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

1059:   Level: intermediate

1061: @*/
1062: PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
1063: {
1065:   MPI_Comm       comm;
1066:   PetscMPIInt    rank, size;
1067:   DM_Network     *network = (DM_Network*)dm->data;

1070:   PetscObjectGetComm((PetscObject)dm,&comm);
1071:   MPI_Comm_rank(comm, &rank);
1072:   MPI_Comm_size(comm, &size);

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

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

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

1085:     PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);
1086:     PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);
1087:     PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);
1088:   } else {
1089:     /* create structures for vertex */
1090:     PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);
1091:     /* create structures for edge */
1092:     PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);
1093:   }

1095:   /* Add viewers */
1096:   PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");
1097:   PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");
1098:   PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");
1099:   PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");
1100:   return(0);
1101: }

1103: /*@
1104:   DMNetworkDistribute - Distributes the network and moves associated component data.

1106:   Collective

1108:   Input Parameter:
1109: + DM - the DMNetwork object
1110: - overlap - The overlap of partitions, 0 is the default

1112:   Notes:
1113:   Distributes the network with <overlap>-overlapping partitioning of the edges.

1115:   Level: intermediate

1117: .seealso: DMNetworkCreate
1118: @*/
1119: PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
1120: {
1121:   MPI_Comm       comm;
1123:   PetscMPIInt    size;
1124:   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
1125:   DM_Network     *newDMnetwork;
1126:   PetscSF        pointsf=NULL;
1127:   DM             newDM;
1128:   PetscInt       j,e,v,offset,*subnetvtx;
1129:   PetscPartitioner         part;
1130:   DMNetworkComponentHeader header;

1133:   PetscObjectGetComm((PetscObject)*dm,&comm);
1134:   MPI_Comm_size(comm, &size);
1135:   if (size == 1) return(0);

1137:   DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);
1138:   newDMnetwork = (DM_Network*)newDM->data;
1139:   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);

1141:   /* Enable runtime options for petscpartitioner */
1142:   DMPlexGetPartitioner(oldDMnetwork->plex,&part);
1143:   PetscPartitionerSetFromOptions(part);

1145:   /* Distribute plex dm and dof section */
1146:   DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);

1148:   /* Distribute dof section */
1149:   PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);
1150:   PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);
1151:   PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);

1153:   /* Distribute data and associated section */
1154:   DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);

1156:   PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);
1157:   DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);
1158:   DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);
1159:   newDMnetwork->nEdges    = newDMnetwork->eEnd - newDMnetwork->eStart;
1160:   newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
1161:   newDMnetwork->NVertices = oldDMnetwork->NVertices;
1162:   newDMnetwork->NEdges    = oldDMnetwork->NEdges;

1164:   /* Set Dof section as the section for dm */
1165:   DMSetLocalSection(newDMnetwork->plex,newDMnetwork->DofSection);
1166:   DMGetGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);

1168:   /* Set up subnetwork info in the newDM */
1169:   newDMnetwork->nsubnet  = oldDMnetwork->nsubnet;
1170:   newDMnetwork->ncsubnet = oldDMnetwork->ncsubnet;
1171:   PetscCalloc1(newDMnetwork->nsubnet,&newDMnetwork->subnet);
1172:   /* Copy over the global number of vertices and edges in each subnetwork. Note that these are already
1173:      calculated in DMNetworkLayoutSetUp()
1174:   */
1175:   for(j=0; j < newDMnetwork->nsubnet; j++) {
1176:     newDMnetwork->subnet[j].Nvtx  = oldDMnetwork->subnet[j].Nvtx;
1177:     newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge;
1178:   }

1180:   for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1181:     PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);
1182:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1183:     newDMnetwork->subnet[header->subnetid].nedge++;
1184:   }

1186:   for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1187:     PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);
1188:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1189:     newDMnetwork->subnet[header->subnetid].nvtx++;
1190:   }

1192:   /* Now create the vertices and edge arrays for the subnetworks */
1193:   PetscCalloc1(newDMnetwork->vEnd-newDMnetwork->vStart,&newDMnetwork->subnetvtx);
1194:   subnetvtx = newDMnetwork->subnetvtx;

1196:   for (j=0; j<newDMnetwork->nsubnet; j++) {
1197:     PetscCalloc1(newDMnetwork->subnet[j].nedge,&newDMnetwork->subnet[j].edges);
1198:     newDMnetwork->subnet[j].vertices = subnetvtx;
1199:     subnetvtx                       += newDMnetwork->subnet[j].nvtx;

1201:     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
1202:        These get updated when the vertices and edges are added. */
1203:     newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0;
1204:   }

1206:   /* Set the vertices and edges in each subnetwork */
1207:   for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1208:     PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);
1209:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1210:     newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++] = e;
1211:   }

1213:   for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1214:     PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);
1215:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1216:     newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++] = v;
1217:   }

1219:   newDM->setupcalled = (*dm)->setupcalled;
1220:   newDMnetwork->distributecalled = PETSC_TRUE;

1222:   /* Destroy point SF */
1223:   PetscSFDestroy(&pointsf);

1225:   DMDestroy(dm);
1226:   *dm  = newDM;
1227:   return(0);
1228: }

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

1233:   Input Parameters:
1234: + masterSF - the original SF structure
1235: - map      - a ISLocalToGlobal mapping that contains the subset of points

1237:   Output Parameters:
1238: . subSF    - a subset of the masterSF for the desired subset.
1239: */

1241: PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) {

1243:   PetscErrorCode        ierr;
1244:   PetscInt              nroots, nleaves, *ilocal_sub;
1245:   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
1246:   PetscInt              *local_points, *remote_points;
1247:   PetscSFNode           *iremote_sub;
1248:   const PetscInt        *ilocal;
1249:   const PetscSFNode     *iremote;

1252:   PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);

1254:   /* Look for leaves that pertain to the subset of points. Get the local ordering */
1255:   PetscMalloc1(nleaves,&ilocal_map);
1256:   ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);
1257:   for (i = 0; i < nleaves; i++) {
1258:     if (ilocal_map[i] != -1) nleaves_sub += 1;
1259:   }
1260:   /* Re-number ilocal with subset numbering. Need information from roots */
1261:   PetscMalloc2(nroots,&local_points,nroots,&remote_points);
1262:   for (i = 0; i < nroots; i++) local_points[i] = i;
1263:   ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);
1264:   PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);
1265:   PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);
1266:   /* Fill up graph using local (that is, local to the subset) numbering. */
1267:   PetscMalloc1(nleaves_sub,&ilocal_sub);
1268:   PetscMalloc1(nleaves_sub,&iremote_sub);
1269:   nleaves_sub = 0;
1270:   for (i = 0; i < nleaves; i++) {
1271:     if (ilocal_map[i] != -1) {
1272:       ilocal_sub[nleaves_sub] = ilocal_map[i];
1273:       iremote_sub[nleaves_sub].rank = iremote[i].rank;
1274:       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
1275:       nleaves_sub += 1;
1276:     }
1277:   }
1278:   PetscFree2(local_points,remote_points);
1279:   ISLocalToGlobalMappingGetSize(map,&nroots_sub);

1281:   /* Create new subSF */
1282:   PetscSFCreate(PETSC_COMM_WORLD,subSF);
1283:   PetscSFSetFromOptions(*subSF);
1284:   PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);
1285:   PetscFree(ilocal_map);
1286:   PetscFree(iremote_sub);
1287:   return(0);
1288: }

1290: /*@C
1291:   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point

1293:   Not Collective

1295:   Input Parameters:
1296: + dm - The DMNetwork object
1297: - p  - the vertex point

1299:   Output Paramters:
1300: + nedges - number of edges connected to this vertex point
1301: - edges  - List of edge points

1303:   Level: intermediate

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

1309: .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices
1310: @*/
1311: PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
1312: {
1314:   DM_Network     *network = (DM_Network*)dm->data;

1317:   DMPlexGetSupportSize(network->plex,vertex,nedges);
1318:   DMPlexGetSupport(network->plex,vertex,edges);
1319:   return(0);
1320: }

1322: /*@C
1323:   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point

1325:   Not Collective

1327:   Input Parameters:
1328: + dm - The DMNetwork object
1329: - p  - the edge point

1331:   Output Paramters:
1332: . vertices  - vertices connected to this edge

1334:   Level: intermediate

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

1340: .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
1341: @*/
1342: PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
1343: {
1345:   DM_Network     *network = (DM_Network*)dm->data;

1348:   DMPlexGetCone(network->plex,edge,vertices);
1349:   return(0);
1350: }

1352: /*@
1353:   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex

1355:   Not Collective

1357:   Input Parameters:
1358: + dm - The DMNetwork object
1359: - p  - the vertex point

1361:   Output Parameter:
1362: . isghost - TRUE if the vertex is a ghost point

1364:   Level: intermediate

1366: .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange
1367: @*/
1368: PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
1369: {
1371:   DM_Network     *network = (DM_Network*)dm->data;
1372:   PetscInt       offsetg;
1373:   PetscSection   sectiong;

1376:   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first");
1377:   *isghost = PETSC_FALSE;
1378:   DMGetGlobalSection(network->plex,&sectiong);
1379:   PetscSectionGetOffset(sectiong,p,&offsetg);
1380:   if (offsetg < 0) *isghost = PETSC_TRUE;
1381:   return(0);
1382: }

1384: PetscErrorCode DMSetUp_Network(DM dm)
1385: {
1387:   DM_Network     *network=(DM_Network*)dm->data;

1390:   DMNetworkComponentSetUp(dm);
1391:   DMNetworkVariablesSetUp(dm);

1393:   DMSetLocalSection(network->plex,network->DofSection);
1394:   DMGetGlobalSection(network->plex,&network->GlobalDofSection);

1396:   dm->setupcalled = PETSC_TRUE;
1397:   DMViewFromOptions(dm,NULL,"-dm_view");
1398:   return(0);
1399: }

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

1405:     Collective

1407:     Input Parameters:
1408: +   dm - The DMNetwork object
1409: .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
1410: -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices

1412:     Level: intermediate

1414: @*/
1415: PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
1416: {
1417:   DM_Network     *network=(DM_Network*)dm->data;
1419:   PetscInt       nVertices = network->nVertices;

1422:   network->userEdgeJacobian   = eflg;
1423:   network->userVertexJacobian = vflg;

1425:   if (eflg && !network->Je) {
1426:     PetscCalloc1(3*network->nEdges,&network->Je);
1427:   }

1429:   if (vflg && !network->Jv && nVertices) {
1430:     PetscInt       i,*vptr,nedges,vStart=network->vStart;
1431:     PetscInt       nedges_total;
1432:     const PetscInt *edges;

1434:     /* count nvertex_total */
1435:     nedges_total = 0;
1436:     PetscMalloc1(nVertices+1,&vptr);

1438:     vptr[0] = 0;
1439:     for (i=0; i<nVertices; i++) {
1440:       DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);
1441:       nedges_total += nedges;
1442:       vptr[i+1] = vptr[i] + 2*nedges + 1;
1443:     }

1445:     PetscCalloc1(2*nedges_total+nVertices,&network->Jv);
1446:     network->Jvptr = vptr;
1447:   }
1448:   return(0);
1449: }

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

1454:     Not Collective

1456:     Input Parameters:
1457: +   dm - The DMNetwork object
1458: .   p  - the edge point
1459: -   J - array (size = 3) of Jacobian submatrices for this edge point:
1460:         J[0]: this edge
1461:         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()

1463:     Level: intermediate

1465: .seealso: DMNetworkVertexSetMatrix
1466: @*/
1467: PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
1468: {
1469:   DM_Network     *network=(DM_Network*)dm->data;

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

1474:   if (J) {
1475:     network->Je[3*p]   = J[0];
1476:     network->Je[3*p+1] = J[1];
1477:     network->Je[3*p+2] = J[2];
1478:   }
1479:   return(0);
1480: }

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

1485:     Not Collective

1487:     Input Parameters:
1488: +   dm - The DMNetwork object
1489: .   p  - the vertex point
1490: -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
1491:         J[0]:       this vertex
1492:         J[1+2*i]:   i-th supporting edge
1493:         J[1+2*i+1]: i-th connected vertex

1495:     Level: intermediate

1497: .seealso: DMNetworkEdgeSetMatrix
1498: @*/
1499: PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
1500: {
1502:   DM_Network     *network=(DM_Network*)dm->data;
1503:   PetscInt       i,*vptr,nedges,vStart=network->vStart;
1504:   const PetscInt *edges;

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

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

1513:     /* Set Jacobian for each supporting edge and connected vertex */
1514:     DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);
1515:     for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
1516:   }
1517:   return(0);
1518: }

1520: PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1521: {
1523:   PetscInt       j;
1524:   PetscScalar    val=(PetscScalar)ncols;

1527:   if (!ghost) {
1528:     for (j=0; j<nrows; j++) {
1529:       VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);
1530:     }
1531:   } else {
1532:     for (j=0; j<nrows; j++) {
1533:       VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);
1534:     }
1535:   }
1536:   return(0);
1537: }

1539: PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1540: {
1542:   PetscInt       j,ncols_u;
1543:   PetscScalar    val;

1546:   if (!ghost) {
1547:     for (j=0; j<nrows; j++) {
1548:       MatGetRow(Ju,j,&ncols_u,NULL,NULL);
1549:       val = (PetscScalar)ncols_u;
1550:       VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);
1551:       MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);
1552:     }
1553:   } else {
1554:     for (j=0; j<nrows; j++) {
1555:       MatGetRow(Ju,j,&ncols_u,NULL,NULL);
1556:       val = (PetscScalar)ncols_u;
1557:       VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);
1558:       MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);
1559:     }
1560:   }
1561:   return(0);
1562: }

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

1569:   if (Ju) {
1570:     MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);
1571:   } else {
1572:     MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);
1573:   }
1574:   return(0);
1575: }

1577: PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1578: {
1580:   PetscInt       j,*cols;
1581:   PetscScalar    *zeros;

1584:   PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);
1585:   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1586:   MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);
1587:   PetscFree2(cols,zeros);
1588:   return(0);
1589: }

1591: PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1592: {
1594:   PetscInt       j,M,N,row,col,ncols_u;
1595:   const PetscInt *cols;
1596:   PetscScalar    zero=0.0;

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

1602:   for (row=0; row<nrows; row++) {
1603:     MatGetRow(Ju,row,&ncols_u,&cols,NULL);
1604:     for (j=0; j<ncols_u; j++) {
1605:       col = cols[j] + cstart;
1606:       MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);
1607:     }
1608:     MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);
1609:   }
1610:   return(0);
1611: }

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

1618:   if (Ju) {
1619:     MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);
1620:   } else {
1621:     MatSetDenseblock_private(nrows,rows,ncols,cstart,J);
1622:   }
1623:   return(0);
1624: }

1626: /* 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.
1627: */
1628: PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
1629: {
1631:   PetscInt       i,size,dof;
1632:   PetscInt       *glob2loc;

1635:   PetscSectionGetStorageSize(localsec,&size);
1636:   PetscMalloc1(size,&glob2loc);

1638:   for (i = 0; i < size; i++) {
1639:     PetscSectionGetOffset(globalsec,i,&dof);
1640:     dof = (dof >= 0) ? dof : -(dof + 1);
1641:     glob2loc[i] = dof;
1642:   }

1644:   ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);
1645: #if 0
1646:   PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);
1647: #endif
1648:   return(0);
1649: }

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

1653: PetscErrorCode DMCreateMatrix_Network_Nest(DM dm,Mat *J)
1654: {
1656:   DM_Network     *network = (DM_Network*)dm->data;
1657:   PetscMPIInt    rank, size;
1658:   PetscInt       eDof,vDof;
1659:   Mat            j11,j12,j21,j22,bA[2][2];
1660:   MPI_Comm       comm;
1661:   ISLocalToGlobalMapping eISMap,vISMap;

1664:   PetscObjectGetComm((PetscObject)dm,&comm);
1665:   MPI_Comm_rank(comm,&rank);
1666:   MPI_Comm_size(comm,&size);

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

1671:   MatCreate(comm, &j11);
1672:   MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);
1673:   MatSetType(j11, MATMPIAIJ);

1675:   MatCreate(comm, &j12);
1676:   MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);
1677:   MatSetType(j12, MATMPIAIJ);

1679:   MatCreate(comm, &j21);
1680:   MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);
1681:   MatSetType(j21, MATMPIAIJ);

1683:   MatCreate(comm, &j22);
1684:   MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);
1685:   MatSetType(j22, MATMPIAIJ);

1687:   bA[0][0] = j11;
1688:   bA[0][1] = j12;
1689:   bA[1][0] = j21;
1690:   bA[1][1] = j22;

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

1695:   MatSetLocalToGlobalMapping(j11,eISMap,eISMap);
1696:   MatSetLocalToGlobalMapping(j12,eISMap,vISMap);
1697:   MatSetLocalToGlobalMapping(j21,vISMap,eISMap);
1698:   MatSetLocalToGlobalMapping(j22,vISMap,vISMap);

1700:   MatSetUp(j11);
1701:   MatSetUp(j12);
1702:   MatSetUp(j21);
1703:   MatSetUp(j22);

1705:   MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);
1706:   MatSetUp(*J);
1707:   MatNestSetVecType(*J,VECNEST);
1708:   MatDestroy(&j11);
1709:   MatDestroy(&j12);
1710:   MatDestroy(&j21);
1711:   MatDestroy(&j22);

1713:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
1714:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
1715:   MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);

1717:   /* Free structures */
1718:   ISLocalToGlobalMappingDestroy(&eISMap);
1719:   ISLocalToGlobalMappingDestroy(&vISMap);
1720:   return(0);
1721: }

1723: PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
1724: {
1726:   DM_Network     *network = (DM_Network*)dm->data;
1727:   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1728:   PetscInt       cstart,ncols,j,e,v;
1729:   PetscBool      ghost,ghost_vc,ghost2,isNest;
1730:   Mat            Juser;
1731:   PetscSection   sectionGlobal;
1732:   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
1733:   const PetscInt *edges,*cone;
1734:   MPI_Comm       comm;
1735:   MatType        mtype;
1736:   Vec            vd_nz,vo_nz;
1737:   PetscInt       *dnnz,*onnz;
1738:   PetscScalar    *vdnz,*vonz;

1741:   mtype = dm->mattype;
1742:   PetscStrcmp(mtype,MATNEST,&isNest);
1743:   if (isNest) {
1744:     DMCreateMatrix_Network_Nest(dm,J);
1745:     return(0);
1746:   }

1748:   if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1749:     /* user does not provide Jacobian blocks */
1750:     DMCreateMatrix_Plex(network->plex,J);
1751:     MatSetDM(*J,dm);
1752:     return(0);
1753:   }

1755:   MatCreate(PetscObjectComm((PetscObject)dm),J);
1756:   DMGetGlobalSection(network->plex,&sectionGlobal);
1757:   PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);
1758:   MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);

1760:   MatSetType(*J,MATAIJ);
1761:   MatSetFromOptions(*J);

1763:   /* (1) Set matrix preallocation */
1764:   /*------------------------------*/
1765:   PetscObjectGetComm((PetscObject)dm,&comm);
1766:   VecCreate(comm,&vd_nz);
1767:   VecSetSizes(vd_nz,localSize,PETSC_DECIDE);
1768:   VecSetFromOptions(vd_nz);
1769:   VecSet(vd_nz,0.0);
1770:   VecDuplicate(vd_nz,&vo_nz);

1772:   /* Set preallocation for edges */
1773:   /*-----------------------------*/
1774:   DMNetworkGetEdgeRange(dm,&eStart,&eEnd);

1776:   PetscMalloc1(localSize,&rows);
1777:   for (e=eStart; e<eEnd; e++) {
1778:     /* Get row indices */
1779:     DMNetworkGetVariableGlobalOffset(dm,e,&rstart);
1780:     DMNetworkGetNumVariables(dm,e,&nrows);
1781:     if (nrows) {
1782:       for (j=0; j<nrows; j++) rows[j] = j + rstart;

1784:       /* Set preallocation for conntected vertices */
1785:       DMNetworkGetConnectedVertices(dm,e,&cone);
1786:       for (v=0; v<2; v++) {
1787:         DMNetworkGetNumVariables(dm,cone[v],&ncols);

1789:         if (network->Je) {
1790:           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1791:         } else Juser = NULL;
1792:         DMNetworkIsGhostVertex(dm,cone[v],&ghost);
1793:         MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);
1794:       }

1796:       /* Set preallocation for edge self */
1797:       cstart = rstart;
1798:       if (network->Je) {
1799:         Juser = network->Je[3*e]; /* Jacobian(e,e) */
1800:       } else Juser = NULL;
1801:       MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);
1802:     }
1803:   }

1805:   /* Set preallocation for vertices */
1806:   /*--------------------------------*/
1807:   DMNetworkGetVertexRange(dm,&vStart,&vEnd);
1808:   if (vEnd - vStart) vptr = network->Jvptr;

1810:   for (v=vStart; v<vEnd; v++) {
1811:     /* Get row indices */
1812:     DMNetworkGetVariableGlobalOffset(dm,v,&rstart);
1813:     DMNetworkGetNumVariables(dm,v,&nrows);
1814:     if (!nrows) continue;

1816:     DMNetworkIsGhostVertex(dm,v,&ghost);
1817:     if (ghost) {
1818:       PetscMalloc1(nrows,&rows_v);
1819:     } else {
1820:       rows_v = rows;
1821:     }

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

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

1828:     for (e=0; e<nedges; e++) {
1829:       /* Supporting edges */
1830:       DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);
1831:       DMNetworkGetNumVariables(dm,edges[e],&ncols);

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

1838:       /* Connected vertices */
1839:       DMNetworkGetConnectedVertices(dm,edges[e],&cone);
1840:       vc = (v == cone[0]) ? cone[1]:cone[0];
1841:       DMNetworkIsGhostVertex(dm,vc,&ghost_vc);

1843:       DMNetworkGetNumVariables(dm,vc,&ncols);

1845:       if (network->Jv) {
1846:         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1847:       } else Juser = NULL;
1848:       if (ghost_vc||ghost) {
1849:         ghost2 = PETSC_TRUE;
1850:       } else {
1851:         ghost2 = PETSC_FALSE;
1852:       }
1853:       MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);
1854:     }

1856:     /* Set preallocation for vertex self */
1857:     DMNetworkIsGhostVertex(dm,v,&ghost);
1858:     if (!ghost) {
1859:       DMNetworkGetVariableGlobalOffset(dm,v,&cstart);
1860:       if (network->Jv) {
1861:         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1862:       } else Juser = NULL;
1863:       MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);
1864:     }
1865:     if (ghost) {
1866:       PetscFree(rows_v);
1867:     }
1868:   }

1870:   VecAssemblyBegin(vd_nz);
1871:   VecAssemblyBegin(vo_nz);

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

1875:   VecAssemblyEnd(vd_nz);
1876:   VecAssemblyEnd(vo_nz);

1878:   VecGetArray(vd_nz,&vdnz);
1879:   VecGetArray(vo_nz,&vonz);
1880:   for (j=0; j<localSize; j++) {
1881:     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1882:     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1883:   }
1884:   VecRestoreArray(vd_nz,&vdnz);
1885:   VecRestoreArray(vo_nz,&vonz);
1886:   VecDestroy(&vd_nz);
1887:   VecDestroy(&vo_nz);

1889:   MatSeqAIJSetPreallocation(*J,0,dnnz);
1890:   MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);
1891:   MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);

1893:   PetscFree2(dnnz,onnz);

1895:   /* (2) Set matrix entries for edges */
1896:   /*----------------------------------*/
1897:   for (e=eStart; e<eEnd; e++) {
1898:     /* Get row indices */
1899:     DMNetworkGetVariableGlobalOffset(dm,e,&rstart);
1900:     DMNetworkGetNumVariables(dm,e,&nrows);
1901:     if (nrows) {
1902:       for (j=0; j<nrows; j++) rows[j] = j + rstart;

1904:       /* Set matrix entries for conntected vertices */
1905:       DMNetworkGetConnectedVertices(dm,e,&cone);
1906:       for (v=0; v<2; v++) {
1907:         DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);
1908:         DMNetworkGetNumVariables(dm,cone[v],&ncols);

1910:         if (network->Je) {
1911:           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1912:         } else Juser = NULL;
1913:         MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);
1914:       }

1916:       /* Set matrix entries for edge self */
1917:       cstart = rstart;
1918:       if (network->Je) {
1919:         Juser = network->Je[3*e]; /* Jacobian(e,e) */
1920:       } else Juser = NULL;
1921:       MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);
1922:     }
1923:   }

1925:   /* Set matrix entries for vertices */
1926:   /*---------------------------------*/
1927:   for (v=vStart; v<vEnd; v++) {
1928:     /* Get row indices */
1929:     DMNetworkGetVariableGlobalOffset(dm,v,&rstart);
1930:     DMNetworkGetNumVariables(dm,v,&nrows);
1931:     if (!nrows) continue;

1933:     DMNetworkIsGhostVertex(dm,v,&ghost);
1934:     if (ghost) {
1935:       PetscMalloc1(nrows,&rows_v);
1936:     } else {
1937:       rows_v = rows;
1938:     }
1939:     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;

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

1944:     for (e=0; e<nedges; e++) {
1945:       /* Supporting edges */
1946:       DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);
1947:       DMNetworkGetNumVariables(dm,edges[e],&ncols);

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

1954:       /* Connected vertices */
1955:       DMNetworkGetConnectedVertices(dm,edges[e],&cone);
1956:       vc = (v == cone[0]) ? cone[1]:cone[0];

1958:       DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);
1959:       DMNetworkGetNumVariables(dm,vc,&ncols);

1961:       if (network->Jv) {
1962:         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1963:       } else Juser = NULL;
1964:       MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);
1965:     }

1967:     /* Set matrix entries for vertex self */
1968:     if (!ghost) {
1969:       DMNetworkGetVariableGlobalOffset(dm,v,&cstart);
1970:       if (network->Jv) {
1971:         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1972:       } else Juser = NULL;
1973:       MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);
1974:     }
1975:     if (ghost) {
1976:       PetscFree(rows_v);
1977:     }
1978:   }
1979:   PetscFree(rows);

1981:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
1982:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);

1984:   MatSetDM(*J,dm);
1985:   return(0);
1986: }

1988: PetscErrorCode DMDestroy_Network(DM dm)
1989: {
1991:   DM_Network     *network = (DM_Network*)dm->data;
1992:   PetscInt       j;

1995:   if (--network->refct > 0) return(0);
1996:   if (network->Je) {
1997:     PetscFree(network->Je);
1998:   }
1999:   if (network->Jv) {
2000:     PetscFree(network->Jvptr);
2001:     PetscFree(network->Jv);
2002:   }

2004:   ISLocalToGlobalMappingDestroy(&network->vertex.mapping);
2005:   PetscSectionDestroy(&network->vertex.DofSection);
2006:   PetscSectionDestroy(&network->vertex.GlobalDofSection);
2007:   if (network->vltog) {
2008:     PetscFree(network->vltog);
2009:   }
2010:   if (network->vertex.sf) {
2011:     PetscSFDestroy(&network->vertex.sf);
2012:   }
2013:   /* edge */
2014:   ISLocalToGlobalMappingDestroy(&network->edge.mapping);
2015:   PetscSectionDestroy(&network->edge.DofSection);
2016:   PetscSectionDestroy(&network->edge.GlobalDofSection);
2017:   if (network->edge.sf) {
2018:     PetscSFDestroy(&network->edge.sf);
2019:   }
2020:   DMDestroy(&network->plex);
2021:   PetscSectionDestroy(&network->DataSection);
2022:   PetscSectionDestroy(&network->DofSection);

2024:   for(j=0; j<network->nsubnet; j++) {
2025:     PetscFree(network->subnet[j].edges);
2026:   }
2027:   PetscFree(network->subnetvtx);

2029:   PetscFree(network->subnet);
2030:   PetscFree(network->componentdataarray);
2031:   PetscFree2(network->header,network->cvalue);
2032:   PetscFree(network);
2033:   return(0);
2034: }

2036: PetscErrorCode DMView_Network(DM dm,PetscViewer viewer)
2037: {
2039:   DM_Network     *network = (DM_Network*)dm->data;
2040:   PetscBool      iascii;
2041:   PetscMPIInt    rank;
2042:   PetscInt       p,nsubnet;

2045:   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first");
2046:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
2049:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
2050:   if (iascii) {
2051:     const PetscInt    *cone,*vtx,*edges;
2052:     PetscInt          vfrom,vto,i,j,nv,ne;

2054:     nsubnet = network->nsubnet - network->ncsubnet; /* num of subnetworks */
2055:     PetscViewerASCIIPushSynchronized(viewer);
2056:     PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] nsubnet: %D; nsubnetCouple: %D; nEdges: %D; nVertices: %D\n",rank,nsubnet,network->ncsubnet,network->nEdges,network->nVertices);

2058:     for (i=0; i<nsubnet; i++) {
2059:       DMNetworkGetSubnetworkInfo(dm,i,&nv,&ne,&vtx,&edges);
2060:       if (ne) {
2061:         PetscViewerASCIISynchronizedPrintf(viewer, "     Subnet %D: nEdges %D, nVertices %D\n",i,ne,nv);
2062:         for (j=0; j<ne; j++) {
2063:           p = edges[j];
2064:           DMNetworkGetConnectedVertices(dm,p,&cone);
2065:           DMNetworkGetGlobalVertexIndex(dm,cone[0],&vfrom);
2066:           DMNetworkGetGlobalVertexIndex(dm,cone[1],&vto);
2067:           DMNetworkGetGlobalEdgeIndex(dm,edges[j],&p);
2068:           PetscViewerASCIISynchronizedPrintf(viewer, "       edge %D: %D----> %D\n",p,vfrom,vto);
2069:         }
2070:       }
2071:     }
2072:     /* Coupling subnets */
2073:     nsubnet = network->nsubnet;
2074:     for (; i<nsubnet; i++) {
2075:       DMNetworkGetSubnetworkInfo(dm,i,&nv,&ne,&vtx,&edges);
2076:       if (ne) {
2077:         PetscViewerASCIISynchronizedPrintf(viewer, "     Subnet %D (couple): nEdges %D, nVertices %D\n",i,ne,nv);
2078:         for (j=0; j<ne; j++) {
2079:           p = edges[j];
2080:           DMNetworkGetConnectedVertices(dm,p,&cone);
2081:           DMNetworkGetGlobalVertexIndex(dm,cone[0],&vfrom);
2082:           DMNetworkGetGlobalVertexIndex(dm,cone[1],&vto);
2083:           PetscViewerASCIISynchronizedPrintf(viewer, "       edge %D: %D----> %D\n",p,vfrom,vto);
2084:         }
2085:       }
2086:     }
2087:     PetscViewerFlush(viewer);
2088:     PetscViewerASCIIPopSynchronized(viewer);
2089:   } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMNetwork writing", ((PetscObject)viewer)->type_name);
2090:   return(0);
2091: }

2093: PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
2094: {
2096:   DM_Network     *network = (DM_Network*)dm->data;

2099:   DMGlobalToLocalBegin(network->plex,g,mode,l);
2100:   return(0);
2101: }

2103: PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
2104: {
2106:   DM_Network     *network = (DM_Network*)dm->data;

2109:   DMGlobalToLocalEnd(network->plex,g,mode,l);
2110:   return(0);
2111: }

2113: PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
2114: {
2116:   DM_Network     *network = (DM_Network*)dm->data;

2119:   DMLocalToGlobalBegin(network->plex,l,mode,g);
2120:   return(0);
2121: }

2123: PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
2124: {
2126:   DM_Network     *network = (DM_Network*)dm->data;

2129:   DMLocalToGlobalEnd(network->plex,l,mode,g);
2130:   return(0);
2131: }

2133: /*@
2134:   DMNetworkGetVertexLocalToGlobalOrdering - Get vertex globle index

2136:   Not collective

2138:   Input Parameters
2139: + dm - the dm object
2140: - vloc - local vertex ordering, start from 0

2142:   Output Parameters
2143: +  vg  - global vertex ordering, start from 0

2145:   Level: Advanced

2147: .seealso: DMNetworkSetVertexLocalToGlobalOrdering()
2148: @*/
2149: PetscErrorCode DMNetworkGetVertexLocalToGlobalOrdering(DM dm,PetscInt vloc,PetscInt *vg)
2150: {
2151:   DM_Network  *network = (DM_Network*)dm->data;
2152:   PetscInt    *vltog = network->vltog;

2155:   if (!vltog) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Must call DMNetworkSetVertexLocalToGlobalOrdering() first");
2156:   *vg = vltog[vloc];
2157:   return(0);
2158: }

2160: /*@
2161:   DMNetworkSetVertexLocalToGlobalOrdering - Create and setup vertex local to globle map

2163:   Collective

2165:   Input Parameters:
2166: + dm - the dm object

2168:   Level: Advanced

2170: .seealso: DMNetworkGetGlobalVertexIndex()
2171: @*/
2172: PetscErrorCode DMNetworkSetVertexLocalToGlobalOrdering(DM dm)
2173: {
2174:   PetscErrorCode    ierr;
2175:   DM_Network        *network=(DM_Network*)dm->data;
2176:   MPI_Comm          comm;
2177:   PetscMPIInt       rank,size,*displs,*recvcounts,remoterank;
2178:   PetscBool         ghost;
2179:   PetscInt          *vltog,nroots,nleaves,i,*vrange,k,N,lidx;
2180:   const PetscSFNode *iremote;
2181:   PetscSF           vsf;
2182:   Vec               Vleaves,Vleaves_seq;
2183:   VecScatter        ctx;
2184:   PetscScalar       *varr,val;
2185:   const PetscScalar *varr_read;

2188:   PetscObjectGetComm((PetscObject)dm,&comm);
2189:   MPI_Comm_size(comm,&size);
2190:   MPI_Comm_rank(comm,&rank);

2192:   if (size == 1) {
2193:     nroots = network->vEnd - network->vStart;
2194:     PetscMalloc1(nroots, &vltog);
2195:     for (i=0; i<nroots; i++) vltog[i] = i;
2196:     network->vltog = vltog;
2197:     return(0);
2198:   }

2200:   if (!network->distributecalled) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE,"Must call DMNetworkDistribute() first");
2201:   if (network->vltog) {
2202:     PetscFree(network->vltog);
2203:   }

2205:   DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);
2206:   PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);
2207:   vsf = network->vertex.sf;

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

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

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

2219:   PetscMalloc1(nroots, &vltog);
2220:   network->vltog = vltog;

2222:   /* Set vltog for non-ghost vertices */
2223:   k = 0;
2224:   for (i=0; i<nroots; i++) {
2225:     DMNetworkIsGhostVertex(dm,i+network->vStart,&ghost);
2226:     if (ghost) continue;
2227:     vltog[i] = vrange[rank] + k++;
2228:   }
2229:   PetscFree3(vrange,displs,recvcounts);

2231:   /* Set vltog for ghost vertices */
2232:   /* (a) create parallel Vleaves and sequential Vleaves_seq to convert local iremote[*].index to global index */
2233:   VecCreate(comm,&Vleaves);
2234:   VecSetSizes(Vleaves,2*nleaves,PETSC_DETERMINE);
2235:   VecSetFromOptions(Vleaves);
2236:   VecGetArray(Vleaves,&varr);
2237:   for (i=0; i<nleaves; i++) {
2238:     varr[2*i]   = (PetscScalar)(iremote[i].rank);  /* rank of remote process */
2239:     varr[2*i+1] = (PetscScalar)(iremote[i].index); /* local index in remote process */
2240:   }
2241:   VecRestoreArray(Vleaves,&varr);

2243:   /* (b) scatter local info to remote processes via VecScatter() */
2244:   VecScatterCreateToAll(Vleaves,&ctx,&Vleaves_seq);
2245:   VecScatterBegin(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD);
2246:   VecScatterEnd(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD);

2248:   /* (c) convert local indices to global indices in parallel vector Vleaves */
2249:   VecGetSize(Vleaves_seq,&N);
2250:   VecGetArrayRead(Vleaves_seq,&varr_read);
2251:   for (i=0; i<N; i+=2) {
2252:     remoterank = (PetscMPIInt)PetscRealPart(varr_read[i]);
2253:     if (remoterank == rank) {
2254:       k = i+1; /* row number */
2255:       lidx = (PetscInt)PetscRealPart(varr_read[i+1]);
2256:       val  = (PetscScalar)vltog[lidx]; /* global index for non-ghost vertex computed above */
2257:       VecSetValues(Vleaves,1,&k,&val,INSERT_VALUES);
2258:     }
2259:   }
2260:   VecRestoreArrayRead(Vleaves_seq,&varr_read);
2261:   VecAssemblyBegin(Vleaves);
2262:   VecAssemblyEnd(Vleaves);

2264:   /* (d) Set vltog for ghost vertices by copying local values of Vleaves */
2265:   VecGetArrayRead(Vleaves,&varr_read);
2266:   k = 0;
2267:   for (i=0; i<nroots; i++) {
2268:     DMNetworkIsGhostVertex(dm,i+network->vStart,&ghost);
2269:     if (!ghost) continue;
2270:     vltog[i] = (PetscInt)PetscRealPart(varr_read[2*k+1]); k++;
2271:   }
2272:   VecRestoreArrayRead(Vleaves,&varr_read);

2274:   VecDestroy(&Vleaves);
2275:   VecDestroy(&Vleaves_seq);
2276:   VecScatterDestroy(&ctx);
2277:   return(0);
2278: }