Actual source code: network.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1:  #include <petsc/private/dmnetworkimpl.h>
  2:  #include <petscdmplex.h>
  3:  #include <petscsf.h>

  5: /*@
  6:   DMNetworkGetPlex - Gets the Plex DM associated with this network DM

  8:   Not collective

 10:   Input Parameters:
 11: + netdm - the dm object
 12: - plexmdm - the plex dm object

 14:   Level: Advanced

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

 23:   *plexdm = network->plex;
 24:   return(0);
 25: }

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

 30:   Collective on DM

 32:   Input Parameters:
 33: + dm - the dm object
 34: . Nsubnet - number of subnetworks
 35: . NsubnetCouple - number of coupling subnetworks
 36: . nV - number of local vertices for each subnetwork and coupling subnetwork
 37: . nE - number of local edges for each subnetwork and coupling subnetwork
 38: . NV - number of global vertices (or PETSC_DETERMINE) for each subnetwork and coupling subnetwork
 39: - NE - number of global edges (or PETSC_DETERMINE) for each subnetwork and coupling subnetwork

 41:    Notes
 42:    If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang.

 44:    You cannot change the sizes once they have been set. nV, nE, NV and NE are arrays of length Nsubnet+NsubnetCouple.

 46:    Level: intermediate

 48: .seealso: DMNetworkCreate()
 49: @*/
 50: PetscErrorCode DMNetworkSetSizes(DM dm,PetscInt Nsubnet,PetscInt NsubnetCouple,PetscInt nV[], PetscInt nE[], PetscInt NV[], PetscInt NE[])
 51: {
 53:   DM_Network     *network = (DM_Network*) dm->data;
 54:   PetscInt       a[2],b[2],i;

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

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

 64:   network->nsubnet  = Nsubnet + NsubnetCouple;
 65:   Nsubnet          += NsubnetCouple;
 66:   network->ncsubnet = NsubnetCouple;
 67:   PetscCalloc1(Nsubnet,&network->subnet);
 68:   for(i=0; i < Nsubnet; i++) {
 69:     if (NV) {
 71:       if (NV[i] > 0 && nV[i] > NV[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local vertex size %D cannot be larger than global vertex size %D",i,nV[i],NV[i]);
 72:     }
 73:     if (NE) {
 75:       if (NE[i] > 0 && nE[i] > NE[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local edge size %D cannot be larger than global edge size %D",i,nE[i],NE[i]);
 76:     }

 78:     a[0] = nV[i]; a[1] = nE[i];
 79:     MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));
 80:     network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1];

 82:     network->subnet[i].id = i;

 84:     network->subnet[i].nvtx = nV[i];
 85:     network->subnet[i].vStart = network->nVertices;
 86:     network->subnet[i].vEnd   = network->subnet[i].vStart + network->subnet[i].Nvtx;
 87:     network->nVertices += network->subnet[i].nvtx;
 88:     network->NVertices += network->subnet[i].Nvtx;

 90:     network->subnet[i].nedge = nE[i];
 91:     network->subnet[i].eStart = network->nEdges;
 92:     network->subnet[i].eEnd = network->subnet[i].eStart + network->subnet[i].Nedge;
 93:     network->nEdges += network->subnet[i].nedge;
 94:     network->NEdges += network->subnet[i].Nedge;
 95:   }
 96:   return(0);
 97: }

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

102:   Logically collective on DM

104:   Input Parameters:
105: + dm - the dm object
106: . edgelist - list of edges for each subnetwork
107: - edgelistCouple - list of edges for each coupling subnetwork

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

113:   Level: intermediate

115:   Example usage:
116:   Consider the following 2 separate networks and a coupling network:

118: .vb
119:  network 0: v0 -> v1 -> v2 -> v3
120:  network 1: v1 -> v2 -> v0
121:  coupling network: network 1: v2 -> network 0: v0
122: .ve

124:  The resulting input
125:    edgelist[0] = [0 1 | 1 2 | 2 3];
126:    edgelist[1] = [1 2 | 2 0]
127:    edgelistCouple[0] = [(network)1 (v)2 (network)0 (v)0].

129: .seealso: DMNetworkCreate, DMNetworkSetSizes
130: @*/
131: PetscErrorCode DMNetworkSetEdgeList(DM dm,PetscInt *edgelist[],PetscInt *edgelistCouple[])
132: {
133:   DM_Network *network = (DM_Network*) dm->data;
134:   PetscInt   i,nsubnet,ncsubnet=network->ncsubnet;

137:   nsubnet = network->nsubnet - ncsubnet;
138:   for(i=0; i < nsubnet; i++) {
139:     network->subnet[i].edgelist = edgelist[i];
140:   }
141:   if (edgelistCouple) {
142:     PetscInt j;
143:     j = 0;
144:     nsubnet = network->nsubnet;
145:     while (i < nsubnet) network->subnet[i++].edgelist = edgelistCouple[j++];
146:   }
147:   return(0);
148: }

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

153:   Collective on DM

155:   Input Parameters
156: . DM - the dmnetwork object

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

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

164:   Level: intermediate

166: .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList
167: @*/
168: PetscErrorCode DMNetworkLayoutSetUp(DM dm)
169: {
171:   DM_Network     *network = (DM_Network*)dm->data;
172:   PetscInt       dim = 1; /* One dimensional network */
173:   PetscInt       numCorners=2,spacedim=2;
174:   double         *vertexcoords=NULL;
175:   PetscInt       i,j,ndata,ctr=0,nsubnet;
176:   PetscInt       k,netid,vid;
177:   PetscInt       *edgelist_couple=NULL;

180:   if (network->nVertices) {
181:     PetscCalloc1(numCorners*network->nVertices,&vertexcoords);
182:   }

184:   /* Create the edgelist for the network by concatenating edgelists of the subnetworks */
185:   nsubnet = network->nsubnet - network->ncsubnet;
186:   PetscCalloc1(2*network->nEdges,&network->edges);
187:   for (i=0; i < nsubnet; i++) {
188:     for (j = 0; j < network->subnet[i].nedge; j++) {
189:       network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
190:       network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
191:       ctr++;
192:     }
193:   }

195:   i       = nsubnet; /* netid of coupling subnet */
196:   nsubnet = network->nsubnet;
197:   while (i < nsubnet) {
198:     edgelist_couple = network->subnet[i].edgelist;
199:     k = 0;
200:     for (j = 0; j < network->subnet[i].nedge; j++) {
201:       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
202:       network->edges[2*ctr] = network->subnet[netid].vStart + vid; k += 2;

204:       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
205:       network->edges[2*ctr+1] = network->subnet[netid].vStart + vid; k+=2;
206:       ctr++;
207:     }
208:     i++;
209:   }

211: #if 0
212:   for(i=0; i < network->nEdges; i++) {
213:     PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);
214:     printf("\n");
215:   }
216: #endif

218: #if defined(PETSC_USE_64BIT_INDICES)
219:   {
220:     int      *edges;
221:     PetscInt ii;
222:     PetscMalloc1(network->nEdges*numCorners,&edges);
223:     for (ii=0; ii<network->nEdges*numCorners; ii++) {
224:       edges[ii] = (int) network->edges[ii];
225:     }
226:     DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,edges,spacedim,vertexcoords,&network->plex);
227:     PetscFree(edges);
228:   }
229: #else
230:     DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);
231:  #endif

