xref: /petsc/src/dm/impls/network/network.c (revision e108cb999e836831742efd3a44f4221fe84ae696)
1af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h>  /*I  "petscdmnetwork.h"  I*/
25f2c45f1SShri Abhyankar #include <petscdmplex.h>
35f2c45f1SShri Abhyankar #include <petscsf.h>
45f2c45f1SShri Abhyankar 
55f2c45f1SShri Abhyankar /*@
6556ed216SShri Abhyankar   DMNetworkGetPlex - Gets the Plex DM associated with this network DM
7556ed216SShri Abhyankar 
8556ed216SShri Abhyankar   Not collective
9556ed216SShri Abhyankar 
10556ed216SShri Abhyankar   Input Parameters:
11556ed216SShri Abhyankar + netdm - the dm object
12556ed216SShri Abhyankar - plexmdm - the plex dm object
13556ed216SShri Abhyankar 
14556ed216SShri Abhyankar   Level: Advanced
15556ed216SShri Abhyankar 
16556ed216SShri Abhyankar .seealso: DMNetworkCreate()
17556ed216SShri Abhyankar @*/
18556ed216SShri Abhyankar PetscErrorCode DMNetworkGetPlex(DM netdm, DM *plexdm)
19556ed216SShri Abhyankar {
20556ed216SShri Abhyankar   DM_Network *network = (DM_Network*) netdm->data;
21556ed216SShri Abhyankar 
22556ed216SShri Abhyankar   PetscFunctionBegin;
23556ed216SShri Abhyankar   *plexdm = network->plex;
24556ed216SShri Abhyankar   PetscFunctionReturn(0);
25556ed216SShri Abhyankar }
26556ed216SShri Abhyankar 
27556ed216SShri Abhyankar /*@
28e2aaf10cSShri Abhyankar   DMNetworkSetSizes - Sets the number of subnetworks, local and global vertices and edges for each subnetwork.
295f2c45f1SShri Abhyankar 
305f2c45f1SShri Abhyankar   Collective on DM
315f2c45f1SShri Abhyankar 
325f2c45f1SShri Abhyankar   Input Parameters:
335f2c45f1SShri Abhyankar + dm - the dm object
34e2aaf10cSShri Abhyankar . Nsubnet - number of subnetworks
35f025b11dSHong Zhang . NsubnetCouple - number of coupling subnetworks
36f025b11dSHong Zhang . nV - number of local vertices for each subnetwork and coupling subnetwork
37f025b11dSHong Zhang . nE - number of local edges for each subnetwork and coupling subnetwork
38f025b11dSHong Zhang . NV - number of global vertices (or PETSC_DETERMINE) for each subnetwork and coupling subnetwork
39f025b11dSHong Zhang - NE - number of global edges (or PETSC_DETERMINE) for each subnetwork and coupling subnetwork
405f2c45f1SShri Abhyankar 
415f2c45f1SShri Abhyankar    Notes
425f2c45f1SShri Abhyankar    If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang.
435f2c45f1SShri Abhyankar 
44f025b11dSHong Zhang    You cannot change the sizes once they have been set. nV, nE, NV and NE are arrays of length Nsubnet+NsubnetCouple.
455f2c45f1SShri Abhyankar 
461b266c99SBarry Smith    Level: intermediate
471b266c99SBarry Smith 
481b266c99SBarry Smith .seealso: DMNetworkCreate()
495f2c45f1SShri Abhyankar @*/
50f025b11dSHong Zhang PetscErrorCode DMNetworkSetSizes(DM dm,PetscInt Nsubnet,PetscInt NsubnetCouple,PetscInt nV[], PetscInt nE[], PetscInt NV[], PetscInt NE[])
515f2c45f1SShri Abhyankar {
525f2c45f1SShri Abhyankar   PetscErrorCode ierr;
535f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
54e2aaf10cSShri Abhyankar   PetscInt       a[2],b[2],i;
555f2c45f1SShri Abhyankar 
565f2c45f1SShri Abhyankar   PetscFunctionBegin;
575f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
58e2aaf10cSShri Abhyankar   if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet);
597765340cSHong Zhang   if (NsubnetCouple < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of coupling subnetworks %D cannot be less than 0",NsubnetCouple);
607765340cSHong Zhang 
617765340cSHong Zhang   if (Nsubnet > 0) PetscValidLogicalCollectiveInt(dm,Nsubnet,2);
627765340cSHong Zhang   if (network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network");
637765340cSHong Zhang 
647765340cSHong Zhang   network->nsubnet  = Nsubnet + NsubnetCouple;
657765340cSHong Zhang   Nsubnet          += NsubnetCouple;
667765340cSHong Zhang   network->ncsubnet = NsubnetCouple;
677765340cSHong Zhang   ierr = PetscCalloc1(Nsubnet,&network->subnet);CHKERRQ(ierr);
687765340cSHong Zhang   for(i=0; i < Nsubnet; i++) {
69f025b11dSHong Zhang     if (NV) {
707765340cSHong Zhang       if (NV[i] > 0) PetscValidLogicalCollectiveInt(dm,NV[i],6);
717765340cSHong Zhang       if (NV[i] > 0 && nV[i] > NV[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local vertex size %D cannot be larger than global vertex size %D",i,nV[i],NV[i]);
72f025b11dSHong Zhang     }
73f025b11dSHong Zhang     if (NE) {
74f025b11dSHong Zhang       if (NE[i] > 0) PetscValidLogicalCollectiveInt(dm,NE[i],7);
757765340cSHong Zhang       if (NE[i] > 0 && nE[i] > NE[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local edge size %D cannot be larger than global edge size %D",i,nE[i],NE[i]);
76f025b11dSHong Zhang     }
77f025b11dSHong Zhang 
787765340cSHong Zhang     a[0] = nV[i]; a[1] = nE[i];
797765340cSHong Zhang     ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
807765340cSHong Zhang     network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1];
817765340cSHong Zhang 
827765340cSHong Zhang     network->subnet[i].id = i;
837765340cSHong Zhang 
847765340cSHong Zhang     network->subnet[i].nvtx = nV[i];
857765340cSHong Zhang     network->subnet[i].vStart = network->nVertices;
867765340cSHong Zhang     network->subnet[i].vEnd   = network->subnet[i].vStart + network->subnet[i].Nvtx;
877765340cSHong Zhang     network->nVertices += network->subnet[i].nvtx;
887765340cSHong Zhang     network->NVertices += network->subnet[i].Nvtx;
897765340cSHong Zhang 
907765340cSHong Zhang     network->subnet[i].nedge = nE[i];
917765340cSHong Zhang     network->subnet[i].eStart = network->nEdges;
927765340cSHong Zhang     network->subnet[i].eEnd = network->subnet[i].eStart + network->subnet[i].Nedge;
937765340cSHong Zhang     network->nEdges += network->subnet[i].nedge;
947765340cSHong Zhang     network->NEdges += network->subnet[i].Nedge;
957765340cSHong Zhang   }
967765340cSHong Zhang   PetscFunctionReturn(0);
977765340cSHong Zhang }
987765340cSHong Zhang 
995f2c45f1SShri Abhyankar /*@
1005f2c45f1SShri Abhyankar   DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network
1015f2c45f1SShri Abhyankar 
1025f2c45f1SShri Abhyankar   Logically collective on DM
1035f2c45f1SShri Abhyankar 
1045f2c45f1SShri Abhyankar   Input Parameters:
105e3e68989SHong Zhang + dm - the dm object
106e3e68989SHong Zhang . edgelist - list of edges for each subnetwork
107e3e68989SHong Zhang - edgelistCouple - list of edges for each coupling subnetwork
1085f2c45f1SShri Abhyankar 
1095f2c45f1SShri Abhyankar   Notes:
1105f2c45f1SShri Abhyankar   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
1115f2c45f1SShri Abhyankar   not be destroyed before the call to DMNetworkLayoutSetUp
1125f2c45f1SShri Abhyankar 
1135f2c45f1SShri Abhyankar   Level: intermediate
1145f2c45f1SShri Abhyankar 
115e3e68989SHong Zhang   Example usage:
116e3e68989SHong Zhang   Consider the following 2 separate networks and a coupling network:
117e3e68989SHong Zhang 
118e3e68989SHong Zhang .vb
119e3e68989SHong Zhang  network 0: v0 -> v1 -> v2 -> v3
120e3e68989SHong Zhang  network 1: v1 -> v2 -> v0
121e3e68989SHong Zhang  coupling network: network 1: v2 -> network 0: v0
122e3e68989SHong Zhang .ve
123e3e68989SHong Zhang 
124e3e68989SHong Zhang  The resulting input
125e3e68989SHong Zhang    edgelist[0] = [0 1 | 1 2 | 2 3];
126e3e68989SHong Zhang    edgelist[1] = [1 2 | 2 0]
127e3e68989SHong Zhang    edgelistCouple[0] = [(network)1 (v)2 (network)0 (v)0].
128e3e68989SHong Zhang 
1295f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes
1305f2c45f1SShri Abhyankar @*/
1314e18019cSBarry Smith PetscErrorCode DMNetworkSetEdgeList(DM dm,PetscInt *edgelist[],PetscInt *edgelistCouple[])
1325f2c45f1SShri Abhyankar {
1335f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network*) dm->data;
134e3e68989SHong Zhang   PetscInt   i,nsubnet,ncsubnet=network->ncsubnet;
1355f2c45f1SShri Abhyankar 
1365f2c45f1SShri Abhyankar   PetscFunctionBegin;
137e3e68989SHong Zhang   nsubnet = network->nsubnet - ncsubnet;
138e3e68989SHong Zhang   for(i=0; i < nsubnet; i++) {
139e2aaf10cSShri Abhyankar     network->subnet[i].edgelist = edgelist[i];
140e2aaf10cSShri Abhyankar   }
141e3e68989SHong Zhang   if (edgelistCouple) {
142e3e68989SHong Zhang     PetscInt j;
143e3e68989SHong Zhang     j = 0;
144e3e68989SHong Zhang     nsubnet = network->nsubnet;
145e3e68989SHong Zhang     while (i < nsubnet) network->subnet[i++].edgelist = edgelistCouple[j++];
146e3e68989SHong Zhang   }
1475f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
1485f2c45f1SShri Abhyankar }
1495f2c45f1SShri Abhyankar 
1505f2c45f1SShri Abhyankar /*@
1515f2c45f1SShri Abhyankar   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
1525f2c45f1SShri Abhyankar 
1535f2c45f1SShri Abhyankar   Collective on DM
1545f2c45f1SShri Abhyankar 
1555f2c45f1SShri Abhyankar   Input Parameters
1565f2c45f1SShri Abhyankar . DM - the dmnetwork object
1575f2c45f1SShri Abhyankar 
1585f2c45f1SShri Abhyankar   Notes:
1595f2c45f1SShri Abhyankar   This routine should be called after the network sizes and edgelists have been provided. It creates
1605f2c45f1SShri Abhyankar   the bare layout of the network and sets up the network to begin insertion of components.
1615f2c45f1SShri Abhyankar 
1625f2c45f1SShri Abhyankar   All the components should be registered before calling this routine.
1635f2c45f1SShri Abhyankar 
1645f2c45f1SShri Abhyankar   Level: intermediate
1655f2c45f1SShri Abhyankar 
1665f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList
1675f2c45f1SShri Abhyankar @*/
1685f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm)
1695f2c45f1SShri Abhyankar {
1705f2c45f1SShri Abhyankar   PetscErrorCode ierr;
1715f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
1725f2c45f1SShri Abhyankar   PetscInt       dim = 1; /* One dimensional network */
173991cf414SHong Zhang   PetscInt       numCorners=2,spacedim=2;
1746500d4abSHong Zhang   double         *vertexcoords=NULL;
1757765340cSHong Zhang   PetscInt       i,j,ndata,ctr=0,nsubnet;
1769ad3091eSHong Zhang   PetscInt       k,netid,vid;
177acdb140fSBarry Smith   PetscInt       *edgelist_couple=NULL;
1786500d4abSHong Zhang 
1796500d4abSHong Zhang   PetscFunctionBegin;
1806500d4abSHong Zhang   if (network->nVertices) {
1816500d4abSHong Zhang     ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr);
1826500d4abSHong Zhang   }
1836500d4abSHong Zhang 
1846500d4abSHong Zhang   /* Create the edgelist for the network by concatenating edgelists of the subnetworks */
1857765340cSHong Zhang   nsubnet = network->nsubnet - network->ncsubnet;
1866500d4abSHong Zhang   ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr);
1877765340cSHong Zhang   for (i=0; i < nsubnet; i++) {
1886500d4abSHong Zhang     for (j = 0; j < network->subnet[i].nedge; j++) {
1896500d4abSHong Zhang       network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
1906500d4abSHong Zhang       network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
1916500d4abSHong Zhang       ctr++;
1926500d4abSHong Zhang     }
1936500d4abSHong Zhang   }
1947765340cSHong Zhang 
1957765340cSHong Zhang   i       = nsubnet; /* netid of coupling subnet */
1967765340cSHong Zhang   nsubnet = network->nsubnet;
1977765340cSHong Zhang   while (i < nsubnet) {
198991cf414SHong Zhang     edgelist_couple = network->subnet[i].edgelist;
1996500d4abSHong Zhang     k = 0;
2006500d4abSHong Zhang     for (j = 0; j < network->subnet[i].nedge; j++) {
2016500d4abSHong Zhang       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
2026500d4abSHong Zhang       network->edges[2*ctr] = network->subnet[netid].vStart + vid; k += 2;
2036500d4abSHong Zhang 
2046500d4abSHong Zhang       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
205991cf414SHong Zhang       network->edges[2*ctr+1] = network->subnet[netid].vStart + vid; k+=2;
2066500d4abSHong Zhang       ctr++;
2076500d4abSHong Zhang     }
2087765340cSHong Zhang     i++;
2097765340cSHong Zhang   }
2106500d4abSHong Zhang 
21134bd1bf5SHong Zhang #if 0
2126500d4abSHong Zhang   for(i=0; i < network->nEdges; i++) {
2136500d4abSHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr);
2146500d4abSHong Zhang     printf("\n");
2156500d4abSHong Zhang   }
2166500d4abSHong Zhang #endif
2176500d4abSHong Zhang 
218acdb140fSBarry Smith #if defined(PETSC_USE_64BIT_INDICES)
219acdb140fSBarry Smith   {
220acdb140fSBarry Smith     int      *edges;
221acdb140fSBarry Smith     PetscInt ii;
222acdb140fSBarry Smith     ierr = PetscMalloc1(network->nEdges*numCorners,&edges);CHKERRQ(ierr);
223acdb140fSBarry Smith     for (ii=0; ii<network->nEdges*numCorners; ii++) {
224acdb140fSBarry Smith       edges[ii] = (int) network->edges[ii];
225acdb140fSBarry Smith     }
226acdb140fSBarry Smith     ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
227acdb140fSBarry Smith     ierr = PetscFree(edges);
228acdb140fSBarry Smith   }
229acdb140fSBarry Smith #else
2306500d4abSHong Zhang     ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
231acdb140fSBarry Smith  #endif
232acdb140fSBarry Smith 
2336500d4abSHong Zhang   if (network->nVertices) {
2346500d4abSHong Zhang     ierr = PetscFree(vertexcoords);CHKERRQ(ierr);
2356500d4abSHong Zhang   }
2366500d4abSHong Zhang   ierr = PetscFree(network->edges);CHKERRQ(ierr);
2376500d4abSHong Zhang 
2386500d4abSHong Zhang   ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr);
2396500d4abSHong Zhang   ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr);
2406500d4abSHong Zhang   ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr);
2416500d4abSHong Zhang 
2426500d4abSHong Zhang   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr);
2436500d4abSHong Zhang   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr);
2446500d4abSHong Zhang   ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr);
2456500d4abSHong Zhang   ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr);
2466500d4abSHong Zhang 
2476500d4abSHong Zhang   /* Create vertices and edges array for the subnetworks */
2486500d4abSHong Zhang   for (j=0; j < network->nsubnet; j++) {
2496500d4abSHong Zhang     ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr);
2506500d4abSHong Zhang     ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr);
2516500d4abSHong Zhang     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
2526500d4abSHong Zhang        These get updated when the vertices and edges are added. */
2536500d4abSHong Zhang     network->subnet[j].nvtx = network->subnet[j].nedge = 0;
2546500d4abSHong Zhang   }
2556500d4abSHong Zhang 
2566500d4abSHong Zhang   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
2576500d4abSHong Zhang   ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr);
2586500d4abSHong Zhang   for (i=network->eStart; i < network->eEnd; i++) {
2596500d4abSHong Zhang     network->header[i].index = i;   /* Global edge number */
2606500d4abSHong Zhang     for (j=0; j < network->nsubnet; j++) {
2616500d4abSHong Zhang       if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) {
2626500d4abSHong Zhang 	network->header[i].subnetid = j; /* Subnetwork id */
2636500d4abSHong Zhang 	network->subnet[j].edges[network->subnet[j].nedge++] = i;
2646500d4abSHong Zhang 	break;
2656500d4abSHong Zhang       }
2666500d4abSHong Zhang     }
2676500d4abSHong Zhang     network->header[i].ndata = 0;
2686500d4abSHong Zhang     ndata = network->header[i].ndata;
2696500d4abSHong Zhang     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
2706500d4abSHong Zhang     network->header[i].offset[ndata] = 0;
2716500d4abSHong Zhang   }
2726500d4abSHong Zhang 
2736500d4abSHong Zhang   for(i=network->vStart; i < network->vEnd; i++) {
2746500d4abSHong Zhang     network->header[i].index = i - network->vStart;
2756500d4abSHong Zhang     for (j=0; j < network->nsubnet; j++) {
2766500d4abSHong Zhang       if ((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) {
2776500d4abSHong Zhang 	network->header[i].subnetid = j;
2786500d4abSHong Zhang 	network->subnet[j].vertices[network->subnet[j].nvtx++] = i;
2796500d4abSHong Zhang 	break;
2806500d4abSHong Zhang       }
2816500d4abSHong Zhang     }
2826500d4abSHong Zhang     network->header[i].ndata = 0;
2836500d4abSHong Zhang     ndata = network->header[i].ndata;
2846500d4abSHong Zhang     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
2856500d4abSHong Zhang     network->header[i].offset[ndata] = 0;
2866500d4abSHong Zhang   }
2876500d4abSHong Zhang 
2886500d4abSHong Zhang   ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr);
2895f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2905f2c45f1SShri Abhyankar }
2915f2c45f1SShri Abhyankar 
29294ef8ddeSSatish Balay /*@C
2932727e31bSShri Abhyankar   DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork
2942727e31bSShri Abhyankar 
2952727e31bSShri Abhyankar   Input Parameters
2962727e31bSShri Abhyankar + dm   - the number object
2972727e31bSShri Abhyankar - id   - the ID (integer) of the subnetwork
2982727e31bSShri Abhyankar 
2992727e31bSShri Abhyankar   Output Parameters
3002727e31bSShri Abhyankar + nv    - number of vertices (local)
3012727e31bSShri Abhyankar . ne    - number of edges (local)
3022727e31bSShri Abhyankar . vtx   - local vertices for this subnetwork
3032727e31bSShri Abhyankar . edge  - local edges for this subnetwork
3042727e31bSShri Abhyankar 
3052727e31bSShri Abhyankar   Notes:
3062727e31bSShri Abhyankar   Cannot call this routine before DMNetworkLayoutSetup()
3072727e31bSShri Abhyankar 
3082727e31bSShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
3092727e31bSShri Abhyankar @*/
3102727e31bSShri Abhyankar PetscErrorCode DMNetworkGetSubnetworkInfo(DM netdm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge)
3112727e31bSShri Abhyankar {
3122727e31bSShri Abhyankar   DM_Network     *network = (DM_Network*) netdm->data;
3132727e31bSShri Abhyankar 
3142727e31bSShri Abhyankar   PetscFunctionBegin;
3152727e31bSShri Abhyankar   *nv = network->subnet[id].nvtx;
3162727e31bSShri Abhyankar   *ne = network->subnet[id].nedge;
3172727e31bSShri Abhyankar   *vtx = network->subnet[id].vertices;
3182727e31bSShri Abhyankar   *edge = network->subnet[id].edges;
3192727e31bSShri Abhyankar   PetscFunctionReturn(0);
3202727e31bSShri Abhyankar }
3212727e31bSShri Abhyankar 
3222727e31bSShri Abhyankar /*@C
3235f2c45f1SShri Abhyankar   DMNetworkRegisterComponent - Registers the network component
3245f2c45f1SShri Abhyankar 
3255f2c45f1SShri Abhyankar   Logically collective on DM
3265f2c45f1SShri Abhyankar 
3275f2c45f1SShri Abhyankar   Input Parameters
3285f2c45f1SShri Abhyankar + dm   - the network object
3295f2c45f1SShri Abhyankar . name - the component name
3305f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data
3315f2c45f1SShri Abhyankar 
3325f2c45f1SShri Abhyankar    Output Parameters
3335f2c45f1SShri Abhyankar .   key - an integer key that defines the component
3345f2c45f1SShri Abhyankar 
3355f2c45f1SShri Abhyankar    Notes
3365f2c45f1SShri Abhyankar    This routine should be called by all processors before calling DMNetworkLayoutSetup().
3375f2c45f1SShri Abhyankar 
3385f2c45f1SShri Abhyankar    Level: intermediate
3395f2c45f1SShri Abhyankar 
3405f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
3415f2c45f1SShri Abhyankar @*/
3425f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key)
3435f2c45f1SShri Abhyankar {
3445f2c45f1SShri Abhyankar   PetscErrorCode        ierr;
3455f2c45f1SShri Abhyankar   DM_Network            *network = (DM_Network*) dm->data;
3465f2c45f1SShri Abhyankar   DMNetworkComponent    *component=&network->component[network->ncomponent];
3475f2c45f1SShri Abhyankar   PetscBool             flg=PETSC_FALSE;
3485f2c45f1SShri Abhyankar   PetscInt              i;
3495f2c45f1SShri Abhyankar 
3505f2c45f1SShri Abhyankar   PetscFunctionBegin;
3515f2c45f1SShri Abhyankar   for (i=0; i < network->ncomponent; i++) {
3525f2c45f1SShri Abhyankar     ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr);
3535f2c45f1SShri Abhyankar     if (flg) {
3545f2c45f1SShri Abhyankar       *key = i;
3555f2c45f1SShri Abhyankar       PetscFunctionReturn(0);
3565f2c45f1SShri Abhyankar     }
3576d64e262SShri Abhyankar   }
3586d64e262SShri Abhyankar   if(network->ncomponent == MAX_COMPONENTS) {
3596d64e262SShri Abhyankar     SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components registered exceeds the max %D",MAX_COMPONENTS);
3605f2c45f1SShri Abhyankar   }
3615f2c45f1SShri Abhyankar 
3625f2c45f1SShri Abhyankar   ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr);
3635f2c45f1SShri Abhyankar   component->size = size/sizeof(DMNetworkComponentGenericDataType);
3645f2c45f1SShri Abhyankar   *key = network->ncomponent;
3655f2c45f1SShri Abhyankar   network->ncomponent++;
3665f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3675f2c45f1SShri Abhyankar }
3685f2c45f1SShri Abhyankar 
3695f2c45f1SShri Abhyankar /*@
3705f2c45f1SShri Abhyankar   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.
3715f2c45f1SShri Abhyankar 
3725f2c45f1SShri Abhyankar   Not Collective
3735f2c45f1SShri Abhyankar 
3745f2c45f1SShri Abhyankar   Input Parameters:
3755f2c45f1SShri Abhyankar + dm - The DMNetwork object
3765f2c45f1SShri Abhyankar 
3775f2c45f1SShri Abhyankar   Output Paramters:
3785f2c45f1SShri Abhyankar + vStart - The first vertex point
3795f2c45f1SShri Abhyankar - vEnd   - One beyond the last vertex point
3805f2c45f1SShri Abhyankar 
3815f2c45f1SShri Abhyankar   Level: intermediate
3825f2c45f1SShri Abhyankar 
3835f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange
3845f2c45f1SShri Abhyankar @*/
3855f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
3865f2c45f1SShri Abhyankar {
3875f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
3885f2c45f1SShri Abhyankar 
3895f2c45f1SShri Abhyankar   PetscFunctionBegin;
3905f2c45f1SShri Abhyankar   if (vStart) *vStart = network->vStart;
3915f2c45f1SShri Abhyankar   if (vEnd) *vEnd = network->vEnd;
3925f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3935f2c45f1SShri Abhyankar }
3945f2c45f1SShri Abhyankar 
3955f2c45f1SShri Abhyankar /*@
3965f2c45f1SShri Abhyankar   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.
3975f2c45f1SShri Abhyankar 
3985f2c45f1SShri Abhyankar   Not Collective
3995f2c45f1SShri Abhyankar 
4005f2c45f1SShri Abhyankar   Input Parameters:
4015f2c45f1SShri Abhyankar + dm - The DMNetwork object
4025f2c45f1SShri Abhyankar 
4035f2c45f1SShri Abhyankar   Output Paramters:
4045f2c45f1SShri Abhyankar + eStart - The first edge point
4055f2c45f1SShri Abhyankar - eEnd   - One beyond the last edge point
4065f2c45f1SShri Abhyankar 
4075f2c45f1SShri Abhyankar   Level: intermediate
4085f2c45f1SShri Abhyankar 
4095f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange
4105f2c45f1SShri Abhyankar @*/
4115f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
4125f2c45f1SShri Abhyankar {
4135f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4145f2c45f1SShri Abhyankar 
4155f2c45f1SShri Abhyankar   PetscFunctionBegin;
4165f2c45f1SShri Abhyankar   if (eStart) *eStart = network->eStart;
4175f2c45f1SShri Abhyankar   if (eEnd) *eEnd = network->eEnd;
4185f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4195f2c45f1SShri Abhyankar }
4205f2c45f1SShri Abhyankar 
4217b6afd5bSHong Zhang /*@
422e85e6aecSHong Zhang   DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge.
4237b6afd5bSHong Zhang 
4247b6afd5bSHong Zhang   Not Collective
4257b6afd5bSHong Zhang 
4267b6afd5bSHong Zhang   Input Parameters:
4277b6afd5bSHong Zhang + dm - DMNetwork object
428e85e6aecSHong Zhang - p  - edge point
4297b6afd5bSHong Zhang 
4307b6afd5bSHong Zhang   Output Paramters:
431e85e6aecSHong Zhang . index - user global numbering for the edge
4327b6afd5bSHong Zhang 
4337b6afd5bSHong Zhang   Level: intermediate
4347b6afd5bSHong Zhang 
435e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalVertexIndex
4367b6afd5bSHong Zhang @*/
437e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
4387b6afd5bSHong Zhang {
4397b6afd5bSHong Zhang   PetscErrorCode    ierr;
4407b6afd5bSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
4417b6afd5bSHong Zhang   PetscInt          offsetp;
4427b6afd5bSHong Zhang   DMNetworkComponentHeader header;
4437b6afd5bSHong Zhang 
4447b6afd5bSHong Zhang   PetscFunctionBegin;
4457b6afd5bSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
4467b6afd5bSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
447e85e6aecSHong Zhang   *index = header->index;
4487b6afd5bSHong Zhang   PetscFunctionReturn(0);
4497b6afd5bSHong Zhang }
4507b6afd5bSHong Zhang 
4515f2c45f1SShri Abhyankar /*@
452e85e6aecSHong Zhang   DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex.
453e85e6aecSHong Zhang 
454e85e6aecSHong Zhang   Not Collective
455e85e6aecSHong Zhang 
456e85e6aecSHong Zhang   Input Parameters:
457e85e6aecSHong Zhang + dm - DMNetwork object
458e85e6aecSHong Zhang - p  - vertex point
459e85e6aecSHong Zhang 
460e85e6aecSHong Zhang   Output Paramters:
461e85e6aecSHong Zhang . index - user global numbering for the vertex
462e85e6aecSHong Zhang 
463e85e6aecSHong Zhang   Level: intermediate
464e85e6aecSHong Zhang 
465e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalEdgeIndex
466e85e6aecSHong Zhang @*/
467e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
468e85e6aecSHong Zhang {
469e85e6aecSHong Zhang   PetscErrorCode    ierr;
470e85e6aecSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
471e85e6aecSHong Zhang   PetscInt          offsetp;
472e85e6aecSHong Zhang   DMNetworkComponentHeader header;
473e85e6aecSHong Zhang 
474e85e6aecSHong Zhang   PetscFunctionBegin;
475e85e6aecSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
476e85e6aecSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
477e85e6aecSHong Zhang   *index = header->index;
478e85e6aecSHong Zhang   PetscFunctionReturn(0);
479e85e6aecSHong Zhang }
480e85e6aecSHong Zhang 
481c3b11c7cSShri Abhyankar /*
482c3b11c7cSShri Abhyankar   DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the
483c3b11c7cSShri Abhyankar                                     component value from the component data array
484c3b11c7cSShri Abhyankar 
485c3b11c7cSShri Abhyankar   Not Collective
486c3b11c7cSShri Abhyankar 
487c3b11c7cSShri Abhyankar   Input Parameters:
488c3b11c7cSShri Abhyankar + dm      - The DMNetwork object
489c3b11c7cSShri Abhyankar . p       - vertex/edge point
490c3b11c7cSShri Abhyankar - compnum - component number
491c3b11c7cSShri Abhyankar 
492c3b11c7cSShri Abhyankar   Output Parameters:
493c3b11c7cSShri Abhyankar + compkey - the key obtained when registering the component
494c3b11c7cSShri Abhyankar - offset  - offset into the component data array associated with the vertex/edge point
495c3b11c7cSShri Abhyankar 
496c3b11c7cSShri Abhyankar   Notes:
497c3b11c7cSShri Abhyankar   Typical usage:
498c3b11c7cSShri Abhyankar 
499c3b11c7cSShri Abhyankar   DMNetworkGetComponentDataArray(dm, &arr);
500c3b11c7cSShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
501c3b11c7cSShri Abhyankar   Loop over vertices or edges
502c3b11c7cSShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
503c3b11c7cSShri Abhyankar     Loop over numcomps
504c3b11c7cSShri Abhyankar       DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset);
505c3b11c7cSShri Abhyankar       compdata = (UserCompDataType)(arr+offset);
506c3b11c7cSShri Abhyankar 
507c3b11c7cSShri Abhyankar   Level: intermediate
508c3b11c7cSShri Abhyankar 
509c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
510c3b11c7cSShri Abhyankar */
511c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
512c3b11c7cSShri Abhyankar {
513c3b11c7cSShri Abhyankar   PetscErrorCode           ierr;
514c3b11c7cSShri Abhyankar   PetscInt                 offsetp;
515c3b11c7cSShri Abhyankar   DMNetworkComponentHeader header;
516c3b11c7cSShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
517c3b11c7cSShri Abhyankar 
518c3b11c7cSShri Abhyankar   PetscFunctionBegin;
519c3b11c7cSShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
520c3b11c7cSShri Abhyankar   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
521c3b11c7cSShri Abhyankar   if (compkey) *compkey = header->key[compnum];
522c3b11c7cSShri Abhyankar   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
523c3b11c7cSShri Abhyankar   PetscFunctionReturn(0);
524c3b11c7cSShri Abhyankar }
525c3b11c7cSShri Abhyankar 
526c3b11c7cSShri Abhyankar /*@
527c3b11c7cSShri Abhyankar   DMNetworkGetComponent - Returns the network component and its key
528c3b11c7cSShri Abhyankar 
529c3b11c7cSShri Abhyankar   Not Collective
530c3b11c7cSShri Abhyankar 
531c3b11c7cSShri Abhyankar   Input Parameters
532c3b11c7cSShri Abhyankar + dm - DMNetwork object
533c3b11c7cSShri Abhyankar . p  - edge or vertex point
534c3b11c7cSShri Abhyankar - compnum - component number
535c3b11c7cSShri Abhyankar 
536c3b11c7cSShri Abhyankar   Output Parameters:
537c3b11c7cSShri Abhyankar + compkey - the key set for this computing during registration
538c3b11c7cSShri Abhyankar - component - the component data
539c3b11c7cSShri Abhyankar 
540c3b11c7cSShri Abhyankar   Notes:
541c3b11c7cSShri Abhyankar   Typical usage:
542c3b11c7cSShri Abhyankar 
543c3b11c7cSShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
544c3b11c7cSShri Abhyankar   Loop over vertices or edges
545c3b11c7cSShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
546c3b11c7cSShri Abhyankar     Loop over numcomps
547c3b11c7cSShri Abhyankar       DMNetworkGetComponent(dm,v,compnum,&key,&component);
548c3b11c7cSShri Abhyankar 
549c3b11c7cSShri Abhyankar   Level: intermediate
550c3b11c7cSShri Abhyankar 
551c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetVariableOffset
552c3b11c7cSShri Abhyankar @*/
553c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *key, void **component)
554c3b11c7cSShri Abhyankar {
555c3b11c7cSShri Abhyankar   PetscErrorCode ierr;
556c3b11c7cSShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
557*e108cb99SStefano Zampini   PetscInt       offsetd = 0;
558c3b11c7cSShri Abhyankar 
559c3b11c7cSShri Abhyankar   PetscFunctionBegin;
560c3b11c7cSShri Abhyankar   ierr = DMNetworkGetComponentKeyOffset(dm,p,compnum,key,&offsetd);CHKERRQ(ierr);
561c3b11c7cSShri Abhyankar   *component = network->componentdataarray+offsetd;
562c3b11c7cSShri Abhyankar   PetscFunctionReturn(0);
563c3b11c7cSShri Abhyankar }
564c3b11c7cSShri Abhyankar 
565e85e6aecSHong Zhang /*@
566325661f6SSatish Balay   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)
5675f2c45f1SShri Abhyankar 
5685f2c45f1SShri Abhyankar   Not Collective
5695f2c45f1SShri Abhyankar 
5705f2c45f1SShri Abhyankar   Input Parameters:
5715f2c45f1SShri Abhyankar + dm           - The DMNetwork object
5725f2c45f1SShri Abhyankar . p            - vertex/edge point
5735f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component
5745f2c45f1SShri Abhyankar - compvalue    - pointer to the data structure for the component
5755f2c45f1SShri Abhyankar 
5765f2c45f1SShri Abhyankar   Level: intermediate
5775f2c45f1SShri Abhyankar 
5785f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
5795f2c45f1SShri Abhyankar @*/
5805f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
5815f2c45f1SShri Abhyankar {
5825f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
58343a39a44SBarry Smith   DMNetworkComponent       *component = &network->component[componentkey];
5845f2c45f1SShri Abhyankar   DMNetworkComponentHeader header = &network->header[p];
5855f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue = &network->cvalue[p];
5865f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
5875f2c45f1SShri Abhyankar 
5885f2c45f1SShri Abhyankar   PetscFunctionBegin;
589fa58f0a9SHong 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);
590fa58f0a9SHong Zhang 
59143a39a44SBarry Smith   header->size[header->ndata] = component->size;
59243a39a44SBarry Smith   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
5935f2c45f1SShri Abhyankar   header->key[header->ndata] = componentkey;
5945f2c45f1SShri Abhyankar   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
5955f2c45f1SShri Abhyankar 
5965f2c45f1SShri Abhyankar   cvalue->data[header->ndata] = (void*)compvalue;
5975f2c45f1SShri Abhyankar   header->ndata++;
5985f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5995f2c45f1SShri Abhyankar }
6005f2c45f1SShri Abhyankar 
6015f2c45f1SShri Abhyankar /*@
6025f2c45f1SShri Abhyankar   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
6035f2c45f1SShri Abhyankar 
6045f2c45f1SShri Abhyankar   Not Collective
6055f2c45f1SShri Abhyankar 
6065f2c45f1SShri Abhyankar   Input Parameters:
6075f2c45f1SShri Abhyankar + dm - The DMNetwork object
6085f2c45f1SShri Abhyankar . p  - vertex/edge point
6095f2c45f1SShri Abhyankar 
6105f2c45f1SShri Abhyankar   Output Parameters:
6115f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge
6125f2c45f1SShri Abhyankar 
6135f2c45f1SShri Abhyankar   Level: intermediate
6145f2c45f1SShri Abhyankar 
6155f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
6165f2c45f1SShri Abhyankar @*/
6175f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
6185f2c45f1SShri Abhyankar {
6195f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6205f2c45f1SShri Abhyankar   PetscInt       offset;
6215f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6225f2c45f1SShri Abhyankar 
6235f2c45f1SShri Abhyankar   PetscFunctionBegin;
6245f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
6255f2c45f1SShri Abhyankar   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
6265f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6275f2c45f1SShri Abhyankar }
6285f2c45f1SShri Abhyankar 
6295f2c45f1SShri Abhyankar /*@
6305f2c45f1SShri Abhyankar   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
6315f2c45f1SShri Abhyankar 
6325f2c45f1SShri Abhyankar   Not Collective
6335f2c45f1SShri Abhyankar 
6345f2c45f1SShri Abhyankar   Input Parameters:
6355f2c45f1SShri Abhyankar + dm     - The DMNetwork object
6365f2c45f1SShri Abhyankar - p      - the edge/vertex point
6375f2c45f1SShri Abhyankar 
6385f2c45f1SShri Abhyankar   Output Parameters:
6395f2c45f1SShri Abhyankar . offset - the offset
6405f2c45f1SShri Abhyankar 
6415f2c45f1SShri Abhyankar   Level: intermediate
6425f2c45f1SShri Abhyankar 
6435f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
6445f2c45f1SShri Abhyankar @*/
6455f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
6465f2c45f1SShri Abhyankar {
6475f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6485f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6495f2c45f1SShri Abhyankar 
6505f2c45f1SShri Abhyankar   PetscFunctionBegin;
6515f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr);
6525f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6535f2c45f1SShri Abhyankar }
6545f2c45f1SShri Abhyankar 
6555f2c45f1SShri Abhyankar /*@
6565f2c45f1SShri Abhyankar   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
6575f2c45f1SShri Abhyankar 
6585f2c45f1SShri Abhyankar   Not Collective
6595f2c45f1SShri Abhyankar 
6605f2c45f1SShri Abhyankar   Input Parameters:
6615f2c45f1SShri Abhyankar + dm      - The DMNetwork object
6625f2c45f1SShri Abhyankar - p       - the edge/vertex point
6635f2c45f1SShri Abhyankar 
6645f2c45f1SShri Abhyankar   Output Parameters:
6655f2c45f1SShri Abhyankar . offsetg - the offset
6665f2c45f1SShri Abhyankar 
6675f2c45f1SShri Abhyankar   Level: intermediate
6685f2c45f1SShri Abhyankar 
6695f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
6705f2c45f1SShri Abhyankar @*/
6715f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
6725f2c45f1SShri Abhyankar {
6735f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6745f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6755f2c45f1SShri Abhyankar 
6765f2c45f1SShri Abhyankar   PetscFunctionBegin;
6775f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr);
6786fefedf4SHong Zhang   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */
6795f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6805f2c45f1SShri Abhyankar }
6815f2c45f1SShri Abhyankar 
68224121865SAdrian Maldonado /*@
68324121865SAdrian Maldonado   DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector.
68424121865SAdrian Maldonado 
68524121865SAdrian Maldonado   Not Collective
68624121865SAdrian Maldonado 
68724121865SAdrian Maldonado   Input Parameters:
68824121865SAdrian Maldonado + dm     - The DMNetwork object
68924121865SAdrian Maldonado - p      - the edge point
69024121865SAdrian Maldonado 
69124121865SAdrian Maldonado   Output Parameters:
69224121865SAdrian Maldonado . offset - the offset
69324121865SAdrian Maldonado 
69424121865SAdrian Maldonado   Level: intermediate
69524121865SAdrian Maldonado 
69624121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
69724121865SAdrian Maldonado @*/
69824121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
69924121865SAdrian Maldonado {
70024121865SAdrian Maldonado   PetscErrorCode ierr;
70124121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
70224121865SAdrian Maldonado 
70324121865SAdrian Maldonado   PetscFunctionBegin;
70424121865SAdrian Maldonado 
70524121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr);
70624121865SAdrian Maldonado   PetscFunctionReturn(0);
70724121865SAdrian Maldonado }
70824121865SAdrian Maldonado 
70924121865SAdrian Maldonado /*@
71024121865SAdrian Maldonado   DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector.
71124121865SAdrian Maldonado 
71224121865SAdrian Maldonado   Not Collective
71324121865SAdrian Maldonado 
71424121865SAdrian Maldonado   Input Parameters:
71524121865SAdrian Maldonado + dm     - The DMNetwork object
71624121865SAdrian Maldonado - p      - the vertex point
71724121865SAdrian Maldonado 
71824121865SAdrian Maldonado   Output Parameters:
71924121865SAdrian Maldonado . offset - the offset
72024121865SAdrian Maldonado 
72124121865SAdrian Maldonado   Level: intermediate
72224121865SAdrian Maldonado 
72324121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
72424121865SAdrian Maldonado @*/
72524121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
72624121865SAdrian Maldonado {
72724121865SAdrian Maldonado   PetscErrorCode ierr;
72824121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
72924121865SAdrian Maldonado 
73024121865SAdrian Maldonado   PetscFunctionBegin;
73124121865SAdrian Maldonado 
73224121865SAdrian Maldonado   p -= network->vStart;
73324121865SAdrian Maldonado 
73424121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr);
73524121865SAdrian Maldonado   PetscFunctionReturn(0);
73624121865SAdrian Maldonado }
7375f2c45f1SShri Abhyankar /*@
7385f2c45f1SShri Abhyankar   DMNetworkAddNumVariables - Add number of variables associated with a given point.
7395f2c45f1SShri Abhyankar 
7405f2c45f1SShri Abhyankar   Not Collective
7415f2c45f1SShri Abhyankar 
7425f2c45f1SShri Abhyankar   Input Parameters:
7435f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
7445f2c45f1SShri Abhyankar . p    - the vertex/edge point
7455f2c45f1SShri Abhyankar - nvar - number of additional variables
7465f2c45f1SShri Abhyankar 
7475f2c45f1SShri Abhyankar   Level: intermediate
7485f2c45f1SShri Abhyankar 
7495f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables
7505f2c45f1SShri Abhyankar @*/
7515f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
7525f2c45f1SShri Abhyankar {
7535f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7545f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7555f2c45f1SShri Abhyankar 
7565f2c45f1SShri Abhyankar   PetscFunctionBegin;
7575f2c45f1SShri Abhyankar   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
7585f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7595f2c45f1SShri Abhyankar }
7605f2c45f1SShri Abhyankar 
76127f51fceSHong Zhang /*@
76227f51fceSHong Zhang   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
76327f51fceSHong Zhang 
76427f51fceSHong Zhang   Not Collective
76527f51fceSHong Zhang 
76627f51fceSHong Zhang   Input Parameters:
76727f51fceSHong Zhang + dm   - The DMNetworkObject
76827f51fceSHong Zhang - p    - the vertex/edge point
76927f51fceSHong Zhang 
77027f51fceSHong Zhang   Output Parameters:
77127f51fceSHong Zhang . nvar - number of variables
77227f51fceSHong Zhang 
77327f51fceSHong Zhang   Level: intermediate
77427f51fceSHong Zhang 
77527f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
77627f51fceSHong Zhang @*/
77727f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
77827f51fceSHong Zhang {
77927f51fceSHong Zhang   PetscErrorCode ierr;
78027f51fceSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
78127f51fceSHong Zhang 
78227f51fceSHong Zhang   PetscFunctionBegin;
78327f51fceSHong Zhang   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
78427f51fceSHong Zhang   PetscFunctionReturn(0);
78527f51fceSHong Zhang }
78627f51fceSHong Zhang 
7875f2c45f1SShri Abhyankar /*@
7885f2c45f1SShri Abhyankar   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
7895f2c45f1SShri Abhyankar 
7905f2c45f1SShri Abhyankar   Not Collective
7915f2c45f1SShri Abhyankar 
7925f2c45f1SShri Abhyankar   Input Parameters:
7935f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
7945f2c45f1SShri Abhyankar . p    - the vertex/edge point
7955f2c45f1SShri Abhyankar - nvar - number of variables
7965f2c45f1SShri Abhyankar 
7975f2c45f1SShri Abhyankar   Level: intermediate
7985f2c45f1SShri Abhyankar 
7995f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables
8005f2c45f1SShri Abhyankar @*/
8015f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
8025f2c45f1SShri Abhyankar {
8035f2c45f1SShri Abhyankar   PetscErrorCode ierr;
8045f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
8055f2c45f1SShri Abhyankar 
8065f2c45f1SShri Abhyankar   PetscFunctionBegin;
8075f2c45f1SShri Abhyankar   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
8085f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
8095f2c45f1SShri Abhyankar }
8105f2c45f1SShri Abhyankar 
8115f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This
8125f2c45f1SShri Abhyankar    function is called during DMSetUp() */
8135f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm)
8145f2c45f1SShri Abhyankar {
8155f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
8165f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
817e53b5ba3SHong Zhang   PetscInt                 arr_size,p,offset,offsetp,ncomp,i;
8185f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
8195f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue;
8205f2c45f1SShri Abhyankar   DMNetworkComponentGenericDataType *componentdataarray;
8215f2c45f1SShri Abhyankar 
8225f2c45f1SShri Abhyankar   PetscFunctionBegin;
8235f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
8245f2c45f1SShri Abhyankar   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
82575b160a0SShri Abhyankar   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
8265f2c45f1SShri Abhyankar   componentdataarray = network->componentdataarray;
8275f2c45f1SShri Abhyankar   for (p = network->pStart; p < network->pEnd; p++) {
8285f2c45f1SShri Abhyankar     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
8295f2c45f1SShri Abhyankar     /* Copy header */
8305f2c45f1SShri Abhyankar     header = &network->header[p];
831302440fdSBarry Smith     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
8325f2c45f1SShri Abhyankar     /* Copy data */
8335f2c45f1SShri Abhyankar     cvalue = &network->cvalue[p];
8345f2c45f1SShri Abhyankar     ncomp = header->ndata;
8355f2c45f1SShri Abhyankar     for (i = 0; i < ncomp; i++) {
8365f2c45f1SShri Abhyankar       offset = offsetp + network->dataheadersize + header->offset[i];
837302440fdSBarry Smith       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
8385f2c45f1SShri Abhyankar     }
8395f2c45f1SShri Abhyankar   }
8405f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
8415f2c45f1SShri Abhyankar }
8425f2c45f1SShri Abhyankar 
8435f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */
8445f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm)
8455f2c45f1SShri Abhyankar {
8465f2c45f1SShri Abhyankar   PetscErrorCode ierr;
8475f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
8485f2c45f1SShri Abhyankar 
8495f2c45f1SShri Abhyankar   PetscFunctionBegin;
8505f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
8515f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
8525f2c45f1SShri Abhyankar }
8535f2c45f1SShri Abhyankar 
8545f2c45f1SShri Abhyankar /*@C
8555f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray - Returns the component data array
8565f2c45f1SShri Abhyankar 
8575f2c45f1SShri Abhyankar   Not Collective
8585f2c45f1SShri Abhyankar 
8595f2c45f1SShri Abhyankar   Input Parameters:
8605f2c45f1SShri Abhyankar . dm - The DMNetwork Object
8615f2c45f1SShri Abhyankar 
8625f2c45f1SShri Abhyankar   Output Parameters:
8635f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components
8645f2c45f1SShri Abhyankar 
8655f2c45f1SShri Abhyankar   Level: intermediate
8665f2c45f1SShri Abhyankar 
867a730d845SHong Zhang .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents
8685f2c45f1SShri Abhyankar @*/
8695f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
8705f2c45f1SShri Abhyankar {
8715f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
8725f2c45f1SShri Abhyankar 
8735f2c45f1SShri Abhyankar   PetscFunctionBegin;
8745f2c45f1SShri Abhyankar   *componentdataarray = network->componentdataarray;
8755f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
8765f2c45f1SShri Abhyankar }
8775f2c45f1SShri Abhyankar 
87824121865SAdrian Maldonado /* Get a subsection from a range of points */
87924121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection)
88024121865SAdrian Maldonado {
88124121865SAdrian Maldonado   PetscErrorCode ierr;
88224121865SAdrian Maldonado   PetscInt       i, nvar;
88324121865SAdrian Maldonado 
88424121865SAdrian Maldonado   PetscFunctionBegin;
88524121865SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr);
88624121865SAdrian Maldonado   ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr);
88724121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
88824121865SAdrian Maldonado     ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr);
88924121865SAdrian Maldonado     ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr);
89024121865SAdrian Maldonado   }
89124121865SAdrian Maldonado 
89224121865SAdrian Maldonado   ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr);
89324121865SAdrian Maldonado   PetscFunctionReturn(0);
89424121865SAdrian Maldonado }
89524121865SAdrian Maldonado 
89624121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */
89724121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
89824121865SAdrian Maldonado {
89924121865SAdrian Maldonado   PetscErrorCode ierr;
90024121865SAdrian Maldonado   PetscInt       i, *subpoints;
90124121865SAdrian Maldonado 
90224121865SAdrian Maldonado   PetscFunctionBegin;
90324121865SAdrian Maldonado   /* Create index sets to map from "points" to "subpoints" */
90424121865SAdrian Maldonado   ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr);
90524121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
90624121865SAdrian Maldonado     subpoints[i - pstart] = i;
90724121865SAdrian Maldonado   }
908459726d8SSatish Balay   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr);
90924121865SAdrian Maldonado   ierr = PetscFree(subpoints);CHKERRQ(ierr);
91024121865SAdrian Maldonado   PetscFunctionReturn(0);
91124121865SAdrian Maldonado }
91224121865SAdrian Maldonado 
91324121865SAdrian Maldonado /*@
91424121865SAdrian Maldonado   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute.
91524121865SAdrian Maldonado 
91624121865SAdrian Maldonado   Collective
91724121865SAdrian Maldonado 
91824121865SAdrian Maldonado   Input Parameters:
91924121865SAdrian Maldonado . dm   - The DMNetworkObject
92024121865SAdrian Maldonado 
92124121865SAdrian Maldonado   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
92224121865SAdrian Maldonado 
92324121865SAdrian Maldonado   points = [0 1 2 3 4 5 6]
92424121865SAdrian Maldonado 
92524121865SAdrian Maldonado   where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]).
92624121865SAdrian Maldonado 
92724121865SAdrian Maldonado   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
92824121865SAdrian Maldonado 
92924121865SAdrian Maldonado   Level: intermediate
93024121865SAdrian Maldonado 
93124121865SAdrian Maldonado @*/
93224121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
93324121865SAdrian Maldonado {
93424121865SAdrian Maldonado   PetscErrorCode ierr;
93524121865SAdrian Maldonado   MPI_Comm       comm;
9369852e123SBarry Smith   PetscMPIInt    rank, size;
93724121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
93824121865SAdrian Maldonado 
939eab1376dSHong Zhang   PetscFunctionBegin;
94024121865SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
94124121865SAdrian Maldonado   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
9429852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
94324121865SAdrian Maldonado 
94424121865SAdrian Maldonado   /* Create maps for vertices and edges */
94524121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr);
94624121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr);
94724121865SAdrian Maldonado 
94824121865SAdrian Maldonado   /* Create local sub-sections */
94924121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr);
95024121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr);
95124121865SAdrian Maldonado 
9529852e123SBarry Smith   if (size > 1) {
95324121865SAdrian Maldonado     ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr);
95424121865SAdrian Maldonado     ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr);
95524121865SAdrian Maldonado   ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr);
95624121865SAdrian Maldonado   ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr);
95724121865SAdrian Maldonado   } else {
95824121865SAdrian Maldonado   /* create structures for vertex */
95924121865SAdrian Maldonado   ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr);
96024121865SAdrian Maldonado   /* create structures for edge */
96124121865SAdrian Maldonado   ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr);
96224121865SAdrian Maldonado   }
96324121865SAdrian Maldonado 
96424121865SAdrian Maldonado 
96524121865SAdrian Maldonado   /* Add viewers */
96624121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr);
96724121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr);
96824121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr);
96924121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr);
97024121865SAdrian Maldonado 
97124121865SAdrian Maldonado   PetscFunctionReturn(0);
97224121865SAdrian Maldonado }
9737b6afd5bSHong Zhang 
9745f2c45f1SShri Abhyankar /*@
9755f2c45f1SShri Abhyankar   DMNetworkDistribute - Distributes the network and moves associated component data.
9765f2c45f1SShri Abhyankar 
9775f2c45f1SShri Abhyankar   Collective
9785f2c45f1SShri Abhyankar 
9795f2c45f1SShri Abhyankar   Input Parameter:
980d3464fd4SAdrian Maldonado + DM - the DMNetwork object
9815f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default
9825f2c45f1SShri Abhyankar 
9835f2c45f1SShri Abhyankar   Notes:
9848b171c8eSHong Zhang   Distributes the network with <overlap>-overlapping partitioning of the edges.
9855f2c45f1SShri Abhyankar 
9865f2c45f1SShri Abhyankar   Level: intermediate
9875f2c45f1SShri Abhyankar 
9885f2c45f1SShri Abhyankar .seealso: DMNetworkCreate
9895f2c45f1SShri Abhyankar @*/
990d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
9915f2c45f1SShri Abhyankar {
992d3464fd4SAdrian Maldonado   MPI_Comm       comm;
9935f2c45f1SShri Abhyankar   PetscErrorCode ierr;
994d3464fd4SAdrian Maldonado   PetscMPIInt    size;
995d3464fd4SAdrian Maldonado   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
996d3464fd4SAdrian Maldonado   DM_Network     *newDMnetwork;
9975f2c45f1SShri Abhyankar   PetscSF        pointsf;
9985f2c45f1SShri Abhyankar   DM             newDM;
99951ac5effSHong Zhang   PetscPartitioner part;
1000b9c6e19dSShri Abhyankar   PetscInt         j,e,v,offset;
1001b9c6e19dSShri Abhyankar   DMNetworkComponentHeader header;
10025f2c45f1SShri Abhyankar 
10035f2c45f1SShri Abhyankar   PetscFunctionBegin;
1004d3464fd4SAdrian Maldonado 
1005d3464fd4SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr);
1006d3464fd4SAdrian Maldonado   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1007d3464fd4SAdrian Maldonado   if (size == 1) PetscFunctionReturn(0);
1008d3464fd4SAdrian Maldonado 
1009d3464fd4SAdrian Maldonado   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr);
10105f2c45f1SShri Abhyankar   newDMnetwork = (DM_Network*)newDM->data;
10115f2c45f1SShri Abhyankar   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
101251ac5effSHong Zhang 
101351ac5effSHong Zhang   /* Enable runtime options for petscpartitioner */
101451ac5effSHong Zhang   ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr);
101551ac5effSHong Zhang   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
101651ac5effSHong Zhang 
10175f2c45f1SShri Abhyankar   /* Distribute plex dm and dof section */
101880cf41d5SMatthew G. Knepley   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
101951ac5effSHong Zhang 
10205f2c45f1SShri Abhyankar   /* Distribute dof section */
1021d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr);
10225f2c45f1SShri Abhyankar   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
1023d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr);
102451ac5effSHong Zhang 
10255f2c45f1SShri Abhyankar   /* Distribute data and associated section */
102631da1fc8SHong Zhang   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
102724121865SAdrian Maldonado 
10285f2c45f1SShri Abhyankar   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
10295f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
10305f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
10315f2c45f1SShri Abhyankar   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
10326fefedf4SHong Zhang   newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
10336fefedf4SHong Zhang   newDMnetwork->NVertices = oldDMnetwork->NVertices;
10345f2c45f1SShri Abhyankar   newDMnetwork->NEdges = oldDMnetwork->NEdges;
103524121865SAdrian Maldonado 
10365f2c45f1SShri Abhyankar   /* Set Dof section as the default section for dm */
10375f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
10385f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
10395f2c45f1SShri Abhyankar 
1040b9c6e19dSShri Abhyankar   /* Set up subnetwork info in the newDM */
1041b9c6e19dSShri Abhyankar   newDMnetwork->nsubnet = oldDMnetwork->nsubnet;
1042b9c6e19dSShri Abhyankar   ierr = PetscCalloc1(newDMnetwork->nsubnet,&newDMnetwork->subnet);CHKERRQ(ierr);
1043b9c6e19dSShri Abhyankar   /* Copy over the global number of vertices and edges in each subnetwork. Note that these are already
1044b9c6e19dSShri Abhyankar      calculated in DMNetworkLayoutSetUp()
1045b9c6e19dSShri Abhyankar   */
1046b9c6e19dSShri Abhyankar   for(j=0; j < newDMnetwork->nsubnet; j++) {
1047b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx;
1048b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge;
1049b9c6e19dSShri Abhyankar   }
1050b9c6e19dSShri Abhyankar 
1051b9c6e19dSShri Abhyankar   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1052b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
1053b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1054b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].nedge++;
1055b9c6e19dSShri Abhyankar   }
1056b9c6e19dSShri Abhyankar 
1057b9c6e19dSShri Abhyankar   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1058b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
1059b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1060b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].nvtx++;
1061b9c6e19dSShri Abhyankar   }
1062b9c6e19dSShri Abhyankar 
1063b9c6e19dSShri Abhyankar   /* Now create the vertices and edge arrays for the subnetworks */
1064b9c6e19dSShri Abhyankar   for(j=0; j < newDMnetwork->nsubnet; j++) {
1065b9c6e19dSShri Abhyankar     ierr = PetscCalloc1(newDMnetwork->subnet[j].nedge,&newDMnetwork->subnet[j].edges);CHKERRQ(ierr);
1066b9c6e19dSShri Abhyankar     ierr = PetscCalloc1(newDMnetwork->subnet[j].nvtx,&newDMnetwork->subnet[j].vertices);CHKERRQ(ierr);
1067b9c6e19dSShri Abhyankar     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
1068b9c6e19dSShri Abhyankar        These get updated when the vertices and edges are added. */
1069b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0;
1070b9c6e19dSShri Abhyankar   }
1071b9c6e19dSShri Abhyankar 
1072b9c6e19dSShri Abhyankar   /* Set the vertices and edges in each subnetwork */
1073b9c6e19dSShri Abhyankar   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1074b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
1075b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1076b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++]  = e;
1077b9c6e19dSShri Abhyankar   }
1078b9c6e19dSShri Abhyankar 
1079b9c6e19dSShri Abhyankar   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1080b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
1081b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1082b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++]  = v;
1083b9c6e19dSShri Abhyankar   }
1084b9c6e19dSShri Abhyankar 
108524121865SAdrian Maldonado   /* Destroy point SF */
108624121865SAdrian Maldonado   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
108724121865SAdrian Maldonado 
1088d3464fd4SAdrian Maldonado   ierr = DMDestroy(dm);CHKERRQ(ierr);
1089d3464fd4SAdrian Maldonado   *dm  = newDM;
10905f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10915f2c45f1SShri Abhyankar }
10925f2c45f1SShri Abhyankar 
109324121865SAdrian Maldonado /*@C
109424121865SAdrian Maldonado   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering.
109524121865SAdrian Maldonado 
109624121865SAdrian Maldonado   Input Parameters:
109724121865SAdrian Maldonado + masterSF - the original SF structure
109824121865SAdrian Maldonado - map      - a ISLocalToGlobal mapping that contains the subset of points
109924121865SAdrian Maldonado 
110024121865SAdrian Maldonado   Output Parameters:
110124121865SAdrian Maldonado . subSF    - a subset of the masterSF for the desired subset.
110224121865SAdrian Maldonado */
110324121865SAdrian Maldonado 
110424121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) {
110524121865SAdrian Maldonado 
110624121865SAdrian Maldonado   PetscErrorCode        ierr;
110724121865SAdrian Maldonado   PetscInt              nroots, nleaves, *ilocal_sub;
110824121865SAdrian Maldonado   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
110924121865SAdrian Maldonado   PetscInt              *local_points, *remote_points;
111024121865SAdrian Maldonado   PetscSFNode           *iremote_sub;
111124121865SAdrian Maldonado   const PetscInt        *ilocal;
111224121865SAdrian Maldonado   const PetscSFNode     *iremote;
111324121865SAdrian Maldonado 
111424121865SAdrian Maldonado   PetscFunctionBegin;
111524121865SAdrian Maldonado   ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
111624121865SAdrian Maldonado 
111724121865SAdrian Maldonado   /* Look for leaves that pertain to the subset of points. Get the local ordering */
111824121865SAdrian Maldonado   ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr);
111924121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr);
112024121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
112124121865SAdrian Maldonado     if (ilocal_map[i] != -1) nleaves_sub += 1;
112224121865SAdrian Maldonado   }
112324121865SAdrian Maldonado   /* Re-number ilocal with subset numbering. Need information from roots */
112424121865SAdrian Maldonado   ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr);
112524121865SAdrian Maldonado   for (i = 0; i < nroots; i++) local_points[i] = i;
112624121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr);
112724121865SAdrian Maldonado   ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
112824121865SAdrian Maldonado   ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
112924121865SAdrian Maldonado   /* Fill up graph using local (that is, local to the subset) numbering. */
11304b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr);
11314b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr);
113224121865SAdrian Maldonado   nleaves_sub = 0;
113324121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
113424121865SAdrian Maldonado     if (ilocal_map[i] != -1) {
113524121865SAdrian Maldonado       ilocal_sub[nleaves_sub] = ilocal_map[i];
11364b70a8deSAdrian Maldonado       iremote_sub[nleaves_sub].rank = iremote[i].rank;
113724121865SAdrian Maldonado       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
113824121865SAdrian Maldonado       nleaves_sub += 1;
113924121865SAdrian Maldonado     }
114024121865SAdrian Maldonado   }
114124121865SAdrian Maldonado   ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr);
114224121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr);
114324121865SAdrian Maldonado 
114424121865SAdrian Maldonado   /* Create new subSF */
114524121865SAdrian Maldonado   ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr);
114624121865SAdrian Maldonado   ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr);
11474b70a8deSAdrian Maldonado   ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr);
114824121865SAdrian Maldonado   ierr = PetscFree(ilocal_map);CHKERRQ(ierr);
11494b70a8deSAdrian Maldonado   ierr = PetscFree(iremote_sub);CHKERRQ(ierr);
115024121865SAdrian Maldonado   PetscFunctionReturn(0);
115124121865SAdrian Maldonado }
115224121865SAdrian Maldonado 
11535f2c45f1SShri Abhyankar /*@C
11545f2c45f1SShri Abhyankar   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
11555f2c45f1SShri Abhyankar 
11565f2c45f1SShri Abhyankar   Not Collective
11575f2c45f1SShri Abhyankar 
11585f2c45f1SShri Abhyankar   Input Parameters:
11595f2c45f1SShri Abhyankar + dm - The DMNetwork object
11605f2c45f1SShri Abhyankar - p  - the vertex point
11615f2c45f1SShri Abhyankar 
11625f2c45f1SShri Abhyankar   Output Paramters:
11635f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point
11645f2c45f1SShri Abhyankar - edges  - List of edge points
11655f2c45f1SShri Abhyankar 
11665f2c45f1SShri Abhyankar   Level: intermediate
11675f2c45f1SShri Abhyankar 
11685f2c45f1SShri Abhyankar   Fortran Notes:
11695f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
11705f2c45f1SShri Abhyankar   include petsc.h90 in your code.
11715f2c45f1SShri Abhyankar 
1172d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices
11735f2c45f1SShri Abhyankar @*/
11745f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
11755f2c45f1SShri Abhyankar {
11765f2c45f1SShri Abhyankar   PetscErrorCode ierr;
11775f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
11785f2c45f1SShri Abhyankar 
11795f2c45f1SShri Abhyankar   PetscFunctionBegin;
11805f2c45f1SShri Abhyankar   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
11815f2c45f1SShri Abhyankar   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
11825f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
11835f2c45f1SShri Abhyankar }
11845f2c45f1SShri Abhyankar 
11855f2c45f1SShri Abhyankar /*@C
1186d842c372SHong Zhang   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
11875f2c45f1SShri Abhyankar 
11885f2c45f1SShri Abhyankar   Not Collective
11895f2c45f1SShri Abhyankar 
11905f2c45f1SShri Abhyankar   Input Parameters:
11915f2c45f1SShri Abhyankar + dm - The DMNetwork object
11925f2c45f1SShri Abhyankar - p  - the edge point
11935f2c45f1SShri Abhyankar 
11945f2c45f1SShri Abhyankar   Output Paramters:
11955f2c45f1SShri Abhyankar . vertices  - vertices connected to this edge
11965f2c45f1SShri Abhyankar 
11975f2c45f1SShri Abhyankar   Level: intermediate
11985f2c45f1SShri Abhyankar 
11995f2c45f1SShri Abhyankar   Fortran Notes:
12005f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
12015f2c45f1SShri Abhyankar   include petsc.h90 in your code.
12025f2c45f1SShri Abhyankar 
12035f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
12045f2c45f1SShri Abhyankar @*/
1205d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
12065f2c45f1SShri Abhyankar {
12075f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12085f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
12095f2c45f1SShri Abhyankar 
12105f2c45f1SShri Abhyankar   PetscFunctionBegin;
12115f2c45f1SShri Abhyankar   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
12125f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12135f2c45f1SShri Abhyankar }
12145f2c45f1SShri Abhyankar 
12155f2c45f1SShri Abhyankar /*@
12165f2c45f1SShri Abhyankar   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
12175f2c45f1SShri Abhyankar 
12185f2c45f1SShri Abhyankar   Not Collective
12195f2c45f1SShri Abhyankar 
12205f2c45f1SShri Abhyankar   Input Parameters:
12215f2c45f1SShri Abhyankar + dm - The DMNetwork object
12225f2c45f1SShri Abhyankar . p  - the vertex point
12235f2c45f1SShri Abhyankar 
12245f2c45f1SShri Abhyankar   Output Parameter:
12255f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point
12265f2c45f1SShri Abhyankar 
12275f2c45f1SShri Abhyankar   Level: intermediate
12285f2c45f1SShri Abhyankar 
1229d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange
12305f2c45f1SShri Abhyankar @*/
12315f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
12325f2c45f1SShri Abhyankar {
12335f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12345f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
12355f2c45f1SShri Abhyankar   PetscInt       offsetg;
12365f2c45f1SShri Abhyankar   PetscSection   sectiong;
12375f2c45f1SShri Abhyankar 
12385f2c45f1SShri Abhyankar   PetscFunctionBegin;
12395f2c45f1SShri Abhyankar   *isghost = PETSC_FALSE;
12405f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
12415f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
12425f2c45f1SShri Abhyankar   if (offsetg < 0) *isghost = PETSC_TRUE;
12435f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12445f2c45f1SShri Abhyankar }
12455f2c45f1SShri Abhyankar 
12465f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm)
12475f2c45f1SShri Abhyankar {
12485f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12495f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
12505f2c45f1SShri Abhyankar 
12515f2c45f1SShri Abhyankar   PetscFunctionBegin;
12525f2c45f1SShri Abhyankar   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
12535f2c45f1SShri Abhyankar   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
12545f2c45f1SShri Abhyankar 
12555f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
12565f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
12575f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12585f2c45f1SShri Abhyankar }
12595f2c45f1SShri Abhyankar 
12601ad426b7SHong Zhang /*@
126117df6e9eSHong Zhang     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
12621ad426b7SHong Zhang                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
12631ad426b7SHong Zhang 
12641ad426b7SHong Zhang     Collective
12651ad426b7SHong Zhang 
12661ad426b7SHong Zhang     Input Parameters:
126783b2e829SHong Zhang +   dm - The DMNetwork object
126883b2e829SHong Zhang .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
126983b2e829SHong Zhang -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
12701ad426b7SHong Zhang 
12711ad426b7SHong Zhang     Level: intermediate
12721ad426b7SHong Zhang 
12731ad426b7SHong Zhang @*/
127483b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
12751ad426b7SHong Zhang {
12761ad426b7SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
12778675203cSHong Zhang   PetscErrorCode ierr;
12781ad426b7SHong Zhang 
12791ad426b7SHong Zhang   PetscFunctionBegin;
128083b2e829SHong Zhang   network->userEdgeJacobian   = eflg;
128183b2e829SHong Zhang   network->userVertexJacobian = vflg;
12828675203cSHong Zhang 
12838675203cSHong Zhang   if (eflg && !network->Je) {
12848675203cSHong Zhang     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
12858675203cSHong Zhang   }
12868675203cSHong Zhang 
12878675203cSHong Zhang   if (vflg && !network->Jv) {
12888675203cSHong Zhang     PetscInt       i,*vptr,nedges,vStart=network->vStart;
12898675203cSHong Zhang     PetscInt       nVertices = network->nVertices,nedges_total;
12908675203cSHong Zhang     const PetscInt *edges;
12918675203cSHong Zhang 
12928675203cSHong Zhang     /* count nvertex_total */
12938675203cSHong Zhang     nedges_total = 0;
12948675203cSHong Zhang     ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr);
12958675203cSHong Zhang 
12968675203cSHong Zhang     vptr[0] = 0;
12978675203cSHong Zhang     for (i=0; i<nVertices; i++) {
12988675203cSHong Zhang       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
12998675203cSHong Zhang       nedges_total += nedges;
13008675203cSHong Zhang       vptr[i+1] = vptr[i] + 2*nedges + 1;
13018675203cSHong Zhang     }
13028675203cSHong Zhang 
13038675203cSHong Zhang     ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr);
13048675203cSHong Zhang     network->Jvptr = vptr;
13058675203cSHong Zhang   }
13061ad426b7SHong Zhang   PetscFunctionReturn(0);
13071ad426b7SHong Zhang }
13081ad426b7SHong Zhang 
13091ad426b7SHong Zhang /*@
131083b2e829SHong Zhang     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
131183b2e829SHong Zhang 
131283b2e829SHong Zhang     Not Collective
131383b2e829SHong Zhang 
131483b2e829SHong Zhang     Input Parameters:
131583b2e829SHong Zhang +   dm - The DMNetwork object
131683b2e829SHong Zhang .   p  - the edge point
13173e97b6e8SHong Zhang -   J - array (size = 3) of Jacobian submatrices for this edge point:
13183e97b6e8SHong Zhang         J[0]: this edge
1319d842c372SHong Zhang         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()
132083b2e829SHong Zhang 
132183b2e829SHong Zhang     Level: intermediate
132283b2e829SHong Zhang 
132383b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix
132483b2e829SHong Zhang @*/
132583b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
132683b2e829SHong Zhang {
132783b2e829SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
132883b2e829SHong Zhang 
132983b2e829SHong Zhang   PetscFunctionBegin;
13308675203cSHong Zhang   if (!network->Je) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
13318675203cSHong Zhang 
13328675203cSHong Zhang   if (J) {
1333883e35e8SHong Zhang     network->Je[3*p]   = J[0];
1334883e35e8SHong Zhang     network->Je[3*p+1] = J[1];
1335883e35e8SHong Zhang     network->Je[3*p+2] = J[2];
13368675203cSHong Zhang   }
133783b2e829SHong Zhang   PetscFunctionReturn(0);
133883b2e829SHong Zhang }
133983b2e829SHong Zhang 
134083b2e829SHong Zhang /*@
134176ddfea5SHong Zhang     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
13421ad426b7SHong Zhang 
13431ad426b7SHong Zhang     Not Collective
13441ad426b7SHong Zhang 
13451ad426b7SHong Zhang     Input Parameters:
13461ad426b7SHong Zhang +   dm - The DMNetwork object
13471ad426b7SHong Zhang .   p  - the vertex point
13483e97b6e8SHong Zhang -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
13493e97b6e8SHong Zhang         J[0]:       this vertex
13503e97b6e8SHong Zhang         J[1+2*i]:   i-th supporting edge
13513e97b6e8SHong Zhang         J[1+2*i+1]: i-th connected vertex
13521ad426b7SHong Zhang 
13531ad426b7SHong Zhang     Level: intermediate
13541ad426b7SHong Zhang 
135583b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix
13561ad426b7SHong Zhang @*/
1357883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
13585f2c45f1SShri Abhyankar {
13595f2c45f1SShri Abhyankar   PetscErrorCode ierr;
13605f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
13618675203cSHong Zhang   PetscInt       i,*vptr,nedges,vStart=network->vStart;
1362883e35e8SHong Zhang   const PetscInt *edges;
13635f2c45f1SShri Abhyankar 
13645f2c45f1SShri Abhyankar   PetscFunctionBegin;
13658675203cSHong Zhang   if (!network->Jv) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
1366883e35e8SHong Zhang 
13678675203cSHong Zhang   if (J) {
1368883e35e8SHong Zhang     vptr = network->Jvptr;
13693e97b6e8SHong Zhang     network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
13703e97b6e8SHong Zhang 
13713e97b6e8SHong Zhang     /* Set Jacobian for each supporting edge and connected vertex */
1372883e35e8SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
1373883e35e8SHong Zhang     for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
13748675203cSHong Zhang   }
1375883e35e8SHong Zhang   PetscFunctionReturn(0);
1376883e35e8SHong Zhang }
1377883e35e8SHong Zhang 
1378e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
13795cf7da58SHong Zhang {
13805cf7da58SHong Zhang   PetscErrorCode ierr;
13815cf7da58SHong Zhang   PetscInt       j;
13825cf7da58SHong Zhang   PetscScalar    val=(PetscScalar)ncols;
13835cf7da58SHong Zhang 
13845cf7da58SHong Zhang   PetscFunctionBegin;
13855cf7da58SHong Zhang   if (!ghost) {
13865cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
13875cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
13885cf7da58SHong Zhang     }
13895cf7da58SHong Zhang   } else {
13905cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
13915cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
13925cf7da58SHong Zhang     }
13935cf7da58SHong Zhang   }
13945cf7da58SHong Zhang   PetscFunctionReturn(0);
13955cf7da58SHong Zhang }
13965cf7da58SHong Zhang 
1397e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
13985cf7da58SHong Zhang {
13995cf7da58SHong Zhang   PetscErrorCode ierr;
14005cf7da58SHong Zhang   PetscInt       j,ncols_u;
14015cf7da58SHong Zhang   PetscScalar    val;
14025cf7da58SHong Zhang 
14035cf7da58SHong Zhang   PetscFunctionBegin;
14045cf7da58SHong Zhang   if (!ghost) {
14055cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
14065cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
14075cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
14085cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
14095cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
14105cf7da58SHong Zhang     }
14115cf7da58SHong Zhang   } else {
14125cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
14135cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
14145cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
14155cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
14165cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
14175cf7da58SHong Zhang     }
14185cf7da58SHong Zhang   }
14195cf7da58SHong Zhang   PetscFunctionReturn(0);
14205cf7da58SHong Zhang }
14215cf7da58SHong Zhang 
1422e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
14235cf7da58SHong Zhang {
14245cf7da58SHong Zhang   PetscErrorCode ierr;
14255cf7da58SHong Zhang 
14265cf7da58SHong Zhang   PetscFunctionBegin;
14275cf7da58SHong Zhang   if (Ju) {
14285cf7da58SHong Zhang     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
14295cf7da58SHong Zhang   } else {
14305cf7da58SHong Zhang     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
14315cf7da58SHong Zhang   }
14325cf7da58SHong Zhang   PetscFunctionReturn(0);
14335cf7da58SHong Zhang }
14345cf7da58SHong Zhang 
1435e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1436883e35e8SHong Zhang {
1437883e35e8SHong Zhang   PetscErrorCode ierr;
1438883e35e8SHong Zhang   PetscInt       j,*cols;
1439883e35e8SHong Zhang   PetscScalar    *zeros;
1440883e35e8SHong Zhang 
1441883e35e8SHong Zhang   PetscFunctionBegin;
1442883e35e8SHong Zhang   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
1443883e35e8SHong Zhang   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1444883e35e8SHong Zhang   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
1445883e35e8SHong Zhang   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
14461ad426b7SHong Zhang   PetscFunctionReturn(0);
14471ad426b7SHong Zhang }
1448a4e85ca8SHong Zhang 
1449e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
14503e97b6e8SHong Zhang {
14513e97b6e8SHong Zhang   PetscErrorCode ierr;
14523e97b6e8SHong Zhang   PetscInt       j,M,N,row,col,ncols_u;
14533e97b6e8SHong Zhang   const PetscInt *cols;
14543e97b6e8SHong Zhang   PetscScalar    zero=0.0;
14553e97b6e8SHong Zhang 
14563e97b6e8SHong Zhang   PetscFunctionBegin;
14573e97b6e8SHong Zhang   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
14583e97b6e8SHong 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);
14593e97b6e8SHong Zhang 
14603e97b6e8SHong Zhang   for (row=0; row<nrows; row++) {
14613e97b6e8SHong Zhang     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
14623e97b6e8SHong Zhang     for (j=0; j<ncols_u; j++) {
14633e97b6e8SHong Zhang       col = cols[j] + cstart;
14643e97b6e8SHong Zhang       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
14653e97b6e8SHong Zhang     }
14663e97b6e8SHong Zhang     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
14673e97b6e8SHong Zhang   }
14683e97b6e8SHong Zhang   PetscFunctionReturn(0);
14693e97b6e8SHong Zhang }
14701ad426b7SHong Zhang 
1471e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1472a4e85ca8SHong Zhang {
1473a4e85ca8SHong Zhang   PetscErrorCode ierr;
1474f4431b8cSHong Zhang 
1475a4e85ca8SHong Zhang   PetscFunctionBegin;
1476a4e85ca8SHong Zhang   if (Ju) {
1477a4e85ca8SHong Zhang     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1478a4e85ca8SHong Zhang   } else {
1479a4e85ca8SHong Zhang     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1480a4e85ca8SHong Zhang   }
1481a4e85ca8SHong Zhang   PetscFunctionReturn(0);
1482a4e85ca8SHong Zhang }
1483a4e85ca8SHong Zhang 
148424121865SAdrian 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.
148524121865SAdrian Maldonado */
148624121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
148724121865SAdrian Maldonado {
148824121865SAdrian Maldonado   PetscErrorCode ierr;
148924121865SAdrian Maldonado   PetscInt       i, size, dof;
149024121865SAdrian Maldonado   PetscInt       *glob2loc;
149124121865SAdrian Maldonado 
149224121865SAdrian Maldonado   PetscFunctionBegin;
149324121865SAdrian Maldonado   ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr);
149424121865SAdrian Maldonado   ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr);
149524121865SAdrian Maldonado 
149624121865SAdrian Maldonado   for (i = 0; i < size; i++) {
149724121865SAdrian Maldonado     ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr);
149824121865SAdrian Maldonado     dof = (dof >= 0) ? dof : -(dof + 1);
149924121865SAdrian Maldonado     glob2loc[i] = dof;
150024121865SAdrian Maldonado   }
150124121865SAdrian Maldonado 
150224121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
150324121865SAdrian Maldonado #if 0
150424121865SAdrian Maldonado   ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
150524121865SAdrian Maldonado #endif
150624121865SAdrian Maldonado   PetscFunctionReturn(0);
150724121865SAdrian Maldonado }
150824121865SAdrian Maldonado 
150901ad2aeeSHong Zhang #include <petsc/private/matimpl.h>
15101ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
15111ad426b7SHong Zhang {
15121ad426b7SHong Zhang   PetscErrorCode ierr;
151324121865SAdrian Maldonado   PetscMPIInt    rank, size;
15141ad426b7SHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
1515a4e85ca8SHong Zhang   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1516840c2264SHong Zhang   PetscInt       cstart,ncols,j,e,v;
151724121865SAdrian Maldonado   PetscBool      ghost,ghost_vc,ghost2,isNest;
1518a4e85ca8SHong Zhang   Mat            Juser;
1519bfbc38dcSHong Zhang   PetscSection   sectionGlobal;
1520447d78afSSatish Balay   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
1521a4e85ca8SHong Zhang   const PetscInt *edges,*cone;
15225cf7da58SHong Zhang   MPI_Comm       comm;
152324121865SAdrian Maldonado   MatType        mtype;
15245cf7da58SHong Zhang   Vec            vd_nz,vo_nz;
15255cf7da58SHong Zhang   PetscInt       *dnnz,*onnz;
15265cf7da58SHong Zhang   PetscScalar    *vdnz,*vonz;
15271ad426b7SHong Zhang 
15281ad426b7SHong Zhang   PetscFunctionBegin;
152924121865SAdrian Maldonado   mtype = dm->mattype;
153024121865SAdrian Maldonado   ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr);
153124121865SAdrian Maldonado 
153224121865SAdrian Maldonado   if (isNest) {
15330731d606SHong Zhang     /* ierr = DMCreateMatrix_Network_Nest(); */
153424121865SAdrian Maldonado     PetscInt   eDof, vDof;
153524121865SAdrian Maldonado     Mat        j11, j12, j21, j22, bA[2][2];
153624121865SAdrian Maldonado     ISLocalToGlobalMapping eISMap, vISMap;
153724121865SAdrian Maldonado 
153824121865SAdrian Maldonado     ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
153924121865SAdrian Maldonado     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
154024121865SAdrian Maldonado     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
154124121865SAdrian Maldonado 
154224121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr);
154324121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr);
154424121865SAdrian Maldonado 
154501ad2aeeSHong Zhang     ierr = MatCreate(comm, &j11);CHKERRQ(ierr);
154624121865SAdrian Maldonado     ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
154724121865SAdrian Maldonado     ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr);
154824121865SAdrian Maldonado 
154901ad2aeeSHong Zhang     ierr = MatCreate(comm, &j12);CHKERRQ(ierr);
155024121865SAdrian Maldonado     ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr);
155124121865SAdrian Maldonado     ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr);
155224121865SAdrian Maldonado 
155301ad2aeeSHong Zhang     ierr = MatCreate(comm, &j21);CHKERRQ(ierr);
155424121865SAdrian Maldonado     ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
155524121865SAdrian Maldonado     ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr);
155624121865SAdrian Maldonado 
155701ad2aeeSHong Zhang     ierr = MatCreate(comm, &j22);CHKERRQ(ierr);
155824121865SAdrian Maldonado     ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
155924121865SAdrian Maldonado     ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr);
156024121865SAdrian Maldonado 
15613f6a6bdaSHong Zhang     bA[0][0] = j11;
15623f6a6bdaSHong Zhang     bA[0][1] = j12;
15633f6a6bdaSHong Zhang     bA[1][0] = j21;
15643f6a6bdaSHong Zhang     bA[1][1] = j22;
156524121865SAdrian Maldonado 
156624121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr);
156724121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr);
156824121865SAdrian Maldonado 
156924121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr);
157024121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr);
157124121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr);
157224121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr);
157324121865SAdrian Maldonado 
157424121865SAdrian Maldonado     ierr = MatSetUp(j11);CHKERRQ(ierr);
157524121865SAdrian Maldonado     ierr = MatSetUp(j12);CHKERRQ(ierr);
157624121865SAdrian Maldonado     ierr = MatSetUp(j21);CHKERRQ(ierr);
157724121865SAdrian Maldonado     ierr = MatSetUp(j22);CHKERRQ(ierr);
157824121865SAdrian Maldonado 
157901ad2aeeSHong Zhang     ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr);
158024121865SAdrian Maldonado     ierr = MatSetUp(*J);CHKERRQ(ierr);
158124121865SAdrian Maldonado     ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr);
158224121865SAdrian Maldonado     ierr = MatDestroy(&j11);CHKERRQ(ierr);
158324121865SAdrian Maldonado     ierr = MatDestroy(&j12);CHKERRQ(ierr);
158424121865SAdrian Maldonado     ierr = MatDestroy(&j21);CHKERRQ(ierr);
158524121865SAdrian Maldonado     ierr = MatDestroy(&j22);CHKERRQ(ierr);
158624121865SAdrian Maldonado 
158724121865SAdrian Maldonado     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158824121865SAdrian Maldonado     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158924121865SAdrian Maldonado 
159024121865SAdrian Maldonado     /* Free structures */
159124121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr);
159224121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr);
159324121865SAdrian Maldonado 
159424121865SAdrian Maldonado     PetscFunctionReturn(0);
159524121865SAdrian Maldonado   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1596a4e85ca8SHong Zhang     /* user does not provide Jacobian blocks */
1597bfbc38dcSHong Zhang     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1598bfbc38dcSHong Zhang     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
15991ad426b7SHong Zhang     PetscFunctionReturn(0);
16001ad426b7SHong Zhang   }
16011ad426b7SHong Zhang 
1602bfbc38dcSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
16032a945128SHong Zhang   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1604bfbc38dcSHong Zhang   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1605bfbc38dcSHong Zhang   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
16062a945128SHong Zhang 
16072a945128SHong Zhang   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
16082a945128SHong Zhang   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
160989898e50SHong Zhang 
161089898e50SHong Zhang   /* (1) Set matrix preallocation */
161189898e50SHong Zhang   /*------------------------------*/
1612840c2264SHong Zhang   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1613840c2264SHong Zhang   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1614840c2264SHong Zhang   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1615840c2264SHong Zhang   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1616840c2264SHong Zhang   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1617840c2264SHong Zhang   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1618840c2264SHong Zhang 
161989898e50SHong Zhang   /* Set preallocation for edges */
162089898e50SHong Zhang   /*-----------------------------*/
1621840c2264SHong Zhang   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1622840c2264SHong Zhang 
1623bdcb62a2SHong Zhang   ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr);
1624840c2264SHong Zhang   for (e=eStart; e<eEnd; e++) {
1625840c2264SHong Zhang     /* Get row indices */
1626840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1627840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1628840c2264SHong Zhang     if (nrows) {
1629840c2264SHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1630840c2264SHong Zhang 
16315cf7da58SHong Zhang       /* Set preallocation for conntected vertices */
1632d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1633840c2264SHong Zhang       for (v=0; v<2; v++) {
1634840c2264SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1635840c2264SHong Zhang 
16368675203cSHong Zhang         if (network->Je) {
1637840c2264SHong Zhang           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
16388675203cSHong Zhang         } else Juser = NULL;
1639840c2264SHong Zhang         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
16405cf7da58SHong Zhang         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1641840c2264SHong Zhang       }
1642840c2264SHong Zhang 
164389898e50SHong Zhang       /* Set preallocation for edge self */
1644840c2264SHong Zhang       cstart = rstart;
16458675203cSHong Zhang       if (network->Je) {
1646840c2264SHong Zhang         Juser = network->Je[3*e]; /* Jacobian(e,e) */
16478675203cSHong Zhang       } else Juser = NULL;
16485cf7da58SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1649840c2264SHong Zhang     }
1650840c2264SHong Zhang   }
1651840c2264SHong Zhang 
165289898e50SHong Zhang   /* Set preallocation for vertices */
165389898e50SHong Zhang   /*--------------------------------*/
1654840c2264SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
16558675203cSHong Zhang   if (vEnd - vStart) vptr = network->Jvptr;
1656840c2264SHong Zhang 
1657840c2264SHong Zhang   for (v=vStart; v<vEnd; v++) {
1658840c2264SHong Zhang     /* Get row indices */
1659840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1660840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1661840c2264SHong Zhang     if (!nrows) continue;
1662840c2264SHong Zhang 
1663bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1664bdcb62a2SHong Zhang     if (ghost) {
1665bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1666bdcb62a2SHong Zhang     } else {
1667bdcb62a2SHong Zhang       rows_v = rows;
1668bdcb62a2SHong Zhang     }
1669bdcb62a2SHong Zhang 
1670bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1671840c2264SHong Zhang 
1672840c2264SHong Zhang     /* Get supporting edges and connected vertices */
1673840c2264SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1674840c2264SHong Zhang 
1675840c2264SHong Zhang     for (e=0; e<nedges; e++) {
1676840c2264SHong Zhang       /* Supporting edges */
1677840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1678840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1679840c2264SHong Zhang 
16808675203cSHong Zhang       if (network->Jv) {
1681840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
16828675203cSHong Zhang       } else Juser = NULL;
1683bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1684840c2264SHong Zhang 
1685840c2264SHong Zhang       /* Connected vertices */
1686d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
1687840c2264SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
1688840c2264SHong Zhang       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1689840c2264SHong Zhang 
1690840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1691840c2264SHong Zhang 
16928675203cSHong Zhang       if (network->Jv) {
1693840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
16948675203cSHong Zhang       } else Juser = NULL;
1695e102a522SHong Zhang       if (ghost_vc||ghost) {
1696e102a522SHong Zhang         ghost2 = PETSC_TRUE;
1697e102a522SHong Zhang       } else {
1698e102a522SHong Zhang         ghost2 = PETSC_FALSE;
1699e102a522SHong Zhang       }
1700e102a522SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr);
1701840c2264SHong Zhang     }
1702840c2264SHong Zhang 
170389898e50SHong Zhang     /* Set preallocation for vertex self */
1704840c2264SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1705840c2264SHong Zhang     if (!ghost) {
1706840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
17078675203cSHong Zhang       if (network->Jv) {
1708840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
17098675203cSHong Zhang       } else Juser = NULL;
1710bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1711840c2264SHong Zhang     }
1712bdcb62a2SHong Zhang     if (ghost) {
1713bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1714bdcb62a2SHong Zhang     }
1715840c2264SHong Zhang   }
1716840c2264SHong Zhang 
1717840c2264SHong Zhang   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1718840c2264SHong Zhang   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
17195cf7da58SHong Zhang 
17205cf7da58SHong Zhang   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
17215cf7da58SHong Zhang 
17225cf7da58SHong Zhang   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1723840c2264SHong Zhang   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1724840c2264SHong Zhang 
1725840c2264SHong Zhang   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1726840c2264SHong Zhang   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1727840c2264SHong Zhang   for (j=0; j<localSize; j++) {
1728e102a522SHong Zhang     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1729e102a522SHong Zhang     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1730840c2264SHong Zhang   }
1731840c2264SHong Zhang   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1732840c2264SHong Zhang   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1733840c2264SHong Zhang   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1734840c2264SHong Zhang   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1735840c2264SHong Zhang 
17365cf7da58SHong Zhang   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
17375cf7da58SHong Zhang   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
17385cf7da58SHong Zhang   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
17395cf7da58SHong Zhang 
17405cf7da58SHong Zhang   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
17415cf7da58SHong Zhang 
174289898e50SHong Zhang   /* (2) Set matrix entries for edges */
174389898e50SHong Zhang   /*----------------------------------*/
17441ad426b7SHong Zhang   for (e=eStart; e<eEnd; e++) {
1745bfbc38dcSHong Zhang     /* Get row indices */
17461ad426b7SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
174717df6e9eSHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
17484b976069SHong Zhang     if (nrows) {
174917df6e9eSHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
17501ad426b7SHong Zhang 
1751bfbc38dcSHong Zhang       /* Set matrix entries for conntected vertices */
1752d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1753bfbc38dcSHong Zhang       for (v=0; v<2; v++) {
1754bfbc38dcSHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1755883e35e8SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
17563e97b6e8SHong Zhang 
17578675203cSHong Zhang         if (network->Je) {
1758a4e85ca8SHong Zhang           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
17598675203cSHong Zhang         } else Juser = NULL;
1760a4e85ca8SHong Zhang         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1761bfbc38dcSHong Zhang       }
176217df6e9eSHong Zhang 
1763bfbc38dcSHong Zhang       /* Set matrix entries for edge self */
17643e97b6e8SHong Zhang       cstart = rstart;
17658675203cSHong Zhang       if (network->Je) {
1766a4e85ca8SHong Zhang         Juser = network->Je[3*e]; /* Jacobian(e,e) */
17678675203cSHong Zhang       } else Juser = NULL;
1768a4e85ca8SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
17691ad426b7SHong Zhang     }
17704b976069SHong Zhang   }
17711ad426b7SHong Zhang 
1772bfbc38dcSHong Zhang   /* Set matrix entries for vertices */
177383b2e829SHong Zhang   /*---------------------------------*/
17741ad426b7SHong Zhang   for (v=vStart; v<vEnd; v++) {
1775bfbc38dcSHong Zhang     /* Get row indices */
1776596e729fSHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1777596e729fSHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
17784b976069SHong Zhang     if (!nrows) continue;
1779596e729fSHong Zhang 
1780bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1781bdcb62a2SHong Zhang     if (ghost) {
1782bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1783bdcb62a2SHong Zhang     } else {
1784bdcb62a2SHong Zhang       rows_v = rows;
1785bdcb62a2SHong Zhang     }
1786bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1787596e729fSHong Zhang 
1788bfbc38dcSHong Zhang     /* Get supporting edges and connected vertices */
1789596e729fSHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1790596e729fSHong Zhang 
1791596e729fSHong Zhang     for (e=0; e<nedges; e++) {
1792bfbc38dcSHong Zhang       /* Supporting edges */
1793596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1794596e729fSHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1795596e729fSHong Zhang 
17968675203cSHong Zhang       if (network->Jv) {
1797a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
17988675203cSHong Zhang       } else Juser = NULL;
1799bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1800596e729fSHong Zhang 
1801bfbc38dcSHong Zhang       /* Connected vertices */
1802d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
18032a945128SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
18042a945128SHong Zhang 
180544aca652SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
180644aca652SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1807a4e85ca8SHong Zhang 
18088675203cSHong Zhang       if (network->Jv) {
1809a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
18108675203cSHong Zhang       } else Juser = NULL;
1811bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1812596e729fSHong Zhang     }
1813596e729fSHong Zhang 
1814bfbc38dcSHong Zhang     /* Set matrix entries for vertex self */
18151ad426b7SHong Zhang     if (!ghost) {
1816596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
18178675203cSHong Zhang       if (network->Jv) {
1818a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
18198675203cSHong Zhang       } else Juser = NULL;
1820bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr);
1821bdcb62a2SHong Zhang     }
1822bdcb62a2SHong Zhang     if (ghost) {
1823bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1824bdcb62a2SHong Zhang     }
18251ad426b7SHong Zhang   }
1826a4e85ca8SHong Zhang   ierr = PetscFree(rows);CHKERRQ(ierr);
1827bdcb62a2SHong Zhang 
18281ad426b7SHong Zhang   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
18291ad426b7SHong Zhang   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1830dd6f46cdSHong Zhang 
18315f2c45f1SShri Abhyankar   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
18325f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
18335f2c45f1SShri Abhyankar }
18345f2c45f1SShri Abhyankar 
18355f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm)
18365f2c45f1SShri Abhyankar {
18375f2c45f1SShri Abhyankar   PetscErrorCode ierr;
18385f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
18392727e31bSShri Abhyankar   PetscInt       j;
18405f2c45f1SShri Abhyankar 
18415f2c45f1SShri Abhyankar   PetscFunctionBegin;
18428415c774SShri Abhyankar   if (--network->refct > 0) PetscFunctionReturn(0);
184383b2e829SHong Zhang   if (network->Je) {
184483b2e829SHong Zhang     ierr = PetscFree(network->Je);CHKERRQ(ierr);
184583b2e829SHong Zhang   }
184683b2e829SHong Zhang   if (network->Jv) {
1847883e35e8SHong Zhang     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
184883b2e829SHong Zhang     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
18491ad426b7SHong Zhang   }
185013c2a604SAdrian Maldonado 
185113c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr);
185213c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr);
185313c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr);
185413c2a604SAdrian Maldonado   if (network->vertex.sf) {
185513c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr);
185613c2a604SAdrian Maldonado   }
185713c2a604SAdrian Maldonado   /* edge */
185813c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr);
185913c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr);
186013c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr);
186113c2a604SAdrian Maldonado   if (network->edge.sf) {
186213c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr);
186313c2a604SAdrian Maldonado   }
18645f2c45f1SShri Abhyankar   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
18655f2c45f1SShri Abhyankar   network->edges = NULL;
18665f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
18675f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
186883b2e829SHong Zhang 
18692727e31bSShri Abhyankar   for(j=0; j < network->nsubnet; j++) {
18702727e31bSShri Abhyankar     ierr = PetscFree(network->subnet[j].edges);CHKERRQ(ierr);
18712727e31bSShri Abhyankar     ierr = PetscFree(network->subnet[j].vertices);CHKERRQ(ierr);
18722727e31bSShri Abhyankar   }
1873e2aaf10cSShri Abhyankar   ierr = PetscFree(network->subnet);CHKERRQ(ierr);
18745f2c45f1SShri Abhyankar   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
18755f2c45f1SShri Abhyankar   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
18765f2c45f1SShri Abhyankar   ierr = PetscFree(network->header);CHKERRQ(ierr);
18775f2c45f1SShri Abhyankar   ierr = PetscFree(network);CHKERRQ(ierr);
18785f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
18795f2c45f1SShri Abhyankar }
18805f2c45f1SShri Abhyankar 
18815f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm,PetscViewer viewer)
18825f2c45f1SShri Abhyankar {
18835f2c45f1SShri Abhyankar   PetscErrorCode ierr;
18845f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
18855f2c45f1SShri Abhyankar 
18865f2c45f1SShri Abhyankar   PetscFunctionBegin;
18875f2c45f1SShri Abhyankar   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
18885f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
18895f2c45f1SShri Abhyankar }
18905f2c45f1SShri Abhyankar 
18915f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
18925f2c45f1SShri Abhyankar {
18935f2c45f1SShri Abhyankar   PetscErrorCode ierr;
18945f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
18955f2c45f1SShri Abhyankar 
18965f2c45f1SShri Abhyankar   PetscFunctionBegin;
18975f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
18985f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
18995f2c45f1SShri Abhyankar }
19005f2c45f1SShri Abhyankar 
19015f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
19025f2c45f1SShri Abhyankar {
19035f2c45f1SShri Abhyankar   PetscErrorCode ierr;
19045f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
19055f2c45f1SShri Abhyankar 
19065f2c45f1SShri Abhyankar   PetscFunctionBegin;
19075f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
19085f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19095f2c45f1SShri Abhyankar }
19105f2c45f1SShri Abhyankar 
19115f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
19125f2c45f1SShri Abhyankar {
19135f2c45f1SShri Abhyankar   PetscErrorCode ierr;
19145f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
19155f2c45f1SShri Abhyankar 
19165f2c45f1SShri Abhyankar   PetscFunctionBegin;
19175f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
19185f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19195f2c45f1SShri Abhyankar }
19205f2c45f1SShri Abhyankar 
19215f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
19225f2c45f1SShri Abhyankar {
19235f2c45f1SShri Abhyankar   PetscErrorCode ierr;
19245f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
19255f2c45f1SShri Abhyankar 
19265f2c45f1SShri Abhyankar   PetscFunctionBegin;
19275f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
19285f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19295f2c45f1SShri Abhyankar }
1930