xref: /petsc/src/dm/impls/network/network.c (revision 3ebf9fb9b401d3c2866effc47ff393849bf95196)
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
34caf410d2SHong Zhang . Nsubnet - global number of subnetworks
35caf410d2SHong Zhang . nV - number of local vertices for each subnetwork
36caf410d2SHong Zhang . nE - number of local edges for each subnetwork
37caf410d2SHong Zhang . NsubnetCouple - global number of coupling subnetworks
38caf410d2SHong Zhang - nec - number of local edges for each coupling subnetwork
395f2c45f1SShri Abhyankar 
40caf410d2SHong Zhang    You cannot change the sizes once they have been set. nV, nE are arrays of length Nsubnet, and nec is array of length NsubnetCouple.
415f2c45f1SShri Abhyankar 
421b266c99SBarry Smith    Level: intermediate
431b266c99SBarry Smith 
441b266c99SBarry Smith .seealso: DMNetworkCreate()
455f2c45f1SShri Abhyankar @*/
46caf410d2SHong Zhang PetscErrorCode DMNetworkSetSizes(DM dm,PetscInt Nsubnet,PetscInt nV[], PetscInt nE[],PetscInt NsubnetCouple,PetscInt nec[])
475f2c45f1SShri Abhyankar {
485f2c45f1SShri Abhyankar   PetscErrorCode ierr;
495f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
50e2aaf10cSShri Abhyankar   PetscInt       a[2],b[2],i;
515f2c45f1SShri Abhyankar 
525f2c45f1SShri Abhyankar   PetscFunctionBegin;
535f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
54e2aaf10cSShri Abhyankar   if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet);
557765340cSHong Zhang   if (NsubnetCouple < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of coupling subnetworks %D cannot be less than 0",NsubnetCouple);
567765340cSHong Zhang 
57caf410d2SHong Zhang   PetscValidLogicalCollectiveInt(dm,Nsubnet,2);
58caf410d2SHong Zhang   if (NsubnetCouple) PetscValidLogicalCollectiveInt(dm,NsubnetCouple,5);
597765340cSHong Zhang   if (network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network");
607765340cSHong Zhang 
61caf410d2SHong Zhang   if (!nV || !nE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Local vertex size or edge size must be provided");
62f025b11dSHong Zhang 
63caf410d2SHong Zhang   network->nsubnet  = Nsubnet + NsubnetCouple;
64caf410d2SHong Zhang   network->ncsubnet = NsubnetCouple;
65caf410d2SHong Zhang   ierr = PetscCalloc1(Nsubnet+NsubnetCouple,&network->subnet);CHKERRQ(ierr);
66caf410d2SHong Zhang 
67caf410d2SHong Zhang   /* ----------------------------------------------------------
68caf410d2SHong Zhang    p=v or e; P=V or E
69caf410d2SHong Zhang    subnet[0].pStart   = 0
70caf410d2SHong Zhang    subnet[i+1].pStart = subnet[i].pEnd = subnet[i].pStart + (nE[i] or NV[i])
71caf410d2SHong Zhang    ----------------------------------------------------------------------- */
72caf410d2SHong Zhang   for (i=0; i < Nsubnet; i++) {
73caf410d2SHong Zhang     /* Get global number of vertices and edges for subnet[i] */
74caf410d2SHong Zhang     a[0] = nV[i]; a[1] = nE[i]; /* local number of vertices (excluding ghost) and edges */
757765340cSHong Zhang     ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
767765340cSHong Zhang     network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1];
777765340cSHong Zhang 
78caf410d2SHong Zhang     network->subnet[i].nvtx   = nV[i]; /* local nvtx, without ghost */
797765340cSHong Zhang 
80caf410d2SHong Zhang     /* global subnet[].vStart and vEnd, used by DMNetworkLayoutSetUp() */
81caf410d2SHong Zhang     network->subnet[i].vStart = network->NVertices;
827765340cSHong Zhang     network->subnet[i].vEnd   = network->subnet[i].vStart + network->subnet[i].Nvtx;
83caf410d2SHong Zhang 
84caf410d2SHong Zhang     network->nVertices += nV[i];
857765340cSHong Zhang     network->NVertices += network->subnet[i].Nvtx;
867765340cSHong Zhang 
877765340cSHong Zhang     network->subnet[i].nedge  = nE[i];
887765340cSHong Zhang     network->subnet[i].eStart = network->nEdges;
89caf410d2SHong Zhang     network->subnet[i].eEnd   = network->subnet[i].eStart + nE[i];
90caf410d2SHong Zhang     network->nEdges += nE[i];
91caf410d2SHong Zhang     network->NEdges += network->subnet[i].Nedge;
92caf410d2SHong Zhang   }
93caf410d2SHong Zhang 
94caf410d2SHong Zhang   /* coupling subnetwork */
95caf410d2SHong Zhang   for (; i < Nsubnet+NsubnetCouple; i++) {
96caf410d2SHong Zhang     /* Get global number of coupling edges for subnet[i] */
97caf410d2SHong Zhang     ierr = MPIU_Allreduce(nec+(i-Nsubnet),&network->subnet[i].Nedge,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
98caf410d2SHong Zhang 
99caf410d2SHong Zhang     network->subnet[i].nvtx   = 0; /* We design coupling subnetwork such that it does not have its own vertices */
100caf410d2SHong Zhang     network->subnet[i].vStart = network->nVertices;
101caf410d2SHong Zhang     network->subnet[i].vEnd   = network->subnet[i].vStart;
102caf410d2SHong Zhang 
103caf410d2SHong Zhang     network->subnet[i].nedge  = nec[i-Nsubnet];
104caf410d2SHong Zhang     network->subnet[i].eStart = network->nEdges;
105caf410d2SHong Zhang     network->subnet[i].eEnd = network->subnet[i].eStart + nec[i-Nsubnet];
106caf410d2SHong Zhang     network->nEdges += nec[i-Nsubnet];
1077765340cSHong Zhang     network->NEdges += network->subnet[i].Nedge;
1087765340cSHong Zhang   }
1097765340cSHong Zhang   PetscFunctionReturn(0);
1107765340cSHong Zhang }
1117765340cSHong Zhang 
1125f2c45f1SShri Abhyankar /*@
1135f2c45f1SShri Abhyankar   DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network
1145f2c45f1SShri Abhyankar 
1155f2c45f1SShri Abhyankar   Logically collective on DM
1165f2c45f1SShri Abhyankar 
1175f2c45f1SShri Abhyankar   Input Parameters:
118e3e68989SHong Zhang + dm - the dm object
119e3e68989SHong Zhang . edgelist - list of edges for each subnetwork
120e3e68989SHong Zhang - edgelistCouple - list of edges for each coupling subnetwork
1215f2c45f1SShri Abhyankar 
1225f2c45f1SShri Abhyankar   Notes:
1235f2c45f1SShri Abhyankar   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
1245f2c45f1SShri Abhyankar   not be destroyed before the call to DMNetworkLayoutSetUp
1255f2c45f1SShri Abhyankar 
1265f2c45f1SShri Abhyankar   Level: intermediate
1275f2c45f1SShri Abhyankar 
128e3e68989SHong Zhang   Example usage:
129e3e68989SHong Zhang   Consider the following 2 separate networks and a coupling network:
130e3e68989SHong Zhang 
131e3e68989SHong Zhang .vb
132e3e68989SHong Zhang  network 0: v0 -> v1 -> v2 -> v3
133e3e68989SHong Zhang  network 1: v1 -> v2 -> v0
134e3e68989SHong Zhang  coupling network: network 1: v2 -> network 0: v0
135e3e68989SHong Zhang .ve
136e3e68989SHong Zhang 
137e3e68989SHong Zhang  The resulting input
138e3e68989SHong Zhang    edgelist[0] = [0 1 | 1 2 | 2 3];
139e3e68989SHong Zhang    edgelist[1] = [1 2 | 2 0]
140e3e68989SHong Zhang    edgelistCouple[0] = [(network)1 (v)2 (network)0 (v)0].
141e3e68989SHong Zhang 
1425f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes
1435f2c45f1SShri Abhyankar @*/
1444e18019cSBarry Smith PetscErrorCode DMNetworkSetEdgeList(DM dm,PetscInt *edgelist[],PetscInt *edgelistCouple[])
1455f2c45f1SShri Abhyankar {
1465f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network*) dm->data;
147e3e68989SHong Zhang   PetscInt   i,nsubnet,ncsubnet=network->ncsubnet;
1485f2c45f1SShri Abhyankar 
1495f2c45f1SShri Abhyankar   PetscFunctionBegin;
150e3e68989SHong Zhang   nsubnet = network->nsubnet - ncsubnet;
151e3e68989SHong Zhang   for(i=0; i < nsubnet; i++) {
152e2aaf10cSShri Abhyankar     network->subnet[i].edgelist = edgelist[i];
153e2aaf10cSShri Abhyankar   }
154e3e68989SHong Zhang   if (edgelistCouple) {
155e3e68989SHong Zhang     PetscInt j;
156e3e68989SHong Zhang     j = 0;
157e3e68989SHong Zhang     nsubnet = network->nsubnet;
158e3e68989SHong Zhang     while (i < nsubnet) network->subnet[i++].edgelist = edgelistCouple[j++];
159e3e68989SHong Zhang   }
1605f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
1615f2c45f1SShri Abhyankar }
1625f2c45f1SShri Abhyankar 
1635f2c45f1SShri Abhyankar /*@
1645f2c45f1SShri Abhyankar   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
1655f2c45f1SShri Abhyankar 
1665f2c45f1SShri Abhyankar   Collective on DM
1675f2c45f1SShri Abhyankar 
1685f2c45f1SShri Abhyankar   Input Parameters
1695f2c45f1SShri Abhyankar . DM - the dmnetwork object
1705f2c45f1SShri Abhyankar 
1715f2c45f1SShri Abhyankar   Notes:
1725f2c45f1SShri Abhyankar   This routine should be called after the network sizes and edgelists have been provided. It creates
1735f2c45f1SShri Abhyankar   the bare layout of the network and sets up the network to begin insertion of components.
1745f2c45f1SShri Abhyankar 
1755f2c45f1SShri Abhyankar   All the components should be registered before calling this routine.
1765f2c45f1SShri Abhyankar 
1775f2c45f1SShri Abhyankar   Level: intermediate
1785f2c45f1SShri Abhyankar 
1795f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList
1805f2c45f1SShri Abhyankar @*/
1815f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm)
1825f2c45f1SShri Abhyankar {
1835f2c45f1SShri Abhyankar   PetscErrorCode ierr;
1845f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
185caf410d2SHong Zhang   PetscInt       numCorners=2,spacedim=2,dim = 1; /* One dimensional network */
186*3ebf9fb9SHong Zhang   PetscReal      *vertexcoords=NULL;
187caf410d2SHong Zhang   PetscInt       i,j,ctr,nsubnet,*eowners,np,*edges,*subnetvtx,vStart;
188caf410d2SHong Zhang   PetscInt       k,netid,vid, *vidxlTog,*edgelist_couple=NULL;
189caf410d2SHong Zhang   const PetscInt *cone;
190caf410d2SHong Zhang   MPI_Comm       comm;
191caf410d2SHong Zhang   PetscMPIInt    size,rank;
1926500d4abSHong Zhang 
1936500d4abSHong Zhang   PetscFunctionBegin;
194caf410d2SHong Zhang   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
195caf410d2SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
196caf410d2SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1976500d4abSHong Zhang 
198caf410d2SHong Zhang   /* Create the local edgelist for the network by concatenating local input edgelists of the subnetworks */
199caf410d2SHong Zhang   ierr = PetscCalloc2(numCorners*network->nVertices,&vertexcoords,2*network->nEdges,&edges);CHKERRQ(ierr);
2007765340cSHong Zhang   nsubnet = network->nsubnet - network->ncsubnet;
201caf410d2SHong Zhang   ctr = 0;
2027765340cSHong Zhang   for (i=0; i < nsubnet; i++) {
2036500d4abSHong Zhang     for (j = 0; j < network->subnet[i].nedge; j++) {
204caf410d2SHong Zhang       edges[2*ctr]   = network->subnet[i].eStart + network->subnet[i].edgelist[2*j];
205caf410d2SHong Zhang       edges[2*ctr+1] = network->subnet[i].eStart + network->subnet[i].edgelist[2*j+1];
2066500d4abSHong Zhang       ctr++;
2076500d4abSHong Zhang     }
2086500d4abSHong Zhang   }
2097765340cSHong Zhang 
210caf410d2SHong Zhang   /* Append local coupling edgelists of the subnetworks */
2117765340cSHong Zhang   i       = nsubnet; /* netid of coupling subnet */
2127765340cSHong Zhang   nsubnet = network->nsubnet;
2137765340cSHong Zhang   while (i < nsubnet) {
214991cf414SHong Zhang     edgelist_couple = network->subnet[i].edgelist;
2156500d4abSHong Zhang     k = 0;
2166500d4abSHong Zhang     for (j = 0; j < network->subnet[i].nedge; j++) {
2176500d4abSHong Zhang       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
218caf410d2SHong Zhang       edges[2*ctr] = network->subnet[netid].vStart + vid; k += 2;
2196500d4abSHong Zhang 
2206500d4abSHong Zhang       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
221caf410d2SHong Zhang       edges[2*ctr+1] = network->subnet[netid].vStart + vid; k+=2;
2226500d4abSHong Zhang       ctr++;
2236500d4abSHong Zhang     }
2247765340cSHong Zhang     i++;
2257765340cSHong Zhang   }
226caf410d2SHong Zhang   /*
227caf410d2SHong Zhang   if (rank == 0) {
228caf410d2SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] edgelist:\n",rank);
2296500d4abSHong Zhang     for(i=0; i < network->nEdges; i++) {
230caf410d2SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",edges[2*i],edges[2*i+1]);CHKERRQ(ierr);
2316500d4abSHong Zhang       printf("\n");
2326500d4abSHong Zhang     }
233caf410d2SHong Zhang   }
234caf410d2SHong Zhang    */
2356500d4abSHong Zhang 
236caf410d2SHong Zhang   /* Create network->plex */
237acdb140fSBarry Smith #if defined(PETSC_USE_64BIT_INDICES)
238acdb140fSBarry Smith   {
239caf410d2SHong Zhang     int *edges64;
240caf410d2SHong Zhang     np = network->nEdges*numCorners;
241caf410d2SHong Zhang     ierr = PetscMalloc1(np,&edges64);CHKERRQ(ierr);
242caf410d2SHong Zhang     for (i=0; i<np; i++) edges64[i] = (int)edges[i];
243caf410d2SHong Zhang 
244caf410d2SHong Zhang     if (size == 1) {
245caf410d2SHong Zhang       ierr = DMPlexCreateFromCellList(comm,dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,(const int*)edges64,spacedim,(const double*)vertexcoords,&network->plex);CHKERRQ(ierr);
246caf410d2SHong Zhang     } else {
247*3ebf9fb9SHong Zhang       ierr = DMPlexCreateFromCellListParallel(comm,dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,(const int*)edges64,spacedim,(const PetscReal*)vertexcoords,NULL,&network->plex);CHKERRQ(ierr);
248acdb140fSBarry Smith     }
249caf410d2SHong Zhang     ierr = PetscFree(edges64);CHKERRQ(ierr);
250acdb140fSBarry Smith   }
251acdb140fSBarry Smith #else
252caf410d2SHong Zhang   if (size == 1) {
253caf410d2SHong Zhang     ierr = DMPlexCreateFromCellList(comm,dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,(const int*)edges,spacedim,(const double*)vertexcoords,&network->plex);CHKERRQ(ierr);
254caf410d2SHong Zhang   } else {
255*3ebf9fb9SHong Zhang     ierr = DMPlexCreateFromCellListParallel(comm,dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,(const int*)edges,spacedim,(const PetscReal*)vertexcoords,NULL,&network->plex);CHKERRQ(ierr);
2566500d4abSHong Zhang   }
257caf410d2SHong Zhang #endif
2586500d4abSHong Zhang 
2596500d4abSHong Zhang   ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr);
2606500d4abSHong Zhang   ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr);
2616500d4abSHong Zhang   ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr);
262caf410d2SHong Zhang   vStart = network->vStart;
2636500d4abSHong Zhang 
264caf410d2SHong Zhang   ierr = PetscSectionCreate(comm,&network->DataSection);CHKERRQ(ierr);
265caf410d2SHong Zhang   ierr = PetscSectionCreate(comm,&network->DofSection);CHKERRQ(ierr);
2666500d4abSHong Zhang   ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr);
2676500d4abSHong Zhang   ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr);
2686500d4abSHong Zhang 
269caf410d2SHong Zhang   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
270caf410d2SHong Zhang   np = network->pEnd - network->pStart;
271caf410d2SHong Zhang   ierr = PetscCalloc2(np,&network->header,np,&network->cvalue);CHKERRQ(ierr);
272caf410d2SHong Zhang 
273caf410d2SHong Zhang   /* Create vidxlTog: maps local vertex index to global index */
274caf410d2SHong Zhang   np = network->vEnd - vStart;
275caf410d2SHong Zhang   ierr = PetscMalloc2(np,&vidxlTog,size+1,&eowners);CHKERRQ(ierr);
276caf410d2SHong Zhang   ctr = 0;
277caf410d2SHong Zhang   for (i=network->eStart; i<network->eEnd; i++) {
278caf410d2SHong Zhang     ierr = DMNetworkGetConnectedVertices(dm,i,&cone);CHKERRQ(ierr);
279caf410d2SHong Zhang     vidxlTog[cone[0] - vStart] = edges[2*ctr];
280caf410d2SHong Zhang     vidxlTog[cone[1] - vStart] = edges[2*ctr+1];
281caf410d2SHong Zhang     ctr++;
282caf410d2SHong Zhang   }
283caf410d2SHong Zhang   ierr = PetscFree2(vertexcoords,edges);CHKERRQ(ierr);
284caf410d2SHong Zhang 
2856500d4abSHong Zhang   /* Create vertices and edges array for the subnetworks */
2866500d4abSHong Zhang   for (j=0; j < network->nsubnet; j++) {
2876500d4abSHong Zhang     ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr);
288caf410d2SHong Zhang 
2896500d4abSHong Zhang     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
2906500d4abSHong Zhang        These get updated when the vertices and edges are added. */
291caf410d2SHong Zhang     network->subnet[j].nvtx  = 0;
292caf410d2SHong Zhang     network->subnet[j].nedge = 0;
2936500d4abSHong Zhang   }
294caf410d2SHong Zhang   ierr = PetscCalloc1(np,&network->subnetvtx);CHKERRQ(ierr);
2956500d4abSHong Zhang 
296caf410d2SHong Zhang 
297caf410d2SHong Zhang   /* Get edge ownership */
298caf410d2SHong Zhang   np = network->eEnd - network->eStart;
299caf410d2SHong Zhang   ierr = MPI_Allgather(&np,1,MPIU_INT,eowners+1,1,MPIU_INT,comm);CHKERRQ(ierr);
300caf410d2SHong Zhang   eowners[0] = 0;
301caf410d2SHong Zhang   for (i=2; i<=size; i++) eowners[i] += eowners[i-1];
302caf410d2SHong Zhang 
303caf410d2SHong Zhang   i = 0; j = 0;
304caf410d2SHong Zhang   while (i < np) { /* local edges, including coupling edges */
305caf410d2SHong Zhang     network->header[i].index = i + eowners[rank];   /* Global edge index */
306caf410d2SHong Zhang 
307caf410d2SHong Zhang     if (j < network->nsubnet && i < network->subnet[j].eEnd) {
3086500d4abSHong Zhang       network->header[i].subnetid = j; /* Subnetwork id */
3096500d4abSHong Zhang       network->subnet[j].edges[network->subnet[j].nedge++] = i;
310caf410d2SHong Zhang 
311caf410d2SHong Zhang       network->header[i].ndata = 0;
312caf410d2SHong Zhang       ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
313caf410d2SHong Zhang       network->header[i].offset[0] = 0;
314caf410d2SHong Zhang       i++;
315caf410d2SHong Zhang     }
316caf410d2SHong Zhang     if (i >= network->subnet[j].eEnd) j++;
317caf410d2SHong Zhang   }
318caf410d2SHong Zhang 
319caf410d2SHong Zhang   /* Count network->subnet[*].nvtx */
320caf410d2SHong Zhang   for (i=vStart; i<network->vEnd; i++) { /* local vertices, including ghosts */
321caf410d2SHong Zhang     k = vidxlTog[i-vStart];
322caf410d2SHong Zhang     for (j=0; j < network->nsubnet; j++) {
323caf410d2SHong Zhang       if (network->subnet[j].vStart <= k && k < network->subnet[j].vEnd) {
324caf410d2SHong Zhang         network->subnet[j].nvtx++;
3256500d4abSHong Zhang         break;
3266500d4abSHong Zhang       }
3276500d4abSHong Zhang     }
3286500d4abSHong Zhang   }
3296500d4abSHong Zhang 
330caf410d2SHong Zhang   /* Set network->subnet[*].vertices on array network->subnetvtx */
331caf410d2SHong Zhang   subnetvtx = network->subnetvtx;
3326500d4abSHong Zhang   for (j=0; j<network->nsubnet; j++) {
333caf410d2SHong Zhang     network->subnet[j].vertices = subnetvtx;
334caf410d2SHong Zhang     subnetvtx                  += network->subnet[j].nvtx;
335caf410d2SHong Zhang     network->subnet[j].nvtx = 0;
336caf410d2SHong Zhang   }
337caf410d2SHong Zhang 
338caf410d2SHong Zhang   /* Set vertex array for the subnetworks */
339caf410d2SHong Zhang   for (i=vStart; i<network->vEnd; i++) { /* local vertices, including ghosts */
340caf410d2SHong Zhang     network->header[i].index = vidxlTog[i-vStart]; /*  Global vertex index */
341caf410d2SHong Zhang 
342caf410d2SHong Zhang     k = vidxlTog[i-vStart];
343caf410d2SHong Zhang     for (j=0; j < network->nsubnet; j++) {
344caf410d2SHong Zhang       if (network->subnet[j].vStart <= k && k < network->subnet[j].vEnd) {
3456500d4abSHong Zhang         network->header[i].subnetid = j;
3466500d4abSHong Zhang         network->subnet[j].vertices[network->subnet[j].nvtx++] = i;
3476500d4abSHong Zhang         break;
3486500d4abSHong Zhang       }
3496500d4abSHong Zhang     }
350caf410d2SHong Zhang 
3516500d4abSHong Zhang     network->header[i].ndata = 0;
3526500d4abSHong Zhang     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
353caf410d2SHong Zhang     network->header[i].offset[0] = 0;
3546500d4abSHong Zhang   }
3556500d4abSHong Zhang 
356caf410d2SHong Zhang   ierr = PetscFree2(vidxlTog,eowners);CHKERRQ(ierr);
3575f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3585f2c45f1SShri Abhyankar }
3595f2c45f1SShri Abhyankar 
36094ef8ddeSSatish Balay /*@C
3612727e31bSShri Abhyankar   DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork
3622727e31bSShri Abhyankar 
3632727e31bSShri Abhyankar   Input Parameters
364caf410d2SHong Zhang + dm - the DM object
3652727e31bSShri Abhyankar - id   - the ID (integer) of the subnetwork
3662727e31bSShri Abhyankar 
3672727e31bSShri Abhyankar   Output Parameters
3682727e31bSShri Abhyankar + nv    - number of vertices (local)
3692727e31bSShri Abhyankar . ne    - number of edges (local)
3702727e31bSShri Abhyankar . vtx   - local vertices for this subnetwork
3712727e31bSShri Abhyankar . edge  - local edges for this subnetwork
3722727e31bSShri Abhyankar 
3732727e31bSShri Abhyankar   Notes:
3742727e31bSShri Abhyankar   Cannot call this routine before DMNetworkLayoutSetup()
3752727e31bSShri Abhyankar 
37606dd6b0eSSatish Balay   Level: intermediate
37706dd6b0eSSatish Balay 
3782727e31bSShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
3792727e31bSShri Abhyankar @*/
380caf410d2SHong Zhang PetscErrorCode DMNetworkGetSubnetworkInfo(DM dm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge)
3812727e31bSShri Abhyankar {
382caf410d2SHong Zhang   DM_Network *network = (DM_Network*)dm->data;
3832727e31bSShri Abhyankar 
3842727e31bSShri Abhyankar   PetscFunctionBegin;
3852727e31bSShri Abhyankar   *nv   = network->subnet[id].nvtx;
3862727e31bSShri Abhyankar   *ne   = network->subnet[id].nedge;
3872727e31bSShri Abhyankar   *vtx  = network->subnet[id].vertices;
3882727e31bSShri Abhyankar   *edge = network->subnet[id].edges;
3892727e31bSShri Abhyankar   PetscFunctionReturn(0);
3902727e31bSShri Abhyankar }
3912727e31bSShri Abhyankar 
3922727e31bSShri Abhyankar /*@C
393caf410d2SHong Zhang   DMNetworkGetSubnetworkCoupleInfo - Returns the info for the coupling subnetwork
394caf410d2SHong Zhang 
395caf410d2SHong Zhang   Input Parameters
396caf410d2SHong Zhang + dm - the DM object
397caf410d2SHong Zhang - id   - the ID (integer) of the coupling subnetwork
398caf410d2SHong Zhang 
399caf410d2SHong Zhang   Output Parameters
400caf410d2SHong Zhang + ne - number of edges (local)
401caf410d2SHong Zhang - edge  - local edges for this coupling subnetwork
402caf410d2SHong Zhang 
403caf410d2SHong Zhang   Notes:
404caf410d2SHong Zhang   Cannot call this routine before DMNetworkLayoutSetup()
405caf410d2SHong Zhang 
406caf410d2SHong Zhang   Level: intermediate
407caf410d2SHong Zhang 
408caf410d2SHong Zhang .seealso: DMNetworkGetSubnetworkInfo, DMNetworkLayoutSetUp, DMNetworkCreate
409caf410d2SHong Zhang @*/
410caf410d2SHong Zhang PetscErrorCode DMNetworkGetSubnetworkCoupleInfo(DM dm,PetscInt id,PetscInt *ne,const PetscInt **edge)
411caf410d2SHong Zhang {
412caf410d2SHong Zhang   DM_Network *net = (DM_Network*)dm->data;
413caf410d2SHong Zhang   PetscInt   id1 = id + net->nsubnet - net->ncsubnet;
414caf410d2SHong Zhang 
415caf410d2SHong Zhang   PetscFunctionBegin;
416caf410d2SHong Zhang   *ne   = net->subnet[id1].nedge;
417caf410d2SHong Zhang   *edge = net->subnet[id1].edges;
418caf410d2SHong Zhang   PetscFunctionReturn(0);
419caf410d2SHong Zhang }
420caf410d2SHong Zhang 
421caf410d2SHong Zhang /*@C
4225f2c45f1SShri Abhyankar   DMNetworkRegisterComponent - Registers the network component
4235f2c45f1SShri Abhyankar 
4245f2c45f1SShri Abhyankar   Logically collective on DM
4255f2c45f1SShri Abhyankar 
4265f2c45f1SShri Abhyankar   Input Parameters
4275f2c45f1SShri Abhyankar + dm   - the network object
4285f2c45f1SShri Abhyankar . name - the component name
4295f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data
4305f2c45f1SShri Abhyankar 
4315f2c45f1SShri Abhyankar    Output Parameters
4325f2c45f1SShri Abhyankar .   key - an integer key that defines the component
4335f2c45f1SShri Abhyankar 
4345f2c45f1SShri Abhyankar    Notes
4355f2c45f1SShri Abhyankar    This routine should be called by all processors before calling DMNetworkLayoutSetup().
4365f2c45f1SShri Abhyankar 
4375f2c45f1SShri Abhyankar    Level: intermediate
4385f2c45f1SShri Abhyankar 
4395f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
4405f2c45f1SShri Abhyankar @*/
441caf410d2SHong Zhang PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,size_t size,PetscInt *key)
4425f2c45f1SShri Abhyankar {
4435f2c45f1SShri Abhyankar   PetscErrorCode        ierr;
4445f2c45f1SShri Abhyankar   DM_Network            *network = (DM_Network*) dm->data;
4455f2c45f1SShri Abhyankar   DMNetworkComponent    *component=&network->component[network->ncomponent];
4465f2c45f1SShri Abhyankar   PetscBool             flg=PETSC_FALSE;
4475f2c45f1SShri Abhyankar   PetscInt              i;
4485f2c45f1SShri Abhyankar 
4495f2c45f1SShri Abhyankar   PetscFunctionBegin;
4505f2c45f1SShri Abhyankar   for (i=0; i < network->ncomponent; i++) {
4515f2c45f1SShri Abhyankar     ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr);
4525f2c45f1SShri Abhyankar     if (flg) {
4535f2c45f1SShri Abhyankar       *key = i;
4545f2c45f1SShri Abhyankar       PetscFunctionReturn(0);
4555f2c45f1SShri Abhyankar     }
4566d64e262SShri Abhyankar   }
4576d64e262SShri Abhyankar   if(network->ncomponent == MAX_COMPONENTS) {
4586d64e262SShri Abhyankar     SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components registered exceeds the max %D",MAX_COMPONENTS);
4595f2c45f1SShri Abhyankar   }
4605f2c45f1SShri Abhyankar 
4615f2c45f1SShri Abhyankar   ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr);
4625f2c45f1SShri Abhyankar   component->size = size/sizeof(DMNetworkComponentGenericDataType);
4635f2c45f1SShri Abhyankar   *key = network->ncomponent;
4645f2c45f1SShri Abhyankar   network->ncomponent++;
4655f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4665f2c45f1SShri Abhyankar }
4675f2c45f1SShri Abhyankar 
4685f2c45f1SShri Abhyankar /*@
4695f2c45f1SShri Abhyankar   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.
4705f2c45f1SShri Abhyankar 
4715f2c45f1SShri Abhyankar   Not Collective
4725f2c45f1SShri Abhyankar 
4735f2c45f1SShri Abhyankar   Input Parameters:
4745f2c45f1SShri Abhyankar + dm - The DMNetwork object
4755f2c45f1SShri Abhyankar 
4765f2c45f1SShri Abhyankar   Output Paramters:
4775f2c45f1SShri Abhyankar + vStart - The first vertex point
4785f2c45f1SShri Abhyankar - vEnd   - One beyond the last vertex point
4795f2c45f1SShri Abhyankar 
4805f2c45f1SShri Abhyankar   Level: intermediate
4815f2c45f1SShri Abhyankar 
4825f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange
4835f2c45f1SShri Abhyankar @*/
4845f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
4855f2c45f1SShri Abhyankar {
4865f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4875f2c45f1SShri Abhyankar 
4885f2c45f1SShri Abhyankar   PetscFunctionBegin;
4895f2c45f1SShri Abhyankar   if (vStart) *vStart = network->vStart;
4905f2c45f1SShri Abhyankar   if (vEnd) *vEnd = network->vEnd;
4915f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4925f2c45f1SShri Abhyankar }
4935f2c45f1SShri Abhyankar 
4945f2c45f1SShri Abhyankar /*@
4955f2c45f1SShri Abhyankar   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.
4965f2c45f1SShri Abhyankar 
4975f2c45f1SShri Abhyankar   Not Collective
4985f2c45f1SShri Abhyankar 
4995f2c45f1SShri Abhyankar   Input Parameters:
5005f2c45f1SShri Abhyankar + dm - The DMNetwork object
5015f2c45f1SShri Abhyankar 
5025f2c45f1SShri Abhyankar   Output Paramters:
5035f2c45f1SShri Abhyankar + eStart - The first edge point
5045f2c45f1SShri Abhyankar - eEnd   - One beyond the last edge point
5055f2c45f1SShri Abhyankar 
5065f2c45f1SShri Abhyankar   Level: intermediate
5075f2c45f1SShri Abhyankar 
5085f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange
5095f2c45f1SShri Abhyankar @*/
5105f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
5115f2c45f1SShri Abhyankar {
5125f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5135f2c45f1SShri Abhyankar 
5145f2c45f1SShri Abhyankar   PetscFunctionBegin;
5155f2c45f1SShri Abhyankar   if (eStart) *eStart = network->eStart;
5165f2c45f1SShri Abhyankar   if (eEnd) *eEnd = network->eEnd;
5175f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5185f2c45f1SShri Abhyankar }
5195f2c45f1SShri Abhyankar 
5207b6afd5bSHong Zhang /*@
521e85e6aecSHong Zhang   DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge.
5227b6afd5bSHong Zhang 
5237b6afd5bSHong Zhang   Not Collective
5247b6afd5bSHong Zhang 
5257b6afd5bSHong Zhang   Input Parameters:
5267b6afd5bSHong Zhang + dm - DMNetwork object
527e85e6aecSHong Zhang - p  - edge point
5287b6afd5bSHong Zhang 
5297b6afd5bSHong Zhang   Output Paramters:
530e85e6aecSHong Zhang . index - user global numbering for the edge
5317b6afd5bSHong Zhang 
5327b6afd5bSHong Zhang   Level: intermediate
5337b6afd5bSHong Zhang 
534e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalVertexIndex
5357b6afd5bSHong Zhang @*/
536e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
5377b6afd5bSHong Zhang {
5387b6afd5bSHong Zhang   PetscErrorCode    ierr;
5397b6afd5bSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
5407b6afd5bSHong Zhang   PetscInt          offsetp;
5417b6afd5bSHong Zhang   DMNetworkComponentHeader header;
5427b6afd5bSHong Zhang 
5437b6afd5bSHong Zhang   PetscFunctionBegin;
544caf410d2SHong Zhang   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first");
5457b6afd5bSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
5467b6afd5bSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
547e85e6aecSHong Zhang   *index = header->index;
5487b6afd5bSHong Zhang   PetscFunctionReturn(0);
5497b6afd5bSHong Zhang }
5507b6afd5bSHong Zhang 
5515f2c45f1SShri Abhyankar /*@
552e85e6aecSHong Zhang   DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex.
553e85e6aecSHong Zhang 
554e85e6aecSHong Zhang   Not Collective
555e85e6aecSHong Zhang 
556e85e6aecSHong Zhang   Input Parameters:
557e85e6aecSHong Zhang + dm - DMNetwork object
558e85e6aecSHong Zhang - p  - vertex point
559e85e6aecSHong Zhang 
560e85e6aecSHong Zhang   Output Paramters:
561e85e6aecSHong Zhang . index - user global numbering for the vertex
562e85e6aecSHong Zhang 
563e85e6aecSHong Zhang   Level: intermediate
564e85e6aecSHong Zhang 
565e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalEdgeIndex
566e85e6aecSHong Zhang @*/
567e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
568e85e6aecSHong Zhang {
569e85e6aecSHong Zhang   PetscErrorCode    ierr;
570e85e6aecSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
571e85e6aecSHong Zhang   PetscInt          offsetp;
572e85e6aecSHong Zhang   DMNetworkComponentHeader header;
573e85e6aecSHong Zhang 
574e85e6aecSHong Zhang   PetscFunctionBegin;
575caf410d2SHong Zhang   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first");
576e85e6aecSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
577e85e6aecSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
578e85e6aecSHong Zhang   *index = header->index;
579e85e6aecSHong Zhang   PetscFunctionReturn(0);
580e85e6aecSHong Zhang }
581e85e6aecSHong Zhang 
582c3b11c7cSShri Abhyankar /*
583c3b11c7cSShri Abhyankar   DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the
584c3b11c7cSShri Abhyankar                                     component value from the component data array
585c3b11c7cSShri Abhyankar 
586c3b11c7cSShri Abhyankar   Not Collective
587c3b11c7cSShri Abhyankar 
588c3b11c7cSShri Abhyankar   Input Parameters:
589c3b11c7cSShri Abhyankar + dm      - The DMNetwork object
590c3b11c7cSShri Abhyankar . p       - vertex/edge point
591c3b11c7cSShri Abhyankar - compnum - component number
592c3b11c7cSShri Abhyankar 
593c3b11c7cSShri Abhyankar   Output Parameters:
594c3b11c7cSShri Abhyankar + compkey - the key obtained when registering the component
595c3b11c7cSShri Abhyankar - offset  - offset into the component data array associated with the vertex/edge point
596c3b11c7cSShri Abhyankar 
597c3b11c7cSShri Abhyankar   Notes:
598c3b11c7cSShri Abhyankar   Typical usage:
599c3b11c7cSShri Abhyankar 
600c3b11c7cSShri Abhyankar   DMNetworkGetComponentDataArray(dm, &arr);
601c3b11c7cSShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
602c3b11c7cSShri Abhyankar   Loop over vertices or edges
603c3b11c7cSShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
604c3b11c7cSShri Abhyankar     Loop over numcomps
605c3b11c7cSShri Abhyankar       DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset);
606c3b11c7cSShri Abhyankar       compdata = (UserCompDataType)(arr+offset);
607c3b11c7cSShri Abhyankar 
608c3b11c7cSShri Abhyankar   Level: intermediate
609c3b11c7cSShri Abhyankar 
610c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
611c3b11c7cSShri Abhyankar */
612c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
613c3b11c7cSShri Abhyankar {
614c3b11c7cSShri Abhyankar   PetscErrorCode           ierr;
615c3b11c7cSShri Abhyankar   PetscInt                 offsetp;
616c3b11c7cSShri Abhyankar   DMNetworkComponentHeader header;
617c3b11c7cSShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
618c3b11c7cSShri Abhyankar 
619c3b11c7cSShri Abhyankar   PetscFunctionBegin;
620c3b11c7cSShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
621c3b11c7cSShri Abhyankar   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
622c3b11c7cSShri Abhyankar   if (compkey) *compkey = header->key[compnum];
623c3b11c7cSShri Abhyankar   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
624c3b11c7cSShri Abhyankar   PetscFunctionReturn(0);
625c3b11c7cSShri Abhyankar }
626c3b11c7cSShri Abhyankar 
627c3b11c7cSShri Abhyankar /*@
628c3b11c7cSShri Abhyankar   DMNetworkGetComponent - Returns the network component and its key
629c3b11c7cSShri Abhyankar 
630c3b11c7cSShri Abhyankar   Not Collective
631c3b11c7cSShri Abhyankar 
632c3b11c7cSShri Abhyankar   Input Parameters
633c3b11c7cSShri Abhyankar + dm - DMNetwork object
634c3b11c7cSShri Abhyankar . p  - edge or vertex point
635c3b11c7cSShri Abhyankar - compnum - component number
636c3b11c7cSShri Abhyankar 
637c3b11c7cSShri Abhyankar   Output Parameters:
638c3b11c7cSShri Abhyankar + compkey - the key set for this computing during registration
639c3b11c7cSShri Abhyankar - component - the component data
640c3b11c7cSShri Abhyankar 
641c3b11c7cSShri Abhyankar   Notes:
642c3b11c7cSShri Abhyankar   Typical usage:
643c3b11c7cSShri Abhyankar 
644c3b11c7cSShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
645c3b11c7cSShri Abhyankar   Loop over vertices or edges
646c3b11c7cSShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
647c3b11c7cSShri Abhyankar     Loop over numcomps
648c3b11c7cSShri Abhyankar       DMNetworkGetComponent(dm,v,compnum,&key,&component);
649c3b11c7cSShri Abhyankar 
650c3b11c7cSShri Abhyankar   Level: intermediate
651c3b11c7cSShri Abhyankar 
652c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetVariableOffset
653c3b11c7cSShri Abhyankar @*/
654c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *key, void **component)
655c3b11c7cSShri Abhyankar {
656c3b11c7cSShri Abhyankar   PetscErrorCode ierr;
657c3b11c7cSShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
658e108cb99SStefano Zampini   PetscInt       offsetd = 0;
659c3b11c7cSShri Abhyankar 
660c3b11c7cSShri Abhyankar   PetscFunctionBegin;
661c3b11c7cSShri Abhyankar   ierr = DMNetworkGetComponentKeyOffset(dm,p,compnum,key,&offsetd);CHKERRQ(ierr);
662c3b11c7cSShri Abhyankar   *component = network->componentdataarray+offsetd;
663c3b11c7cSShri Abhyankar   PetscFunctionReturn(0);
664c3b11c7cSShri Abhyankar }
665c3b11c7cSShri Abhyankar 
666e85e6aecSHong Zhang /*@
667325661f6SSatish Balay   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)
6685f2c45f1SShri Abhyankar 
6695f2c45f1SShri Abhyankar   Not Collective
6705f2c45f1SShri Abhyankar 
6715f2c45f1SShri Abhyankar   Input Parameters:
6725f2c45f1SShri Abhyankar + dm           - The DMNetwork object
6735f2c45f1SShri Abhyankar . p            - vertex/edge point
6745f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component
6755f2c45f1SShri Abhyankar - compvalue    - pointer to the data structure for the component
6765f2c45f1SShri Abhyankar 
6775f2c45f1SShri Abhyankar   Level: intermediate
6785f2c45f1SShri Abhyankar 
6795f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
6805f2c45f1SShri Abhyankar @*/
6815f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
6825f2c45f1SShri Abhyankar {
6835f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
68443a39a44SBarry Smith   DMNetworkComponent       *component = &network->component[componentkey];
6855f2c45f1SShri Abhyankar   DMNetworkComponentHeader header = &network->header[p];
6865f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue = &network->cvalue[p];
6875f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
6885f2c45f1SShri Abhyankar 
6895f2c45f1SShri Abhyankar   PetscFunctionBegin;
690fa58f0a9SHong 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);
691fa58f0a9SHong Zhang 
69243a39a44SBarry Smith   header->size[header->ndata] = component->size;
69343a39a44SBarry Smith   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
6945f2c45f1SShri Abhyankar   header->key[header->ndata] = componentkey;
6955f2c45f1SShri Abhyankar   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
6965f2c45f1SShri Abhyankar 
6975f2c45f1SShri Abhyankar   cvalue->data[header->ndata] = (void*)compvalue;
6985f2c45f1SShri Abhyankar   header->ndata++;
6995f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7005f2c45f1SShri Abhyankar }
7015f2c45f1SShri Abhyankar 
7025f2c45f1SShri Abhyankar /*@
7035f2c45f1SShri Abhyankar   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
7045f2c45f1SShri Abhyankar 
7055f2c45f1SShri Abhyankar   Not Collective
7065f2c45f1SShri Abhyankar 
7075f2c45f1SShri Abhyankar   Input Parameters:
7085f2c45f1SShri Abhyankar + dm - The DMNetwork object
7095f2c45f1SShri Abhyankar . p  - vertex/edge point
7105f2c45f1SShri Abhyankar 
7115f2c45f1SShri Abhyankar   Output Parameters:
7125f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge
7135f2c45f1SShri Abhyankar 
7145f2c45f1SShri Abhyankar   Level: intermediate
7155f2c45f1SShri Abhyankar 
7165f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
7175f2c45f1SShri Abhyankar @*/
7185f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
7195f2c45f1SShri Abhyankar {
7205f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7215f2c45f1SShri Abhyankar   PetscInt       offset;
7225f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7235f2c45f1SShri Abhyankar 
7245f2c45f1SShri Abhyankar   PetscFunctionBegin;
7255f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
7265f2c45f1SShri Abhyankar   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
7275f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7285f2c45f1SShri Abhyankar }
7295f2c45f1SShri Abhyankar 
7305f2c45f1SShri Abhyankar /*@
7315f2c45f1SShri Abhyankar   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
7325f2c45f1SShri Abhyankar 
7335f2c45f1SShri Abhyankar   Not Collective
7345f2c45f1SShri Abhyankar 
7355f2c45f1SShri Abhyankar   Input Parameters:
7365f2c45f1SShri Abhyankar + dm     - The DMNetwork object
7375f2c45f1SShri Abhyankar - p      - the edge/vertex point
7385f2c45f1SShri Abhyankar 
7395f2c45f1SShri Abhyankar   Output Parameters:
7405f2c45f1SShri Abhyankar . offset - the offset
7415f2c45f1SShri Abhyankar 
7425f2c45f1SShri Abhyankar   Level: intermediate
7435f2c45f1SShri Abhyankar 
7445f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
7455f2c45f1SShri Abhyankar @*/
7465f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
7475f2c45f1SShri Abhyankar {
7485f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7495f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7505f2c45f1SShri Abhyankar 
7515f2c45f1SShri Abhyankar   PetscFunctionBegin;
7525f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr);
7535f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7545f2c45f1SShri Abhyankar }
7555f2c45f1SShri Abhyankar 
7565f2c45f1SShri Abhyankar /*@
7575f2c45f1SShri Abhyankar   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
7585f2c45f1SShri Abhyankar 
7595f2c45f1SShri Abhyankar   Not Collective
7605f2c45f1SShri Abhyankar 
7615f2c45f1SShri Abhyankar   Input Parameters:
7625f2c45f1SShri Abhyankar + dm      - The DMNetwork object
7635f2c45f1SShri Abhyankar - p       - the edge/vertex point
7645f2c45f1SShri Abhyankar 
7655f2c45f1SShri Abhyankar   Output Parameters:
7665f2c45f1SShri Abhyankar . offsetg - the offset
7675f2c45f1SShri Abhyankar 
7685f2c45f1SShri Abhyankar   Level: intermediate
7695f2c45f1SShri Abhyankar 
7705f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
7715f2c45f1SShri Abhyankar @*/
7725f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
7735f2c45f1SShri Abhyankar {
7745f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7755f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7765f2c45f1SShri Abhyankar 
7775f2c45f1SShri Abhyankar   PetscFunctionBegin;
7785f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr);
7796fefedf4SHong Zhang   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */
7805f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7815f2c45f1SShri Abhyankar }
7825f2c45f1SShri Abhyankar 
78324121865SAdrian Maldonado /*@
78424121865SAdrian Maldonado   DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector.
78524121865SAdrian Maldonado 
78624121865SAdrian Maldonado   Not Collective
78724121865SAdrian Maldonado 
78824121865SAdrian Maldonado   Input Parameters:
78924121865SAdrian Maldonado + dm     - The DMNetwork object
79024121865SAdrian Maldonado - p      - the edge point
79124121865SAdrian Maldonado 
79224121865SAdrian Maldonado   Output Parameters:
79324121865SAdrian Maldonado . offset - the offset
79424121865SAdrian Maldonado 
79524121865SAdrian Maldonado   Level: intermediate
79624121865SAdrian Maldonado 
79724121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
79824121865SAdrian Maldonado @*/
79924121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
80024121865SAdrian Maldonado {
80124121865SAdrian Maldonado   PetscErrorCode ierr;
80224121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
80324121865SAdrian Maldonado 
80424121865SAdrian Maldonado   PetscFunctionBegin;
80524121865SAdrian Maldonado 
80624121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr);
80724121865SAdrian Maldonado   PetscFunctionReturn(0);
80824121865SAdrian Maldonado }
80924121865SAdrian Maldonado 
81024121865SAdrian Maldonado /*@
81124121865SAdrian Maldonado   DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector.
81224121865SAdrian Maldonado 
81324121865SAdrian Maldonado   Not Collective
81424121865SAdrian Maldonado 
81524121865SAdrian Maldonado   Input Parameters:
81624121865SAdrian Maldonado + dm     - The DMNetwork object
81724121865SAdrian Maldonado - p      - the vertex point
81824121865SAdrian Maldonado 
81924121865SAdrian Maldonado   Output Parameters:
82024121865SAdrian Maldonado . offset - the offset
82124121865SAdrian Maldonado 
82224121865SAdrian Maldonado   Level: intermediate
82324121865SAdrian Maldonado 
82424121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
82524121865SAdrian Maldonado @*/
82624121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
82724121865SAdrian Maldonado {
82824121865SAdrian Maldonado   PetscErrorCode ierr;
82924121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
83024121865SAdrian Maldonado 
83124121865SAdrian Maldonado   PetscFunctionBegin;
83224121865SAdrian Maldonado 
83324121865SAdrian Maldonado   p -= network->vStart;
83424121865SAdrian Maldonado 
83524121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr);
83624121865SAdrian Maldonado   PetscFunctionReturn(0);
83724121865SAdrian Maldonado }
8385f2c45f1SShri Abhyankar /*@
8395f2c45f1SShri Abhyankar   DMNetworkAddNumVariables - Add number of variables associated with a given point.
8405f2c45f1SShri Abhyankar 
8415f2c45f1SShri Abhyankar   Not Collective
8425f2c45f1SShri Abhyankar 
8435f2c45f1SShri Abhyankar   Input Parameters:
8445f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
8455f2c45f1SShri Abhyankar . p    - the vertex/edge point
8465f2c45f1SShri Abhyankar - nvar - number of additional variables
8475f2c45f1SShri Abhyankar 
8485f2c45f1SShri Abhyankar   Level: intermediate
8495f2c45f1SShri Abhyankar 
8505f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables
8515f2c45f1SShri Abhyankar @*/
8525f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
8535f2c45f1SShri Abhyankar {
8545f2c45f1SShri Abhyankar   PetscErrorCode ierr;
8555f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
8565f2c45f1SShri Abhyankar 
8575f2c45f1SShri Abhyankar   PetscFunctionBegin;
8585f2c45f1SShri Abhyankar   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
8595f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
8605f2c45f1SShri Abhyankar }
8615f2c45f1SShri Abhyankar 
86227f51fceSHong Zhang /*@
86327f51fceSHong Zhang   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
86427f51fceSHong Zhang 
86527f51fceSHong Zhang   Not Collective
86627f51fceSHong Zhang 
86727f51fceSHong Zhang   Input Parameters:
86827f51fceSHong Zhang + dm   - The DMNetworkObject
86927f51fceSHong Zhang - p    - the vertex/edge point
87027f51fceSHong Zhang 
87127f51fceSHong Zhang   Output Parameters:
87227f51fceSHong Zhang . nvar - number of variables
87327f51fceSHong Zhang 
87427f51fceSHong Zhang   Level: intermediate
87527f51fceSHong Zhang 
87627f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
87727f51fceSHong Zhang @*/
87827f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
87927f51fceSHong Zhang {
88027f51fceSHong Zhang   PetscErrorCode ierr;
88127f51fceSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
88227f51fceSHong Zhang 
88327f51fceSHong Zhang   PetscFunctionBegin;
88427f51fceSHong Zhang   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
88527f51fceSHong Zhang   PetscFunctionReturn(0);
88627f51fceSHong Zhang }
88727f51fceSHong Zhang 
8885f2c45f1SShri Abhyankar /*@
8895f2c45f1SShri Abhyankar   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
8905f2c45f1SShri Abhyankar 
8915f2c45f1SShri Abhyankar   Not Collective
8925f2c45f1SShri Abhyankar 
8935f2c45f1SShri Abhyankar   Input Parameters:
8945f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
8955f2c45f1SShri Abhyankar . p    - the vertex/edge point
8965f2c45f1SShri Abhyankar - nvar - number of variables
8975f2c45f1SShri Abhyankar 
8985f2c45f1SShri Abhyankar   Level: intermediate
8995f2c45f1SShri Abhyankar 
9005f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables
9015f2c45f1SShri Abhyankar @*/
9025f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
9035f2c45f1SShri Abhyankar {
9045f2c45f1SShri Abhyankar   PetscErrorCode ierr;
9055f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
9065f2c45f1SShri Abhyankar 
9075f2c45f1SShri Abhyankar   PetscFunctionBegin;
9085f2c45f1SShri Abhyankar   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
9095f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9105f2c45f1SShri Abhyankar }
9115f2c45f1SShri Abhyankar 
9125f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This
9135f2c45f1SShri Abhyankar    function is called during DMSetUp() */
9145f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm)
9155f2c45f1SShri Abhyankar {
9165f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
9175f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
918e53b5ba3SHong Zhang   PetscInt                 arr_size,p,offset,offsetp,ncomp,i;
9195f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
9205f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue;
9215f2c45f1SShri Abhyankar   DMNetworkComponentGenericDataType *componentdataarray;
9225f2c45f1SShri Abhyankar 
9235f2c45f1SShri Abhyankar   PetscFunctionBegin;
9245f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
9255f2c45f1SShri Abhyankar   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
92675b160a0SShri Abhyankar   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
9275f2c45f1SShri Abhyankar   componentdataarray = network->componentdataarray;
9285f2c45f1SShri Abhyankar   for (p = network->pStart; p < network->pEnd; p++) {
9295f2c45f1SShri Abhyankar     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
9305f2c45f1SShri Abhyankar     /* Copy header */
9315f2c45f1SShri Abhyankar     header = &network->header[p];
932302440fdSBarry Smith     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
9335f2c45f1SShri Abhyankar     /* Copy data */
9345f2c45f1SShri Abhyankar     cvalue = &network->cvalue[p];
9355f2c45f1SShri Abhyankar     ncomp = header->ndata;
9365f2c45f1SShri Abhyankar     for (i = 0; i < ncomp; i++) {
9375f2c45f1SShri Abhyankar       offset = offsetp + network->dataheadersize + header->offset[i];
938302440fdSBarry Smith       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
9395f2c45f1SShri Abhyankar     }
9405f2c45f1SShri Abhyankar   }
9415f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9425f2c45f1SShri Abhyankar }
9435f2c45f1SShri Abhyankar 
9445f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */
9455f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm)
9465f2c45f1SShri Abhyankar {
9475f2c45f1SShri Abhyankar   PetscErrorCode ierr;
9485f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
9495f2c45f1SShri Abhyankar 
9505f2c45f1SShri Abhyankar   PetscFunctionBegin;
9515f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
9525f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9535f2c45f1SShri Abhyankar }
9545f2c45f1SShri Abhyankar 
9555f2c45f1SShri Abhyankar /*@C
9565f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray - Returns the component data array
9575f2c45f1SShri Abhyankar 
9585f2c45f1SShri Abhyankar   Not Collective
9595f2c45f1SShri Abhyankar 
9605f2c45f1SShri Abhyankar   Input Parameters:
9615f2c45f1SShri Abhyankar . dm - The DMNetwork Object
9625f2c45f1SShri Abhyankar 
9635f2c45f1SShri Abhyankar   Output Parameters:
9645f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components
9655f2c45f1SShri Abhyankar 
9665f2c45f1SShri Abhyankar   Level: intermediate
9675f2c45f1SShri Abhyankar 
968a730d845SHong Zhang .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents
9695f2c45f1SShri Abhyankar @*/
9705f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
9715f2c45f1SShri Abhyankar {
9725f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
9735f2c45f1SShri Abhyankar 
9745f2c45f1SShri Abhyankar   PetscFunctionBegin;
9755f2c45f1SShri Abhyankar   *componentdataarray = network->componentdataarray;
9765f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9775f2c45f1SShri Abhyankar }
9785f2c45f1SShri Abhyankar 
97924121865SAdrian Maldonado /* Get a subsection from a range of points */
98024121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection)
98124121865SAdrian Maldonado {
98224121865SAdrian Maldonado   PetscErrorCode ierr;
98324121865SAdrian Maldonado   PetscInt       i, nvar;
98424121865SAdrian Maldonado 
98524121865SAdrian Maldonado   PetscFunctionBegin;
98624121865SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr);
98724121865SAdrian Maldonado   ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr);
98824121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
98924121865SAdrian Maldonado     ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr);
99024121865SAdrian Maldonado     ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr);
99124121865SAdrian Maldonado   }
99224121865SAdrian Maldonado 
99324121865SAdrian Maldonado   ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr);
99424121865SAdrian Maldonado   PetscFunctionReturn(0);
99524121865SAdrian Maldonado }
99624121865SAdrian Maldonado 
99724121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */
99824121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
99924121865SAdrian Maldonado {
100024121865SAdrian Maldonado   PetscErrorCode ierr;
100124121865SAdrian Maldonado   PetscInt       i, *subpoints;
100224121865SAdrian Maldonado 
100324121865SAdrian Maldonado   PetscFunctionBegin;
100424121865SAdrian Maldonado   /* Create index sets to map from "points" to "subpoints" */
100524121865SAdrian Maldonado   ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr);
100624121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
100724121865SAdrian Maldonado     subpoints[i - pstart] = i;
100824121865SAdrian Maldonado   }
1009459726d8SSatish Balay   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr);
101024121865SAdrian Maldonado   ierr = PetscFree(subpoints);CHKERRQ(ierr);
101124121865SAdrian Maldonado   PetscFunctionReturn(0);
101224121865SAdrian Maldonado }
101324121865SAdrian Maldonado 
101424121865SAdrian Maldonado /*@
101524121865SAdrian Maldonado   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute.
101624121865SAdrian Maldonado 
101724121865SAdrian Maldonado   Collective
101824121865SAdrian Maldonado 
101924121865SAdrian Maldonado   Input Parameters:
102024121865SAdrian Maldonado . dm   - The DMNetworkObject
102124121865SAdrian Maldonado 
102224121865SAdrian Maldonado   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
102324121865SAdrian Maldonado 
102424121865SAdrian Maldonado   points = [0 1 2 3 4 5 6]
102524121865SAdrian Maldonado 
1026caf410d2SHong Zhang   where edges = [0,1,2,3] and vertices = [4,5,6]. The new orderings will be specific to the subset (i.e vertices = [0,1,2] <- [4,5,6]).
102724121865SAdrian Maldonado 
102824121865SAdrian Maldonado   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
102924121865SAdrian Maldonado 
103024121865SAdrian Maldonado   Level: intermediate
103124121865SAdrian Maldonado 
103224121865SAdrian Maldonado @*/
103324121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
103424121865SAdrian Maldonado {
103524121865SAdrian Maldonado   PetscErrorCode ierr;
103624121865SAdrian Maldonado   MPI_Comm       comm;
10379852e123SBarry Smith   PetscMPIInt    rank, size;
103824121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
103924121865SAdrian Maldonado 
1040eab1376dSHong Zhang   PetscFunctionBegin;
104124121865SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
104224121865SAdrian Maldonado   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
10439852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
104424121865SAdrian Maldonado 
104524121865SAdrian Maldonado   /* Create maps for vertices and edges */
104624121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr);
104724121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr);
104824121865SAdrian Maldonado 
104924121865SAdrian Maldonado   /* Create local sub-sections */
105024121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr);
105124121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr);
105224121865SAdrian Maldonado 
10539852e123SBarry Smith   if (size > 1) {
105424121865SAdrian Maldonado     ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr);
105524121865SAdrian Maldonado     ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr);
105624121865SAdrian Maldonado   ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr);
105724121865SAdrian Maldonado   ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr);
105824121865SAdrian Maldonado   } else {
105924121865SAdrian Maldonado   /* create structures for vertex */
106024121865SAdrian Maldonado   ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr);
106124121865SAdrian Maldonado   /* create structures for edge */
106224121865SAdrian Maldonado   ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr);
106324121865SAdrian Maldonado   }
106424121865SAdrian Maldonado 
106524121865SAdrian Maldonado 
106624121865SAdrian Maldonado   /* Add viewers */
106724121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr);
106824121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr);
106924121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr);
107024121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr);
107124121865SAdrian Maldonado 
107224121865SAdrian Maldonado   PetscFunctionReturn(0);
107324121865SAdrian Maldonado }
10747b6afd5bSHong Zhang 
10755f2c45f1SShri Abhyankar /*@
10765f2c45f1SShri Abhyankar   DMNetworkDistribute - Distributes the network and moves associated component data.
10775f2c45f1SShri Abhyankar 
10785f2c45f1SShri Abhyankar   Collective
10795f2c45f1SShri Abhyankar 
10805f2c45f1SShri Abhyankar   Input Parameter:
1081d3464fd4SAdrian Maldonado + DM - the DMNetwork object
10825f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default
10835f2c45f1SShri Abhyankar 
10845f2c45f1SShri Abhyankar   Notes:
10858b171c8eSHong Zhang   Distributes the network with <overlap>-overlapping partitioning of the edges.
10865f2c45f1SShri Abhyankar 
10875f2c45f1SShri Abhyankar   Level: intermediate
10885f2c45f1SShri Abhyankar 
10895f2c45f1SShri Abhyankar .seealso: DMNetworkCreate
10905f2c45f1SShri Abhyankar @*/
1091d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
10925f2c45f1SShri Abhyankar {
1093d3464fd4SAdrian Maldonado   MPI_Comm       comm;
10945f2c45f1SShri Abhyankar   PetscErrorCode ierr;
1095d3464fd4SAdrian Maldonado   PetscMPIInt    size;
1096d3464fd4SAdrian Maldonado   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
1097d3464fd4SAdrian Maldonado   DM_Network     *newDMnetwork;
1098caf410d2SHong Zhang   PetscSF        pointsf=NULL;
10995f2c45f1SShri Abhyankar   DM             newDM;
1100caf410d2SHong Zhang   PetscInt       j,e,v,offset,*subnetvtx;
110151ac5effSHong Zhang   PetscPartitioner         part;
1102b9c6e19dSShri Abhyankar   DMNetworkComponentHeader header;
11035f2c45f1SShri Abhyankar 
11045f2c45f1SShri Abhyankar   PetscFunctionBegin;
1105d3464fd4SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr);
1106d3464fd4SAdrian Maldonado   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1107d3464fd4SAdrian Maldonado   if (size == 1) PetscFunctionReturn(0);
1108d3464fd4SAdrian Maldonado 
1109d3464fd4SAdrian Maldonado   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr);
11105f2c45f1SShri Abhyankar   newDMnetwork = (DM_Network*)newDM->data;
11115f2c45f1SShri Abhyankar   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
111251ac5effSHong Zhang 
111351ac5effSHong Zhang   /* Enable runtime options for petscpartitioner */
111451ac5effSHong Zhang   ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr);
111551ac5effSHong Zhang   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
111651ac5effSHong Zhang 
11175f2c45f1SShri Abhyankar   /* Distribute plex dm and dof section */
111880cf41d5SMatthew G. Knepley   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
111951ac5effSHong Zhang 
11205f2c45f1SShri Abhyankar   /* Distribute dof section */
1121d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr);
11225f2c45f1SShri Abhyankar   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
1123d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr);
112451ac5effSHong Zhang 
11255f2c45f1SShri Abhyankar   /* Distribute data and associated section */
112631da1fc8SHong Zhang   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
112724121865SAdrian Maldonado 
11285f2c45f1SShri Abhyankar   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
11295f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
11305f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
11315f2c45f1SShri Abhyankar   newDMnetwork->nEdges    = newDMnetwork->eEnd - newDMnetwork->eStart;
11326fefedf4SHong Zhang   newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
11336fefedf4SHong Zhang   newDMnetwork->NVertices = oldDMnetwork->NVertices;
11345f2c45f1SShri Abhyankar   newDMnetwork->NEdges    = oldDMnetwork->NEdges;
113524121865SAdrian Maldonado 
11365f2c45f1SShri Abhyankar   /* Set Dof section as the default section for dm */
1137e87a4003SBarry Smith   ierr = DMSetSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
1138e87a4003SBarry Smith   ierr = DMGetGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
11395f2c45f1SShri Abhyankar 
1140b9c6e19dSShri Abhyankar   /* Set up subnetwork info in the newDM */
1141b9c6e19dSShri Abhyankar   newDMnetwork->nsubnet  = oldDMnetwork->nsubnet;
1142caf410d2SHong Zhang   newDMnetwork->ncsubnet = oldDMnetwork->ncsubnet;
1143b9c6e19dSShri Abhyankar   ierr = PetscCalloc1(newDMnetwork->nsubnet,&newDMnetwork->subnet);CHKERRQ(ierr);
1144b9c6e19dSShri Abhyankar   /* Copy over the global number of vertices and edges in each subnetwork. Note that these are already
1145b9c6e19dSShri Abhyankar      calculated in DMNetworkLayoutSetUp()
1146b9c6e19dSShri Abhyankar   */
1147b9c6e19dSShri Abhyankar   for(j=0; j < newDMnetwork->nsubnet; j++) {
1148b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].Nvtx  = oldDMnetwork->subnet[j].Nvtx;
1149b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge;
1150b9c6e19dSShri Abhyankar   }
1151b9c6e19dSShri Abhyankar 
1152b9c6e19dSShri Abhyankar   for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1153b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
1154b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1155b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].nedge++;
1156b9c6e19dSShri Abhyankar   }
1157b9c6e19dSShri Abhyankar 
1158b9c6e19dSShri Abhyankar   for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1159b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
1160b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1161b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].nvtx++;
1162b9c6e19dSShri Abhyankar   }
1163b9c6e19dSShri Abhyankar 
1164b9c6e19dSShri Abhyankar   /* Now create the vertices and edge arrays for the subnetworks */
1165caf410d2SHong Zhang   ierr = PetscCalloc1(newDMnetwork->vEnd-newDMnetwork->vStart,&newDMnetwork->subnetvtx);CHKERRQ(ierr);
1166caf410d2SHong Zhang   subnetvtx = newDMnetwork->subnetvtx;
1167caf410d2SHong Zhang 
1168b9c6e19dSShri Abhyankar   for (j=0; j<newDMnetwork->nsubnet; j++) {
1169b9c6e19dSShri Abhyankar     ierr = PetscCalloc1(newDMnetwork->subnet[j].nedge,&newDMnetwork->subnet[j].edges);CHKERRQ(ierr);
1170caf410d2SHong Zhang     newDMnetwork->subnet[j].vertices = subnetvtx;
1171caf410d2SHong Zhang     subnetvtx                       += newDMnetwork->subnet[j].nvtx;
1172caf410d2SHong Zhang 
1173b9c6e19dSShri Abhyankar     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
1174b9c6e19dSShri Abhyankar        These get updated when the vertices and edges are added. */
1175b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0;
1176b9c6e19dSShri Abhyankar   }
1177b9c6e19dSShri Abhyankar 
1178b9c6e19dSShri Abhyankar   /* Set the vertices and edges in each subnetwork */
1179b9c6e19dSShri Abhyankar   for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1180b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
1181b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1182b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++] = e;
1183b9c6e19dSShri Abhyankar   }
1184b9c6e19dSShri Abhyankar 
1185b9c6e19dSShri Abhyankar   for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1186b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
1187b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1188b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++] = v;
1189b9c6e19dSShri Abhyankar   }
1190b9c6e19dSShri Abhyankar 
1191caf410d2SHong Zhang   newDM->setupcalled = (*dm)->setupcalled;
1192caf410d2SHong Zhang 
119324121865SAdrian Maldonado   /* Destroy point SF */
119424121865SAdrian Maldonado   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
119524121865SAdrian Maldonado 
1196d3464fd4SAdrian Maldonado   ierr = DMDestroy(dm);CHKERRQ(ierr);
1197d3464fd4SAdrian Maldonado   *dm  = newDM;
11985f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
11995f2c45f1SShri Abhyankar }
12005f2c45f1SShri Abhyankar 
120124121865SAdrian Maldonado /*@C
120224121865SAdrian Maldonado   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering.
120324121865SAdrian Maldonado 
120424121865SAdrian Maldonado   Input Parameters:
120524121865SAdrian Maldonado + masterSF - the original SF structure
120624121865SAdrian Maldonado - map      - a ISLocalToGlobal mapping that contains the subset of points
120724121865SAdrian Maldonado 
120824121865SAdrian Maldonado   Output Parameters:
120924121865SAdrian Maldonado . subSF    - a subset of the masterSF for the desired subset.
121024121865SAdrian Maldonado */
121124121865SAdrian Maldonado 
121224121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) {
121324121865SAdrian Maldonado 
121424121865SAdrian Maldonado   PetscErrorCode        ierr;
121524121865SAdrian Maldonado   PetscInt              nroots, nleaves, *ilocal_sub;
121624121865SAdrian Maldonado   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
121724121865SAdrian Maldonado   PetscInt              *local_points, *remote_points;
121824121865SAdrian Maldonado   PetscSFNode           *iremote_sub;
121924121865SAdrian Maldonado   const PetscInt        *ilocal;
122024121865SAdrian Maldonado   const PetscSFNode     *iremote;
122124121865SAdrian Maldonado 
122224121865SAdrian Maldonado   PetscFunctionBegin;
122324121865SAdrian Maldonado   ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
122424121865SAdrian Maldonado 
122524121865SAdrian Maldonado   /* Look for leaves that pertain to the subset of points. Get the local ordering */
122624121865SAdrian Maldonado   ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr);
122724121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr);
122824121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
122924121865SAdrian Maldonado     if (ilocal_map[i] != -1) nleaves_sub += 1;
123024121865SAdrian Maldonado   }
123124121865SAdrian Maldonado   /* Re-number ilocal with subset numbering. Need information from roots */
123224121865SAdrian Maldonado   ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr);
123324121865SAdrian Maldonado   for (i = 0; i < nroots; i++) local_points[i] = i;
123424121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr);
123524121865SAdrian Maldonado   ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
123624121865SAdrian Maldonado   ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
123724121865SAdrian Maldonado   /* Fill up graph using local (that is, local to the subset) numbering. */
12384b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr);
12394b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr);
124024121865SAdrian Maldonado   nleaves_sub = 0;
124124121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
124224121865SAdrian Maldonado     if (ilocal_map[i] != -1) {
124324121865SAdrian Maldonado       ilocal_sub[nleaves_sub] = ilocal_map[i];
12444b70a8deSAdrian Maldonado       iremote_sub[nleaves_sub].rank = iremote[i].rank;
124524121865SAdrian Maldonado       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
124624121865SAdrian Maldonado       nleaves_sub += 1;
124724121865SAdrian Maldonado     }
124824121865SAdrian Maldonado   }
124924121865SAdrian Maldonado   ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr);
125024121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr);
125124121865SAdrian Maldonado 
125224121865SAdrian Maldonado   /* Create new subSF */
125324121865SAdrian Maldonado   ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr);
125424121865SAdrian Maldonado   ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr);
12554b70a8deSAdrian Maldonado   ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr);
125624121865SAdrian Maldonado   ierr = PetscFree(ilocal_map);CHKERRQ(ierr);
12574b70a8deSAdrian Maldonado   ierr = PetscFree(iremote_sub);CHKERRQ(ierr);
125824121865SAdrian Maldonado   PetscFunctionReturn(0);
125924121865SAdrian Maldonado }
126024121865SAdrian Maldonado 
12615f2c45f1SShri Abhyankar /*@C
12625f2c45f1SShri Abhyankar   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
12635f2c45f1SShri Abhyankar 
12645f2c45f1SShri Abhyankar   Not Collective
12655f2c45f1SShri Abhyankar 
12665f2c45f1SShri Abhyankar   Input Parameters:
12675f2c45f1SShri Abhyankar + dm - The DMNetwork object
12685f2c45f1SShri Abhyankar - p  - the vertex point
12695f2c45f1SShri Abhyankar 
12705f2c45f1SShri Abhyankar   Output Paramters:
12715f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point
12725f2c45f1SShri Abhyankar - edges  - List of edge points
12735f2c45f1SShri Abhyankar 
12745f2c45f1SShri Abhyankar   Level: intermediate
12755f2c45f1SShri Abhyankar 
12765f2c45f1SShri Abhyankar   Fortran Notes:
12775f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
12785f2c45f1SShri Abhyankar   include petsc.h90 in your code.
12795f2c45f1SShri Abhyankar 
1280d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices
12815f2c45f1SShri Abhyankar @*/
12825f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
12835f2c45f1SShri Abhyankar {
12845f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12855f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
12865f2c45f1SShri Abhyankar 
12875f2c45f1SShri Abhyankar   PetscFunctionBegin;
12885f2c45f1SShri Abhyankar   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
12895f2c45f1SShri Abhyankar   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
12905f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12915f2c45f1SShri Abhyankar }
12925f2c45f1SShri Abhyankar 
12935f2c45f1SShri Abhyankar /*@C
1294d842c372SHong Zhang   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
12955f2c45f1SShri Abhyankar 
12965f2c45f1SShri Abhyankar   Not Collective
12975f2c45f1SShri Abhyankar 
12985f2c45f1SShri Abhyankar   Input Parameters:
12995f2c45f1SShri Abhyankar + dm - The DMNetwork object
13005f2c45f1SShri Abhyankar - p  - the edge point
13015f2c45f1SShri Abhyankar 
13025f2c45f1SShri Abhyankar   Output Paramters:
13035f2c45f1SShri Abhyankar . vertices  - vertices connected to this edge
13045f2c45f1SShri Abhyankar 
13055f2c45f1SShri Abhyankar   Level: intermediate
13065f2c45f1SShri Abhyankar 
13075f2c45f1SShri Abhyankar   Fortran Notes:
13085f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
13095f2c45f1SShri Abhyankar   include petsc.h90 in your code.
13105f2c45f1SShri Abhyankar 
13115f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
13125f2c45f1SShri Abhyankar @*/
1313d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
13145f2c45f1SShri Abhyankar {
13155f2c45f1SShri Abhyankar   PetscErrorCode ierr;
13165f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
13175f2c45f1SShri Abhyankar 
13185f2c45f1SShri Abhyankar   PetscFunctionBegin;
13195f2c45f1SShri Abhyankar   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
13205f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
13215f2c45f1SShri Abhyankar }
13225f2c45f1SShri Abhyankar 
13235f2c45f1SShri Abhyankar /*@
13245f2c45f1SShri Abhyankar   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
13255f2c45f1SShri Abhyankar 
13265f2c45f1SShri Abhyankar   Not Collective
13275f2c45f1SShri Abhyankar 
13285f2c45f1SShri Abhyankar   Input Parameters:
13295f2c45f1SShri Abhyankar + dm - The DMNetwork object
13305f2c45f1SShri Abhyankar . p  - the vertex point
13315f2c45f1SShri Abhyankar 
13325f2c45f1SShri Abhyankar   Output Parameter:
13335f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point
13345f2c45f1SShri Abhyankar 
13355f2c45f1SShri Abhyankar   Level: intermediate
13365f2c45f1SShri Abhyankar 
1337d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange
13385f2c45f1SShri Abhyankar @*/
13395f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
13405f2c45f1SShri Abhyankar {
13415f2c45f1SShri Abhyankar   PetscErrorCode ierr;
13425f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
13435f2c45f1SShri Abhyankar   PetscInt       offsetg;
13445f2c45f1SShri Abhyankar   PetscSection   sectiong;
13455f2c45f1SShri Abhyankar 
13465f2c45f1SShri Abhyankar   PetscFunctionBegin;
1347caf410d2SHong Zhang   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first");
13485f2c45f1SShri Abhyankar   *isghost = PETSC_FALSE;
1349e87a4003SBarry Smith   ierr = DMGetGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
13505f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
13515f2c45f1SShri Abhyankar   if (offsetg < 0) *isghost = PETSC_TRUE;
13525f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
13535f2c45f1SShri Abhyankar }
13545f2c45f1SShri Abhyankar 
13555f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm)
13565f2c45f1SShri Abhyankar {
13575f2c45f1SShri Abhyankar   PetscErrorCode ierr;
13585f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
13595f2c45f1SShri Abhyankar 
13605f2c45f1SShri Abhyankar   PetscFunctionBegin;
13615f2c45f1SShri Abhyankar   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
13625f2c45f1SShri Abhyankar   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
13635f2c45f1SShri Abhyankar 
1364e87a4003SBarry Smith   ierr = DMSetSection(network->plex,network->DofSection);CHKERRQ(ierr);
1365e87a4003SBarry Smith   ierr = DMGetGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
13665f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
13675f2c45f1SShri Abhyankar }
13685f2c45f1SShri Abhyankar 
13691ad426b7SHong Zhang /*@
137017df6e9eSHong Zhang     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
13711ad426b7SHong Zhang                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
13721ad426b7SHong Zhang 
13731ad426b7SHong Zhang     Collective
13741ad426b7SHong Zhang 
13751ad426b7SHong Zhang     Input Parameters:
137683b2e829SHong Zhang +   dm - The DMNetwork object
137783b2e829SHong Zhang .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
137883b2e829SHong Zhang -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
13791ad426b7SHong Zhang 
13801ad426b7SHong Zhang     Level: intermediate
13811ad426b7SHong Zhang 
13821ad426b7SHong Zhang @*/
138383b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
13841ad426b7SHong Zhang {
13851ad426b7SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
13868675203cSHong Zhang   PetscErrorCode ierr;
138766b4e0ffSHong Zhang   PetscInt       nVertices = network->nVertices;
13881ad426b7SHong Zhang 
13891ad426b7SHong Zhang   PetscFunctionBegin;
139083b2e829SHong Zhang   network->userEdgeJacobian   = eflg;
139183b2e829SHong Zhang   network->userVertexJacobian = vflg;
13928675203cSHong Zhang 
13938675203cSHong Zhang   if (eflg && !network->Je) {
13948675203cSHong Zhang     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
13958675203cSHong Zhang   }
13968675203cSHong Zhang 
139766b4e0ffSHong Zhang   if (vflg && !network->Jv && nVertices) {
13988675203cSHong Zhang     PetscInt       i,*vptr,nedges,vStart=network->vStart;
139966b4e0ffSHong Zhang     PetscInt       nedges_total;
14008675203cSHong Zhang     const PetscInt *edges;
14018675203cSHong Zhang 
14028675203cSHong Zhang     /* count nvertex_total */
14038675203cSHong Zhang     nedges_total = 0;
14048675203cSHong Zhang     ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr);
14058675203cSHong Zhang 
14068675203cSHong Zhang     vptr[0] = 0;
14078675203cSHong Zhang     for (i=0; i<nVertices; i++) {
14088675203cSHong Zhang       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
14098675203cSHong Zhang       nedges_total += nedges;
14108675203cSHong Zhang       vptr[i+1] = vptr[i] + 2*nedges + 1;
14118675203cSHong Zhang     }
14128675203cSHong Zhang 
14138675203cSHong Zhang     ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr);
14148675203cSHong Zhang     network->Jvptr = vptr;
14158675203cSHong Zhang   }
14161ad426b7SHong Zhang   PetscFunctionReturn(0);
14171ad426b7SHong Zhang }
14181ad426b7SHong Zhang 
14191ad426b7SHong Zhang /*@
142083b2e829SHong Zhang     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
142183b2e829SHong Zhang 
142283b2e829SHong Zhang     Not Collective
142383b2e829SHong Zhang 
142483b2e829SHong Zhang     Input Parameters:
142583b2e829SHong Zhang +   dm - The DMNetwork object
142683b2e829SHong Zhang .   p  - the edge point
14273e97b6e8SHong Zhang -   J - array (size = 3) of Jacobian submatrices for this edge point:
14283e97b6e8SHong Zhang         J[0]: this edge
1429d842c372SHong Zhang         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()
143083b2e829SHong Zhang 
143183b2e829SHong Zhang     Level: intermediate
143283b2e829SHong Zhang 
143383b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix
143483b2e829SHong Zhang @*/
143583b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
143683b2e829SHong Zhang {
143783b2e829SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
143883b2e829SHong Zhang 
143983b2e829SHong Zhang   PetscFunctionBegin;
14408675203cSHong Zhang   if (!network->Je) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
14418675203cSHong Zhang 
14428675203cSHong Zhang   if (J) {
1443883e35e8SHong Zhang     network->Je[3*p]   = J[0];
1444883e35e8SHong Zhang     network->Je[3*p+1] = J[1];
1445883e35e8SHong Zhang     network->Je[3*p+2] = J[2];
14468675203cSHong Zhang   }
144783b2e829SHong Zhang   PetscFunctionReturn(0);
144883b2e829SHong Zhang }
144983b2e829SHong Zhang 
145083b2e829SHong Zhang /*@
145176ddfea5SHong Zhang     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
14521ad426b7SHong Zhang 
14531ad426b7SHong Zhang     Not Collective
14541ad426b7SHong Zhang 
14551ad426b7SHong Zhang     Input Parameters:
14561ad426b7SHong Zhang +   dm - The DMNetwork object
14571ad426b7SHong Zhang .   p  - the vertex point
14583e97b6e8SHong Zhang -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
14593e97b6e8SHong Zhang         J[0]:       this vertex
14603e97b6e8SHong Zhang         J[1+2*i]:   i-th supporting edge
14613e97b6e8SHong Zhang         J[1+2*i+1]: i-th connected vertex
14621ad426b7SHong Zhang 
14631ad426b7SHong Zhang     Level: intermediate
14641ad426b7SHong Zhang 
146583b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix
14661ad426b7SHong Zhang @*/
1467883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
14685f2c45f1SShri Abhyankar {
14695f2c45f1SShri Abhyankar   PetscErrorCode ierr;
14705f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
14718675203cSHong Zhang   PetscInt       i,*vptr,nedges,vStart=network->vStart;
1472883e35e8SHong Zhang   const PetscInt *edges;
14735f2c45f1SShri Abhyankar 
14745f2c45f1SShri Abhyankar   PetscFunctionBegin;
14758675203cSHong Zhang   if (!network->Jv) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
1476883e35e8SHong Zhang 
14778675203cSHong Zhang   if (J) {
1478883e35e8SHong Zhang     vptr = network->Jvptr;
14793e97b6e8SHong Zhang     network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
14803e97b6e8SHong Zhang 
14813e97b6e8SHong Zhang     /* Set Jacobian for each supporting edge and connected vertex */
1482883e35e8SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
1483883e35e8SHong Zhang     for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
14848675203cSHong Zhang   }
1485883e35e8SHong Zhang   PetscFunctionReturn(0);
1486883e35e8SHong Zhang }
1487883e35e8SHong Zhang 
1488e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
14895cf7da58SHong Zhang {
14905cf7da58SHong Zhang   PetscErrorCode ierr;
14915cf7da58SHong Zhang   PetscInt       j;
14925cf7da58SHong Zhang   PetscScalar    val=(PetscScalar)ncols;
14935cf7da58SHong Zhang 
14945cf7da58SHong Zhang   PetscFunctionBegin;
14955cf7da58SHong Zhang   if (!ghost) {
14965cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
14975cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
14985cf7da58SHong Zhang     }
14995cf7da58SHong Zhang   } else {
15005cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
15015cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
15025cf7da58SHong Zhang     }
15035cf7da58SHong Zhang   }
15045cf7da58SHong Zhang   PetscFunctionReturn(0);
15055cf7da58SHong Zhang }
15065cf7da58SHong Zhang 
1507e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
15085cf7da58SHong Zhang {
15095cf7da58SHong Zhang   PetscErrorCode ierr;
15105cf7da58SHong Zhang   PetscInt       j,ncols_u;
15115cf7da58SHong Zhang   PetscScalar    val;
15125cf7da58SHong Zhang 
15135cf7da58SHong Zhang   PetscFunctionBegin;
15145cf7da58SHong Zhang   if (!ghost) {
15155cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
15165cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
15175cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
15185cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
15195cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
15205cf7da58SHong Zhang     }
15215cf7da58SHong Zhang   } else {
15225cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
15235cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
15245cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
15255cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
15265cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
15275cf7da58SHong Zhang     }
15285cf7da58SHong Zhang   }
15295cf7da58SHong Zhang   PetscFunctionReturn(0);
15305cf7da58SHong Zhang }
15315cf7da58SHong Zhang 
1532e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
15335cf7da58SHong Zhang {
15345cf7da58SHong Zhang   PetscErrorCode ierr;
15355cf7da58SHong Zhang 
15365cf7da58SHong Zhang   PetscFunctionBegin;
15375cf7da58SHong Zhang   if (Ju) {
15385cf7da58SHong Zhang     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
15395cf7da58SHong Zhang   } else {
15405cf7da58SHong Zhang     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
15415cf7da58SHong Zhang   }
15425cf7da58SHong Zhang   PetscFunctionReturn(0);
15435cf7da58SHong Zhang }
15445cf7da58SHong Zhang 
1545e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1546883e35e8SHong Zhang {
1547883e35e8SHong Zhang   PetscErrorCode ierr;
1548883e35e8SHong Zhang   PetscInt       j,*cols;
1549883e35e8SHong Zhang   PetscScalar    *zeros;
1550883e35e8SHong Zhang 
1551883e35e8SHong Zhang   PetscFunctionBegin;
1552883e35e8SHong Zhang   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
1553883e35e8SHong Zhang   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1554883e35e8SHong Zhang   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
1555883e35e8SHong Zhang   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
15561ad426b7SHong Zhang   PetscFunctionReturn(0);
15571ad426b7SHong Zhang }
1558a4e85ca8SHong Zhang 
1559e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
15603e97b6e8SHong Zhang {
15613e97b6e8SHong Zhang   PetscErrorCode ierr;
15623e97b6e8SHong Zhang   PetscInt       j,M,N,row,col,ncols_u;
15633e97b6e8SHong Zhang   const PetscInt *cols;
15643e97b6e8SHong Zhang   PetscScalar    zero=0.0;
15653e97b6e8SHong Zhang 
15663e97b6e8SHong Zhang   PetscFunctionBegin;
15673e97b6e8SHong Zhang   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
15683e97b6e8SHong 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);
15693e97b6e8SHong Zhang 
15703e97b6e8SHong Zhang   for (row=0; row<nrows; row++) {
15713e97b6e8SHong Zhang     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
15723e97b6e8SHong Zhang     for (j=0; j<ncols_u; j++) {
15733e97b6e8SHong Zhang       col = cols[j] + cstart;
15743e97b6e8SHong Zhang       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
15753e97b6e8SHong Zhang     }
15763e97b6e8SHong Zhang     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
15773e97b6e8SHong Zhang   }
15783e97b6e8SHong Zhang   PetscFunctionReturn(0);
15793e97b6e8SHong Zhang }
15801ad426b7SHong Zhang 
1581e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1582a4e85ca8SHong Zhang {
1583a4e85ca8SHong Zhang   PetscErrorCode ierr;
1584f4431b8cSHong Zhang 
1585a4e85ca8SHong Zhang   PetscFunctionBegin;
1586a4e85ca8SHong Zhang   if (Ju) {
1587a4e85ca8SHong Zhang     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1588a4e85ca8SHong Zhang   } else {
1589a4e85ca8SHong Zhang     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1590a4e85ca8SHong Zhang   }
1591a4e85ca8SHong Zhang   PetscFunctionReturn(0);
1592a4e85ca8SHong Zhang }
1593a4e85ca8SHong Zhang 
159424121865SAdrian 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.
159524121865SAdrian Maldonado */
159624121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
159724121865SAdrian Maldonado {
159824121865SAdrian Maldonado   PetscErrorCode ierr;
159924121865SAdrian Maldonado   PetscInt       i, size, dof;
160024121865SAdrian Maldonado   PetscInt       *glob2loc;
160124121865SAdrian Maldonado 
160224121865SAdrian Maldonado   PetscFunctionBegin;
160324121865SAdrian Maldonado   ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr);
160424121865SAdrian Maldonado   ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr);
160524121865SAdrian Maldonado 
160624121865SAdrian Maldonado   for (i = 0; i < size; i++) {
160724121865SAdrian Maldonado     ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr);
160824121865SAdrian Maldonado     dof = (dof >= 0) ? dof : -(dof + 1);
160924121865SAdrian Maldonado     glob2loc[i] = dof;
161024121865SAdrian Maldonado   }
161124121865SAdrian Maldonado 
161224121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
161324121865SAdrian Maldonado #if 0
161424121865SAdrian Maldonado   ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
161524121865SAdrian Maldonado #endif
161624121865SAdrian Maldonado   PetscFunctionReturn(0);
161724121865SAdrian Maldonado }
161824121865SAdrian Maldonado 
161901ad2aeeSHong Zhang #include <petsc/private/matimpl.h>
16201ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
16211ad426b7SHong Zhang {
16221ad426b7SHong Zhang   PetscErrorCode ierr;
162324121865SAdrian Maldonado   PetscMPIInt    rank, size;
16241ad426b7SHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
1625a4e85ca8SHong Zhang   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1626840c2264SHong Zhang   PetscInt       cstart,ncols,j,e,v;
162724121865SAdrian Maldonado   PetscBool      ghost,ghost_vc,ghost2,isNest;
1628a4e85ca8SHong Zhang   Mat            Juser;
1629bfbc38dcSHong Zhang   PetscSection   sectionGlobal;
1630447d78afSSatish Balay   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
1631a4e85ca8SHong Zhang   const PetscInt *edges,*cone;
16325cf7da58SHong Zhang   MPI_Comm       comm;
163324121865SAdrian Maldonado   MatType        mtype;
16345cf7da58SHong Zhang   Vec            vd_nz,vo_nz;
16355cf7da58SHong Zhang   PetscInt       *dnnz,*onnz;
16365cf7da58SHong Zhang   PetscScalar    *vdnz,*vonz;
16371ad426b7SHong Zhang 
16381ad426b7SHong Zhang   PetscFunctionBegin;
163924121865SAdrian Maldonado   mtype = dm->mattype;
164024121865SAdrian Maldonado   ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr);
164124121865SAdrian Maldonado 
164224121865SAdrian Maldonado   if (isNest) {
16430731d606SHong Zhang     /* ierr = DMCreateMatrix_Network_Nest(); */
164424121865SAdrian Maldonado     PetscInt   eDof, vDof;
164524121865SAdrian Maldonado     Mat        j11, j12, j21, j22, bA[2][2];
164624121865SAdrian Maldonado     ISLocalToGlobalMapping eISMap, vISMap;
164724121865SAdrian Maldonado 
164824121865SAdrian Maldonado     ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
164924121865SAdrian Maldonado     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
165024121865SAdrian Maldonado     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
165124121865SAdrian Maldonado 
165224121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr);
165324121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr);
165424121865SAdrian Maldonado 
165501ad2aeeSHong Zhang     ierr = MatCreate(comm, &j11);CHKERRQ(ierr);
165624121865SAdrian Maldonado     ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
165724121865SAdrian Maldonado     ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr);
165824121865SAdrian Maldonado 
165901ad2aeeSHong Zhang     ierr = MatCreate(comm, &j12);CHKERRQ(ierr);
166024121865SAdrian Maldonado     ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr);
166124121865SAdrian Maldonado     ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr);
166224121865SAdrian Maldonado 
166301ad2aeeSHong Zhang     ierr = MatCreate(comm, &j21);CHKERRQ(ierr);
166424121865SAdrian Maldonado     ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
166524121865SAdrian Maldonado     ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr);
166624121865SAdrian Maldonado 
166701ad2aeeSHong Zhang     ierr = MatCreate(comm, &j22);CHKERRQ(ierr);
166824121865SAdrian Maldonado     ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
166924121865SAdrian Maldonado     ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr);
167024121865SAdrian Maldonado 
16713f6a6bdaSHong Zhang     bA[0][0] = j11;
16723f6a6bdaSHong Zhang     bA[0][1] = j12;
16733f6a6bdaSHong Zhang     bA[1][0] = j21;
16743f6a6bdaSHong Zhang     bA[1][1] = j22;
167524121865SAdrian Maldonado 
167624121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr);
167724121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr);
167824121865SAdrian Maldonado 
167924121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr);
168024121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr);
168124121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr);
168224121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr);
168324121865SAdrian Maldonado 
168424121865SAdrian Maldonado     ierr = MatSetUp(j11);CHKERRQ(ierr);
168524121865SAdrian Maldonado     ierr = MatSetUp(j12);CHKERRQ(ierr);
168624121865SAdrian Maldonado     ierr = MatSetUp(j21);CHKERRQ(ierr);
168724121865SAdrian Maldonado     ierr = MatSetUp(j22);CHKERRQ(ierr);
168824121865SAdrian Maldonado 
168901ad2aeeSHong Zhang     ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr);
169024121865SAdrian Maldonado     ierr = MatSetUp(*J);CHKERRQ(ierr);
169124121865SAdrian Maldonado     ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr);
169224121865SAdrian Maldonado     ierr = MatDestroy(&j11);CHKERRQ(ierr);
169324121865SAdrian Maldonado     ierr = MatDestroy(&j12);CHKERRQ(ierr);
169424121865SAdrian Maldonado     ierr = MatDestroy(&j21);CHKERRQ(ierr);
169524121865SAdrian Maldonado     ierr = MatDestroy(&j22);CHKERRQ(ierr);
169624121865SAdrian Maldonado 
169724121865SAdrian Maldonado     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169824121865SAdrian Maldonado     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169924121865SAdrian Maldonado 
170024121865SAdrian Maldonado     /* Free structures */
170124121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr);
170224121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr);
170324121865SAdrian Maldonado 
170424121865SAdrian Maldonado     PetscFunctionReturn(0);
170524121865SAdrian Maldonado   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1706a4e85ca8SHong Zhang     /* user does not provide Jacobian blocks */
1707bfbc38dcSHong Zhang     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1708bfbc38dcSHong Zhang     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
17091ad426b7SHong Zhang     PetscFunctionReturn(0);
17101ad426b7SHong Zhang   }
17111ad426b7SHong Zhang 
1712bfbc38dcSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
1713e87a4003SBarry Smith   ierr = DMGetGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1714bfbc38dcSHong Zhang   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1715bfbc38dcSHong Zhang   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
17162a945128SHong Zhang 
17172a945128SHong Zhang   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
17182a945128SHong Zhang   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
171989898e50SHong Zhang 
172089898e50SHong Zhang   /* (1) Set matrix preallocation */
172189898e50SHong Zhang   /*------------------------------*/
1722840c2264SHong Zhang   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1723840c2264SHong Zhang   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1724840c2264SHong Zhang   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1725840c2264SHong Zhang   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1726840c2264SHong Zhang   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1727840c2264SHong Zhang   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1728840c2264SHong Zhang 
172989898e50SHong Zhang   /* Set preallocation for edges */
173089898e50SHong Zhang   /*-----------------------------*/
1731840c2264SHong Zhang   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1732840c2264SHong Zhang 
1733bdcb62a2SHong Zhang   ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr);
1734840c2264SHong Zhang   for (e=eStart; e<eEnd; e++) {
1735840c2264SHong Zhang     /* Get row indices */
1736840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1737840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1738840c2264SHong Zhang     if (nrows) {
1739840c2264SHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1740840c2264SHong Zhang 
17415cf7da58SHong Zhang       /* Set preallocation for conntected vertices */
1742d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1743840c2264SHong Zhang       for (v=0; v<2; v++) {
1744840c2264SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1745840c2264SHong Zhang 
17468675203cSHong Zhang         if (network->Je) {
1747840c2264SHong Zhang           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
17488675203cSHong Zhang         } else Juser = NULL;
1749840c2264SHong Zhang         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
17505cf7da58SHong Zhang         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1751840c2264SHong Zhang       }
1752840c2264SHong Zhang 
175389898e50SHong Zhang       /* Set preallocation for edge self */
1754840c2264SHong Zhang       cstart = rstart;
17558675203cSHong Zhang       if (network->Je) {
1756840c2264SHong Zhang         Juser = network->Je[3*e]; /* Jacobian(e,e) */
17578675203cSHong Zhang       } else Juser = NULL;
17585cf7da58SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1759840c2264SHong Zhang     }
1760840c2264SHong Zhang   }
1761840c2264SHong Zhang 
176289898e50SHong Zhang   /* Set preallocation for vertices */
176389898e50SHong Zhang   /*--------------------------------*/
1764840c2264SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
17658675203cSHong Zhang   if (vEnd - vStart) vptr = network->Jvptr;
1766840c2264SHong Zhang 
1767840c2264SHong Zhang   for (v=vStart; v<vEnd; v++) {
1768840c2264SHong Zhang     /* Get row indices */
1769840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1770840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1771840c2264SHong Zhang     if (!nrows) continue;
1772840c2264SHong Zhang 
1773bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1774bdcb62a2SHong Zhang     if (ghost) {
1775bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1776bdcb62a2SHong Zhang     } else {
1777bdcb62a2SHong Zhang       rows_v = rows;
1778bdcb62a2SHong Zhang     }
1779bdcb62a2SHong Zhang 
1780bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1781840c2264SHong Zhang 
1782840c2264SHong Zhang     /* Get supporting edges and connected vertices */
1783840c2264SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1784840c2264SHong Zhang 
1785840c2264SHong Zhang     for (e=0; e<nedges; e++) {
1786840c2264SHong Zhang       /* Supporting edges */
1787840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1788840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1789840c2264SHong Zhang 
17908675203cSHong Zhang       if (network->Jv) {
1791840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
17928675203cSHong Zhang       } else Juser = NULL;
1793bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1794840c2264SHong Zhang 
1795840c2264SHong Zhang       /* Connected vertices */
1796d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
1797840c2264SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
1798840c2264SHong Zhang       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1799840c2264SHong Zhang 
1800840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1801840c2264SHong Zhang 
18028675203cSHong Zhang       if (network->Jv) {
1803840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
18048675203cSHong Zhang       } else Juser = NULL;
1805e102a522SHong Zhang       if (ghost_vc||ghost) {
1806e102a522SHong Zhang         ghost2 = PETSC_TRUE;
1807e102a522SHong Zhang       } else {
1808e102a522SHong Zhang         ghost2 = PETSC_FALSE;
1809e102a522SHong Zhang       }
1810e102a522SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr);
1811840c2264SHong Zhang     }
1812840c2264SHong Zhang 
181389898e50SHong Zhang     /* Set preallocation for vertex self */
1814840c2264SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1815840c2264SHong Zhang     if (!ghost) {
1816840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
18178675203cSHong Zhang       if (network->Jv) {
1818840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
18198675203cSHong Zhang       } else Juser = NULL;
1820bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1821840c2264SHong Zhang     }
1822bdcb62a2SHong Zhang     if (ghost) {
1823bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1824bdcb62a2SHong Zhang     }
1825840c2264SHong Zhang   }
1826840c2264SHong Zhang 
1827840c2264SHong Zhang   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1828840c2264SHong Zhang   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
18295cf7da58SHong Zhang 
18305cf7da58SHong Zhang   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
18315cf7da58SHong Zhang 
18325cf7da58SHong Zhang   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1833840c2264SHong Zhang   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1834840c2264SHong Zhang 
1835840c2264SHong Zhang   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1836840c2264SHong Zhang   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1837840c2264SHong Zhang   for (j=0; j<localSize; j++) {
1838e102a522SHong Zhang     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1839e102a522SHong Zhang     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1840840c2264SHong Zhang   }
1841840c2264SHong Zhang   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1842840c2264SHong Zhang   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1843840c2264SHong Zhang   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1844840c2264SHong Zhang   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1845840c2264SHong Zhang 
18465cf7da58SHong Zhang   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
18475cf7da58SHong Zhang   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
18485cf7da58SHong Zhang   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
18495cf7da58SHong Zhang 
18505cf7da58SHong Zhang   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
18515cf7da58SHong Zhang 
185289898e50SHong Zhang   /* (2) Set matrix entries for edges */
185389898e50SHong Zhang   /*----------------------------------*/
18541ad426b7SHong Zhang   for (e=eStart; e<eEnd; e++) {
1855bfbc38dcSHong Zhang     /* Get row indices */
18561ad426b7SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
185717df6e9eSHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
18584b976069SHong Zhang     if (nrows) {
185917df6e9eSHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
18601ad426b7SHong Zhang 
1861bfbc38dcSHong Zhang       /* Set matrix entries for conntected vertices */
1862d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1863bfbc38dcSHong Zhang       for (v=0; v<2; v++) {
1864bfbc38dcSHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1865883e35e8SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
18663e97b6e8SHong Zhang 
18678675203cSHong Zhang         if (network->Je) {
1868a4e85ca8SHong Zhang           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
18698675203cSHong Zhang         } else Juser = NULL;
1870a4e85ca8SHong Zhang         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1871bfbc38dcSHong Zhang       }
187217df6e9eSHong Zhang 
1873bfbc38dcSHong Zhang       /* Set matrix entries for edge self */
18743e97b6e8SHong Zhang       cstart = rstart;
18758675203cSHong Zhang       if (network->Je) {
1876a4e85ca8SHong Zhang         Juser = network->Je[3*e]; /* Jacobian(e,e) */
18778675203cSHong Zhang       } else Juser = NULL;
1878a4e85ca8SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
18791ad426b7SHong Zhang     }
18804b976069SHong Zhang   }
18811ad426b7SHong Zhang 
1882bfbc38dcSHong Zhang   /* Set matrix entries for vertices */
188383b2e829SHong Zhang   /*---------------------------------*/
18841ad426b7SHong Zhang   for (v=vStart; v<vEnd; v++) {
1885bfbc38dcSHong Zhang     /* Get row indices */
1886596e729fSHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1887596e729fSHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
18884b976069SHong Zhang     if (!nrows) continue;
1889596e729fSHong Zhang 
1890bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1891bdcb62a2SHong Zhang     if (ghost) {
1892bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1893bdcb62a2SHong Zhang     } else {
1894bdcb62a2SHong Zhang       rows_v = rows;
1895bdcb62a2SHong Zhang     }
1896bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1897596e729fSHong Zhang 
1898bfbc38dcSHong Zhang     /* Get supporting edges and connected vertices */
1899596e729fSHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1900596e729fSHong Zhang 
1901596e729fSHong Zhang     for (e=0; e<nedges; e++) {
1902bfbc38dcSHong Zhang       /* Supporting edges */
1903596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1904596e729fSHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1905596e729fSHong Zhang 
19068675203cSHong Zhang       if (network->Jv) {
1907a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
19088675203cSHong Zhang       } else Juser = NULL;
1909bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1910596e729fSHong Zhang 
1911bfbc38dcSHong Zhang       /* Connected vertices */
1912d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
19132a945128SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
19142a945128SHong Zhang 
191544aca652SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
191644aca652SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1917a4e85ca8SHong Zhang 
19188675203cSHong Zhang       if (network->Jv) {
1919a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
19208675203cSHong Zhang       } else Juser = NULL;
1921bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1922596e729fSHong Zhang     }
1923596e729fSHong Zhang 
1924bfbc38dcSHong Zhang     /* Set matrix entries for vertex self */
19251ad426b7SHong Zhang     if (!ghost) {
1926596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
19278675203cSHong Zhang       if (network->Jv) {
1928a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
19298675203cSHong Zhang       } else Juser = NULL;
1930bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr);
1931bdcb62a2SHong Zhang     }
1932bdcb62a2SHong Zhang     if (ghost) {
1933bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1934bdcb62a2SHong Zhang     }
19351ad426b7SHong Zhang   }
1936a4e85ca8SHong Zhang   ierr = PetscFree(rows);CHKERRQ(ierr);
1937bdcb62a2SHong Zhang 
19381ad426b7SHong Zhang   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19391ad426b7SHong Zhang   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1940dd6f46cdSHong Zhang 
19415f2c45f1SShri Abhyankar   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
19425f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19435f2c45f1SShri Abhyankar }
19445f2c45f1SShri Abhyankar 
19455f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm)
19465f2c45f1SShri Abhyankar {
19475f2c45f1SShri Abhyankar   PetscErrorCode ierr;
19485f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
19492727e31bSShri Abhyankar   PetscInt       j;
19505f2c45f1SShri Abhyankar 
19515f2c45f1SShri Abhyankar   PetscFunctionBegin;
19528415c774SShri Abhyankar   if (--network->refct > 0) PetscFunctionReturn(0);
195383b2e829SHong Zhang   if (network->Je) {
195483b2e829SHong Zhang     ierr = PetscFree(network->Je);CHKERRQ(ierr);
195583b2e829SHong Zhang   }
195683b2e829SHong Zhang   if (network->Jv) {
1957883e35e8SHong Zhang     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
195883b2e829SHong Zhang     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
19591ad426b7SHong Zhang   }
196013c2a604SAdrian Maldonado 
196113c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr);
196213c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr);
196313c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr);
196413c2a604SAdrian Maldonado   if (network->vertex.sf) {
196513c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr);
196613c2a604SAdrian Maldonado   }
196713c2a604SAdrian Maldonado   /* edge */
196813c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr);
196913c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr);
197013c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr);
197113c2a604SAdrian Maldonado   if (network->edge.sf) {
197213c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr);
197313c2a604SAdrian Maldonado   }
19745f2c45f1SShri Abhyankar   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
19755f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
19765f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
197783b2e829SHong Zhang 
19782727e31bSShri Abhyankar   for(j=0; j<network->nsubnet; j++) {
19792727e31bSShri Abhyankar     ierr = PetscFree(network->subnet[j].edges);CHKERRQ(ierr);
19802727e31bSShri Abhyankar   }
1981caf410d2SHong Zhang   ierr = PetscFree(network->subnetvtx);CHKERRQ(ierr);
1982caf410d2SHong Zhang 
1983e2aaf10cSShri Abhyankar   ierr = PetscFree(network->subnet);CHKERRQ(ierr);
19845f2c45f1SShri Abhyankar   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
1985caf410d2SHong Zhang   ierr = PetscFree2(network->header,network->cvalue);CHKERRQ(ierr);
19865f2c45f1SShri Abhyankar   ierr = PetscFree(network);CHKERRQ(ierr);
19875f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19885f2c45f1SShri Abhyankar }
19895f2c45f1SShri Abhyankar 
19905f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm,PetscViewer viewer)
19915f2c45f1SShri Abhyankar {
19925f2c45f1SShri Abhyankar   PetscErrorCode ierr;
19935f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
1994caf410d2SHong Zhang   PetscBool      iascii;
1995caf410d2SHong Zhang   PetscMPIInt    rank;
1996caf410d2SHong Zhang   PetscInt       p,nsubnet;
19975f2c45f1SShri Abhyankar 
19985f2c45f1SShri Abhyankar   PetscFunctionBegin;
1999caf410d2SHong Zhang   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first");
2000caf410d2SHong Zhang   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
2001caf410d2SHong Zhang   PetscValidHeaderSpecific(dm,DM_CLASSID, 1);
2002caf410d2SHong Zhang   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
2003caf410d2SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
2004caf410d2SHong Zhang   if (iascii) {
2005caf410d2SHong Zhang     const PetscInt    *cone,*vtx,*edges;
2006caf410d2SHong Zhang     PetscInt          vfrom,vto,i,j,nv,ne;
2007caf410d2SHong Zhang 
2008caf410d2SHong Zhang     nsubnet = network->nsubnet - network->ncsubnet; /* num of subnetworks */
2009caf410d2SHong Zhang     ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
2010caf410d2SHong Zhang     ierr = PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] nsubnet: %D; nsubnetCouple: %D; nEdges: %D; nVertices: %D\n",rank,nsubnet,network->ncsubnet,network->nEdges,network->nVertices);CHKERRQ(ierr);
2011caf410d2SHong Zhang 
2012caf410d2SHong Zhang     for (i=0; i<nsubnet; i++) {
2013caf410d2SHong Zhang       ierr = DMNetworkGetSubnetworkInfo(dm,i,&nv,&ne,&vtx,&edges);CHKERRQ(ierr);
2014caf410d2SHong Zhang       if (ne) {
2015caf410d2SHong Zhang         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "     Subnet %D: nEdges %D, nVertices %D\n",i,ne,nv);CHKERRQ(ierr);
2016caf410d2SHong Zhang         for (j=0; j<ne; j++) {
2017caf410d2SHong Zhang           p = edges[j];
2018caf410d2SHong Zhang           ierr = DMNetworkGetConnectedVertices(dm,p,&cone);CHKERRQ(ierr);
2019caf410d2SHong Zhang           ierr = DMNetworkGetGlobalVertexIndex(dm,cone[0],&vfrom);CHKERRQ(ierr);
2020caf410d2SHong Zhang           ierr = DMNetworkGetGlobalVertexIndex(dm,cone[1],&vto);CHKERRQ(ierr);
2021caf410d2SHong Zhang           ierr = DMNetworkGetGlobalEdgeIndex(dm,edges[j],&p);CHKERRQ(ierr);
2022caf410d2SHong Zhang           ierr = PetscViewerASCIISynchronizedPrintf(viewer, "       edge %D: %D----> %D\n",p,vfrom,vto);CHKERRQ(ierr);
2023caf410d2SHong Zhang         }
2024caf410d2SHong Zhang       }
2025caf410d2SHong Zhang     }
2026caf410d2SHong Zhang     /* Coupling subnets */
2027caf410d2SHong Zhang     nsubnet = network->nsubnet;
2028caf410d2SHong Zhang     for (; i<nsubnet; i++) {
2029caf410d2SHong Zhang       ierr = DMNetworkGetSubnetworkInfo(dm,i,&nv,&ne,&vtx,&edges);CHKERRQ(ierr);
2030caf410d2SHong Zhang       if (ne) {
2031caf410d2SHong Zhang         ierr = PetscViewerASCIISynchronizedPrintf(viewer, "     Subnet %D (couple): nEdges %D, nVertices %D\n",i,ne,nv);CHKERRQ(ierr);
2032caf410d2SHong Zhang         for (j=0; j<ne; j++) {
2033caf410d2SHong Zhang           p = edges[j];
2034caf410d2SHong Zhang           ierr = DMNetworkGetConnectedVertices(dm,p,&cone);CHKERRQ(ierr);
2035caf410d2SHong Zhang           ierr = DMNetworkGetGlobalVertexIndex(dm,cone[0],&vfrom);CHKERRQ(ierr);
2036caf410d2SHong Zhang           ierr = DMNetworkGetGlobalVertexIndex(dm,cone[1],&vto);CHKERRQ(ierr);
2037caf410d2SHong Zhang           ierr = PetscViewerASCIISynchronizedPrintf(viewer, "       edge %D: %D----> %D\n",p,vfrom,vto);CHKERRQ(ierr);
2038caf410d2SHong Zhang         }
2039caf410d2SHong Zhang       }
2040caf410d2SHong Zhang     }
2041caf410d2SHong Zhang     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
2042caf410d2SHong Zhang     ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
2043caf410d2SHong Zhang   } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMNetwork writing", ((PetscObject)viewer)->type_name);
20445f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
20455f2c45f1SShri Abhyankar }
20465f2c45f1SShri Abhyankar 
20475f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
20485f2c45f1SShri Abhyankar {
20495f2c45f1SShri Abhyankar   PetscErrorCode ierr;
20505f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
20515f2c45f1SShri Abhyankar 
20525f2c45f1SShri Abhyankar   PetscFunctionBegin;
20535f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
20545f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
20555f2c45f1SShri Abhyankar }
20565f2c45f1SShri Abhyankar 
20575f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
20585f2c45f1SShri Abhyankar {
20595f2c45f1SShri Abhyankar   PetscErrorCode ierr;
20605f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
20615f2c45f1SShri Abhyankar 
20625f2c45f1SShri Abhyankar   PetscFunctionBegin;
20635f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
20645f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
20655f2c45f1SShri Abhyankar }
20665f2c45f1SShri Abhyankar 
20675f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
20685f2c45f1SShri Abhyankar {
20695f2c45f1SShri Abhyankar   PetscErrorCode ierr;
20705f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
20715f2c45f1SShri Abhyankar 
20725f2c45f1SShri Abhyankar   PetscFunctionBegin;
20735f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
20745f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
20755f2c45f1SShri Abhyankar }
20765f2c45f1SShri Abhyankar 
20775f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
20785f2c45f1SShri Abhyankar {
20795f2c45f1SShri Abhyankar   PetscErrorCode ierr;
20805f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
20815f2c45f1SShri Abhyankar 
20825f2c45f1SShri Abhyankar   PetscFunctionBegin;
20835f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
20845f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
20855f2c45f1SShri Abhyankar }
2086