233:   if (network->nVertices) {
234:     PetscFree(vertexcoords);
235:   }
236:   PetscFree(network->edges);

238:   DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);
239:   DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);
240:   DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);

242:   PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);
243:   PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);
244:   PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);
245:   PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);

247:   /* Create vertices and edges array for the subnetworks */
248:   for (j=0; j < network->nsubnet; j++) {
249:     PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);
250:     PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);
251:     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
252:        These get updated when the vertices and edges are added. */
253:     network->subnet[j].nvtx = network->subnet[j].nedge = 0;
254:   }

256:   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
257:   PetscCalloc1(network->pEnd-network->pStart,&network->header);
258:   for (i=network->eStart; i < network->eEnd; i++) {
259:     network->header[i].index = i;   /* Global edge number */
260:     for (j=0; j < network->nsubnet; j++) {
261:       if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) {
262:         network->header[i].subnetid = j; /* Subnetwork id */
263:         network->subnet[j].edges[network->subnet[j].nedge++] = i;
264:         break;
265:       }
266:     }
267:     network->header[i].ndata = 0;
268:     ndata = network->header[i].ndata;
269:     PetscSectionAddDof(network->DataSection,i,network->dataheadersize);
270:     network->header[i].offset[ndata] = 0;
271:   }

273:   for(i=network->vStart; i < network->vEnd; i++) {
274:     network->header[i].index = i - network->vStart;
275:     for (j=0; j < network->nsubnet; j++) {
276:       if ((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) {
277:         network->header[i].subnetid = j;
278:         network->subnet[j].vertices[network->subnet[j].nvtx++] = i;
279:         break;
280:       }
281:     }
282:     network->header[i].ndata = 0;
283:     ndata = network->header[i].ndata;
284:     PetscSectionAddDof(network->DataSection,i,network->dataheadersize);
285:     network->header[i].offset[ndata] = 0;
286:   }

288:   PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);
289:   return(0);
290: }

292: /*@C
293:   DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork

295:   Input Parameters
296: + dm   - the number object
297: - id   - the ID (integer) of the subnetwork

299:   Output Parameters
300: + nv    - number of vertices (local)
301: . ne    - number of edges (local)
302: . vtx   - local vertices for this subnetwork
303: . edge  - local edges for this subnetwork

305:   Notes:
306:   Cannot call this routine before DMNetworkLayoutSetup()

308:   Level: intermediate

310: .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
311: @*/
312: PetscErrorCode DMNetworkGetSubnetworkInfo(DM netdm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge)
313: {
314:   DM_Network     *network = (DM_Network*) netdm->data;

317:   *nv = network->subnet[id].nvtx;
318:   *ne = network->subnet[id].nedge;
319:   *vtx = network->subnet[id].vertices;
320:   *edge = network->subnet[id].edges;
321:   return(0);
322: }

324: /*@C
325:   DMNetworkRegisterComponent - Registers the network component

327:   Logically collective on DM

329:   Input Parameters
330: + dm   - the network object
331: . name - the component name
332: - size - the storage size in bytes for this component data

334:    Output Parameters
335: .   key - an integer key that defines the component

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

340:    Level: intermediate

342: .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
343: @*/
344: PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key)
345: {
346:   PetscErrorCode        ierr;
347:   DM_Network            *network = (DM_Network*) dm->data;
348:   DMNetworkComponent    *component=&network->component[network->ncomponent];
349:   PetscBool             flg=PETSC_FALSE;
350:   PetscInt              i;

353:   for (i=0; i < network->ncomponent; i++) {
354:     PetscStrcmp(component->name,name,&flg);
355:     if (flg) {
356:       *key = i;
357:       return(0);
358:     }
359:   }
360:   if(network->ncomponent == MAX_COMPONENTS) {
361:     SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components registered exceeds the max %D",MAX_COMPONENTS);
362:   }

364:   PetscStrcpy(component->name,name);
365:   component->size = size/sizeof(DMNetworkComponentGenericDataType);
366:   *key = network->ncomponent;
367:   network->ncomponent++;
368:   return(0);
369: }

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

374:   Not Collective

376:   Input Parameters:
377: + dm - The DMNetwork object

379:   Output Paramters:
380: + vStart - The first vertex point
381: - vEnd   - One beyond the last vertex point

383:   Level: intermediate

385: .seealso: DMNetworkGetEdgeRange
386: @*/
387: PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
388: {
389:   DM_Network     *network = (DM_Network*)dm->data;

392:   if (vStart) *vStart = network->vStart;
393:   if (vEnd) *vEnd = network->vEnd;
394:   return(0);
395: }

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

400:   Not Collective

402:   Input Parameters:
403: + dm - The DMNetwork object

405:   Output Paramters:
406: + eStart - The first edge point
407: - eEnd   - One beyond the last edge point

409:   Level: intermediate

411: .seealso: DMNetworkGetVertexRange
412: @*/
413: PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
414: {
415:   DM_Network     *network = (DM_Network*)dm->data;

418:   if (eStart) *eStart = network->eStart;
419:   if (eEnd) *eEnd = network->eEnd;
420:   return(0);
421: }

423: /*@
424:   DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge.

426:   Not Collective

428:   Input Parameters:
429: + dm - DMNetwork object
430: - p  - edge point

432:   Output Paramters:
433: . index - user global numbering for the edge

435:   Level: intermediate

437: .seealso: DMNetworkGetGlobalVertexIndex
438: @*/
439: PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
440: {
441:   PetscErrorCode    ierr;
442:   DM_Network        *network = (DM_Network*)dm->data;
443:   PetscInt          offsetp;
444:   DMNetworkComponentHeader header;

447:   PetscSectionGetOffset(network->DataSection,p,&offsetp);
448:   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
449:   *index = header->index;
450:   return(0);
451: }

