xref: /petsc/src/dm/impls/network/network.c (revision b9c6e19d9dc64277a806394b87c516adb6f1d372)
1af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h>  /*I  "petscdmnetwork.h"  I*/
25f2c45f1SShri Abhyankar #include <petscdmplex.h>
35f2c45f1SShri Abhyankar #include <petscsf.h>
45f2c45f1SShri Abhyankar 
55f2c45f1SShri Abhyankar /*@
6556ed216SShri Abhyankar   DMNetworkGetPlex - Gets the Plex DM associated with this network DM
7556ed216SShri Abhyankar 
8556ed216SShri Abhyankar   Not collective
9556ed216SShri Abhyankar 
10556ed216SShri Abhyankar   Input Parameters:
11556ed216SShri Abhyankar + netdm - the dm object
12556ed216SShri Abhyankar - plexmdm - the plex dm object
13556ed216SShri Abhyankar 
14556ed216SShri Abhyankar   Level: Advanced
15556ed216SShri Abhyankar 
16556ed216SShri Abhyankar .seealso: DMNetworkCreate()
17556ed216SShri Abhyankar @*/
18556ed216SShri Abhyankar PetscErrorCode DMNetworkGetPlex(DM netdm, DM *plexdm)
19556ed216SShri Abhyankar {
20556ed216SShri Abhyankar   DM_Network     *network = (DM_Network*) netdm->data;
21556ed216SShri Abhyankar 
22556ed216SShri Abhyankar   PetscFunctionBegin;
23556ed216SShri Abhyankar   *plexdm = network->plex;
24556ed216SShri Abhyankar   PetscFunctionReturn(0);
25556ed216SShri Abhyankar }
26556ed216SShri Abhyankar 
27556ed216SShri Abhyankar /*@
28e2aaf10cSShri Abhyankar   DMNetworkSetSizes - Sets the number of subnetworks,local and global vertices and edges for each subnetwork.
295f2c45f1SShri Abhyankar 
305f2c45f1SShri Abhyankar   Collective on DM
315f2c45f1SShri Abhyankar 
325f2c45f1SShri Abhyankar   Input Parameters:
335f2c45f1SShri Abhyankar + dm - the dm object
34e2aaf10cSShri Abhyankar . Nsubnet - number of subnetworks
35e2aaf10cSShri Abhyankar . nV - number of local vertices for each subnetwork
36e2aaf10cSShri Abhyankar . nE - number of local edges for each subnetwork
37e2aaf10cSShri Abhyankar . NV - number of global vertices (or PETSC_DETERMINE) for each subnetwork
38e2aaf10cSShri Abhyankar - NE - number of global edges (or PETSC_DETERMINE) for each subnetwork
395f2c45f1SShri Abhyankar 
405f2c45f1SShri Abhyankar    Notes
415f2c45f1SShri Abhyankar    If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang.
425f2c45f1SShri Abhyankar 
435f2c45f1SShri Abhyankar    You cannot change the sizes once they have been set
445f2c45f1SShri Abhyankar 
451b266c99SBarry Smith    Level: intermediate
461b266c99SBarry Smith 
471b266c99SBarry Smith .seealso: DMNetworkCreate()
485f2c45f1SShri Abhyankar @*/
49e2aaf10cSShri Abhyankar PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt Nsubnet, PetscInt nV[], PetscInt nE[], PetscInt NV[], PetscInt NE[])
505f2c45f1SShri Abhyankar {
515f2c45f1SShri Abhyankar   PetscErrorCode ierr;
525f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
53e2aaf10cSShri Abhyankar   PetscInt       a[2],b[2],i;
545f2c45f1SShri Abhyankar 
555f2c45f1SShri Abhyankar   PetscFunctionBegin;
565f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
57e2aaf10cSShri Abhyankar   if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet);
58e2aaf10cSShri Abhyankar 
59e2aaf10cSShri Abhyankar   if(Nsubnet > 0) PetscValidLogicalCollectiveInt(dm,Nsubnet,2);
60e2aaf10cSShri Abhyankar   if(network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network");
61e2aaf10cSShri Abhyankar 
62e2aaf10cSShri Abhyankar   network->nsubnet = Nsubnet;
63e2aaf10cSShri Abhyankar   ierr = PetscCalloc1(Nsubnet,&network->subnet);CHKERRQ(ierr);
64e2aaf10cSShri Abhyankar   for(i=0; i < network->nsubnet; i++) {
65e2aaf10cSShri Abhyankar     if (NV[i] > 0) PetscValidLogicalCollectiveInt(dm,NV[i],5);
66e2aaf10cSShri Abhyankar     if (NE[i] > 0) PetscValidLogicalCollectiveInt(dm,NE[i],6);
67e2aaf10cSShri Abhyankar     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]);
68e2aaf10cSShri Abhyankar     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]);
69e2aaf10cSShri Abhyankar     a[0] = nV[i]; a[1] = nE[i];
70b2566f29SBarry Smith     ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
71e2aaf10cSShri Abhyankar     network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1];
72e2aaf10cSShri Abhyankar 
73e2aaf10cSShri Abhyankar     network->subnet[i].id = i;
74e2aaf10cSShri Abhyankar 
75e2aaf10cSShri Abhyankar     network->subnet[i].nvtx = nV[i];
76e2aaf10cSShri Abhyankar     network->subnet[i].vStart = network->nVertices;
77e2aaf10cSShri Abhyankar     network->subnet[i].vEnd   = network->subnet[i].vStart + network->subnet[i].Nvtx;
78e2aaf10cSShri Abhyankar     network->nVertices += network->subnet[i].nvtx;
79e2aaf10cSShri Abhyankar     network->NVertices += network->subnet[i].Nvtx;
80e2aaf10cSShri Abhyankar 
81e2aaf10cSShri Abhyankar     network->subnet[i].nedge = nE[i];
82e2aaf10cSShri Abhyankar     network->subnet[i].eStart = network->nEdges;
83e2aaf10cSShri Abhyankar     network->subnet[i].eEnd = network->subnet[i].eStart + network->subnet[i].Nedge;
84e2aaf10cSShri Abhyankar     network->nEdges += network->subnet[i].nedge;
85e2aaf10cSShri Abhyankar     network->NEdges += network->subnet[i].Nedge;
865f2c45f1SShri Abhyankar   }
875f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
885f2c45f1SShri Abhyankar }
895f2c45f1SShri Abhyankar 
905f2c45f1SShri Abhyankar /*@
915f2c45f1SShri Abhyankar   DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network
925f2c45f1SShri Abhyankar 
935f2c45f1SShri Abhyankar   Logically collective on DM
945f2c45f1SShri Abhyankar 
955f2c45f1SShri Abhyankar   Input Parameters:
96e2aaf10cSShri Abhyankar . edges - list of edges for each subnetwork
975f2c45f1SShri Abhyankar 
985f2c45f1SShri Abhyankar   Notes:
995f2c45f1SShri Abhyankar   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
1005f2c45f1SShri Abhyankar   not be destroyed before the call to DMNetworkLayoutSetUp
1015f2c45f1SShri Abhyankar 
1025f2c45f1SShri Abhyankar   Level: intermediate
1035f2c45f1SShri Abhyankar 
1045f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes
1055f2c45f1SShri Abhyankar @*/
106e2aaf10cSShri Abhyankar PetscErrorCode DMNetworkSetEdgeList(DM dm, int *edgelist[])
1075f2c45f1SShri Abhyankar {
1085f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network*) dm->data;
109e2aaf10cSShri Abhyankar   PetscInt   i;
1105f2c45f1SShri Abhyankar 
1115f2c45f1SShri Abhyankar   PetscFunctionBegin;
112e2aaf10cSShri Abhyankar   for(i=0; i < network->nsubnet; i++) {
113e2aaf10cSShri Abhyankar     network->subnet[i].edgelist = edgelist[i];
114e2aaf10cSShri Abhyankar   }
1155f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
1165f2c45f1SShri Abhyankar }
1175f2c45f1SShri Abhyankar 
1185f2c45f1SShri Abhyankar /*@
1195f2c45f1SShri Abhyankar   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
1205f2c45f1SShri Abhyankar 
1215f2c45f1SShri Abhyankar   Collective on DM
1225f2c45f1SShri Abhyankar 
1235f2c45f1SShri Abhyankar   Input Parameters
1245f2c45f1SShri Abhyankar . DM - the dmnetwork object
1255f2c45f1SShri Abhyankar 
1265f2c45f1SShri Abhyankar   Notes:
1275f2c45f1SShri Abhyankar   This routine should be called after the network sizes and edgelists have been provided. It creates
1285f2c45f1SShri Abhyankar   the bare layout of the network and sets up the network to begin insertion of components.
1295f2c45f1SShri Abhyankar 
1305f2c45f1SShri Abhyankar   All the components should be registered before calling this routine.
1315f2c45f1SShri Abhyankar 
1325f2c45f1SShri Abhyankar   Level: intermediate
1335f2c45f1SShri Abhyankar 
1345f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList
1355f2c45f1SShri Abhyankar @*/
1365f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm)
1375f2c45f1SShri Abhyankar {
1385f2c45f1SShri Abhyankar   PetscErrorCode ierr;
1395f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
1405f2c45f1SShri Abhyankar   PetscInt       dim = 1; /* One dimensional network */
1415f2c45f1SShri Abhyankar   PetscInt       numCorners=2;
1425f2c45f1SShri Abhyankar   PetscInt       spacedim=2;
1435f2c45f1SShri Abhyankar   double         *vertexcoords=NULL;
144e2aaf10cSShri Abhyankar   PetscInt       i,j;
1455f2c45f1SShri Abhyankar   PetscInt       ndata;
146e2aaf10cSShri Abhyankar   PetscInt       ctr=0;
1475f2c45f1SShri Abhyankar 
1485f2c45f1SShri Abhyankar   PetscFunctionBegin;
149e2aaf10cSShri Abhyankar 
1506fefedf4SHong Zhang   if (network->nVertices) {
1516fefedf4SHong Zhang     ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr);
1525f2c45f1SShri Abhyankar   }
153e2aaf10cSShri Abhyankar 
154e2aaf10cSShri Abhyankar   /* Create the edgelist for the network by concatenating edgelists of the subnetworks */
155e2aaf10cSShri Abhyankar   ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr);
156e2aaf10cSShri Abhyankar   for(i=0; i < network->nsubnet; i++) {
157e2aaf10cSShri Abhyankar     for(j = 0; j < network->subnet[i].nedge; j++) {
158e2aaf10cSShri Abhyankar       network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
159e2aaf10cSShri Abhyankar       network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
160e2aaf10cSShri Abhyankar       ctr++;
161e2aaf10cSShri Abhyankar     }
162e2aaf10cSShri Abhyankar   }
163e2aaf10cSShri Abhyankar 
164e2aaf10cSShri Abhyankar #if 0
165e2aaf10cSShri Abhyankar   for(i=0; i < network->nEdges; i++) {
166e2aaf10cSShri Abhyankar     ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr);
167e2aaf10cSShri Abhyankar   }
168e2aaf10cSShri Abhyankar #endif
169e2aaf10cSShri Abhyankar 
1706fefedf4SHong Zhang   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
1716fefedf4SHong Zhang   if (network->nVertices) {
1725f2c45f1SShri Abhyankar     ierr = PetscFree(vertexcoords);CHKERRQ(ierr);
1735f2c45f1SShri Abhyankar   }
174e2aaf10cSShri Abhyankar   ierr = PetscFree(network->edges);CHKERRQ(ierr);
175e2aaf10cSShri Abhyankar 
1765f2c45f1SShri Abhyankar   ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr);
1775f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr);
1785f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr);
1795f2c45f1SShri Abhyankar 
1805f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr);
1815f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr);
1825f2c45f1SShri Abhyankar   ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr);
1835f2c45f1SShri Abhyankar   ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr);
1845f2c45f1SShri Abhyankar 
1852727e31bSShri Abhyankar   /* Create vertices and edges array for the subnetworks */
1862727e31bSShri Abhyankar   for(j=0; j < network->nsubnet; j++) {
1872727e31bSShri Abhyankar     ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr);
1882727e31bSShri Abhyankar     ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr);
1892727e31bSShri Abhyankar     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
1902727e31bSShri Abhyankar        These get updated when the vertices and edges are added. */
1912727e31bSShri Abhyankar     network->subnet[j].nvtx = network->subnet[j].nedge = 0;
1922727e31bSShri Abhyankar   }
1932727e31bSShri Abhyankar 
1945f2c45f1SShri Abhyankar   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
1956caa05f4SBarry Smith   ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr);
196e2aaf10cSShri Abhyankar   for(i=network->eStart; i < network->eEnd; i++) {
197e2aaf10cSShri Abhyankar     network->header[i].index = i;   /* Global edge number */
198e2aaf10cSShri Abhyankar     for(j=0; j < network->nsubnet; j++) {
199e2aaf10cSShri Abhyankar       if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) {
200e2aaf10cSShri Abhyankar 	network->header[i].subnetid = j; /* Subnetwork id */
2012727e31bSShri Abhyankar 	network->subnet[j].edges[network->subnet[j].nedge++] = i;
202e2aaf10cSShri Abhyankar 	break;
203e2aaf10cSShri Abhyankar       }
2047b6afd5bSHong Zhang     }
2055f2c45f1SShri Abhyankar     network->header[i].ndata = 0;
2065f2c45f1SShri Abhyankar     ndata = network->header[i].ndata;
2075f2c45f1SShri Abhyankar     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
2085f2c45f1SShri Abhyankar     network->header[i].offset[ndata] = 0;
2095f2c45f1SShri Abhyankar   }
210e2aaf10cSShri Abhyankar 
211e2aaf10cSShri Abhyankar   for(i=network->vStart; i < network->vEnd; i++) {
212e2aaf10cSShri Abhyankar     network->header[i].index = i - network->vStart;
213e2aaf10cSShri Abhyankar     for(j=0; j < network->nsubnet; j++) {
214e2aaf10cSShri Abhyankar       if((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) {
215e2aaf10cSShri Abhyankar 	network->header[i].subnetid = j;
2162727e31bSShri Abhyankar 	network->subnet[j].vertices[network->subnet[j].nvtx++] = i;
217e2aaf10cSShri Abhyankar 	break;
218e2aaf10cSShri Abhyankar       }
219e2aaf10cSShri Abhyankar     }
220e2aaf10cSShri Abhyankar     network->header[i].ndata = 0;
221e2aaf10cSShri Abhyankar     ndata = network->header[i].ndata;
222e2aaf10cSShri Abhyankar     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
223e2aaf10cSShri Abhyankar     network->header[i].offset[ndata] = 0;
224e2aaf10cSShri Abhyankar   }
225e2aaf10cSShri Abhyankar 
226854ce69bSBarry Smith   ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr);
227e2aaf10cSShri Abhyankar 
2285f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2295f2c45f1SShri Abhyankar }
2305f2c45f1SShri Abhyankar 
23194ef8ddeSSatish Balay /*@C
2322727e31bSShri Abhyankar   DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork
2332727e31bSShri Abhyankar 
2342727e31bSShri Abhyankar   Input Parameters
2352727e31bSShri Abhyankar + dm   - the number object
2362727e31bSShri Abhyankar - id   - the ID (integer) of the subnetwork
2372727e31bSShri Abhyankar 
2382727e31bSShri Abhyankar   Output Parameters
2392727e31bSShri Abhyankar + nv    - number of vertices (local)
2402727e31bSShri Abhyankar . ne    - number of edges (local)
2412727e31bSShri Abhyankar . vtx   - local vertices for this subnetwork
2422727e31bSShri Abhyankar . edge  - local edges for this subnetwork
2432727e31bSShri Abhyankar 
2442727e31bSShri Abhyankar   Notes:
2452727e31bSShri Abhyankar   Cannot call this routine before DMNetworkLayoutSetup()
2462727e31bSShri Abhyankar 
2472727e31bSShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
2482727e31bSShri Abhyankar @*/
2492727e31bSShri Abhyankar PetscErrorCode DMNetworkGetSubnetworkInfo(DM netdm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge)
2502727e31bSShri Abhyankar {
2512727e31bSShri Abhyankar   DM_Network     *network = (DM_Network*) netdm->data;
2522727e31bSShri Abhyankar 
2532727e31bSShri Abhyankar   PetscFunctionBegin;
2542727e31bSShri Abhyankar   *nv = network->subnet[id].nvtx;
2552727e31bSShri Abhyankar   *ne = network->subnet[id].nedge;
2562727e31bSShri Abhyankar   *vtx = network->subnet[id].vertices;
2572727e31bSShri Abhyankar   *edge = network->subnet[id].edges;
2582727e31bSShri Abhyankar   PetscFunctionReturn(0);
2592727e31bSShri Abhyankar }
2602727e31bSShri Abhyankar 
2612727e31bSShri Abhyankar /*@C
2625f2c45f1SShri Abhyankar   DMNetworkRegisterComponent - Registers the network component
2635f2c45f1SShri Abhyankar 
2645f2c45f1SShri Abhyankar   Logically collective on DM
2655f2c45f1SShri Abhyankar 
2665f2c45f1SShri Abhyankar   Input Parameters
2675f2c45f1SShri Abhyankar + dm   - the network object
2685f2c45f1SShri Abhyankar . name - the component name
2695f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data
2705f2c45f1SShri Abhyankar 
2715f2c45f1SShri Abhyankar    Output Parameters
2725f2c45f1SShri Abhyankar .   key - an integer key that defines the component
2735f2c45f1SShri Abhyankar 
2745f2c45f1SShri Abhyankar    Notes
2755f2c45f1SShri Abhyankar    This routine should be called by all processors before calling DMNetworkLayoutSetup().
2765f2c45f1SShri Abhyankar 
2775f2c45f1SShri Abhyankar    Level: intermediate
2785f2c45f1SShri Abhyankar 
2795f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
2805f2c45f1SShri Abhyankar @*/
2815f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key)
2825f2c45f1SShri Abhyankar {
2835f2c45f1SShri Abhyankar   PetscErrorCode        ierr;
2845f2c45f1SShri Abhyankar   DM_Network            *network = (DM_Network*) dm->data;
2855f2c45f1SShri Abhyankar   DMNetworkComponent    *component=&network->component[network->ncomponent];
2865f2c45f1SShri Abhyankar   PetscBool             flg=PETSC_FALSE;
2875f2c45f1SShri Abhyankar   PetscInt              i;
2885f2c45f1SShri Abhyankar 
2895f2c45f1SShri Abhyankar   PetscFunctionBegin;
2905f2c45f1SShri Abhyankar 
2915f2c45f1SShri Abhyankar   for (i=0; i < network->ncomponent; i++) {
2925f2c45f1SShri Abhyankar     ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr);
2935f2c45f1SShri Abhyankar     if (flg) {
2945f2c45f1SShri Abhyankar       *key = i;
2955f2c45f1SShri Abhyankar       PetscFunctionReturn(0);
2965f2c45f1SShri Abhyankar     }
2975f2c45f1SShri Abhyankar   }
2985f2c45f1SShri Abhyankar 
2995f2c45f1SShri Abhyankar   ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr);
3005f2c45f1SShri Abhyankar   component->size = size/sizeof(DMNetworkComponentGenericDataType);
3015f2c45f1SShri Abhyankar   *key = network->ncomponent;
3025f2c45f1SShri Abhyankar   network->ncomponent++;
3035f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3045f2c45f1SShri Abhyankar }
3055f2c45f1SShri Abhyankar 
3065f2c45f1SShri Abhyankar /*@
3075f2c45f1SShri Abhyankar   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.
3085f2c45f1SShri Abhyankar 
3095f2c45f1SShri Abhyankar   Not Collective
3105f2c45f1SShri Abhyankar 
3115f2c45f1SShri Abhyankar   Input Parameters:
3125f2c45f1SShri Abhyankar + dm - The DMNetwork object
3135f2c45f1SShri Abhyankar 
3145f2c45f1SShri Abhyankar   Output Paramters:
3155f2c45f1SShri Abhyankar + vStart - The first vertex point
3165f2c45f1SShri Abhyankar - vEnd   - One beyond the last vertex point
3175f2c45f1SShri Abhyankar 
3185f2c45f1SShri Abhyankar   Level: intermediate
3195f2c45f1SShri Abhyankar 
3205f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange
3215f2c45f1SShri Abhyankar @*/
3225f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
3235f2c45f1SShri Abhyankar {
3245f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
3255f2c45f1SShri Abhyankar 
3265f2c45f1SShri Abhyankar   PetscFunctionBegin;
3275f2c45f1SShri Abhyankar   if (vStart) *vStart = network->vStart;
3285f2c45f1SShri Abhyankar   if (vEnd) *vEnd = network->vEnd;
3295f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3305f2c45f1SShri Abhyankar }
3315f2c45f1SShri Abhyankar 
3325f2c45f1SShri Abhyankar /*@
3335f2c45f1SShri Abhyankar   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.
3345f2c45f1SShri Abhyankar 
3355f2c45f1SShri Abhyankar   Not Collective
3365f2c45f1SShri Abhyankar 
3375f2c45f1SShri Abhyankar   Input Parameters:
3385f2c45f1SShri Abhyankar + dm - The DMNetwork object
3395f2c45f1SShri Abhyankar 
3405f2c45f1SShri Abhyankar   Output Paramters:
3415f2c45f1SShri Abhyankar + eStart - The first edge point
3425f2c45f1SShri Abhyankar - eEnd   - One beyond the last edge point
3435f2c45f1SShri Abhyankar 
3445f2c45f1SShri Abhyankar   Level: intermediate
3455f2c45f1SShri Abhyankar 
3465f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange
3475f2c45f1SShri Abhyankar @*/
3485f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
3495f2c45f1SShri Abhyankar {
3505f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
3515f2c45f1SShri Abhyankar 
3525f2c45f1SShri Abhyankar   PetscFunctionBegin;
3535f2c45f1SShri Abhyankar   if (eStart) *eStart = network->eStart;
3545f2c45f1SShri Abhyankar   if (eEnd) *eEnd = network->eEnd;
3555f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3565f2c45f1SShri Abhyankar }
3575f2c45f1SShri Abhyankar 
3587b6afd5bSHong Zhang /*@
359e85e6aecSHong Zhang   DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge.
3607b6afd5bSHong Zhang 
3617b6afd5bSHong Zhang   Not Collective
3627b6afd5bSHong Zhang 
3637b6afd5bSHong Zhang   Input Parameters:
3647b6afd5bSHong Zhang + dm - DMNetwork object
365e85e6aecSHong Zhang - p  - edge point
3667b6afd5bSHong Zhang 
3677b6afd5bSHong Zhang   Output Paramters:
368e85e6aecSHong Zhang . index - user global numbering for the edge
3697b6afd5bSHong Zhang 
3707b6afd5bSHong Zhang   Level: intermediate
3717b6afd5bSHong Zhang 
372e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalVertexIndex
3737b6afd5bSHong Zhang @*/
374e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
3757b6afd5bSHong Zhang {
3767b6afd5bSHong Zhang   PetscErrorCode    ierr;
3777b6afd5bSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
3787b6afd5bSHong Zhang   PetscInt          offsetp;
3797b6afd5bSHong Zhang   DMNetworkComponentHeader header;
3807b6afd5bSHong Zhang 
3817b6afd5bSHong Zhang   PetscFunctionBegin;
3827b6afd5bSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
3837b6afd5bSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
384e85e6aecSHong Zhang   *index = header->index;
3857b6afd5bSHong Zhang   PetscFunctionReturn(0);
3867b6afd5bSHong Zhang }
3877b6afd5bSHong Zhang 
3885f2c45f1SShri Abhyankar /*@
389e85e6aecSHong Zhang   DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex.
390e85e6aecSHong Zhang 
391e85e6aecSHong Zhang   Not Collective
392e85e6aecSHong Zhang 
393e85e6aecSHong Zhang   Input Parameters:
394e85e6aecSHong Zhang + dm - DMNetwork object
395e85e6aecSHong Zhang - p  - vertex point
396e85e6aecSHong Zhang 
397e85e6aecSHong Zhang   Output Paramters:
398e85e6aecSHong Zhang . index - user global numbering for the vertex
399e85e6aecSHong Zhang 
400e85e6aecSHong Zhang   Level: intermediate
401e85e6aecSHong Zhang 
402e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalEdgeIndex
403e85e6aecSHong Zhang @*/
404e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
405e85e6aecSHong Zhang {
406e85e6aecSHong Zhang   PetscErrorCode    ierr;
407e85e6aecSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
408e85e6aecSHong Zhang   PetscInt          offsetp;
409e85e6aecSHong Zhang   DMNetworkComponentHeader header;
410e85e6aecSHong Zhang 
411e85e6aecSHong Zhang   PetscFunctionBegin;
412e85e6aecSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
413e85e6aecSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
414e85e6aecSHong Zhang   *index = header->index;
415e85e6aecSHong Zhang   PetscFunctionReturn(0);
416e85e6aecSHong Zhang }
417e85e6aecSHong Zhang 
418e85e6aecSHong Zhang /*@
419325661f6SSatish Balay   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)
4205f2c45f1SShri Abhyankar 
4215f2c45f1SShri Abhyankar   Not Collective
4225f2c45f1SShri Abhyankar 
4235f2c45f1SShri Abhyankar   Input Parameters:
4245f2c45f1SShri Abhyankar + dm           - The DMNetwork object
4255f2c45f1SShri Abhyankar . p            - vertex/edge point
4265f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component
4275f2c45f1SShri Abhyankar - compvalue    - pointer to the data structure for the component
4285f2c45f1SShri Abhyankar 
4295f2c45f1SShri Abhyankar   Level: intermediate
4305f2c45f1SShri Abhyankar 
4315f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
4325f2c45f1SShri Abhyankar @*/
4335f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
4345f2c45f1SShri Abhyankar {
4355f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
43643a39a44SBarry Smith   DMNetworkComponent       *component = &network->component[componentkey];
4375f2c45f1SShri Abhyankar   DMNetworkComponentHeader header = &network->header[p];
4385f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue = &network->cvalue[p];
4395f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
4405f2c45f1SShri Abhyankar 
4415f2c45f1SShri Abhyankar   PetscFunctionBegin;
442fa58f0a9SHong Zhang   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);
443fa58f0a9SHong Zhang 
44443a39a44SBarry Smith   header->size[header->ndata] = component->size;
44543a39a44SBarry Smith   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
4465f2c45f1SShri Abhyankar   header->key[header->ndata] = componentkey;
4475f2c45f1SShri Abhyankar   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
4485f2c45f1SShri Abhyankar 
4495f2c45f1SShri Abhyankar   cvalue->data[header->ndata] = (void*)compvalue;
4505f2c45f1SShri Abhyankar   header->ndata++;
4515f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4525f2c45f1SShri Abhyankar }
4535f2c45f1SShri Abhyankar 
4545f2c45f1SShri Abhyankar /*@
4555f2c45f1SShri Abhyankar   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
4565f2c45f1SShri Abhyankar 
4575f2c45f1SShri Abhyankar   Not Collective
4585f2c45f1SShri Abhyankar 
4595f2c45f1SShri Abhyankar   Input Parameters:
4605f2c45f1SShri Abhyankar + dm - The DMNetwork object
4615f2c45f1SShri Abhyankar . p  - vertex/edge point
4625f2c45f1SShri Abhyankar 
4635f2c45f1SShri Abhyankar   Output Parameters:
4645f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge
4655f2c45f1SShri Abhyankar 
4665f2c45f1SShri Abhyankar   Level: intermediate
4675f2c45f1SShri Abhyankar 
4685f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
4695f2c45f1SShri Abhyankar @*/
4705f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
4715f2c45f1SShri Abhyankar {
4725f2c45f1SShri Abhyankar   PetscErrorCode ierr;
4735f2c45f1SShri Abhyankar   PetscInt       offset;
4745f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4755f2c45f1SShri Abhyankar 
4765f2c45f1SShri Abhyankar   PetscFunctionBegin;
4775f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
4785f2c45f1SShri Abhyankar   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
4795f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4805f2c45f1SShri Abhyankar }
4815f2c45f1SShri Abhyankar 
4825f2c45f1SShri Abhyankar /*@
483a730d845SHong Zhang   DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the
4845f2c45f1SShri Abhyankar                                     component value from the component data array
4855f2c45f1SShri Abhyankar 
4865f2c45f1SShri Abhyankar   Not Collective
4875f2c45f1SShri Abhyankar 
4885f2c45f1SShri Abhyankar   Input Parameters:
4895f2c45f1SShri Abhyankar + dm      - The DMNetwork object
4905f2c45f1SShri Abhyankar . p       - vertex/edge point
4915f2c45f1SShri Abhyankar - compnum - component number
4925f2c45f1SShri Abhyankar 
4935f2c45f1SShri Abhyankar   Output Parameters:
4945f2c45f1SShri Abhyankar + compkey - the key obtained when registering the component
4955f2c45f1SShri Abhyankar - offset  - offset into the component data array associated with the vertex/edge point
4965f2c45f1SShri Abhyankar 
4975f2c45f1SShri Abhyankar   Notes:
4985f2c45f1SShri Abhyankar   Typical usage:
4995f2c45f1SShri Abhyankar 
5005f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray(dm, &arr);
5015f2c45f1SShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
5025f2c45f1SShri Abhyankar   Loop over vertices or edges
5035f2c45f1SShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
5045f2c45f1SShri Abhyankar     Loop over numcomps
505a730d845SHong Zhang       DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset);
5065f2c45f1SShri Abhyankar       compdata = (UserCompDataType)(arr+offset);
5075f2c45f1SShri Abhyankar 
5085f2c45f1SShri Abhyankar   Level: intermediate
5095f2c45f1SShri Abhyankar 
5105f2c45f1SShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
5115f2c45f1SShri Abhyankar @*/
512a730d845SHong Zhang PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
5135f2c45f1SShri Abhyankar {
5145f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
5155f2c45f1SShri Abhyankar   PetscInt                 offsetp;
5165f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
5175f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
5185f2c45f1SShri Abhyankar 
5195f2c45f1SShri Abhyankar   PetscFunctionBegin;
5205f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
5215f2c45f1SShri Abhyankar   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
522d36f4e81SHong Zhang   if (compkey) *compkey = header->key[compnum];
523d36f4e81SHong Zhang   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
5245f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5255f2c45f1SShri Abhyankar }
5265f2c45f1SShri Abhyankar 
5275f2c45f1SShri Abhyankar /*@
5285f2c45f1SShri Abhyankar   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
5295f2c45f1SShri Abhyankar 
5305f2c45f1SShri Abhyankar   Not Collective
5315f2c45f1SShri Abhyankar 
5325f2c45f1SShri Abhyankar   Input Parameters:
5335f2c45f1SShri Abhyankar + dm     - The DMNetwork object
5345f2c45f1SShri Abhyankar - p      - the edge/vertex point
5355f2c45f1SShri Abhyankar 
5365f2c45f1SShri Abhyankar   Output Parameters:
5375f2c45f1SShri Abhyankar . offset - the offset
5385f2c45f1SShri Abhyankar 
5395f2c45f1SShri Abhyankar   Level: intermediate
5405f2c45f1SShri Abhyankar 
5415f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
5425f2c45f1SShri Abhyankar @*/
5435f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
5445f2c45f1SShri Abhyankar {
5455f2c45f1SShri Abhyankar   PetscErrorCode ierr;
5465f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5475f2c45f1SShri Abhyankar 
5485f2c45f1SShri Abhyankar   PetscFunctionBegin;
5495f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr);
5505f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5515f2c45f1SShri Abhyankar }
5525f2c45f1SShri Abhyankar 
5535f2c45f1SShri Abhyankar /*@
5545f2c45f1SShri Abhyankar   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
5555f2c45f1SShri Abhyankar 
5565f2c45f1SShri Abhyankar   Not Collective
5575f2c45f1SShri Abhyankar 
5585f2c45f1SShri Abhyankar   Input Parameters:
5595f2c45f1SShri Abhyankar + dm      - The DMNetwork object
5605f2c45f1SShri Abhyankar - p       - the edge/vertex point
5615f2c45f1SShri Abhyankar 
5625f2c45f1SShri Abhyankar   Output Parameters:
5635f2c45f1SShri Abhyankar . offsetg - the offset
5645f2c45f1SShri Abhyankar 
5655f2c45f1SShri Abhyankar   Level: intermediate
5665f2c45f1SShri Abhyankar 
5675f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
5685f2c45f1SShri Abhyankar @*/
5695f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
5705f2c45f1SShri Abhyankar {
5715f2c45f1SShri Abhyankar   PetscErrorCode ierr;
5725f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5735f2c45f1SShri Abhyankar 
5745f2c45f1SShri Abhyankar   PetscFunctionBegin;
5755f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr);
5766fefedf4SHong Zhang   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */
5775f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5785f2c45f1SShri Abhyankar }
5795f2c45f1SShri Abhyankar 
58024121865SAdrian Maldonado /*@
58124121865SAdrian Maldonado   DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector.
58224121865SAdrian Maldonado 
58324121865SAdrian Maldonado   Not Collective
58424121865SAdrian Maldonado 
58524121865SAdrian Maldonado   Input Parameters:
58624121865SAdrian Maldonado + dm     - The DMNetwork object
58724121865SAdrian Maldonado - p      - the edge point
58824121865SAdrian Maldonado 
58924121865SAdrian Maldonado   Output Parameters:
59024121865SAdrian Maldonado . offset - the offset
59124121865SAdrian Maldonado 
59224121865SAdrian Maldonado   Level: intermediate
59324121865SAdrian Maldonado 
59424121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
59524121865SAdrian Maldonado @*/
59624121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
59724121865SAdrian Maldonado {
59824121865SAdrian Maldonado   PetscErrorCode ierr;
59924121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
60024121865SAdrian Maldonado 
60124121865SAdrian Maldonado   PetscFunctionBegin;
60224121865SAdrian Maldonado 
60324121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr);
60424121865SAdrian Maldonado   PetscFunctionReturn(0);
60524121865SAdrian Maldonado }
60624121865SAdrian Maldonado 
60724121865SAdrian Maldonado /*@
60824121865SAdrian Maldonado   DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector.
60924121865SAdrian Maldonado 
61024121865SAdrian Maldonado   Not Collective
61124121865SAdrian Maldonado 
61224121865SAdrian Maldonado   Input Parameters:
61324121865SAdrian Maldonado + dm     - The DMNetwork object
61424121865SAdrian Maldonado - p      - the vertex point
61524121865SAdrian Maldonado 
61624121865SAdrian Maldonado   Output Parameters:
61724121865SAdrian Maldonado . offset - the offset
61824121865SAdrian Maldonado 
61924121865SAdrian Maldonado   Level: intermediate
62024121865SAdrian Maldonado 
62124121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
62224121865SAdrian Maldonado @*/
62324121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
62424121865SAdrian Maldonado {
62524121865SAdrian Maldonado   PetscErrorCode ierr;
62624121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
62724121865SAdrian Maldonado 
62824121865SAdrian Maldonado   PetscFunctionBegin;
62924121865SAdrian Maldonado 
63024121865SAdrian Maldonado   p -= network->vStart;
63124121865SAdrian Maldonado 
63224121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr);
63324121865SAdrian Maldonado   PetscFunctionReturn(0);
63424121865SAdrian Maldonado }
6355f2c45f1SShri Abhyankar /*@
6365f2c45f1SShri Abhyankar   DMNetworkAddNumVariables - Add number of variables associated with a given point.
6375f2c45f1SShri Abhyankar 
6385f2c45f1SShri Abhyankar   Not Collective
6395f2c45f1SShri Abhyankar 
6405f2c45f1SShri Abhyankar   Input Parameters:
6415f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
6425f2c45f1SShri Abhyankar . p    - the vertex/edge point
6435f2c45f1SShri Abhyankar - nvar - number of additional variables
6445f2c45f1SShri Abhyankar 
6455f2c45f1SShri Abhyankar   Level: intermediate
6465f2c45f1SShri Abhyankar 
6475f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables
6485f2c45f1SShri Abhyankar @*/
6495f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
6505f2c45f1SShri Abhyankar {
6515f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6525f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6535f2c45f1SShri Abhyankar 
6545f2c45f1SShri Abhyankar   PetscFunctionBegin;
6555f2c45f1SShri Abhyankar   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
6565f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6575f2c45f1SShri Abhyankar }
6585f2c45f1SShri Abhyankar 
65927f51fceSHong Zhang /*@
66027f51fceSHong Zhang   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
66127f51fceSHong Zhang 
66227f51fceSHong Zhang   Not Collective
66327f51fceSHong Zhang 
66427f51fceSHong Zhang   Input Parameters:
66527f51fceSHong Zhang + dm   - The DMNetworkObject
66627f51fceSHong Zhang - p    - the vertex/edge point
66727f51fceSHong Zhang 
66827f51fceSHong Zhang   Output Parameters:
66927f51fceSHong Zhang . nvar - number of variables
67027f51fceSHong Zhang 
67127f51fceSHong Zhang   Level: intermediate
67227f51fceSHong Zhang 
67327f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
67427f51fceSHong Zhang @*/
67527f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
67627f51fceSHong Zhang {
67727f51fceSHong Zhang   PetscErrorCode ierr;
67827f51fceSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
67927f51fceSHong Zhang 
68027f51fceSHong Zhang   PetscFunctionBegin;
68127f51fceSHong Zhang   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
68227f51fceSHong Zhang   PetscFunctionReturn(0);
68327f51fceSHong Zhang }
68427f51fceSHong Zhang 
6855f2c45f1SShri Abhyankar /*@
6865f2c45f1SShri Abhyankar   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
6875f2c45f1SShri Abhyankar 
6885f2c45f1SShri Abhyankar   Not Collective
6895f2c45f1SShri Abhyankar 
6905f2c45f1SShri Abhyankar   Input Parameters:
6915f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
6925f2c45f1SShri Abhyankar . p    - the vertex/edge point
6935f2c45f1SShri Abhyankar - nvar - number of variables
6945f2c45f1SShri Abhyankar 
6955f2c45f1SShri Abhyankar   Level: intermediate
6965f2c45f1SShri Abhyankar 
6975f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables
6985f2c45f1SShri Abhyankar @*/
6995f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
7005f2c45f1SShri Abhyankar {
7015f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7025f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7035f2c45f1SShri Abhyankar 
7045f2c45f1SShri Abhyankar   PetscFunctionBegin;
7055f2c45f1SShri Abhyankar   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
7065f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7075f2c45f1SShri Abhyankar }
7085f2c45f1SShri Abhyankar 
7095f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This
7105f2c45f1SShri Abhyankar    function is called during DMSetUp() */
7115f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm)
7125f2c45f1SShri Abhyankar {
7135f2c45f1SShri Abhyankar   PetscErrorCode              ierr;
7145f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7155f2c45f1SShri Abhyankar   PetscInt                    arr_size;
7165f2c45f1SShri Abhyankar   PetscInt                    p,offset,offsetp;
7175f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
7185f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue;
7195f2c45f1SShri Abhyankar   DMNetworkComponentGenericDataType      *componentdataarray;
7205f2c45f1SShri Abhyankar   PetscInt ncomp, i;
7215f2c45f1SShri Abhyankar 
7225f2c45f1SShri Abhyankar   PetscFunctionBegin;
7235f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
7245f2c45f1SShri Abhyankar   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
72575b160a0SShri Abhyankar   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
7265f2c45f1SShri Abhyankar   componentdataarray = network->componentdataarray;
7275f2c45f1SShri Abhyankar   for (p = network->pStart; p < network->pEnd; p++) {
7285f2c45f1SShri Abhyankar     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
7295f2c45f1SShri Abhyankar     /* Copy header */
7305f2c45f1SShri Abhyankar     header = &network->header[p];
731302440fdSBarry Smith     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
7325f2c45f1SShri Abhyankar     /* Copy data */
7335f2c45f1SShri Abhyankar     cvalue = &network->cvalue[p];
7345f2c45f1SShri Abhyankar     ncomp = header->ndata;
7355f2c45f1SShri Abhyankar     for (i = 0; i < ncomp; i++) {
7365f2c45f1SShri Abhyankar       offset = offsetp + network->dataheadersize + header->offset[i];
737302440fdSBarry Smith       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
7385f2c45f1SShri Abhyankar     }
7395f2c45f1SShri Abhyankar   }
7405f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7415f2c45f1SShri Abhyankar }
7425f2c45f1SShri Abhyankar 
7435f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */
7445f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm)
7455f2c45f1SShri Abhyankar {
7465f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7475f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7485f2c45f1SShri Abhyankar 
7495f2c45f1SShri Abhyankar   PetscFunctionBegin;
7505f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
7515f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7525f2c45f1SShri Abhyankar }
7535f2c45f1SShri Abhyankar 
7545f2c45f1SShri Abhyankar /*@C
7555f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray - Returns the component data array
7565f2c45f1SShri Abhyankar 
7575f2c45f1SShri Abhyankar   Not Collective
7585f2c45f1SShri Abhyankar 
7595f2c45f1SShri Abhyankar   Input Parameters:
7605f2c45f1SShri Abhyankar . dm - The DMNetwork Object
7615f2c45f1SShri Abhyankar 
7625f2c45f1SShri Abhyankar   Output Parameters:
7635f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components
7645f2c45f1SShri Abhyankar 
7655f2c45f1SShri Abhyankar   Level: intermediate
7665f2c45f1SShri Abhyankar 
767a730d845SHong Zhang .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents
7685f2c45f1SShri Abhyankar @*/
7695f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
7705f2c45f1SShri Abhyankar {
7715f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7725f2c45f1SShri Abhyankar 
7735f2c45f1SShri Abhyankar   PetscFunctionBegin;
7745f2c45f1SShri Abhyankar   *componentdataarray = network->componentdataarray;
7755f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7765f2c45f1SShri Abhyankar }
7775f2c45f1SShri Abhyankar 
77824121865SAdrian Maldonado /* Get a subsection from a range of points */
77924121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection)
78024121865SAdrian Maldonado {
78124121865SAdrian Maldonado   PetscErrorCode ierr;
78224121865SAdrian Maldonado   PetscInt       i, nvar;
78324121865SAdrian Maldonado 
78424121865SAdrian Maldonado   PetscFunctionBegin;
78524121865SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr);
78624121865SAdrian Maldonado   ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr);
78724121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
78824121865SAdrian Maldonado     ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr);
78924121865SAdrian Maldonado     ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr);
79024121865SAdrian Maldonado   }
79124121865SAdrian Maldonado 
79224121865SAdrian Maldonado   ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr);
79324121865SAdrian Maldonado   PetscFunctionReturn(0);
79424121865SAdrian Maldonado }
79524121865SAdrian Maldonado 
79624121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */
79724121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
79824121865SAdrian Maldonado {
79924121865SAdrian Maldonado   PetscErrorCode ierr;
80024121865SAdrian Maldonado   PetscInt       i, *subpoints;
80124121865SAdrian Maldonado 
80224121865SAdrian Maldonado   PetscFunctionBegin;
80324121865SAdrian Maldonado   /* Create index sets to map from "points" to "subpoints" */
80424121865SAdrian Maldonado   ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr);
80524121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
80624121865SAdrian Maldonado     subpoints[i - pstart] = i;
80724121865SAdrian Maldonado   }
808459726d8SSatish Balay   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr);
80924121865SAdrian Maldonado   ierr = PetscFree(subpoints);CHKERRQ(ierr);
81024121865SAdrian Maldonado   PetscFunctionReturn(0);
81124121865SAdrian Maldonado }
81224121865SAdrian Maldonado 
81324121865SAdrian Maldonado /*@
81424121865SAdrian Maldonado   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute.
81524121865SAdrian Maldonado 
81624121865SAdrian Maldonado   Collective
81724121865SAdrian Maldonado 
81824121865SAdrian Maldonado   Input Parameters:
81924121865SAdrian Maldonado . dm   - The DMNetworkObject
82024121865SAdrian Maldonado 
82124121865SAdrian Maldonado   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
82224121865SAdrian Maldonado 
82324121865SAdrian Maldonado   points = [0 1 2 3 4 5 6]
82424121865SAdrian Maldonado 
82524121865SAdrian Maldonado   where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]).
82624121865SAdrian Maldonado 
82724121865SAdrian Maldonado   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
82824121865SAdrian Maldonado 
82924121865SAdrian Maldonado   Level: intermediate
83024121865SAdrian Maldonado 
83124121865SAdrian Maldonado @*/
83224121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
83324121865SAdrian Maldonado {
83424121865SAdrian Maldonado   PetscErrorCode ierr;
83524121865SAdrian Maldonado   MPI_Comm       comm;
8369852e123SBarry Smith   PetscMPIInt    rank, size;
83724121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
83824121865SAdrian Maldonado 
839eab1376dSHong Zhang   PetscFunctionBegin;
84024121865SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
84124121865SAdrian Maldonado   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
8429852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
84324121865SAdrian Maldonado 
84424121865SAdrian Maldonado   /* Create maps for vertices and edges */
84524121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr);
84624121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr);
84724121865SAdrian Maldonado 
84824121865SAdrian Maldonado   /* Create local sub-sections */
84924121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr);
85024121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr);
85124121865SAdrian Maldonado 
8529852e123SBarry Smith   if (size > 1) {
85324121865SAdrian Maldonado     ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr);
85424121865SAdrian Maldonado     ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr);
85524121865SAdrian Maldonado   ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr);
85624121865SAdrian Maldonado   ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr);
85724121865SAdrian Maldonado   } else {
85824121865SAdrian Maldonado   /* create structures for vertex */
85924121865SAdrian Maldonado   ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr);
86024121865SAdrian Maldonado   /* create structures for edge */
86124121865SAdrian Maldonado   ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr);
86224121865SAdrian Maldonado   }
86324121865SAdrian Maldonado 
86424121865SAdrian Maldonado 
86524121865SAdrian Maldonado   /* Add viewers */
86624121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr);
86724121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr);
86824121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr);
86924121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr);
87024121865SAdrian Maldonado 
87124121865SAdrian Maldonado   PetscFunctionReturn(0);
87224121865SAdrian Maldonado }
8737b6afd5bSHong Zhang 
8745f2c45f1SShri Abhyankar /*@
8755f2c45f1SShri Abhyankar   DMNetworkDistribute - Distributes the network and moves associated component data.
8765f2c45f1SShri Abhyankar 
8775f2c45f1SShri Abhyankar   Collective
8785f2c45f1SShri Abhyankar 
8795f2c45f1SShri Abhyankar   Input Parameter:
880d3464fd4SAdrian Maldonado + DM - the DMNetwork object
8815f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default
8825f2c45f1SShri Abhyankar 
8835f2c45f1SShri Abhyankar   Notes:
8848b171c8eSHong Zhang   Distributes the network with <overlap>-overlapping partitioning of the edges.
8855f2c45f1SShri Abhyankar 
8865f2c45f1SShri Abhyankar   Level: intermediate
8875f2c45f1SShri Abhyankar 
8885f2c45f1SShri Abhyankar .seealso: DMNetworkCreate
8895f2c45f1SShri Abhyankar @*/
890d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
8915f2c45f1SShri Abhyankar {
892d3464fd4SAdrian Maldonado   MPI_Comm       comm;
8935f2c45f1SShri Abhyankar   PetscErrorCode ierr;
894d3464fd4SAdrian Maldonado   PetscMPIInt    size;
895d3464fd4SAdrian Maldonado   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
896d3464fd4SAdrian Maldonado   DM_Network     *newDMnetwork;
8975f2c45f1SShri Abhyankar   PetscSF        pointsf;
8985f2c45f1SShri Abhyankar   DM             newDM;
89951ac5effSHong Zhang   PetscPartitioner part;
900*b9c6e19dSShri Abhyankar   PetscInt         j,e,v,offset;
901*b9c6e19dSShri Abhyankar   DMNetworkComponentHeader header;
9025f2c45f1SShri Abhyankar 
9035f2c45f1SShri Abhyankar   PetscFunctionBegin;
904d3464fd4SAdrian Maldonado 
905d3464fd4SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr);
906d3464fd4SAdrian Maldonado   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
907d3464fd4SAdrian Maldonado   if (size == 1) PetscFunctionReturn(0);
908d3464fd4SAdrian Maldonado 
909d3464fd4SAdrian Maldonado   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr);
9105f2c45f1SShri Abhyankar   newDMnetwork = (DM_Network*)newDM->data;
9115f2c45f1SShri Abhyankar   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
91251ac5effSHong Zhang 
91351ac5effSHong Zhang   /* Enable runtime options for petscpartitioner */
91451ac5effSHong Zhang   ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr);
91551ac5effSHong Zhang   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
91651ac5effSHong Zhang 
9175f2c45f1SShri Abhyankar   /* Distribute plex dm and dof section */
91880cf41d5SMatthew G. Knepley   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
91951ac5effSHong Zhang 
9205f2c45f1SShri Abhyankar   /* Distribute dof section */
921d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr);
9225f2c45f1SShri Abhyankar   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
923d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr);
92451ac5effSHong Zhang 
9255f2c45f1SShri Abhyankar   /* Distribute data and associated section */
92631da1fc8SHong Zhang   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
92724121865SAdrian Maldonado 
9285f2c45f1SShri Abhyankar   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
9295f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
9305f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
9315f2c45f1SShri Abhyankar   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
9326fefedf4SHong Zhang   newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
9336fefedf4SHong Zhang   newDMnetwork->NVertices = oldDMnetwork->NVertices;
9345f2c45f1SShri Abhyankar   newDMnetwork->NEdges = oldDMnetwork->NEdges;
93524121865SAdrian Maldonado 
9365f2c45f1SShri Abhyankar   /* Set Dof section as the default section for dm */
9375f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
9385f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
9395f2c45f1SShri Abhyankar 
940*b9c6e19dSShri Abhyankar   /* Set up subnetwork info in the newDM */
941*b9c6e19dSShri Abhyankar   newDMnetwork->nsubnet = oldDMnetwork->nsubnet;
942*b9c6e19dSShri Abhyankar   ierr = PetscCalloc1(newDMnetwork->nsubnet,&newDMnetwork->subnet);CHKERRQ(ierr);
943*b9c6e19dSShri Abhyankar   /* Copy over the global number of vertices and edges in each subnetwork. Note that these are already
944*b9c6e19dSShri Abhyankar      calculated in DMNetworkLayoutSetUp()
945*b9c6e19dSShri Abhyankar   */
946*b9c6e19dSShri Abhyankar   for(j=0; j < newDMnetwork->nsubnet; j++) {
947*b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx;
948*b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge;
949*b9c6e19dSShri Abhyankar   }
950*b9c6e19dSShri Abhyankar 
951*b9c6e19dSShri Abhyankar   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
952*b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
953*b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
954*b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].nedge++;
955*b9c6e19dSShri Abhyankar   }
956*b9c6e19dSShri Abhyankar 
957*b9c6e19dSShri Abhyankar   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
958*b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
959*b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
960*b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].nvtx++;
961*b9c6e19dSShri Abhyankar   }
962*b9c6e19dSShri Abhyankar 
963*b9c6e19dSShri Abhyankar   /* Now create the vertices and edge arrays for the subnetworks */
964*b9c6e19dSShri Abhyankar   for(j=0; j < newDMnetwork->nsubnet; j++) {
965*b9c6e19dSShri Abhyankar     ierr = PetscCalloc1(newDMnetwork->subnet[j].nedge,&newDMnetwork->subnet[j].edges);CHKERRQ(ierr);
966*b9c6e19dSShri Abhyankar     ierr = PetscCalloc1(newDMnetwork->subnet[j].nvtx,&newDMnetwork->subnet[j].vertices);CHKERRQ(ierr);
967*b9c6e19dSShri Abhyankar     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
968*b9c6e19dSShri Abhyankar        These get updated when the vertices and edges are added. */
969*b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0;
970*b9c6e19dSShri Abhyankar   }
971*b9c6e19dSShri Abhyankar 
972*b9c6e19dSShri Abhyankar   /* Set the vertices and edges in each subnetwork */
973*b9c6e19dSShri Abhyankar   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
974*b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
975*b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
976*b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++]  = e;
977*b9c6e19dSShri Abhyankar   }
978*b9c6e19dSShri Abhyankar 
979*b9c6e19dSShri Abhyankar   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
980*b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
981*b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
982*b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++]  = v;
983*b9c6e19dSShri Abhyankar   }
984*b9c6e19dSShri Abhyankar 
98524121865SAdrian Maldonado   /* Destroy point SF */
98624121865SAdrian Maldonado   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
98724121865SAdrian Maldonado 
988d3464fd4SAdrian Maldonado   ierr = DMDestroy(dm);CHKERRQ(ierr);
989d3464fd4SAdrian Maldonado   *dm  = newDM;
9905f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9915f2c45f1SShri Abhyankar }
9925f2c45f1SShri Abhyankar 
99324121865SAdrian Maldonado /*@C
99424121865SAdrian Maldonado   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering.
99524121865SAdrian Maldonado 
99624121865SAdrian Maldonado   Input Parameters:
99724121865SAdrian Maldonado + masterSF - the original SF structure
99824121865SAdrian Maldonado - map      - a ISLocalToGlobal mapping that contains the subset of points
99924121865SAdrian Maldonado 
100024121865SAdrian Maldonado   Output Parameters:
100124121865SAdrian Maldonado . subSF    - a subset of the masterSF for the desired subset.
100224121865SAdrian Maldonado */
100324121865SAdrian Maldonado 
100424121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) {
100524121865SAdrian Maldonado 
100624121865SAdrian Maldonado   PetscErrorCode        ierr;
100724121865SAdrian Maldonado   PetscInt              nroots, nleaves, *ilocal_sub;
100824121865SAdrian Maldonado   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
100924121865SAdrian Maldonado   PetscInt              *local_points, *remote_points;
101024121865SAdrian Maldonado   PetscSFNode           *iremote_sub;
101124121865SAdrian Maldonado   const PetscInt        *ilocal;
101224121865SAdrian Maldonado   const PetscSFNode     *iremote;
101324121865SAdrian Maldonado 
101424121865SAdrian Maldonado   PetscFunctionBegin;
101524121865SAdrian Maldonado   ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
101624121865SAdrian Maldonado 
101724121865SAdrian Maldonado   /* Look for leaves that pertain to the subset of points. Get the local ordering */
101824121865SAdrian Maldonado   ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr);
101924121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr);
102024121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
102124121865SAdrian Maldonado     if (ilocal_map[i] != -1) nleaves_sub += 1;
102224121865SAdrian Maldonado   }
102324121865SAdrian Maldonado   /* Re-number ilocal with subset numbering. Need information from roots */
102424121865SAdrian Maldonado   ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr);
102524121865SAdrian Maldonado   for (i = 0; i < nroots; i++) local_points[i] = i;
102624121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr);
102724121865SAdrian Maldonado   ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
102824121865SAdrian Maldonado   ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
102924121865SAdrian Maldonado   /* Fill up graph using local (that is, local to the subset) numbering. */
10304b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr);
10314b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr);
103224121865SAdrian Maldonado   nleaves_sub = 0;
103324121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
103424121865SAdrian Maldonado     if (ilocal_map[i] != -1) {
103524121865SAdrian Maldonado       ilocal_sub[nleaves_sub] = ilocal_map[i];
10364b70a8deSAdrian Maldonado       iremote_sub[nleaves_sub].rank = iremote[i].rank;
103724121865SAdrian Maldonado       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
103824121865SAdrian Maldonado       nleaves_sub += 1;
103924121865SAdrian Maldonado     }
104024121865SAdrian Maldonado   }
104124121865SAdrian Maldonado   ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr);
104224121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr);
104324121865SAdrian Maldonado 
104424121865SAdrian Maldonado   /* Create new subSF */
104524121865SAdrian Maldonado   ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr);
104624121865SAdrian Maldonado   ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr);
10474b70a8deSAdrian Maldonado   ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr);
104824121865SAdrian Maldonado   ierr = PetscFree(ilocal_map);CHKERRQ(ierr);
10494b70a8deSAdrian Maldonado   ierr = PetscFree(iremote_sub);CHKERRQ(ierr);
105024121865SAdrian Maldonado   PetscFunctionReturn(0);
105124121865SAdrian Maldonado }
105224121865SAdrian Maldonado 
10535f2c45f1SShri Abhyankar /*@C
10545f2c45f1SShri Abhyankar   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
10555f2c45f1SShri Abhyankar 
10565f2c45f1SShri Abhyankar   Not Collective
10575f2c45f1SShri Abhyankar 
10585f2c45f1SShri Abhyankar   Input Parameters:
10595f2c45f1SShri Abhyankar + dm - The DMNetwork object
10605f2c45f1SShri Abhyankar - p  - the vertex point
10615f2c45f1SShri Abhyankar 
10625f2c45f1SShri Abhyankar   Output Paramters:
10635f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point
10645f2c45f1SShri Abhyankar - edges  - List of edge points
10655f2c45f1SShri Abhyankar 
10665f2c45f1SShri Abhyankar   Level: intermediate
10675f2c45f1SShri Abhyankar 
10685f2c45f1SShri Abhyankar   Fortran Notes:
10695f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
10705f2c45f1SShri Abhyankar   include petsc.h90 in your code.
10715f2c45f1SShri Abhyankar 
1072d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices
10735f2c45f1SShri Abhyankar @*/
10745f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
10755f2c45f1SShri Abhyankar {
10765f2c45f1SShri Abhyankar   PetscErrorCode ierr;
10775f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
10785f2c45f1SShri Abhyankar 
10795f2c45f1SShri Abhyankar   PetscFunctionBegin;
10805f2c45f1SShri Abhyankar   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
10815f2c45f1SShri Abhyankar   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
10825f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10835f2c45f1SShri Abhyankar }
10845f2c45f1SShri Abhyankar 
10855f2c45f1SShri Abhyankar /*@C
1086d842c372SHong Zhang   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
10875f2c45f1SShri Abhyankar 
10885f2c45f1SShri Abhyankar   Not Collective
10895f2c45f1SShri Abhyankar 
10905f2c45f1SShri Abhyankar   Input Parameters:
10915f2c45f1SShri Abhyankar + dm - The DMNetwork object
10925f2c45f1SShri Abhyankar - p  - the edge point
10935f2c45f1SShri Abhyankar 
10945f2c45f1SShri Abhyankar   Output Paramters:
10955f2c45f1SShri Abhyankar . vertices  - vertices connected to this edge
10965f2c45f1SShri Abhyankar 
10975f2c45f1SShri Abhyankar   Level: intermediate
10985f2c45f1SShri Abhyankar 
10995f2c45f1SShri Abhyankar   Fortran Notes:
11005f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
11015f2c45f1SShri Abhyankar   include petsc.h90 in your code.
11025f2c45f1SShri Abhyankar 
11035f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
11045f2c45f1SShri Abhyankar @*/
1105d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
11065f2c45f1SShri Abhyankar {
11075f2c45f1SShri Abhyankar   PetscErrorCode ierr;
11085f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
11095f2c45f1SShri Abhyankar 
11105f2c45f1SShri Abhyankar   PetscFunctionBegin;
11115f2c45f1SShri Abhyankar   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
11125f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
11135f2c45f1SShri Abhyankar }
11145f2c45f1SShri Abhyankar 
11155f2c45f1SShri Abhyankar /*@
11165f2c45f1SShri Abhyankar   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
11175f2c45f1SShri Abhyankar 
11185f2c45f1SShri Abhyankar   Not Collective
11195f2c45f1SShri Abhyankar 
11205f2c45f1SShri Abhyankar   Input Parameters:
11215f2c45f1SShri Abhyankar + dm - The DMNetwork object
11225f2c45f1SShri Abhyankar . p  - the vertex point
11235f2c45f1SShri Abhyankar 
11245f2c45f1SShri Abhyankar   Output Parameter:
11255f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point
11265f2c45f1SShri Abhyankar 
11275f2c45f1SShri Abhyankar   Level: intermediate
11285f2c45f1SShri Abhyankar 
1129d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange
11305f2c45f1SShri Abhyankar @*/
11315f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
11325f2c45f1SShri Abhyankar {
11335f2c45f1SShri Abhyankar   PetscErrorCode ierr;
11345f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
11355f2c45f1SShri Abhyankar   PetscInt       offsetg;
11365f2c45f1SShri Abhyankar   PetscSection   sectiong;
11375f2c45f1SShri Abhyankar 
11385f2c45f1SShri Abhyankar   PetscFunctionBegin;
11395f2c45f1SShri Abhyankar   *isghost = PETSC_FALSE;
11405f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
11415f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
11425f2c45f1SShri Abhyankar   if (offsetg < 0) *isghost = PETSC_TRUE;
11435f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
11445f2c45f1SShri Abhyankar }
11455f2c45f1SShri Abhyankar 
11465f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm)
11475f2c45f1SShri Abhyankar {
11485f2c45f1SShri Abhyankar   PetscErrorCode ierr;
11495f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
11505f2c45f1SShri Abhyankar 
11515f2c45f1SShri Abhyankar   PetscFunctionBegin;
11525f2c45f1SShri Abhyankar   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
11535f2c45f1SShri Abhyankar   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
11545f2c45f1SShri Abhyankar 
11555f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
11565f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
11575f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
11585f2c45f1SShri Abhyankar }
11595f2c45f1SShri Abhyankar 
11601ad426b7SHong Zhang /*@
116117df6e9eSHong Zhang     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
11621ad426b7SHong Zhang                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
11631ad426b7SHong Zhang 
11641ad426b7SHong Zhang     Collective
11651ad426b7SHong Zhang 
11661ad426b7SHong Zhang     Input Parameters:
116783b2e829SHong Zhang +   dm - The DMNetwork object
116883b2e829SHong Zhang .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
116983b2e829SHong Zhang -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
11701ad426b7SHong Zhang 
11711ad426b7SHong Zhang     Level: intermediate
11721ad426b7SHong Zhang 
11731ad426b7SHong Zhang @*/
117483b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
11751ad426b7SHong Zhang {
11761ad426b7SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
11771ad426b7SHong Zhang 
11781ad426b7SHong Zhang   PetscFunctionBegin;
117983b2e829SHong Zhang   network->userEdgeJacobian   = eflg;
118083b2e829SHong Zhang   network->userVertexJacobian = vflg;
11811ad426b7SHong Zhang   PetscFunctionReturn(0);
11821ad426b7SHong Zhang }
11831ad426b7SHong Zhang 
11841ad426b7SHong Zhang /*@
118583b2e829SHong Zhang     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
118683b2e829SHong Zhang 
118783b2e829SHong Zhang     Not Collective
118883b2e829SHong Zhang 
118983b2e829SHong Zhang     Input Parameters:
119083b2e829SHong Zhang +   dm - The DMNetwork object
119183b2e829SHong Zhang .   p  - the edge point
11923e97b6e8SHong Zhang -   J - array (size = 3) of Jacobian submatrices for this edge point:
11933e97b6e8SHong Zhang         J[0]: this edge
1194d842c372SHong Zhang         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()
119583b2e829SHong Zhang 
119683b2e829SHong Zhang     Level: intermediate
119783b2e829SHong Zhang 
119883b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix
119983b2e829SHong Zhang @*/
120083b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
120183b2e829SHong Zhang {
120283b2e829SHong Zhang   PetscErrorCode ierr;
120383b2e829SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
120483b2e829SHong Zhang 
120583b2e829SHong Zhang   PetscFunctionBegin;
1206883e35e8SHong Zhang   if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
120783b2e829SHong Zhang   if (!network->Je) {
1208883e35e8SHong Zhang     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
120983b2e829SHong Zhang   }
1210883e35e8SHong Zhang   network->Je[3*p]   = J[0];
1211883e35e8SHong Zhang   network->Je[3*p+1] = J[1];
1212883e35e8SHong Zhang   network->Je[3*p+2] = J[2];
121383b2e829SHong Zhang   PetscFunctionReturn(0);
121483b2e829SHong Zhang }
121583b2e829SHong Zhang 
121683b2e829SHong Zhang /*@
121776ddfea5SHong Zhang     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
12181ad426b7SHong Zhang 
12191ad426b7SHong Zhang     Not Collective
12201ad426b7SHong Zhang 
12211ad426b7SHong Zhang     Input Parameters:
12221ad426b7SHong Zhang +   dm - The DMNetwork object
12231ad426b7SHong Zhang .   p  - the vertex point
12243e97b6e8SHong Zhang -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
12253e97b6e8SHong Zhang         J[0]:       this vertex
12263e97b6e8SHong Zhang         J[1+2*i]:   i-th supporting edge
12273e97b6e8SHong Zhang         J[1+2*i+1]: i-th connected vertex
12281ad426b7SHong Zhang 
12291ad426b7SHong Zhang     Level: intermediate
12301ad426b7SHong Zhang 
123183b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix
12321ad426b7SHong Zhang @*/
1233883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
12345f2c45f1SShri Abhyankar {
12355f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12365f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
1237f4431b8cSHong Zhang   PetscInt       i,*vptr,nedges,vStart,vEnd;
1238883e35e8SHong Zhang   const PetscInt *edges;
12395f2c45f1SShri Abhyankar 
12405f2c45f1SShri Abhyankar   PetscFunctionBegin;
1241883e35e8SHong Zhang   if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
1242883e35e8SHong Zhang 
1243883e35e8SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1244f4431b8cSHong Zhang 
124583b2e829SHong Zhang   if (!network->Jv) {
12466fefedf4SHong Zhang     PetscInt nVertices = network->nVertices,nedges_total;
1247883e35e8SHong Zhang     /* count nvertex_total */
12483e97b6e8SHong Zhang     nedges_total = 0;
1249883e35e8SHong Zhang     ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
12506fefedf4SHong Zhang     ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr);
1251f4431b8cSHong Zhang 
1252883e35e8SHong Zhang     vptr[0] = 0;
12536fefedf4SHong Zhang     for (i=0; i<nVertices; i++) {
1254f4431b8cSHong Zhang       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
1255883e35e8SHong Zhang       nedges_total += nedges;
1256f4431b8cSHong Zhang       vptr[i+1] = vptr[i] + 2*nedges + 1;
12571ad426b7SHong Zhang     }
12583e97b6e8SHong Zhang 
12596fefedf4SHong Zhang     ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr);
1260883e35e8SHong Zhang     network->Jvptr = vptr;
1261883e35e8SHong Zhang   }
1262883e35e8SHong Zhang 
1263883e35e8SHong Zhang   vptr = network->Jvptr;
12643e97b6e8SHong Zhang   network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
12653e97b6e8SHong Zhang 
12663e97b6e8SHong Zhang   /* Set Jacobian for each supporting edge and connected vertex */
1267883e35e8SHong Zhang   ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
1268883e35e8SHong Zhang   for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
1269883e35e8SHong Zhang   PetscFunctionReturn(0);
1270883e35e8SHong Zhang }
1271883e35e8SHong Zhang 
1272e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
12735cf7da58SHong Zhang {
12745cf7da58SHong Zhang   PetscErrorCode ierr;
12755cf7da58SHong Zhang   PetscInt       j;
12765cf7da58SHong Zhang   PetscScalar    val=(PetscScalar)ncols;
12775cf7da58SHong Zhang 
12785cf7da58SHong Zhang   PetscFunctionBegin;
12795cf7da58SHong Zhang   if (!ghost) {
12805cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
12815cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
12825cf7da58SHong Zhang     }
12835cf7da58SHong Zhang   } else {
12845cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
12855cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
12865cf7da58SHong Zhang     }
12875cf7da58SHong Zhang   }
12885cf7da58SHong Zhang   PetscFunctionReturn(0);
12895cf7da58SHong Zhang }
12905cf7da58SHong Zhang 
1291e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
12925cf7da58SHong Zhang {
12935cf7da58SHong Zhang   PetscErrorCode ierr;
12945cf7da58SHong Zhang   PetscInt       j,ncols_u;
12955cf7da58SHong Zhang   PetscScalar    val;
12965cf7da58SHong Zhang 
12975cf7da58SHong Zhang   PetscFunctionBegin;
12985cf7da58SHong Zhang   if (!ghost) {
12995cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
13005cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
13015cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
13025cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
13035cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
13045cf7da58SHong Zhang     }
13055cf7da58SHong Zhang   } else {
13065cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
13075cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
13085cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
13095cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
13105cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
13115cf7da58SHong Zhang     }
13125cf7da58SHong Zhang   }
13135cf7da58SHong Zhang   PetscFunctionReturn(0);
13145cf7da58SHong Zhang }
13155cf7da58SHong Zhang 
1316e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
13175cf7da58SHong Zhang {
13185cf7da58SHong Zhang   PetscErrorCode ierr;
13195cf7da58SHong Zhang 
13205cf7da58SHong Zhang   PetscFunctionBegin;
13215cf7da58SHong Zhang   if (Ju) {
13225cf7da58SHong Zhang     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
13235cf7da58SHong Zhang   } else {
13245cf7da58SHong Zhang     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
13255cf7da58SHong Zhang   }
13265cf7da58SHong Zhang   PetscFunctionReturn(0);
13275cf7da58SHong Zhang }
13285cf7da58SHong Zhang 
1329e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1330883e35e8SHong Zhang {
1331883e35e8SHong Zhang   PetscErrorCode ierr;
1332883e35e8SHong Zhang   PetscInt       j,*cols;
1333883e35e8SHong Zhang   PetscScalar    *zeros;
1334883e35e8SHong Zhang 
1335883e35e8SHong Zhang   PetscFunctionBegin;
1336883e35e8SHong Zhang   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
1337883e35e8SHong Zhang   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1338883e35e8SHong Zhang   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
1339883e35e8SHong Zhang   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
13401ad426b7SHong Zhang   PetscFunctionReturn(0);
13411ad426b7SHong Zhang }
1342a4e85ca8SHong Zhang 
1343e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
13443e97b6e8SHong Zhang {
13453e97b6e8SHong Zhang   PetscErrorCode ierr;
13463e97b6e8SHong Zhang   PetscInt       j,M,N,row,col,ncols_u;
13473e97b6e8SHong Zhang   const PetscInt *cols;
13483e97b6e8SHong Zhang   PetscScalar    zero=0.0;
13493e97b6e8SHong Zhang 
13503e97b6e8SHong Zhang   PetscFunctionBegin;
13513e97b6e8SHong Zhang   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
13523e97b6e8SHong Zhang   if (nrows != M || ncols != N) SETERRQ4(PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%D by %D must equal %D by %D",nrows,ncols,M,N);
13533e97b6e8SHong Zhang 
13543e97b6e8SHong Zhang   for (row=0; row<nrows; row++) {
13553e97b6e8SHong Zhang     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
13563e97b6e8SHong Zhang     for (j=0; j<ncols_u; j++) {
13573e97b6e8SHong Zhang       col = cols[j] + cstart;
13583e97b6e8SHong Zhang       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
13593e97b6e8SHong Zhang     }
13603e97b6e8SHong Zhang     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
13613e97b6e8SHong Zhang   }
13623e97b6e8SHong Zhang   PetscFunctionReturn(0);
13633e97b6e8SHong Zhang }
13641ad426b7SHong Zhang 
1365e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1366a4e85ca8SHong Zhang {
1367a4e85ca8SHong Zhang   PetscErrorCode ierr;
1368f4431b8cSHong Zhang 
1369a4e85ca8SHong Zhang   PetscFunctionBegin;
1370a4e85ca8SHong Zhang   if (Ju) {
1371a4e85ca8SHong Zhang     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1372a4e85ca8SHong Zhang   } else {
1373a4e85ca8SHong Zhang     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1374a4e85ca8SHong Zhang   }
1375a4e85ca8SHong Zhang   PetscFunctionReturn(0);
1376a4e85ca8SHong Zhang }
1377a4e85ca8SHong Zhang 
137824121865SAdrian Maldonado /* 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.
137924121865SAdrian Maldonado */
138024121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
138124121865SAdrian Maldonado {
138224121865SAdrian Maldonado   PetscErrorCode ierr;
138324121865SAdrian Maldonado   PetscInt       i, size, dof;
138424121865SAdrian Maldonado   PetscInt       *glob2loc;
138524121865SAdrian Maldonado 
138624121865SAdrian Maldonado   PetscFunctionBegin;
138724121865SAdrian Maldonado   ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr);
138824121865SAdrian Maldonado   ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr);
138924121865SAdrian Maldonado 
139024121865SAdrian Maldonado   for (i = 0; i < size; i++) {
139124121865SAdrian Maldonado     ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr);
139224121865SAdrian Maldonado     dof = (dof >= 0) ? dof : -(dof + 1);
139324121865SAdrian Maldonado     glob2loc[i] = dof;
139424121865SAdrian Maldonado   }
139524121865SAdrian Maldonado 
139624121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
139724121865SAdrian Maldonado #if 0
139824121865SAdrian Maldonado   ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
139924121865SAdrian Maldonado #endif
140024121865SAdrian Maldonado   PetscFunctionReturn(0);
140124121865SAdrian Maldonado }
140224121865SAdrian Maldonado 
140301ad2aeeSHong Zhang #include <petsc/private/matimpl.h>
14041ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
14051ad426b7SHong Zhang {
14061ad426b7SHong Zhang   PetscErrorCode ierr;
140724121865SAdrian Maldonado   PetscMPIInt    rank, size;
14081ad426b7SHong Zhang   DM_Network     *network = (DM_Network*) dm->data;
1409a4e85ca8SHong Zhang   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1410840c2264SHong Zhang   PetscInt       cstart,ncols,j,e,v;
141124121865SAdrian Maldonado   PetscBool      ghost,ghost_vc,ghost2,isNest;
1412a4e85ca8SHong Zhang   Mat            Juser;
1413bfbc38dcSHong Zhang   PetscSection   sectionGlobal;
1414447d78afSSatish Balay   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
1415a4e85ca8SHong Zhang   const PetscInt *edges,*cone;
14165cf7da58SHong Zhang   MPI_Comm       comm;
141724121865SAdrian Maldonado   MatType        mtype;
14185cf7da58SHong Zhang   Vec            vd_nz,vo_nz;
14195cf7da58SHong Zhang   PetscInt       *dnnz,*onnz;
14205cf7da58SHong Zhang   PetscScalar    *vdnz,*vonz;
14211ad426b7SHong Zhang 
14221ad426b7SHong Zhang   PetscFunctionBegin;
142324121865SAdrian Maldonado   mtype = dm->mattype;
142424121865SAdrian Maldonado   ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr);
142524121865SAdrian Maldonado 
142624121865SAdrian Maldonado   if (isNest) {
14270731d606SHong Zhang     /* ierr = DMCreateMatrix_Network_Nest(); */
142824121865SAdrian Maldonado     PetscInt   eDof, vDof;
142924121865SAdrian Maldonado     Mat        j11, j12, j21, j22, bA[2][2];
143024121865SAdrian Maldonado     ISLocalToGlobalMapping eISMap, vISMap;
143124121865SAdrian Maldonado 
143224121865SAdrian Maldonado     ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
143324121865SAdrian Maldonado     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
143424121865SAdrian Maldonado     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
143524121865SAdrian Maldonado 
143624121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr);
143724121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr);
143824121865SAdrian Maldonado 
143901ad2aeeSHong Zhang     ierr = MatCreate(comm, &j11);CHKERRQ(ierr);
144024121865SAdrian Maldonado     ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
144124121865SAdrian Maldonado     ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr);
144224121865SAdrian Maldonado 
144301ad2aeeSHong Zhang     ierr = MatCreate(comm, &j12);CHKERRQ(ierr);
144424121865SAdrian Maldonado     ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr);
144524121865SAdrian Maldonado     ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr);
144624121865SAdrian Maldonado 
144701ad2aeeSHong Zhang     ierr = MatCreate(comm, &j21);CHKERRQ(ierr);
144824121865SAdrian Maldonado     ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
144924121865SAdrian Maldonado     ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr);
145024121865SAdrian Maldonado 
145101ad2aeeSHong Zhang     ierr = MatCreate(comm, &j22);CHKERRQ(ierr);
145224121865SAdrian Maldonado     ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
145324121865SAdrian Maldonado     ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr);
145424121865SAdrian Maldonado 
14553f6a6bdaSHong Zhang     bA[0][0] = j11;
14563f6a6bdaSHong Zhang     bA[0][1] = j12;
14573f6a6bdaSHong Zhang     bA[1][0] = j21;
14583f6a6bdaSHong Zhang     bA[1][1] = j22;
145924121865SAdrian Maldonado 
146024121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr);
146124121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr);
146224121865SAdrian Maldonado 
146324121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr);
146424121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr);
146524121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr);
146624121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr);
146724121865SAdrian Maldonado 
146824121865SAdrian Maldonado     ierr = MatSetUp(j11);CHKERRQ(ierr);
146924121865SAdrian Maldonado     ierr = MatSetUp(j12);CHKERRQ(ierr);
147024121865SAdrian Maldonado     ierr = MatSetUp(j21);CHKERRQ(ierr);
147124121865SAdrian Maldonado     ierr = MatSetUp(j22);CHKERRQ(ierr);
147224121865SAdrian Maldonado 
147301ad2aeeSHong Zhang     ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr);
147424121865SAdrian Maldonado     ierr = MatSetUp(*J);CHKERRQ(ierr);
147524121865SAdrian Maldonado     ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr);
147624121865SAdrian Maldonado     ierr = MatDestroy(&j11);CHKERRQ(ierr);
147724121865SAdrian Maldonado     ierr = MatDestroy(&j12);CHKERRQ(ierr);
147824121865SAdrian Maldonado     ierr = MatDestroy(&j21);CHKERRQ(ierr);
147924121865SAdrian Maldonado     ierr = MatDestroy(&j22);CHKERRQ(ierr);
148024121865SAdrian Maldonado 
148124121865SAdrian Maldonado     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
148224121865SAdrian Maldonado     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
148324121865SAdrian Maldonado 
148424121865SAdrian Maldonado     /* Free structures */
148524121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr);
148624121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr);
148724121865SAdrian Maldonado 
148824121865SAdrian Maldonado     PetscFunctionReturn(0);
148924121865SAdrian Maldonado   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1490a4e85ca8SHong Zhang     /* user does not provide Jacobian blocks */
1491bfbc38dcSHong Zhang     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1492bfbc38dcSHong Zhang     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
14931ad426b7SHong Zhang     PetscFunctionReturn(0);
14941ad426b7SHong Zhang   }
14951ad426b7SHong Zhang 
1496bfbc38dcSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
14972a945128SHong Zhang   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1498bfbc38dcSHong Zhang   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1499bfbc38dcSHong Zhang   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
15002a945128SHong Zhang 
15012a945128SHong Zhang   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
15022a945128SHong Zhang   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
150389898e50SHong Zhang 
150489898e50SHong Zhang   /* (1) Set matrix preallocation */
150589898e50SHong Zhang   /*------------------------------*/
1506840c2264SHong Zhang   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1507840c2264SHong Zhang   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1508840c2264SHong Zhang   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1509840c2264SHong Zhang   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1510840c2264SHong Zhang   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1511840c2264SHong Zhang   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1512840c2264SHong Zhang 
151389898e50SHong Zhang   /* Set preallocation for edges */
151489898e50SHong Zhang   /*-----------------------------*/
1515840c2264SHong Zhang   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1516840c2264SHong Zhang 
1517bdcb62a2SHong Zhang   ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr);
1518840c2264SHong Zhang   for (e=eStart; e<eEnd; e++) {
1519840c2264SHong Zhang     /* Get row indices */
1520840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1521840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1522840c2264SHong Zhang     if (nrows) {
1523840c2264SHong Zhang       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1524840c2264SHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1525840c2264SHong Zhang 
15265cf7da58SHong Zhang       /* Set preallocation for conntected vertices */
1527d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1528840c2264SHong Zhang       for (v=0; v<2; v++) {
1529840c2264SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1530840c2264SHong Zhang 
1531840c2264SHong Zhang         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1532840c2264SHong Zhang         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
15335cf7da58SHong Zhang         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1534840c2264SHong Zhang       }
1535840c2264SHong Zhang 
153689898e50SHong Zhang       /* Set preallocation for edge self */
1537840c2264SHong Zhang       cstart = rstart;
1538840c2264SHong Zhang       Juser = network->Je[3*e]; /* Jacobian(e,e) */
15395cf7da58SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1540840c2264SHong Zhang     }
1541840c2264SHong Zhang   }
1542840c2264SHong Zhang 
154389898e50SHong Zhang   /* Set preallocation for vertices */
154489898e50SHong Zhang   /*--------------------------------*/
1545840c2264SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1546840c2264SHong Zhang   if (vEnd - vStart) {
1547840c2264SHong Zhang     if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Jv");
1548840c2264SHong Zhang     vptr = network->Jvptr;
1549840c2264SHong Zhang     if (!vptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide vptr");
1550840c2264SHong Zhang   }
1551840c2264SHong Zhang 
1552840c2264SHong Zhang   for (v=vStart; v<vEnd; v++) {
1553840c2264SHong Zhang     /* Get row indices */
1554840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1555840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1556840c2264SHong Zhang     if (!nrows) continue;
1557840c2264SHong Zhang 
1558bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1559bdcb62a2SHong Zhang     if (ghost) {
1560bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1561bdcb62a2SHong Zhang     } else {
1562bdcb62a2SHong Zhang       rows_v = rows;
1563bdcb62a2SHong Zhang     }
1564bdcb62a2SHong Zhang 
1565bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1566840c2264SHong Zhang 
1567840c2264SHong Zhang     /* Get supporting edges and connected vertices */
1568840c2264SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1569840c2264SHong Zhang 
1570840c2264SHong Zhang     for (e=0; e<nedges; e++) {
1571840c2264SHong Zhang       /* Supporting edges */
1572840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1573840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1574840c2264SHong Zhang 
1575840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1576bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1577840c2264SHong Zhang 
1578840c2264SHong Zhang       /* Connected vertices */
1579d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
1580840c2264SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
1581840c2264SHong Zhang       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1582840c2264SHong Zhang 
1583840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1584840c2264SHong Zhang 
1585840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1586e102a522SHong Zhang       if (ghost_vc||ghost) {
1587e102a522SHong Zhang         ghost2 = PETSC_TRUE;
1588e102a522SHong Zhang       } else {
1589e102a522SHong Zhang         ghost2 = PETSC_FALSE;
1590e102a522SHong Zhang       }
1591e102a522SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr);
1592840c2264SHong Zhang     }
1593840c2264SHong Zhang 
159489898e50SHong Zhang     /* Set preallocation for vertex self */
1595840c2264SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1596840c2264SHong Zhang     if (!ghost) {
1597840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1598840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1599bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1600840c2264SHong Zhang     }
1601bdcb62a2SHong Zhang     if (ghost) {
1602bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1603bdcb62a2SHong Zhang     }
1604840c2264SHong Zhang   }
1605840c2264SHong Zhang 
1606840c2264SHong Zhang   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1607840c2264SHong Zhang   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
16085cf7da58SHong Zhang 
16095cf7da58SHong Zhang   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
16105cf7da58SHong Zhang 
16115cf7da58SHong Zhang   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1612840c2264SHong Zhang   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1613840c2264SHong Zhang 
1614840c2264SHong Zhang   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1615840c2264SHong Zhang   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1616840c2264SHong Zhang   for (j=0; j<localSize; j++) {
1617e102a522SHong Zhang     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1618e102a522SHong Zhang     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1619840c2264SHong Zhang   }
1620840c2264SHong Zhang   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1621840c2264SHong Zhang   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1622840c2264SHong Zhang   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1623840c2264SHong Zhang   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1624840c2264SHong Zhang 
16255cf7da58SHong Zhang   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
16265cf7da58SHong Zhang   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
16275cf7da58SHong Zhang   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
16285cf7da58SHong Zhang 
16295cf7da58SHong Zhang   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
16305cf7da58SHong Zhang 
163189898e50SHong Zhang   /* (2) Set matrix entries for edges */
163289898e50SHong Zhang   /*----------------------------------*/
16331ad426b7SHong Zhang   for (e=eStart; e<eEnd; e++) {
1634bfbc38dcSHong Zhang     /* Get row indices */
16351ad426b7SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
163617df6e9eSHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
16374b976069SHong Zhang     if (nrows) {
16384b976069SHong Zhang       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1639bdcb62a2SHong Zhang 
164017df6e9eSHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
16411ad426b7SHong Zhang 
1642bfbc38dcSHong Zhang       /* Set matrix entries for conntected vertices */
1643d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1644bfbc38dcSHong Zhang       for (v=0; v<2; v++) {
1645bfbc38dcSHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1646883e35e8SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
16473e97b6e8SHong Zhang 
1648a4e85ca8SHong Zhang         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1649a4e85ca8SHong Zhang         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1650bfbc38dcSHong Zhang       }
165117df6e9eSHong Zhang 
1652bfbc38dcSHong Zhang       /* Set matrix entries for edge self */
16533e97b6e8SHong Zhang       cstart = rstart;
1654a4e85ca8SHong Zhang       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1655a4e85ca8SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
16561ad426b7SHong Zhang     }
16574b976069SHong Zhang   }
16581ad426b7SHong Zhang 
1659bfbc38dcSHong Zhang   /* Set matrix entries for vertices */
166083b2e829SHong Zhang   /*---------------------------------*/
16611ad426b7SHong Zhang   for (v=vStart; v<vEnd; v++) {
1662bfbc38dcSHong Zhang     /* Get row indices */
1663596e729fSHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1664596e729fSHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
16654b976069SHong Zhang     if (!nrows) continue;
1666596e729fSHong Zhang 
1667bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1668bdcb62a2SHong Zhang     if (ghost) {
1669bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1670bdcb62a2SHong Zhang     } else {
1671bdcb62a2SHong Zhang       rows_v = rows;
1672bdcb62a2SHong Zhang     }
1673bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1674596e729fSHong Zhang 
1675bfbc38dcSHong Zhang     /* Get supporting edges and connected vertices */
1676596e729fSHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1677596e729fSHong Zhang 
1678596e729fSHong Zhang     for (e=0; e<nedges; e++) {
1679bfbc38dcSHong Zhang       /* Supporting edges */
1680596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1681596e729fSHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1682596e729fSHong Zhang 
1683a4e85ca8SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1684bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1685596e729fSHong Zhang 
1686bfbc38dcSHong Zhang       /* Connected vertices */
1687d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
16882a945128SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
16892a945128SHong Zhang 
169044aca652SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
169144aca652SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1692a4e85ca8SHong Zhang 
1693a4e85ca8SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1694bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1695596e729fSHong Zhang     }
1696596e729fSHong Zhang 
1697bfbc38dcSHong Zhang     /* Set matrix entries for vertex self */
16981ad426b7SHong Zhang     if (!ghost) {
1699596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1700a4e85ca8SHong Zhang       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1701bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr);
1702bdcb62a2SHong Zhang     }
1703bdcb62a2SHong Zhang     if (ghost) {
1704bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1705bdcb62a2SHong Zhang     }
17061ad426b7SHong Zhang   }
1707a4e85ca8SHong Zhang   ierr = PetscFree(rows);CHKERRQ(ierr);
1708bdcb62a2SHong Zhang 
17091ad426b7SHong Zhang   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17101ad426b7SHong Zhang   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1711dd6f46cdSHong Zhang 
17125f2c45f1SShri Abhyankar   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
17135f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
17145f2c45f1SShri Abhyankar }
17155f2c45f1SShri Abhyankar 
17165f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm)
17175f2c45f1SShri Abhyankar {
17185f2c45f1SShri Abhyankar   PetscErrorCode ierr;
17195f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
17202727e31bSShri Abhyankar   PetscInt       j;
17215f2c45f1SShri Abhyankar 
17225f2c45f1SShri Abhyankar   PetscFunctionBegin;
17238415c774SShri Abhyankar   if (--network->refct > 0) PetscFunctionReturn(0);
172483b2e829SHong Zhang   if (network->Je) {
172583b2e829SHong Zhang     ierr = PetscFree(network->Je);CHKERRQ(ierr);
172683b2e829SHong Zhang   }
172783b2e829SHong Zhang   if (network->Jv) {
1728883e35e8SHong Zhang     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
172983b2e829SHong Zhang     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
17301ad426b7SHong Zhang   }
173113c2a604SAdrian Maldonado 
173213c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr);
173313c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr);
173413c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr);
173513c2a604SAdrian Maldonado   if (network->vertex.sf) {
173613c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr);
173713c2a604SAdrian Maldonado   }
173813c2a604SAdrian Maldonado   /* edge */
173913c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr);
174013c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr);
174113c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr);
174213c2a604SAdrian Maldonado   if (network->edge.sf) {
174313c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr);
174413c2a604SAdrian Maldonado   }
17455f2c45f1SShri Abhyankar   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
17465f2c45f1SShri Abhyankar   network->edges = NULL;
17475f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
17485f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
174983b2e829SHong Zhang 
17502727e31bSShri Abhyankar   for(j=0; j < network->nsubnet; j++) {
17512727e31bSShri Abhyankar     ierr = PetscFree(network->subnet[j].edges);CHKERRQ(ierr);
17522727e31bSShri Abhyankar     ierr = PetscFree(network->subnet[j].vertices);CHKERRQ(ierr);
17532727e31bSShri Abhyankar   }
1754e2aaf10cSShri Abhyankar   ierr = PetscFree(network->subnet);CHKERRQ(ierr);
17555f2c45f1SShri Abhyankar   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
17565f2c45f1SShri Abhyankar   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
17575f2c45f1SShri Abhyankar   ierr = PetscFree(network->header);CHKERRQ(ierr);
17585f2c45f1SShri Abhyankar   ierr = PetscFree(network);CHKERRQ(ierr);
17595f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
17605f2c45f1SShri Abhyankar }
17615f2c45f1SShri Abhyankar 
17625f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
17635f2c45f1SShri Abhyankar {
17645f2c45f1SShri Abhyankar   PetscErrorCode ierr;
17655f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
17665f2c45f1SShri Abhyankar 
17675f2c45f1SShri Abhyankar   PetscFunctionBegin;
17685f2c45f1SShri Abhyankar   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
17695f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
17705f2c45f1SShri Abhyankar }
17715f2c45f1SShri Abhyankar 
17725f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
17735f2c45f1SShri Abhyankar {
17745f2c45f1SShri Abhyankar   PetscErrorCode ierr;
17755f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
17765f2c45f1SShri Abhyankar 
17775f2c45f1SShri Abhyankar   PetscFunctionBegin;
17785f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
17795f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
17805f2c45f1SShri Abhyankar }
17815f2c45f1SShri Abhyankar 
17825f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
17835f2c45f1SShri Abhyankar {
17845f2c45f1SShri Abhyankar   PetscErrorCode ierr;
17855f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
17865f2c45f1SShri Abhyankar 
17875f2c45f1SShri Abhyankar   PetscFunctionBegin;
17885f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
17895f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
17905f2c45f1SShri Abhyankar }
17915f2c45f1SShri Abhyankar 
17925f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
17935f2c45f1SShri Abhyankar {
17945f2c45f1SShri Abhyankar   PetscErrorCode ierr;
17955f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
17965f2c45f1SShri Abhyankar 
17975f2c45f1SShri Abhyankar   PetscFunctionBegin;
17985f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
17995f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
18005f2c45f1SShri Abhyankar }
18015f2c45f1SShri Abhyankar 
18025f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
18035f2c45f1SShri Abhyankar {
18045f2c45f1SShri Abhyankar   PetscErrorCode ierr;
18055f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
18065f2c45f1SShri Abhyankar 
18075f2c45f1SShri Abhyankar   PetscFunctionBegin;
18085f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
18095f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
18105f2c45f1SShri Abhyankar }
1811