xref: /petsc/src/dm/impls/network/network.c (revision a2b725a8db0d6bf6cc2a1c6df7dd8029aadfff6e)
1af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h>  /*I  "petscdmnetwork.h"  I*/
25f2c45f1SShri Abhyankar #include <petscdmplex.h>
35f2c45f1SShri Abhyankar #include <petscsf.h>
45f2c45f1SShri Abhyankar 
55f2c45f1SShri Abhyankar /*@
6556ed216SShri Abhyankar   DMNetworkGetPlex - Gets the Plex DM associated with this network DM
7556ed216SShri Abhyankar 
8556ed216SShri Abhyankar   Not collective
9556ed216SShri Abhyankar 
10556ed216SShri Abhyankar   Input Parameters:
11556ed216SShri Abhyankar + netdm - the dm object
12556ed216SShri Abhyankar - plexmdm - the plex dm object
13556ed216SShri Abhyankar 
14556ed216SShri Abhyankar   Level: Advanced
15556ed216SShri Abhyankar 
16556ed216SShri Abhyankar .seealso: DMNetworkCreate()
17556ed216SShri Abhyankar @*/
18556ed216SShri Abhyankar PetscErrorCode DMNetworkGetPlex(DM netdm, DM *plexdm)
19556ed216SShri Abhyankar {
20556ed216SShri Abhyankar   DM_Network *network = (DM_Network*) netdm->data;
21556ed216SShri Abhyankar 
22556ed216SShri Abhyankar   PetscFunctionBegin;
23556ed216SShri Abhyankar   *plexdm = network->plex;
24556ed216SShri Abhyankar   PetscFunctionReturn(0);
25556ed216SShri Abhyankar }
26556ed216SShri Abhyankar 
27556ed216SShri Abhyankar /*@
28e2aaf10cSShri Abhyankar   DMNetworkSetSizes - Sets the number of subnetworks, local and global vertices and edges for each subnetwork.
295f2c45f1SShri Abhyankar 
305f2c45f1SShri Abhyankar   Collective on DM
315f2c45f1SShri Abhyankar 
325f2c45f1SShri Abhyankar   Input Parameters:
335f2c45f1SShri Abhyankar + dm - the dm object
34e2aaf10cSShri Abhyankar . Nsubnet - number of subnetworks
35f025b11dSHong Zhang . NsubnetCouple - number of coupling subnetworks
36f025b11dSHong Zhang . nV - number of local vertices for each subnetwork and coupling subnetwork
37f025b11dSHong Zhang . nE - number of local edges for each subnetwork and coupling subnetwork
38f025b11dSHong Zhang . NV - number of global vertices (or PETSC_DETERMINE) for each subnetwork and coupling subnetwork
39f025b11dSHong Zhang - NE - number of global edges (or PETSC_DETERMINE) for each subnetwork and coupling subnetwork
405f2c45f1SShri Abhyankar 
415f2c45f1SShri Abhyankar    Notes
425f2c45f1SShri Abhyankar    If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang.
435f2c45f1SShri Abhyankar 
44f025b11dSHong Zhang    You cannot change the sizes once they have been set. nV, nE, NV and NE are arrays of length Nsubnet+NsubnetCouple.
455f2c45f1SShri Abhyankar 
461b266c99SBarry Smith    Level: intermediate
471b266c99SBarry Smith 
481b266c99SBarry Smith .seealso: DMNetworkCreate()
495f2c45f1SShri Abhyankar @*/
50f025b11dSHong Zhang PetscErrorCode DMNetworkSetSizes(DM dm,PetscInt Nsubnet,PetscInt NsubnetCouple,PetscInt nV[], PetscInt nE[], PetscInt NV[], PetscInt NE[])
515f2c45f1SShri Abhyankar {
525f2c45f1SShri Abhyankar   PetscErrorCode ierr;
535f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
54e2aaf10cSShri Abhyankar   PetscInt       a[2],b[2],i;
555f2c45f1SShri Abhyankar 
565f2c45f1SShri Abhyankar   PetscFunctionBegin;
575f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
58e2aaf10cSShri Abhyankar   if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet);
597765340cSHong Zhang   if (NsubnetCouple < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of coupling subnetworks %D cannot be less than 0",NsubnetCouple);
607765340cSHong Zhang 
617765340cSHong Zhang   if (Nsubnet > 0) PetscValidLogicalCollectiveInt(dm,Nsubnet,2);
627765340cSHong Zhang   if (network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network");
637765340cSHong Zhang 
647765340cSHong Zhang   network->nsubnet  = Nsubnet + NsubnetCouple;
657765340cSHong Zhang   Nsubnet          += NsubnetCouple;
667765340cSHong Zhang   network->ncsubnet = NsubnetCouple;
677765340cSHong Zhang   ierr = PetscCalloc1(Nsubnet,&network->subnet);CHKERRQ(ierr);
687765340cSHong Zhang   for(i=0; i < Nsubnet; i++) {
69f025b11dSHong Zhang     if (NV) {
707765340cSHong Zhang       if (NV[i] > 0) PetscValidLogicalCollectiveInt(dm,NV[i],6);
717765340cSHong Zhang       if (NV[i] > 0 && nV[i] > NV[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local vertex size %D cannot be larger than global vertex size %D",i,nV[i],NV[i]);
72f025b11dSHong Zhang     }
73f025b11dSHong Zhang     if (NE) {
74f025b11dSHong Zhang       if (NE[i] > 0) PetscValidLogicalCollectiveInt(dm,NE[i],7);
757765340cSHong Zhang       if (NE[i] > 0 && nE[i] > NE[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local edge size %D cannot be larger than global edge size %D",i,nE[i],NE[i]);
76f025b11dSHong Zhang     }
77f025b11dSHong Zhang 
787765340cSHong Zhang     a[0] = nV[i]; a[1] = nE[i];
797765340cSHong Zhang     ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
807765340cSHong Zhang     network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1];
817765340cSHong Zhang 
827765340cSHong Zhang     network->subnet[i].id = i;
837765340cSHong Zhang 
847765340cSHong Zhang     network->subnet[i].nvtx = nV[i];
857765340cSHong Zhang     network->subnet[i].vStart = network->nVertices;
867765340cSHong Zhang     network->subnet[i].vEnd   = network->subnet[i].vStart + network->subnet[i].Nvtx;
877765340cSHong Zhang     network->nVertices += network->subnet[i].nvtx;
887765340cSHong Zhang     network->NVertices += network->subnet[i].Nvtx;
897765340cSHong Zhang 
907765340cSHong Zhang     network->subnet[i].nedge = nE[i];
917765340cSHong Zhang     network->subnet[i].eStart = network->nEdges;
927765340cSHong Zhang     network->subnet[i].eEnd = network->subnet[i].eStart + network->subnet[i].Nedge;
937765340cSHong Zhang     network->nEdges += network->subnet[i].nedge;
947765340cSHong Zhang     network->NEdges += network->subnet[i].Nedge;
957765340cSHong Zhang   }
967765340cSHong Zhang   PetscFunctionReturn(0);
977765340cSHong Zhang }
987765340cSHong Zhang 
995f2c45f1SShri Abhyankar /*@
1005f2c45f1SShri Abhyankar   DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network
1015f2c45f1SShri Abhyankar 
1025f2c45f1SShri Abhyankar   Logically collective on DM
1035f2c45f1SShri Abhyankar 
1045f2c45f1SShri Abhyankar   Input Parameters:
105e3e68989SHong Zhang + dm - the dm object
106e3e68989SHong Zhang . edgelist - list of edges for each subnetwork
107e3e68989SHong Zhang - edgelistCouple - list of edges for each coupling subnetwork
1085f2c45f1SShri Abhyankar 
1095f2c45f1SShri Abhyankar   Notes:
1105f2c45f1SShri Abhyankar   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
1115f2c45f1SShri Abhyankar   not be destroyed before the call to DMNetworkLayoutSetUp
1125f2c45f1SShri Abhyankar 
1135f2c45f1SShri Abhyankar   Level: intermediate
1145f2c45f1SShri Abhyankar 
115e3e68989SHong Zhang   Example usage:
116e3e68989SHong Zhang   Consider the following 2 separate networks and a coupling network:
117e3e68989SHong Zhang 
118e3e68989SHong Zhang .vb
119e3e68989SHong Zhang  network 0: v0 -> v1 -> v2 -> v3
120e3e68989SHong Zhang  network 1: v1 -> v2 -> v0
121e3e68989SHong Zhang  coupling network: network 1: v2 -> network 0: v0
122e3e68989SHong Zhang .ve
123e3e68989SHong Zhang 
124e3e68989SHong Zhang  The resulting input
125e3e68989SHong Zhang    edgelist[0] = [0 1 | 1 2 | 2 3];
126e3e68989SHong Zhang    edgelist[1] = [1 2 | 2 0]
127e3e68989SHong Zhang    edgelistCouple[0] = [(network)1 (v)2 (network)0 (v)0].
128e3e68989SHong Zhang 
1295f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes
1305f2c45f1SShri Abhyankar @*/
1314e18019cSBarry Smith PetscErrorCode DMNetworkSetEdgeList(DM dm,PetscInt *edgelist[],PetscInt *edgelistCouple[])
1325f2c45f1SShri Abhyankar {
1335f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network*) dm->data;
134e3e68989SHong Zhang   PetscInt   i,nsubnet,ncsubnet=network->ncsubnet;
1355f2c45f1SShri Abhyankar 
1365f2c45f1SShri Abhyankar   PetscFunctionBegin;
137e3e68989SHong Zhang   nsubnet = network->nsubnet - ncsubnet;
138e3e68989SHong Zhang   for(i=0; i < nsubnet; i++) {
139e2aaf10cSShri Abhyankar     network->subnet[i].edgelist = edgelist[i];
140e2aaf10cSShri Abhyankar   }
141e3e68989SHong Zhang   if (edgelistCouple) {
142e3e68989SHong Zhang     PetscInt j;
143e3e68989SHong Zhang     j = 0;
144e3e68989SHong Zhang     nsubnet = network->nsubnet;
145e3e68989SHong Zhang     while (i < nsubnet) network->subnet[i++].edgelist = edgelistCouple[j++];
146e3e68989SHong Zhang   }
1475f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
1485f2c45f1SShri Abhyankar }
1495f2c45f1SShri Abhyankar 
1505f2c45f1SShri Abhyankar /*@
1515f2c45f1SShri Abhyankar   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
1525f2c45f1SShri Abhyankar 
1535f2c45f1SShri Abhyankar   Collective on DM
1545f2c45f1SShri Abhyankar 
1555f2c45f1SShri Abhyankar   Input Parameters
1565f2c45f1SShri Abhyankar . DM - the dmnetwork object
1575f2c45f1SShri Abhyankar 
1585f2c45f1SShri Abhyankar   Notes:
1595f2c45f1SShri Abhyankar   This routine should be called after the network sizes and edgelists have been provided. It creates
1605f2c45f1SShri Abhyankar   the bare layout of the network and sets up the network to begin insertion of components.
1615f2c45f1SShri Abhyankar 
1625f2c45f1SShri Abhyankar   All the components should be registered before calling this routine.
1635f2c45f1SShri Abhyankar 
1645f2c45f1SShri Abhyankar   Level: intermediate
1655f2c45f1SShri Abhyankar 
1665f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList
1675f2c45f1SShri Abhyankar @*/
1685f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm)
1695f2c45f1SShri Abhyankar {
1705f2c45f1SShri Abhyankar   PetscErrorCode ierr;
1715f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
1725f2c45f1SShri Abhyankar   PetscInt       dim = 1; /* One dimensional network */
173991cf414SHong Zhang   PetscInt       numCorners=2,spacedim=2;
1746500d4abSHong Zhang   double         *vertexcoords=NULL;
1757765340cSHong Zhang   PetscInt       i,j,ndata,ctr=0,nsubnet;
1769ad3091eSHong Zhang   PetscInt       k,netid,vid;
177acdb140fSBarry Smith   PetscInt       *edgelist_couple=NULL;
1786500d4abSHong Zhang 
1796500d4abSHong Zhang   PetscFunctionBegin;
1806500d4abSHong Zhang   if (network->nVertices) {
1816500d4abSHong Zhang     ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr);
1826500d4abSHong Zhang   }
1836500d4abSHong Zhang 
1846500d4abSHong Zhang   /* Create the edgelist for the network by concatenating edgelists of the subnetworks */
1857765340cSHong Zhang   nsubnet = network->nsubnet - network->ncsubnet;
1866500d4abSHong Zhang   ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr);
1877765340cSHong Zhang   for (i=0; i < nsubnet; i++) {
1886500d4abSHong Zhang     for (j = 0; j < network->subnet[i].nedge; j++) {
1896500d4abSHong Zhang       network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
1906500d4abSHong Zhang       network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
1916500d4abSHong Zhang       ctr++;
1926500d4abSHong Zhang     }
1936500d4abSHong Zhang   }
1947765340cSHong Zhang 
1957765340cSHong Zhang   i       = nsubnet; /* netid of coupling subnet */
1967765340cSHong Zhang   nsubnet = network->nsubnet;
1977765340cSHong Zhang   while (i < nsubnet) {
198991cf414SHong Zhang     edgelist_couple = network->subnet[i].edgelist;
1996500d4abSHong Zhang     k = 0;
2006500d4abSHong Zhang     for (j = 0; j < network->subnet[i].nedge; j++) {
2016500d4abSHong Zhang       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
2026500d4abSHong Zhang       network->edges[2*ctr] = network->subnet[netid].vStart + vid; k += 2;
2036500d4abSHong Zhang 
2046500d4abSHong Zhang       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
205991cf414SHong Zhang       network->edges[2*ctr+1] = network->subnet[netid].vStart + vid; k+=2;
2066500d4abSHong Zhang       ctr++;
2076500d4abSHong Zhang     }
2087765340cSHong Zhang     i++;
2097765340cSHong Zhang   }
2106500d4abSHong Zhang 
21134bd1bf5SHong Zhang #if 0
2126500d4abSHong Zhang   for(i=0; i < network->nEdges; i++) {
2136500d4abSHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr);
2146500d4abSHong Zhang     printf("\n");
2156500d4abSHong Zhang   }
2166500d4abSHong Zhang #endif
2176500d4abSHong Zhang 
218acdb140fSBarry Smith #if defined(PETSC_USE_64BIT_INDICES)
219acdb140fSBarry Smith   {
220acdb140fSBarry Smith     int      *edges;
221acdb140fSBarry Smith     PetscInt ii;
222acdb140fSBarry Smith     ierr = PetscMalloc1(network->nEdges*numCorners,&edges);CHKERRQ(ierr);
223acdb140fSBarry Smith     for (ii=0; ii<network->nEdges*numCorners; ii++) {
224acdb140fSBarry Smith       edges[ii] = (int) network->edges[ii];
225acdb140fSBarry Smith     }
226acdb140fSBarry Smith     ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
227acdb140fSBarry Smith     ierr = PetscFree(edges);
228acdb140fSBarry Smith   }
229acdb140fSBarry Smith #else
2306500d4abSHong Zhang     ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
231acdb140fSBarry Smith  #endif
232acdb140fSBarry Smith 
2336500d4abSHong Zhang   if (network->nVertices) {
2346500d4abSHong Zhang     ierr = PetscFree(vertexcoords);CHKERRQ(ierr);
2356500d4abSHong Zhang   }
2366500d4abSHong Zhang   ierr = PetscFree(network->edges);CHKERRQ(ierr);
2376500d4abSHong Zhang 
2386500d4abSHong Zhang   ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr);
2396500d4abSHong Zhang   ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr);
2406500d4abSHong Zhang   ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr);
2416500d4abSHong Zhang 
2426500d4abSHong Zhang   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr);
2436500d4abSHong Zhang   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr);
2446500d4abSHong Zhang   ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr);
2456500d4abSHong Zhang   ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr);
2466500d4abSHong Zhang 
2476500d4abSHong Zhang   /* Create vertices and edges array for the subnetworks */
2486500d4abSHong Zhang   for (j=0; j < network->nsubnet; j++) {
2496500d4abSHong Zhang     ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr);
2506500d4abSHong Zhang     ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr);
2516500d4abSHong Zhang     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
2526500d4abSHong Zhang        These get updated when the vertices and edges are added. */
2536500d4abSHong Zhang     network->subnet[j].nvtx = network->subnet[j].nedge = 0;
2546500d4abSHong Zhang   }
2556500d4abSHong Zhang 
2566500d4abSHong Zhang   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
2576500d4abSHong Zhang   ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr);
2586500d4abSHong Zhang   for (i=network->eStart; i < network->eEnd; i++) {
2596500d4abSHong Zhang     network->header[i].index = i;   /* Global edge number */
2606500d4abSHong Zhang     for (j=0; j < network->nsubnet; j++) {
2616500d4abSHong Zhang       if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) {
2626500d4abSHong Zhang 	network->header[i].subnetid = j; /* Subnetwork id */
2636500d4abSHong Zhang 	network->subnet[j].edges[network->subnet[j].nedge++] = i;
2646500d4abSHong Zhang 	break;
2656500d4abSHong Zhang       }
2666500d4abSHong Zhang     }
2676500d4abSHong Zhang     network->header[i].ndata = 0;
2686500d4abSHong Zhang     ndata = network->header[i].ndata;
2696500d4abSHong Zhang     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
2706500d4abSHong Zhang     network->header[i].offset[ndata] = 0;
2716500d4abSHong Zhang   }
2726500d4abSHong Zhang 
2736500d4abSHong Zhang   for(i=network->vStart; i < network->vEnd; i++) {
2746500d4abSHong Zhang     network->header[i].index = i - network->vStart;
2756500d4abSHong Zhang     for (j=0; j < network->nsubnet; j++) {
2766500d4abSHong Zhang       if ((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) {
2776500d4abSHong Zhang 	network->header[i].subnetid = j;
2786500d4abSHong Zhang 	network->subnet[j].vertices[network->subnet[j].nvtx++] = i;
2796500d4abSHong Zhang 	break;
2806500d4abSHong Zhang       }
2816500d4abSHong Zhang     }
2826500d4abSHong Zhang     network->header[i].ndata = 0;
2836500d4abSHong Zhang     ndata = network->header[i].ndata;
2846500d4abSHong Zhang     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
2856500d4abSHong Zhang     network->header[i].offset[ndata] = 0;
2866500d4abSHong Zhang   }
2876500d4abSHong Zhang 
2886500d4abSHong Zhang   ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr);
2895f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2905f2c45f1SShri Abhyankar }
2915f2c45f1SShri Abhyankar 
29294ef8ddeSSatish Balay /*@C
2932727e31bSShri Abhyankar   DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork
2942727e31bSShri Abhyankar 
2952727e31bSShri Abhyankar   Input Parameters
2962727e31bSShri Abhyankar + dm   - the number object
2972727e31bSShri Abhyankar - id   - the ID (integer) of the subnetwork
2982727e31bSShri Abhyankar 
2992727e31bSShri Abhyankar   Output Parameters
3002727e31bSShri Abhyankar + nv    - number of vertices (local)
3012727e31bSShri Abhyankar . ne    - number of edges (local)
3022727e31bSShri Abhyankar . vtx   - local vertices for this subnetwork
303*a2b725a8SWilliam Gropp - edge  - local edges for this subnetwork
3042727e31bSShri Abhyankar 
3052727e31bSShri Abhyankar   Notes:
3062727e31bSShri Abhyankar   Cannot call this routine before DMNetworkLayoutSetup()
3072727e31bSShri Abhyankar 
30806dd6b0eSSatish Balay   Level: intermediate
30906dd6b0eSSatish Balay 
3102727e31bSShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
3112727e31bSShri Abhyankar @*/
3122727e31bSShri Abhyankar PetscErrorCode DMNetworkGetSubnetworkInfo(DM netdm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge)
3132727e31bSShri Abhyankar {
3142727e31bSShri Abhyankar   DM_Network     *network = (DM_Network*) netdm->data;
3152727e31bSShri Abhyankar 
3162727e31bSShri Abhyankar   PetscFunctionBegin;
3172727e31bSShri Abhyankar   *nv = network->subnet[id].nvtx;
3182727e31bSShri Abhyankar   *ne = network->subnet[id].nedge;
3192727e31bSShri Abhyankar   *vtx = network->subnet[id].vertices;
3202727e31bSShri Abhyankar   *edge = network->subnet[id].edges;
3212727e31bSShri Abhyankar   PetscFunctionReturn(0);
3222727e31bSShri Abhyankar }
3232727e31bSShri Abhyankar 
3242727e31bSShri Abhyankar /*@C
3255f2c45f1SShri Abhyankar   DMNetworkRegisterComponent - Registers the network component
3265f2c45f1SShri Abhyankar 
3275f2c45f1SShri Abhyankar   Logically collective on DM
3285f2c45f1SShri Abhyankar 
3295f2c45f1SShri Abhyankar   Input Parameters
3305f2c45f1SShri Abhyankar + dm   - the network object
3315f2c45f1SShri Abhyankar . name - the component name
3325f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data
3335f2c45f1SShri Abhyankar 
3345f2c45f1SShri Abhyankar    Output Parameters
3355f2c45f1SShri Abhyankar .   key - an integer key that defines the component
3365f2c45f1SShri Abhyankar 
3375f2c45f1SShri Abhyankar    Notes
3385f2c45f1SShri Abhyankar    This routine should be called by all processors before calling DMNetworkLayoutSetup().
3395f2c45f1SShri Abhyankar 
3405f2c45f1SShri Abhyankar    Level: intermediate
3415f2c45f1SShri Abhyankar 
3425f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
3435f2c45f1SShri Abhyankar @*/
3445f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key)
3455f2c45f1SShri Abhyankar {
3465f2c45f1SShri Abhyankar   PetscErrorCode        ierr;
3475f2c45f1SShri Abhyankar   DM_Network            *network = (DM_Network*) dm->data;
3485f2c45f1SShri Abhyankar   DMNetworkComponent    *component=&network->component[network->ncomponent];
3495f2c45f1SShri Abhyankar   PetscBool             flg=PETSC_FALSE;
3505f2c45f1SShri Abhyankar   PetscInt              i;
3515f2c45f1SShri Abhyankar 
3525f2c45f1SShri Abhyankar   PetscFunctionBegin;
3535f2c45f1SShri Abhyankar   for (i=0; i < network->ncomponent; i++) {
3545f2c45f1SShri Abhyankar     ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr);
3555f2c45f1SShri Abhyankar     if (flg) {
3565f2c45f1SShri Abhyankar       *key = i;
3575f2c45f1SShri Abhyankar       PetscFunctionReturn(0);
3585f2c45f1SShri Abhyankar     }
3596d64e262SShri Abhyankar   }
3606d64e262SShri Abhyankar   if(network->ncomponent == MAX_COMPONENTS) {
3616d64e262SShri Abhyankar     SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components registered exceeds the max %D",MAX_COMPONENTS);
3625f2c45f1SShri Abhyankar   }
3635f2c45f1SShri Abhyankar 
3645f2c45f1SShri Abhyankar   ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr);
3655f2c45f1SShri Abhyankar   component->size = size/sizeof(DMNetworkComponentGenericDataType);
3665f2c45f1SShri Abhyankar   *key = network->ncomponent;
3675f2c45f1SShri Abhyankar   network->ncomponent++;
3685f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3695f2c45f1SShri Abhyankar }
3705f2c45f1SShri Abhyankar 
3715f2c45f1SShri Abhyankar /*@
3725f2c45f1SShri Abhyankar   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.
3735f2c45f1SShri Abhyankar 
3745f2c45f1SShri Abhyankar   Not Collective
3755f2c45f1SShri Abhyankar 
3765f2c45f1SShri Abhyankar   Input Parameters:
377*a2b725a8SWilliam Gropp . dm - The DMNetwork object
3785f2c45f1SShri Abhyankar 
3795f2c45f1SShri Abhyankar   Output Paramters:
3805f2c45f1SShri Abhyankar + vStart - The first vertex point
3815f2c45f1SShri Abhyankar - vEnd   - One beyond the last vertex point
3825f2c45f1SShri Abhyankar 
3835f2c45f1SShri Abhyankar   Level: intermediate
3845f2c45f1SShri Abhyankar 
3855f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange
3865f2c45f1SShri Abhyankar @*/
3875f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
3885f2c45f1SShri Abhyankar {
3895f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
3905f2c45f1SShri Abhyankar 
3915f2c45f1SShri Abhyankar   PetscFunctionBegin;
3925f2c45f1SShri Abhyankar   if (vStart) *vStart = network->vStart;
3935f2c45f1SShri Abhyankar   if (vEnd) *vEnd = network->vEnd;
3945f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3955f2c45f1SShri Abhyankar }
3965f2c45f1SShri Abhyankar 
3975f2c45f1SShri Abhyankar /*@
3985f2c45f1SShri Abhyankar   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.
3995f2c45f1SShri Abhyankar 
4005f2c45f1SShri Abhyankar   Not Collective
4015f2c45f1SShri Abhyankar 
4025f2c45f1SShri Abhyankar   Input Parameters:
403*a2b725a8SWilliam Gropp . dm - The DMNetwork object
4045f2c45f1SShri Abhyankar 
4055f2c45f1SShri Abhyankar   Output Paramters:
4065f2c45f1SShri Abhyankar + eStart - The first edge point
4075f2c45f1SShri Abhyankar - eEnd   - One beyond the last edge point
4085f2c45f1SShri Abhyankar 
4095f2c45f1SShri Abhyankar   Level: intermediate
4105f2c45f1SShri Abhyankar 
4115f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange
4125f2c45f1SShri Abhyankar @*/
4135f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
4145f2c45f1SShri Abhyankar {
4155f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4165f2c45f1SShri Abhyankar 
4175f2c45f1SShri Abhyankar   PetscFunctionBegin;
4185f2c45f1SShri Abhyankar   if (eStart) *eStart = network->eStart;
4195f2c45f1SShri Abhyankar   if (eEnd) *eEnd = network->eEnd;
4205f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4215f2c45f1SShri Abhyankar }
4225f2c45f1SShri Abhyankar 
4237b6afd5bSHong Zhang /*@
424e85e6aecSHong Zhang   DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge.
4257b6afd5bSHong Zhang 
4267b6afd5bSHong Zhang   Not Collective
4277b6afd5bSHong Zhang 
4287b6afd5bSHong Zhang   Input Parameters:
4297b6afd5bSHong Zhang + dm - DMNetwork object
430e85e6aecSHong Zhang - p  - edge point
4317b6afd5bSHong Zhang 
4327b6afd5bSHong Zhang   Output Paramters:
433e85e6aecSHong Zhang . index - user global numbering for the edge
4347b6afd5bSHong Zhang 
4357b6afd5bSHong Zhang   Level: intermediate
4367b6afd5bSHong Zhang 
437e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalVertexIndex
4387b6afd5bSHong Zhang @*/
439e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
4407b6afd5bSHong Zhang {
4417b6afd5bSHong Zhang   PetscErrorCode    ierr;
4427b6afd5bSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
4437b6afd5bSHong Zhang   PetscInt          offsetp;
4447b6afd5bSHong Zhang   DMNetworkComponentHeader header;
4457b6afd5bSHong Zhang 
4467b6afd5bSHong Zhang   PetscFunctionBegin;
4477b6afd5bSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
4487b6afd5bSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
449e85e6aecSHong Zhang   *index = header->index;
4507b6afd5bSHong Zhang   PetscFunctionReturn(0);
4517b6afd5bSHong Zhang }
4527b6afd5bSHong Zhang 
4535f2c45f1SShri Abhyankar /*@
454e85e6aecSHong Zhang   DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex.
455e85e6aecSHong Zhang 
456e85e6aecSHong Zhang   Not Collective
457e85e6aecSHong Zhang 
458e85e6aecSHong Zhang   Input Parameters:
459e85e6aecSHong Zhang + dm - DMNetwork object
460e85e6aecSHong Zhang - p  - vertex point
461e85e6aecSHong Zhang 
462e85e6aecSHong Zhang   Output Paramters:
463e85e6aecSHong Zhang . index - user global numbering for the vertex
464e85e6aecSHong Zhang 
465e85e6aecSHong Zhang   Level: intermediate
466e85e6aecSHong Zhang 
467e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalEdgeIndex
468e85e6aecSHong Zhang @*/
469e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
470e85e6aecSHong Zhang {
471e85e6aecSHong Zhang   PetscErrorCode    ierr;
472e85e6aecSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
473e85e6aecSHong Zhang   PetscInt          offsetp;
474e85e6aecSHong Zhang   DMNetworkComponentHeader header;
475e85e6aecSHong Zhang 
476e85e6aecSHong Zhang   PetscFunctionBegin;
477e85e6aecSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
478e85e6aecSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
479e85e6aecSHong Zhang   *index = header->index;
480e85e6aecSHong Zhang   PetscFunctionReturn(0);
481e85e6aecSHong Zhang }
482e85e6aecSHong Zhang 
483c3b11c7cSShri Abhyankar /*
484c3b11c7cSShri Abhyankar   DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the
485c3b11c7cSShri Abhyankar                                     component value from the component data array
486c3b11c7cSShri Abhyankar 
487c3b11c7cSShri Abhyankar   Not Collective
488c3b11c7cSShri Abhyankar 
489c3b11c7cSShri Abhyankar   Input Parameters:
490c3b11c7cSShri Abhyankar + dm      - The DMNetwork object
491c3b11c7cSShri Abhyankar . p       - vertex/edge point
492c3b11c7cSShri Abhyankar - compnum - component number
493c3b11c7cSShri Abhyankar 
494c3b11c7cSShri Abhyankar   Output Parameters:
495c3b11c7cSShri Abhyankar + compkey - the key obtained when registering the component
496c3b11c7cSShri Abhyankar - offset  - offset into the component data array associated with the vertex/edge point
497c3b11c7cSShri Abhyankar 
498c3b11c7cSShri Abhyankar   Notes:
499c3b11c7cSShri Abhyankar   Typical usage:
500c3b11c7cSShri Abhyankar 
501c3b11c7cSShri Abhyankar   DMNetworkGetComponentDataArray(dm, &arr);
502c3b11c7cSShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
503c3b11c7cSShri Abhyankar   Loop over vertices or edges
504c3b11c7cSShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
505c3b11c7cSShri Abhyankar     Loop over numcomps
506c3b11c7cSShri Abhyankar       DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset);
507c3b11c7cSShri Abhyankar       compdata = (UserCompDataType)(arr+offset);
508c3b11c7cSShri Abhyankar 
509c3b11c7cSShri Abhyankar   Level: intermediate
510c3b11c7cSShri Abhyankar 
511c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
512c3b11c7cSShri Abhyankar */
513c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
514c3b11c7cSShri Abhyankar {
515c3b11c7cSShri Abhyankar   PetscErrorCode           ierr;
516c3b11c7cSShri Abhyankar   PetscInt                 offsetp;
517c3b11c7cSShri Abhyankar   DMNetworkComponentHeader header;
518c3b11c7cSShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
519c3b11c7cSShri Abhyankar 
520c3b11c7cSShri Abhyankar   PetscFunctionBegin;
521c3b11c7cSShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
522c3b11c7cSShri Abhyankar   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
523c3b11c7cSShri Abhyankar   if (compkey) *compkey = header->key[compnum];
524c3b11c7cSShri Abhyankar   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
525c3b11c7cSShri Abhyankar   PetscFunctionReturn(0);
526c3b11c7cSShri Abhyankar }
527c3b11c7cSShri Abhyankar 
528c3b11c7cSShri Abhyankar /*@
529c3b11c7cSShri Abhyankar   DMNetworkGetComponent - Returns the network component and its key
530c3b11c7cSShri Abhyankar 
531c3b11c7cSShri Abhyankar   Not Collective
532c3b11c7cSShri Abhyankar 
533c3b11c7cSShri Abhyankar   Input Parameters
534c3b11c7cSShri Abhyankar + dm - DMNetwork object
535c3b11c7cSShri Abhyankar . p  - edge or vertex point
536c3b11c7cSShri Abhyankar - compnum - component number
537c3b11c7cSShri Abhyankar 
538c3b11c7cSShri Abhyankar   Output Parameters:
539c3b11c7cSShri Abhyankar + compkey - the key set for this computing during registration
540c3b11c7cSShri Abhyankar - component - the component data
541c3b11c7cSShri Abhyankar 
542c3b11c7cSShri Abhyankar   Notes:
543c3b11c7cSShri Abhyankar   Typical usage:
544c3b11c7cSShri Abhyankar 
545c3b11c7cSShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
546c3b11c7cSShri Abhyankar   Loop over vertices or edges
547c3b11c7cSShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
548c3b11c7cSShri Abhyankar     Loop over numcomps
549c3b11c7cSShri Abhyankar       DMNetworkGetComponent(dm,v,compnum,&key,&component);
550c3b11c7cSShri Abhyankar 
551c3b11c7cSShri Abhyankar   Level: intermediate
552c3b11c7cSShri Abhyankar 
553c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetVariableOffset
554c3b11c7cSShri Abhyankar @*/
555c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *key, void **component)
556c3b11c7cSShri Abhyankar {
557c3b11c7cSShri Abhyankar   PetscErrorCode ierr;
558c3b11c7cSShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
559e108cb99SStefano Zampini   PetscInt       offsetd = 0;
560c3b11c7cSShri Abhyankar 
561c3b11c7cSShri Abhyankar   PetscFunctionBegin;
562c3b11c7cSShri Abhyankar   ierr = DMNetworkGetComponentKeyOffset(dm,p,compnum,key,&offsetd);CHKERRQ(ierr);
563c3b11c7cSShri Abhyankar   *component = network->componentdataarray+offsetd;
564c3b11c7cSShri Abhyankar   PetscFunctionReturn(0);
565c3b11c7cSShri Abhyankar }
566c3b11c7cSShri Abhyankar 
567e85e6aecSHong Zhang /*@
568325661f6SSatish Balay   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)
5695f2c45f1SShri Abhyankar 
5705f2c45f1SShri Abhyankar   Not Collective
5715f2c45f1SShri Abhyankar 
5725f2c45f1SShri Abhyankar   Input Parameters:
5735f2c45f1SShri Abhyankar + dm           - The DMNetwork object
5745f2c45f1SShri Abhyankar . p            - vertex/edge point
5755f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component
5765f2c45f1SShri Abhyankar - compvalue    - pointer to the data structure for the component
5775f2c45f1SShri Abhyankar 
5785f2c45f1SShri Abhyankar   Level: intermediate
5795f2c45f1SShri Abhyankar 
5805f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
5815f2c45f1SShri Abhyankar @*/
5825f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
5835f2c45f1SShri Abhyankar {
5845f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
58543a39a44SBarry Smith   DMNetworkComponent       *component = &network->component[componentkey];
5865f2c45f1SShri Abhyankar   DMNetworkComponentHeader header = &network->header[p];
5875f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue = &network->cvalue[p];
5885f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
5895f2c45f1SShri Abhyankar 
5905f2c45f1SShri Abhyankar   PetscFunctionBegin;
591fa58f0a9SHong 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);
592fa58f0a9SHong Zhang 
59343a39a44SBarry Smith   header->size[header->ndata] = component->size;
59443a39a44SBarry Smith   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
5955f2c45f1SShri Abhyankar   header->key[header->ndata] = componentkey;
5965f2c45f1SShri Abhyankar   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
5975f2c45f1SShri Abhyankar 
5985f2c45f1SShri Abhyankar   cvalue->data[header->ndata] = (void*)compvalue;
5995f2c45f1SShri Abhyankar   header->ndata++;
6005f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6015f2c45f1SShri Abhyankar }
6025f2c45f1SShri Abhyankar 
6035f2c45f1SShri Abhyankar /*@
6045f2c45f1SShri Abhyankar   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
6055f2c45f1SShri Abhyankar 
6065f2c45f1SShri Abhyankar   Not Collective
6075f2c45f1SShri Abhyankar 
6085f2c45f1SShri Abhyankar   Input Parameters:
6095f2c45f1SShri Abhyankar + dm - The DMNetwork object
610*a2b725a8SWilliam Gropp - p  - vertex/edge point
6115f2c45f1SShri Abhyankar 
6125f2c45f1SShri Abhyankar   Output Parameters:
6135f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge
6145f2c45f1SShri Abhyankar 
6155f2c45f1SShri Abhyankar   Level: intermediate
6165f2c45f1SShri Abhyankar 
6175f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
6185f2c45f1SShri Abhyankar @*/
6195f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
6205f2c45f1SShri Abhyankar {
6215f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6225f2c45f1SShri Abhyankar   PetscInt       offset;
6235f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6245f2c45f1SShri Abhyankar 
6255f2c45f1SShri Abhyankar   PetscFunctionBegin;
6265f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
6275f2c45f1SShri Abhyankar   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
6285f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6295f2c45f1SShri Abhyankar }
6305f2c45f1SShri Abhyankar 
6315f2c45f1SShri Abhyankar /*@
6325f2c45f1SShri Abhyankar   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
6335f2c45f1SShri Abhyankar 
6345f2c45f1SShri Abhyankar   Not Collective
6355f2c45f1SShri Abhyankar 
6365f2c45f1SShri Abhyankar   Input Parameters:
6375f2c45f1SShri Abhyankar + dm     - The DMNetwork object
6385f2c45f1SShri Abhyankar - p      - the edge/vertex point
6395f2c45f1SShri Abhyankar 
6405f2c45f1SShri Abhyankar   Output Parameters:
6415f2c45f1SShri Abhyankar . offset - the offset
6425f2c45f1SShri Abhyankar 
6435f2c45f1SShri Abhyankar   Level: intermediate
6445f2c45f1SShri Abhyankar 
6455f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
6465f2c45f1SShri Abhyankar @*/
6475f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
6485f2c45f1SShri Abhyankar {
6495f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6505f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6515f2c45f1SShri Abhyankar 
6525f2c45f1SShri Abhyankar   PetscFunctionBegin;
6535f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr);
6545f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6555f2c45f1SShri Abhyankar }
6565f2c45f1SShri Abhyankar 
6575f2c45f1SShri Abhyankar /*@
6585f2c45f1SShri Abhyankar   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
6595f2c45f1SShri Abhyankar 
6605f2c45f1SShri Abhyankar   Not Collective
6615f2c45f1SShri Abhyankar 
6625f2c45f1SShri Abhyankar   Input Parameters:
6635f2c45f1SShri Abhyankar + dm      - The DMNetwork object
6645f2c45f1SShri Abhyankar - p       - the edge/vertex point
6655f2c45f1SShri Abhyankar 
6665f2c45f1SShri Abhyankar   Output Parameters:
6675f2c45f1SShri Abhyankar . offsetg - the offset
6685f2c45f1SShri Abhyankar 
6695f2c45f1SShri Abhyankar   Level: intermediate
6705f2c45f1SShri Abhyankar 
6715f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
6725f2c45f1SShri Abhyankar @*/
6735f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
6745f2c45f1SShri Abhyankar {
6755f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6765f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6775f2c45f1SShri Abhyankar 
6785f2c45f1SShri Abhyankar   PetscFunctionBegin;
6795f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr);
6806fefedf4SHong Zhang   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */
6815f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6825f2c45f1SShri Abhyankar }
6835f2c45f1SShri Abhyankar 
68424121865SAdrian Maldonado /*@
68524121865SAdrian Maldonado   DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector.
68624121865SAdrian Maldonado 
68724121865SAdrian Maldonado   Not Collective
68824121865SAdrian Maldonado 
68924121865SAdrian Maldonado   Input Parameters:
69024121865SAdrian Maldonado + dm     - The DMNetwork object
69124121865SAdrian Maldonado - p      - the edge point
69224121865SAdrian Maldonado 
69324121865SAdrian Maldonado   Output Parameters:
69424121865SAdrian Maldonado . offset - the offset
69524121865SAdrian Maldonado 
69624121865SAdrian Maldonado   Level: intermediate
69724121865SAdrian Maldonado 
69824121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
69924121865SAdrian Maldonado @*/
70024121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
70124121865SAdrian Maldonado {
70224121865SAdrian Maldonado   PetscErrorCode ierr;
70324121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
70424121865SAdrian Maldonado 
70524121865SAdrian Maldonado   PetscFunctionBegin;
70624121865SAdrian Maldonado 
70724121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr);
70824121865SAdrian Maldonado   PetscFunctionReturn(0);
70924121865SAdrian Maldonado }
71024121865SAdrian Maldonado 
71124121865SAdrian Maldonado /*@
71224121865SAdrian Maldonado   DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector.
71324121865SAdrian Maldonado 
71424121865SAdrian Maldonado   Not Collective
71524121865SAdrian Maldonado 
71624121865SAdrian Maldonado   Input Parameters:
71724121865SAdrian Maldonado + dm     - The DMNetwork object
71824121865SAdrian Maldonado - p      - the vertex point
71924121865SAdrian Maldonado 
72024121865SAdrian Maldonado   Output Parameters:
72124121865SAdrian Maldonado . offset - the offset
72224121865SAdrian Maldonado 
72324121865SAdrian Maldonado   Level: intermediate
72424121865SAdrian Maldonado 
72524121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
72624121865SAdrian Maldonado @*/
72724121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
72824121865SAdrian Maldonado {
72924121865SAdrian Maldonado   PetscErrorCode ierr;
73024121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
73124121865SAdrian Maldonado 
73224121865SAdrian Maldonado   PetscFunctionBegin;
73324121865SAdrian Maldonado 
73424121865SAdrian Maldonado   p -= network->vStart;
73524121865SAdrian Maldonado 
73624121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr);
73724121865SAdrian Maldonado   PetscFunctionReturn(0);
73824121865SAdrian Maldonado }
7395f2c45f1SShri Abhyankar /*@
7405f2c45f1SShri Abhyankar   DMNetworkAddNumVariables - Add number of variables associated with a given point.
7415f2c45f1SShri Abhyankar 
7425f2c45f1SShri Abhyankar   Not Collective
7435f2c45f1SShri Abhyankar 
7445f2c45f1SShri Abhyankar   Input Parameters:
7455f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
7465f2c45f1SShri Abhyankar . p    - the vertex/edge point
7475f2c45f1SShri Abhyankar - nvar - number of additional variables
7485f2c45f1SShri Abhyankar 
7495f2c45f1SShri Abhyankar   Level: intermediate
7505f2c45f1SShri Abhyankar 
7515f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables
7525f2c45f1SShri Abhyankar @*/
7535f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
7545f2c45f1SShri Abhyankar {
7555f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7565f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7575f2c45f1SShri Abhyankar 
7585f2c45f1SShri Abhyankar   PetscFunctionBegin;
7595f2c45f1SShri Abhyankar   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
7605f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7615f2c45f1SShri Abhyankar }
7625f2c45f1SShri Abhyankar 
76327f51fceSHong Zhang /*@
76427f51fceSHong Zhang   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
76527f51fceSHong Zhang 
76627f51fceSHong Zhang   Not Collective
76727f51fceSHong Zhang 
76827f51fceSHong Zhang   Input Parameters:
76927f51fceSHong Zhang + dm   - The DMNetworkObject
77027f51fceSHong Zhang - p    - the vertex/edge point
77127f51fceSHong Zhang 
77227f51fceSHong Zhang   Output Parameters:
77327f51fceSHong Zhang . nvar - number of variables
77427f51fceSHong Zhang 
77527f51fceSHong Zhang   Level: intermediate
77627f51fceSHong Zhang 
77727f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
77827f51fceSHong Zhang @*/
77927f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
78027f51fceSHong Zhang {
78127f51fceSHong Zhang   PetscErrorCode ierr;
78227f51fceSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
78327f51fceSHong Zhang 
78427f51fceSHong Zhang   PetscFunctionBegin;
78527f51fceSHong Zhang   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
78627f51fceSHong Zhang   PetscFunctionReturn(0);
78727f51fceSHong Zhang }
78827f51fceSHong Zhang 
7895f2c45f1SShri Abhyankar /*@
7905f2c45f1SShri Abhyankar   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
7915f2c45f1SShri Abhyankar 
7925f2c45f1SShri Abhyankar   Not Collective
7935f2c45f1SShri Abhyankar 
7945f2c45f1SShri Abhyankar   Input Parameters:
7955f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
7965f2c45f1SShri Abhyankar . p    - the vertex/edge point
7975f2c45f1SShri Abhyankar - nvar - number of variables
7985f2c45f1SShri Abhyankar 
7995f2c45f1SShri Abhyankar   Level: intermediate
8005f2c45f1SShri Abhyankar 
8015f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables
8025f2c45f1SShri Abhyankar @*/
8035f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
8045f2c45f1SShri Abhyankar {
8055f2c45f1SShri Abhyankar   PetscErrorCode ierr;
8065f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
8075f2c45f1SShri Abhyankar 
8085f2c45f1SShri Abhyankar   PetscFunctionBegin;
8095f2c45f1SShri Abhyankar   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
8105f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
8115f2c45f1SShri Abhyankar }
8125f2c45f1SShri Abhyankar 
8135f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This
8145f2c45f1SShri Abhyankar    function is called during DMSetUp() */
8155f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm)
8165f2c45f1SShri Abhyankar {
8175f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
8185f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
819e53b5ba3SHong Zhang   PetscInt                 arr_size,p,offset,offsetp,ncomp,i;
8205f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
8215f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue;
8225f2c45f1SShri Abhyankar   DMNetworkComponentGenericDataType *componentdataarray;
8235f2c45f1SShri Abhyankar 
8245f2c45f1SShri Abhyankar   PetscFunctionBegin;
8255f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
8265f2c45f1SShri Abhyankar   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
82775b160a0SShri Abhyankar   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
8285f2c45f1SShri Abhyankar   componentdataarray = network->componentdataarray;
8295f2c45f1SShri Abhyankar   for (p = network->pStart; p < network->pEnd; p++) {
8305f2c45f1SShri Abhyankar     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
8315f2c45f1SShri Abhyankar     /* Copy header */
8325f2c45f1SShri Abhyankar     header = &network->header[p];
833302440fdSBarry Smith     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
8345f2c45f1SShri Abhyankar     /* Copy data */
8355f2c45f1SShri Abhyankar     cvalue = &network->cvalue[p];
8365f2c45f1SShri Abhyankar     ncomp = header->ndata;
8375f2c45f1SShri Abhyankar     for (i = 0; i < ncomp; i++) {
8385f2c45f1SShri Abhyankar       offset = offsetp + network->dataheadersize + header->offset[i];
839302440fdSBarry Smith       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
8405f2c45f1SShri Abhyankar     }
8415f2c45f1SShri Abhyankar   }
8425f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
8435f2c45f1SShri Abhyankar }
8445f2c45f1SShri Abhyankar 
8455f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */
8465f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm)
8475f2c45f1SShri Abhyankar {
8485f2c45f1SShri Abhyankar   PetscErrorCode ierr;
8495f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
8505f2c45f1SShri Abhyankar 
8515f2c45f1SShri Abhyankar   PetscFunctionBegin;
8525f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
8535f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
8545f2c45f1SShri Abhyankar }
8555f2c45f1SShri Abhyankar 
8565f2c45f1SShri Abhyankar /*@C
8575f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray - Returns the component data array
8585f2c45f1SShri Abhyankar 
8595f2c45f1SShri Abhyankar   Not Collective
8605f2c45f1SShri Abhyankar 
8615f2c45f1SShri Abhyankar   Input Parameters:
8625f2c45f1SShri Abhyankar . dm - The DMNetwork Object
8635f2c45f1SShri Abhyankar 
8645f2c45f1SShri Abhyankar   Output Parameters:
8655f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components
8665f2c45f1SShri Abhyankar 
8675f2c45f1SShri Abhyankar   Level: intermediate
8685f2c45f1SShri Abhyankar 
869a730d845SHong Zhang .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents
8705f2c45f1SShri Abhyankar @*/
8715f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
8725f2c45f1SShri Abhyankar {
8735f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
8745f2c45f1SShri Abhyankar 
8755f2c45f1SShri Abhyankar   PetscFunctionBegin;
8765f2c45f1SShri Abhyankar   *componentdataarray = network->componentdataarray;
8775f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
8785f2c45f1SShri Abhyankar }
8795f2c45f1SShri Abhyankar 
88024121865SAdrian Maldonado /* Get a subsection from a range of points */
88124121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection)
88224121865SAdrian Maldonado {
88324121865SAdrian Maldonado   PetscErrorCode ierr;
88424121865SAdrian Maldonado   PetscInt       i, nvar;
88524121865SAdrian Maldonado 
88624121865SAdrian Maldonado   PetscFunctionBegin;
88724121865SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr);
88824121865SAdrian Maldonado   ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr);
88924121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
89024121865SAdrian Maldonado     ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr);
89124121865SAdrian Maldonado     ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr);
89224121865SAdrian Maldonado   }
89324121865SAdrian Maldonado 
89424121865SAdrian Maldonado   ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr);
89524121865SAdrian Maldonado   PetscFunctionReturn(0);
89624121865SAdrian Maldonado }
89724121865SAdrian Maldonado 
89824121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */
89924121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
90024121865SAdrian Maldonado {
90124121865SAdrian Maldonado   PetscErrorCode ierr;
90224121865SAdrian Maldonado   PetscInt       i, *subpoints;
90324121865SAdrian Maldonado 
90424121865SAdrian Maldonado   PetscFunctionBegin;
90524121865SAdrian Maldonado   /* Create index sets to map from "points" to "subpoints" */
90624121865SAdrian Maldonado   ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr);
90724121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
90824121865SAdrian Maldonado     subpoints[i - pstart] = i;
90924121865SAdrian Maldonado   }
910459726d8SSatish Balay   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr);
91124121865SAdrian Maldonado   ierr = PetscFree(subpoints);CHKERRQ(ierr);
91224121865SAdrian Maldonado   PetscFunctionReturn(0);
91324121865SAdrian Maldonado }
91424121865SAdrian Maldonado 
91524121865SAdrian Maldonado /*@
91624121865SAdrian Maldonado   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute.
91724121865SAdrian Maldonado 
91824121865SAdrian Maldonado   Collective
91924121865SAdrian Maldonado 
92024121865SAdrian Maldonado   Input Parameters:
92124121865SAdrian Maldonado . dm   - The DMNetworkObject
92224121865SAdrian Maldonado 
92324121865SAdrian Maldonado   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
92424121865SAdrian Maldonado 
92524121865SAdrian Maldonado   points = [0 1 2 3 4 5 6]
92624121865SAdrian Maldonado 
92724121865SAdrian Maldonado   where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]).
92824121865SAdrian Maldonado 
92924121865SAdrian Maldonado   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
93024121865SAdrian Maldonado 
93124121865SAdrian Maldonado   Level: intermediate
93224121865SAdrian Maldonado 
93324121865SAdrian Maldonado @*/
93424121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
93524121865SAdrian Maldonado {
93624121865SAdrian Maldonado   PetscErrorCode ierr;
93724121865SAdrian Maldonado   MPI_Comm       comm;
9389852e123SBarry Smith   PetscMPIInt    rank, size;
93924121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
94024121865SAdrian Maldonado 
941eab1376dSHong Zhang   PetscFunctionBegin;
94224121865SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
94324121865SAdrian Maldonado   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
9449852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
94524121865SAdrian Maldonado 
94624121865SAdrian Maldonado   /* Create maps for vertices and edges */
94724121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr);
94824121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr);
94924121865SAdrian Maldonado 
95024121865SAdrian Maldonado   /* Create local sub-sections */
95124121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr);
95224121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr);
95324121865SAdrian Maldonado 
9549852e123SBarry Smith   if (size > 1) {
95524121865SAdrian Maldonado     ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr);
95624121865SAdrian Maldonado     ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr);
95724121865SAdrian Maldonado   ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr);
95824121865SAdrian Maldonado   ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr);
95924121865SAdrian Maldonado   } else {
96024121865SAdrian Maldonado   /* create structures for vertex */
96124121865SAdrian Maldonado   ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr);
96224121865SAdrian Maldonado   /* create structures for edge */
96324121865SAdrian Maldonado   ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr);
96424121865SAdrian Maldonado   }
96524121865SAdrian Maldonado 
96624121865SAdrian Maldonado 
96724121865SAdrian Maldonado   /* Add viewers */
96824121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr);
96924121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr);
97024121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr);
97124121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr);
97224121865SAdrian Maldonado 
97324121865SAdrian Maldonado   PetscFunctionReturn(0);
97424121865SAdrian Maldonado }
9757b6afd5bSHong Zhang 
9765f2c45f1SShri Abhyankar /*@
9775f2c45f1SShri Abhyankar   DMNetworkDistribute - Distributes the network and moves associated component data.
9785f2c45f1SShri Abhyankar 
9795f2c45f1SShri Abhyankar   Collective
9805f2c45f1SShri Abhyankar 
9815f2c45f1SShri Abhyankar   Input Parameter:
982d3464fd4SAdrian Maldonado + DM - the DMNetwork object
9835f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default
9845f2c45f1SShri Abhyankar 
9855f2c45f1SShri Abhyankar   Notes:
9868b171c8eSHong Zhang   Distributes the network with <overlap>-overlapping partitioning of the edges.
9875f2c45f1SShri Abhyankar 
9885f2c45f1SShri Abhyankar   Level: intermediate
9895f2c45f1SShri Abhyankar 
9905f2c45f1SShri Abhyankar .seealso: DMNetworkCreate
9915f2c45f1SShri Abhyankar @*/
992d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
9935f2c45f1SShri Abhyankar {
994d3464fd4SAdrian Maldonado   MPI_Comm       comm;
9955f2c45f1SShri Abhyankar   PetscErrorCode ierr;
996d3464fd4SAdrian Maldonado   PetscMPIInt    size;
997d3464fd4SAdrian Maldonado   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
998d3464fd4SAdrian Maldonado   DM_Network     *newDMnetwork;
9995f2c45f1SShri Abhyankar   PetscSF        pointsf;
10005f2c45f1SShri Abhyankar   DM             newDM;
100151ac5effSHong Zhang   PetscPartitioner part;
1002b9c6e19dSShri Abhyankar   PetscInt         j,e,v,offset;
1003b9c6e19dSShri Abhyankar   DMNetworkComponentHeader header;
10045f2c45f1SShri Abhyankar 
10055f2c45f1SShri Abhyankar   PetscFunctionBegin;
1006d3464fd4SAdrian Maldonado 
1007d3464fd4SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr);
1008d3464fd4SAdrian Maldonado   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1009d3464fd4SAdrian Maldonado   if (size == 1) PetscFunctionReturn(0);
1010d3464fd4SAdrian Maldonado 
1011d3464fd4SAdrian Maldonado   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr);
10125f2c45f1SShri Abhyankar   newDMnetwork = (DM_Network*)newDM->data;
10135f2c45f1SShri Abhyankar   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
101451ac5effSHong Zhang 
101551ac5effSHong Zhang   /* Enable runtime options for petscpartitioner */
101651ac5effSHong Zhang   ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr);
101751ac5effSHong Zhang   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
101851ac5effSHong Zhang 
10195f2c45f1SShri Abhyankar   /* Distribute plex dm and dof section */
102080cf41d5SMatthew G. Knepley   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
102151ac5effSHong Zhang 
10225f2c45f1SShri Abhyankar   /* Distribute dof section */
1023d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr);
10245f2c45f1SShri Abhyankar   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
1025d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr);
102651ac5effSHong Zhang 
10275f2c45f1SShri Abhyankar   /* Distribute data and associated section */
102831da1fc8SHong Zhang   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
102924121865SAdrian Maldonado 
10305f2c45f1SShri Abhyankar   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
10315f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
10325f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
10335f2c45f1SShri Abhyankar   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
10346fefedf4SHong Zhang   newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
10356fefedf4SHong Zhang   newDMnetwork->NVertices = oldDMnetwork->NVertices;
10365f2c45f1SShri Abhyankar   newDMnetwork->NEdges = oldDMnetwork->NEdges;
103724121865SAdrian Maldonado 
10385f2c45f1SShri Abhyankar   /* Set Dof section as the default section for dm */
1039e87a4003SBarry Smith   ierr = DMSetSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
1040e87a4003SBarry Smith   ierr = DMGetGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
10415f2c45f1SShri Abhyankar 
1042b9c6e19dSShri Abhyankar   /* Set up subnetwork info in the newDM */
1043b9c6e19dSShri Abhyankar   newDMnetwork->nsubnet = oldDMnetwork->nsubnet;
1044b9c6e19dSShri Abhyankar   ierr = PetscCalloc1(newDMnetwork->nsubnet,&newDMnetwork->subnet);CHKERRQ(ierr);
1045b9c6e19dSShri Abhyankar   /* Copy over the global number of vertices and edges in each subnetwork. Note that these are already
1046b9c6e19dSShri Abhyankar      calculated in DMNetworkLayoutSetUp()
1047b9c6e19dSShri Abhyankar   */
1048b9c6e19dSShri Abhyankar   for(j=0; j < newDMnetwork->nsubnet; j++) {
1049b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx;
1050b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge;
1051b9c6e19dSShri Abhyankar   }
1052b9c6e19dSShri Abhyankar 
1053b9c6e19dSShri Abhyankar   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1054b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
1055b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1056b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].nedge++;
1057b9c6e19dSShri Abhyankar   }
1058b9c6e19dSShri Abhyankar 
1059b9c6e19dSShri Abhyankar   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1060b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
1061b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1062b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].nvtx++;
1063b9c6e19dSShri Abhyankar   }
1064b9c6e19dSShri Abhyankar 
1065b9c6e19dSShri Abhyankar   /* Now create the vertices and edge arrays for the subnetworks */
1066b9c6e19dSShri Abhyankar   for(j=0; j < newDMnetwork->nsubnet; j++) {
1067b9c6e19dSShri Abhyankar     ierr = PetscCalloc1(newDMnetwork->subnet[j].nedge,&newDMnetwork->subnet[j].edges);CHKERRQ(ierr);
1068b9c6e19dSShri Abhyankar     ierr = PetscCalloc1(newDMnetwork->subnet[j].nvtx,&newDMnetwork->subnet[j].vertices);CHKERRQ(ierr);
1069b9c6e19dSShri Abhyankar     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
1070b9c6e19dSShri Abhyankar        These get updated when the vertices and edges are added. */
1071b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0;
1072b9c6e19dSShri Abhyankar   }
1073b9c6e19dSShri Abhyankar 
1074b9c6e19dSShri Abhyankar   /* Set the vertices and edges in each subnetwork */
1075b9c6e19dSShri Abhyankar   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1076b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
1077b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1078b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++]  = e;
1079b9c6e19dSShri Abhyankar   }
1080b9c6e19dSShri Abhyankar 
1081b9c6e19dSShri Abhyankar   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1082b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
1083b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1084b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++]  = v;
1085b9c6e19dSShri Abhyankar   }
1086b9c6e19dSShri Abhyankar 
108724121865SAdrian Maldonado   /* Destroy point SF */
108824121865SAdrian Maldonado   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
108924121865SAdrian Maldonado 
1090d3464fd4SAdrian Maldonado   ierr = DMDestroy(dm);CHKERRQ(ierr);
1091d3464fd4SAdrian Maldonado   *dm  = newDM;
10925f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10935f2c45f1SShri Abhyankar }
10945f2c45f1SShri Abhyankar 
109524121865SAdrian Maldonado /*@C
109624121865SAdrian Maldonado   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering.
109724121865SAdrian Maldonado 
109824121865SAdrian Maldonado   Input Parameters:
109924121865SAdrian Maldonado + masterSF - the original SF structure
110024121865SAdrian Maldonado - map      - a ISLocalToGlobal mapping that contains the subset of points
110124121865SAdrian Maldonado 
110224121865SAdrian Maldonado   Output Parameters:
110324121865SAdrian Maldonado . subSF    - a subset of the masterSF for the desired subset.
110424121865SAdrian Maldonado */
110524121865SAdrian Maldonado 
110624121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) {
110724121865SAdrian Maldonado 
110824121865SAdrian Maldonado   PetscErrorCode        ierr;
110924121865SAdrian Maldonado   PetscInt              nroots, nleaves, *ilocal_sub;
111024121865SAdrian Maldonado   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
111124121865SAdrian Maldonado   PetscInt              *local_points, *remote_points;
111224121865SAdrian Maldonado   PetscSFNode           *iremote_sub;
111324121865SAdrian Maldonado   const PetscInt        *ilocal;
111424121865SAdrian Maldonado   const PetscSFNode     *iremote;
111524121865SAdrian Maldonado 
111624121865SAdrian Maldonado   PetscFunctionBegin;
111724121865SAdrian Maldonado   ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
111824121865SAdrian Maldonado 
111924121865SAdrian Maldonado   /* Look for leaves that pertain to the subset of points. Get the local ordering */
112024121865SAdrian Maldonado   ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr);
112124121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr);
112224121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
112324121865SAdrian Maldonado     if (ilocal_map[i] != -1) nleaves_sub += 1;
112424121865SAdrian Maldonado   }
112524121865SAdrian Maldonado   /* Re-number ilocal with subset numbering. Need information from roots */
112624121865SAdrian Maldonado   ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr);
112724121865SAdrian Maldonado   for (i = 0; i < nroots; i++) local_points[i] = i;
112824121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr);
112924121865SAdrian Maldonado   ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
113024121865SAdrian Maldonado   ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
113124121865SAdrian Maldonado   /* Fill up graph using local (that is, local to the subset) numbering. */
11324b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr);
11334b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr);
113424121865SAdrian Maldonado   nleaves_sub = 0;
113524121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
113624121865SAdrian Maldonado     if (ilocal_map[i] != -1) {
113724121865SAdrian Maldonado       ilocal_sub[nleaves_sub] = ilocal_map[i];
11384b70a8deSAdrian Maldonado       iremote_sub[nleaves_sub].rank = iremote[i].rank;
113924121865SAdrian Maldonado       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
114024121865SAdrian Maldonado       nleaves_sub += 1;
114124121865SAdrian Maldonado     }
114224121865SAdrian Maldonado   }
114324121865SAdrian Maldonado   ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr);
114424121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr);
114524121865SAdrian Maldonado 
114624121865SAdrian Maldonado   /* Create new subSF */
114724121865SAdrian Maldonado   ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr);
114824121865SAdrian Maldonado   ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr);
11494b70a8deSAdrian Maldonado   ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr);
115024121865SAdrian Maldonado   ierr = PetscFree(ilocal_map);CHKERRQ(ierr);
11514b70a8deSAdrian Maldonado   ierr = PetscFree(iremote_sub);CHKERRQ(ierr);
115224121865SAdrian Maldonado   PetscFunctionReturn(0);
115324121865SAdrian Maldonado }
115424121865SAdrian Maldonado 
11555f2c45f1SShri Abhyankar /*@C
11565f2c45f1SShri Abhyankar   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
11575f2c45f1SShri Abhyankar 
11585f2c45f1SShri Abhyankar   Not Collective
11595f2c45f1SShri Abhyankar 
11605f2c45f1SShri Abhyankar   Input Parameters:
11615f2c45f1SShri Abhyankar + dm - The DMNetwork object
11625f2c45f1SShri Abhyankar - p  - the vertex point
11635f2c45f1SShri Abhyankar 
11645f2c45f1SShri Abhyankar   Output Paramters:
11655f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point
11665f2c45f1SShri Abhyankar - edges  - List of edge points
11675f2c45f1SShri Abhyankar 
11685f2c45f1SShri Abhyankar   Level: intermediate
11695f2c45f1SShri Abhyankar 
11705f2c45f1SShri Abhyankar   Fortran Notes:
11715f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
11725f2c45f1SShri Abhyankar   include petsc.h90 in your code.
11735f2c45f1SShri Abhyankar 
1174d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices
11755f2c45f1SShri Abhyankar @*/
11765f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
11775f2c45f1SShri Abhyankar {
11785f2c45f1SShri Abhyankar   PetscErrorCode ierr;
11795f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
11805f2c45f1SShri Abhyankar 
11815f2c45f1SShri Abhyankar   PetscFunctionBegin;
11825f2c45f1SShri Abhyankar   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
11835f2c45f1SShri Abhyankar   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
11845f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
11855f2c45f1SShri Abhyankar }
11865f2c45f1SShri Abhyankar 
11875f2c45f1SShri Abhyankar /*@C
1188d842c372SHong Zhang   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
11895f2c45f1SShri Abhyankar 
11905f2c45f1SShri Abhyankar   Not Collective
11915f2c45f1SShri Abhyankar 
11925f2c45f1SShri Abhyankar   Input Parameters:
11935f2c45f1SShri Abhyankar + dm - The DMNetwork object
11945f2c45f1SShri Abhyankar - p  - the edge point
11955f2c45f1SShri Abhyankar 
11965f2c45f1SShri Abhyankar   Output Paramters:
11975f2c45f1SShri Abhyankar . vertices  - vertices connected to this edge
11985f2c45f1SShri Abhyankar 
11995f2c45f1SShri Abhyankar   Level: intermediate
12005f2c45f1SShri Abhyankar 
12015f2c45f1SShri Abhyankar   Fortran Notes:
12025f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
12035f2c45f1SShri Abhyankar   include petsc.h90 in your code.
12045f2c45f1SShri Abhyankar 
12055f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
12065f2c45f1SShri Abhyankar @*/
1207d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
12085f2c45f1SShri Abhyankar {
12095f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12105f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
12115f2c45f1SShri Abhyankar 
12125f2c45f1SShri Abhyankar   PetscFunctionBegin;
12135f2c45f1SShri Abhyankar   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
12145f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12155f2c45f1SShri Abhyankar }
12165f2c45f1SShri Abhyankar 
12175f2c45f1SShri Abhyankar /*@
12185f2c45f1SShri Abhyankar   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
12195f2c45f1SShri Abhyankar 
12205f2c45f1SShri Abhyankar   Not Collective
12215f2c45f1SShri Abhyankar 
12225f2c45f1SShri Abhyankar   Input Parameters:
12235f2c45f1SShri Abhyankar + dm - The DMNetwork object
1224*a2b725a8SWilliam Gropp - p  - the vertex point
12255f2c45f1SShri Abhyankar 
12265f2c45f1SShri Abhyankar   Output Parameter:
12275f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point
12285f2c45f1SShri Abhyankar 
12295f2c45f1SShri Abhyankar   Level: intermediate
12305f2c45f1SShri Abhyankar 
1231d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange
12325f2c45f1SShri Abhyankar @*/
12335f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
12345f2c45f1SShri Abhyankar {
12355f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12365f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
12375f2c45f1SShri Abhyankar   PetscInt       offsetg;
12385f2c45f1SShri Abhyankar   PetscSection   sectiong;
12395f2c45f1SShri Abhyankar 
12405f2c45f1SShri Abhyankar   PetscFunctionBegin;
12415f2c45f1SShri Abhyankar   *isghost = PETSC_FALSE;
1242e87a4003SBarry Smith   ierr = DMGetGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
12435f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
12445f2c45f1SShri Abhyankar   if (offsetg < 0) *isghost = PETSC_TRUE;
12455f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12465f2c45f1SShri Abhyankar }
12475f2c45f1SShri Abhyankar 
12485f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm)
12495f2c45f1SShri Abhyankar {
12505f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12515f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
12525f2c45f1SShri Abhyankar 
12535f2c45f1SShri Abhyankar   PetscFunctionBegin;
12545f2c45f1SShri Abhyankar   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
12555f2c45f1SShri Abhyankar   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
12565f2c45f1SShri Abhyankar 
1257e87a4003SBarry Smith   ierr = DMSetSection(network->plex,network->DofSection);CHKERRQ(ierr);
1258e87a4003SBarry Smith   ierr = DMGetGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
12595f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12605f2c45f1SShri Abhyankar }
12615f2c45f1SShri Abhyankar 
12621ad426b7SHong Zhang /*@
126317df6e9eSHong Zhang     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
12641ad426b7SHong Zhang                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
12651ad426b7SHong Zhang 
12661ad426b7SHong Zhang     Collective
12671ad426b7SHong Zhang 
12681ad426b7SHong Zhang     Input Parameters:
126983b2e829SHong Zhang +   dm - The DMNetwork object
127083b2e829SHong Zhang .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
127183b2e829SHong Zhang -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
12721ad426b7SHong Zhang 
12731ad426b7SHong Zhang     Level: intermediate
12741ad426b7SHong Zhang 
12751ad426b7SHong Zhang @*/
127683b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
12771ad426b7SHong Zhang {
12781ad426b7SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
12798675203cSHong Zhang   PetscErrorCode ierr;
128066b4e0ffSHong Zhang   PetscInt       nVertices = network->nVertices;
12811ad426b7SHong Zhang 
12821ad426b7SHong Zhang   PetscFunctionBegin;
128383b2e829SHong Zhang   network->userEdgeJacobian   = eflg;
128483b2e829SHong Zhang   network->userVertexJacobian = vflg;
12858675203cSHong Zhang 
12868675203cSHong Zhang   if (eflg && !network->Je) {
12878675203cSHong Zhang     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
12888675203cSHong Zhang   }
12898675203cSHong Zhang 
129066b4e0ffSHong Zhang   if (vflg && !network->Jv && nVertices) {
12918675203cSHong Zhang     PetscInt       i,*vptr,nedges,vStart=network->vStart;
129266b4e0ffSHong Zhang     PetscInt       nedges_total;
12938675203cSHong Zhang     const PetscInt *edges;
12948675203cSHong Zhang 
12958675203cSHong Zhang     /* count nvertex_total */
12968675203cSHong Zhang     nedges_total = 0;
12978675203cSHong Zhang     ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr);
12988675203cSHong Zhang 
12998675203cSHong Zhang     vptr[0] = 0;
13008675203cSHong Zhang     for (i=0; i<nVertices; i++) {
13018675203cSHong Zhang       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
13028675203cSHong Zhang       nedges_total += nedges;
13038675203cSHong Zhang       vptr[i+1] = vptr[i] + 2*nedges + 1;
13048675203cSHong Zhang     }
13058675203cSHong Zhang 
13068675203cSHong Zhang     ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr);
13078675203cSHong Zhang     network->Jvptr = vptr;
13088675203cSHong Zhang   }
13091ad426b7SHong Zhang   PetscFunctionReturn(0);
13101ad426b7SHong Zhang }
13111ad426b7SHong Zhang 
13121ad426b7SHong Zhang /*@
131383b2e829SHong Zhang     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
131483b2e829SHong Zhang 
131583b2e829SHong Zhang     Not Collective
131683b2e829SHong Zhang 
131783b2e829SHong Zhang     Input Parameters:
131883b2e829SHong Zhang +   dm - The DMNetwork object
131983b2e829SHong Zhang .   p  - the edge point
13203e97b6e8SHong Zhang -   J - array (size = 3) of Jacobian submatrices for this edge point:
13213e97b6e8SHong Zhang         J[0]: this edge
1322d842c372SHong Zhang         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()
132383b2e829SHong Zhang 
132483b2e829SHong Zhang     Level: intermediate
132583b2e829SHong Zhang 
132683b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix
132783b2e829SHong Zhang @*/
132883b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
132983b2e829SHong Zhang {
133083b2e829SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
133183b2e829SHong Zhang 
133283b2e829SHong Zhang   PetscFunctionBegin;
13338675203cSHong Zhang   if (!network->Je) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
13348675203cSHong Zhang 
13358675203cSHong Zhang   if (J) {
1336883e35e8SHong Zhang     network->Je[3*p]   = J[0];
1337883e35e8SHong Zhang     network->Je[3*p+1] = J[1];
1338883e35e8SHong Zhang     network->Je[3*p+2] = J[2];
13398675203cSHong Zhang   }
134083b2e829SHong Zhang   PetscFunctionReturn(0);
134183b2e829SHong Zhang }
134283b2e829SHong Zhang 
134383b2e829SHong Zhang /*@
134476ddfea5SHong Zhang     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
13451ad426b7SHong Zhang 
13461ad426b7SHong Zhang     Not Collective
13471ad426b7SHong Zhang 
13481ad426b7SHong Zhang     Input Parameters:
13491ad426b7SHong Zhang +   dm - The DMNetwork object
13501ad426b7SHong Zhang .   p  - the vertex point
13513e97b6e8SHong Zhang -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
13523e97b6e8SHong Zhang         J[0]:       this vertex
13533e97b6e8SHong Zhang         J[1+2*i]:   i-th supporting edge
13543e97b6e8SHong Zhang         J[1+2*i+1]: i-th connected vertex
13551ad426b7SHong Zhang 
13561ad426b7SHong Zhang     Level: intermediate
13571ad426b7SHong Zhang 
135883b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix
13591ad426b7SHong Zhang @*/
1360883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
13615f2c45f1SShri Abhyankar {
13625f2c45f1SShri Abhyankar   PetscErrorCode ierr;
13635f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
13648675203cSHong Zhang   PetscInt       i,*vptr,nedges,vStart=network->vStart;
1365883e35e8SHong Zhang   const PetscInt *edges;
13665f2c45f1SShri Abhyankar 
13675f2c45f1SShri Abhyankar   PetscFunctionBegin;
13688675203cSHong Zhang   if (!network->Jv) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
1369883e35e8SHong Zhang 
13708675203cSHong Zhang   if (J) {
1371883e35e8SHong Zhang     vptr = network->Jvptr;
13723e97b6e8SHong Zhang     network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
13733e97b6e8SHong Zhang 
13743e97b6e8SHong Zhang     /* Set Jacobian for each supporting edge and connected vertex */
1375883e35e8SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
1376883e35e8SHong Zhang     for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
13778675203cSHong Zhang   }
1378883e35e8SHong Zhang   PetscFunctionReturn(0);
1379883e35e8SHong Zhang }
1380883e35e8SHong Zhang 
1381e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
13825cf7da58SHong Zhang {
13835cf7da58SHong Zhang   PetscErrorCode ierr;
13845cf7da58SHong Zhang   PetscInt       j;
13855cf7da58SHong Zhang   PetscScalar    val=(PetscScalar)ncols;
13865cf7da58SHong Zhang 
13875cf7da58SHong Zhang   PetscFunctionBegin;
13885cf7da58SHong Zhang   if (!ghost) {
13895cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
13905cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
13915cf7da58SHong Zhang     }
13925cf7da58SHong Zhang   } else {
13935cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
13945cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
13955cf7da58SHong Zhang     }
13965cf7da58SHong Zhang   }
13975cf7da58SHong Zhang   PetscFunctionReturn(0);
13985cf7da58SHong Zhang }
13995cf7da58SHong Zhang 
1400e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
14015cf7da58SHong Zhang {
14025cf7da58SHong Zhang   PetscErrorCode ierr;
14035cf7da58SHong Zhang   PetscInt       j,ncols_u;
14045cf7da58SHong Zhang   PetscScalar    val;
14055cf7da58SHong Zhang 
14065cf7da58SHong Zhang   PetscFunctionBegin;
14075cf7da58SHong Zhang   if (!ghost) {
14085cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
14095cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
14105cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
14115cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
14125cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
14135cf7da58SHong Zhang     }
14145cf7da58SHong Zhang   } else {
14155cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
14165cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
14175cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
14185cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
14195cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
14205cf7da58SHong Zhang     }
14215cf7da58SHong Zhang   }
14225cf7da58SHong Zhang   PetscFunctionReturn(0);
14235cf7da58SHong Zhang }
14245cf7da58SHong Zhang 
1425e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
14265cf7da58SHong Zhang {
14275cf7da58SHong Zhang   PetscErrorCode ierr;
14285cf7da58SHong Zhang 
14295cf7da58SHong Zhang   PetscFunctionBegin;
14305cf7da58SHong Zhang   if (Ju) {
14315cf7da58SHong Zhang     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
14325cf7da58SHong Zhang   } else {
14335cf7da58SHong Zhang     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
14345cf7da58SHong Zhang   }
14355cf7da58SHong Zhang   PetscFunctionReturn(0);
14365cf7da58SHong Zhang }
14375cf7da58SHong Zhang 
1438e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1439883e35e8SHong Zhang {
1440883e35e8SHong Zhang   PetscErrorCode ierr;
1441883e35e8SHong Zhang   PetscInt       j,*cols;
1442883e35e8SHong Zhang   PetscScalar    *zeros;
1443883e35e8SHong Zhang 
1444883e35e8SHong Zhang   PetscFunctionBegin;
1445883e35e8SHong Zhang   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
1446883e35e8SHong Zhang   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1447883e35e8SHong Zhang   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
1448883e35e8SHong Zhang   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
14491ad426b7SHong Zhang   PetscFunctionReturn(0);
14501ad426b7SHong Zhang }
1451a4e85ca8SHong Zhang 
1452e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
14533e97b6e8SHong Zhang {
14543e97b6e8SHong Zhang   PetscErrorCode ierr;
14553e97b6e8SHong Zhang   PetscInt       j,M,N,row,col,ncols_u;
14563e97b6e8SHong Zhang   const PetscInt *cols;
14573e97b6e8SHong Zhang   PetscScalar    zero=0.0;
14583e97b6e8SHong Zhang 
14593e97b6e8SHong Zhang   PetscFunctionBegin;
14603e97b6e8SHong Zhang   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
14613e97b6e8SHong 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);
14623e97b6e8SHong Zhang 
14633e97b6e8SHong Zhang   for (row=0; row<nrows; row++) {
14643e97b6e8SHong Zhang     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
14653e97b6e8SHong Zhang     for (j=0; j<ncols_u; j++) {
14663e97b6e8SHong Zhang       col = cols[j] + cstart;
14673e97b6e8SHong Zhang       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
14683e97b6e8SHong Zhang     }
14693e97b6e8SHong Zhang     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
14703e97b6e8SHong Zhang   }
14713e97b6e8SHong Zhang   PetscFunctionReturn(0);
14723e97b6e8SHong Zhang }
14731ad426b7SHong Zhang 
1474e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1475a4e85ca8SHong Zhang {
1476a4e85ca8SHong Zhang   PetscErrorCode ierr;
1477f4431b8cSHong Zhang 
1478a4e85ca8SHong Zhang   PetscFunctionBegin;
1479a4e85ca8SHong Zhang   if (Ju) {
1480a4e85ca8SHong Zhang     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1481a4e85ca8SHong Zhang   } else {
1482a4e85ca8SHong Zhang     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1483a4e85ca8SHong Zhang   }
1484a4e85ca8SHong Zhang   PetscFunctionReturn(0);
1485a4e85ca8SHong Zhang }
1486a4e85ca8SHong Zhang 
148724121865SAdrian 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.
148824121865SAdrian Maldonado */
148924121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
149024121865SAdrian Maldonado {
149124121865SAdrian Maldonado   PetscErrorCode ierr;
149224121865SAdrian Maldonado   PetscInt       i, size, dof;
149324121865SAdrian Maldonado   PetscInt       *glob2loc;
149424121865SAdrian Maldonado 
149524121865SAdrian Maldonado   PetscFunctionBegin;
149624121865SAdrian Maldonado   ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr);
149724121865SAdrian Maldonado   ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr);
149824121865SAdrian Maldonado 
149924121865SAdrian Maldonado   for (i = 0; i < size; i++) {
150024121865SAdrian Maldonado     ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr);
150124121865SAdrian Maldonado     dof = (dof >= 0) ? dof : -(dof + 1);
150224121865SAdrian Maldonado     glob2loc[i] = dof;
150324121865SAdrian Maldonado   }
150424121865SAdrian Maldonado 
150524121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
150624121865SAdrian Maldonado #if 0
150724121865SAdrian Maldonado   ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
150824121865SAdrian Maldonado #endif
150924121865SAdrian Maldonado   PetscFunctionReturn(0);
151024121865SAdrian Maldonado }
151124121865SAdrian Maldonado 
151201ad2aeeSHong Zhang #include <petsc/private/matimpl.h>
15131ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
15141ad426b7SHong Zhang {
15151ad426b7SHong Zhang   PetscErrorCode ierr;
151624121865SAdrian Maldonado   PetscMPIInt    rank, size;
15171ad426b7SHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
1518a4e85ca8SHong Zhang   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1519840c2264SHong Zhang   PetscInt       cstart,ncols,j,e,v;
152024121865SAdrian Maldonado   PetscBool      ghost,ghost_vc,ghost2,isNest;
1521a4e85ca8SHong Zhang   Mat            Juser;
1522bfbc38dcSHong Zhang   PetscSection   sectionGlobal;
1523447d78afSSatish Balay   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
1524a4e85ca8SHong Zhang   const PetscInt *edges,*cone;
15255cf7da58SHong Zhang   MPI_Comm       comm;
152624121865SAdrian Maldonado   MatType        mtype;
15275cf7da58SHong Zhang   Vec            vd_nz,vo_nz;
15285cf7da58SHong Zhang   PetscInt       *dnnz,*onnz;
15295cf7da58SHong Zhang   PetscScalar    *vdnz,*vonz;
15301ad426b7SHong Zhang 
15311ad426b7SHong Zhang   PetscFunctionBegin;
153224121865SAdrian Maldonado   mtype = dm->mattype;
153324121865SAdrian Maldonado   ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr);
153424121865SAdrian Maldonado 
153524121865SAdrian Maldonado   if (isNest) {
15360731d606SHong Zhang     /* ierr = DMCreateMatrix_Network_Nest(); */
153724121865SAdrian Maldonado     PetscInt   eDof, vDof;
153824121865SAdrian Maldonado     Mat        j11, j12, j21, j22, bA[2][2];
153924121865SAdrian Maldonado     ISLocalToGlobalMapping eISMap, vISMap;
154024121865SAdrian Maldonado 
154124121865SAdrian Maldonado     ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
154224121865SAdrian Maldonado     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
154324121865SAdrian Maldonado     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
154424121865SAdrian Maldonado 
154524121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr);
154624121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr);
154724121865SAdrian Maldonado 
154801ad2aeeSHong Zhang     ierr = MatCreate(comm, &j11);CHKERRQ(ierr);
154924121865SAdrian Maldonado     ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
155024121865SAdrian Maldonado     ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr);
155124121865SAdrian Maldonado 
155201ad2aeeSHong Zhang     ierr = MatCreate(comm, &j12);CHKERRQ(ierr);
155324121865SAdrian Maldonado     ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr);
155424121865SAdrian Maldonado     ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr);
155524121865SAdrian Maldonado 
155601ad2aeeSHong Zhang     ierr = MatCreate(comm, &j21);CHKERRQ(ierr);
155724121865SAdrian Maldonado     ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
155824121865SAdrian Maldonado     ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr);
155924121865SAdrian Maldonado 
156001ad2aeeSHong Zhang     ierr = MatCreate(comm, &j22);CHKERRQ(ierr);
156124121865SAdrian Maldonado     ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
156224121865SAdrian Maldonado     ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr);
156324121865SAdrian Maldonado 
15643f6a6bdaSHong Zhang     bA[0][0] = j11;
15653f6a6bdaSHong Zhang     bA[0][1] = j12;
15663f6a6bdaSHong Zhang     bA[1][0] = j21;
15673f6a6bdaSHong Zhang     bA[1][1] = j22;
156824121865SAdrian Maldonado 
156924121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr);
157024121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr);
157124121865SAdrian Maldonado 
157224121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr);
157324121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr);
157424121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr);
157524121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr);
157624121865SAdrian Maldonado 
157724121865SAdrian Maldonado     ierr = MatSetUp(j11);CHKERRQ(ierr);
157824121865SAdrian Maldonado     ierr = MatSetUp(j12);CHKERRQ(ierr);
157924121865SAdrian Maldonado     ierr = MatSetUp(j21);CHKERRQ(ierr);
158024121865SAdrian Maldonado     ierr = MatSetUp(j22);CHKERRQ(ierr);
158124121865SAdrian Maldonado 
158201ad2aeeSHong Zhang     ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr);
158324121865SAdrian Maldonado     ierr = MatSetUp(*J);CHKERRQ(ierr);
158424121865SAdrian Maldonado     ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr);
158524121865SAdrian Maldonado     ierr = MatDestroy(&j11);CHKERRQ(ierr);
158624121865SAdrian Maldonado     ierr = MatDestroy(&j12);CHKERRQ(ierr);
158724121865SAdrian Maldonado     ierr = MatDestroy(&j21);CHKERRQ(ierr);
158824121865SAdrian Maldonado     ierr = MatDestroy(&j22);CHKERRQ(ierr);
158924121865SAdrian Maldonado 
159024121865SAdrian Maldonado     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
159124121865SAdrian Maldonado     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
159224121865SAdrian Maldonado 
159324121865SAdrian Maldonado     /* Free structures */
159424121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr);
159524121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr);
159624121865SAdrian Maldonado 
159724121865SAdrian Maldonado     PetscFunctionReturn(0);
159824121865SAdrian Maldonado   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1599a4e85ca8SHong Zhang     /* user does not provide Jacobian blocks */
1600bfbc38dcSHong Zhang     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1601bfbc38dcSHong Zhang     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
16021ad426b7SHong Zhang     PetscFunctionReturn(0);
16031ad426b7SHong Zhang   }
16041ad426b7SHong Zhang 
1605bfbc38dcSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
1606e87a4003SBarry Smith   ierr = DMGetGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1607bfbc38dcSHong Zhang   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1608bfbc38dcSHong Zhang   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
16092a945128SHong Zhang 
16102a945128SHong Zhang   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
16112a945128SHong Zhang   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
161289898e50SHong Zhang 
161389898e50SHong Zhang   /* (1) Set matrix preallocation */
161489898e50SHong Zhang   /*------------------------------*/
1615840c2264SHong Zhang   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1616840c2264SHong Zhang   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1617840c2264SHong Zhang   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1618840c2264SHong Zhang   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1619840c2264SHong Zhang   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1620840c2264SHong Zhang   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1621840c2264SHong Zhang 
162289898e50SHong Zhang   /* Set preallocation for edges */
162389898e50SHong Zhang   /*-----------------------------*/
1624840c2264SHong Zhang   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1625840c2264SHong Zhang 
1626bdcb62a2SHong Zhang   ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr);
1627840c2264SHong Zhang   for (e=eStart; e<eEnd; e++) {
1628840c2264SHong Zhang     /* Get row indices */
1629840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1630840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1631840c2264SHong Zhang     if (nrows) {
1632840c2264SHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1633840c2264SHong Zhang 
16345cf7da58SHong Zhang       /* Set preallocation for conntected vertices */
1635d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1636840c2264SHong Zhang       for (v=0; v<2; v++) {
1637840c2264SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1638840c2264SHong Zhang 
16398675203cSHong Zhang         if (network->Je) {
1640840c2264SHong Zhang           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
16418675203cSHong Zhang         } else Juser = NULL;
1642840c2264SHong Zhang         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
16435cf7da58SHong Zhang         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1644840c2264SHong Zhang       }
1645840c2264SHong Zhang 
164689898e50SHong Zhang       /* Set preallocation for edge self */
1647840c2264SHong Zhang       cstart = rstart;
16488675203cSHong Zhang       if (network->Je) {
1649840c2264SHong Zhang         Juser = network->Je[3*e]; /* Jacobian(e,e) */
16508675203cSHong Zhang       } else Juser = NULL;
16515cf7da58SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1652840c2264SHong Zhang     }
1653840c2264SHong Zhang   }
1654840c2264SHong Zhang 
165589898e50SHong Zhang   /* Set preallocation for vertices */
165689898e50SHong Zhang   /*--------------------------------*/
1657840c2264SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
16588675203cSHong Zhang   if (vEnd - vStart) vptr = network->Jvptr;
1659840c2264SHong Zhang 
1660840c2264SHong Zhang   for (v=vStart; v<vEnd; v++) {
1661840c2264SHong Zhang     /* Get row indices */
1662840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1663840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1664840c2264SHong Zhang     if (!nrows) continue;
1665840c2264SHong Zhang 
1666bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1667bdcb62a2SHong Zhang     if (ghost) {
1668bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1669bdcb62a2SHong Zhang     } else {
1670bdcb62a2SHong Zhang       rows_v = rows;
1671bdcb62a2SHong Zhang     }
1672bdcb62a2SHong Zhang 
1673bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1674840c2264SHong Zhang 
1675840c2264SHong Zhang     /* Get supporting edges and connected vertices */
1676840c2264SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1677840c2264SHong Zhang 
1678840c2264SHong Zhang     for (e=0; e<nedges; e++) {
1679840c2264SHong Zhang       /* Supporting edges */
1680840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1681840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1682840c2264SHong Zhang 
16838675203cSHong Zhang       if (network->Jv) {
1684840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
16858675203cSHong Zhang       } else Juser = NULL;
1686bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1687840c2264SHong Zhang 
1688840c2264SHong Zhang       /* Connected vertices */
1689d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
1690840c2264SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
1691840c2264SHong Zhang       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1692840c2264SHong Zhang 
1693840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1694840c2264SHong Zhang 
16958675203cSHong Zhang       if (network->Jv) {
1696840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
16978675203cSHong Zhang       } else Juser = NULL;
1698e102a522SHong Zhang       if (ghost_vc||ghost) {
1699e102a522SHong Zhang         ghost2 = PETSC_TRUE;
1700e102a522SHong Zhang       } else {
1701e102a522SHong Zhang         ghost2 = PETSC_FALSE;
1702e102a522SHong Zhang       }
1703e102a522SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr);
1704840c2264SHong Zhang     }
1705840c2264SHong Zhang 
170689898e50SHong Zhang     /* Set preallocation for vertex self */
1707840c2264SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1708840c2264SHong Zhang     if (!ghost) {
1709840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
17108675203cSHong Zhang       if (network->Jv) {
1711840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
17128675203cSHong Zhang       } else Juser = NULL;
1713bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1714840c2264SHong Zhang     }
1715bdcb62a2SHong Zhang     if (ghost) {
1716bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1717bdcb62a2SHong Zhang     }
1718840c2264SHong Zhang   }
1719840c2264SHong Zhang 
1720840c2264SHong Zhang   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1721840c2264SHong Zhang   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
17225cf7da58SHong Zhang 
17235cf7da58SHong Zhang   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
17245cf7da58SHong Zhang 
17255cf7da58SHong Zhang   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1726840c2264SHong Zhang   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1727840c2264SHong Zhang 
1728840c2264SHong Zhang   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1729840c2264SHong Zhang   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1730840c2264SHong Zhang   for (j=0; j<localSize; j++) {
1731e102a522SHong Zhang     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1732e102a522SHong Zhang     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1733840c2264SHong Zhang   }
1734840c2264SHong Zhang   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1735840c2264SHong Zhang   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1736840c2264SHong Zhang   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1737840c2264SHong Zhang   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1738840c2264SHong Zhang 
17395cf7da58SHong Zhang   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
17405cf7da58SHong Zhang   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
17415cf7da58SHong Zhang   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
17425cf7da58SHong Zhang 
17435cf7da58SHong Zhang   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
17445cf7da58SHong Zhang 
174589898e50SHong Zhang   /* (2) Set matrix entries for edges */
174689898e50SHong Zhang   /*----------------------------------*/
17471ad426b7SHong Zhang   for (e=eStart; e<eEnd; e++) {
1748bfbc38dcSHong Zhang     /* Get row indices */
17491ad426b7SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
175017df6e9eSHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
17514b976069SHong Zhang     if (nrows) {
175217df6e9eSHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
17531ad426b7SHong Zhang 
1754bfbc38dcSHong Zhang       /* Set matrix entries for conntected vertices */
1755d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1756bfbc38dcSHong Zhang       for (v=0; v<2; v++) {
1757bfbc38dcSHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1758883e35e8SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
17593e97b6e8SHong Zhang 
17608675203cSHong Zhang         if (network->Je) {
1761a4e85ca8SHong Zhang           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
17628675203cSHong Zhang         } else Juser = NULL;
1763a4e85ca8SHong Zhang         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1764bfbc38dcSHong Zhang       }
176517df6e9eSHong Zhang 
1766bfbc38dcSHong Zhang       /* Set matrix entries for edge self */
17673e97b6e8SHong Zhang       cstart = rstart;
17688675203cSHong Zhang       if (network->Je) {
1769a4e85ca8SHong Zhang         Juser = network->Je[3*e]; /* Jacobian(e,e) */
17708675203cSHong Zhang       } else Juser = NULL;
1771a4e85ca8SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
17721ad426b7SHong Zhang     }
17734b976069SHong Zhang   }
17741ad426b7SHong Zhang 
1775bfbc38dcSHong Zhang   /* Set matrix entries for vertices */
177683b2e829SHong Zhang   /*---------------------------------*/
17771ad426b7SHong Zhang   for (v=vStart; v<vEnd; v++) {
1778bfbc38dcSHong Zhang     /* Get row indices */
1779596e729fSHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1780596e729fSHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
17814b976069SHong Zhang     if (!nrows) continue;
1782596e729fSHong Zhang 
1783bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1784bdcb62a2SHong Zhang     if (ghost) {
1785bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1786bdcb62a2SHong Zhang     } else {
1787bdcb62a2SHong Zhang       rows_v = rows;
1788bdcb62a2SHong Zhang     }
1789bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1790596e729fSHong Zhang 
1791bfbc38dcSHong Zhang     /* Get supporting edges and connected vertices */
1792596e729fSHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1793596e729fSHong Zhang 
1794596e729fSHong Zhang     for (e=0; e<nedges; e++) {
1795bfbc38dcSHong Zhang       /* Supporting edges */
1796596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1797596e729fSHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1798596e729fSHong Zhang 
17998675203cSHong Zhang       if (network->Jv) {
1800a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
18018675203cSHong Zhang       } else Juser = NULL;
1802bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1803596e729fSHong Zhang 
1804bfbc38dcSHong Zhang       /* Connected vertices */
1805d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
18062a945128SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
18072a945128SHong Zhang 
180844aca652SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
180944aca652SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1810a4e85ca8SHong Zhang 
18118675203cSHong Zhang       if (network->Jv) {
1812a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
18138675203cSHong Zhang       } else Juser = NULL;
1814bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1815596e729fSHong Zhang     }
1816596e729fSHong Zhang 
1817bfbc38dcSHong Zhang     /* Set matrix entries for vertex self */
18181ad426b7SHong Zhang     if (!ghost) {
1819596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
18208675203cSHong Zhang       if (network->Jv) {
1821a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
18228675203cSHong Zhang       } else Juser = NULL;
1823bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr);
1824bdcb62a2SHong Zhang     }
1825bdcb62a2SHong Zhang     if (ghost) {
1826bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1827bdcb62a2SHong Zhang     }
18281ad426b7SHong Zhang   }
1829a4e85ca8SHong Zhang   ierr = PetscFree(rows);CHKERRQ(ierr);
1830bdcb62a2SHong Zhang 
18311ad426b7SHong Zhang   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
18321ad426b7SHong Zhang   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1833dd6f46cdSHong Zhang 
18345f2c45f1SShri Abhyankar   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
18355f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
18365f2c45f1SShri Abhyankar }
18375f2c45f1SShri Abhyankar 
18385f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm)
18395f2c45f1SShri Abhyankar {
18405f2c45f1SShri Abhyankar   PetscErrorCode ierr;
18415f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
18422727e31bSShri Abhyankar   PetscInt       j;
18435f2c45f1SShri Abhyankar 
18445f2c45f1SShri Abhyankar   PetscFunctionBegin;
18458415c774SShri Abhyankar   if (--network->refct > 0) PetscFunctionReturn(0);
184683b2e829SHong Zhang   if (network->Je) {
184783b2e829SHong Zhang     ierr = PetscFree(network->Je);CHKERRQ(ierr);
184883b2e829SHong Zhang   }
184983b2e829SHong Zhang   if (network->Jv) {
1850883e35e8SHong Zhang     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
185183b2e829SHong Zhang     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
18521ad426b7SHong Zhang   }
185313c2a604SAdrian Maldonado 
185413c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr);
185513c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr);
185613c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr);
185713c2a604SAdrian Maldonado   if (network->vertex.sf) {
185813c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr);
185913c2a604SAdrian Maldonado   }
186013c2a604SAdrian Maldonado   /* edge */
186113c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr);
186213c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr);
186313c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr);
186413c2a604SAdrian Maldonado   if (network->edge.sf) {
186513c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr);
186613c2a604SAdrian Maldonado   }
18675f2c45f1SShri Abhyankar   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
18685f2c45f1SShri Abhyankar   network->edges = NULL;
18695f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
18705f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
187183b2e829SHong Zhang 
18722727e31bSShri Abhyankar   for(j=0; j < network->nsubnet; j++) {
18732727e31bSShri Abhyankar     ierr = PetscFree(network->subnet[j].edges);CHKERRQ(ierr);
18742727e31bSShri Abhyankar     ierr = PetscFree(network->subnet[j].vertices);CHKERRQ(ierr);
18752727e31bSShri Abhyankar   }
1876e2aaf10cSShri Abhyankar   ierr = PetscFree(network->subnet);CHKERRQ(ierr);
18775f2c45f1SShri Abhyankar   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
18785f2c45f1SShri Abhyankar   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
18795f2c45f1SShri Abhyankar   ierr = PetscFree(network->header);CHKERRQ(ierr);
18805f2c45f1SShri Abhyankar   ierr = PetscFree(network);CHKERRQ(ierr);
18815f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
18825f2c45f1SShri Abhyankar }
18835f2c45f1SShri Abhyankar 
18845f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm,PetscViewer viewer)
18855f2c45f1SShri Abhyankar {
18865f2c45f1SShri Abhyankar   PetscErrorCode ierr;
18875f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
18885f2c45f1SShri Abhyankar 
18895f2c45f1SShri Abhyankar   PetscFunctionBegin;
18905f2c45f1SShri Abhyankar   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
18915f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
18925f2c45f1SShri Abhyankar }
18935f2c45f1SShri Abhyankar 
18945f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
18955f2c45f1SShri Abhyankar {
18965f2c45f1SShri Abhyankar   PetscErrorCode ierr;
18975f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
18985f2c45f1SShri Abhyankar 
18995f2c45f1SShri Abhyankar   PetscFunctionBegin;
19005f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
19015f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19025f2c45f1SShri Abhyankar }
19035f2c45f1SShri Abhyankar 
19045f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
19055f2c45f1SShri Abhyankar {
19065f2c45f1SShri Abhyankar   PetscErrorCode ierr;
19075f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
19085f2c45f1SShri Abhyankar 
19095f2c45f1SShri Abhyankar   PetscFunctionBegin;
19105f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
19115f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19125f2c45f1SShri Abhyankar }
19135f2c45f1SShri Abhyankar 
19145f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
19155f2c45f1SShri Abhyankar {
19165f2c45f1SShri Abhyankar   PetscErrorCode ierr;
19175f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
19185f2c45f1SShri Abhyankar 
19195f2c45f1SShri Abhyankar   PetscFunctionBegin;
19205f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
19215f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19225f2c45f1SShri Abhyankar }
19235f2c45f1SShri Abhyankar 
19245f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
19255f2c45f1SShri Abhyankar {
19265f2c45f1SShri Abhyankar   PetscErrorCode ierr;
19275f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
19285f2c45f1SShri Abhyankar 
19295f2c45f1SShri Abhyankar   PetscFunctionBegin;
19305f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
19315f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19325f2c45f1SShri Abhyankar }
1933