453: /*@
454:   DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex.

456:   Not Collective

458:   Input Parameters:
459: + dm - DMNetwork object
460: - p  - vertex point

462:   Output Paramters:
463: . index - user global numbering for the vertex

465:   Level: intermediate

467: .seealso: DMNetworkGetGlobalEdgeIndex
468: @*/
469: PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
470: {
471:   PetscErrorCode    ierr;
472:   DM_Network        *network = (DM_Network*)dm->data;
473:   PetscInt          offsetp;
474:   DMNetworkComponentHeader header;

477:   PetscSectionGetOffset(network->DataSection,p,&offsetp);
478:   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
479:   *index = header->index;
480:   return(0);
481: }

483: /*
484:   DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the
485:                                     component value from the component data array

487:   Not Collective

489:   Input Parameters:
490: + dm      - The DMNetwork object
491: . p       - vertex/edge point
492: - compnum - component number

494:   Output Parameters:
495: + compkey - the key obtained when registering the component
496: - offset  - offset into the component data array associated with the vertex/edge point

498:   Notes:
499:   Typical usage:

501:   DMNetworkGetComponentDataArray(dm, &arr);
502:   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
503:   Loop over vertices or edges
504:     DMNetworkGetNumComponents(dm,v,&numcomps);
505:     Loop over numcomps
506:       DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset);
507:       compdata = (UserCompDataType)(arr+offset);

509:   Level: intermediate

511: .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
512: */
513: PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
514: {
515:   PetscErrorCode           ierr;
516:   PetscInt                 offsetp;
517:   DMNetworkComponentHeader header;
518:   DM_Network               *network = (DM_Network*)dm->data;

521:   PetscSectionGetOffset(network->DataSection,p,&offsetp);
522:   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
523:   if (compkey) *compkey = header->key[compnum];
524:   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
525:   return(0);
526: }

528: /*@
529:   DMNetworkGetComponent - Returns the network component and its key

531:   Not Collective

533:   Input Parameters
534: + dm - DMNetwork object
535: . p  - edge or vertex point
536: - compnum - component number

538:   Output Parameters:
539: + compkey - the key set for this computing during registration
540: - component - the component data

542:   Notes:
543:   Typical usage:

545:   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
546:   Loop over vertices or edges
547:     DMNetworkGetNumComponents(dm,v,&numcomps);
548:     Loop over numcomps
549:       DMNetworkGetComponent(dm,v,compnum,&key,&component);

551:   Level: intermediate

553: .seealso: DMNetworkGetNumComponents, DMNetworkGetVariableOffset
554: @*/
555: PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *key, void **component)
556: {
558:   DM_Network     *network = (DM_Network*)dm->data;
559:   PetscInt       offsetd = 0;

562:   DMNetworkGetComponentKeyOffset(dm,p,compnum,key,&offsetd);
563:   *component = network->componentdataarray+offsetd;
564:   return(0);
565: }

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

570:   Not Collective

572:   Input Parameters:
573: + dm           - The DMNetwork object
574: . p            - vertex/edge point
575: . componentkey - component key returned while registering the component
576: - compvalue    - pointer to the data structure for the component

578:   Level: intermediate

580: .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
581: @*/
582: PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
583: {
584:   DM_Network               *network = (DM_Network*)dm->data;
585:   DMNetworkComponent       *component = &network->component[componentkey];
586:   DMNetworkComponentHeader header = &network->header[p];
587:   DMNetworkComponentValue  cvalue = &network->cvalue[p];
588:   PetscErrorCode           ierr;

591:   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);

593:   header->size[header->ndata] = component->size;
594:   PetscSectionAddDof(network->DataSection,p,component->size);
595:   header->key[header->ndata] = componentkey;
596:   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];

598:   cvalue->data[header->ndata] = (void*)compvalue;
599:   header->ndata++;
600:   return(0);
601: }

603: /*@
604:   DMNetworkGetNumComponents - Get the number of components at a vertex/edge

606:   Not Collective

608:   Input Parameters:
609: + dm - The DMNetwork object
610: . p  - vertex/edge point

612:   Output Parameters:
613: . numcomponents - Number of components at the vertex/edge

615:   Level: intermediate

617: .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
618: @*/
619: PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
620: {
622:   PetscInt       offset;
623:   DM_Network     *network = (DM_Network*)dm->data;

626:   PetscSectionGetOffset(network->DataSection,p,&offset);
627:   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
628:   return(0);
629: }

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

634:   Not Collective

636:   Input Parameters:
637: + dm     - The DMNetwork object
638: - p      - the edge/vertex point

640:   Output Parameters:
641: . offset - the offset

643:   Level: intermediate

645: .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
646: @*/
647: PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
648: {
650:   DM_Network     *network = (DM_Network*)dm->data;

653:   PetscSectionGetOffset(network->plex->defaultSection,p,offset);
654:   return(0);
655: }

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

660:   Not Collective

662:   Input Parameters:
663: + dm      - The DMNetwork object
664: - p       - the edge/vertex point

666:   Output Parameters:
667: . offsetg - the offset

669:   Level: intermediate

671: .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
672: @*/
673: PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
674: {
676:   DM_Network     *network = (DM_Network*)dm->data;

679:   PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);
680:   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */
681:   return(0);
682: }

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

687:   Not Collective

689:   Input Parameters:
690: + dm     - The DMNetwork object
691: - p      - the edge point

693:   Output Parameters:
694: . offset - the offset

696:   Level: intermediate

698: .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
699: @*/
700: PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
701: {
703:   DM_Network     *network = (DM_Network*)dm->data;


707:   PetscSectionGetOffset(network->edge.DofSection,p,offset);
708:   return(0);
709: }

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

714:   Not Collective

716:   Input Parameters:
717: + dm     - The DMNetwork object
718: - p      - the vertex point

720:   Output Parameters:
721: . offset - the offset

723:   Level: intermediate

725: .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
726: @*/
727: PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
728: {
730:   DM_Network     *network = (DM_Network*)dm->data;


734:   p -= network->vStart;

736:   PetscSectionGetOffset(network->vertex.DofSection,p,offset);
737:   return(0);
738: }
739: /*@
740:   DMNetworkAddNumVariables - Add number of variables associated with a given point.

742:   Not Collective

744:   Input Parameters:
745: + dm   - The DMNetworkObject
746: . p    - the vertex/edge point
747: - nvar - number of additional variables

749:   Level: intermediate

751: .seealso: DMNetworkSetNumVariables
752: @*/
753: PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
754: {
756:   DM_Network     *network = (DM_Network*)dm->data;

759:   PetscSectionAddDof(network->DofSection,p,nvar);
760:   return(0);
761: }

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

