xref: /petsc/src/dm/impls/network/network.c (revision e2aaf10c0406ebc0742d360f5ef07fdd39120198)
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 /*@
28*e2aaf10cSShri 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
34*e2aaf10cSShri Abhyankar . Nsubnet - number of subnetworks
35*e2aaf10cSShri Abhyankar . nV - number of local vertices for each subnetwork
36*e2aaf10cSShri Abhyankar . nE - number of local edges for each subnetwork
37*e2aaf10cSShri Abhyankar . NV - number of global vertices (or PETSC_DETERMINE) for each subnetwork
38*e2aaf10cSShri 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 @*/
49*e2aaf10cSShri 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;
53*e2aaf10cSShri Abhyankar   PetscInt       a[2],b[2],i;
545f2c45f1SShri Abhyankar 
555f2c45f1SShri Abhyankar   PetscFunctionBegin;
565f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
57*e2aaf10cSShri Abhyankar   if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet);
58*e2aaf10cSShri Abhyankar 
59*e2aaf10cSShri Abhyankar   if(Nsubnet > 0) PetscValidLogicalCollectiveInt(dm,Nsubnet,2);
60*e2aaf10cSShri Abhyankar   if(network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network");
61*e2aaf10cSShri Abhyankar 
62*e2aaf10cSShri Abhyankar   network->nsubnet = Nsubnet;
63*e2aaf10cSShri Abhyankar   ierr = PetscCalloc1(Nsubnet,&network->subnet);CHKERRQ(ierr);
64*e2aaf10cSShri Abhyankar   for(i=0; i < network->nsubnet; i++) {
65*e2aaf10cSShri Abhyankar     if (NV[i] > 0) PetscValidLogicalCollectiveInt(dm,NV[i],5);
66*e2aaf10cSShri Abhyankar     if (NE[i] > 0) PetscValidLogicalCollectiveInt(dm,NE[i],6);
67*e2aaf10cSShri 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]);
68*e2aaf10cSShri 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]);
69*e2aaf10cSShri 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);
71*e2aaf10cSShri Abhyankar     network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1];
72*e2aaf10cSShri Abhyankar 
73*e2aaf10cSShri Abhyankar     network->subnet[i].id = i;
74*e2aaf10cSShri Abhyankar 
75*e2aaf10cSShri Abhyankar     network->subnet[i].nvtx = nV[i];
76*e2aaf10cSShri Abhyankar     network->subnet[i].vStart = network->nVertices;
77*e2aaf10cSShri Abhyankar     network->subnet[i].vEnd   = network->subnet[i].vStart + network->subnet[i].Nvtx;
78*e2aaf10cSShri Abhyankar     network->nVertices += network->subnet[i].nvtx;
79*e2aaf10cSShri Abhyankar     network->NVertices += network->subnet[i].Nvtx;
80*e2aaf10cSShri Abhyankar 
81*e2aaf10cSShri Abhyankar     network->subnet[i].nedge = nE[i];
82*e2aaf10cSShri Abhyankar     network->subnet[i].eStart = network->nEdges;
83*e2aaf10cSShri Abhyankar     network->subnet[i].eEnd = network->subnet[i].eStart + network->subnet[i].Nedge;
84*e2aaf10cSShri Abhyankar     network->nEdges += network->subnet[i].nedge;
85*e2aaf10cSShri 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:
96*e2aaf10cSShri 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 @*/
106*e2aaf10cSShri Abhyankar PetscErrorCode DMNetworkSetEdgeList(DM dm, int *edgelist[])
1075f2c45f1SShri Abhyankar {
1085f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network*) dm->data;
109*e2aaf10cSShri Abhyankar   PetscInt   i;
1105f2c45f1SShri Abhyankar 
1115f2c45f1SShri Abhyankar   PetscFunctionBegin;
112*e2aaf10cSShri Abhyankar   for(i=0; i < network->nsubnet; i++) {
113*e2aaf10cSShri Abhyankar     network->subnet[i].edgelist = edgelist[i];
114*e2aaf10cSShri 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;
144*e2aaf10cSShri Abhyankar   PetscInt       i,j;
1455f2c45f1SShri Abhyankar   PetscInt       ndata;
146*e2aaf10cSShri Abhyankar   PetscInt       ctr=0;
1475f2c45f1SShri Abhyankar 
1485f2c45f1SShri Abhyankar   PetscFunctionBegin;
149*e2aaf10cSShri Abhyankar 
1506fefedf4SHong Zhang   if (network->nVertices) {
1516fefedf4SHong Zhang     ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr);
1525f2c45f1SShri Abhyankar   }
153*e2aaf10cSShri Abhyankar 
154*e2aaf10cSShri Abhyankar   /* Create the edgelist for the network by concatenating edgelists of the subnetworks */
155*e2aaf10cSShri Abhyankar   ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr);
156*e2aaf10cSShri Abhyankar   for(i=0; i < network->nsubnet; i++) {
157*e2aaf10cSShri Abhyankar     for(j = 0; j < network->subnet[i].nedge; j++) {
158*e2aaf10cSShri Abhyankar       network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
159*e2aaf10cSShri Abhyankar       network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
160*e2aaf10cSShri Abhyankar       ctr++;
161*e2aaf10cSShri Abhyankar     }
162*e2aaf10cSShri Abhyankar   }
163*e2aaf10cSShri Abhyankar 
164*e2aaf10cSShri Abhyankar #if 0
165*e2aaf10cSShri Abhyankar   for(i=0; i < network->nEdges; i++) {
166*e2aaf10cSShri Abhyankar     ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr);
167*e2aaf10cSShri Abhyankar   }
168*e2aaf10cSShri Abhyankar #endif
169*e2aaf10cSShri 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   }
174*e2aaf10cSShri Abhyankar   ierr = PetscFree(network->edges);CHKERRQ(ierr);
175*e2aaf10cSShri 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 
1855f2c45f1SShri Abhyankar   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
1866caa05f4SBarry Smith   ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr);
187*e2aaf10cSShri Abhyankar   for(i=network->eStart; i < network->eEnd; i++) {
188*e2aaf10cSShri Abhyankar     network->header[i].index = i;   /* Global edge number */
189*e2aaf10cSShri Abhyankar     for(j=0; j < network->nsubnet; j++) {
190*e2aaf10cSShri Abhyankar       if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) {
191*e2aaf10cSShri Abhyankar 	network->header[i].subnetid = j; /* Subnetwork id */
192*e2aaf10cSShri Abhyankar 	break;
193*e2aaf10cSShri Abhyankar       }
1947b6afd5bSHong Zhang     }
1955f2c45f1SShri Abhyankar     network->header[i].ndata = 0;
1965f2c45f1SShri Abhyankar     ndata = network->header[i].ndata;
1975f2c45f1SShri Abhyankar     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
1985f2c45f1SShri Abhyankar     network->header[i].offset[ndata] = 0;
1995f2c45f1SShri Abhyankar   }
200*e2aaf10cSShri Abhyankar 
201*e2aaf10cSShri Abhyankar   for(i=network->vStart; i < network->vEnd; i++) {
202*e2aaf10cSShri Abhyankar     network->header[i].index = i - network->vStart;
203*e2aaf10cSShri Abhyankar     for(j=0; j < network->nsubnet; j++) {
204*e2aaf10cSShri Abhyankar       if((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) {
205*e2aaf10cSShri Abhyankar 	network->header[i].subnetid = j;
206*e2aaf10cSShri Abhyankar 	break;
207*e2aaf10cSShri Abhyankar       }
208*e2aaf10cSShri Abhyankar     }
209*e2aaf10cSShri Abhyankar     network->header[i].ndata = 0;
210*e2aaf10cSShri Abhyankar     ndata = network->header[i].ndata;
211*e2aaf10cSShri Abhyankar     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
212*e2aaf10cSShri Abhyankar     network->header[i].offset[ndata] = 0;
213*e2aaf10cSShri Abhyankar   }
214*e2aaf10cSShri Abhyankar 
215854ce69bSBarry Smith   ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr);
216*e2aaf10cSShri Abhyankar 
2175f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2185f2c45f1SShri Abhyankar }
2195f2c45f1SShri Abhyankar 
22094ef8ddeSSatish Balay /*@C
2215f2c45f1SShri Abhyankar   DMNetworkRegisterComponent - Registers the network component
2225f2c45f1SShri Abhyankar 
2235f2c45f1SShri Abhyankar   Logically collective on DM
2245f2c45f1SShri Abhyankar 
2255f2c45f1SShri Abhyankar   Input Parameters
2265f2c45f1SShri Abhyankar + dm   - the network object
2275f2c45f1SShri Abhyankar . name - the component name
2285f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data
2295f2c45f1SShri Abhyankar 
2305f2c45f1SShri Abhyankar    Output Parameters
2315f2c45f1SShri Abhyankar .   key - an integer key that defines the component
2325f2c45f1SShri Abhyankar 
2335f2c45f1SShri Abhyankar    Notes
2345f2c45f1SShri Abhyankar    This routine should be called by all processors before calling DMNetworkLayoutSetup().
2355f2c45f1SShri Abhyankar 
2365f2c45f1SShri Abhyankar    Level: intermediate
2375f2c45f1SShri Abhyankar 
2385f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
2395f2c45f1SShri Abhyankar @*/
2405f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key)
2415f2c45f1SShri Abhyankar {
2425f2c45f1SShri Abhyankar   PetscErrorCode        ierr;
2435f2c45f1SShri Abhyankar   DM_Network            *network = (DM_Network*) dm->data;
2445f2c45f1SShri Abhyankar   DMNetworkComponent    *component=&network->component[network->ncomponent];
2455f2c45f1SShri Abhyankar   PetscBool             flg=PETSC_FALSE;
2465f2c45f1SShri Abhyankar   PetscInt              i;
2475f2c45f1SShri Abhyankar 
2485f2c45f1SShri Abhyankar   PetscFunctionBegin;
2495f2c45f1SShri Abhyankar 
2505f2c45f1SShri Abhyankar   for (i=0; i < network->ncomponent; i++) {
2515f2c45f1SShri Abhyankar     ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr);
2525f2c45f1SShri Abhyankar     if (flg) {
2535f2c45f1SShri Abhyankar       *key = i;
2545f2c45f1SShri Abhyankar       PetscFunctionReturn(0);
2555f2c45f1SShri Abhyankar     }
2565f2c45f1SShri Abhyankar   }
2575f2c45f1SShri Abhyankar 
2585f2c45f1SShri Abhyankar   ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr);
2595f2c45f1SShri Abhyankar   component->size = size/sizeof(DMNetworkComponentGenericDataType);
2605f2c45f1SShri Abhyankar   *key = network->ncomponent;
2615f2c45f1SShri Abhyankar   network->ncomponent++;
2625f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2635f2c45f1SShri Abhyankar }
2645f2c45f1SShri Abhyankar 
2655f2c45f1SShri Abhyankar /*@
2665f2c45f1SShri Abhyankar   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.
2675f2c45f1SShri Abhyankar 
2685f2c45f1SShri Abhyankar   Not Collective
2695f2c45f1SShri Abhyankar 
2705f2c45f1SShri Abhyankar   Input Parameters:
2715f2c45f1SShri Abhyankar + dm - The DMNetwork object
2725f2c45f1SShri Abhyankar 
2735f2c45f1SShri Abhyankar   Output Paramters:
2745f2c45f1SShri Abhyankar + vStart - The first vertex point
2755f2c45f1SShri Abhyankar - vEnd   - One beyond the last vertex point
2765f2c45f1SShri Abhyankar 
2775f2c45f1SShri Abhyankar   Level: intermediate
2785f2c45f1SShri Abhyankar 
2795f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange
2805f2c45f1SShri Abhyankar @*/
2815f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
2825f2c45f1SShri Abhyankar {
2835f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
2845f2c45f1SShri Abhyankar 
2855f2c45f1SShri Abhyankar   PetscFunctionBegin;
2865f2c45f1SShri Abhyankar   if (vStart) *vStart = network->vStart;
2875f2c45f1SShri Abhyankar   if (vEnd) *vEnd = network->vEnd;
2885f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2895f2c45f1SShri Abhyankar }
2905f2c45f1SShri Abhyankar 
2915f2c45f1SShri Abhyankar /*@
2925f2c45f1SShri Abhyankar   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.
2935f2c45f1SShri Abhyankar 
2945f2c45f1SShri Abhyankar   Not Collective
2955f2c45f1SShri Abhyankar 
2965f2c45f1SShri Abhyankar   Input Parameters:
2975f2c45f1SShri Abhyankar + dm - The DMNetwork object
2985f2c45f1SShri Abhyankar 
2995f2c45f1SShri Abhyankar   Output Paramters:
3005f2c45f1SShri Abhyankar + eStart - The first edge point
3015f2c45f1SShri Abhyankar - eEnd   - One beyond the last edge point
3025f2c45f1SShri Abhyankar 
3035f2c45f1SShri Abhyankar   Level: intermediate
3045f2c45f1SShri Abhyankar 
3055f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange
3065f2c45f1SShri Abhyankar @*/
3075f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
3085f2c45f1SShri Abhyankar {
3095f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
3105f2c45f1SShri Abhyankar 
3115f2c45f1SShri Abhyankar   PetscFunctionBegin;
3125f2c45f1SShri Abhyankar   if (eStart) *eStart = network->eStart;
3135f2c45f1SShri Abhyankar   if (eEnd) *eEnd = network->eEnd;
3145f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3155f2c45f1SShri Abhyankar }
3165f2c45f1SShri Abhyankar 
3177b6afd5bSHong Zhang /*@
318e85e6aecSHong Zhang   DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge.
3197b6afd5bSHong Zhang 
3207b6afd5bSHong Zhang   Not Collective
3217b6afd5bSHong Zhang 
3227b6afd5bSHong Zhang   Input Parameters:
3237b6afd5bSHong Zhang + dm - DMNetwork object
324e85e6aecSHong Zhang - p  - edge point
3257b6afd5bSHong Zhang 
3267b6afd5bSHong Zhang   Output Paramters:
327e85e6aecSHong Zhang . index - user global numbering for the edge
3287b6afd5bSHong Zhang 
3297b6afd5bSHong Zhang   Level: intermediate
3307b6afd5bSHong Zhang 
331e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalVertexIndex
3327b6afd5bSHong Zhang @*/
333e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
3347b6afd5bSHong Zhang {
3357b6afd5bSHong Zhang   PetscErrorCode    ierr;
3367b6afd5bSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
3377b6afd5bSHong Zhang   PetscInt          offsetp;
3387b6afd5bSHong Zhang   DMNetworkComponentHeader header;
3397b6afd5bSHong Zhang 
3407b6afd5bSHong Zhang   PetscFunctionBegin;
3417b6afd5bSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
3427b6afd5bSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
343e85e6aecSHong Zhang   *index = header->index;
3447b6afd5bSHong Zhang   PetscFunctionReturn(0);
3457b6afd5bSHong Zhang }
3467b6afd5bSHong Zhang 
3475f2c45f1SShri Abhyankar /*@
348e85e6aecSHong Zhang   DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex.
349e85e6aecSHong Zhang 
350e85e6aecSHong Zhang   Not Collective
351e85e6aecSHong Zhang 
352e85e6aecSHong Zhang   Input Parameters:
353e85e6aecSHong Zhang + dm - DMNetwork object
354e85e6aecSHong Zhang - p  - vertex point
355e85e6aecSHong Zhang 
356e85e6aecSHong Zhang   Output Paramters:
357e85e6aecSHong Zhang . index - user global numbering for the vertex
358e85e6aecSHong Zhang 
359e85e6aecSHong Zhang   Level: intermediate
360e85e6aecSHong Zhang 
361e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalEdgeIndex
362e85e6aecSHong Zhang @*/
363e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
364e85e6aecSHong Zhang {
365e85e6aecSHong Zhang   PetscErrorCode    ierr;
366e85e6aecSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
367e85e6aecSHong Zhang   PetscInt          offsetp;
368e85e6aecSHong Zhang   DMNetworkComponentHeader header;
369e85e6aecSHong Zhang 
370e85e6aecSHong Zhang   PetscFunctionBegin;
371e85e6aecSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
372e85e6aecSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
373e85e6aecSHong Zhang   *index = header->index;
374e85e6aecSHong Zhang   PetscFunctionReturn(0);
375e85e6aecSHong Zhang }
376e85e6aecSHong Zhang 
377e85e6aecSHong Zhang /*@
378325661f6SSatish Balay   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)
3795f2c45f1SShri Abhyankar 
3805f2c45f1SShri Abhyankar   Not Collective
3815f2c45f1SShri Abhyankar 
3825f2c45f1SShri Abhyankar   Input Parameters:
3835f2c45f1SShri Abhyankar + dm           - The DMNetwork object
3845f2c45f1SShri Abhyankar . p            - vertex/edge point
3855f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component
3865f2c45f1SShri Abhyankar - compvalue    - pointer to the data structure for the component
3875f2c45f1SShri Abhyankar 
3885f2c45f1SShri Abhyankar   Level: intermediate
3895f2c45f1SShri Abhyankar 
3905f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
3915f2c45f1SShri Abhyankar @*/
3925f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
3935f2c45f1SShri Abhyankar {
3945f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
39543a39a44SBarry Smith   DMNetworkComponent       *component = &network->component[componentkey];
3965f2c45f1SShri Abhyankar   DMNetworkComponentHeader header = &network->header[p];
3975f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue = &network->cvalue[p];
3985f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
3995f2c45f1SShri Abhyankar 
4005f2c45f1SShri Abhyankar   PetscFunctionBegin;
401fa58f0a9SHong 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);
402fa58f0a9SHong Zhang 
40343a39a44SBarry Smith   header->size[header->ndata] = component->size;
40443a39a44SBarry Smith   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
4055f2c45f1SShri Abhyankar   header->key[header->ndata] = componentkey;
4065f2c45f1SShri Abhyankar   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
4075f2c45f1SShri Abhyankar 
4085f2c45f1SShri Abhyankar   cvalue->data[header->ndata] = (void*)compvalue;
4095f2c45f1SShri Abhyankar   header->ndata++;
4105f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4115f2c45f1SShri Abhyankar }
4125f2c45f1SShri Abhyankar 
4135f2c45f1SShri Abhyankar /*@
4145f2c45f1SShri Abhyankar   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
4155f2c45f1SShri Abhyankar 
4165f2c45f1SShri Abhyankar   Not Collective
4175f2c45f1SShri Abhyankar 
4185f2c45f1SShri Abhyankar   Input Parameters:
4195f2c45f1SShri Abhyankar + dm - The DMNetwork object
4205f2c45f1SShri Abhyankar . p  - vertex/edge point
4215f2c45f1SShri Abhyankar 
4225f2c45f1SShri Abhyankar   Output Parameters:
4235f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge
4245f2c45f1SShri Abhyankar 
4255f2c45f1SShri Abhyankar   Level: intermediate
4265f2c45f1SShri Abhyankar 
4275f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
4285f2c45f1SShri Abhyankar @*/
4295f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
4305f2c45f1SShri Abhyankar {
4315f2c45f1SShri Abhyankar   PetscErrorCode ierr;
4325f2c45f1SShri Abhyankar   PetscInt       offset;
4335f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4345f2c45f1SShri Abhyankar 
4355f2c45f1SShri Abhyankar   PetscFunctionBegin;
4365f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
4375f2c45f1SShri Abhyankar   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
4385f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4395f2c45f1SShri Abhyankar }
4405f2c45f1SShri Abhyankar 
4415f2c45f1SShri Abhyankar /*@
442a730d845SHong Zhang   DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the
4435f2c45f1SShri Abhyankar                                     component value from the component data array
4445f2c45f1SShri Abhyankar 
4455f2c45f1SShri Abhyankar   Not Collective
4465f2c45f1SShri Abhyankar 
4475f2c45f1SShri Abhyankar   Input Parameters:
4485f2c45f1SShri Abhyankar + dm      - The DMNetwork object
4495f2c45f1SShri Abhyankar . p       - vertex/edge point
4505f2c45f1SShri Abhyankar - compnum - component number
4515f2c45f1SShri Abhyankar 
4525f2c45f1SShri Abhyankar   Output Parameters:
4535f2c45f1SShri Abhyankar + compkey - the key obtained when registering the component
4545f2c45f1SShri Abhyankar - offset  - offset into the component data array associated with the vertex/edge point
4555f2c45f1SShri Abhyankar 
4565f2c45f1SShri Abhyankar   Notes:
4575f2c45f1SShri Abhyankar   Typical usage:
4585f2c45f1SShri Abhyankar 
4595f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray(dm, &arr);
4605f2c45f1SShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
4615f2c45f1SShri Abhyankar   Loop over vertices or edges
4625f2c45f1SShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
4635f2c45f1SShri Abhyankar     Loop over numcomps
464a730d845SHong Zhang       DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset);
4655f2c45f1SShri Abhyankar       compdata = (UserCompDataType)(arr+offset);
4665f2c45f1SShri Abhyankar 
4675f2c45f1SShri Abhyankar   Level: intermediate
4685f2c45f1SShri Abhyankar 
4695f2c45f1SShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
4705f2c45f1SShri Abhyankar @*/
471a730d845SHong Zhang PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
4725f2c45f1SShri Abhyankar {
4735f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
4745f2c45f1SShri Abhyankar   PetscInt                 offsetp;
4755f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
4765f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
4775f2c45f1SShri Abhyankar 
4785f2c45f1SShri Abhyankar   PetscFunctionBegin;
4795f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
4805f2c45f1SShri Abhyankar   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
481d36f4e81SHong Zhang   if (compkey) *compkey = header->key[compnum];
482d36f4e81SHong Zhang   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
4835f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4845f2c45f1SShri Abhyankar }
4855f2c45f1SShri Abhyankar 
4865f2c45f1SShri Abhyankar /*@
4875f2c45f1SShri Abhyankar   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
4885f2c45f1SShri Abhyankar 
4895f2c45f1SShri Abhyankar   Not Collective
4905f2c45f1SShri Abhyankar 
4915f2c45f1SShri Abhyankar   Input Parameters:
4925f2c45f1SShri Abhyankar + dm     - The DMNetwork object
4935f2c45f1SShri Abhyankar - p      - the edge/vertex point
4945f2c45f1SShri Abhyankar 
4955f2c45f1SShri Abhyankar   Output Parameters:
4965f2c45f1SShri Abhyankar . offset - the offset
4975f2c45f1SShri Abhyankar 
4985f2c45f1SShri Abhyankar   Level: intermediate
4995f2c45f1SShri Abhyankar 
5005f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
5015f2c45f1SShri Abhyankar @*/
5025f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
5035f2c45f1SShri Abhyankar {
5045f2c45f1SShri Abhyankar   PetscErrorCode ierr;
5055f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5065f2c45f1SShri Abhyankar 
5075f2c45f1SShri Abhyankar   PetscFunctionBegin;
5085f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr);
5095f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5105f2c45f1SShri Abhyankar }
5115f2c45f1SShri Abhyankar 
5125f2c45f1SShri Abhyankar /*@
5135f2c45f1SShri Abhyankar   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
5145f2c45f1SShri Abhyankar 
5155f2c45f1SShri Abhyankar   Not Collective
5165f2c45f1SShri Abhyankar 
5175f2c45f1SShri Abhyankar   Input Parameters:
5185f2c45f1SShri Abhyankar + dm      - The DMNetwork object
5195f2c45f1SShri Abhyankar - p       - the edge/vertex point
5205f2c45f1SShri Abhyankar 
5215f2c45f1SShri Abhyankar   Output Parameters:
5225f2c45f1SShri Abhyankar . offsetg - the offset
5235f2c45f1SShri Abhyankar 
5245f2c45f1SShri Abhyankar   Level: intermediate
5255f2c45f1SShri Abhyankar 
5265f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
5275f2c45f1SShri Abhyankar @*/
5285f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
5295f2c45f1SShri Abhyankar {
5305f2c45f1SShri Abhyankar   PetscErrorCode ierr;
5315f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5325f2c45f1SShri Abhyankar 
5335f2c45f1SShri Abhyankar   PetscFunctionBegin;
5345f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr);
5356fefedf4SHong Zhang   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */
5365f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5375f2c45f1SShri Abhyankar }
5385f2c45f1SShri Abhyankar 
53924121865SAdrian Maldonado /*@
54024121865SAdrian Maldonado   DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector.
54124121865SAdrian Maldonado 
54224121865SAdrian Maldonado   Not Collective
54324121865SAdrian Maldonado 
54424121865SAdrian Maldonado   Input Parameters:
54524121865SAdrian Maldonado + dm     - The DMNetwork object
54624121865SAdrian Maldonado - p      - the edge point
54724121865SAdrian Maldonado 
54824121865SAdrian Maldonado   Output Parameters:
54924121865SAdrian Maldonado . offset - the offset
55024121865SAdrian Maldonado 
55124121865SAdrian Maldonado   Level: intermediate
55224121865SAdrian Maldonado 
55324121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
55424121865SAdrian Maldonado @*/
55524121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
55624121865SAdrian Maldonado {
55724121865SAdrian Maldonado   PetscErrorCode ierr;
55824121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
55924121865SAdrian Maldonado 
56024121865SAdrian Maldonado   PetscFunctionBegin;
56124121865SAdrian Maldonado 
56224121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr);
56324121865SAdrian Maldonado   PetscFunctionReturn(0);
56424121865SAdrian Maldonado }
56524121865SAdrian Maldonado 
56624121865SAdrian Maldonado /*@
56724121865SAdrian Maldonado   DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector.
56824121865SAdrian Maldonado 
56924121865SAdrian Maldonado   Not Collective
57024121865SAdrian Maldonado 
57124121865SAdrian Maldonado   Input Parameters:
57224121865SAdrian Maldonado + dm     - The DMNetwork object
57324121865SAdrian Maldonado - p      - the vertex point
57424121865SAdrian Maldonado 
57524121865SAdrian Maldonado   Output Parameters:
57624121865SAdrian Maldonado . offset - the offset
57724121865SAdrian Maldonado 
57824121865SAdrian Maldonado   Level: intermediate
57924121865SAdrian Maldonado 
58024121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
58124121865SAdrian Maldonado @*/
58224121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
58324121865SAdrian Maldonado {
58424121865SAdrian Maldonado   PetscErrorCode ierr;
58524121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
58624121865SAdrian Maldonado 
58724121865SAdrian Maldonado   PetscFunctionBegin;
58824121865SAdrian Maldonado 
58924121865SAdrian Maldonado   p -= network->vStart;
59024121865SAdrian Maldonado 
59124121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr);
59224121865SAdrian Maldonado   PetscFunctionReturn(0);
59324121865SAdrian Maldonado }
5945f2c45f1SShri Abhyankar /*@
5955f2c45f1SShri Abhyankar   DMNetworkAddNumVariables - Add number of variables associated with a given point.
5965f2c45f1SShri Abhyankar 
5975f2c45f1SShri Abhyankar   Not Collective
5985f2c45f1SShri Abhyankar 
5995f2c45f1SShri Abhyankar   Input Parameters:
6005f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
6015f2c45f1SShri Abhyankar . p    - the vertex/edge point
6025f2c45f1SShri Abhyankar - nvar - number of additional variables
6035f2c45f1SShri Abhyankar 
6045f2c45f1SShri Abhyankar   Level: intermediate
6055f2c45f1SShri Abhyankar 
6065f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables
6075f2c45f1SShri Abhyankar @*/
6085f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
6095f2c45f1SShri Abhyankar {
6105f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6115f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6125f2c45f1SShri Abhyankar 
6135f2c45f1SShri Abhyankar   PetscFunctionBegin;
6145f2c45f1SShri Abhyankar   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
6155f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6165f2c45f1SShri Abhyankar }
6175f2c45f1SShri Abhyankar 
61827f51fceSHong Zhang /*@
61927f51fceSHong Zhang   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
62027f51fceSHong Zhang 
62127f51fceSHong Zhang   Not Collective
62227f51fceSHong Zhang 
62327f51fceSHong Zhang   Input Parameters:
62427f51fceSHong Zhang + dm   - The DMNetworkObject
62527f51fceSHong Zhang - p    - the vertex/edge point
62627f51fceSHong Zhang 
62727f51fceSHong Zhang   Output Parameters:
62827f51fceSHong Zhang . nvar - number of variables
62927f51fceSHong Zhang 
63027f51fceSHong Zhang   Level: intermediate
63127f51fceSHong Zhang 
63227f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
63327f51fceSHong Zhang @*/
63427f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
63527f51fceSHong Zhang {
63627f51fceSHong Zhang   PetscErrorCode ierr;
63727f51fceSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
63827f51fceSHong Zhang 
63927f51fceSHong Zhang   PetscFunctionBegin;
64027f51fceSHong Zhang   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
64127f51fceSHong Zhang   PetscFunctionReturn(0);
64227f51fceSHong Zhang }
64327f51fceSHong Zhang 
6445f2c45f1SShri Abhyankar /*@
6455f2c45f1SShri Abhyankar   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
6465f2c45f1SShri Abhyankar 
6475f2c45f1SShri Abhyankar   Not Collective
6485f2c45f1SShri Abhyankar 
6495f2c45f1SShri Abhyankar   Input Parameters:
6505f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
6515f2c45f1SShri Abhyankar . p    - the vertex/edge point
6525f2c45f1SShri Abhyankar - nvar - number of variables
6535f2c45f1SShri Abhyankar 
6545f2c45f1SShri Abhyankar   Level: intermediate
6555f2c45f1SShri Abhyankar 
6565f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables
6575f2c45f1SShri Abhyankar @*/
6585f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
6595f2c45f1SShri Abhyankar {
6605f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6615f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6625f2c45f1SShri Abhyankar 
6635f2c45f1SShri Abhyankar   PetscFunctionBegin;
6645f2c45f1SShri Abhyankar   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
6655f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6665f2c45f1SShri Abhyankar }
6675f2c45f1SShri Abhyankar 
6685f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This
6695f2c45f1SShri Abhyankar    function is called during DMSetUp() */
6705f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm)
6715f2c45f1SShri Abhyankar {
6725f2c45f1SShri Abhyankar   PetscErrorCode              ierr;
6735f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6745f2c45f1SShri Abhyankar   PetscInt                    arr_size;
6755f2c45f1SShri Abhyankar   PetscInt                    p,offset,offsetp;
6765f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
6775f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue;
6785f2c45f1SShri Abhyankar   DMNetworkComponentGenericDataType      *componentdataarray;
6795f2c45f1SShri Abhyankar   PetscInt ncomp, i;
6805f2c45f1SShri Abhyankar 
6815f2c45f1SShri Abhyankar   PetscFunctionBegin;
6825f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
6835f2c45f1SShri Abhyankar   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
68475b160a0SShri Abhyankar   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
6855f2c45f1SShri Abhyankar   componentdataarray = network->componentdataarray;
6865f2c45f1SShri Abhyankar   for (p = network->pStart; p < network->pEnd; p++) {
6875f2c45f1SShri Abhyankar     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
6885f2c45f1SShri Abhyankar     /* Copy header */
6895f2c45f1SShri Abhyankar     header = &network->header[p];
690302440fdSBarry Smith     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
6915f2c45f1SShri Abhyankar     /* Copy data */
6925f2c45f1SShri Abhyankar     cvalue = &network->cvalue[p];
6935f2c45f1SShri Abhyankar     ncomp = header->ndata;
6945f2c45f1SShri Abhyankar     for (i = 0; i < ncomp; i++) {
6955f2c45f1SShri Abhyankar       offset = offsetp + network->dataheadersize + header->offset[i];
696302440fdSBarry Smith       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
6975f2c45f1SShri Abhyankar     }
6985f2c45f1SShri Abhyankar   }
6995f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7005f2c45f1SShri Abhyankar }
7015f2c45f1SShri Abhyankar 
7025f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */
7035f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm)
7045f2c45f1SShri Abhyankar {
7055f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7065f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7075f2c45f1SShri Abhyankar 
7085f2c45f1SShri Abhyankar   PetscFunctionBegin;
7095f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
7105f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7115f2c45f1SShri Abhyankar }
7125f2c45f1SShri Abhyankar 
7135f2c45f1SShri Abhyankar /*@C
7145f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray - Returns the component data array
7155f2c45f1SShri Abhyankar 
7165f2c45f1SShri Abhyankar   Not Collective
7175f2c45f1SShri Abhyankar 
7185f2c45f1SShri Abhyankar   Input Parameters:
7195f2c45f1SShri Abhyankar . dm - The DMNetwork Object
7205f2c45f1SShri Abhyankar 
7215f2c45f1SShri Abhyankar   Output Parameters:
7225f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components
7235f2c45f1SShri Abhyankar 
7245f2c45f1SShri Abhyankar   Level: intermediate
7255f2c45f1SShri Abhyankar 
726a730d845SHong Zhang .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents
7275f2c45f1SShri Abhyankar @*/
7285f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
7295f2c45f1SShri Abhyankar {
7305f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7315f2c45f1SShri Abhyankar 
7325f2c45f1SShri Abhyankar   PetscFunctionBegin;
7335f2c45f1SShri Abhyankar   *componentdataarray = network->componentdataarray;
7345f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7355f2c45f1SShri Abhyankar }
7365f2c45f1SShri Abhyankar 
73724121865SAdrian Maldonado /* Get a subsection from a range of points */
73824121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection)
73924121865SAdrian Maldonado {
74024121865SAdrian Maldonado   PetscErrorCode ierr;
74124121865SAdrian Maldonado   PetscInt       i, nvar;
74224121865SAdrian Maldonado 
74324121865SAdrian Maldonado   PetscFunctionBegin;
74424121865SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr);
74524121865SAdrian Maldonado   ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr);
74624121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
74724121865SAdrian Maldonado     ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr);
74824121865SAdrian Maldonado     ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr);
74924121865SAdrian Maldonado   }
75024121865SAdrian Maldonado 
75124121865SAdrian Maldonado   ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr);
75224121865SAdrian Maldonado   PetscFunctionReturn(0);
75324121865SAdrian Maldonado }
75424121865SAdrian Maldonado 
75524121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */
75624121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
75724121865SAdrian Maldonado {
75824121865SAdrian Maldonado   PetscErrorCode ierr;
75924121865SAdrian Maldonado   PetscInt       i, *subpoints;
76024121865SAdrian Maldonado 
76124121865SAdrian Maldonado   PetscFunctionBegin;
76224121865SAdrian Maldonado   /* Create index sets to map from "points" to "subpoints" */
76324121865SAdrian Maldonado   ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr);
76424121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
76524121865SAdrian Maldonado     subpoints[i - pstart] = i;
76624121865SAdrian Maldonado   }
767459726d8SSatish Balay   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr);
76824121865SAdrian Maldonado   ierr = PetscFree(subpoints);CHKERRQ(ierr);
76924121865SAdrian Maldonado   PetscFunctionReturn(0);
77024121865SAdrian Maldonado }
77124121865SAdrian Maldonado 
77224121865SAdrian Maldonado /*@
77324121865SAdrian Maldonado   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute.
77424121865SAdrian Maldonado 
77524121865SAdrian Maldonado   Collective
77624121865SAdrian Maldonado 
77724121865SAdrian Maldonado   Input Parameters:
77824121865SAdrian Maldonado . dm   - The DMNetworkObject
77924121865SAdrian Maldonado 
78024121865SAdrian Maldonado   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
78124121865SAdrian Maldonado 
78224121865SAdrian Maldonado   points = [0 1 2 3 4 5 6]
78324121865SAdrian Maldonado 
78424121865SAdrian Maldonado   where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]).
78524121865SAdrian Maldonado 
78624121865SAdrian Maldonado   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
78724121865SAdrian Maldonado 
78824121865SAdrian Maldonado   Level: intermediate
78924121865SAdrian Maldonado 
79024121865SAdrian Maldonado @*/
79124121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
79224121865SAdrian Maldonado {
79324121865SAdrian Maldonado   PetscErrorCode ierr;
79424121865SAdrian Maldonado   MPI_Comm       comm;
7959852e123SBarry Smith   PetscMPIInt    rank, size;
79624121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
79724121865SAdrian Maldonado 
798eab1376dSHong Zhang   PetscFunctionBegin;
79924121865SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
80024121865SAdrian Maldonado   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
8019852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
80224121865SAdrian Maldonado 
80324121865SAdrian Maldonado   /* Create maps for vertices and edges */
80424121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr);
80524121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr);
80624121865SAdrian Maldonado 
80724121865SAdrian Maldonado   /* Create local sub-sections */
80824121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr);
80924121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr);
81024121865SAdrian Maldonado 
8119852e123SBarry Smith   if (size > 1) {
81224121865SAdrian Maldonado     ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr);
81324121865SAdrian Maldonado     ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr);
81424121865SAdrian Maldonado   ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr);
81524121865SAdrian Maldonado   ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr);
81624121865SAdrian Maldonado   } else {
81724121865SAdrian Maldonado   /* create structures for vertex */
81824121865SAdrian Maldonado   ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr);
81924121865SAdrian Maldonado   /* create structures for edge */
82024121865SAdrian Maldonado   ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr);
82124121865SAdrian Maldonado   }
82224121865SAdrian Maldonado 
82324121865SAdrian Maldonado 
82424121865SAdrian Maldonado   /* Add viewers */
82524121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr);
82624121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr);
82724121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr);
82824121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr);
82924121865SAdrian Maldonado 
83024121865SAdrian Maldonado   PetscFunctionReturn(0);
83124121865SAdrian Maldonado }
8327b6afd5bSHong Zhang 
8335f2c45f1SShri Abhyankar /*@
8345f2c45f1SShri Abhyankar   DMNetworkDistribute - Distributes the network and moves associated component data.
8355f2c45f1SShri Abhyankar 
8365f2c45f1SShri Abhyankar   Collective
8375f2c45f1SShri Abhyankar 
8385f2c45f1SShri Abhyankar   Input Parameter:
839d3464fd4SAdrian Maldonado + DM - the DMNetwork object
8405f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default
8415f2c45f1SShri Abhyankar 
8425f2c45f1SShri Abhyankar   Notes:
8438b171c8eSHong Zhang   Distributes the network with <overlap>-overlapping partitioning of the edges.
8445f2c45f1SShri Abhyankar 
8455f2c45f1SShri Abhyankar   Level: intermediate
8465f2c45f1SShri Abhyankar 
8475f2c45f1SShri Abhyankar .seealso: DMNetworkCreate
8485f2c45f1SShri Abhyankar @*/
849d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
8505f2c45f1SShri Abhyankar {
851d3464fd4SAdrian Maldonado   MPI_Comm       comm;
8525f2c45f1SShri Abhyankar   PetscErrorCode ierr;
853d3464fd4SAdrian Maldonado   PetscMPIInt    size;
854d3464fd4SAdrian Maldonado   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
855d3464fd4SAdrian Maldonado   DM_Network     *newDMnetwork;
8565f2c45f1SShri Abhyankar   PetscSF        pointsf;
8575f2c45f1SShri Abhyankar   DM             newDM;
85851ac5effSHong Zhang   PetscPartitioner part;
8595f2c45f1SShri Abhyankar 
8605f2c45f1SShri Abhyankar   PetscFunctionBegin;
861d3464fd4SAdrian Maldonado 
862d3464fd4SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr);
863d3464fd4SAdrian Maldonado   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
864d3464fd4SAdrian Maldonado   if (size == 1) PetscFunctionReturn(0);
865d3464fd4SAdrian Maldonado 
866d3464fd4SAdrian Maldonado   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr);
8675f2c45f1SShri Abhyankar   newDMnetwork = (DM_Network*)newDM->data;
8685f2c45f1SShri Abhyankar   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
86951ac5effSHong Zhang 
87051ac5effSHong Zhang   /* Enable runtime options for petscpartitioner */
87151ac5effSHong Zhang   ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr);
87251ac5effSHong Zhang   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
87351ac5effSHong Zhang 
8745f2c45f1SShri Abhyankar   /* Distribute plex dm and dof section */
87580cf41d5SMatthew G. Knepley   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
87651ac5effSHong Zhang 
8775f2c45f1SShri Abhyankar   /* Distribute dof section */
878d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr);
8795f2c45f1SShri Abhyankar   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
880d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr);
88151ac5effSHong Zhang 
8825f2c45f1SShri Abhyankar   /* Distribute data and associated section */
88331da1fc8SHong Zhang   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
88424121865SAdrian Maldonado 
8855f2c45f1SShri Abhyankar   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
8865f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
8875f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
8885f2c45f1SShri Abhyankar   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
8896fefedf4SHong Zhang   newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
8906fefedf4SHong Zhang   newDMnetwork->NVertices = oldDMnetwork->NVertices;
8915f2c45f1SShri Abhyankar   newDMnetwork->NEdges = oldDMnetwork->NEdges;
89224121865SAdrian Maldonado 
8935f2c45f1SShri Abhyankar   /* Set Dof section as the default section for dm */
8945f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
8955f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
8965f2c45f1SShri Abhyankar 
89724121865SAdrian Maldonado   /* Destroy point SF */
89824121865SAdrian Maldonado   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
89924121865SAdrian Maldonado 
900d3464fd4SAdrian Maldonado   ierr = DMDestroy(dm);CHKERRQ(ierr);
901d3464fd4SAdrian Maldonado   *dm  = newDM;
9025f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9035f2c45f1SShri Abhyankar }
9045f2c45f1SShri Abhyankar 
90524121865SAdrian Maldonado /*@C
90624121865SAdrian Maldonado   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering.
90724121865SAdrian Maldonado 
90824121865SAdrian Maldonado   Input Parameters:
90924121865SAdrian Maldonado + masterSF - the original SF structure
91024121865SAdrian Maldonado - map      - a ISLocalToGlobal mapping that contains the subset of points
91124121865SAdrian Maldonado 
91224121865SAdrian Maldonado   Output Parameters:
91324121865SAdrian Maldonado . subSF    - a subset of the masterSF for the desired subset.
91424121865SAdrian Maldonado */
91524121865SAdrian Maldonado 
91624121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) {
91724121865SAdrian Maldonado 
91824121865SAdrian Maldonado   PetscErrorCode        ierr;
91924121865SAdrian Maldonado   PetscInt              nroots, nleaves, *ilocal_sub;
92024121865SAdrian Maldonado   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
92124121865SAdrian Maldonado   PetscInt              *local_points, *remote_points;
92224121865SAdrian Maldonado   PetscSFNode           *iremote_sub;
92324121865SAdrian Maldonado   const PetscInt        *ilocal;
92424121865SAdrian Maldonado   const PetscSFNode     *iremote;
92524121865SAdrian Maldonado 
92624121865SAdrian Maldonado   PetscFunctionBegin;
92724121865SAdrian Maldonado   ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
92824121865SAdrian Maldonado 
92924121865SAdrian Maldonado   /* Look for leaves that pertain to the subset of points. Get the local ordering */
93024121865SAdrian Maldonado   ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr);
93124121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr);
93224121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
93324121865SAdrian Maldonado     if (ilocal_map[i] != -1) nleaves_sub += 1;
93424121865SAdrian Maldonado   }
93524121865SAdrian Maldonado   /* Re-number ilocal with subset numbering. Need information from roots */
93624121865SAdrian Maldonado   ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr);
93724121865SAdrian Maldonado   for (i = 0; i < nroots; i++) local_points[i] = i;
93824121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr);
93924121865SAdrian Maldonado   ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
94024121865SAdrian Maldonado   ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
94124121865SAdrian Maldonado   /* Fill up graph using local (that is, local to the subset) numbering. */
9424b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr);
9434b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr);
94424121865SAdrian Maldonado   nleaves_sub = 0;
94524121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
94624121865SAdrian Maldonado     if (ilocal_map[i] != -1) {
94724121865SAdrian Maldonado       ilocal_sub[nleaves_sub] = ilocal_map[i];
9484b70a8deSAdrian Maldonado       iremote_sub[nleaves_sub].rank = iremote[i].rank;
94924121865SAdrian Maldonado       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
95024121865SAdrian Maldonado       nleaves_sub += 1;
95124121865SAdrian Maldonado     }
95224121865SAdrian Maldonado   }
95324121865SAdrian Maldonado   ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr);
95424121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr);
95524121865SAdrian Maldonado 
95624121865SAdrian Maldonado   /* Create new subSF */
95724121865SAdrian Maldonado   ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr);
95824121865SAdrian Maldonado   ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr);
9594b70a8deSAdrian Maldonado   ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr);
96024121865SAdrian Maldonado   ierr = PetscFree(ilocal_map);CHKERRQ(ierr);
9614b70a8deSAdrian Maldonado   ierr = PetscFree(iremote_sub);CHKERRQ(ierr);
96224121865SAdrian Maldonado   PetscFunctionReturn(0);
96324121865SAdrian Maldonado }
96424121865SAdrian Maldonado 
9655f2c45f1SShri Abhyankar /*@C
9665f2c45f1SShri Abhyankar   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
9675f2c45f1SShri Abhyankar 
9685f2c45f1SShri Abhyankar   Not Collective
9695f2c45f1SShri Abhyankar 
9705f2c45f1SShri Abhyankar   Input Parameters:
9715f2c45f1SShri Abhyankar + dm - The DMNetwork object
9725f2c45f1SShri Abhyankar - p  - the vertex point
9735f2c45f1SShri Abhyankar 
9745f2c45f1SShri Abhyankar   Output Paramters:
9755f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point
9765f2c45f1SShri Abhyankar - edges  - List of edge points
9775f2c45f1SShri Abhyankar 
9785f2c45f1SShri Abhyankar   Level: intermediate
9795f2c45f1SShri Abhyankar 
9805f2c45f1SShri Abhyankar   Fortran Notes:
9815f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
9825f2c45f1SShri Abhyankar   include petsc.h90 in your code.
9835f2c45f1SShri Abhyankar 
984d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices
9855f2c45f1SShri Abhyankar @*/
9865f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
9875f2c45f1SShri Abhyankar {
9885f2c45f1SShri Abhyankar   PetscErrorCode ierr;
9895f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
9905f2c45f1SShri Abhyankar 
9915f2c45f1SShri Abhyankar   PetscFunctionBegin;
9925f2c45f1SShri Abhyankar   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
9935f2c45f1SShri Abhyankar   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
9945f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9955f2c45f1SShri Abhyankar }
9965f2c45f1SShri Abhyankar 
9975f2c45f1SShri Abhyankar /*@C
998d842c372SHong Zhang   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
9995f2c45f1SShri Abhyankar 
10005f2c45f1SShri Abhyankar   Not Collective
10015f2c45f1SShri Abhyankar 
10025f2c45f1SShri Abhyankar   Input Parameters:
10035f2c45f1SShri Abhyankar + dm - The DMNetwork object
10045f2c45f1SShri Abhyankar - p  - the edge point
10055f2c45f1SShri Abhyankar 
10065f2c45f1SShri Abhyankar   Output Paramters:
10075f2c45f1SShri Abhyankar . vertices  - vertices connected to this edge
10085f2c45f1SShri Abhyankar 
10095f2c45f1SShri Abhyankar   Level: intermediate
10105f2c45f1SShri Abhyankar 
10115f2c45f1SShri Abhyankar   Fortran Notes:
10125f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
10135f2c45f1SShri Abhyankar   include petsc.h90 in your code.
10145f2c45f1SShri Abhyankar 
10155f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
10165f2c45f1SShri Abhyankar @*/
1017d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
10185f2c45f1SShri Abhyankar {
10195f2c45f1SShri Abhyankar   PetscErrorCode ierr;
10205f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
10215f2c45f1SShri Abhyankar 
10225f2c45f1SShri Abhyankar   PetscFunctionBegin;
10235f2c45f1SShri Abhyankar   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
10245f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10255f2c45f1SShri Abhyankar }
10265f2c45f1SShri Abhyankar 
10275f2c45f1SShri Abhyankar /*@
10285f2c45f1SShri Abhyankar   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
10295f2c45f1SShri Abhyankar 
10305f2c45f1SShri Abhyankar   Not Collective
10315f2c45f1SShri Abhyankar 
10325f2c45f1SShri Abhyankar   Input Parameters:
10335f2c45f1SShri Abhyankar + dm - The DMNetwork object
10345f2c45f1SShri Abhyankar . p  - the vertex point
10355f2c45f1SShri Abhyankar 
10365f2c45f1SShri Abhyankar   Output Parameter:
10375f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point
10385f2c45f1SShri Abhyankar 
10395f2c45f1SShri Abhyankar   Level: intermediate
10405f2c45f1SShri Abhyankar 
1041d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange
10425f2c45f1SShri Abhyankar @*/
10435f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
10445f2c45f1SShri Abhyankar {
10455f2c45f1SShri Abhyankar   PetscErrorCode ierr;
10465f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
10475f2c45f1SShri Abhyankar   PetscInt       offsetg;
10485f2c45f1SShri Abhyankar   PetscSection   sectiong;
10495f2c45f1SShri Abhyankar 
10505f2c45f1SShri Abhyankar   PetscFunctionBegin;
10515f2c45f1SShri Abhyankar   *isghost = PETSC_FALSE;
10525f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
10535f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
10545f2c45f1SShri Abhyankar   if (offsetg < 0) *isghost = PETSC_TRUE;
10555f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10565f2c45f1SShri Abhyankar }
10575f2c45f1SShri Abhyankar 
10585f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm)
10595f2c45f1SShri Abhyankar {
10605f2c45f1SShri Abhyankar   PetscErrorCode ierr;
10615f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
10625f2c45f1SShri Abhyankar 
10635f2c45f1SShri Abhyankar   PetscFunctionBegin;
10645f2c45f1SShri Abhyankar   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
10655f2c45f1SShri Abhyankar   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
10665f2c45f1SShri Abhyankar 
10675f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
10685f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
10695f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10705f2c45f1SShri Abhyankar }
10715f2c45f1SShri Abhyankar 
10721ad426b7SHong Zhang /*@
107317df6e9eSHong Zhang     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
10741ad426b7SHong Zhang                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
10751ad426b7SHong Zhang 
10761ad426b7SHong Zhang     Collective
10771ad426b7SHong Zhang 
10781ad426b7SHong Zhang     Input Parameters:
107983b2e829SHong Zhang +   dm - The DMNetwork object
108083b2e829SHong Zhang .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
108183b2e829SHong Zhang -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
10821ad426b7SHong Zhang 
10831ad426b7SHong Zhang     Level: intermediate
10841ad426b7SHong Zhang 
10851ad426b7SHong Zhang @*/
108683b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
10871ad426b7SHong Zhang {
10881ad426b7SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
10891ad426b7SHong Zhang 
10901ad426b7SHong Zhang   PetscFunctionBegin;
109183b2e829SHong Zhang   network->userEdgeJacobian   = eflg;
109283b2e829SHong Zhang   network->userVertexJacobian = vflg;
10931ad426b7SHong Zhang   PetscFunctionReturn(0);
10941ad426b7SHong Zhang }
10951ad426b7SHong Zhang 
10961ad426b7SHong Zhang /*@
109783b2e829SHong Zhang     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
109883b2e829SHong Zhang 
109983b2e829SHong Zhang     Not Collective
110083b2e829SHong Zhang 
110183b2e829SHong Zhang     Input Parameters:
110283b2e829SHong Zhang +   dm - The DMNetwork object
110383b2e829SHong Zhang .   p  - the edge point
11043e97b6e8SHong Zhang -   J - array (size = 3) of Jacobian submatrices for this edge point:
11053e97b6e8SHong Zhang         J[0]: this edge
1106d842c372SHong Zhang         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()
110783b2e829SHong Zhang 
110883b2e829SHong Zhang     Level: intermediate
110983b2e829SHong Zhang 
111083b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix
111183b2e829SHong Zhang @*/
111283b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
111383b2e829SHong Zhang {
111483b2e829SHong Zhang   PetscErrorCode ierr;
111583b2e829SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
111683b2e829SHong Zhang 
111783b2e829SHong Zhang   PetscFunctionBegin;
1118883e35e8SHong Zhang   if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
111983b2e829SHong Zhang   if (!network->Je) {
1120883e35e8SHong Zhang     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
112183b2e829SHong Zhang   }
1122883e35e8SHong Zhang   network->Je[3*p]   = J[0];
1123883e35e8SHong Zhang   network->Je[3*p+1] = J[1];
1124883e35e8SHong Zhang   network->Je[3*p+2] = J[2];
112583b2e829SHong Zhang   PetscFunctionReturn(0);
112683b2e829SHong Zhang }
112783b2e829SHong Zhang 
112883b2e829SHong Zhang /*@
112976ddfea5SHong Zhang     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
11301ad426b7SHong Zhang 
11311ad426b7SHong Zhang     Not Collective
11321ad426b7SHong Zhang 
11331ad426b7SHong Zhang     Input Parameters:
11341ad426b7SHong Zhang +   dm - The DMNetwork object
11351ad426b7SHong Zhang .   p  - the vertex point
11363e97b6e8SHong Zhang -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
11373e97b6e8SHong Zhang         J[0]:       this vertex
11383e97b6e8SHong Zhang         J[1+2*i]:   i-th supporting edge
11393e97b6e8SHong Zhang         J[1+2*i+1]: i-th connected vertex
11401ad426b7SHong Zhang 
11411ad426b7SHong Zhang     Level: intermediate
11421ad426b7SHong Zhang 
114383b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix
11441ad426b7SHong Zhang @*/
1145883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
11465f2c45f1SShri Abhyankar {
11475f2c45f1SShri Abhyankar   PetscErrorCode ierr;
11485f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
1149f4431b8cSHong Zhang   PetscInt       i,*vptr,nedges,vStart,vEnd;
1150883e35e8SHong Zhang   const PetscInt *edges;
11515f2c45f1SShri Abhyankar 
11525f2c45f1SShri Abhyankar   PetscFunctionBegin;
1153883e35e8SHong Zhang   if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
1154883e35e8SHong Zhang 
1155883e35e8SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1156f4431b8cSHong Zhang 
115783b2e829SHong Zhang   if (!network->Jv) {
11586fefedf4SHong Zhang     PetscInt nVertices = network->nVertices,nedges_total;
1159883e35e8SHong Zhang     /* count nvertex_total */
11603e97b6e8SHong Zhang     nedges_total = 0;
1161883e35e8SHong Zhang     ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
11626fefedf4SHong Zhang     ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr);
1163f4431b8cSHong Zhang 
1164883e35e8SHong Zhang     vptr[0] = 0;
11656fefedf4SHong Zhang     for (i=0; i<nVertices; i++) {
1166f4431b8cSHong Zhang       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
1167883e35e8SHong Zhang       nedges_total += nedges;
1168f4431b8cSHong Zhang       vptr[i+1] = vptr[i] + 2*nedges + 1;
11691ad426b7SHong Zhang     }
11703e97b6e8SHong Zhang 
11716fefedf4SHong Zhang     ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr);
1172883e35e8SHong Zhang     network->Jvptr = vptr;
1173883e35e8SHong Zhang   }
1174883e35e8SHong Zhang 
1175883e35e8SHong Zhang   vptr = network->Jvptr;
11763e97b6e8SHong Zhang   network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
11773e97b6e8SHong Zhang 
11783e97b6e8SHong Zhang   /* Set Jacobian for each supporting edge and connected vertex */
1179883e35e8SHong Zhang   ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
1180883e35e8SHong Zhang   for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
1181883e35e8SHong Zhang   PetscFunctionReturn(0);
1182883e35e8SHong Zhang }
1183883e35e8SHong Zhang 
1184e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
11855cf7da58SHong Zhang {
11865cf7da58SHong Zhang   PetscErrorCode ierr;
11875cf7da58SHong Zhang   PetscInt       j;
11885cf7da58SHong Zhang   PetscScalar    val=(PetscScalar)ncols;
11895cf7da58SHong Zhang 
11905cf7da58SHong Zhang   PetscFunctionBegin;
11915cf7da58SHong Zhang   if (!ghost) {
11925cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
11935cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
11945cf7da58SHong Zhang     }
11955cf7da58SHong Zhang   } else {
11965cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
11975cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
11985cf7da58SHong Zhang     }
11995cf7da58SHong Zhang   }
12005cf7da58SHong Zhang   PetscFunctionReturn(0);
12015cf7da58SHong Zhang }
12025cf7da58SHong Zhang 
1203e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
12045cf7da58SHong Zhang {
12055cf7da58SHong Zhang   PetscErrorCode ierr;
12065cf7da58SHong Zhang   PetscInt       j,ncols_u;
12075cf7da58SHong Zhang   PetscScalar    val;
12085cf7da58SHong Zhang 
12095cf7da58SHong Zhang   PetscFunctionBegin;
12105cf7da58SHong Zhang   if (!ghost) {
12115cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
12125cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
12135cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
12145cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
12155cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
12165cf7da58SHong Zhang     }
12175cf7da58SHong Zhang   } else {
12185cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
12195cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
12205cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
12215cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
12225cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
12235cf7da58SHong Zhang     }
12245cf7da58SHong Zhang   }
12255cf7da58SHong Zhang   PetscFunctionReturn(0);
12265cf7da58SHong Zhang }
12275cf7da58SHong Zhang 
1228e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
12295cf7da58SHong Zhang {
12305cf7da58SHong Zhang   PetscErrorCode ierr;
12315cf7da58SHong Zhang 
12325cf7da58SHong Zhang   PetscFunctionBegin;
12335cf7da58SHong Zhang   if (Ju) {
12345cf7da58SHong Zhang     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
12355cf7da58SHong Zhang   } else {
12365cf7da58SHong Zhang     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
12375cf7da58SHong Zhang   }
12385cf7da58SHong Zhang   PetscFunctionReturn(0);
12395cf7da58SHong Zhang }
12405cf7da58SHong Zhang 
1241e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1242883e35e8SHong Zhang {
1243883e35e8SHong Zhang   PetscErrorCode ierr;
1244883e35e8SHong Zhang   PetscInt       j,*cols;
1245883e35e8SHong Zhang   PetscScalar    *zeros;
1246883e35e8SHong Zhang 
1247883e35e8SHong Zhang   PetscFunctionBegin;
1248883e35e8SHong Zhang   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
1249883e35e8SHong Zhang   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1250883e35e8SHong Zhang   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
1251883e35e8SHong Zhang   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
12521ad426b7SHong Zhang   PetscFunctionReturn(0);
12531ad426b7SHong Zhang }
1254a4e85ca8SHong Zhang 
1255e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
12563e97b6e8SHong Zhang {
12573e97b6e8SHong Zhang   PetscErrorCode ierr;
12583e97b6e8SHong Zhang   PetscInt       j,M,N,row,col,ncols_u;
12593e97b6e8SHong Zhang   const PetscInt *cols;
12603e97b6e8SHong Zhang   PetscScalar    zero=0.0;
12613e97b6e8SHong Zhang 
12623e97b6e8SHong Zhang   PetscFunctionBegin;
12633e97b6e8SHong Zhang   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
12643e97b6e8SHong 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);
12653e97b6e8SHong Zhang 
12663e97b6e8SHong Zhang   for (row=0; row<nrows; row++) {
12673e97b6e8SHong Zhang     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
12683e97b6e8SHong Zhang     for (j=0; j<ncols_u; j++) {
12693e97b6e8SHong Zhang       col = cols[j] + cstart;
12703e97b6e8SHong Zhang       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
12713e97b6e8SHong Zhang     }
12723e97b6e8SHong Zhang     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
12733e97b6e8SHong Zhang   }
12743e97b6e8SHong Zhang   PetscFunctionReturn(0);
12753e97b6e8SHong Zhang }
12761ad426b7SHong Zhang 
1277e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1278a4e85ca8SHong Zhang {
1279a4e85ca8SHong Zhang   PetscErrorCode ierr;
1280f4431b8cSHong Zhang 
1281a4e85ca8SHong Zhang   PetscFunctionBegin;
1282a4e85ca8SHong Zhang   if (Ju) {
1283a4e85ca8SHong Zhang     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1284a4e85ca8SHong Zhang   } else {
1285a4e85ca8SHong Zhang     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1286a4e85ca8SHong Zhang   }
1287a4e85ca8SHong Zhang   PetscFunctionReturn(0);
1288a4e85ca8SHong Zhang }
1289a4e85ca8SHong Zhang 
129024121865SAdrian 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.
129124121865SAdrian Maldonado */
129224121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
129324121865SAdrian Maldonado {
129424121865SAdrian Maldonado   PetscErrorCode ierr;
129524121865SAdrian Maldonado   PetscInt       i, size, dof;
129624121865SAdrian Maldonado   PetscInt       *glob2loc;
129724121865SAdrian Maldonado 
129824121865SAdrian Maldonado   PetscFunctionBegin;
129924121865SAdrian Maldonado   ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr);
130024121865SAdrian Maldonado   ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr);
130124121865SAdrian Maldonado 
130224121865SAdrian Maldonado   for (i = 0; i < size; i++) {
130324121865SAdrian Maldonado     ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr);
130424121865SAdrian Maldonado     dof = (dof >= 0) ? dof : -(dof + 1);
130524121865SAdrian Maldonado     glob2loc[i] = dof;
130624121865SAdrian Maldonado   }
130724121865SAdrian Maldonado 
130824121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
130924121865SAdrian Maldonado #if 0
131024121865SAdrian Maldonado   ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
131124121865SAdrian Maldonado #endif
131224121865SAdrian Maldonado   PetscFunctionReturn(0);
131324121865SAdrian Maldonado }
131424121865SAdrian Maldonado 
131501ad2aeeSHong Zhang #include <petsc/private/matimpl.h>
13161ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
13171ad426b7SHong Zhang {
13181ad426b7SHong Zhang   PetscErrorCode ierr;
131924121865SAdrian Maldonado   PetscMPIInt    rank, size;
13201ad426b7SHong Zhang   DM_Network     *network = (DM_Network*) dm->data;
1321a4e85ca8SHong Zhang   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1322840c2264SHong Zhang   PetscInt       cstart,ncols,j,e,v;
132324121865SAdrian Maldonado   PetscBool      ghost,ghost_vc,ghost2,isNest;
1324a4e85ca8SHong Zhang   Mat            Juser;
1325bfbc38dcSHong Zhang   PetscSection   sectionGlobal;
1326447d78afSSatish Balay   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
1327a4e85ca8SHong Zhang   const PetscInt *edges,*cone;
13285cf7da58SHong Zhang   MPI_Comm       comm;
132924121865SAdrian Maldonado   MatType        mtype;
13305cf7da58SHong Zhang   Vec            vd_nz,vo_nz;
13315cf7da58SHong Zhang   PetscInt       *dnnz,*onnz;
13325cf7da58SHong Zhang   PetscScalar    *vdnz,*vonz;
13331ad426b7SHong Zhang 
13341ad426b7SHong Zhang   PetscFunctionBegin;
133524121865SAdrian Maldonado   mtype = dm->mattype;
133624121865SAdrian Maldonado   ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr);
133724121865SAdrian Maldonado 
133824121865SAdrian Maldonado   if (isNest) {
13390731d606SHong Zhang     /* ierr = DMCreateMatrix_Network_Nest(); */
134024121865SAdrian Maldonado     PetscInt   eDof, vDof;
134124121865SAdrian Maldonado     Mat        j11, j12, j21, j22, bA[2][2];
134224121865SAdrian Maldonado     ISLocalToGlobalMapping eISMap, vISMap;
134324121865SAdrian Maldonado 
134424121865SAdrian Maldonado     ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
134524121865SAdrian Maldonado     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
134624121865SAdrian Maldonado     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
134724121865SAdrian Maldonado 
134824121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr);
134924121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr);
135024121865SAdrian Maldonado 
135101ad2aeeSHong Zhang     ierr = MatCreate(comm, &j11);CHKERRQ(ierr);
135224121865SAdrian Maldonado     ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
135324121865SAdrian Maldonado     ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr);
135424121865SAdrian Maldonado 
135501ad2aeeSHong Zhang     ierr = MatCreate(comm, &j12);CHKERRQ(ierr);
135624121865SAdrian Maldonado     ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr);
135724121865SAdrian Maldonado     ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr);
135824121865SAdrian Maldonado 
135901ad2aeeSHong Zhang     ierr = MatCreate(comm, &j21);CHKERRQ(ierr);
136024121865SAdrian Maldonado     ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
136124121865SAdrian Maldonado     ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr);
136224121865SAdrian Maldonado 
136301ad2aeeSHong Zhang     ierr = MatCreate(comm, &j22);CHKERRQ(ierr);
136424121865SAdrian Maldonado     ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
136524121865SAdrian Maldonado     ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr);
136624121865SAdrian Maldonado 
13673f6a6bdaSHong Zhang     bA[0][0] = j11;
13683f6a6bdaSHong Zhang     bA[0][1] = j12;
13693f6a6bdaSHong Zhang     bA[1][0] = j21;
13703f6a6bdaSHong Zhang     bA[1][1] = j22;
137124121865SAdrian Maldonado 
137224121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr);
137324121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr);
137424121865SAdrian Maldonado 
137524121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr);
137624121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr);
137724121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr);
137824121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr);
137924121865SAdrian Maldonado 
138024121865SAdrian Maldonado     ierr = MatSetUp(j11);CHKERRQ(ierr);
138124121865SAdrian Maldonado     ierr = MatSetUp(j12);CHKERRQ(ierr);
138224121865SAdrian Maldonado     ierr = MatSetUp(j21);CHKERRQ(ierr);
138324121865SAdrian Maldonado     ierr = MatSetUp(j22);CHKERRQ(ierr);
138424121865SAdrian Maldonado 
138501ad2aeeSHong Zhang     ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr);
138624121865SAdrian Maldonado     ierr = MatSetUp(*J);CHKERRQ(ierr);
138724121865SAdrian Maldonado     ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr);
138824121865SAdrian Maldonado     ierr = MatDestroy(&j11);CHKERRQ(ierr);
138924121865SAdrian Maldonado     ierr = MatDestroy(&j12);CHKERRQ(ierr);
139024121865SAdrian Maldonado     ierr = MatDestroy(&j21);CHKERRQ(ierr);
139124121865SAdrian Maldonado     ierr = MatDestroy(&j22);CHKERRQ(ierr);
139224121865SAdrian Maldonado 
139324121865SAdrian Maldonado     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
139424121865SAdrian Maldonado     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
139524121865SAdrian Maldonado 
139624121865SAdrian Maldonado     /* Free structures */
139724121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr);
139824121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr);
139924121865SAdrian Maldonado 
140024121865SAdrian Maldonado     PetscFunctionReturn(0);
140124121865SAdrian Maldonado   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1402a4e85ca8SHong Zhang     /* user does not provide Jacobian blocks */
1403bfbc38dcSHong Zhang     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1404bfbc38dcSHong Zhang     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
14051ad426b7SHong Zhang     PetscFunctionReturn(0);
14061ad426b7SHong Zhang   }
14071ad426b7SHong Zhang 
1408bfbc38dcSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
14092a945128SHong Zhang   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1410bfbc38dcSHong Zhang   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1411bfbc38dcSHong Zhang   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
14122a945128SHong Zhang 
14132a945128SHong Zhang   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
14142a945128SHong Zhang   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
141589898e50SHong Zhang 
141689898e50SHong Zhang   /* (1) Set matrix preallocation */
141789898e50SHong Zhang   /*------------------------------*/
1418840c2264SHong Zhang   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1419840c2264SHong Zhang   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1420840c2264SHong Zhang   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1421840c2264SHong Zhang   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1422840c2264SHong Zhang   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1423840c2264SHong Zhang   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1424840c2264SHong Zhang 
142589898e50SHong Zhang   /* Set preallocation for edges */
142689898e50SHong Zhang   /*-----------------------------*/
1427840c2264SHong Zhang   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1428840c2264SHong Zhang 
1429bdcb62a2SHong Zhang   ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr);
1430840c2264SHong Zhang   for (e=eStart; e<eEnd; e++) {
1431840c2264SHong Zhang     /* Get row indices */
1432840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1433840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1434840c2264SHong Zhang     if (nrows) {
1435840c2264SHong Zhang       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1436840c2264SHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1437840c2264SHong Zhang 
14385cf7da58SHong Zhang       /* Set preallocation for conntected vertices */
1439d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1440840c2264SHong Zhang       for (v=0; v<2; v++) {
1441840c2264SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1442840c2264SHong Zhang 
1443840c2264SHong Zhang         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1444840c2264SHong Zhang         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
14455cf7da58SHong Zhang         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1446840c2264SHong Zhang       }
1447840c2264SHong Zhang 
144889898e50SHong Zhang       /* Set preallocation for edge self */
1449840c2264SHong Zhang       cstart = rstart;
1450840c2264SHong Zhang       Juser = network->Je[3*e]; /* Jacobian(e,e) */
14515cf7da58SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1452840c2264SHong Zhang     }
1453840c2264SHong Zhang   }
1454840c2264SHong Zhang 
145589898e50SHong Zhang   /* Set preallocation for vertices */
145689898e50SHong Zhang   /*--------------------------------*/
1457840c2264SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1458840c2264SHong Zhang   if (vEnd - vStart) {
1459840c2264SHong Zhang     if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Jv");
1460840c2264SHong Zhang     vptr = network->Jvptr;
1461840c2264SHong Zhang     if (!vptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide vptr");
1462840c2264SHong Zhang   }
1463840c2264SHong Zhang 
1464840c2264SHong Zhang   for (v=vStart; v<vEnd; v++) {
1465840c2264SHong Zhang     /* Get row indices */
1466840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1467840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1468840c2264SHong Zhang     if (!nrows) continue;
1469840c2264SHong Zhang 
1470bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1471bdcb62a2SHong Zhang     if (ghost) {
1472bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1473bdcb62a2SHong Zhang     } else {
1474bdcb62a2SHong Zhang       rows_v = rows;
1475bdcb62a2SHong Zhang     }
1476bdcb62a2SHong Zhang 
1477bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1478840c2264SHong Zhang 
1479840c2264SHong Zhang     /* Get supporting edges and connected vertices */
1480840c2264SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1481840c2264SHong Zhang 
1482840c2264SHong Zhang     for (e=0; e<nedges; e++) {
1483840c2264SHong Zhang       /* Supporting edges */
1484840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1485840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1486840c2264SHong Zhang 
1487840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1488bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1489840c2264SHong Zhang 
1490840c2264SHong Zhang       /* Connected vertices */
1491d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
1492840c2264SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
1493840c2264SHong Zhang       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1494840c2264SHong Zhang 
1495840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1496840c2264SHong Zhang 
1497840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1498e102a522SHong Zhang       if (ghost_vc||ghost) {
1499e102a522SHong Zhang         ghost2 = PETSC_TRUE;
1500e102a522SHong Zhang       } else {
1501e102a522SHong Zhang         ghost2 = PETSC_FALSE;
1502e102a522SHong Zhang       }
1503e102a522SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr);
1504840c2264SHong Zhang     }
1505840c2264SHong Zhang 
150689898e50SHong Zhang     /* Set preallocation for vertex self */
1507840c2264SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1508840c2264SHong Zhang     if (!ghost) {
1509840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1510840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1511bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1512840c2264SHong Zhang     }
1513bdcb62a2SHong Zhang     if (ghost) {
1514bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1515bdcb62a2SHong Zhang     }
1516840c2264SHong Zhang   }
1517840c2264SHong Zhang 
1518840c2264SHong Zhang   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1519840c2264SHong Zhang   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
15205cf7da58SHong Zhang 
15215cf7da58SHong Zhang   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
15225cf7da58SHong Zhang 
15235cf7da58SHong Zhang   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1524840c2264SHong Zhang   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1525840c2264SHong Zhang 
1526840c2264SHong Zhang   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1527840c2264SHong Zhang   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1528840c2264SHong Zhang   for (j=0; j<localSize; j++) {
1529e102a522SHong Zhang     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1530e102a522SHong Zhang     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1531840c2264SHong Zhang   }
1532840c2264SHong Zhang   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1533840c2264SHong Zhang   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1534840c2264SHong Zhang   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1535840c2264SHong Zhang   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1536840c2264SHong Zhang 
15375cf7da58SHong Zhang   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
15385cf7da58SHong Zhang   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
15395cf7da58SHong Zhang   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
15405cf7da58SHong Zhang 
15415cf7da58SHong Zhang   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
15425cf7da58SHong Zhang 
154389898e50SHong Zhang   /* (2) Set matrix entries for edges */
154489898e50SHong Zhang   /*----------------------------------*/
15451ad426b7SHong Zhang   for (e=eStart; e<eEnd; e++) {
1546bfbc38dcSHong Zhang     /* Get row indices */
15471ad426b7SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
154817df6e9eSHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
15494b976069SHong Zhang     if (nrows) {
15504b976069SHong Zhang       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1551bdcb62a2SHong Zhang 
155217df6e9eSHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
15531ad426b7SHong Zhang 
1554bfbc38dcSHong Zhang       /* Set matrix entries for conntected vertices */
1555d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1556bfbc38dcSHong Zhang       for (v=0; v<2; v++) {
1557bfbc38dcSHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1558883e35e8SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
15593e97b6e8SHong Zhang 
1560a4e85ca8SHong Zhang         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1561a4e85ca8SHong Zhang         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1562bfbc38dcSHong Zhang       }
156317df6e9eSHong Zhang 
1564bfbc38dcSHong Zhang       /* Set matrix entries for edge self */
15653e97b6e8SHong Zhang       cstart = rstart;
1566a4e85ca8SHong Zhang       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1567a4e85ca8SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
15681ad426b7SHong Zhang     }
15694b976069SHong Zhang   }
15701ad426b7SHong Zhang 
1571bfbc38dcSHong Zhang   /* Set matrix entries for vertices */
157283b2e829SHong Zhang   /*---------------------------------*/
15731ad426b7SHong Zhang   for (v=vStart; v<vEnd; v++) {
1574bfbc38dcSHong Zhang     /* Get row indices */
1575596e729fSHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1576596e729fSHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
15774b976069SHong Zhang     if (!nrows) continue;
1578596e729fSHong Zhang 
1579bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1580bdcb62a2SHong Zhang     if (ghost) {
1581bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1582bdcb62a2SHong Zhang     } else {
1583bdcb62a2SHong Zhang       rows_v = rows;
1584bdcb62a2SHong Zhang     }
1585bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1586596e729fSHong Zhang 
1587bfbc38dcSHong Zhang     /* Get supporting edges and connected vertices */
1588596e729fSHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1589596e729fSHong Zhang 
1590596e729fSHong Zhang     for (e=0; e<nedges; e++) {
1591bfbc38dcSHong Zhang       /* Supporting edges */
1592596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1593596e729fSHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1594596e729fSHong Zhang 
1595a4e85ca8SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1596bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1597596e729fSHong Zhang 
1598bfbc38dcSHong Zhang       /* Connected vertices */
1599d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
16002a945128SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
16012a945128SHong Zhang 
160244aca652SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
160344aca652SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1604a4e85ca8SHong Zhang 
1605a4e85ca8SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1606bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1607596e729fSHong Zhang     }
1608596e729fSHong Zhang 
1609bfbc38dcSHong Zhang     /* Set matrix entries for vertex self */
16101ad426b7SHong Zhang     if (!ghost) {
1611596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1612a4e85ca8SHong Zhang       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1613bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr);
1614bdcb62a2SHong Zhang     }
1615bdcb62a2SHong Zhang     if (ghost) {
1616bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1617bdcb62a2SHong Zhang     }
16181ad426b7SHong Zhang   }
1619a4e85ca8SHong Zhang   ierr = PetscFree(rows);CHKERRQ(ierr);
1620bdcb62a2SHong Zhang 
16211ad426b7SHong Zhang   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
16221ad426b7SHong Zhang   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1623dd6f46cdSHong Zhang 
16245f2c45f1SShri Abhyankar   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
16255f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16265f2c45f1SShri Abhyankar }
16275f2c45f1SShri Abhyankar 
16285f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm)
16295f2c45f1SShri Abhyankar {
16305f2c45f1SShri Abhyankar   PetscErrorCode ierr;
16315f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
16325f2c45f1SShri Abhyankar 
16335f2c45f1SShri Abhyankar   PetscFunctionBegin;
16348415c774SShri Abhyankar   if (--network->refct > 0) PetscFunctionReturn(0);
163583b2e829SHong Zhang   if (network->Je) {
163683b2e829SHong Zhang     ierr = PetscFree(network->Je);CHKERRQ(ierr);
163783b2e829SHong Zhang   }
163883b2e829SHong Zhang   if (network->Jv) {
1639883e35e8SHong Zhang     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
164083b2e829SHong Zhang     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
16411ad426b7SHong Zhang   }
164213c2a604SAdrian Maldonado 
164313c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr);
164413c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr);
164513c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr);
164613c2a604SAdrian Maldonado   if (network->vertex.sf) {
164713c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr);
164813c2a604SAdrian Maldonado   }
164913c2a604SAdrian Maldonado   /* edge */
165013c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr);
165113c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr);
165213c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr);
165313c2a604SAdrian Maldonado   if (network->edge.sf) {
165413c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr);
165513c2a604SAdrian Maldonado   }
16565f2c45f1SShri Abhyankar   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
16575f2c45f1SShri Abhyankar   network->edges = NULL;
16585f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
16595f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
166083b2e829SHong Zhang 
1661*e2aaf10cSShri Abhyankar   ierr = PetscFree(network->subnet);CHKERRQ(ierr);
16625f2c45f1SShri Abhyankar   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
16635f2c45f1SShri Abhyankar   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
16645f2c45f1SShri Abhyankar   ierr = PetscFree(network->header);CHKERRQ(ierr);
16655f2c45f1SShri Abhyankar   ierr = PetscFree(network);CHKERRQ(ierr);
16665f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16675f2c45f1SShri Abhyankar }
16685f2c45f1SShri Abhyankar 
16695f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
16705f2c45f1SShri Abhyankar {
16715f2c45f1SShri Abhyankar   PetscErrorCode ierr;
16725f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
16735f2c45f1SShri Abhyankar 
16745f2c45f1SShri Abhyankar   PetscFunctionBegin;
16755f2c45f1SShri Abhyankar   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
16765f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16775f2c45f1SShri Abhyankar }
16785f2c45f1SShri Abhyankar 
16795f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
16805f2c45f1SShri Abhyankar {
16815f2c45f1SShri Abhyankar   PetscErrorCode ierr;
16825f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
16835f2c45f1SShri Abhyankar 
16845f2c45f1SShri Abhyankar   PetscFunctionBegin;
16855f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
16865f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16875f2c45f1SShri Abhyankar }
16885f2c45f1SShri Abhyankar 
16895f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
16905f2c45f1SShri Abhyankar {
16915f2c45f1SShri Abhyankar   PetscErrorCode ierr;
16925f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
16935f2c45f1SShri Abhyankar 
16945f2c45f1SShri Abhyankar   PetscFunctionBegin;
16955f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
16965f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16975f2c45f1SShri Abhyankar }
16985f2c45f1SShri Abhyankar 
16995f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
17005f2c45f1SShri Abhyankar {
17015f2c45f1SShri Abhyankar   PetscErrorCode ierr;
17025f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
17035f2c45f1SShri Abhyankar 
17045f2c45f1SShri Abhyankar   PetscFunctionBegin;
17055f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
17065f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
17075f2c45f1SShri Abhyankar }
17085f2c45f1SShri Abhyankar 
17095f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
17105f2c45f1SShri Abhyankar {
17115f2c45f1SShri Abhyankar   PetscErrorCode ierr;
17125f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
17135f2c45f1SShri Abhyankar 
17145f2c45f1SShri Abhyankar   PetscFunctionBegin;
17155f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
17165f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
17175f2c45f1SShri Abhyankar }
1718