xref: /petsc/src/dm/impls/network/network.c (revision 8675203cc7b097cc684504c0b41aa1e4ff7cb748)
1af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h>  /*I  "petscdmnetwork.h"  I*/
25f2c45f1SShri Abhyankar #include <petscdmplex.h>
35f2c45f1SShri Abhyankar #include <petscsf.h>
45f2c45f1SShri Abhyankar 
55f2c45f1SShri Abhyankar /*@
6556ed216SShri Abhyankar   DMNetworkGetPlex - Gets the Plex DM associated with this network DM
7556ed216SShri Abhyankar 
8556ed216SShri Abhyankar   Not collective
9556ed216SShri Abhyankar 
10556ed216SShri Abhyankar   Input Parameters:
11556ed216SShri Abhyankar + netdm - the dm object
12556ed216SShri Abhyankar - plexmdm - the plex dm object
13556ed216SShri Abhyankar 
14556ed216SShri Abhyankar   Level: Advanced
15556ed216SShri Abhyankar 
16556ed216SShri Abhyankar .seealso: DMNetworkCreate()
17556ed216SShri Abhyankar @*/
18556ed216SShri Abhyankar PetscErrorCode DMNetworkGetPlex(DM netdm, DM *plexdm)
19556ed216SShri Abhyankar {
20556ed216SShri Abhyankar   DM_Network     *network = (DM_Network*) netdm->data;
21556ed216SShri Abhyankar 
22556ed216SShri Abhyankar   PetscFunctionBegin;
23556ed216SShri Abhyankar   *plexdm = network->plex;
24556ed216SShri Abhyankar   PetscFunctionReturn(0);
25556ed216SShri Abhyankar }
26556ed216SShri Abhyankar 
27556ed216SShri Abhyankar /*@
28e2aaf10cSShri Abhyankar   DMNetworkSetSizes - Sets the number of subnetworks,local and global vertices and edges for each subnetwork.
295f2c45f1SShri Abhyankar 
305f2c45f1SShri Abhyankar   Collective on DM
315f2c45f1SShri Abhyankar 
325f2c45f1SShri Abhyankar   Input Parameters:
335f2c45f1SShri Abhyankar + dm - the dm object
34e2aaf10cSShri Abhyankar . Nsubnet - number of subnetworks
35e2aaf10cSShri Abhyankar . nV - number of local vertices for each subnetwork
36e2aaf10cSShri Abhyankar . nE - number of local edges for each subnetwork
37e2aaf10cSShri Abhyankar . NV - number of global vertices (or PETSC_DETERMINE) for each subnetwork
38e2aaf10cSShri Abhyankar - NE - number of global edges (or PETSC_DETERMINE) for each subnetwork
395f2c45f1SShri Abhyankar 
405f2c45f1SShri Abhyankar    Notes
415f2c45f1SShri Abhyankar    If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang.
425f2c45f1SShri Abhyankar 
435f2c45f1SShri Abhyankar    You cannot change the sizes once they have been set
445f2c45f1SShri Abhyankar 
451b266c99SBarry Smith    Level: intermediate
461b266c99SBarry Smith 
471b266c99SBarry Smith .seealso: DMNetworkCreate()
485f2c45f1SShri Abhyankar @*/
49e2aaf10cSShri Abhyankar PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt Nsubnet, PetscInt nV[], PetscInt nE[], PetscInt NV[], PetscInt NE[])
505f2c45f1SShri Abhyankar {
515f2c45f1SShri Abhyankar   PetscErrorCode ierr;
525f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
53e2aaf10cSShri Abhyankar   PetscInt       a[2],b[2],i;
545f2c45f1SShri Abhyankar 
555f2c45f1SShri Abhyankar   PetscFunctionBegin;
565f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
57e2aaf10cSShri Abhyankar   if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet);
58e2aaf10cSShri Abhyankar 
59e2aaf10cSShri Abhyankar   if(Nsubnet > 0) PetscValidLogicalCollectiveInt(dm,Nsubnet,2);
60e2aaf10cSShri Abhyankar   if(network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network");
61e2aaf10cSShri Abhyankar 
62e2aaf10cSShri Abhyankar   network->nsubnet = Nsubnet;
63e2aaf10cSShri Abhyankar   ierr = PetscCalloc1(Nsubnet,&network->subnet);CHKERRQ(ierr);
64e2aaf10cSShri Abhyankar   for(i=0; i < network->nsubnet; i++) {
65e2aaf10cSShri Abhyankar     if (NV[i] > 0) PetscValidLogicalCollectiveInt(dm,NV[i],5);
66e2aaf10cSShri Abhyankar     if (NE[i] > 0) PetscValidLogicalCollectiveInt(dm,NE[i],6);
67e2aaf10cSShri Abhyankar     if (NV[i] > 0 && nV[i] > NV[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local vertex size %D cannot be larger than global vertex size %D",i,nV[i],NV[i]);
68e2aaf10cSShri Abhyankar     if (NE[i] > 0 && nE[i] > NE[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local edge size %D cannot be larger than global edge size %D",i,nE[i],NE[i]);
69e2aaf10cSShri Abhyankar     a[0] = nV[i]; a[1] = nE[i];
70b2566f29SBarry Smith     ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
71e2aaf10cSShri Abhyankar     network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1];
72e2aaf10cSShri Abhyankar 
73e2aaf10cSShri Abhyankar     network->subnet[i].id = i;
74e2aaf10cSShri Abhyankar 
75e2aaf10cSShri Abhyankar     network->subnet[i].nvtx = nV[i];
76e2aaf10cSShri Abhyankar     network->subnet[i].vStart = network->nVertices;
77e2aaf10cSShri Abhyankar     network->subnet[i].vEnd   = network->subnet[i].vStart + network->subnet[i].Nvtx;
78e2aaf10cSShri Abhyankar     network->nVertices += network->subnet[i].nvtx;
79e2aaf10cSShri Abhyankar     network->NVertices += network->subnet[i].Nvtx;
80e2aaf10cSShri Abhyankar 
81e2aaf10cSShri Abhyankar     network->subnet[i].nedge = nE[i];
82e2aaf10cSShri Abhyankar     network->subnet[i].eStart = network->nEdges;
83e2aaf10cSShri Abhyankar     network->subnet[i].eEnd = network->subnet[i].eStart + network->subnet[i].Nedge;
84e2aaf10cSShri Abhyankar     network->nEdges += network->subnet[i].nedge;
85e2aaf10cSShri Abhyankar     network->NEdges += network->subnet[i].Nedge;
865f2c45f1SShri Abhyankar   }
875f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
885f2c45f1SShri Abhyankar }
895f2c45f1SShri Abhyankar 
905f2c45f1SShri Abhyankar /*@
915f2c45f1SShri Abhyankar   DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network
925f2c45f1SShri Abhyankar 
935f2c45f1SShri Abhyankar   Logically collective on DM
945f2c45f1SShri Abhyankar 
955f2c45f1SShri Abhyankar   Input Parameters:
96e2aaf10cSShri Abhyankar . edges - list of edges for each subnetwork
975f2c45f1SShri Abhyankar 
985f2c45f1SShri Abhyankar   Notes:
995f2c45f1SShri Abhyankar   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
1005f2c45f1SShri Abhyankar   not be destroyed before the call to DMNetworkLayoutSetUp
1015f2c45f1SShri Abhyankar 
1025f2c45f1SShri Abhyankar   Level: intermediate
1035f2c45f1SShri Abhyankar 
1045f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes
1055f2c45f1SShri Abhyankar @*/
106e2aaf10cSShri Abhyankar PetscErrorCode DMNetworkSetEdgeList(DM dm, int *edgelist[])
1075f2c45f1SShri Abhyankar {
1085f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network*) dm->data;
109e2aaf10cSShri Abhyankar   PetscInt   i;
1105f2c45f1SShri Abhyankar 
1115f2c45f1SShri Abhyankar   PetscFunctionBegin;
112e2aaf10cSShri Abhyankar   for(i=0; i < network->nsubnet; i++) {
113e2aaf10cSShri Abhyankar     network->subnet[i].edgelist = edgelist[i];
114e2aaf10cSShri Abhyankar   }
1155f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
1165f2c45f1SShri Abhyankar }
1175f2c45f1SShri Abhyankar 
1185f2c45f1SShri Abhyankar /*@
1195f2c45f1SShri Abhyankar   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
1205f2c45f1SShri Abhyankar 
1215f2c45f1SShri Abhyankar   Collective on DM
1225f2c45f1SShri Abhyankar 
1235f2c45f1SShri Abhyankar   Input Parameters
1245f2c45f1SShri Abhyankar . DM - the dmnetwork object
1255f2c45f1SShri Abhyankar 
1265f2c45f1SShri Abhyankar   Notes:
1275f2c45f1SShri Abhyankar   This routine should be called after the network sizes and edgelists have been provided. It creates
1285f2c45f1SShri Abhyankar   the bare layout of the network and sets up the network to begin insertion of components.
1295f2c45f1SShri Abhyankar 
1305f2c45f1SShri Abhyankar   All the components should be registered before calling this routine.
1315f2c45f1SShri Abhyankar 
1325f2c45f1SShri Abhyankar   Level: intermediate
1335f2c45f1SShri Abhyankar 
1345f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList
1355f2c45f1SShri Abhyankar @*/
1365f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm)
1375f2c45f1SShri Abhyankar {
1385f2c45f1SShri Abhyankar   PetscErrorCode ierr;
1395f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
1405f2c45f1SShri Abhyankar   PetscInt       dim = 1; /* One dimensional network */
1415f2c45f1SShri Abhyankar   PetscInt       numCorners=2;
1425f2c45f1SShri Abhyankar   PetscInt       spacedim=2;
1435f2c45f1SShri Abhyankar   double         *vertexcoords=NULL;
144e2aaf10cSShri Abhyankar   PetscInt       i,j;
1455f2c45f1SShri Abhyankar   PetscInt       ndata;
146e2aaf10cSShri Abhyankar   PetscInt       ctr=0;
1475f2c45f1SShri Abhyankar 
1485f2c45f1SShri Abhyankar   PetscFunctionBegin;
149e2aaf10cSShri Abhyankar 
1506fefedf4SHong Zhang   if (network->nVertices) {
1516fefedf4SHong Zhang     ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr);
1525f2c45f1SShri Abhyankar   }
153e2aaf10cSShri Abhyankar 
154e2aaf10cSShri Abhyankar   /* Create the edgelist for the network by concatenating edgelists of the subnetworks */
155e2aaf10cSShri Abhyankar   ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr);
156e2aaf10cSShri Abhyankar   for(i=0; i < network->nsubnet; i++) {
157e2aaf10cSShri Abhyankar     for(j = 0; j < network->subnet[i].nedge; j++) {
158e2aaf10cSShri Abhyankar       network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
159e2aaf10cSShri Abhyankar       network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
160e2aaf10cSShri Abhyankar       ctr++;
161e2aaf10cSShri Abhyankar     }
162e2aaf10cSShri Abhyankar   }
163e2aaf10cSShri Abhyankar 
164e2aaf10cSShri Abhyankar #if 0
165e2aaf10cSShri Abhyankar   for(i=0; i < network->nEdges; i++) {
166e2aaf10cSShri Abhyankar     ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr);
167e2aaf10cSShri Abhyankar   }
168e2aaf10cSShri Abhyankar #endif
169e2aaf10cSShri Abhyankar 
1706fefedf4SHong Zhang   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
1716fefedf4SHong Zhang   if (network->nVertices) {
1725f2c45f1SShri Abhyankar     ierr = PetscFree(vertexcoords);CHKERRQ(ierr);
1735f2c45f1SShri Abhyankar   }
174e2aaf10cSShri Abhyankar   ierr = PetscFree(network->edges);CHKERRQ(ierr);
175e2aaf10cSShri Abhyankar 
1765f2c45f1SShri Abhyankar   ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr);
1775f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr);
1785f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr);
1795f2c45f1SShri Abhyankar 
1805f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr);
1815f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr);
1825f2c45f1SShri Abhyankar   ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr);
1835f2c45f1SShri Abhyankar   ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr);
1845f2c45f1SShri Abhyankar 
1852727e31bSShri Abhyankar   /* Create vertices and edges array for the subnetworks */
1862727e31bSShri Abhyankar   for(j=0; j < network->nsubnet; j++) {
1872727e31bSShri Abhyankar     ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr);
1882727e31bSShri Abhyankar     ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr);
1892727e31bSShri Abhyankar     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
1902727e31bSShri Abhyankar        These get updated when the vertices and edges are added. */
1912727e31bSShri Abhyankar     network->subnet[j].nvtx = network->subnet[j].nedge = 0;
1922727e31bSShri Abhyankar   }
1932727e31bSShri Abhyankar 
1945f2c45f1SShri Abhyankar   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
1956caa05f4SBarry Smith   ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr);
196e2aaf10cSShri Abhyankar   for(i=network->eStart; i < network->eEnd; i++) {
197e2aaf10cSShri Abhyankar     network->header[i].index = i;   /* Global edge number */
198e2aaf10cSShri Abhyankar     for(j=0; j < network->nsubnet; j++) {
199e2aaf10cSShri Abhyankar       if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) {
200e2aaf10cSShri Abhyankar 	network->header[i].subnetid = j; /* Subnetwork id */
2012727e31bSShri Abhyankar 	network->subnet[j].edges[network->subnet[j].nedge++] = i;
202e2aaf10cSShri Abhyankar 	break;
203e2aaf10cSShri Abhyankar       }
2047b6afd5bSHong Zhang     }
2055f2c45f1SShri Abhyankar     network->header[i].ndata = 0;
2065f2c45f1SShri Abhyankar     ndata = network->header[i].ndata;
2075f2c45f1SShri Abhyankar     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
2085f2c45f1SShri Abhyankar     network->header[i].offset[ndata] = 0;
2095f2c45f1SShri Abhyankar   }
210e2aaf10cSShri Abhyankar 
211e2aaf10cSShri Abhyankar   for(i=network->vStart; i < network->vEnd; i++) {
212e2aaf10cSShri Abhyankar     network->header[i].index = i - network->vStart;
213e2aaf10cSShri Abhyankar     for(j=0; j < network->nsubnet; j++) {
214e2aaf10cSShri Abhyankar       if((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) {
215e2aaf10cSShri Abhyankar 	network->header[i].subnetid = j;
2162727e31bSShri Abhyankar 	network->subnet[j].vertices[network->subnet[j].nvtx++] = i;
217e2aaf10cSShri Abhyankar 	break;
218e2aaf10cSShri Abhyankar       }
219e2aaf10cSShri Abhyankar     }
220e2aaf10cSShri Abhyankar     network->header[i].ndata = 0;
221e2aaf10cSShri Abhyankar     ndata = network->header[i].ndata;
222e2aaf10cSShri Abhyankar     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
223e2aaf10cSShri Abhyankar     network->header[i].offset[ndata] = 0;
224e2aaf10cSShri Abhyankar   }
225e2aaf10cSShri Abhyankar 
226854ce69bSBarry Smith   ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr);
227e2aaf10cSShri Abhyankar 
2285f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2295f2c45f1SShri Abhyankar }
2305f2c45f1SShri Abhyankar 
23194ef8ddeSSatish Balay /*@C
2322727e31bSShri Abhyankar   DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork
2332727e31bSShri Abhyankar 
2342727e31bSShri Abhyankar   Input Parameters
2352727e31bSShri Abhyankar + dm   - the number object
2362727e31bSShri Abhyankar - id   - the ID (integer) of the subnetwork
2372727e31bSShri Abhyankar 
2382727e31bSShri Abhyankar   Output Parameters
2392727e31bSShri Abhyankar + nv    - number of vertices (local)
2402727e31bSShri Abhyankar . ne    - number of edges (local)
2412727e31bSShri Abhyankar . vtx   - local vertices for this subnetwork
2422727e31bSShri Abhyankar . edge  - local edges for this subnetwork
2432727e31bSShri Abhyankar 
2442727e31bSShri Abhyankar   Notes:
2452727e31bSShri Abhyankar   Cannot call this routine before DMNetworkLayoutSetup()
2462727e31bSShri Abhyankar 
2472727e31bSShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
2482727e31bSShri Abhyankar @*/
2492727e31bSShri Abhyankar PetscErrorCode DMNetworkGetSubnetworkInfo(DM netdm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge)
2502727e31bSShri Abhyankar {
2512727e31bSShri Abhyankar   DM_Network     *network = (DM_Network*) netdm->data;
2522727e31bSShri Abhyankar 
2532727e31bSShri Abhyankar   PetscFunctionBegin;
2542727e31bSShri Abhyankar   *nv = network->subnet[id].nvtx;
2552727e31bSShri Abhyankar   *ne = network->subnet[id].nedge;
2562727e31bSShri Abhyankar   *vtx = network->subnet[id].vertices;
2572727e31bSShri Abhyankar   *edge = network->subnet[id].edges;
2582727e31bSShri Abhyankar   PetscFunctionReturn(0);
2592727e31bSShri Abhyankar }
2602727e31bSShri Abhyankar 
2612727e31bSShri Abhyankar /*@C
2625f2c45f1SShri Abhyankar   DMNetworkRegisterComponent - Registers the network component
2635f2c45f1SShri Abhyankar 
2645f2c45f1SShri Abhyankar   Logically collective on DM
2655f2c45f1SShri Abhyankar 
2665f2c45f1SShri Abhyankar   Input Parameters
2675f2c45f1SShri Abhyankar + dm   - the network object
2685f2c45f1SShri Abhyankar . name - the component name
2695f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data
2705f2c45f1SShri Abhyankar 
2715f2c45f1SShri Abhyankar    Output Parameters
2725f2c45f1SShri Abhyankar .   key - an integer key that defines the component
2735f2c45f1SShri Abhyankar 
2745f2c45f1SShri Abhyankar    Notes
2755f2c45f1SShri Abhyankar    This routine should be called by all processors before calling DMNetworkLayoutSetup().
2765f2c45f1SShri Abhyankar 
2775f2c45f1SShri Abhyankar    Level: intermediate
2785f2c45f1SShri Abhyankar 
2795f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
2805f2c45f1SShri Abhyankar @*/
2815f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key)
2825f2c45f1SShri Abhyankar {
2835f2c45f1SShri Abhyankar   PetscErrorCode        ierr;
2845f2c45f1SShri Abhyankar   DM_Network            *network = (DM_Network*) dm->data;
2855f2c45f1SShri Abhyankar   DMNetworkComponent    *component=&network->component[network->ncomponent];
2865f2c45f1SShri Abhyankar   PetscBool             flg=PETSC_FALSE;
2875f2c45f1SShri Abhyankar   PetscInt              i;
2885f2c45f1SShri Abhyankar 
2895f2c45f1SShri Abhyankar   PetscFunctionBegin;
2905f2c45f1SShri Abhyankar 
2915f2c45f1SShri Abhyankar   for (i=0; i < network->ncomponent; i++) {
2925f2c45f1SShri Abhyankar     ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr);
2935f2c45f1SShri Abhyankar     if (flg) {
2945f2c45f1SShri Abhyankar       *key = i;
2955f2c45f1SShri Abhyankar       PetscFunctionReturn(0);
2965f2c45f1SShri Abhyankar     }
2975f2c45f1SShri Abhyankar   }
2985f2c45f1SShri Abhyankar 
2995f2c45f1SShri Abhyankar   ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr);
3005f2c45f1SShri Abhyankar   component->size = size/sizeof(DMNetworkComponentGenericDataType);
3015f2c45f1SShri Abhyankar   *key = network->ncomponent;
3025f2c45f1SShri Abhyankar   network->ncomponent++;
3035f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3045f2c45f1SShri Abhyankar }
3055f2c45f1SShri Abhyankar 
3065f2c45f1SShri Abhyankar /*@
3075f2c45f1SShri Abhyankar   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.
3085f2c45f1SShri Abhyankar 
3095f2c45f1SShri Abhyankar   Not Collective
3105f2c45f1SShri Abhyankar 
3115f2c45f1SShri Abhyankar   Input Parameters:
3125f2c45f1SShri Abhyankar + dm - The DMNetwork object
3135f2c45f1SShri Abhyankar 
3145f2c45f1SShri Abhyankar   Output Paramters:
3155f2c45f1SShri Abhyankar + vStart - The first vertex point
3165f2c45f1SShri Abhyankar - vEnd   - One beyond the last vertex point
3175f2c45f1SShri Abhyankar 
3185f2c45f1SShri Abhyankar   Level: intermediate
3195f2c45f1SShri Abhyankar 
3205f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange
3215f2c45f1SShri Abhyankar @*/
3225f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
3235f2c45f1SShri Abhyankar {
3245f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
3255f2c45f1SShri Abhyankar 
3265f2c45f1SShri Abhyankar   PetscFunctionBegin;
3275f2c45f1SShri Abhyankar   if (vStart) *vStart = network->vStart;
3285f2c45f1SShri Abhyankar   if (vEnd) *vEnd = network->vEnd;
3295f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3305f2c45f1SShri Abhyankar }
3315f2c45f1SShri Abhyankar 
3325f2c45f1SShri Abhyankar /*@
3335f2c45f1SShri Abhyankar   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.
3345f2c45f1SShri Abhyankar 
3355f2c45f1SShri Abhyankar   Not Collective
3365f2c45f1SShri Abhyankar 
3375f2c45f1SShri Abhyankar   Input Parameters:
3385f2c45f1SShri Abhyankar + dm - The DMNetwork object
3395f2c45f1SShri Abhyankar 
3405f2c45f1SShri Abhyankar   Output Paramters:
3415f2c45f1SShri Abhyankar + eStart - The first edge point
3425f2c45f1SShri Abhyankar - eEnd   - One beyond the last edge point
3435f2c45f1SShri Abhyankar 
3445f2c45f1SShri Abhyankar   Level: intermediate
3455f2c45f1SShri Abhyankar 
3465f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange
3475f2c45f1SShri Abhyankar @*/
3485f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
3495f2c45f1SShri Abhyankar {
3505f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
3515f2c45f1SShri Abhyankar 
3525f2c45f1SShri Abhyankar   PetscFunctionBegin;
3535f2c45f1SShri Abhyankar   if (eStart) *eStart = network->eStart;
3545f2c45f1SShri Abhyankar   if (eEnd) *eEnd = network->eEnd;
3555f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3565f2c45f1SShri Abhyankar }
3575f2c45f1SShri Abhyankar 
3587b6afd5bSHong Zhang /*@
359e85e6aecSHong Zhang   DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge.
3607b6afd5bSHong Zhang 
3617b6afd5bSHong Zhang   Not Collective
3627b6afd5bSHong Zhang 
3637b6afd5bSHong Zhang   Input Parameters:
3647b6afd5bSHong Zhang + dm - DMNetwork object
365e85e6aecSHong Zhang - p  - edge point
3667b6afd5bSHong Zhang 
3677b6afd5bSHong Zhang   Output Paramters:
368e85e6aecSHong Zhang . index - user global numbering for the edge
3697b6afd5bSHong Zhang 
3707b6afd5bSHong Zhang   Level: intermediate
3717b6afd5bSHong Zhang 
372e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalVertexIndex
3737b6afd5bSHong Zhang @*/
374e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
3757b6afd5bSHong Zhang {
3767b6afd5bSHong Zhang   PetscErrorCode    ierr;
3777b6afd5bSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
3787b6afd5bSHong Zhang   PetscInt          offsetp;
3797b6afd5bSHong Zhang   DMNetworkComponentHeader header;
3807b6afd5bSHong Zhang 
3817b6afd5bSHong Zhang   PetscFunctionBegin;
3827b6afd5bSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
3837b6afd5bSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
384e85e6aecSHong Zhang   *index = header->index;
3857b6afd5bSHong Zhang   PetscFunctionReturn(0);
3867b6afd5bSHong Zhang }
3877b6afd5bSHong Zhang 
3885f2c45f1SShri Abhyankar /*@
389e85e6aecSHong Zhang   DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex.
390e85e6aecSHong Zhang 
391e85e6aecSHong Zhang   Not Collective
392e85e6aecSHong Zhang 
393e85e6aecSHong Zhang   Input Parameters:
394e85e6aecSHong Zhang + dm - DMNetwork object
395e85e6aecSHong Zhang - p  - vertex point
396e85e6aecSHong Zhang 
397e85e6aecSHong Zhang   Output Paramters:
398e85e6aecSHong Zhang . index - user global numbering for the vertex
399e85e6aecSHong Zhang 
400e85e6aecSHong Zhang   Level: intermediate
401e85e6aecSHong Zhang 
402e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalEdgeIndex
403e85e6aecSHong Zhang @*/
404e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
405e85e6aecSHong Zhang {
406e85e6aecSHong Zhang   PetscErrorCode    ierr;
407e85e6aecSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
408e85e6aecSHong Zhang   PetscInt          offsetp;
409e85e6aecSHong Zhang   DMNetworkComponentHeader header;
410e85e6aecSHong Zhang 
411e85e6aecSHong Zhang   PetscFunctionBegin;
412e85e6aecSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
413e85e6aecSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
414e85e6aecSHong Zhang   *index = header->index;
415e85e6aecSHong Zhang   PetscFunctionReturn(0);
416e85e6aecSHong Zhang }
417e85e6aecSHong Zhang 
418e85e6aecSHong Zhang /*@
419325661f6SSatish Balay   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)
4205f2c45f1SShri Abhyankar 
4215f2c45f1SShri Abhyankar   Not Collective
4225f2c45f1SShri Abhyankar 
4235f2c45f1SShri Abhyankar   Input Parameters:
4245f2c45f1SShri Abhyankar + dm           - The DMNetwork object
4255f2c45f1SShri Abhyankar . p            - vertex/edge point
4265f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component
4275f2c45f1SShri Abhyankar - compvalue    - pointer to the data structure for the component
4285f2c45f1SShri Abhyankar 
4295f2c45f1SShri Abhyankar   Level: intermediate
4305f2c45f1SShri Abhyankar 
4315f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
4325f2c45f1SShri Abhyankar @*/
4335f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
4345f2c45f1SShri Abhyankar {
4355f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
43643a39a44SBarry Smith   DMNetworkComponent       *component = &network->component[componentkey];
4375f2c45f1SShri Abhyankar   DMNetworkComponentHeader header = &network->header[p];
4385f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue = &network->cvalue[p];
4395f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
4405f2c45f1SShri Abhyankar 
4415f2c45f1SShri Abhyankar   PetscFunctionBegin;
442fa58f0a9SHong Zhang   if (header->ndata == MAX_DATA_AT_POINT) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components at a point exceeds the max %D",MAX_DATA_AT_POINT);
443fa58f0a9SHong Zhang 
44443a39a44SBarry Smith   header->size[header->ndata] = component->size;
44543a39a44SBarry Smith   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
4465f2c45f1SShri Abhyankar   header->key[header->ndata] = componentkey;
4475f2c45f1SShri Abhyankar   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
4485f2c45f1SShri Abhyankar 
4495f2c45f1SShri Abhyankar   cvalue->data[header->ndata] = (void*)compvalue;
4505f2c45f1SShri Abhyankar   header->ndata++;
4515f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4525f2c45f1SShri Abhyankar }
4535f2c45f1SShri Abhyankar 
4545f2c45f1SShri Abhyankar /*@
4555f2c45f1SShri Abhyankar   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
4565f2c45f1SShri Abhyankar 
4575f2c45f1SShri Abhyankar   Not Collective
4585f2c45f1SShri Abhyankar 
4595f2c45f1SShri Abhyankar   Input Parameters:
4605f2c45f1SShri Abhyankar + dm - The DMNetwork object
4615f2c45f1SShri Abhyankar . p  - vertex/edge point
4625f2c45f1SShri Abhyankar 
4635f2c45f1SShri Abhyankar   Output Parameters:
4645f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge
4655f2c45f1SShri Abhyankar 
4665f2c45f1SShri Abhyankar   Level: intermediate
4675f2c45f1SShri Abhyankar 
4685f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
4695f2c45f1SShri Abhyankar @*/
4705f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
4715f2c45f1SShri Abhyankar {
4725f2c45f1SShri Abhyankar   PetscErrorCode ierr;
4735f2c45f1SShri Abhyankar   PetscInt       offset;
4745f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4755f2c45f1SShri Abhyankar 
4765f2c45f1SShri Abhyankar   PetscFunctionBegin;
4775f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
4785f2c45f1SShri Abhyankar   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
4795f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4805f2c45f1SShri Abhyankar }
4815f2c45f1SShri Abhyankar 
4825f2c45f1SShri Abhyankar /*@
483a730d845SHong Zhang   DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the
4845f2c45f1SShri Abhyankar                                     component value from the component data array
4855f2c45f1SShri Abhyankar 
4865f2c45f1SShri Abhyankar   Not Collective
4875f2c45f1SShri Abhyankar 
4885f2c45f1SShri Abhyankar   Input Parameters:
4895f2c45f1SShri Abhyankar + dm      - The DMNetwork object
4905f2c45f1SShri Abhyankar . p       - vertex/edge point
4915f2c45f1SShri Abhyankar - compnum - component number
4925f2c45f1SShri Abhyankar 
4935f2c45f1SShri Abhyankar   Output Parameters:
4945f2c45f1SShri Abhyankar + compkey - the key obtained when registering the component
4955f2c45f1SShri Abhyankar - offset  - offset into the component data array associated with the vertex/edge point
4965f2c45f1SShri Abhyankar 
4975f2c45f1SShri Abhyankar   Notes:
4985f2c45f1SShri Abhyankar   Typical usage:
4995f2c45f1SShri Abhyankar 
5005f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray(dm, &arr);
5015f2c45f1SShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
5025f2c45f1SShri Abhyankar   Loop over vertices or edges
5035f2c45f1SShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
5045f2c45f1SShri Abhyankar     Loop over numcomps
505a730d845SHong Zhang       DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset);
5065f2c45f1SShri Abhyankar       compdata = (UserCompDataType)(arr+offset);
5075f2c45f1SShri Abhyankar 
5085f2c45f1SShri Abhyankar   Level: intermediate
5095f2c45f1SShri Abhyankar 
5105f2c45f1SShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
5115f2c45f1SShri Abhyankar @*/
512a730d845SHong Zhang PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
5135f2c45f1SShri Abhyankar {
5145f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
5155f2c45f1SShri Abhyankar   PetscInt                 offsetp;
5165f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
5175f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
5185f2c45f1SShri Abhyankar 
5195f2c45f1SShri Abhyankar   PetscFunctionBegin;
5205f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
5215f2c45f1SShri Abhyankar   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
522d36f4e81SHong Zhang   if (compkey) *compkey = header->key[compnum];
523d36f4e81SHong Zhang   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
5245f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5255f2c45f1SShri Abhyankar }
5265f2c45f1SShri Abhyankar 
5275f2c45f1SShri Abhyankar /*@
5285f2c45f1SShri Abhyankar   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
5295f2c45f1SShri Abhyankar 
5305f2c45f1SShri Abhyankar   Not Collective
5315f2c45f1SShri Abhyankar 
5325f2c45f1SShri Abhyankar   Input Parameters:
5335f2c45f1SShri Abhyankar + dm     - The DMNetwork object
5345f2c45f1SShri Abhyankar - p      - the edge/vertex point
5355f2c45f1SShri Abhyankar 
5365f2c45f1SShri Abhyankar   Output Parameters:
5375f2c45f1SShri Abhyankar . offset - the offset
5385f2c45f1SShri Abhyankar 
5395f2c45f1SShri Abhyankar   Level: intermediate
5405f2c45f1SShri Abhyankar 
5415f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
5425f2c45f1SShri Abhyankar @*/
5435f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
5445f2c45f1SShri Abhyankar {
5455f2c45f1SShri Abhyankar   PetscErrorCode ierr;
5465f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5475f2c45f1SShri Abhyankar 
5485f2c45f1SShri Abhyankar   PetscFunctionBegin;
5495f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr);
5505f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5515f2c45f1SShri Abhyankar }
5525f2c45f1SShri Abhyankar 
5535f2c45f1SShri Abhyankar /*@
5545f2c45f1SShri Abhyankar   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
5555f2c45f1SShri Abhyankar 
5565f2c45f1SShri Abhyankar   Not Collective
5575f2c45f1SShri Abhyankar 
5585f2c45f1SShri Abhyankar   Input Parameters:
5595f2c45f1SShri Abhyankar + dm      - The DMNetwork object
5605f2c45f1SShri Abhyankar - p       - the edge/vertex point
5615f2c45f1SShri Abhyankar 
5625f2c45f1SShri Abhyankar   Output Parameters:
5635f2c45f1SShri Abhyankar . offsetg - the offset
5645f2c45f1SShri Abhyankar 
5655f2c45f1SShri Abhyankar   Level: intermediate
5665f2c45f1SShri Abhyankar 
5675f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
5685f2c45f1SShri Abhyankar @*/
5695f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
5705f2c45f1SShri Abhyankar {
5715f2c45f1SShri Abhyankar   PetscErrorCode ierr;
5725f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5735f2c45f1SShri Abhyankar 
5745f2c45f1SShri Abhyankar   PetscFunctionBegin;
5755f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr);
5766fefedf4SHong Zhang   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */
5775f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5785f2c45f1SShri Abhyankar }
5795f2c45f1SShri Abhyankar 
58024121865SAdrian Maldonado /*@
58124121865SAdrian Maldonado   DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector.
58224121865SAdrian Maldonado 
58324121865SAdrian Maldonado   Not Collective
58424121865SAdrian Maldonado 
58524121865SAdrian Maldonado   Input Parameters:
58624121865SAdrian Maldonado + dm     - The DMNetwork object
58724121865SAdrian Maldonado - p      - the edge point
58824121865SAdrian Maldonado 
58924121865SAdrian Maldonado   Output Parameters:
59024121865SAdrian Maldonado . offset - the offset
59124121865SAdrian Maldonado 
59224121865SAdrian Maldonado   Level: intermediate
59324121865SAdrian Maldonado 
59424121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
59524121865SAdrian Maldonado @*/
59624121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
59724121865SAdrian Maldonado {
59824121865SAdrian Maldonado   PetscErrorCode ierr;
59924121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
60024121865SAdrian Maldonado 
60124121865SAdrian Maldonado   PetscFunctionBegin;
60224121865SAdrian Maldonado 
60324121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr);
60424121865SAdrian Maldonado   PetscFunctionReturn(0);
60524121865SAdrian Maldonado }
60624121865SAdrian Maldonado 
60724121865SAdrian Maldonado /*@
60824121865SAdrian Maldonado   DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector.
60924121865SAdrian Maldonado 
61024121865SAdrian Maldonado   Not Collective
61124121865SAdrian Maldonado 
61224121865SAdrian Maldonado   Input Parameters:
61324121865SAdrian Maldonado + dm     - The DMNetwork object
61424121865SAdrian Maldonado - p      - the vertex point
61524121865SAdrian Maldonado 
61624121865SAdrian Maldonado   Output Parameters:
61724121865SAdrian Maldonado . offset - the offset
61824121865SAdrian Maldonado 
61924121865SAdrian Maldonado   Level: intermediate
62024121865SAdrian Maldonado 
62124121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
62224121865SAdrian Maldonado @*/
62324121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
62424121865SAdrian Maldonado {
62524121865SAdrian Maldonado   PetscErrorCode ierr;
62624121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
62724121865SAdrian Maldonado 
62824121865SAdrian Maldonado   PetscFunctionBegin;
62924121865SAdrian Maldonado 
63024121865SAdrian Maldonado   p -= network->vStart;
63124121865SAdrian Maldonado 
63224121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr);
63324121865SAdrian Maldonado   PetscFunctionReturn(0);
63424121865SAdrian Maldonado }
6355f2c45f1SShri Abhyankar /*@
6365f2c45f1SShri Abhyankar   DMNetworkAddNumVariables - Add number of variables associated with a given point.
6375f2c45f1SShri Abhyankar 
6385f2c45f1SShri Abhyankar   Not Collective
6395f2c45f1SShri Abhyankar 
6405f2c45f1SShri Abhyankar   Input Parameters:
6415f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
6425f2c45f1SShri Abhyankar . p    - the vertex/edge point
6435f2c45f1SShri Abhyankar - nvar - number of additional variables
6445f2c45f1SShri Abhyankar 
6455f2c45f1SShri Abhyankar   Level: intermediate
6465f2c45f1SShri Abhyankar 
6475f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables
6485f2c45f1SShri Abhyankar @*/
6495f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
6505f2c45f1SShri Abhyankar {
6515f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6525f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6535f2c45f1SShri Abhyankar 
6545f2c45f1SShri Abhyankar   PetscFunctionBegin;
6555f2c45f1SShri Abhyankar   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
6565f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6575f2c45f1SShri Abhyankar }
6585f2c45f1SShri Abhyankar 
65927f51fceSHong Zhang /*@
66027f51fceSHong Zhang   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
66127f51fceSHong Zhang 
66227f51fceSHong Zhang   Not Collective
66327f51fceSHong Zhang 
66427f51fceSHong Zhang   Input Parameters:
66527f51fceSHong Zhang + dm   - The DMNetworkObject
66627f51fceSHong Zhang - p    - the vertex/edge point
66727f51fceSHong Zhang 
66827f51fceSHong Zhang   Output Parameters:
66927f51fceSHong Zhang . nvar - number of variables
67027f51fceSHong Zhang 
67127f51fceSHong Zhang   Level: intermediate
67227f51fceSHong Zhang 
67327f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
67427f51fceSHong Zhang @*/
67527f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
67627f51fceSHong Zhang {
67727f51fceSHong Zhang   PetscErrorCode ierr;
67827f51fceSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
67927f51fceSHong Zhang 
68027f51fceSHong Zhang   PetscFunctionBegin;
68127f51fceSHong Zhang   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
68227f51fceSHong Zhang   PetscFunctionReturn(0);
68327f51fceSHong Zhang }
68427f51fceSHong Zhang 
6855f2c45f1SShri Abhyankar /*@
6865f2c45f1SShri Abhyankar   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
6875f2c45f1SShri Abhyankar 
6885f2c45f1SShri Abhyankar   Not Collective
6895f2c45f1SShri Abhyankar 
6905f2c45f1SShri Abhyankar   Input Parameters:
6915f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
6925f2c45f1SShri Abhyankar . p    - the vertex/edge point
6935f2c45f1SShri Abhyankar - nvar - number of variables
6945f2c45f1SShri Abhyankar 
6955f2c45f1SShri Abhyankar   Level: intermediate
6965f2c45f1SShri Abhyankar 
6975f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables
6985f2c45f1SShri Abhyankar @*/
6995f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
7005f2c45f1SShri Abhyankar {
7015f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7025f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7035f2c45f1SShri Abhyankar 
7045f2c45f1SShri Abhyankar   PetscFunctionBegin;
7055f2c45f1SShri Abhyankar   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
7065f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7075f2c45f1SShri Abhyankar }
7085f2c45f1SShri Abhyankar 
7095f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This
7105f2c45f1SShri Abhyankar    function is called during DMSetUp() */
7115f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm)
7125f2c45f1SShri Abhyankar {
7135f2c45f1SShri Abhyankar   PetscErrorCode              ierr;
7145f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7155f2c45f1SShri Abhyankar   PetscInt                    arr_size;
7165f2c45f1SShri Abhyankar   PetscInt                    p,offset,offsetp;
7175f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
7185f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue;
7195f2c45f1SShri Abhyankar   DMNetworkComponentGenericDataType      *componentdataarray;
7205f2c45f1SShri Abhyankar   PetscInt ncomp, i;
7215f2c45f1SShri Abhyankar 
7225f2c45f1SShri Abhyankar   PetscFunctionBegin;
7235f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
7245f2c45f1SShri Abhyankar   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
72575b160a0SShri Abhyankar   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
7265f2c45f1SShri Abhyankar   componentdataarray = network->componentdataarray;
7275f2c45f1SShri Abhyankar   for (p = network->pStart; p < network->pEnd; p++) {
7285f2c45f1SShri Abhyankar     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
7295f2c45f1SShri Abhyankar     /* Copy header */
7305f2c45f1SShri Abhyankar     header = &network->header[p];
731302440fdSBarry Smith     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
7325f2c45f1SShri Abhyankar     /* Copy data */
7335f2c45f1SShri Abhyankar     cvalue = &network->cvalue[p];
7345f2c45f1SShri Abhyankar     ncomp = header->ndata;
7355f2c45f1SShri Abhyankar     for (i = 0; i < ncomp; i++) {
7365f2c45f1SShri Abhyankar       offset = offsetp + network->dataheadersize + header->offset[i];
737302440fdSBarry Smith       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
7385f2c45f1SShri Abhyankar     }
7395f2c45f1SShri Abhyankar   }
7405f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7415f2c45f1SShri Abhyankar }
7425f2c45f1SShri Abhyankar 
7435f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */
7445f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm)
7455f2c45f1SShri Abhyankar {
7465f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7475f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7485f2c45f1SShri Abhyankar 
7495f2c45f1SShri Abhyankar   PetscFunctionBegin;
7505f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
7515f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7525f2c45f1SShri Abhyankar }
7535f2c45f1SShri Abhyankar 
7545f2c45f1SShri Abhyankar /*@C
7555f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray - Returns the component data array
7565f2c45f1SShri Abhyankar 
7575f2c45f1SShri Abhyankar   Not Collective
7585f2c45f1SShri Abhyankar 
7595f2c45f1SShri Abhyankar   Input Parameters:
7605f2c45f1SShri Abhyankar . dm - The DMNetwork Object
7615f2c45f1SShri Abhyankar 
7625f2c45f1SShri Abhyankar   Output Parameters:
7635f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components
7645f2c45f1SShri Abhyankar 
7655f2c45f1SShri Abhyankar   Level: intermediate
7665f2c45f1SShri Abhyankar 
767a730d845SHong Zhang .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents
7685f2c45f1SShri Abhyankar @*/
7695f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
7705f2c45f1SShri Abhyankar {
7715f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7725f2c45f1SShri Abhyankar 
7735f2c45f1SShri Abhyankar   PetscFunctionBegin;
7745f2c45f1SShri Abhyankar   *componentdataarray = network->componentdataarray;
7755f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7765f2c45f1SShri Abhyankar }
7775f2c45f1SShri Abhyankar 
77824121865SAdrian Maldonado /* Get a subsection from a range of points */
77924121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection)
78024121865SAdrian Maldonado {
78124121865SAdrian Maldonado   PetscErrorCode ierr;
78224121865SAdrian Maldonado   PetscInt       i, nvar;
78324121865SAdrian Maldonado 
78424121865SAdrian Maldonado   PetscFunctionBegin;
78524121865SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr);
78624121865SAdrian Maldonado   ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr);
78724121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
78824121865SAdrian Maldonado     ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr);
78924121865SAdrian Maldonado     ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr);
79024121865SAdrian Maldonado   }
79124121865SAdrian Maldonado 
79224121865SAdrian Maldonado   ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr);
79324121865SAdrian Maldonado   PetscFunctionReturn(0);
79424121865SAdrian Maldonado }
79524121865SAdrian Maldonado 
79624121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */
79724121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
79824121865SAdrian Maldonado {
79924121865SAdrian Maldonado   PetscErrorCode ierr;
80024121865SAdrian Maldonado   PetscInt       i, *subpoints;
80124121865SAdrian Maldonado 
80224121865SAdrian Maldonado   PetscFunctionBegin;
80324121865SAdrian Maldonado   /* Create index sets to map from "points" to "subpoints" */
80424121865SAdrian Maldonado   ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr);
80524121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
80624121865SAdrian Maldonado     subpoints[i - pstart] = i;
80724121865SAdrian Maldonado   }
808459726d8SSatish Balay   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr);
80924121865SAdrian Maldonado   ierr = PetscFree(subpoints);CHKERRQ(ierr);
81024121865SAdrian Maldonado   PetscFunctionReturn(0);
81124121865SAdrian Maldonado }
81224121865SAdrian Maldonado 
81324121865SAdrian Maldonado /*@
81424121865SAdrian Maldonado   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute.
81524121865SAdrian Maldonado 
81624121865SAdrian Maldonado   Collective
81724121865SAdrian Maldonado 
81824121865SAdrian Maldonado   Input Parameters:
81924121865SAdrian Maldonado . dm   - The DMNetworkObject
82024121865SAdrian Maldonado 
82124121865SAdrian Maldonado   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
82224121865SAdrian Maldonado 
82324121865SAdrian Maldonado   points = [0 1 2 3 4 5 6]
82424121865SAdrian Maldonado 
82524121865SAdrian Maldonado   where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]).
82624121865SAdrian Maldonado 
82724121865SAdrian Maldonado   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
82824121865SAdrian Maldonado 
82924121865SAdrian Maldonado   Level: intermediate
83024121865SAdrian Maldonado 
83124121865SAdrian Maldonado @*/
83224121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
83324121865SAdrian Maldonado {
83424121865SAdrian Maldonado   PetscErrorCode ierr;
83524121865SAdrian Maldonado   MPI_Comm       comm;
8369852e123SBarry Smith   PetscMPIInt    rank, size;
83724121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
83824121865SAdrian Maldonado 
839eab1376dSHong Zhang   PetscFunctionBegin;
84024121865SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
84124121865SAdrian Maldonado   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
8429852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
84324121865SAdrian Maldonado 
84424121865SAdrian Maldonado   /* Create maps for vertices and edges */
84524121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr);
84624121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr);
84724121865SAdrian Maldonado 
84824121865SAdrian Maldonado   /* Create local sub-sections */
84924121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr);
85024121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr);
85124121865SAdrian Maldonado 
8529852e123SBarry Smith   if (size > 1) {
85324121865SAdrian Maldonado     ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr);
85424121865SAdrian Maldonado     ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr);
85524121865SAdrian Maldonado   ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr);
85624121865SAdrian Maldonado   ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr);
85724121865SAdrian Maldonado   } else {
85824121865SAdrian Maldonado   /* create structures for vertex */
85924121865SAdrian Maldonado   ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr);
86024121865SAdrian Maldonado   /* create structures for edge */
86124121865SAdrian Maldonado   ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr);
86224121865SAdrian Maldonado   }
86324121865SAdrian Maldonado 
86424121865SAdrian Maldonado 
86524121865SAdrian Maldonado   /* Add viewers */
86624121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr);
86724121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr);
86824121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr);
86924121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr);
87024121865SAdrian Maldonado 
87124121865SAdrian Maldonado   PetscFunctionReturn(0);
87224121865SAdrian Maldonado }
8737b6afd5bSHong Zhang 
8745f2c45f1SShri Abhyankar /*@
8755f2c45f1SShri Abhyankar   DMNetworkDistribute - Distributes the network and moves associated component data.
8765f2c45f1SShri Abhyankar 
8775f2c45f1SShri Abhyankar   Collective
8785f2c45f1SShri Abhyankar 
8795f2c45f1SShri Abhyankar   Input Parameter:
880d3464fd4SAdrian Maldonado + DM - the DMNetwork object
8815f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default
8825f2c45f1SShri Abhyankar 
8835f2c45f1SShri Abhyankar   Notes:
8848b171c8eSHong Zhang   Distributes the network with <overlap>-overlapping partitioning of the edges.
8855f2c45f1SShri Abhyankar 
8865f2c45f1SShri Abhyankar   Level: intermediate
8875f2c45f1SShri Abhyankar 
8885f2c45f1SShri Abhyankar .seealso: DMNetworkCreate
8895f2c45f1SShri Abhyankar @*/
890d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
8915f2c45f1SShri Abhyankar {
892d3464fd4SAdrian Maldonado   MPI_Comm       comm;
8935f2c45f1SShri Abhyankar   PetscErrorCode ierr;
894d3464fd4SAdrian Maldonado   PetscMPIInt    size;
895d3464fd4SAdrian Maldonado   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
896d3464fd4SAdrian Maldonado   DM_Network     *newDMnetwork;
8975f2c45f1SShri Abhyankar   PetscSF        pointsf;
8985f2c45f1SShri Abhyankar   DM             newDM;
89951ac5effSHong Zhang   PetscPartitioner part;
900b9c6e19dSShri Abhyankar   PetscInt         j,e,v,offset;
901b9c6e19dSShri Abhyankar   DMNetworkComponentHeader header;
9025f2c45f1SShri Abhyankar 
9035f2c45f1SShri Abhyankar   PetscFunctionBegin;
904d3464fd4SAdrian Maldonado 
905d3464fd4SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr);
906d3464fd4SAdrian Maldonado   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
907d3464fd4SAdrian Maldonado   if (size == 1) PetscFunctionReturn(0);
908d3464fd4SAdrian Maldonado 
909d3464fd4SAdrian Maldonado   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr);
9105f2c45f1SShri Abhyankar   newDMnetwork = (DM_Network*)newDM->data;
9115f2c45f1SShri Abhyankar   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
91251ac5effSHong Zhang 
91351ac5effSHong Zhang   /* Enable runtime options for petscpartitioner */
91451ac5effSHong Zhang   ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr);
91551ac5effSHong Zhang   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
91651ac5effSHong Zhang 
9175f2c45f1SShri Abhyankar   /* Distribute plex dm and dof section */
91880cf41d5SMatthew G. Knepley   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
91951ac5effSHong Zhang 
9205f2c45f1SShri Abhyankar   /* Distribute dof section */
921d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr);
9225f2c45f1SShri Abhyankar   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
923d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr);
92451ac5effSHong Zhang 
9255f2c45f1SShri Abhyankar   /* Distribute data and associated section */
92631da1fc8SHong Zhang   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
92724121865SAdrian Maldonado 
9285f2c45f1SShri Abhyankar   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
9295f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
9305f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
9315f2c45f1SShri Abhyankar   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
9326fefedf4SHong Zhang   newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
9336fefedf4SHong Zhang   newDMnetwork->NVertices = oldDMnetwork->NVertices;
9345f2c45f1SShri Abhyankar   newDMnetwork->NEdges = oldDMnetwork->NEdges;
93524121865SAdrian Maldonado 
9365f2c45f1SShri Abhyankar   /* Set Dof section as the default section for dm */
9375f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
9385f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
9395f2c45f1SShri Abhyankar 
940b9c6e19dSShri Abhyankar   /* Set up subnetwork info in the newDM */
941b9c6e19dSShri Abhyankar   newDMnetwork->nsubnet = oldDMnetwork->nsubnet;
942b9c6e19dSShri Abhyankar   ierr = PetscCalloc1(newDMnetwork->nsubnet,&newDMnetwork->subnet);CHKERRQ(ierr);
943b9c6e19dSShri Abhyankar   /* Copy over the global number of vertices and edges in each subnetwork. Note that these are already
944b9c6e19dSShri Abhyankar      calculated in DMNetworkLayoutSetUp()
945b9c6e19dSShri Abhyankar   */
946b9c6e19dSShri Abhyankar   for(j=0; j < newDMnetwork->nsubnet; j++) {
947b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx;
948b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge;
949b9c6e19dSShri Abhyankar   }
950b9c6e19dSShri Abhyankar 
951b9c6e19dSShri Abhyankar   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
952b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
953b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
954b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].nedge++;
955b9c6e19dSShri Abhyankar   }
956b9c6e19dSShri Abhyankar 
957b9c6e19dSShri Abhyankar   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
958b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
959b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
960b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].nvtx++;
961b9c6e19dSShri Abhyankar   }
962b9c6e19dSShri Abhyankar 
963b9c6e19dSShri Abhyankar   /* Now create the vertices and edge arrays for the subnetworks */
964b9c6e19dSShri Abhyankar   for(j=0; j < newDMnetwork->nsubnet; j++) {
965b9c6e19dSShri Abhyankar     ierr = PetscCalloc1(newDMnetwork->subnet[j].nedge,&newDMnetwork->subnet[j].edges);CHKERRQ(ierr);
966b9c6e19dSShri Abhyankar     ierr = PetscCalloc1(newDMnetwork->subnet[j].nvtx,&newDMnetwork->subnet[j].vertices);CHKERRQ(ierr);
967b9c6e19dSShri Abhyankar     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
968b9c6e19dSShri Abhyankar        These get updated when the vertices and edges are added. */
969b9c6e19dSShri Abhyankar     newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0;
970b9c6e19dSShri Abhyankar   }
971b9c6e19dSShri Abhyankar 
972b9c6e19dSShri Abhyankar   /* Set the vertices and edges in each subnetwork */
973b9c6e19dSShri Abhyankar   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
974b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
975b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
976b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++]  = e;
977b9c6e19dSShri Abhyankar   }
978b9c6e19dSShri Abhyankar 
979b9c6e19dSShri Abhyankar   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
980b9c6e19dSShri Abhyankar     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
981b9c6e19dSShri Abhyankar     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
982b9c6e19dSShri Abhyankar     newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++]  = v;
983b9c6e19dSShri Abhyankar   }
984b9c6e19dSShri Abhyankar 
98524121865SAdrian Maldonado   /* Destroy point SF */
98624121865SAdrian Maldonado   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
98724121865SAdrian Maldonado 
988d3464fd4SAdrian Maldonado   ierr = DMDestroy(dm);CHKERRQ(ierr);
989d3464fd4SAdrian Maldonado   *dm  = newDM;
9905f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9915f2c45f1SShri Abhyankar }
9925f2c45f1SShri Abhyankar 
99324121865SAdrian Maldonado /*@C
99424121865SAdrian Maldonado   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering.
99524121865SAdrian Maldonado 
99624121865SAdrian Maldonado   Input Parameters:
99724121865SAdrian Maldonado + masterSF - the original SF structure
99824121865SAdrian Maldonado - map      - a ISLocalToGlobal mapping that contains the subset of points
99924121865SAdrian Maldonado 
100024121865SAdrian Maldonado   Output Parameters:
100124121865SAdrian Maldonado . subSF    - a subset of the masterSF for the desired subset.
100224121865SAdrian Maldonado */
100324121865SAdrian Maldonado 
100424121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) {
100524121865SAdrian Maldonado 
100624121865SAdrian Maldonado   PetscErrorCode        ierr;
100724121865SAdrian Maldonado   PetscInt              nroots, nleaves, *ilocal_sub;
100824121865SAdrian Maldonado   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
100924121865SAdrian Maldonado   PetscInt              *local_points, *remote_points;
101024121865SAdrian Maldonado   PetscSFNode           *iremote_sub;
101124121865SAdrian Maldonado   const PetscInt        *ilocal;
101224121865SAdrian Maldonado   const PetscSFNode     *iremote;
101324121865SAdrian Maldonado 
101424121865SAdrian Maldonado   PetscFunctionBegin;
101524121865SAdrian Maldonado   ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
101624121865SAdrian Maldonado 
101724121865SAdrian Maldonado   /* Look for leaves that pertain to the subset of points. Get the local ordering */
101824121865SAdrian Maldonado   ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr);
101924121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr);
102024121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
102124121865SAdrian Maldonado     if (ilocal_map[i] != -1) nleaves_sub += 1;
102224121865SAdrian Maldonado   }
102324121865SAdrian Maldonado   /* Re-number ilocal with subset numbering. Need information from roots */
102424121865SAdrian Maldonado   ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr);
102524121865SAdrian Maldonado   for (i = 0; i < nroots; i++) local_points[i] = i;
102624121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr);
102724121865SAdrian Maldonado   ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
102824121865SAdrian Maldonado   ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
102924121865SAdrian Maldonado   /* Fill up graph using local (that is, local to the subset) numbering. */
10304b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr);
10314b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr);
103224121865SAdrian Maldonado   nleaves_sub = 0;
103324121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
103424121865SAdrian Maldonado     if (ilocal_map[i] != -1) {
103524121865SAdrian Maldonado       ilocal_sub[nleaves_sub] = ilocal_map[i];
10364b70a8deSAdrian Maldonado       iremote_sub[nleaves_sub].rank = iremote[i].rank;
103724121865SAdrian Maldonado       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
103824121865SAdrian Maldonado       nleaves_sub += 1;
103924121865SAdrian Maldonado     }
104024121865SAdrian Maldonado   }
104124121865SAdrian Maldonado   ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr);
104224121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr);
104324121865SAdrian Maldonado 
104424121865SAdrian Maldonado   /* Create new subSF */
104524121865SAdrian Maldonado   ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr);
104624121865SAdrian Maldonado   ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr);
10474b70a8deSAdrian Maldonado   ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr);
104824121865SAdrian Maldonado   ierr = PetscFree(ilocal_map);CHKERRQ(ierr);
10494b70a8deSAdrian Maldonado   ierr = PetscFree(iremote_sub);CHKERRQ(ierr);
105024121865SAdrian Maldonado   PetscFunctionReturn(0);
105124121865SAdrian Maldonado }
105224121865SAdrian Maldonado 
10535f2c45f1SShri Abhyankar /*@C
10545f2c45f1SShri Abhyankar   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
10555f2c45f1SShri Abhyankar 
10565f2c45f1SShri Abhyankar   Not Collective
10575f2c45f1SShri Abhyankar 
10585f2c45f1SShri Abhyankar   Input Parameters:
10595f2c45f1SShri Abhyankar + dm - The DMNetwork object
10605f2c45f1SShri Abhyankar - p  - the vertex point
10615f2c45f1SShri Abhyankar 
10625f2c45f1SShri Abhyankar   Output Paramters:
10635f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point
10645f2c45f1SShri Abhyankar - edges  - List of edge points
10655f2c45f1SShri Abhyankar 
10665f2c45f1SShri Abhyankar   Level: intermediate
10675f2c45f1SShri Abhyankar 
10685f2c45f1SShri Abhyankar   Fortran Notes:
10695f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
10705f2c45f1SShri Abhyankar   include petsc.h90 in your code.
10715f2c45f1SShri Abhyankar 
1072d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices
10735f2c45f1SShri Abhyankar @*/
10745f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
10755f2c45f1SShri Abhyankar {
10765f2c45f1SShri Abhyankar   PetscErrorCode ierr;
10775f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
10785f2c45f1SShri Abhyankar 
10795f2c45f1SShri Abhyankar   PetscFunctionBegin;
10805f2c45f1SShri Abhyankar   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
10815f2c45f1SShri Abhyankar   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
10825f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10835f2c45f1SShri Abhyankar }
10845f2c45f1SShri Abhyankar 
10855f2c45f1SShri Abhyankar /*@C
1086d842c372SHong Zhang   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
10875f2c45f1SShri Abhyankar 
10885f2c45f1SShri Abhyankar   Not Collective
10895f2c45f1SShri Abhyankar 
10905f2c45f1SShri Abhyankar   Input Parameters:
10915f2c45f1SShri Abhyankar + dm - The DMNetwork object
10925f2c45f1SShri Abhyankar - p  - the edge point
10935f2c45f1SShri Abhyankar 
10945f2c45f1SShri Abhyankar   Output Paramters:
10955f2c45f1SShri Abhyankar . vertices  - vertices connected to this edge
10965f2c45f1SShri Abhyankar 
10975f2c45f1SShri Abhyankar   Level: intermediate
10985f2c45f1SShri Abhyankar 
10995f2c45f1SShri Abhyankar   Fortran Notes:
11005f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
11015f2c45f1SShri Abhyankar   include petsc.h90 in your code.
11025f2c45f1SShri Abhyankar 
11035f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
11045f2c45f1SShri Abhyankar @*/
1105d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
11065f2c45f1SShri Abhyankar {
11075f2c45f1SShri Abhyankar   PetscErrorCode ierr;
11085f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
11095f2c45f1SShri Abhyankar 
11105f2c45f1SShri Abhyankar   PetscFunctionBegin;
11115f2c45f1SShri Abhyankar   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
11125f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
11135f2c45f1SShri Abhyankar }
11145f2c45f1SShri Abhyankar 
11155f2c45f1SShri Abhyankar /*@
11165f2c45f1SShri Abhyankar   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
11175f2c45f1SShri Abhyankar 
11185f2c45f1SShri Abhyankar   Not Collective
11195f2c45f1SShri Abhyankar 
11205f2c45f1SShri Abhyankar   Input Parameters:
11215f2c45f1SShri Abhyankar + dm - The DMNetwork object
11225f2c45f1SShri Abhyankar . p  - the vertex point
11235f2c45f1SShri Abhyankar 
11245f2c45f1SShri Abhyankar   Output Parameter:
11255f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point
11265f2c45f1SShri Abhyankar 
11275f2c45f1SShri Abhyankar   Level: intermediate
11285f2c45f1SShri Abhyankar 
1129d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange
11305f2c45f1SShri Abhyankar @*/
11315f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
11325f2c45f1SShri Abhyankar {
11335f2c45f1SShri Abhyankar   PetscErrorCode ierr;
11345f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
11355f2c45f1SShri Abhyankar   PetscInt       offsetg;
11365f2c45f1SShri Abhyankar   PetscSection   sectiong;
11375f2c45f1SShri Abhyankar 
11385f2c45f1SShri Abhyankar   PetscFunctionBegin;
11395f2c45f1SShri Abhyankar   *isghost = PETSC_FALSE;
11405f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
11415f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
11425f2c45f1SShri Abhyankar   if (offsetg < 0) *isghost = PETSC_TRUE;
11435f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
11445f2c45f1SShri Abhyankar }
11455f2c45f1SShri Abhyankar 
11465f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm)
11475f2c45f1SShri Abhyankar {
11485f2c45f1SShri Abhyankar   PetscErrorCode ierr;
11495f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
11505f2c45f1SShri Abhyankar 
11515f2c45f1SShri Abhyankar   PetscFunctionBegin;
11525f2c45f1SShri Abhyankar   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
11535f2c45f1SShri Abhyankar   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
11545f2c45f1SShri Abhyankar 
11555f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
11565f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
11575f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
11585f2c45f1SShri Abhyankar }
11595f2c45f1SShri Abhyankar 
11601ad426b7SHong Zhang /*@
116117df6e9eSHong Zhang     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
11621ad426b7SHong Zhang                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
11631ad426b7SHong Zhang 
11641ad426b7SHong Zhang     Collective
11651ad426b7SHong Zhang 
11661ad426b7SHong Zhang     Input Parameters:
116783b2e829SHong Zhang +   dm - The DMNetwork object
116883b2e829SHong Zhang .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
116983b2e829SHong Zhang -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
11701ad426b7SHong Zhang 
11711ad426b7SHong Zhang     Level: intermediate
11721ad426b7SHong Zhang 
11731ad426b7SHong Zhang @*/
117483b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
11751ad426b7SHong Zhang {
11761ad426b7SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
1177*8675203cSHong Zhang   PetscErrorCode ierr;
11781ad426b7SHong Zhang 
11791ad426b7SHong Zhang   PetscFunctionBegin;
118083b2e829SHong Zhang   network->userEdgeJacobian   = eflg;
118183b2e829SHong Zhang   network->userVertexJacobian = vflg;
1182*8675203cSHong Zhang 
1183*8675203cSHong Zhang   if (eflg && !network->Je) {
1184*8675203cSHong Zhang     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
1185*8675203cSHong Zhang   }
1186*8675203cSHong Zhang 
1187*8675203cSHong Zhang   if (vflg && !network->Jv) {
1188*8675203cSHong Zhang     PetscInt       i,*vptr,nedges,vStart=network->vStart;
1189*8675203cSHong Zhang     PetscInt       nVertices = network->nVertices,nedges_total;
1190*8675203cSHong Zhang     const PetscInt *edges;
1191*8675203cSHong Zhang 
1192*8675203cSHong Zhang     /* count nvertex_total */
1193*8675203cSHong Zhang     nedges_total = 0;
1194*8675203cSHong Zhang     ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr);
1195*8675203cSHong Zhang 
1196*8675203cSHong Zhang     vptr[0] = 0;
1197*8675203cSHong Zhang     for (i=0; i<nVertices; i++) {
1198*8675203cSHong Zhang       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
1199*8675203cSHong Zhang       nedges_total += nedges;
1200*8675203cSHong Zhang       vptr[i+1] = vptr[i] + 2*nedges + 1;
1201*8675203cSHong Zhang     }
1202*8675203cSHong Zhang 
1203*8675203cSHong Zhang     ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr);
1204*8675203cSHong Zhang     network->Jvptr = vptr;
1205*8675203cSHong Zhang   }
12061ad426b7SHong Zhang   PetscFunctionReturn(0);
12071ad426b7SHong Zhang }
12081ad426b7SHong Zhang 
12091ad426b7SHong Zhang /*@
121083b2e829SHong Zhang     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
121183b2e829SHong Zhang 
121283b2e829SHong Zhang     Not Collective
121383b2e829SHong Zhang 
121483b2e829SHong Zhang     Input Parameters:
121583b2e829SHong Zhang +   dm - The DMNetwork object
121683b2e829SHong Zhang .   p  - the edge point
12173e97b6e8SHong Zhang -   J - array (size = 3) of Jacobian submatrices for this edge point:
12183e97b6e8SHong Zhang         J[0]: this edge
1219d842c372SHong Zhang         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()
122083b2e829SHong Zhang 
122183b2e829SHong Zhang     Level: intermediate
122283b2e829SHong Zhang 
122383b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix
122483b2e829SHong Zhang @*/
122583b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
122683b2e829SHong Zhang {
122783b2e829SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
122883b2e829SHong Zhang 
122983b2e829SHong Zhang   PetscFunctionBegin;
1230*8675203cSHong Zhang   if (!network->Je) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
1231*8675203cSHong Zhang 
1232*8675203cSHong Zhang   if (J) {
1233883e35e8SHong Zhang     network->Je[3*p]   = J[0];
1234883e35e8SHong Zhang     network->Je[3*p+1] = J[1];
1235883e35e8SHong Zhang     network->Je[3*p+2] = J[2];
1236*8675203cSHong Zhang   }
123783b2e829SHong Zhang   PetscFunctionReturn(0);
123883b2e829SHong Zhang }
123983b2e829SHong Zhang 
124083b2e829SHong Zhang /*@
124176ddfea5SHong Zhang     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
12421ad426b7SHong Zhang 
12431ad426b7SHong Zhang     Not Collective
12441ad426b7SHong Zhang 
12451ad426b7SHong Zhang     Input Parameters:
12461ad426b7SHong Zhang +   dm - The DMNetwork object
12471ad426b7SHong Zhang .   p  - the vertex point
12483e97b6e8SHong Zhang -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
12493e97b6e8SHong Zhang         J[0]:       this vertex
12503e97b6e8SHong Zhang         J[1+2*i]:   i-th supporting edge
12513e97b6e8SHong Zhang         J[1+2*i+1]: i-th connected vertex
12521ad426b7SHong Zhang 
12531ad426b7SHong Zhang     Level: intermediate
12541ad426b7SHong Zhang 
125583b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix
12561ad426b7SHong Zhang @*/
1257883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
12585f2c45f1SShri Abhyankar {
12595f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12605f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
1261*8675203cSHong Zhang   PetscInt       i,*vptr,nedges,vStart=network->vStart;
1262883e35e8SHong Zhang   const PetscInt *edges;
12635f2c45f1SShri Abhyankar 
12645f2c45f1SShri Abhyankar   PetscFunctionBegin;
1265*8675203cSHong Zhang   if (!network->Jv) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
1266883e35e8SHong Zhang 
1267*8675203cSHong Zhang   if (J) {
1268883e35e8SHong Zhang     vptr = network->Jvptr;
12693e97b6e8SHong Zhang     network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
12703e97b6e8SHong Zhang 
12713e97b6e8SHong Zhang     /* Set Jacobian for each supporting edge and connected vertex */
1272883e35e8SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
1273883e35e8SHong Zhang     for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
1274*8675203cSHong Zhang   }
1275883e35e8SHong Zhang   PetscFunctionReturn(0);
1276883e35e8SHong Zhang }
1277883e35e8SHong Zhang 
1278e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
12795cf7da58SHong Zhang {
12805cf7da58SHong Zhang   PetscErrorCode ierr;
12815cf7da58SHong Zhang   PetscInt       j;
12825cf7da58SHong Zhang   PetscScalar    val=(PetscScalar)ncols;
12835cf7da58SHong Zhang 
12845cf7da58SHong Zhang   PetscFunctionBegin;
12855cf7da58SHong Zhang   if (!ghost) {
12865cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
12875cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
12885cf7da58SHong Zhang     }
12895cf7da58SHong Zhang   } else {
12905cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
12915cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
12925cf7da58SHong Zhang     }
12935cf7da58SHong Zhang   }
12945cf7da58SHong Zhang   PetscFunctionReturn(0);
12955cf7da58SHong Zhang }
12965cf7da58SHong Zhang 
1297e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
12985cf7da58SHong Zhang {
12995cf7da58SHong Zhang   PetscErrorCode ierr;
13005cf7da58SHong Zhang   PetscInt       j,ncols_u;
13015cf7da58SHong Zhang   PetscScalar    val;
13025cf7da58SHong Zhang 
13035cf7da58SHong Zhang   PetscFunctionBegin;
13045cf7da58SHong Zhang   if (!ghost) {
13055cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
13065cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
13075cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
13085cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
13095cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
13105cf7da58SHong Zhang     }
13115cf7da58SHong Zhang   } else {
13125cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
13135cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
13145cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
13155cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
13165cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
13175cf7da58SHong Zhang     }
13185cf7da58SHong Zhang   }
13195cf7da58SHong Zhang   PetscFunctionReturn(0);
13205cf7da58SHong Zhang }
13215cf7da58SHong Zhang 
1322e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
13235cf7da58SHong Zhang {
13245cf7da58SHong Zhang   PetscErrorCode ierr;
13255cf7da58SHong Zhang 
13265cf7da58SHong Zhang   PetscFunctionBegin;
13275cf7da58SHong Zhang   if (Ju) {
13285cf7da58SHong Zhang     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
13295cf7da58SHong Zhang   } else {
13305cf7da58SHong Zhang     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
13315cf7da58SHong Zhang   }
13325cf7da58SHong Zhang   PetscFunctionReturn(0);
13335cf7da58SHong Zhang }
13345cf7da58SHong Zhang 
1335e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1336883e35e8SHong Zhang {
1337883e35e8SHong Zhang   PetscErrorCode ierr;
1338883e35e8SHong Zhang   PetscInt       j,*cols;
1339883e35e8SHong Zhang   PetscScalar    *zeros;
1340883e35e8SHong Zhang 
1341883e35e8SHong Zhang   PetscFunctionBegin;
1342883e35e8SHong Zhang   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
1343883e35e8SHong Zhang   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1344883e35e8SHong Zhang   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
1345883e35e8SHong Zhang   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
13461ad426b7SHong Zhang   PetscFunctionReturn(0);
13471ad426b7SHong Zhang }
1348a4e85ca8SHong Zhang 
1349e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
13503e97b6e8SHong Zhang {
13513e97b6e8SHong Zhang   PetscErrorCode ierr;
13523e97b6e8SHong Zhang   PetscInt       j,M,N,row,col,ncols_u;
13533e97b6e8SHong Zhang   const PetscInt *cols;
13543e97b6e8SHong Zhang   PetscScalar    zero=0.0;
13553e97b6e8SHong Zhang 
13563e97b6e8SHong Zhang   PetscFunctionBegin;
13573e97b6e8SHong Zhang   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
13583e97b6e8SHong 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);
13593e97b6e8SHong Zhang 
13603e97b6e8SHong Zhang   for (row=0; row<nrows; row++) {
13613e97b6e8SHong Zhang     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
13623e97b6e8SHong Zhang     for (j=0; j<ncols_u; j++) {
13633e97b6e8SHong Zhang       col = cols[j] + cstart;
13643e97b6e8SHong Zhang       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
13653e97b6e8SHong Zhang     }
13663e97b6e8SHong Zhang     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
13673e97b6e8SHong Zhang   }
13683e97b6e8SHong Zhang   PetscFunctionReturn(0);
13693e97b6e8SHong Zhang }
13701ad426b7SHong Zhang 
1371e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1372a4e85ca8SHong Zhang {
1373a4e85ca8SHong Zhang   PetscErrorCode ierr;
1374f4431b8cSHong Zhang 
1375a4e85ca8SHong Zhang   PetscFunctionBegin;
1376a4e85ca8SHong Zhang   if (Ju) {
1377a4e85ca8SHong Zhang     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1378a4e85ca8SHong Zhang   } else {
1379a4e85ca8SHong Zhang     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1380a4e85ca8SHong Zhang   }
1381a4e85ca8SHong Zhang   PetscFunctionReturn(0);
1382a4e85ca8SHong Zhang }
1383a4e85ca8SHong Zhang 
138424121865SAdrian 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.
138524121865SAdrian Maldonado */
138624121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
138724121865SAdrian Maldonado {
138824121865SAdrian Maldonado   PetscErrorCode ierr;
138924121865SAdrian Maldonado   PetscInt       i, size, dof;
139024121865SAdrian Maldonado   PetscInt       *glob2loc;
139124121865SAdrian Maldonado 
139224121865SAdrian Maldonado   PetscFunctionBegin;
139324121865SAdrian Maldonado   ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr);
139424121865SAdrian Maldonado   ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr);
139524121865SAdrian Maldonado 
139624121865SAdrian Maldonado   for (i = 0; i < size; i++) {
139724121865SAdrian Maldonado     ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr);
139824121865SAdrian Maldonado     dof = (dof >= 0) ? dof : -(dof + 1);
139924121865SAdrian Maldonado     glob2loc[i] = dof;
140024121865SAdrian Maldonado   }
140124121865SAdrian Maldonado 
140224121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
140324121865SAdrian Maldonado #if 0
140424121865SAdrian Maldonado   ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
140524121865SAdrian Maldonado #endif
140624121865SAdrian Maldonado   PetscFunctionReturn(0);
140724121865SAdrian Maldonado }
140824121865SAdrian Maldonado 
140901ad2aeeSHong Zhang #include <petsc/private/matimpl.h>
14101ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
14111ad426b7SHong Zhang {
14121ad426b7SHong Zhang   PetscErrorCode ierr;
141324121865SAdrian Maldonado   PetscMPIInt    rank, size;
14141ad426b7SHong Zhang   DM_Network     *network = (DM_Network*) dm->data;
1415a4e85ca8SHong Zhang   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1416840c2264SHong Zhang   PetscInt       cstart,ncols,j,e,v;
141724121865SAdrian Maldonado   PetscBool      ghost,ghost_vc,ghost2,isNest;
1418a4e85ca8SHong Zhang   Mat            Juser;
1419bfbc38dcSHong Zhang   PetscSection   sectionGlobal;
1420447d78afSSatish Balay   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
1421a4e85ca8SHong Zhang   const PetscInt *edges,*cone;
14225cf7da58SHong Zhang   MPI_Comm       comm;
142324121865SAdrian Maldonado   MatType        mtype;
14245cf7da58SHong Zhang   Vec            vd_nz,vo_nz;
14255cf7da58SHong Zhang   PetscInt       *dnnz,*onnz;
14265cf7da58SHong Zhang   PetscScalar    *vdnz,*vonz;
14271ad426b7SHong Zhang 
14281ad426b7SHong Zhang   PetscFunctionBegin;
142924121865SAdrian Maldonado   mtype = dm->mattype;
143024121865SAdrian Maldonado   ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr);
143124121865SAdrian Maldonado 
143224121865SAdrian Maldonado   if (isNest) {
14330731d606SHong Zhang     /* ierr = DMCreateMatrix_Network_Nest(); */
143424121865SAdrian Maldonado     PetscInt   eDof, vDof;
143524121865SAdrian Maldonado     Mat        j11, j12, j21, j22, bA[2][2];
143624121865SAdrian Maldonado     ISLocalToGlobalMapping eISMap, vISMap;
143724121865SAdrian Maldonado 
143824121865SAdrian Maldonado     ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
143924121865SAdrian Maldonado     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
144024121865SAdrian Maldonado     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
144124121865SAdrian Maldonado 
144224121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr);
144324121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr);
144424121865SAdrian Maldonado 
144501ad2aeeSHong Zhang     ierr = MatCreate(comm, &j11);CHKERRQ(ierr);
144624121865SAdrian Maldonado     ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
144724121865SAdrian Maldonado     ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr);
144824121865SAdrian Maldonado 
144901ad2aeeSHong Zhang     ierr = MatCreate(comm, &j12);CHKERRQ(ierr);
145024121865SAdrian Maldonado     ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr);
145124121865SAdrian Maldonado     ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr);
145224121865SAdrian Maldonado 
145301ad2aeeSHong Zhang     ierr = MatCreate(comm, &j21);CHKERRQ(ierr);
145424121865SAdrian Maldonado     ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
145524121865SAdrian Maldonado     ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr);
145624121865SAdrian Maldonado 
145701ad2aeeSHong Zhang     ierr = MatCreate(comm, &j22);CHKERRQ(ierr);
145824121865SAdrian Maldonado     ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
145924121865SAdrian Maldonado     ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr);
146024121865SAdrian Maldonado 
14613f6a6bdaSHong Zhang     bA[0][0] = j11;
14623f6a6bdaSHong Zhang     bA[0][1] = j12;
14633f6a6bdaSHong Zhang     bA[1][0] = j21;
14643f6a6bdaSHong Zhang     bA[1][1] = j22;
146524121865SAdrian Maldonado 
146624121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr);
146724121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr);
146824121865SAdrian Maldonado 
146924121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr);
147024121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr);
147124121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr);
147224121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr);
147324121865SAdrian Maldonado 
147424121865SAdrian Maldonado     ierr = MatSetUp(j11);CHKERRQ(ierr);
147524121865SAdrian Maldonado     ierr = MatSetUp(j12);CHKERRQ(ierr);
147624121865SAdrian Maldonado     ierr = MatSetUp(j21);CHKERRQ(ierr);
147724121865SAdrian Maldonado     ierr = MatSetUp(j22);CHKERRQ(ierr);
147824121865SAdrian Maldonado 
147901ad2aeeSHong Zhang     ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr);
148024121865SAdrian Maldonado     ierr = MatSetUp(*J);CHKERRQ(ierr);
148124121865SAdrian Maldonado     ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr);
148224121865SAdrian Maldonado     ierr = MatDestroy(&j11);CHKERRQ(ierr);
148324121865SAdrian Maldonado     ierr = MatDestroy(&j12);CHKERRQ(ierr);
148424121865SAdrian Maldonado     ierr = MatDestroy(&j21);CHKERRQ(ierr);
148524121865SAdrian Maldonado     ierr = MatDestroy(&j22);CHKERRQ(ierr);
148624121865SAdrian Maldonado 
148724121865SAdrian Maldonado     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
148824121865SAdrian Maldonado     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
148924121865SAdrian Maldonado 
149024121865SAdrian Maldonado     /* Free structures */
149124121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr);
149224121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr);
149324121865SAdrian Maldonado 
149424121865SAdrian Maldonado     PetscFunctionReturn(0);
149524121865SAdrian Maldonado   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1496a4e85ca8SHong Zhang     /* user does not provide Jacobian blocks */
1497bfbc38dcSHong Zhang     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1498bfbc38dcSHong Zhang     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
14991ad426b7SHong Zhang     PetscFunctionReturn(0);
15001ad426b7SHong Zhang   }
15011ad426b7SHong Zhang 
1502bfbc38dcSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
15032a945128SHong Zhang   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1504bfbc38dcSHong Zhang   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1505bfbc38dcSHong Zhang   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
15062a945128SHong Zhang 
15072a945128SHong Zhang   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
15082a945128SHong Zhang   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
150989898e50SHong Zhang 
151089898e50SHong Zhang   /* (1) Set matrix preallocation */
151189898e50SHong Zhang   /*------------------------------*/
1512840c2264SHong Zhang   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1513840c2264SHong Zhang   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1514840c2264SHong Zhang   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1515840c2264SHong Zhang   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1516840c2264SHong Zhang   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1517840c2264SHong Zhang   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1518840c2264SHong Zhang 
151989898e50SHong Zhang   /* Set preallocation for edges */
152089898e50SHong Zhang   /*-----------------------------*/
1521840c2264SHong Zhang   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1522840c2264SHong Zhang 
1523bdcb62a2SHong Zhang   ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr);
1524840c2264SHong Zhang   for (e=eStart; e<eEnd; e++) {
1525840c2264SHong Zhang     /* Get row indices */
1526840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1527840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1528840c2264SHong Zhang     if (nrows) {
1529840c2264SHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1530840c2264SHong Zhang 
15315cf7da58SHong Zhang       /* Set preallocation for conntected vertices */
1532d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1533840c2264SHong Zhang       for (v=0; v<2; v++) {
1534840c2264SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1535840c2264SHong Zhang 
1536*8675203cSHong Zhang         if (network->Je) {
1537840c2264SHong Zhang           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1538*8675203cSHong Zhang         } else Juser = NULL;
1539840c2264SHong Zhang         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
15405cf7da58SHong Zhang         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1541840c2264SHong Zhang       }
1542840c2264SHong Zhang 
154389898e50SHong Zhang       /* Set preallocation for edge self */
1544840c2264SHong Zhang       cstart = rstart;
1545*8675203cSHong Zhang       if (network->Je) {
1546840c2264SHong Zhang         Juser = network->Je[3*e]; /* Jacobian(e,e) */
1547*8675203cSHong Zhang       } else Juser = NULL;
15485cf7da58SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1549840c2264SHong Zhang     }
1550840c2264SHong Zhang   }
1551840c2264SHong Zhang 
155289898e50SHong Zhang   /* Set preallocation for vertices */
155389898e50SHong Zhang   /*--------------------------------*/
1554840c2264SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1555*8675203cSHong Zhang   if (vEnd - vStart) vptr = network->Jvptr;
1556840c2264SHong Zhang 
1557840c2264SHong Zhang   for (v=vStart; v<vEnd; v++) {
1558840c2264SHong Zhang     /* Get row indices */
1559840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1560840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1561840c2264SHong Zhang     if (!nrows) continue;
1562840c2264SHong Zhang 
1563bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1564bdcb62a2SHong Zhang     if (ghost) {
1565bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1566bdcb62a2SHong Zhang     } else {
1567bdcb62a2SHong Zhang       rows_v = rows;
1568bdcb62a2SHong Zhang     }
1569bdcb62a2SHong Zhang 
1570bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1571840c2264SHong Zhang 
1572840c2264SHong Zhang     /* Get supporting edges and connected vertices */
1573840c2264SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1574840c2264SHong Zhang 
1575840c2264SHong Zhang     for (e=0; e<nedges; e++) {
1576840c2264SHong Zhang       /* Supporting edges */
1577840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1578840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1579840c2264SHong Zhang 
1580*8675203cSHong Zhang       if (network->Jv) {
1581840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1582*8675203cSHong Zhang       } else Juser = NULL;
1583bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1584840c2264SHong Zhang 
1585840c2264SHong Zhang       /* Connected vertices */
1586d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
1587840c2264SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
1588840c2264SHong Zhang       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1589840c2264SHong Zhang 
1590840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1591840c2264SHong Zhang 
1592*8675203cSHong Zhang       if (network->Jv) {
1593840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1594*8675203cSHong Zhang       } else Juser = NULL;
1595e102a522SHong Zhang       if (ghost_vc||ghost) {
1596e102a522SHong Zhang         ghost2 = PETSC_TRUE;
1597e102a522SHong Zhang       } else {
1598e102a522SHong Zhang         ghost2 = PETSC_FALSE;
1599e102a522SHong Zhang       }
1600e102a522SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr);
1601840c2264SHong Zhang     }
1602840c2264SHong Zhang 
160389898e50SHong Zhang     /* Set preallocation for vertex self */
1604840c2264SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1605840c2264SHong Zhang     if (!ghost) {
1606840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1607*8675203cSHong Zhang       if (network->Jv) {
1608840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1609*8675203cSHong Zhang       } else Juser = NULL;
1610bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1611840c2264SHong Zhang     }
1612bdcb62a2SHong Zhang     if (ghost) {
1613bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1614bdcb62a2SHong Zhang     }
1615840c2264SHong Zhang   }
1616840c2264SHong Zhang 
1617840c2264SHong Zhang   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1618840c2264SHong Zhang   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
16195cf7da58SHong Zhang 
16205cf7da58SHong Zhang   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
16215cf7da58SHong Zhang 
16225cf7da58SHong Zhang   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1623840c2264SHong Zhang   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1624840c2264SHong Zhang 
1625840c2264SHong Zhang   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1626840c2264SHong Zhang   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1627840c2264SHong Zhang   for (j=0; j<localSize; j++) {
1628e102a522SHong Zhang     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1629e102a522SHong Zhang     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1630840c2264SHong Zhang   }
1631840c2264SHong Zhang   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1632840c2264SHong Zhang   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1633840c2264SHong Zhang   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1634840c2264SHong Zhang   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1635840c2264SHong Zhang 
16365cf7da58SHong Zhang   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
16375cf7da58SHong Zhang   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
16385cf7da58SHong Zhang   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
16395cf7da58SHong Zhang 
16405cf7da58SHong Zhang   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
16415cf7da58SHong Zhang 
164289898e50SHong Zhang   /* (2) Set matrix entries for edges */
164389898e50SHong Zhang   /*----------------------------------*/
16441ad426b7SHong Zhang   for (e=eStart; e<eEnd; e++) {
1645bfbc38dcSHong Zhang     /* Get row indices */
16461ad426b7SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
164717df6e9eSHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
16484b976069SHong Zhang     if (nrows) {
164917df6e9eSHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
16501ad426b7SHong Zhang 
1651bfbc38dcSHong Zhang       /* Set matrix entries for conntected vertices */
1652d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1653bfbc38dcSHong Zhang       for (v=0; v<2; v++) {
1654bfbc38dcSHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1655883e35e8SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
16563e97b6e8SHong Zhang 
1657*8675203cSHong Zhang         if (network->Je) {
1658a4e85ca8SHong Zhang           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1659*8675203cSHong Zhang         } else Juser = NULL;
1660a4e85ca8SHong Zhang         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1661bfbc38dcSHong Zhang       }
166217df6e9eSHong Zhang 
1663bfbc38dcSHong Zhang       /* Set matrix entries for edge self */
16643e97b6e8SHong Zhang       cstart = rstart;
1665*8675203cSHong Zhang       if (network->Je) {
1666a4e85ca8SHong Zhang         Juser = network->Je[3*e]; /* Jacobian(e,e) */
1667*8675203cSHong Zhang       } else Juser = NULL;
1668a4e85ca8SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
16691ad426b7SHong Zhang     }
16704b976069SHong Zhang   }
16711ad426b7SHong Zhang 
1672bfbc38dcSHong Zhang   /* Set matrix entries for vertices */
167383b2e829SHong Zhang   /*---------------------------------*/
16741ad426b7SHong Zhang   for (v=vStart; v<vEnd; v++) {
1675bfbc38dcSHong Zhang     /* Get row indices */
1676596e729fSHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1677596e729fSHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
16784b976069SHong Zhang     if (!nrows) continue;
1679596e729fSHong Zhang 
1680bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1681bdcb62a2SHong Zhang     if (ghost) {
1682bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1683bdcb62a2SHong Zhang     } else {
1684bdcb62a2SHong Zhang       rows_v = rows;
1685bdcb62a2SHong Zhang     }
1686bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1687596e729fSHong Zhang 
1688bfbc38dcSHong Zhang     /* Get supporting edges and connected vertices */
1689596e729fSHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1690596e729fSHong Zhang 
1691596e729fSHong Zhang     for (e=0; e<nedges; e++) {
1692bfbc38dcSHong Zhang       /* Supporting edges */
1693596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1694596e729fSHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1695596e729fSHong Zhang 
1696*8675203cSHong Zhang       if (network->Jv) {
1697a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1698*8675203cSHong Zhang       } else Juser = NULL;
1699bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1700596e729fSHong Zhang 
1701bfbc38dcSHong Zhang       /* Connected vertices */
1702d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
17032a945128SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
17042a945128SHong Zhang 
170544aca652SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
170644aca652SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1707a4e85ca8SHong Zhang 
1708*8675203cSHong Zhang       if (network->Jv) {
1709a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1710*8675203cSHong Zhang       } else Juser = NULL;
1711bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1712596e729fSHong Zhang     }
1713596e729fSHong Zhang 
1714bfbc38dcSHong Zhang     /* Set matrix entries for vertex self */
17151ad426b7SHong Zhang     if (!ghost) {
1716596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1717*8675203cSHong Zhang       if (network->Jv) {
1718a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1719*8675203cSHong Zhang       } else Juser = NULL;
1720bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr);
1721bdcb62a2SHong Zhang     }
1722bdcb62a2SHong Zhang     if (ghost) {
1723bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1724bdcb62a2SHong Zhang     }
17251ad426b7SHong Zhang   }
1726a4e85ca8SHong Zhang   ierr = PetscFree(rows);CHKERRQ(ierr);
1727bdcb62a2SHong Zhang 
17281ad426b7SHong Zhang   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17291ad426b7SHong Zhang   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1730dd6f46cdSHong Zhang 
17315f2c45f1SShri Abhyankar   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
17325f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
17335f2c45f1SShri Abhyankar }
17345f2c45f1SShri Abhyankar 
17355f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm)
17365f2c45f1SShri Abhyankar {
17375f2c45f1SShri Abhyankar   PetscErrorCode ierr;
17385f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
17392727e31bSShri Abhyankar   PetscInt       j;
17405f2c45f1SShri Abhyankar 
17415f2c45f1SShri Abhyankar   PetscFunctionBegin;
17428415c774SShri Abhyankar   if (--network->refct > 0) PetscFunctionReturn(0);
174383b2e829SHong Zhang   if (network->Je) {
174483b2e829SHong Zhang     ierr = PetscFree(network->Je);CHKERRQ(ierr);
174583b2e829SHong Zhang   }
174683b2e829SHong Zhang   if (network->Jv) {
1747883e35e8SHong Zhang     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
174883b2e829SHong Zhang     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
17491ad426b7SHong Zhang   }
175013c2a604SAdrian Maldonado 
175113c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr);
175213c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr);
175313c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr);
175413c2a604SAdrian Maldonado   if (network->vertex.sf) {
175513c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr);
175613c2a604SAdrian Maldonado   }
175713c2a604SAdrian Maldonado   /* edge */
175813c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr);
175913c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr);
176013c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr);
176113c2a604SAdrian Maldonado   if (network->edge.sf) {
176213c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr);
176313c2a604SAdrian Maldonado   }
17645f2c45f1SShri Abhyankar   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
17655f2c45f1SShri Abhyankar   network->edges = NULL;
17665f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
17675f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
176883b2e829SHong Zhang 
17692727e31bSShri Abhyankar   for(j=0; j < network->nsubnet; j++) {
17702727e31bSShri Abhyankar     ierr = PetscFree(network->subnet[j].edges);CHKERRQ(ierr);
17712727e31bSShri Abhyankar     ierr = PetscFree(network->subnet[j].vertices);CHKERRQ(ierr);
17722727e31bSShri Abhyankar   }
1773e2aaf10cSShri Abhyankar   ierr = PetscFree(network->subnet);CHKERRQ(ierr);
17745f2c45f1SShri Abhyankar   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
17755f2c45f1SShri Abhyankar   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
17765f2c45f1SShri Abhyankar   ierr = PetscFree(network->header);CHKERRQ(ierr);
17775f2c45f1SShri Abhyankar   ierr = PetscFree(network);CHKERRQ(ierr);
17785f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
17795f2c45f1SShri Abhyankar }
17805f2c45f1SShri Abhyankar 
17815f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
17825f2c45f1SShri Abhyankar {
17835f2c45f1SShri Abhyankar   PetscErrorCode ierr;
17845f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
17855f2c45f1SShri Abhyankar 
17865f2c45f1SShri Abhyankar   PetscFunctionBegin;
17875f2c45f1SShri Abhyankar   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
17885f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
17895f2c45f1SShri Abhyankar }
17905f2c45f1SShri Abhyankar 
17915f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
17925f2c45f1SShri Abhyankar {
17935f2c45f1SShri Abhyankar   PetscErrorCode ierr;
17945f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
17955f2c45f1SShri Abhyankar 
17965f2c45f1SShri Abhyankar   PetscFunctionBegin;
17975f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
17985f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
17995f2c45f1SShri Abhyankar }
18005f2c45f1SShri Abhyankar 
18015f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
18025f2c45f1SShri Abhyankar {
18035f2c45f1SShri Abhyankar   PetscErrorCode ierr;
18045f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
18055f2c45f1SShri Abhyankar 
18065f2c45f1SShri Abhyankar   PetscFunctionBegin;
18075f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
18085f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
18095f2c45f1SShri Abhyankar }
18105f2c45f1SShri Abhyankar 
18115f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
18125f2c45f1SShri Abhyankar {
18135f2c45f1SShri Abhyankar   PetscErrorCode ierr;
18145f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
18155f2c45f1SShri Abhyankar 
18165f2c45f1SShri Abhyankar   PetscFunctionBegin;
18175f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
18185f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
18195f2c45f1SShri Abhyankar }
18205f2c45f1SShri Abhyankar 
18215f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
18225f2c45f1SShri Abhyankar {
18235f2c45f1SShri Abhyankar   PetscErrorCode ierr;
18245f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
18255f2c45f1SShri Abhyankar 
18265f2c45f1SShri Abhyankar   PetscFunctionBegin;
18275f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
18285f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
18295f2c45f1SShri Abhyankar }
1830