766:   Not Collective

768:   Input Parameters:
769: + dm   - The DMNetworkObject
770: - p    - the vertex/edge point

772:   Output Parameters:
773: . nvar - number of variables

775:   Level: intermediate

777: .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
778: @*/
779: PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
780: {
782:   DM_Network     *network = (DM_Network*)dm->data;

785:   PetscSectionGetDof(network->DofSection,p,nvar);
786:   return(0);
787: }

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

792:   Not Collective

794:   Input Parameters:
795: + dm   - The DMNetworkObject
796: . p    - the vertex/edge point
797: - nvar - number of variables

799:   Level: intermediate

801: .seealso: DMNetworkAddNumVariables
802: @*/
803: PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
804: {
806:   DM_Network     *network = (DM_Network*)dm->data;

809:   PetscSectionSetDof(network->DofSection,p,nvar);
810:   return(0);
811: }

813: /* Sets up the array that holds the data for all components and its associated section. This
814:    function is called during DMSetUp() */
815: PetscErrorCode DMNetworkComponentSetUp(DM dm)
816: {
817:   PetscErrorCode           ierr;
818:   DM_Network               *network = (DM_Network*)dm->data;
819:   PetscInt                 arr_size,p,offset,offsetp,ncomp,i;
820:   DMNetworkComponentHeader header;
821:   DMNetworkComponentValue  cvalue;
822:   DMNetworkComponentGenericDataType *componentdataarray;

825:   PetscSectionSetUp(network->DataSection);
826:   PetscSectionGetStorageSize(network->DataSection,&arr_size);
827:   PetscMalloc1(arr_size,&network->componentdataarray);
828:   componentdataarray = network->componentdataarray;
829:   for (p = network->pStart; p < network->pEnd; p++) {
830:     PetscSectionGetOffset(network->DataSection,p,&offsetp);
831:     /* Copy header */
832:     header = &network->header[p];
833:     PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));
834:     /* Copy data */
835:     cvalue = &network->cvalue[p];
836:     ncomp = header->ndata;
837:     for (i = 0; i < ncomp; i++) {
838:       offset = offsetp + network->dataheadersize + header->offset[i];
839:       PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));
840:     }
841:   }
842:   return(0);
843: }

845: /* Sets up the section for dofs. This routine is called during DMSetUp() */
846: PetscErrorCode DMNetworkVariablesSetUp(DM dm)
847: {
849:   DM_Network     *network = (DM_Network*)dm->data;

852:   PetscSectionSetUp(network->DofSection);
853:   return(0);
854: }

856: /*@C
857:   DMNetworkGetComponentDataArray - Returns the component data array

859:   Not Collective

861:   Input Parameters:
862: . dm - The DMNetwork Object

864:   Output Parameters:
865: . componentdataarray - array that holds data for all components

867:   Level: intermediate

869: .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents
870: @*/
871: PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
872: {
873:   DM_Network     *network = (DM_Network*)dm->data;

876:   *componentdataarray = network->componentdataarray;
877:   return(0);
878: }

880: /* Get a subsection from a range of points */
881: PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection)
882: {
884:   PetscInt       i, nvar;

887:   PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);
888:   PetscSectionSetChart(*subsection, 0, pend - pstart);
889:   for (i = pstart; i < pend; i++) {
890:     PetscSectionGetDof(master,i,&nvar);
891:     PetscSectionSetDof(*subsection, i - pstart, nvar);
892:   }

894:   PetscSectionSetUp(*subsection);
895:   return(0);
896: }

898: /* Create a submap of points with a GlobalToLocal structure */
899: PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
900: {
902:   PetscInt       i, *subpoints;

905:   /* Create index sets to map from "points" to "subpoints" */
906:   PetscMalloc1(pend - pstart, &subpoints);
907:   for (i = pstart; i < pend; i++) {
908:     subpoints[i - pstart] = i;
909:   }
910:   ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);
911:   PetscFree(subpoints);
912:   return(0);
913: }

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

918:   Collective

920:   Input Parameters:
921: . dm   - The DMNetworkObject

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

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

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

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

931:   Level: intermediate

933: @*/
934: PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
935: {
937:   MPI_Comm       comm;
938:   PetscMPIInt    rank, size;
939:   DM_Network     *network = (DM_Network*)dm->data;

942:   PetscObjectGetComm((PetscObject)dm,&comm);
943:   MPI_Comm_rank(comm, &rank);
944:   MPI_Comm_size(comm, &size);

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

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

954:   if (size > 1) {
955:     PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);
956:     PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);
957:   PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);
958:   PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);
959:   } else {
960:   /* create structures for vertex */
961:   PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);
962:   /* create structures for edge */
963:   PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);
964:   }


967:   /* Add viewers */
968:   PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");
969:   PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");
970:   PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");
971:   PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");

973:   return(0);
974: }

976: /*@
977:   DMNetworkDistribute - Distributes the network and moves associated component data.

979:   Collective

981:   Input Parameter:
982: + DM - the DMNetwork object
983: - overlap - The overlap of partitions, 0 is the default

985:   Notes:
986:   Distributes the network with <overlap>-overlapping partitioning of the edges.

988:   Level: intermediate

990: .seealso: DMNetworkCreate
991: @*/
992: PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
993: {
994:   MPI_Comm       comm;
996:   PetscMPIInt    size;
997:   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
998:   DM_Network     *newDMnetwork;
999:   PetscSF        pointsf;
1000:   DM             newDM;
1001:   PetscPartitioner part;
1002:   PetscInt         j,e,v,offset;
1003:   DMNetworkComponentHeader header;


1007:   PetscObjectGetComm((PetscObject)*dm,&comm);
1008:   MPI_Comm_size(comm, &size);
1009:   if (size == 1) return(0);

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

1015:   /* Enable runtime options for petscpartitioner */
1016:   DMPlexGetPartitioner(oldDMnetwork->plex,&part);
1017:   PetscPartitionerSetFromOptions(part);

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

1022:   /* Distribute dof section */
1023:   PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);
1024:   PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);
1025:   PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);

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

1030:   PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);
1031:   DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);
1032:   DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);
1033:   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
1034:   newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
1035:   newDMnetwork->NVertices = oldDMnetwork->NVertices;
1036:   newDMnetwork->NEdges = oldDMnetwork->NEdges;

