xref: /petsc/src/dm/impls/network/network.c (revision e3e6898962c06b8733ca06fc25d789138d5927ee)
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:
105*e3e68989SHong Zhang + dm - the dm object
106*e3e68989SHong Zhang . edgelist - list of edges for each subnetwork
107*e3e68989SHong 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 
115*e3e68989SHong Zhang   Example usage:
116*e3e68989SHong Zhang   Consider the following 2 separate networks and a coupling network:
117*e3e68989SHong Zhang 
118*e3e68989SHong Zhang .vb
119*e3e68989SHong Zhang  network 0: v0 -> v1 -> v2 -> v3
120*e3e68989SHong Zhang  network 1: v1 -> v2 -> v0
121*e3e68989SHong Zhang  coupling network: network 1: v2 -> network 0: v0
122*e3e68989SHong Zhang .ve
123*e3e68989SHong Zhang 
124*e3e68989SHong Zhang  The resulting input
125*e3e68989SHong Zhang    edgelist[0] = [0 1 | 1 2 | 2 3];
126*e3e68989SHong Zhang    edgelist[1] = [1 2 | 2 0]
127*e3e68989SHong Zhang    edgelistCouple[0] = [(network)1 (v)2 (network)0 (v)0].
128*e3e68989SHong Zhang 
1295f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes
1305f2c45f1SShri Abhyankar @*/
131*e3e68989SHong Zhang PetscErrorCode DMNetworkSetEdgeList(DM dm,int *edgelist[],int *edgelistCouple[])
1325f2c45f1SShri Abhyankar {
1335f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network*) dm->data;
134*e3e68989SHong Zhang   PetscInt   i,nsubnet,ncsubnet=network->ncsubnet;
1355f2c45f1SShri Abhyankar 
1365f2c45f1SShri Abhyankar   PetscFunctionBegin;
137*e3e68989SHong Zhang   nsubnet = network->nsubnet - ncsubnet;
138*e3e68989SHong Zhang   for(i=0; i < nsubnet; i++) {
139e2aaf10cSShri Abhyankar     network->subnet[i].edgelist = edgelist[i];
140e2aaf10cSShri Abhyankar   }
141*e3e68989SHong Zhang   if (edgelistCouple) {
142*e3e68989SHong Zhang     PetscInt j;
143*e3e68989SHong Zhang     j = 0;
144*e3e68989SHong Zhang     nsubnet = network->nsubnet;
145*e3e68989SHong Zhang     while (i < nsubnet) network->subnet[i++].edgelist = edgelistCouple[j++];
146*e3e68989SHong 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 */
1735f2c45f1SShri Abhyankar   PetscInt       numCorners=2;
1745f2c45f1SShri Abhyankar   PetscInt       spacedim=2;
1755f2c45f1SShri Abhyankar   double         *vertexcoords=NULL;
176e2aaf10cSShri Abhyankar   PetscInt       i,j;
1775f2c45f1SShri Abhyankar   PetscInt       ndata;
178e2aaf10cSShri Abhyankar   PetscInt       ctr=0;
1795f2c45f1SShri Abhyankar 
1805f2c45f1SShri Abhyankar   PetscFunctionBegin;
1816fefedf4SHong Zhang   if (network->nVertices) {
1826fefedf4SHong Zhang     ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr);
1835f2c45f1SShri Abhyankar   }
184e2aaf10cSShri Abhyankar 
185e2aaf10cSShri Abhyankar   /* Create the edgelist for the network by concatenating edgelists of the subnetworks */
186e2aaf10cSShri Abhyankar   ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr);
187e2aaf10cSShri Abhyankar   for(i=0; i < network->nsubnet; i++) {
188e2aaf10cSShri Abhyankar     for(j = 0; j < network->subnet[i].nedge; j++) {
189e2aaf10cSShri Abhyankar       network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
190e2aaf10cSShri Abhyankar       network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
191e2aaf10cSShri Abhyankar       ctr++;
192e2aaf10cSShri Abhyankar     }
193e2aaf10cSShri Abhyankar   }
194e2aaf10cSShri Abhyankar 
19561de3474SHong Zhang #if 0
196e2aaf10cSShri Abhyankar   for(i=0; i < network->nEdges; i++) {
197e2aaf10cSShri Abhyankar     ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr);
198e2aaf10cSShri Abhyankar   }
199e2aaf10cSShri Abhyankar #endif
200e2aaf10cSShri Abhyankar 
2016fefedf4SHong Zhang   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
2026fefedf4SHong Zhang   if (network->nVertices) {
2035f2c45f1SShri Abhyankar     ierr = PetscFree(vertexcoords);CHKERRQ(ierr);
2045f2c45f1SShri Abhyankar   }
205e2aaf10cSShri Abhyankar   ierr = PetscFree(network->edges);CHKERRQ(ierr);
206e2aaf10cSShri Abhyankar 
2075f2c45f1SShri Abhyankar   ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr);
2085f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr);
2095f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr);
2105f2c45f1SShri Abhyankar 
2115f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr);
2125f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr);
2135f2c45f1SShri Abhyankar   ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr);
2145f2c45f1SShri Abhyankar   ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr);
2155f2c45f1SShri Abhyankar 
2162727e31bSShri Abhyankar   /* Create vertices and edges array for the subnetworks */
2172727e31bSShri Abhyankar   for(j=0; j < network->nsubnet; j++) {
2182727e31bSShri Abhyankar     ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr);
2192727e31bSShri Abhyankar     ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr);
2202727e31bSShri Abhyankar     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
2212727e31bSShri Abhyankar        These get updated when the vertices and edges are added. */
2222727e31bSShri Abhyankar     network->subnet[j].nvtx = network->subnet[j].nedge = 0;
2232727e31bSShri Abhyankar   }
2242727e31bSShri Abhyankar 
2255f2c45f1SShri Abhyankar   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
2266caa05f4SBarry Smith   ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr);
227e2aaf10cSShri Abhyankar   for(i=network->eStart; i < network->eEnd; i++) {
228e2aaf10cSShri Abhyankar     network->header[i].index = i;   /* Global edge number */
229e2aaf10cSShri Abhyankar     for(j=0; j < network->nsubnet; j++) {
230e2aaf10cSShri Abhyankar       if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) {
231e2aaf10cSShri Abhyankar 	network->header[i].subnetid = j; /* Subnetwork id */
2322727e31bSShri Abhyankar 	network->subnet[j].edges[network->subnet[j].nedge++] = i;
233e2aaf10cSShri Abhyankar 	break;
234e2aaf10cSShri Abhyankar       }
2357b6afd5bSHong Zhang     }
2365f2c45f1SShri Abhyankar     network->header[i].ndata = 0;
2375f2c45f1SShri Abhyankar     ndata = network->header[i].ndata;
2385f2c45f1SShri Abhyankar     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
2395f2c45f1SShri Abhyankar     network->header[i].offset[ndata] = 0;
2405f2c45f1SShri Abhyankar   }
241e2aaf10cSShri Abhyankar 
242e2aaf10cSShri Abhyankar   for(i=network->vStart; i < network->vEnd; i++) {
243e2aaf10cSShri Abhyankar     network->header[i].index = i - network->vStart;
244e2aaf10cSShri Abhyankar     for(j=0; j < network->nsubnet; j++) {
245e2aaf10cSShri Abhyankar       if((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) {
246e2aaf10cSShri Abhyankar 	network->header[i].subnetid = j;
2472727e31bSShri Abhyankar 	network->subnet[j].vertices[network->subnet[j].nvtx++] = i;
248e2aaf10cSShri Abhyankar 	break;
249e2aaf10cSShri Abhyankar       }
250e2aaf10cSShri Abhyankar     }
251e2aaf10cSShri Abhyankar     network->header[i].ndata = 0;
252e2aaf10cSShri Abhyankar     ndata = network->header[i].ndata;
253e2aaf10cSShri Abhyankar     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
254e2aaf10cSShri Abhyankar     network->header[i].offset[ndata] = 0;
255e2aaf10cSShri Abhyankar   }
256e2aaf10cSShri Abhyankar 
257854ce69bSBarry Smith   ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr);
2586500d4abSHong Zhang   PetscFunctionReturn(0);
2596500d4abSHong Zhang }
260e2aaf10cSShri Abhyankar 
2616500d4abSHong Zhang PetscErrorCode DMNetworkLayoutSetUpCoupled(DM dm)
2626500d4abSHong Zhang {
2636500d4abSHong Zhang   PetscErrorCode ierr;
2646500d4abSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
2656500d4abSHong Zhang   PetscInt       dim = 1; /* One dimensional network */
266991cf414SHong Zhang   PetscInt       numCorners=2,spacedim=2;
2676500d4abSHong Zhang   double         *vertexcoords=NULL;
2687765340cSHong Zhang   PetscInt       i,j,ndata,ctr=0,nsubnet;
269991cf414SHong Zhang   PetscInt       *edgelist_couple=NULL,k,netid,vid;
2706500d4abSHong Zhang 
2716500d4abSHong Zhang   PetscFunctionBegin;
2726500d4abSHong Zhang   if (network->nVertices) {
2736500d4abSHong Zhang     ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr);
2746500d4abSHong Zhang   }
2756500d4abSHong Zhang 
2766500d4abSHong Zhang   /* Create the edgelist for the network by concatenating edgelists of the subnetworks */
2777765340cSHong Zhang   nsubnet = network->nsubnet - network->ncsubnet;
2786500d4abSHong Zhang   ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr);
2797765340cSHong Zhang   for (i=0; i < nsubnet; i++) {
2806500d4abSHong Zhang     for (j = 0; j < network->subnet[i].nedge; j++) {
2816500d4abSHong Zhang       network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
2826500d4abSHong Zhang       network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
2836500d4abSHong Zhang       ctr++;
2846500d4abSHong Zhang     }
2856500d4abSHong Zhang   }
2867765340cSHong Zhang 
2877765340cSHong Zhang   i       = nsubnet; /* netid of coupling subnet */
2887765340cSHong Zhang   nsubnet = network->nsubnet;
2897765340cSHong Zhang   while (i < nsubnet) {
290991cf414SHong Zhang     edgelist_couple = network->subnet[i].edgelist;
2916500d4abSHong Zhang     k = 0;
2926500d4abSHong Zhang     for (j = 0; j < network->subnet[i].nedge; j++) {
2936500d4abSHong Zhang       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
2946500d4abSHong Zhang       network->edges[2*ctr] = network->subnet[netid].vStart + vid; k += 2;
2956500d4abSHong Zhang 
2966500d4abSHong Zhang       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
297991cf414SHong Zhang       network->edges[2*ctr+1] = network->subnet[netid].vStart + vid; k+=2;
2986500d4abSHong Zhang       ctr++;
2996500d4abSHong Zhang     }
3007765340cSHong Zhang     i++;
3017765340cSHong Zhang   }
3026500d4abSHong Zhang 
30334bd1bf5SHong Zhang #if 0
3046500d4abSHong Zhang   for(i=0; i < network->nEdges; i++) {
3056500d4abSHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr);
3066500d4abSHong Zhang     printf("\n");
3076500d4abSHong Zhang   }
3086500d4abSHong Zhang #endif
3096500d4abSHong Zhang 
3106500d4abSHong Zhang   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
3116500d4abSHong Zhang   if (network->nVertices) {
3126500d4abSHong Zhang     ierr = PetscFree(vertexcoords);CHKERRQ(ierr);
3136500d4abSHong Zhang   }
3146500d4abSHong Zhang   ierr = PetscFree(network->edges);CHKERRQ(ierr);
3156500d4abSHong Zhang 
3166500d4abSHong Zhang   ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr);
3176500d4abSHong Zhang   ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr);
3186500d4abSHong Zhang   ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr);
3196500d4abSHong Zhang 
3206500d4abSHong Zhang   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr);
3216500d4abSHong Zhang   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr);
3226500d4abSHong Zhang   ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr);
3236500d4abSHong Zhang   ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr);
3246500d4abSHong Zhang 
3256500d4abSHong Zhang   /* Create vertices and edges array for the subnetworks */
3266500d4abSHong Zhang   for (j=0; j < network->nsubnet; j++) {
3276500d4abSHong Zhang     ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr);
3286500d4abSHong Zhang     ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr);
3296500d4abSHong Zhang     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
3306500d4abSHong Zhang        These get updated when the vertices and edges are added. */
3316500d4abSHong Zhang     network->subnet[j].nvtx = network->subnet[j].nedge = 0;
3326500d4abSHong Zhang   }
3336500d4abSHong Zhang 
3346500d4abSHong Zhang   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
3356500d4abSHong Zhang   ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr);
3366500d4abSHong Zhang   for (i=network->eStart; i < network->eEnd; i++) {
3376500d4abSHong Zhang     network->header[i].index = i;   /* Global edge number */
3386500d4abSHong Zhang     for (j=0; j < network->nsubnet; j++) {
3396500d4abSHong Zhang       if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) {
3406500d4abSHong Zhang 	network->header[i].subnetid = j; /* Subnetwork id */
3416500d4abSHong Zhang 	network->subnet[j].edges[network->subnet[j].nedge++] = i;
3426500d4abSHong Zhang 	break;
3436500d4abSHong Zhang       }
3446500d4abSHong Zhang     }
3456500d4abSHong Zhang     network->header[i].ndata = 0;
3466500d4abSHong Zhang     ndata = network->header[i].ndata;
3476500d4abSHong Zhang     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
3486500d4abSHong Zhang     network->header[i].offset[ndata] = 0;
3496500d4abSHong Zhang   }
3506500d4abSHong Zhang 
3516500d4abSHong Zhang   for(i=network->vStart; i < network->vEnd; i++) {
3526500d4abSHong Zhang     network->header[i].index = i - network->vStart;
3536500d4abSHong Zhang     for (j=0; j < network->nsubnet; j++) {
3546500d4abSHong Zhang       if ((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) {
3556500d4abSHong Zhang 	network->header[i].subnetid = j;
3566500d4abSHong Zhang 	network->subnet[j].vertices[network->subnet[j].nvtx++] = i;
3576500d4abSHong Zhang 	break;
3586500d4abSHong Zhang       }
3596500d4abSHong Zhang     }
3606500d4abSHong Zhang     network->header[i].ndata = 0;
3616500d4abSHong Zhang     ndata = network->header[i].ndata;
3626500d4abSHong Zhang     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
3636500d4abSHong Zhang     network->header[i].offset[ndata] = 0;
3646500d4abSHong Zhang   }
3656500d4abSHong Zhang 
3666500d4abSHong Zhang   ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr);
3675f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3685f2c45f1SShri Abhyankar }
3695f2c45f1SShri Abhyankar 
37094ef8ddeSSatish Balay /*@C
3712727e31bSShri Abhyankar   DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork
3722727e31bSShri Abhyankar 
3732727e31bSShri Abhyankar   Input Parameters
3742727e31bSShri Abhyankar + dm   - the number object
3752727e31bSShri Abhyankar - id   - the ID (integer) of the subnetwork
3762727e31bSShri Abhyankar 
3772727e31bSShri Abhyankar   Output Parameters
3782727e31bSShri Abhyankar + nv    - number of vertices (local)
3792727e31bSShri Abhyankar . ne    - number of edges (local)
3802727e31bSShri Abhyankar . vtx   - local vertices for this subnetwork
3812727e31bSShri Abhyankar . edge  - local edges for this subnetwork
3822727e31bSShri Abhyankar 
3832727e31bSShri Abhyankar   Notes:
3842727e31bSShri Abhyankar   Cannot call this routine before DMNetworkLayoutSetup()
3852727e31bSShri Abhyankar 
3862727e31bSShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
3872727e31bSShri Abhyankar @*/
3882727e31bSShri Abhyankar PetscErrorCode DMNetworkGetSubnetworkInfo(DM netdm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge)
3892727e31bSShri Abhyankar {
3902727e31bSShri Abhyankar   DM_Network     *network = (DM_Network*) netdm->data;
3912727e31bSShri Abhyankar 
3922727e31bSShri Abhyankar   PetscFunctionBegin;
3932727e31bSShri Abhyankar   *nv = network->subnet[id].nvtx;
3942727e31bSShri Abhyankar   *ne = network->subnet[id].nedge;
3952727e31bSShri Abhyankar   *vtx = network->subnet[id].vertices;
3962727e31bSShri Abhyankar   *edge = network->subnet[id].edges;
3972727e31bSShri Abhyankar   PetscFunctionReturn(0);
3982727e31bSShri Abhyankar }
3992727e31bSShri Abhyankar 
4002727e31bSShri Abhyankar /*@C
4015f2c45f1SShri Abhyankar   DMNetworkRegisterComponent - Registers the network component
4025f2c45f1SShri Abhyankar 
4035f2c45f1SShri Abhyankar   Logically collective on DM
4045f2c45f1SShri Abhyankar 
4055f2c45f1SShri Abhyankar   Input Parameters
4065f2c45f1SShri Abhyankar + dm   - the network object
4075f2c45f1SShri Abhyankar . name - the component name
4085f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data
4095f2c45f1SShri Abhyankar 
4105f2c45f1SShri Abhyankar    Output Parameters
4115f2c45f1SShri Abhyankar .   key - an integer key that defines the component
4125f2c45f1SShri Abhyankar 
4135f2c45f1SShri Abhyankar    Notes
4145f2c45f1SShri Abhyankar    This routine should be called by all processors before calling DMNetworkLayoutSetup().
4155f2c45f1SShri Abhyankar 
4165f2c45f1SShri Abhyankar    Level: intermediate
4175f2c45f1SShri Abhyankar 
4185f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
4195f2c45f1SShri Abhyankar @*/
4205f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key)
4215f2c45f1SShri Abhyankar {
4225f2c45f1SShri Abhyankar   PetscErrorCode        ierr;
4235f2c45f1SShri Abhyankar   DM_Network            *network = (DM_Network*) dm->data;
4245f2c45f1SShri Abhyankar   DMNetworkComponent    *component=&network->component[network->ncomponent];
4255f2c45f1SShri Abhyankar   PetscBool             flg=PETSC_FALSE;
4265f2c45f1SShri Abhyankar   PetscInt              i;
4275f2c45f1SShri Abhyankar 
4285f2c45f1SShri Abhyankar   PetscFunctionBegin;
4295f2c45f1SShri Abhyankar   for (i=0; i < network->ncomponent; i++) {
4305f2c45f1SShri Abhyankar     ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr);
4315f2c45f1SShri Abhyankar     if (flg) {
4325f2c45f1SShri Abhyankar       *key = i;
4335f2c45f1SShri Abhyankar       PetscFunctionReturn(0);
4345f2c45f1SShri Abhyankar     }
4356d64e262SShri Abhyankar   }
4366d64e262SShri Abhyankar   if(network->ncomponent == MAX_COMPONENTS) {
4376d64e262SShri Abhyankar     SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components registered exceeds the max %D",MAX_COMPONENTS);
4385f2c45f1SShri Abhyankar   }
4395f2c45f1SShri Abhyankar 
4405f2c45f1SShri Abhyankar   ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr);
4415f2c45f1SShri Abhyankar   component->size = size/sizeof(DMNetworkComponentGenericDataType);
4425f2c45f1SShri Abhyankar   *key = network->ncomponent;
4435f2c45f1SShri Abhyankar   network->ncomponent++;
4445f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4455f2c45f1SShri Abhyankar }
4465f2c45f1SShri Abhyankar 
4475f2c45f1SShri Abhyankar /*@
4485f2c45f1SShri Abhyankar   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.
4495f2c45f1SShri Abhyankar 
4505f2c45f1SShri Abhyankar   Not Collective
4515f2c45f1SShri Abhyankar 
4525f2c45f1SShri Abhyankar   Input Parameters:
4535f2c45f1SShri Abhyankar + dm - The DMNetwork object
4545f2c45f1SShri Abhyankar 
4555f2c45f1SShri Abhyankar   Output Paramters:
4565f2c45f1SShri Abhyankar + vStart - The first vertex point
4575f2c45f1SShri Abhyankar - vEnd   - One beyond the last vertex point
4585f2c45f1SShri Abhyankar 
4595f2c45f1SShri Abhyankar   Level: intermediate
4605f2c45f1SShri Abhyankar 
4615f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange
4625f2c45f1SShri Abhyankar @*/
4635f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
4645f2c45f1SShri Abhyankar {
4655f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4665f2c45f1SShri Abhyankar 
4675f2c45f1SShri Abhyankar   PetscFunctionBegin;
4685f2c45f1SShri Abhyankar   if (vStart) *vStart = network->vStart;
4695f2c45f1SShri Abhyankar   if (vEnd) *vEnd = network->vEnd;
4705f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4715f2c45f1SShri Abhyankar }
4725f2c45f1SShri Abhyankar 
4735f2c45f1SShri Abhyankar /*@
4745f2c45f1SShri Abhyankar   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.
4755f2c45f1SShri Abhyankar 
4765f2c45f1SShri Abhyankar   Not Collective
4775f2c45f1SShri Abhyankar 
4785f2c45f1SShri Abhyankar   Input Parameters:
4795f2c45f1SShri Abhyankar + dm - The DMNetwork object
4805f2c45f1SShri Abhyankar 
4815f2c45f1SShri Abhyankar   Output Paramters:
4825f2c45f1SShri Abhyankar + eStart - The first edge point
4835f2c45f1SShri Abhyankar - eEnd   - One beyond the last edge point
4845f2c45f1SShri Abhyankar 
4855f2c45f1SShri Abhyankar   Level: intermediate
4865f2c45f1SShri Abhyankar 
4875f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange
4885f2c45f1SShri Abhyankar @*/
4895f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
4905f2c45f1SShri Abhyankar {
4915f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4925f2c45f1SShri Abhyankar 
4935f2c45f1SShri Abhyankar   PetscFunctionBegin;
4945f2c45f1SShri Abhyankar   if (eStart) *eStart = network->eStart;
4955f2c45f1SShri Abhyankar   if (eEnd) *eEnd = network->eEnd;
4965f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4975f2c45f1SShri Abhyankar }
4985f2c45f1SShri Abhyankar 
4997b6afd5bSHong Zhang /*@
500e85e6aecSHong Zhang   DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge.
5017b6afd5bSHong Zhang 
5027b6afd5bSHong Zhang   Not Collective
5037b6afd5bSHong Zhang 
5047b6afd5bSHong Zhang   Input Parameters:
5057b6afd5bSHong Zhang + dm - DMNetwork object
506e85e6aecSHong Zhang - p  - edge point
5077b6afd5bSHong Zhang 
5087b6afd5bSHong Zhang   Output Paramters:
509e85e6aecSHong Zhang . index - user global numbering for the edge
5107b6afd5bSHong Zhang 
5117b6afd5bSHong Zhang   Level: intermediate
5127b6afd5bSHong Zhang 
513e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalVertexIndex
5147b6afd5bSHong Zhang @*/
515e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
5167b6afd5bSHong Zhang {
5177b6afd5bSHong Zhang   PetscErrorCode    ierr;
5187b6afd5bSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
5197b6afd5bSHong Zhang   PetscInt          offsetp;
5207b6afd5bSHong Zhang   DMNetworkComponentHeader header;
5217b6afd5bSHong Zhang 
5227b6afd5bSHong Zhang   PetscFunctionBegin;
5237b6afd5bSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
5247b6afd5bSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
525e85e6aecSHong Zhang   *index = header->index;
5267b6afd5bSHong Zhang   PetscFunctionReturn(0);
5277b6afd5bSHong Zhang }
5287b6afd5bSHong Zhang 
5295f2c45f1SShri Abhyankar /*@
530e85e6aecSHong Zhang   DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex.
531e85e6aecSHong Zhang 
532e85e6aecSHong Zhang   Not Collective
533e85e6aecSHong Zhang 
534e85e6aecSHong Zhang   Input Parameters:
535e85e6aecSHong Zhang + dm - DMNetwork object
536e85e6aecSHong Zhang - p  - vertex point
537e85e6aecSHong Zhang 
538e85e6aecSHong Zhang   Output Paramters:
539e85e6aecSHong Zhang . index - user global numbering for the vertex
540e85e6aecSHong Zhang 
541e85e6aecSHong Zhang   Level: intermediate
542e85e6aecSHong Zhang 
543e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalEdgeIndex
544e85e6aecSHong Zhang @*/
545e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
546e85e6aecSHong Zhang {
547e85e6aecSHong Zhang   PetscErrorCode    ierr;
548e85e6aecSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
549e85e6aecSHong Zhang   PetscInt          offsetp;
550e85e6aecSHong Zhang   DMNetworkComponentHeader header;
551e85e6aecSHong Zhang 
552e85e6aecSHong Zhang   PetscFunctionBegin;
553e85e6aecSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
554e85e6aecSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
555e85e6aecSHong Zhang   *index = header->index;
556e85e6aecSHong Zhang   PetscFunctionReturn(0);
557e85e6aecSHong Zhang }
558e85e6aecSHong Zhang 
559c3b11c7cSShri Abhyankar /*
560c3b11c7cSShri Abhyankar   DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the
561c3b11c7cSShri Abhyankar                                     component value from the component data array
562c3b11c7cSShri Abhyankar 
563c3b11c7cSShri Abhyankar   Not Collective
564c3b11c7cSShri Abhyankar 
565c3b11c7cSShri Abhyankar   Input Parameters:
566c3b11c7cSShri Abhyankar + dm      - The DMNetwork object
567c3b11c7cSShri Abhyankar . p       - vertex/edge point
568c3b11c7cSShri Abhyankar - compnum - component number
569c3b11c7cSShri Abhyankar 
570c3b11c7cSShri Abhyankar   Output Parameters:
571c3b11c7cSShri Abhyankar + compkey - the key obtained when registering the component
572c3b11c7cSShri Abhyankar - offset  - offset into the component data array associated with the vertex/edge point
573c3b11c7cSShri Abhyankar 
574c3b11c7cSShri Abhyankar   Notes:
575c3b11c7cSShri Abhyankar   Typical usage:
576c3b11c7cSShri Abhyankar 
577c3b11c7cSShri Abhyankar   DMNetworkGetComponentDataArray(dm, &arr);
578c3b11c7cSShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
579c3b11c7cSShri Abhyankar   Loop over vertices or edges
580c3b11c7cSShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
581c3b11c7cSShri Abhyankar     Loop over numcomps
582c3b11c7cSShri Abhyankar       DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset);
583c3b11c7cSShri Abhyankar       compdata = (UserCompDataType)(arr+offset);
584c3b11c7cSShri Abhyankar 
585c3b11c7cSShri Abhyankar   Level: intermediate
586c3b11c7cSShri Abhyankar 
587c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
588c3b11c7cSShri Abhyankar */
589c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
590c3b11c7cSShri Abhyankar {
591c3b11c7cSShri Abhyankar   PetscErrorCode           ierr;
592c3b11c7cSShri Abhyankar   PetscInt                 offsetp;
593c3b11c7cSShri Abhyankar   DMNetworkComponentHeader header;
594c3b11c7cSShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
595c3b11c7cSShri Abhyankar 
596c3b11c7cSShri Abhyankar   PetscFunctionBegin;
597c3b11c7cSShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
598c3b11c7cSShri Abhyankar   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
599c3b11c7cSShri Abhyankar   if (compkey) *compkey = header->key[compnum];
600c3b11c7cSShri Abhyankar   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
601c3b11c7cSShri Abhyankar   PetscFunctionReturn(0);
602c3b11c7cSShri Abhyankar }
603c3b11c7cSShri Abhyankar 
604c3b11c7cSShri Abhyankar /*@
605c3b11c7cSShri Abhyankar   DMNetworkGetComponent - Returns the network component and its key
606c3b11c7cSShri Abhyankar 
607c3b11c7cSShri Abhyankar   Not Collective
608c3b11c7cSShri Abhyankar 
609c3b11c7cSShri Abhyankar   Input Parameters
610c3b11c7cSShri Abhyankar + dm - DMNetwork object
611c3b11c7cSShri Abhyankar . p  - edge or vertex point
612c3b11c7cSShri Abhyankar - compnum - component number
613c3b11c7cSShri Abhyankar 
614c3b11c7cSShri Abhyankar   Output Parameters:
615c3b11c7cSShri Abhyankar + compkey - the key set for this computing during registration
616c3b11c7cSShri Abhyankar - component - the component data
617c3b11c7cSShri Abhyankar 
618c3b11c7cSShri Abhyankar   Notes:
619c3b11c7cSShri Abhyankar   Typical usage:
620c3b11c7cSShri Abhyankar 
621c3b11c7cSShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
622c3b11c7cSShri Abhyankar   Loop over vertices or edges
623c3b11c7cSShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
624c3b11c7cSShri Abhyankar     Loop over numcomps
625c3b11c7cSShri Abhyankar       DMNetworkGetComponent(dm,v,compnum,&key,&component);
626c3b11c7cSShri Abhyankar 
627c3b11c7cSShri Abhyankar   Level: intermediate
628c3b11c7cSShri Abhyankar 
629c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetVariableOffset
630c3b11c7cSShri Abhyankar @*/
631c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *key, void **component)
632c3b11c7cSShri Abhyankar {
633c3b11c7cSShri Abhyankar   PetscErrorCode ierr;
634c3b11c7cSShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
635c3b11c7cSShri Abhyankar   PetscInt       offsetd;
636c3b11c7cSShri Abhyankar 
637c3b11c7cSShri Abhyankar   PetscFunctionBegin;
638c3b11c7cSShri Abhyankar 
639c3b11c7cSShri Abhyankar   ierr = DMNetworkGetComponentKeyOffset(dm,p,compnum,key,&offsetd);CHKERRQ(ierr);
640c3b11c7cSShri Abhyankar   *component = network->componentdataarray+offsetd;
641c3b11c7cSShri Abhyankar 
642c3b11c7cSShri Abhyankar   PetscFunctionReturn(0);
643c3b11c7cSShri Abhyankar }
644c3b11c7cSShri Abhyankar 
645e85e6aecSHong Zhang /*@
646325661f6SSatish Balay   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)
6475f2c45f1SShri Abhyankar 
6485f2c45f1SShri Abhyankar   Not Collective
6495f2c45f1SShri Abhyankar 
6505f2c45f1SShri Abhyankar   Input Parameters:
6515f2c45f1SShri Abhyankar + dm           - The DMNetwork object
6525f2c45f1SShri Abhyankar . p            - vertex/edge point
6535f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component
6545f2c45f1SShri Abhyankar - compvalue    - pointer to the data structure for the component
6555f2c45f1SShri Abhyankar 
6565f2c45f1SShri Abhyankar   Level: intermediate
6575f2c45f1SShri Abhyankar 
6585f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
6595f2c45f1SShri Abhyankar @*/
6605f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
6615f2c45f1SShri Abhyankar {
6625f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
66343a39a44SBarry Smith   DMNetworkComponent       *component = &network->component[componentkey];
6645f2c45f1SShri Abhyankar   DMNetworkComponentHeader header = &network->header[p];
6655f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue = &network->cvalue[p];
6665f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
6675f2c45f1SShri Abhyankar 
6685f2c45f1SShri Abhyankar   PetscFunctionBegin;
669fa58f0a9SHong 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);
670fa58f0a9SHong Zhang 
67143a39a44SBarry Smith   header->size[header->ndata] = component->size;
67243a39a44SBarry Smith   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
6735f2c45f1SShri Abhyankar   header->key[header->ndata] = componentkey;
6745f2c45f1SShri Abhyankar   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
6755f2c45f1SShri Abhyankar 
6765f2c45f1SShri Abhyankar   cvalue->data[header->ndata] = (void*)compvalue;
6775f2c45f1SShri Abhyankar   header->ndata++;
6785f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6795f2c45f1SShri Abhyankar }
6805f2c45f1SShri Abhyankar 
6815f2c45f1SShri Abhyankar /*@
6825f2c45f1SShri Abhyankar   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
6835f2c45f1SShri Abhyankar 
6845f2c45f1SShri Abhyankar   Not Collective
6855f2c45f1SShri Abhyankar 
6865f2c45f1SShri Abhyankar   Input Parameters:
6875f2c45f1SShri Abhyankar + dm - The DMNetwork object
6885f2c45f1SShri Abhyankar . p  - vertex/edge point
6895f2c45f1SShri Abhyankar 
6905f2c45f1SShri Abhyankar   Output Parameters:
6915f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge
6925f2c45f1SShri Abhyankar 
6935f2c45f1SShri Abhyankar   Level: intermediate
6945f2c45f1SShri Abhyankar 
6955f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
6965f2c45f1SShri Abhyankar @*/
6975f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
6985f2c45f1SShri Abhyankar {
6995f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7005f2c45f1SShri Abhyankar   PetscInt       offset;
7015f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7025f2c45f1SShri Abhyankar 
7035f2c45f1SShri Abhyankar   PetscFunctionBegin;
7045f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
7055f2c45f1SShri Abhyankar   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
7065f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7075f2c45f1SShri Abhyankar }
7085f2c45f1SShri Abhyankar 
7095f2c45f1SShri Abhyankar /*@
7105f2c45f1SShri Abhyankar   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
7115f2c45f1SShri Abhyankar 
7125f2c45f1SShri Abhyankar   Not Collective
7135f2c45f1SShri Abhyankar 
7145f2c45f1SShri Abhyankar   Input Parameters:
7155f2c45f1SShri Abhyankar + dm     - The DMNetwork object
7165f2c45f1SShri Abhyankar - p      - the edge/vertex point
7175f2c45f1SShri Abhyankar 
7185f2c45f1SShri Abhyankar   Output Parameters:
7195f2c45f1SShri Abhyankar . offset - the offset
7205f2c45f1SShri Abhyankar 
7215f2c45f1SShri Abhyankar   Level: intermediate
7225f2c45f1SShri Abhyankar 
7235f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
7245f2c45f1SShri Abhyankar @*/
7255f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
7265f2c45f1SShri Abhyankar {
7275f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7285f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7295f2c45f1SShri Abhyankar 
7305f2c45f1SShri Abhyankar   PetscFunctionBegin;
7315f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr);
7325f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7335f2c45f1SShri Abhyankar }
7345f2c45f1SShri Abhyankar 
7355f2c45f1SShri Abhyankar /*@
7365f2c45f1SShri Abhyankar   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
7375f2c45f1SShri Abhyankar 
7385f2c45f1SShri Abhyankar   Not Collective
7395f2c45f1SShri Abhyankar 
7405f2c45f1SShri Abhyankar   Input Parameters:
7415f2c45f1SShri Abhyankar + dm      - The DMNetwork object
7425f2c45f1SShri Abhyankar - p       - the edge/vertex point
7435f2c45f1SShri Abhyankar 
7445f2c45f1SShri Abhyankar   Output Parameters:
7455f2c45f1SShri Abhyankar . offsetg - the offset
7465f2c45f1SShri Abhyankar 
7475f2c45f1SShri Abhyankar   Level: intermediate
7485f2c45f1SShri Abhyankar 
7495f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
7505f2c45f1SShri Abhyankar @*/
7515f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
7525f2c45f1SShri Abhyankar {
7535f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7545f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7555f2c45f1SShri Abhyankar 
7565f2c45f1SShri Abhyankar   PetscFunctionBegin;
7575f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr);
7586fefedf4SHong Zhang   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */
7595f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7605f2c45f1SShri Abhyankar }
7615f2c45f1SShri Abhyankar 
76224121865SAdrian Maldonado /*@
76324121865SAdrian Maldonado   DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector.
76424121865SAdrian Maldonado 
76524121865SAdrian Maldonado   Not Collective
76624121865SAdrian Maldonado 
76724121865SAdrian Maldonado   Input Parameters:
76824121865SAdrian Maldonado + dm     - The DMNetwork object
76924121865SAdrian Maldonado - p      - the edge point
77024121865SAdrian Maldonado 
77124121865SAdrian Maldonado   Output Parameters:
77224121865SAdrian Maldonado . offset - the offset
77324121865SAdrian Maldonado 
77424121865SAdrian Maldonado   Level: intermediate
77524121865SAdrian Maldonado 
77624121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
77724121865SAdrian Maldonado @*/
77824121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
77924121865SAdrian Maldonado {
78024121865SAdrian Maldonado   PetscErrorCode ierr;
78124121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
78224121865SAdrian Maldonado 
78324121865SAdrian Maldonado   PetscFunctionBegin;
78424121865SAdrian Maldonado 
78524121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr);
78624121865SAdrian Maldonado   PetscFunctionReturn(0);
78724121865SAdrian Maldonado }
78824121865SAdrian Maldonado 
78924121865SAdrian Maldonado /*@
79024121865SAdrian Maldonado   DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector.
79124121865SAdrian Maldonado 
79224121865SAdrian Maldonado   Not Collective
79324121865SAdrian Maldonado 
79424121865SAdrian Maldonado   Input Parameters:
79524121865SAdrian Maldonado + dm     - The DMNetwork object
79624121865SAdrian Maldonado - p      - the vertex point
79724121865SAdrian Maldonado 
79824121865SAdrian Maldonado   Output Parameters:
79924121865SAdrian Maldonado . offset - the offset
80024121865SAdrian Maldonado 
80124121865SAdrian Maldonado   Level: intermediate
80224121865SAdrian Maldonado 
80324121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
80424121865SAdrian Maldonado @*/
80524121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
80624121865SAdrian Maldonado {
80724121865SAdrian Maldonado   PetscErrorCode ierr;
80824121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
80924121865SAdrian Maldonado 
81024121865SAdrian Maldonado   PetscFunctionBegin;
81124121865SAdrian Maldonado 
81224121865SAdrian Maldonado   p -= network->vStart;
81324121865SAdrian Maldonado 
81424121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr);
81524121865SAdrian Maldonado   PetscFunctionReturn(0);
81624121865SAdrian Maldonado }
8175f2c45f1SShri Abhyankar /*@
8185f2c45f1SShri Abhyankar   DMNetworkAddNumVariables - Add number of variables associated with a given point.
8195f2c45f1SShri Abhyankar 
8205f2c45f1SShri Abhyankar   Not Collective
8215f2c45f1SShri Abhyankar 
8225f2c45f1SShri Abhyankar   Input Parameters:
8235f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
8245f2c45f1SShri Abhyankar . p    - the vertex/edge point
8255f2c45f1SShri Abhyankar - nvar - number of additional variables
8265f2c45f1SShri Abhyankar 
8275f2c45f1SShri Abhyankar   Level: intermediate
8285f2c45f1SShri Abhyankar 
8295f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables
8305f2c45f1SShri Abhyankar @*/
8315f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
8325f2c45f1SShri Abhyankar {
8335f2c45f1SShri Abhyankar   PetscErrorCode ierr;
8345f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
8355f2c45f1SShri Abhyankar 
8365f2c45f1SShri Abhyankar   PetscFunctionBegin;
8375f2c45f1SShri Abhyankar   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
8385f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
8395f2c45f1SShri Abhyankar }
8405f2c45f1SShri Abhyankar 
84127f51fceSHong Zhang /*@
84227f51fceSHong Zhang   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
84327f51fceSHong Zhang 
84427f51fceSHong Zhang   Not Collective
84527f51fceSHong Zhang 
84627f51fceSHong Zhang   Input Parameters:
84727f51fceSHong Zhang + dm   - The DMNetworkObject
84827f51fceSHong Zhang - p    - the vertex/edge point
84927f51fceSHong Zhang 
85027f51fceSHong Zhang   Output Parameters:
85127f51fceSHong Zhang . nvar - number of variables
85227f51fceSHong Zhang 
85327f51fceSHong Zhang   Level: intermediate
85427f51fceSHong Zhang 
85527f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
85627f51fceSHong Zhang @*/
85727f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
85827f51fceSHong Zhang {
85927f51fceSHong Zhang   PetscErrorCode ierr;
86027f51fceSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
86127f51fceSHong Zhang 
86227f51fceSHong Zhang   PetscFunctionBegin;
86327f51fceSHong Zhang   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
86427f51fceSHong Zhang   PetscFunctionReturn(0);
86527f51fceSHong Zhang }
86627f51fceSHong Zhang 
8675f2c45f1SShri Abhyankar /*@
8685f2c45f1SShri Abhyankar   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
8695f2c45f1SShri Abhyankar 
8705f2c45f1SShri Abhyankar   Not Collective
8715f2c45f1SShri Abhyankar 
8725f2c45f1SShri Abhyankar   Input Parameters:
8735f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
8745f2c45f1SShri Abhyankar . p    - the vertex/edge point
8755f2c45f1SShri Abhyankar - nvar - number of variables
8765f2c45f1SShri Abhyankar 
8775f2c45f1SShri Abhyankar   Level: intermediate
8785f2c45f1SShri Abhyankar 
8795f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables
8805f2c45f1SShri Abhyankar @*/
8815f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
8825f2c45f1SShri Abhyankar {
8835f2c45f1SShri Abhyankar   PetscErrorCode ierr;
8845f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
8855f2c45f1SShri Abhyankar 
8865f2c45f1SShri Abhyankar   PetscFunctionBegin;
8875f2c45f1SShri Abhyankar   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
8885f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
8895f2c45f1SShri Abhyankar }
8905f2c45f1SShri Abhyankar 
8915f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This
8925f2c45f1SShri Abhyankar    function is called during DMSetUp() */
8935f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm)
8945f2c45f1SShri Abhyankar {
8955f2c45f1SShri Abhyankar   PetscErrorCode              ierr;
8965f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
8975f2c45f1SShri Abhyankar   PetscInt                    arr_size;
8985f2c45f1SShri Abhyankar   PetscInt                    p,offset,offsetp;
8995f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
9005f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue;
9015f2c45f1SShri Abhyankar   DMNetworkComponentGenericDataType      *componentdataarray;
9025f2c45f1SShri Abhyankar   PetscInt ncomp, i;
9035f2c45f1SShri Abhyankar 
9045f2c45f1SShri Abhyankar   PetscFunctionBegin;
9055f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
9065f2c45f1SShri Abhyankar   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
90775b160a0SShri Abhyankar   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
9085f2c45f1SShri Abhyankar   componentdataarray = network->componentdataarray;
9095f2c45f1SShri Abhyankar   for (p = network->pStart; p < network->pEnd; p++) {
9105f2c45f1SShri Abhyankar     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
9115f2c45f1SShri Abhyankar     /* Copy header */
9125f2c45f1SShri Abhyankar     header = &network->header[p];
913302440fdSBarry Smith     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
9145f2c45f1SShri Abhyankar     /* Copy data */
9155f2c45f1SShri Abhyankar     cvalue = &network->cvalue[p];
9165f2c45f1SShri Abhyankar     ncomp = header->ndata;
9175f2c45f1SShri Abhyankar     for (i = 0; i < ncomp; i++) {
9185f2c45f1SShri Abhyankar       offset = offsetp + network->dataheadersize + header->offset[i];
919302440fdSBarry Smith       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
9205f2c45f1SShri Abhyankar     }
9215f2c45f1SShri Abhyankar   }
9225f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9235f2c45f1SShri Abhyankar }
9245f2c45f1SShri Abhyankar 
9255f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */
9265f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm)
9275f2c45f1SShri Abhyankar {
9285f2c45f1SShri Abhyankar   PetscErrorCode ierr;
9295f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
9305f2c45f1SShri Abhyankar 
9315f2c45f1SShri Abhyankar   PetscFunctionBegin;
9325f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
9335f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9345f2c45f1SShri Abhyankar }
9355f2c45f1SShri Abhyankar 
9365f2c45f1SShri Abhyankar /*@C
9375f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray - Returns the component data array
9385f2c45f1SShri Abhyankar 
9395f2c45f1SShri Abhyankar   Not Collective
9405f2c45f1SShri Abhyankar 
9415f2c45f1SShri Abhyankar   Input Parameters:
9425f2c45f1SShri Abhyankar . dm - The DMNetwork Object
9435f2c45f1SShri Abhyankar 
9445f2c45f1SShri Abhyankar   Output Parameters:
9455f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components
9465f2c45f1SShri Abhyankar 
9475f2c45f1SShri Abhyankar   Level: intermediate
9485f2c45f1SShri Abhyankar 
949a730d845SHong Zhang .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents
9505f2c45f1SShri Abhyankar @*/
9515f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
9525f2c45f1SShri Abhyankar {
9535f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
9545f2c45f1SShri Abhyankar 
9555f2c45f1SShri Abhyankar   PetscFunctionBegin;
9565f2c45f1SShri Abhyankar   *componentdataarray = network->componentdataarray;
9575f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9585f2c45f1SShri Abhyankar }
9595f2c45f1SShri Abhyankar 
96024121865SAdrian Maldonado /* Get a subsection from a range of points */
96124121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection)
96224121865SAdrian Maldonado {
96324121865SAdrian Maldonado   PetscErrorCode ierr;
96424121865SAdrian Maldonado   PetscInt       i, nvar;
96524121865SAdrian Maldonado 
96624121865SAdrian Maldonado   PetscFunctionBegin;
96724121865SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr);
96824121865SAdrian Maldonado   ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr);
96924121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
97024121865SAdrian Maldonado     ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr);
97124121865SAdrian Maldonado     ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr);
97224121865SAdrian Maldonado   }
97324121865SAdrian Maldonado 
97424121865SAdrian Maldonado   ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr);
97524121865SAdrian Maldonado   PetscFunctionReturn(0);
97624121865SAdrian Maldonado }
97724121865SAdrian Maldonado 
97824121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */
97924121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
98024121865SAdrian Maldonado {
98124121865SAdrian Maldonado   PetscErrorCode ierr;
98224121865SAdrian Maldonado   PetscInt       i, *subpoints;
98324121865SAdrian Maldonado 
98424121865SAdrian Maldonado   PetscFunctionBegin;
98524121865SAdrian Maldonado   /* Create index sets to map from "points" to "subpoints" */
98624121865SAdrian Maldonado   ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr);
98724121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
98824121865SAdrian Maldonado     subpoints[i - pstart] = i;
98924121865SAdrian Maldonado   }
990459726d8SSatish Balay   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr);
99124121865SAdrian Maldonado   ierr = PetscFree(subpoints);CHKERRQ(ierr);
99224121865SAdrian Maldonado   PetscFunctionReturn(0);
99324121865SAdrian Maldonado }
99424121865SAdrian Maldonado 
99524121865SAdrian Maldonado /*@
99624121865SAdrian Maldonado   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute.
99724121865SAdrian Maldonado 
99824121865SAdrian Maldonado   Collective
99924121865SAdrian Maldonado 
100024121865SAdrian Maldonado   Input Parameters:
100124121865SAdrian Maldonado . dm   - The DMNetworkObject
100224121865SAdrian Maldonado 
100324121865SAdrian Maldonado   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
100424121865SAdrian Maldonado 
100524121865SAdrian Maldonado   points = [0 1 2 3 4 5 6]
100624121865SAdrian Maldonado 
100724121865SAdrian Maldonado   where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]).
100824121865SAdrian Maldonado 
100924121865SAdrian Maldonado   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
101024121865SAdrian Maldonado 
101124121865SAdrian Maldonado   Level: intermediate
101224121865SAdrian Maldonado 
101324121865SAdrian Maldonado @*/
101424121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
101524121865SAdrian Maldonado {
101624121865SAdrian Maldonado   PetscErrorCode ierr;
101724121865SAdrian Maldonado   MPI_Comm       comm;
10189852e123SBarry Smith   PetscMPIInt    rank, size;
101924121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
102024121865SAdrian Maldonado 
1021eab1376dSHong Zhang   PetscFunctionBegin;
102224121865SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
102324121865SAdrian Maldonado   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
10249852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
102524121865SAdrian Maldonado 
102624121865SAdrian Maldonado   /* Create maps for vertices and edges */
102724121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr);
102824121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr);
102924121865SAdrian Maldonado 
103024121865SAdrian Maldonado   /* Create local sub-sections */
103124121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr);
103224121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr);
103324121865SAdrian Maldonado 
10349852e123SBarry Smith   if (size > 1) {
103524121865SAdrian Maldonado     ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr);
103624121865SAdrian Maldonado     ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr);
103724121865SAdrian Maldonado   ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr);
103824121865SAdrian Maldonado   ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr);
103924121865SAdrian Maldonado   } else {
104024121865SAdrian Maldonado   /* create structures for vertex */
104124121865SAdrian Maldonado   ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr);
104224121865SAdrian Maldonado   /* create structures for edge */
104324121865SAdrian Maldonado   ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr);
104424121865SAdrian Maldonado   }
104524121865SAdrian Maldonado 
104624121865SAdrian Maldonado 
104724121865SAdrian Maldonado   /* Add viewers */
104824121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr);
104924121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr);
105024121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr);
105124121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr);
105224121865SAdrian Maldonado 
105324121865SAdrian Maldonado   PetscFunctionReturn(0);
105424121865SAdrian Maldonado }
10557b6afd5bSHong Zhang 
10565f2c45f1SShri Abhyankar /*@
10575f2c45f1SShri Abhyankar   DMNetworkDistribute - Distributes the network and moves associated component data.
10585f2c45f1SShri Abhyankar 
10595f2c45f1SShri Abhyankar   Collective
10605f2c45f1SShri Abhyankar 
10615f2c45f1SShri Abhyankar   Input Parameter:
1062d3464fd4SAdrian Maldonado + DM - the DMNetwork object
10635f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default
10645f2c45f1SShri Abhyankar 
10655f2c45f1SShri Abhyankar   Notes:
10668b171c8eSHong Zhang   Distributes the network with <overlap>-overlapping partitioning of the edges.
10675f2c45f1SShri Abhyankar 
10685f2c45f1SShri Abhyankar   Level: intermediate
10695f2c45f1SShri Abhyankar 
10705f2c45f1SShri Abhyankar .seealso: DMNetworkCreate
10715f2c45f1SShri Abhyankar @*/
1072d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
10735f2c45f1SShri Abhyankar {
1074d3464fd4SAdrian Maldonado   MPI_Comm       comm;
10755f2c45f1SShri Abhyankar   PetscErrorCode ierr;
1076d3464fd4SAdrian Maldonado   PetscMPIInt    size;
1077d3464fd4SAdrian Maldonado   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
1078d3464fd4SAdrian Maldonado   DM_Network     *newDMnetwork;
10795f2c45f1SShri Abhyankar   PetscSF        pointsf;
10805f2c45f1SShri Abhyankar   DM             newDM;
108151ac5effSHong Zhang   PetscPartitioner part;
1082b9c6e19dSShri Abhyankar   PetscInt         j,e,v,offset;
1083b9c6e19dSShri Abhyankar   DMNetworkComponentHeader header;
10845f2c45f1SShri Abhyankar 
10855f2c45f1SShri Abhyankar   PetscFunctionBegin;
1086d3464fd4SAdrian Maldonado 
1087d3464fd4SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr);
1088d3464fd4SAdrian Maldonado   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1089d3464fd4SAdrian Maldonado   if (size == 1) PetscFunctionReturn(0);
1090d3464fd4SAdrian Maldonado 
1091d3464fd4SAdrian Maldonado   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr);
10925f2c45f1SShri Abhyankar   newDMnetwork = (DM_Network*)newDM->data;
10935f2c45f1SShri Abhyankar   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
109451ac5effSHong Zhang 
109551ac5effSHong Zhang   /* Enable runtime options for petscpartitioner */
109651ac5effSHong Zhang   ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr);
109751ac5effSHong Zhang   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
109851ac5effSHong Zhang 
10995f2c45f1SShri Abhyankar   /* Distribute plex dm and dof section */
110080cf41d5SMatthew G. Knepley   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
110151ac5effSHong Zhang 
11025f2c45f1SShri Abhyankar   /* Distribute dof section */
1103d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr);
11045f2c45f1SShri Abhyankar   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
1105d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr);
110651ac5effSHong Zhang 
11075f2c45f1SShri Abhyankar   /* Distribute data and associated section */
110831da1fc8SHong Zhang   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
110924121865SAdrian Maldonado 
11105f2c45f1SShri Abhyankar   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
11115f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
11125f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
11135f2c45f1SShri Abhyankar   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
11146fefedf4SHong Zhang   newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
11156fefedf4SHong Zhang   newDMnetwork->NVertices = oldDMnetwork->NVertices;
11165f2c45f1SShri Abhyankar   newDMnetwork->NEdges = oldDMnetwork->NEdges;
111724121865SAdrian Maldonado 
11185f2c45f1SShri Abhyankar   /* Set Dof section as the default section for dm */
11195f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
11205f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
11215f2c45f1SShri Abhyankar 
1122b9c6e19dSShri Abhyankar   /* Set up subnetwork info in the newDM */
1123b9c6e19dSShri Abhyankar   newDMnetwork->nsubnet = oldDMnetwork->nsubnet;
1124b9c6e19dSShri Abhyankar   ierr = PetscCalloc1(newDMnetwork->nsubnet,&newDMnetwork->subnet);CHKERRQ(ierr);
1125b9c6e19dSShri Abhyankar   /* Copy over the global number of vertices and edges in each subnetwork. Note that these are already
1126b9c6e19dSShri Abhyankar      calculated in DMNetworkLayoutSetUp()
1127b9c6e19dSShri Abhyankar   */
1128b9c6e19dSShri Abhyankar   for(j=0; j < newDMnetwork->nsubnet; j++) {
1129b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx;
1130b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge;
1131b9c6e19dSShri Abhyankar   }
1132b9c6e19dSShri Abhyankar 
1133b9c6e19dSShri Abhyankar   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1134b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
1135b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1136b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].nedge++;
1137b9c6e19dSShri Abhyankar   }
1138b9c6e19dSShri Abhyankar 
1139b9c6e19dSShri Abhyankar   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1140b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
1141b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1142b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].nvtx++;
1143b9c6e19dSShri Abhyankar   }
1144b9c6e19dSShri Abhyankar 
1145b9c6e19dSShri Abhyankar   /* Now create the vertices and edge arrays for the subnetworks */
1146b9c6e19dSShri Abhyankar   for(j=0; j < newDMnetwork->nsubnet; j++) {
1147b9c6e19dSShri Abhyankar     ierr = PetscCalloc1(newDMnetwork->subnet[j].nedge,&newDMnetwork->subnet[j].edges);CHKERRQ(ierr);
1148b9c6e19dSShri Abhyankar     ierr = PetscCalloc1(newDMnetwork->subnet[j].nvtx,&newDMnetwork->subnet[j].vertices);CHKERRQ(ierr);
1149b9c6e19dSShri Abhyankar     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
1150b9c6e19dSShri Abhyankar        These get updated when the vertices and edges are added. */
1151b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0;
1152b9c6e19dSShri Abhyankar   }
1153b9c6e19dSShri Abhyankar 
1154b9c6e19dSShri Abhyankar   /* Set the vertices and edges in each subnetwork */
1155b9c6e19dSShri Abhyankar   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1156b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
1157b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1158b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++]  = e;
1159b9c6e19dSShri Abhyankar   }
1160b9c6e19dSShri Abhyankar 
1161b9c6e19dSShri Abhyankar   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1162b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
1163b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1164b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++]  = v;
1165b9c6e19dSShri Abhyankar   }
1166b9c6e19dSShri Abhyankar 
116724121865SAdrian Maldonado   /* Destroy point SF */
116824121865SAdrian Maldonado   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
116924121865SAdrian Maldonado 
1170d3464fd4SAdrian Maldonado   ierr = DMDestroy(dm);CHKERRQ(ierr);
1171d3464fd4SAdrian Maldonado   *dm  = newDM;
11725f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
11735f2c45f1SShri Abhyankar }
11745f2c45f1SShri Abhyankar 
117524121865SAdrian Maldonado /*@C
117624121865SAdrian Maldonado   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering.
117724121865SAdrian Maldonado 
117824121865SAdrian Maldonado   Input Parameters:
117924121865SAdrian Maldonado + masterSF - the original SF structure
118024121865SAdrian Maldonado - map      - a ISLocalToGlobal mapping that contains the subset of points
118124121865SAdrian Maldonado 
118224121865SAdrian Maldonado   Output Parameters:
118324121865SAdrian Maldonado . subSF    - a subset of the masterSF for the desired subset.
118424121865SAdrian Maldonado */
118524121865SAdrian Maldonado 
118624121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) {
118724121865SAdrian Maldonado 
118824121865SAdrian Maldonado   PetscErrorCode        ierr;
118924121865SAdrian Maldonado   PetscInt              nroots, nleaves, *ilocal_sub;
119024121865SAdrian Maldonado   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
119124121865SAdrian Maldonado   PetscInt              *local_points, *remote_points;
119224121865SAdrian Maldonado   PetscSFNode           *iremote_sub;
119324121865SAdrian Maldonado   const PetscInt        *ilocal;
119424121865SAdrian Maldonado   const PetscSFNode     *iremote;
119524121865SAdrian Maldonado 
119624121865SAdrian Maldonado   PetscFunctionBegin;
119724121865SAdrian Maldonado   ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
119824121865SAdrian Maldonado 
119924121865SAdrian Maldonado   /* Look for leaves that pertain to the subset of points. Get the local ordering */
120024121865SAdrian Maldonado   ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr);
120124121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr);
120224121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
120324121865SAdrian Maldonado     if (ilocal_map[i] != -1) nleaves_sub += 1;
120424121865SAdrian Maldonado   }
120524121865SAdrian Maldonado   /* Re-number ilocal with subset numbering. Need information from roots */
120624121865SAdrian Maldonado   ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr);
120724121865SAdrian Maldonado   for (i = 0; i < nroots; i++) local_points[i] = i;
120824121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr);
120924121865SAdrian Maldonado   ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
121024121865SAdrian Maldonado   ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
121124121865SAdrian Maldonado   /* Fill up graph using local (that is, local to the subset) numbering. */
12124b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr);
12134b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr);
121424121865SAdrian Maldonado   nleaves_sub = 0;
121524121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
121624121865SAdrian Maldonado     if (ilocal_map[i] != -1) {
121724121865SAdrian Maldonado       ilocal_sub[nleaves_sub] = ilocal_map[i];
12184b70a8deSAdrian Maldonado       iremote_sub[nleaves_sub].rank = iremote[i].rank;
121924121865SAdrian Maldonado       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
122024121865SAdrian Maldonado       nleaves_sub += 1;
122124121865SAdrian Maldonado     }
122224121865SAdrian Maldonado   }
122324121865SAdrian Maldonado   ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr);
122424121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr);
122524121865SAdrian Maldonado 
122624121865SAdrian Maldonado   /* Create new subSF */
122724121865SAdrian Maldonado   ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr);
122824121865SAdrian Maldonado   ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr);
12294b70a8deSAdrian Maldonado   ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr);
123024121865SAdrian Maldonado   ierr = PetscFree(ilocal_map);CHKERRQ(ierr);
12314b70a8deSAdrian Maldonado   ierr = PetscFree(iremote_sub);CHKERRQ(ierr);
123224121865SAdrian Maldonado   PetscFunctionReturn(0);
123324121865SAdrian Maldonado }
123424121865SAdrian Maldonado 
12355f2c45f1SShri Abhyankar /*@C
12365f2c45f1SShri Abhyankar   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
12375f2c45f1SShri Abhyankar 
12385f2c45f1SShri Abhyankar   Not Collective
12395f2c45f1SShri Abhyankar 
12405f2c45f1SShri Abhyankar   Input Parameters:
12415f2c45f1SShri Abhyankar + dm - The DMNetwork object
12425f2c45f1SShri Abhyankar - p  - the vertex point
12435f2c45f1SShri Abhyankar 
12445f2c45f1SShri Abhyankar   Output Paramters:
12455f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point
12465f2c45f1SShri Abhyankar - edges  - List of edge points
12475f2c45f1SShri Abhyankar 
12485f2c45f1SShri Abhyankar   Level: intermediate
12495f2c45f1SShri Abhyankar 
12505f2c45f1SShri Abhyankar   Fortran Notes:
12515f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
12525f2c45f1SShri Abhyankar   include petsc.h90 in your code.
12535f2c45f1SShri Abhyankar 
1254d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices
12555f2c45f1SShri Abhyankar @*/
12565f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
12575f2c45f1SShri Abhyankar {
12585f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12595f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
12605f2c45f1SShri Abhyankar 
12615f2c45f1SShri Abhyankar   PetscFunctionBegin;
12625f2c45f1SShri Abhyankar   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
12635f2c45f1SShri Abhyankar   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
12645f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12655f2c45f1SShri Abhyankar }
12665f2c45f1SShri Abhyankar 
12675f2c45f1SShri Abhyankar /*@C
1268d842c372SHong Zhang   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
12695f2c45f1SShri Abhyankar 
12705f2c45f1SShri Abhyankar   Not Collective
12715f2c45f1SShri Abhyankar 
12725f2c45f1SShri Abhyankar   Input Parameters:
12735f2c45f1SShri Abhyankar + dm - The DMNetwork object
12745f2c45f1SShri Abhyankar - p  - the edge point
12755f2c45f1SShri Abhyankar 
12765f2c45f1SShri Abhyankar   Output Paramters:
12775f2c45f1SShri Abhyankar . vertices  - vertices connected to this edge
12785f2c45f1SShri Abhyankar 
12795f2c45f1SShri Abhyankar   Level: intermediate
12805f2c45f1SShri Abhyankar 
12815f2c45f1SShri Abhyankar   Fortran Notes:
12825f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
12835f2c45f1SShri Abhyankar   include petsc.h90 in your code.
12845f2c45f1SShri Abhyankar 
12855f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
12865f2c45f1SShri Abhyankar @*/
1287d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
12885f2c45f1SShri Abhyankar {
12895f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12905f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
12915f2c45f1SShri Abhyankar 
12925f2c45f1SShri Abhyankar   PetscFunctionBegin;
12935f2c45f1SShri Abhyankar   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
12945f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12955f2c45f1SShri Abhyankar }
12965f2c45f1SShri Abhyankar 
12975f2c45f1SShri Abhyankar /*@
12985f2c45f1SShri Abhyankar   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
12995f2c45f1SShri Abhyankar 
13005f2c45f1SShri Abhyankar   Not Collective
13015f2c45f1SShri Abhyankar 
13025f2c45f1SShri Abhyankar   Input Parameters:
13035f2c45f1SShri Abhyankar + dm - The DMNetwork object
13045f2c45f1SShri Abhyankar . p  - the vertex point
13055f2c45f1SShri Abhyankar 
13065f2c45f1SShri Abhyankar   Output Parameter:
13075f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point
13085f2c45f1SShri Abhyankar 
13095f2c45f1SShri Abhyankar   Level: intermediate
13105f2c45f1SShri Abhyankar 
1311d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange
13125f2c45f1SShri Abhyankar @*/
13135f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
13145f2c45f1SShri Abhyankar {
13155f2c45f1SShri Abhyankar   PetscErrorCode ierr;
13165f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
13175f2c45f1SShri Abhyankar   PetscInt       offsetg;
13185f2c45f1SShri Abhyankar   PetscSection   sectiong;
13195f2c45f1SShri Abhyankar 
13205f2c45f1SShri Abhyankar   PetscFunctionBegin;
13215f2c45f1SShri Abhyankar   *isghost = PETSC_FALSE;
13225f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
13235f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
13245f2c45f1SShri Abhyankar   if (offsetg < 0) *isghost = PETSC_TRUE;
13255f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
13265f2c45f1SShri Abhyankar }
13275f2c45f1SShri Abhyankar 
13285f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm)
13295f2c45f1SShri Abhyankar {
13305f2c45f1SShri Abhyankar   PetscErrorCode ierr;
13315f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
13325f2c45f1SShri Abhyankar 
13335f2c45f1SShri Abhyankar   PetscFunctionBegin;
13345f2c45f1SShri Abhyankar   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
13355f2c45f1SShri Abhyankar   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
13365f2c45f1SShri Abhyankar 
13375f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
13385f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
13395f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
13405f2c45f1SShri Abhyankar }
13415f2c45f1SShri Abhyankar 
13421ad426b7SHong Zhang /*@
134317df6e9eSHong Zhang     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
13441ad426b7SHong Zhang                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
13451ad426b7SHong Zhang 
13461ad426b7SHong Zhang     Collective
13471ad426b7SHong Zhang 
13481ad426b7SHong Zhang     Input Parameters:
134983b2e829SHong Zhang +   dm - The DMNetwork object
135083b2e829SHong Zhang .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
135183b2e829SHong Zhang -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
13521ad426b7SHong Zhang 
13531ad426b7SHong Zhang     Level: intermediate
13541ad426b7SHong Zhang 
13551ad426b7SHong Zhang @*/
135683b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
13571ad426b7SHong Zhang {
13581ad426b7SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
13598675203cSHong Zhang   PetscErrorCode ierr;
13601ad426b7SHong Zhang 
13611ad426b7SHong Zhang   PetscFunctionBegin;
136283b2e829SHong Zhang   network->userEdgeJacobian   = eflg;
136383b2e829SHong Zhang   network->userVertexJacobian = vflg;
13648675203cSHong Zhang 
13658675203cSHong Zhang   if (eflg && !network->Je) {
13668675203cSHong Zhang     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
13678675203cSHong Zhang   }
13688675203cSHong Zhang 
13698675203cSHong Zhang   if (vflg && !network->Jv) {
13708675203cSHong Zhang     PetscInt       i,*vptr,nedges,vStart=network->vStart;
13718675203cSHong Zhang     PetscInt       nVertices = network->nVertices,nedges_total;
13728675203cSHong Zhang     const PetscInt *edges;
13738675203cSHong Zhang 
13748675203cSHong Zhang     /* count nvertex_total */
13758675203cSHong Zhang     nedges_total = 0;
13768675203cSHong Zhang     ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr);
13778675203cSHong Zhang 
13788675203cSHong Zhang     vptr[0] = 0;
13798675203cSHong Zhang     for (i=0; i<nVertices; i++) {
13808675203cSHong Zhang       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
13818675203cSHong Zhang       nedges_total += nedges;
13828675203cSHong Zhang       vptr[i+1] = vptr[i] + 2*nedges + 1;
13838675203cSHong Zhang     }
13848675203cSHong Zhang 
13858675203cSHong Zhang     ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr);
13868675203cSHong Zhang     network->Jvptr = vptr;
13878675203cSHong Zhang   }
13881ad426b7SHong Zhang   PetscFunctionReturn(0);
13891ad426b7SHong Zhang }
13901ad426b7SHong Zhang 
13911ad426b7SHong Zhang /*@
139283b2e829SHong Zhang     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
139383b2e829SHong Zhang 
139483b2e829SHong Zhang     Not Collective
139583b2e829SHong Zhang 
139683b2e829SHong Zhang     Input Parameters:
139783b2e829SHong Zhang +   dm - The DMNetwork object
139883b2e829SHong Zhang .   p  - the edge point
13993e97b6e8SHong Zhang -   J - array (size = 3) of Jacobian submatrices for this edge point:
14003e97b6e8SHong Zhang         J[0]: this edge
1401d842c372SHong Zhang         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()
140283b2e829SHong Zhang 
140383b2e829SHong Zhang     Level: intermediate
140483b2e829SHong Zhang 
140583b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix
140683b2e829SHong Zhang @*/
140783b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
140883b2e829SHong Zhang {
140983b2e829SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
141083b2e829SHong Zhang 
141183b2e829SHong Zhang   PetscFunctionBegin;
14128675203cSHong Zhang   if (!network->Je) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
14138675203cSHong Zhang 
14148675203cSHong Zhang   if (J) {
1415883e35e8SHong Zhang     network->Je[3*p]   = J[0];
1416883e35e8SHong Zhang     network->Je[3*p+1] = J[1];
1417883e35e8SHong Zhang     network->Je[3*p+2] = J[2];
14188675203cSHong Zhang   }
141983b2e829SHong Zhang   PetscFunctionReturn(0);
142083b2e829SHong Zhang }
142183b2e829SHong Zhang 
142283b2e829SHong Zhang /*@
142376ddfea5SHong Zhang     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
14241ad426b7SHong Zhang 
14251ad426b7SHong Zhang     Not Collective
14261ad426b7SHong Zhang 
14271ad426b7SHong Zhang     Input Parameters:
14281ad426b7SHong Zhang +   dm - The DMNetwork object
14291ad426b7SHong Zhang .   p  - the vertex point
14303e97b6e8SHong Zhang -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
14313e97b6e8SHong Zhang         J[0]:       this vertex
14323e97b6e8SHong Zhang         J[1+2*i]:   i-th supporting edge
14333e97b6e8SHong Zhang         J[1+2*i+1]: i-th connected vertex
14341ad426b7SHong Zhang 
14351ad426b7SHong Zhang     Level: intermediate
14361ad426b7SHong Zhang 
143783b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix
14381ad426b7SHong Zhang @*/
1439883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
14405f2c45f1SShri Abhyankar {
14415f2c45f1SShri Abhyankar   PetscErrorCode ierr;
14425f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
14438675203cSHong Zhang   PetscInt       i,*vptr,nedges,vStart=network->vStart;
1444883e35e8SHong Zhang   const PetscInt *edges;
14455f2c45f1SShri Abhyankar 
14465f2c45f1SShri Abhyankar   PetscFunctionBegin;
14478675203cSHong Zhang   if (!network->Jv) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
1448883e35e8SHong Zhang 
14498675203cSHong Zhang   if (J) {
1450883e35e8SHong Zhang     vptr = network->Jvptr;
14513e97b6e8SHong Zhang     network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
14523e97b6e8SHong Zhang 
14533e97b6e8SHong Zhang     /* Set Jacobian for each supporting edge and connected vertex */
1454883e35e8SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
1455883e35e8SHong Zhang     for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
14568675203cSHong Zhang   }
1457883e35e8SHong Zhang   PetscFunctionReturn(0);
1458883e35e8SHong Zhang }
1459883e35e8SHong Zhang 
1460e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
14615cf7da58SHong Zhang {
14625cf7da58SHong Zhang   PetscErrorCode ierr;
14635cf7da58SHong Zhang   PetscInt       j;
14645cf7da58SHong Zhang   PetscScalar    val=(PetscScalar)ncols;
14655cf7da58SHong Zhang 
14665cf7da58SHong Zhang   PetscFunctionBegin;
14675cf7da58SHong Zhang   if (!ghost) {
14685cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
14695cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
14705cf7da58SHong Zhang     }
14715cf7da58SHong Zhang   } else {
14725cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
14735cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
14745cf7da58SHong Zhang     }
14755cf7da58SHong Zhang   }
14765cf7da58SHong Zhang   PetscFunctionReturn(0);
14775cf7da58SHong Zhang }
14785cf7da58SHong Zhang 
1479e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
14805cf7da58SHong Zhang {
14815cf7da58SHong Zhang   PetscErrorCode ierr;
14825cf7da58SHong Zhang   PetscInt       j,ncols_u;
14835cf7da58SHong Zhang   PetscScalar    val;
14845cf7da58SHong Zhang 
14855cf7da58SHong Zhang   PetscFunctionBegin;
14865cf7da58SHong Zhang   if (!ghost) {
14875cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
14885cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
14895cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
14905cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
14915cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
14925cf7da58SHong Zhang     }
14935cf7da58SHong Zhang   } else {
14945cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
14955cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
14965cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
14975cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
14985cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
14995cf7da58SHong Zhang     }
15005cf7da58SHong Zhang   }
15015cf7da58SHong Zhang   PetscFunctionReturn(0);
15025cf7da58SHong Zhang }
15035cf7da58SHong Zhang 
1504e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
15055cf7da58SHong Zhang {
15065cf7da58SHong Zhang   PetscErrorCode ierr;
15075cf7da58SHong Zhang 
15085cf7da58SHong Zhang   PetscFunctionBegin;
15095cf7da58SHong Zhang   if (Ju) {
15105cf7da58SHong Zhang     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
15115cf7da58SHong Zhang   } else {
15125cf7da58SHong Zhang     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
15135cf7da58SHong Zhang   }
15145cf7da58SHong Zhang   PetscFunctionReturn(0);
15155cf7da58SHong Zhang }
15165cf7da58SHong Zhang 
1517e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1518883e35e8SHong Zhang {
1519883e35e8SHong Zhang   PetscErrorCode ierr;
1520883e35e8SHong Zhang   PetscInt       j,*cols;
1521883e35e8SHong Zhang   PetscScalar    *zeros;
1522883e35e8SHong Zhang 
1523883e35e8SHong Zhang   PetscFunctionBegin;
1524883e35e8SHong Zhang   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
1525883e35e8SHong Zhang   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1526883e35e8SHong Zhang   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
1527883e35e8SHong Zhang   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
15281ad426b7SHong Zhang   PetscFunctionReturn(0);
15291ad426b7SHong Zhang }
1530a4e85ca8SHong Zhang 
1531e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
15323e97b6e8SHong Zhang {
15333e97b6e8SHong Zhang   PetscErrorCode ierr;
15343e97b6e8SHong Zhang   PetscInt       j,M,N,row,col,ncols_u;
15353e97b6e8SHong Zhang   const PetscInt *cols;
15363e97b6e8SHong Zhang   PetscScalar    zero=0.0;
15373e97b6e8SHong Zhang 
15383e97b6e8SHong Zhang   PetscFunctionBegin;
15393e97b6e8SHong Zhang   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
15403e97b6e8SHong 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);
15413e97b6e8SHong Zhang 
15423e97b6e8SHong Zhang   for (row=0; row<nrows; row++) {
15433e97b6e8SHong Zhang     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
15443e97b6e8SHong Zhang     for (j=0; j<ncols_u; j++) {
15453e97b6e8SHong Zhang       col = cols[j] + cstart;
15463e97b6e8SHong Zhang       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
15473e97b6e8SHong Zhang     }
15483e97b6e8SHong Zhang     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
15493e97b6e8SHong Zhang   }
15503e97b6e8SHong Zhang   PetscFunctionReturn(0);
15513e97b6e8SHong Zhang }
15521ad426b7SHong Zhang 
1553e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1554a4e85ca8SHong Zhang {
1555a4e85ca8SHong Zhang   PetscErrorCode ierr;
1556f4431b8cSHong Zhang 
1557a4e85ca8SHong Zhang   PetscFunctionBegin;
1558a4e85ca8SHong Zhang   if (Ju) {
1559a4e85ca8SHong Zhang     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1560a4e85ca8SHong Zhang   } else {
1561a4e85ca8SHong Zhang     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1562a4e85ca8SHong Zhang   }
1563a4e85ca8SHong Zhang   PetscFunctionReturn(0);
1564a4e85ca8SHong Zhang }
1565a4e85ca8SHong Zhang 
156624121865SAdrian 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.
156724121865SAdrian Maldonado */
156824121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
156924121865SAdrian Maldonado {
157024121865SAdrian Maldonado   PetscErrorCode ierr;
157124121865SAdrian Maldonado   PetscInt       i, size, dof;
157224121865SAdrian Maldonado   PetscInt       *glob2loc;
157324121865SAdrian Maldonado 
157424121865SAdrian Maldonado   PetscFunctionBegin;
157524121865SAdrian Maldonado   ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr);
157624121865SAdrian Maldonado   ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr);
157724121865SAdrian Maldonado 
157824121865SAdrian Maldonado   for (i = 0; i < size; i++) {
157924121865SAdrian Maldonado     ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr);
158024121865SAdrian Maldonado     dof = (dof >= 0) ? dof : -(dof + 1);
158124121865SAdrian Maldonado     glob2loc[i] = dof;
158224121865SAdrian Maldonado   }
158324121865SAdrian Maldonado 
158424121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
158524121865SAdrian Maldonado #if 0
158624121865SAdrian Maldonado   ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
158724121865SAdrian Maldonado #endif
158824121865SAdrian Maldonado   PetscFunctionReturn(0);
158924121865SAdrian Maldonado }
159024121865SAdrian Maldonado 
159101ad2aeeSHong Zhang #include <petsc/private/matimpl.h>
15921ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
15931ad426b7SHong Zhang {
15941ad426b7SHong Zhang   PetscErrorCode ierr;
159524121865SAdrian Maldonado   PetscMPIInt    rank, size;
15961ad426b7SHong Zhang   DM_Network     *network = (DM_Network*) dm->data;
1597a4e85ca8SHong Zhang   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1598840c2264SHong Zhang   PetscInt       cstart,ncols,j,e,v;
159924121865SAdrian Maldonado   PetscBool      ghost,ghost_vc,ghost2,isNest;
1600a4e85ca8SHong Zhang   Mat            Juser;
1601bfbc38dcSHong Zhang   PetscSection   sectionGlobal;
1602447d78afSSatish Balay   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
1603a4e85ca8SHong Zhang   const PetscInt *edges,*cone;
16045cf7da58SHong Zhang   MPI_Comm       comm;
160524121865SAdrian Maldonado   MatType        mtype;
16065cf7da58SHong Zhang   Vec            vd_nz,vo_nz;
16075cf7da58SHong Zhang   PetscInt       *dnnz,*onnz;
16085cf7da58SHong Zhang   PetscScalar    *vdnz,*vonz;
16091ad426b7SHong Zhang 
16101ad426b7SHong Zhang   PetscFunctionBegin;
161124121865SAdrian Maldonado   mtype = dm->mattype;
161224121865SAdrian Maldonado   ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr);
161324121865SAdrian Maldonado 
161424121865SAdrian Maldonado   if (isNest) {
16150731d606SHong Zhang     /* ierr = DMCreateMatrix_Network_Nest(); */
161624121865SAdrian Maldonado     PetscInt   eDof, vDof;
161724121865SAdrian Maldonado     Mat        j11, j12, j21, j22, bA[2][2];
161824121865SAdrian Maldonado     ISLocalToGlobalMapping eISMap, vISMap;
161924121865SAdrian Maldonado 
162024121865SAdrian Maldonado     ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
162124121865SAdrian Maldonado     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
162224121865SAdrian Maldonado     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
162324121865SAdrian Maldonado 
162424121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr);
162524121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr);
162624121865SAdrian Maldonado 
162701ad2aeeSHong Zhang     ierr = MatCreate(comm, &j11);CHKERRQ(ierr);
162824121865SAdrian Maldonado     ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
162924121865SAdrian Maldonado     ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr);
163024121865SAdrian Maldonado 
163101ad2aeeSHong Zhang     ierr = MatCreate(comm, &j12);CHKERRQ(ierr);
163224121865SAdrian Maldonado     ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr);
163324121865SAdrian Maldonado     ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr);
163424121865SAdrian Maldonado 
163501ad2aeeSHong Zhang     ierr = MatCreate(comm, &j21);CHKERRQ(ierr);
163624121865SAdrian Maldonado     ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
163724121865SAdrian Maldonado     ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr);
163824121865SAdrian Maldonado 
163901ad2aeeSHong Zhang     ierr = MatCreate(comm, &j22);CHKERRQ(ierr);
164024121865SAdrian Maldonado     ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
164124121865SAdrian Maldonado     ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr);
164224121865SAdrian Maldonado 
16433f6a6bdaSHong Zhang     bA[0][0] = j11;
16443f6a6bdaSHong Zhang     bA[0][1] = j12;
16453f6a6bdaSHong Zhang     bA[1][0] = j21;
16463f6a6bdaSHong Zhang     bA[1][1] = j22;
164724121865SAdrian Maldonado 
164824121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr);
164924121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr);
165024121865SAdrian Maldonado 
165124121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr);
165224121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr);
165324121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr);
165424121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr);
165524121865SAdrian Maldonado 
165624121865SAdrian Maldonado     ierr = MatSetUp(j11);CHKERRQ(ierr);
165724121865SAdrian Maldonado     ierr = MatSetUp(j12);CHKERRQ(ierr);
165824121865SAdrian Maldonado     ierr = MatSetUp(j21);CHKERRQ(ierr);
165924121865SAdrian Maldonado     ierr = MatSetUp(j22);CHKERRQ(ierr);
166024121865SAdrian Maldonado 
166101ad2aeeSHong Zhang     ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr);
166224121865SAdrian Maldonado     ierr = MatSetUp(*J);CHKERRQ(ierr);
166324121865SAdrian Maldonado     ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr);
166424121865SAdrian Maldonado     ierr = MatDestroy(&j11);CHKERRQ(ierr);
166524121865SAdrian Maldonado     ierr = MatDestroy(&j12);CHKERRQ(ierr);
166624121865SAdrian Maldonado     ierr = MatDestroy(&j21);CHKERRQ(ierr);
166724121865SAdrian Maldonado     ierr = MatDestroy(&j22);CHKERRQ(ierr);
166824121865SAdrian Maldonado 
166924121865SAdrian Maldonado     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
167024121865SAdrian Maldonado     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
167124121865SAdrian Maldonado 
167224121865SAdrian Maldonado     /* Free structures */
167324121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr);
167424121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr);
167524121865SAdrian Maldonado 
167624121865SAdrian Maldonado     PetscFunctionReturn(0);
167724121865SAdrian Maldonado   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1678a4e85ca8SHong Zhang     /* user does not provide Jacobian blocks */
1679bfbc38dcSHong Zhang     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1680bfbc38dcSHong Zhang     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
16811ad426b7SHong Zhang     PetscFunctionReturn(0);
16821ad426b7SHong Zhang   }
16831ad426b7SHong Zhang 
1684bfbc38dcSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
16852a945128SHong Zhang   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1686bfbc38dcSHong Zhang   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1687bfbc38dcSHong Zhang   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
16882a945128SHong Zhang 
16892a945128SHong Zhang   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
16902a945128SHong Zhang   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
169189898e50SHong Zhang 
169289898e50SHong Zhang   /* (1) Set matrix preallocation */
169389898e50SHong Zhang   /*------------------------------*/
1694840c2264SHong Zhang   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1695840c2264SHong Zhang   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1696840c2264SHong Zhang   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1697840c2264SHong Zhang   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1698840c2264SHong Zhang   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1699840c2264SHong Zhang   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1700840c2264SHong Zhang 
170189898e50SHong Zhang   /* Set preallocation for edges */
170289898e50SHong Zhang   /*-----------------------------*/
1703840c2264SHong Zhang   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1704840c2264SHong Zhang 
1705bdcb62a2SHong Zhang   ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr);
1706840c2264SHong Zhang   for (e=eStart; e<eEnd; e++) {
1707840c2264SHong Zhang     /* Get row indices */
1708840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1709840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1710840c2264SHong Zhang     if (nrows) {
1711840c2264SHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1712840c2264SHong Zhang 
17135cf7da58SHong Zhang       /* Set preallocation for conntected vertices */
1714d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1715840c2264SHong Zhang       for (v=0; v<2; v++) {
1716840c2264SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1717840c2264SHong Zhang 
17188675203cSHong Zhang         if (network->Je) {
1719840c2264SHong Zhang           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
17208675203cSHong Zhang         } else Juser = NULL;
1721840c2264SHong Zhang         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
17225cf7da58SHong Zhang         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1723840c2264SHong Zhang       }
1724840c2264SHong Zhang 
172589898e50SHong Zhang       /* Set preallocation for edge self */
1726840c2264SHong Zhang       cstart = rstart;
17278675203cSHong Zhang       if (network->Je) {
1728840c2264SHong Zhang         Juser = network->Je[3*e]; /* Jacobian(e,e) */
17298675203cSHong Zhang       } else Juser = NULL;
17305cf7da58SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1731840c2264SHong Zhang     }
1732840c2264SHong Zhang   }
1733840c2264SHong Zhang 
173489898e50SHong Zhang   /* Set preallocation for vertices */
173589898e50SHong Zhang   /*--------------------------------*/
1736840c2264SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
17378675203cSHong Zhang   if (vEnd - vStart) vptr = network->Jvptr;
1738840c2264SHong Zhang 
1739840c2264SHong Zhang   for (v=vStart; v<vEnd; v++) {
1740840c2264SHong Zhang     /* Get row indices */
1741840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1742840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1743840c2264SHong Zhang     if (!nrows) continue;
1744840c2264SHong Zhang 
1745bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1746bdcb62a2SHong Zhang     if (ghost) {
1747bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1748bdcb62a2SHong Zhang     } else {
1749bdcb62a2SHong Zhang       rows_v = rows;
1750bdcb62a2SHong Zhang     }
1751bdcb62a2SHong Zhang 
1752bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1753840c2264SHong Zhang 
1754840c2264SHong Zhang     /* Get supporting edges and connected vertices */
1755840c2264SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1756840c2264SHong Zhang 
1757840c2264SHong Zhang     for (e=0; e<nedges; e++) {
1758840c2264SHong Zhang       /* Supporting edges */
1759840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1760840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1761840c2264SHong Zhang 
17628675203cSHong Zhang       if (network->Jv) {
1763840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
17648675203cSHong Zhang       } else Juser = NULL;
1765bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1766840c2264SHong Zhang 
1767840c2264SHong Zhang       /* Connected vertices */
1768d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
1769840c2264SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
1770840c2264SHong Zhang       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1771840c2264SHong Zhang 
1772840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1773840c2264SHong Zhang 
17748675203cSHong Zhang       if (network->Jv) {
1775840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
17768675203cSHong Zhang       } else Juser = NULL;
1777e102a522SHong Zhang       if (ghost_vc||ghost) {
1778e102a522SHong Zhang         ghost2 = PETSC_TRUE;
1779e102a522SHong Zhang       } else {
1780e102a522SHong Zhang         ghost2 = PETSC_FALSE;
1781e102a522SHong Zhang       }
1782e102a522SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr);
1783840c2264SHong Zhang     }
1784840c2264SHong Zhang 
178589898e50SHong Zhang     /* Set preallocation for vertex self */
1786840c2264SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1787840c2264SHong Zhang     if (!ghost) {
1788840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
17898675203cSHong Zhang       if (network->Jv) {
1790840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
17918675203cSHong Zhang       } else Juser = NULL;
1792bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1793840c2264SHong Zhang     }
1794bdcb62a2SHong Zhang     if (ghost) {
1795bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1796bdcb62a2SHong Zhang     }
1797840c2264SHong Zhang   }
1798840c2264SHong Zhang 
1799840c2264SHong Zhang   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1800840c2264SHong Zhang   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
18015cf7da58SHong Zhang 
18025cf7da58SHong Zhang   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
18035cf7da58SHong Zhang 
18045cf7da58SHong Zhang   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1805840c2264SHong Zhang   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1806840c2264SHong Zhang 
1807840c2264SHong Zhang   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1808840c2264SHong Zhang   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1809840c2264SHong Zhang   for (j=0; j<localSize; j++) {
1810e102a522SHong Zhang     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1811e102a522SHong Zhang     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1812840c2264SHong Zhang   }
1813840c2264SHong Zhang   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1814840c2264SHong Zhang   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1815840c2264SHong Zhang   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1816840c2264SHong Zhang   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1817840c2264SHong Zhang 
18185cf7da58SHong Zhang   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
18195cf7da58SHong Zhang   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
18205cf7da58SHong Zhang   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
18215cf7da58SHong Zhang 
18225cf7da58SHong Zhang   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
18235cf7da58SHong Zhang 
182489898e50SHong Zhang   /* (2) Set matrix entries for edges */
182589898e50SHong Zhang   /*----------------------------------*/
18261ad426b7SHong Zhang   for (e=eStart; e<eEnd; e++) {
1827bfbc38dcSHong Zhang     /* Get row indices */
18281ad426b7SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
182917df6e9eSHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
18304b976069SHong Zhang     if (nrows) {
183117df6e9eSHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
18321ad426b7SHong Zhang 
1833bfbc38dcSHong Zhang       /* Set matrix entries for conntected vertices */
1834d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1835bfbc38dcSHong Zhang       for (v=0; v<2; v++) {
1836bfbc38dcSHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1837883e35e8SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
18383e97b6e8SHong Zhang 
18398675203cSHong Zhang         if (network->Je) {
1840a4e85ca8SHong Zhang           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
18418675203cSHong Zhang         } else Juser = NULL;
1842a4e85ca8SHong Zhang         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1843bfbc38dcSHong Zhang       }
184417df6e9eSHong Zhang 
1845bfbc38dcSHong Zhang       /* Set matrix entries for edge self */
18463e97b6e8SHong Zhang       cstart = rstart;
18478675203cSHong Zhang       if (network->Je) {
1848a4e85ca8SHong Zhang         Juser = network->Je[3*e]; /* Jacobian(e,e) */
18498675203cSHong Zhang       } else Juser = NULL;
1850a4e85ca8SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
18511ad426b7SHong Zhang     }
18524b976069SHong Zhang   }
18531ad426b7SHong Zhang 
1854bfbc38dcSHong Zhang   /* Set matrix entries for vertices */
185583b2e829SHong Zhang   /*---------------------------------*/
18561ad426b7SHong Zhang   for (v=vStart; v<vEnd; v++) {
1857bfbc38dcSHong Zhang     /* Get row indices */
1858596e729fSHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1859596e729fSHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
18604b976069SHong Zhang     if (!nrows) continue;
1861596e729fSHong Zhang 
1862bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1863bdcb62a2SHong Zhang     if (ghost) {
1864bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1865bdcb62a2SHong Zhang     } else {
1866bdcb62a2SHong Zhang       rows_v = rows;
1867bdcb62a2SHong Zhang     }
1868bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1869596e729fSHong Zhang 
1870bfbc38dcSHong Zhang     /* Get supporting edges and connected vertices */
1871596e729fSHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1872596e729fSHong Zhang 
1873596e729fSHong Zhang     for (e=0; e<nedges; e++) {
1874bfbc38dcSHong Zhang       /* Supporting edges */
1875596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1876596e729fSHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1877596e729fSHong Zhang 
18788675203cSHong Zhang       if (network->Jv) {
1879a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
18808675203cSHong Zhang       } else Juser = NULL;
1881bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1882596e729fSHong Zhang 
1883bfbc38dcSHong Zhang       /* Connected vertices */
1884d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
18852a945128SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
18862a945128SHong Zhang 
188744aca652SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
188844aca652SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1889a4e85ca8SHong Zhang 
18908675203cSHong Zhang       if (network->Jv) {
1891a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
18928675203cSHong Zhang       } else Juser = NULL;
1893bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1894596e729fSHong Zhang     }
1895596e729fSHong Zhang 
1896bfbc38dcSHong Zhang     /* Set matrix entries for vertex self */
18971ad426b7SHong Zhang     if (!ghost) {
1898596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
18998675203cSHong Zhang       if (network->Jv) {
1900a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
19018675203cSHong Zhang       } else Juser = NULL;
1902bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr);
1903bdcb62a2SHong Zhang     }
1904bdcb62a2SHong Zhang     if (ghost) {
1905bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1906bdcb62a2SHong Zhang     }
19071ad426b7SHong Zhang   }
1908a4e85ca8SHong Zhang   ierr = PetscFree(rows);CHKERRQ(ierr);
1909bdcb62a2SHong Zhang 
19101ad426b7SHong Zhang   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19111ad426b7SHong Zhang   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1912dd6f46cdSHong Zhang 
19135f2c45f1SShri Abhyankar   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
19145f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19155f2c45f1SShri Abhyankar }
19165f2c45f1SShri Abhyankar 
19175f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm)
19185f2c45f1SShri Abhyankar {
19195f2c45f1SShri Abhyankar   PetscErrorCode ierr;
19205f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
19212727e31bSShri Abhyankar   PetscInt       j;
19225f2c45f1SShri Abhyankar 
19235f2c45f1SShri Abhyankar   PetscFunctionBegin;
19248415c774SShri Abhyankar   if (--network->refct > 0) PetscFunctionReturn(0);
192583b2e829SHong Zhang   if (network->Je) {
192683b2e829SHong Zhang     ierr = PetscFree(network->Je);CHKERRQ(ierr);
192783b2e829SHong Zhang   }
192883b2e829SHong Zhang   if (network->Jv) {
1929883e35e8SHong Zhang     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
193083b2e829SHong Zhang     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
19311ad426b7SHong Zhang   }
193213c2a604SAdrian Maldonado 
193313c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr);
193413c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr);
193513c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr);
193613c2a604SAdrian Maldonado   if (network->vertex.sf) {
193713c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr);
193813c2a604SAdrian Maldonado   }
193913c2a604SAdrian Maldonado   /* edge */
194013c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr);
194113c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr);
194213c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr);
194313c2a604SAdrian Maldonado   if (network->edge.sf) {
194413c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr);
194513c2a604SAdrian Maldonado   }
19465f2c45f1SShri Abhyankar   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
19475f2c45f1SShri Abhyankar   network->edges = NULL;
19485f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
19495f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
195083b2e829SHong Zhang 
19512727e31bSShri Abhyankar   for(j=0; j < network->nsubnet; j++) {
19522727e31bSShri Abhyankar     ierr = PetscFree(network->subnet[j].edges);CHKERRQ(ierr);
19532727e31bSShri Abhyankar     ierr = PetscFree(network->subnet[j].vertices);CHKERRQ(ierr);
19542727e31bSShri Abhyankar   }
1955e2aaf10cSShri Abhyankar   ierr = PetscFree(network->subnet);CHKERRQ(ierr);
19565f2c45f1SShri Abhyankar   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
19575f2c45f1SShri Abhyankar   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
19585f2c45f1SShri Abhyankar   ierr = PetscFree(network->header);CHKERRQ(ierr);
19595f2c45f1SShri Abhyankar   ierr = PetscFree(network);CHKERRQ(ierr);
19605f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19615f2c45f1SShri Abhyankar }
19625f2c45f1SShri Abhyankar 
19635f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
19645f2c45f1SShri Abhyankar {
19655f2c45f1SShri Abhyankar   PetscErrorCode ierr;
19665f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
19675f2c45f1SShri Abhyankar 
19685f2c45f1SShri Abhyankar   PetscFunctionBegin;
19695f2c45f1SShri Abhyankar   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
19705f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19715f2c45f1SShri Abhyankar }
19725f2c45f1SShri Abhyankar 
19735f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
19745f2c45f1SShri Abhyankar {
19755f2c45f1SShri Abhyankar   PetscErrorCode ierr;
19765f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
19775f2c45f1SShri Abhyankar 
19785f2c45f1SShri Abhyankar   PetscFunctionBegin;
19795f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
19805f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19815f2c45f1SShri Abhyankar }
19825f2c45f1SShri Abhyankar 
19835f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
19845f2c45f1SShri Abhyankar {
19855f2c45f1SShri Abhyankar   PetscErrorCode ierr;
19865f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
19875f2c45f1SShri Abhyankar 
19885f2c45f1SShri Abhyankar   PetscFunctionBegin;
19895f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
19905f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19915f2c45f1SShri Abhyankar }
19925f2c45f1SShri Abhyankar 
19935f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
19945f2c45f1SShri Abhyankar {
19955f2c45f1SShri Abhyankar   PetscErrorCode ierr;
19965f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
19975f2c45f1SShri Abhyankar 
19985f2c45f1SShri Abhyankar   PetscFunctionBegin;
19995f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
20005f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
20015f2c45f1SShri Abhyankar }
20025f2c45f1SShri Abhyankar 
20035f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
20045f2c45f1SShri Abhyankar {
20055f2c45f1SShri Abhyankar   PetscErrorCode ierr;
20065f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
20075f2c45f1SShri Abhyankar 
20085f2c45f1SShri Abhyankar   PetscFunctionBegin;
20095f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
20105f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
20115f2c45f1SShri Abhyankar }
2012