1038:   /* Set Dof section as the default section for dm */
1039:   DMSetSection(newDMnetwork->plex,newDMnetwork->DofSection);
1040:   DMGetGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);

1042:   /* Set up subnetwork info in the newDM */
1043:   newDMnetwork->nsubnet = oldDMnetwork->nsubnet;
1044:   PetscCalloc1(newDMnetwork->nsubnet,&newDMnetwork->subnet);
1045:   /* Copy over the global number of vertices and edges in each subnetwork. Note that these are already
1046:      calculated in DMNetworkLayoutSetUp()
1047:   */
1048:   for(j=0; j < newDMnetwork->nsubnet; j++) {
1049:     newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx;
1050:     newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge;
1051:   }

1053:   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1054:     PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);
1055:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1056:     newDMnetwork->subnet[header->subnetid].nedge++;
1057:   }

1059:   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1060:     PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);
1061:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1062:     newDMnetwork->subnet[header->subnetid].nvtx++;
1063:   }

1065:   /* Now create the vertices and edge arrays for the subnetworks */
1066:   for(j=0; j < newDMnetwork->nsubnet; j++) {
1067:     PetscCalloc1(newDMnetwork->subnet[j].nedge,&newDMnetwork->subnet[j].edges);
1068:     PetscCalloc1(newDMnetwork->subnet[j].nvtx,&newDMnetwork->subnet[j].vertices);
1069:     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
1070:        These get updated when the vertices and edges are added. */
1071:     newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0;
1072:   }

1074:   /* Set the vertices and edges in each subnetwork */
1075:   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1076:     PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);
1077:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1078:     newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++]  = e;
1079:   }

1081:   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1082:     PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);
1083:     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1084:     newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++]  = v;
1085:   }

1087:   /* Destroy point SF */
1088:   PetscSFDestroy(&pointsf);

1090:   DMDestroy(dm);
1091:   *dm  = newDM;
1092:   return(0);
1093: }

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

1098:   Input Parameters:
1099: + masterSF - the original SF structure
1100: - map      - a ISLocalToGlobal mapping that contains the subset of points

1102:   Output Parameters:
1103: . subSF    - a subset of the masterSF for the desired subset.
1104: */

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

1108:   PetscErrorCode        ierr;
1109:   PetscInt              nroots, nleaves, *ilocal_sub;
1110:   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
1111:   PetscInt              *local_points, *remote_points;
1112:   PetscSFNode           *iremote_sub;
1113:   const PetscInt        *ilocal;
1114:   const PetscSFNode     *iremote;

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

1119:   /* Look for leaves that pertain to the subset of points. Get the local ordering */
1120:   PetscMalloc1(nleaves,&ilocal_map);
1121:   ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);
1122:   for (i = 0; i < nleaves; i++) {
1123:     if (ilocal_map[i] != -1) nleaves_sub += 1;
1124:   }
1125:   /* Re-number ilocal with subset numbering. Need information from roots */
1126:   PetscMalloc2(nroots,&local_points,nroots,&remote_points);
1127:   for (i = 0; i < nroots; i++) local_points[i] = i;
1128:   ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);
1129:   PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);
1130:   PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);
1131:   /* Fill up graph using local (that is, local to the subset) numbering. */
1132:   PetscMalloc1(nleaves_sub,&ilocal_sub);
1133:   PetscMalloc1(nleaves_sub,&iremote_sub);
1134:   nleaves_sub = 0;
1135:   for (i = 0; i < nleaves; i++) {
1136:     if (ilocal_map[i] != -1) {
1137:       ilocal_sub[nleaves_sub] = ilocal_map[i];
1138:       iremote_sub[nleaves_sub].rank = iremote[i].rank;
1139:       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
1140:       nleaves_sub += 1;
1141:     }
1142:   }
1143:   PetscFree2(local_points,remote_points);
1144:   ISLocalToGlobalMappingGetSize(map,&nroots_sub);

1146:   /* Create new subSF */
1147:   PetscSFCreate(PETSC_COMM_WORLD,subSF);
1148:   PetscSFSetFromOptions(*subSF);
1149:   PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);
1150:   PetscFree(ilocal_map);
1151:   PetscFree(iremote_sub);
1152:   return(0);
1153: }

1155: /*@C
1156:   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point

1158:   Not Collective

1160:   Input Parameters:
1161: + dm - The DMNetwork object
1162: - p  - the vertex point

1164:   Output Paramters:
1165: + nedges - number of edges connected to this vertex point
1166: - edges  - List of edge points

1168:   Level: intermediate

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

1174: .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices
1175: @*/
1176: PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
1177: {
1179:   DM_Network     *network = (DM_Network*)dm->data;

1182:   DMPlexGetSupportSize(network->plex,vertex,nedges);
1183:   DMPlexGetSupport(network->plex,vertex,edges);
1184:   return(0);
1185: }

1187: /*@C
1188:   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point

1190:   Not Collective

1192:   Input Parameters:
1193: + dm - The DMNetwork object
1194: - p  - the edge point

1196:   Output Paramters:
1197: . vertices  - vertices connected to this edge

1199:   Level: intermediate

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

1205: .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
1206: @*/
1207: PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
1208: {
1210:   DM_Network     *network = (DM_Network*)dm->data;

1213:   DMPlexGetCone(network->plex,edge,vertices);
1214:   return(0);
1215: }

1217: /*@
1218:   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex

1220:   Not Collective

1222:   Input Parameters:
1223: + dm - The DMNetwork object
1224: . p  - the vertex point

1226:   Output Parameter:
1227: . isghost - TRUE if the vertex is a ghost point

1229:   Level: intermediate

1231: .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange
1232: @*/
1233: PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
1234: {
1236:   DM_Network     *network = (DM_Network*)dm->data;
1237:   PetscInt       offsetg;
1238:   PetscSection   sectiong;

1241:   *isghost = PETSC_FALSE;
1242:   DMGetGlobalSection(network->plex,&sectiong);
1243:   PetscSectionGetOffset(sectiong,p,&offsetg);
1244:   if (offsetg < 0) *isghost = PETSC_TRUE;
1245:   return(0);
1246: }

1248: PetscErrorCode DMSetUp_Network(DM dm)
1249: {
1251:   DM_Network     *network=(DM_Network*)dm->data;

1254:   DMNetworkComponentSetUp(dm);
1255:   DMNetworkVariablesSetUp(dm);

1257:   DMSetSection(network->plex,network->DofSection);
1258:   DMGetGlobalSection(network->plex,&network->GlobalDofSection);
1259:   return(0);
1260: }

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

1266:     Collective

1268:     Input Parameters:
1269: +   dm - The DMNetwork object
1270: .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
1271: -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices

1273:     Level: intermediate

1275: @*/
1276: PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
1277: {
1278:   DM_Network     *network=(DM_Network*)dm->data;
1280:   PetscInt       nVertices = network->nVertices;

1283:   network->userEdgeJacobian   = eflg;
1284:   network->userVertexJacobian = vflg;

1286:   if (eflg && !network->Je) {
1287:     PetscCalloc1(3*network->nEdges,&network->Je);
1288:   }

1290:   if (vflg && !network->Jv && nVertices) {
1291:     PetscInt       i,*vptr,nedges,vStart=network->vStart;
1292:     PetscInt       nedges_total;
1293:     const PetscInt *edges;

1295:     /* count nvertex_total */
1296:     nedges_total = 0;
1297:     PetscMalloc1(nVertices+1,&vptr);

1299:     vptr[0] = 0;
1300:     for (i=0; i<nVertices; i++) {
1301:       DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);
1302:       nedges_total += nedges;
1303:       vptr[i+1] = vptr[i] + 2*nedges + 1;
1304:     }

1306:     PetscCalloc1(2*nedges_total+nVertices,&network->Jv);
1307:     network->Jvptr = vptr;
1308:   }
1309:   return(0);
1310: }

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

1315:     Not Collective

1317:     Input Parameters:
1318: +   dm - The DMNetwork object
1319: .   p  - the edge point
1320: -   J - array (size = 3) of Jacobian submatrices for this edge point:
1321:         J[0]: this edge
1322:         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()

1324:     Level: intermediate

1326: .seealso: DMNetworkVertexSetMatrix
1327: @*/
1328: PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
1329: {
1330:   DM_Network     *network=(DM_Network*)dm->data;

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

1335:   if (J) {
1336:     network->Je[3*p]   = J[0];
1337:     network->Je[3*p+1] = J[1];
1338:     network->Je[3*p+2] = J[2];
1339:   }
1340:   return(0);
1341: }

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

1346:     Not Collective

1348:     Input Parameters:
1349: +   dm - The DMNetwork object
1350: .   p  - the vertex point
1351: -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
1352:         J[0]:       this vertex
1353:         J[1+2*i]:   i-th supporting edge
1354:         J[1+2*i+1]: i-th connected vertex

1356:     Level: intermediate

1358: .seealso: DMNetworkEdgeSetMatrix
1359: @*/
1360: PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
1361: {
1363:   DM_Network     *network=(DM_Network*)dm->data;
1364:   PetscInt       i,*vptr,nedges,vStart=network->vStart;
1365:   const PetscInt *edges;

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

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

1374:     /* Set Jacobian for each supporting edge and connected vertex */
1375:     DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);
1376:     for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
1377:   }
1378:   return(0);
1379: }

1381: PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1382: {
1384:   PetscInt       j;
1385:   PetscScalar    val=(PetscScalar)ncols;

1388:   if (!ghost) {
1389:     for (j=0; j<nrows; j++) {
1390:       VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);
1391:     }
1392:   } else {
1393:     for (j=0; j<nrows; j++) {
1394:       VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);
1395:     }
1396:   }
1397:   return(0);
1398: }

1400: PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1401: {
1403:   PetscInt       j,ncols_u;
1404:   PetscScalar    val;

1407:   if (!ghost) {
1408:     for (j=0; j<nrows; j++) {
1409:       MatGetRow(Ju,j,&ncols_u,NULL,NULL);
1410:       val = (PetscScalar)ncols_u;
1411:       VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);
1412:       MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);
1413:     }
1414:   } else {
1415:     for (j=0; j<nrows; j++) {
1416:       MatGetRow(Ju,j,&ncols_u,NULL,NULL);
1417:       val = (PetscScalar)ncols_u;
1418:       VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);
1419:       MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);
1420:     }
1421:   }
1422:   return(0);
1423: }

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

1430:   if (Ju) {
1431:     MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);
1432:   } else {
1433:     MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);
1434:   }
1435:   return(0);
1436: }

1438: PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1439: {
1441:   PetscInt       j,*cols;
1442:   PetscScalar    *zeros;

1445:   PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);
1446:   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1447:   MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);
1448:   PetscFree2(cols,zeros);
1449:   return(0);
1450: }

1452: PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1453: {
1455:   PetscInt       j,M,N,row,col,ncols_u;
1456:   const PetscInt *cols;
1457:   PetscScalar    zero=0.0;

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

1463:   for (row=0; row<nrows; row++) {
1464:     MatGetRow(Ju,row,&ncols_u,&cols,NULL);
1465:     for (j=0; j<ncols_u; j++) {
1466:       col = cols[j] + cstart;
1467:       MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);
1468:     }
1469:     MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);
1470:   }
1471:   return(0);
1472: }

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

1479:   if (Ju) {
1480:     MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);
1481:   } else {
1482:     MatSetDenseblock_private(nrows,rows,ncols,cstart,J);
1483:   }
1484:   return(0);
1485: }

1487: /* 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.
1488: */
1489: PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
1490: {
1492:   PetscInt       i, size, dof;
1493:   PetscInt       *glob2loc;

1496:   PetscSectionGetStorageSize(localsec,&size);
1497:   PetscMalloc1(size,&glob2loc);

1499:   for (i = 0; i < size; i++) {
1500:     PetscSectionGetOffset(globalsec,i,&dof);
1501:     dof = (dof >= 0) ? dof : -(dof + 1);
1502:     glob2loc[i] = dof;
1503:   }

1505:   ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);
1506: #if 0
1507:   PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);
1508: #endif
1509:   return(0);
1510: }

1512:  #include <petsc/private/matimpl.h>
1513: PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
1514: {
1516:   PetscMPIInt    rank, size;
1517:   DM_Network     *network = (DM_Network*)dm->data;
1518:   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1519:   PetscInt       cstart,ncols,j,e,v;
1520:   PetscBool      ghost,ghost_vc,ghost2,isNest;
1521:   Mat            Juser;
1522:   PetscSection   sectionGlobal;
1523:   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
1524:   const PetscInt *edges,*cone;
1525:   MPI_Comm       comm;
1526:   MatType        mtype;
1527:   Vec            vd_nz,vo_nz;
1528:   PetscInt       *dnnz,*onnz;
1529:   PetscScalar    *vdnz,*vonz;

1532:   mtype = dm->mattype;
1533:   PetscStrcmp(mtype, MATNEST, &isNest);

1535:   if (isNest) {
1536:     /* DMCreateMatrix_Network_Nest(); */
1537:     PetscInt   eDof, vDof;
1538:     Mat        j11, j12, j21, j22, bA[2][2];
1539:     ISLocalToGlobalMapping eISMap, vISMap;

1541:     PetscObjectGetComm((PetscObject)dm,&comm);
1542:     MPI_Comm_rank(comm,&rank);
1543:     MPI_Comm_size(comm,&size);

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

1548:     MatCreate(comm, &j11);
1549:     MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);
1550:     MatSetType(j11, MATMPIAIJ);

1552:     MatCreate(comm, &j12);
1553:     MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);
1554:     MatSetType(j12, MATMPIAIJ);

1556:     MatCreate(comm, &j21);
1557:     MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);
1558:     MatSetType(j21, MATMPIAIJ);

1560:     MatCreate(comm, &j22);
1561:     MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);
1562:     MatSetType(j22, MATMPIAIJ);

1564:     bA[0][0] = j11;
1565:     bA[0][1] = j12;
1566:     bA[1][0] = j21;
1567:     bA[1][1] = j22;

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

1572:     MatSetLocalToGlobalMapping(j11,eISMap,eISMap);
1573:     MatSetLocalToGlobalMapping(j12,eISMap,vISMap);
1574:     MatSetLocalToGlobalMapping(j21,vISMap,eISMap);
1575:     MatSetLocalToGlobalMapping(j22,vISMap,vISMap);

1577:     MatSetUp(j11);
1578:     MatSetUp(j12);
1579:     MatSetUp(j21);
1580:     MatSetUp(j22);

1582:     MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);
1583:     MatSetUp(*J);
1584:     MatNestSetVecType(*J,VECNEST);
1585:     MatDestroy(&j11);
1586:     MatDestroy(&j12);
1587:     MatDestroy(&j21);
1588:     MatDestroy(&j22);

1590:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
1591:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);

1593:     /* Free structures */
1594:     ISLocalToGlobalMappingDestroy(&eISMap);
1595:     ISLocalToGlobalMappingDestroy(&vISMap);

1597:     return(0);
1598:   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1599:     /* user does not provide Jacobian blocks */
1600:     DMCreateMatrix(network->plex,J);
1601:     MatSetDM(*J,dm);
1602:     return(0);
1603:   }

1605:   MatCreate(PetscObjectComm((PetscObject)dm),J);
1606:   DMGetGlobalSection(network->plex,&sectionGlobal);
1607:   PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);
1608:   MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);

1610:   MatSetType(*J,MATAIJ);
1611:   MatSetFromOptions(*J);

1613:   /* (1) Set matrix preallocation */
1614:   /*------------------------------*/
1615:   PetscObjectGetComm((PetscObject)dm,&comm);
1616:   VecCreate(comm,&vd_nz);
1617:   VecSetSizes(vd_nz,localSize,PETSC_DECIDE);
1618:   VecSetFromOptions(vd_nz);
1619:   VecSet(vd_nz,0.0);
1620:   VecDuplicate(vd_nz,&vo_nz);

1622:   /* Set preallocation for edges */
1623:   /*-----------------------------*/
1624:   DMNetworkGetEdgeRange(dm,&eStart,&eEnd);

1626:   PetscMalloc1(localSize,&rows);
1627:   for (e=eStart; e<eEnd; e++) {
1628:     /* Get row indices */
1629:     DMNetworkGetVariableGlobalOffset(dm,e,&rstart);
1630:     DMNetworkGetNumVariables(dm,e,&nrows);
1631:     if (nrows) {
1632:       for (j=0; j<nrows; j++) rows[j] = j + rstart;

1634:       /* Set preallocation for conntected vertices */
1635:       DMNetworkGetConnectedVertices(dm,e,&cone);
1636:       for (v=0; v<2; v++) {
1637:         DMNetworkGetNumVariables(dm,cone[v],&ncols);

1639:         if (network->Je) {
1640:           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1641:         } else Juser = NULL;
1642:         DMNetworkIsGhostVertex(dm,cone[v],&ghost);
1643:         MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);
1644:       }

1646:       /* Set preallocation for edge self */
1647:       cstart = rstart;
1648:       if (network->Je) {
1649:         Juser = network->Je[3*e]; /* Jacobian(e,e) */
1650:       } else Juser = NULL;
1651:       MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);
1652:     }
1653:   }

1655:   /* Set preallocation for vertices */
1656:   /*--------------------------------*/
1657:   DMNetworkGetVertexRange(dm,&vStart,&vEnd);
1658:   if (vEnd - vStart) vptr = network->Jvptr;

1660:   for (v=vStart; v<vEnd; v++) {
1661:     /* Get row indices */
1662:     DMNetworkGetVariableGlobalOffset(dm,v,&rstart);
1663:     DMNetworkGetNumVariables(dm,v,&nrows);
1664:     if (!nrows) continue;

1666:     DMNetworkIsGhostVertex(dm,v,&ghost);
1667:     if (ghost) {
1668:       PetscMalloc1(nrows,&rows_v);
1669:     } else {
1670:       rows_v = rows;
1671:     }

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

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

1678:     for (e=0; e<nedges; e++) {
1679:       /* Supporting edges */
1680:       DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);
1681:       DMNetworkGetNumVariables(dm,edges[e],&ncols);

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

1688:       /* Connected vertices */
1689:       DMNetworkGetConnectedVertices(dm,edges[e],&cone);
1690:       vc = (v == cone[0]) ? cone[1]:cone[0];
1691:       DMNetworkIsGhostVertex(dm,vc,&ghost_vc);

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

1695:       if (network->Jv) {
1696:         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1697:       } else Juser = NULL;
1698:       if (ghost_vc||ghost) {
1699:         ghost2 = PETSC_TRUE;
1700:       } else {
1701:         ghost2 = PETSC_FALSE;
1702:       }
1703:       MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);
1704:     }

1706:     /* Set preallocation for vertex self */
1707:     DMNetworkIsGhostVertex(dm,v,&ghost);
1708:     if (!ghost) {
1709:       DMNetworkGetVariableGlobalOffset(dm,v,&cstart);
1710:       if (network->Jv) {
1711:         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1712:       } else Juser = NULL;
1713:       MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);
1714:     }
1715:     if (ghost) {
1716:       PetscFree(rows_v);
1717:     }
1718:   }

1720:   VecAssemblyBegin(vd_nz);
1721:   VecAssemblyBegin(vo_nz);

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

1725:   VecAssemblyEnd(vd_nz);
1726:   VecAssemblyEnd(vo_nz);

1728:   VecGetArray(vd_nz,&vdnz);
1729:   VecGetArray(vo_nz,&vonz);
1730:   for (j=0; j<localSize; j++) {
1731:     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1732:     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1733:   }
1734:   VecRestoreArray(vd_nz,&vdnz);
1735:   VecRestoreArray(vo_nz,&vonz);
1736:   VecDestroy(&vd_nz);
1737:   VecDestroy(&vo_nz);

1739:   MatSeqAIJSetPreallocation(*J,0,dnnz);
1740:   MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);
1741:   MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);

1743:   PetscFree2(dnnz,onnz);

1745:   /* (2) Set matrix entries for edges */
1746:   /*----------------------------------*/
1747:   for (e=eStart; e<eEnd; e++) {
1748:     /* Get row indices */
1749:     DMNetworkGetVariableGlobalOffset(dm,e,&rstart);
1750:     DMNetworkGetNumVariables(dm,e,&nrows);
1751:     if (nrows) {
1752:       for (j=0; j<nrows; j++) rows[j] = j + rstart;

1754:       /* Set matrix entries for conntected vertices */
1755:       DMNetworkGetConnectedVertices(dm,e,&cone);
1756:       for (v=0; v<2; v++) {
1757:         DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);
1758:         DMNetworkGetNumVariables(dm,cone[v],&ncols);

1760:         if (network->Je) {
1761:           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1762:         } else Juser = NULL;
1763:         MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);
1764:       }

1766:       /* Set matrix entries for edge self */
1767:       cstart = rstart;
1768:       if (network->Je) {
1769:         Juser = network->Je[3*e]; /* Jacobian(e,e) */
1770:       } else Juser = NULL;
1771:       MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);
1772:     }
1773:   }

1775:   /* Set matrix entries for vertices */
1776:   /*---------------------------------*/
1777:   for (v=vStart; v<vEnd; v++) {
1778:     /* Get row indices */
1779:     DMNetworkGetVariableGlobalOffset(dm,v,&rstart);
1780:     DMNetworkGetNumVariables(dm,v,&nrows);
1781:     if (!nrows) continue;

1783:     DMNetworkIsGhostVertex(dm,v,&ghost);
1784:     if (ghost) {
1785:       PetscMalloc1(nrows,&rows_v);
1786:     } else {
1787:       rows_v = rows;
1788:     }
1789:     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;

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

1794:     for (e=0; e<nedges; e++) {
1795:       /* Supporting edges */
1796:       DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);
1797:       DMNetworkGetNumVariables(dm,edges[e],&ncols);

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

1804:       /* Connected vertices */
1805:       DMNetworkGetConnectedVertices(dm,edges[e],&cone);
1806:       vc = (v == cone[0]) ? cone[1]:cone[0];

1808:       DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);
1809:       DMNetworkGetNumVariables(dm,vc,&ncols);

1811:       if (network->Jv) {
1812:         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1813:       } else Juser = NULL;
1814:       MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);
1815:     }

1817:     /* Set matrix entries for vertex self */
1818:     if (!ghost) {
1819:       DMNetworkGetVariableGlobalOffset(dm,v,&cstart);
1820:       if (network->Jv) {
1821:         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1822:       } else Juser = NULL;
1823:       MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);
1824:     }
1825:     if (ghost) {
1826:       PetscFree(rows_v);
1827:     }
1828:   }
1829:   PetscFree(rows);

1831:   MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
1832:   MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);

1834:   MatSetDM(*J,dm);
1835:   return(0);
1836: }

1838: PetscErrorCode DMDestroy_Network(DM dm)
1839: {
1841:   DM_Network     *network = (DM_Network*)dm->data;
1842:   PetscInt       j;

1845:   if (--network->refct > 0) return(0);
1846:   if (network->Je) {
1847:     PetscFree(network->Je);
1848:   }
1849:   if (network->Jv) {
1850:     PetscFree(network->Jvptr);
1851:     PetscFree(network->Jv);
1852:   }

1854:   ISLocalToGlobalMappingDestroy(&network->vertex.mapping);
1855:   PetscSectionDestroy(&network->vertex.DofSection);
1856:   PetscSectionDestroy(&network->vertex.GlobalDofSection);
1857:   if (network->vertex.sf) {
1858:     PetscSFDestroy(&network->vertex.sf);
1859:   }
1860:   /* edge */
1861:   ISLocalToGlobalMappingDestroy(&network->edge.mapping);
1862:   PetscSectionDestroy(&network->edge.DofSection);
1863:   PetscSectionDestroy(&network->edge.GlobalDofSection);
1864:   if (network->edge.sf) {
1865:     PetscSFDestroy(&network->edge.sf);
1866:   }
1867:   DMDestroy(&network->plex);
1868:   network->edges = NULL;
1869:   PetscSectionDestroy(&network->DataSection);
1870:   PetscSectionDestroy(&network->DofSection);

1872:   for(j=0; j < network->nsubnet; j++) {
1873:     PetscFree(network->subnet[j].edges);
1874:     PetscFree(network->subnet[j].vertices);
1875:   }
1876:   PetscFree(network->subnet);
1877:   PetscFree(network->componentdataarray);
1878:   PetscFree(network->cvalue);
1879:   PetscFree(network->header);
1880:   PetscFree(network);
1881:   return(0);
1882: }

1884: PetscErrorCode DMView_Network(DM dm,PetscViewer viewer)
1885: {
1887:   DM_Network     *network = (DM_Network*)dm->data;

1890:   DMView(network->plex,viewer);
1891:   return(0);
1892: }

1894: PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
1895: {
1897:   DM_Network     *network = (DM_Network*)dm->data;

1900:   DMGlobalToLocalBegin(network->plex,g,mode,l);
1901:   return(0);
1902: }

1904: PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
1905: {
1907:   DM_Network     *network = (DM_Network*)dm->data;

1910:   DMGlobalToLocalEnd(network->plex,g,mode,l);
1911:   return(0);
1912: }

1914: PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
1915: {
1917:   DM_Network     *network = (DM_Network*)dm->data;

1920:   DMLocalToGlobalBegin(network->plex,l,mode,g);
1921:   return(0);
1922: }

1924: PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
1925: {
1927:   DM_Network     *network = (DM_Network*)dm->data;

1930:   DMLocalToGlobalEnd(network->plex,l,mode,g);
1931:   return(0);
1932: }