xref: /petsc/src/dm/impls/network/network.c (revision daad07d386296cdcbb87925ef5f1432ee4a24ec4)
1af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h>  /*I  "petscdmnetwork.h"  I*/
25f2c45f1SShri Abhyankar 
3e3b1a053SGetnet Betrie PetscLogEvent DMNetwork_LayoutSetUp;
4e3b1a053SGetnet Betrie PetscLogEvent DMNetwork_SetUpNetwork;
5e3b1a053SGetnet Betrie PetscLogEvent DMNetwork_Distribute;
6e3b1a053SGetnet Betrie 
754dfd506SHong Zhang /*
854dfd506SHong Zhang   Creates the component header and value objects for a network point
954dfd506SHong Zhang */
1054dfd506SHong Zhang static PetscErrorCode SetUpNetworkHeaderComponentValue(DM dm,DMNetworkComponentHeader header,DMNetworkComponentValue cvalue)
1154dfd506SHong Zhang {
1254dfd506SHong Zhang   PetscFunctionBegin;
1354dfd506SHong Zhang   /* Allocate arrays for component information */
149566063dSJacob Faibussowitsch   PetscCall(PetscCalloc5(header->maxcomps,&header->size,header->maxcomps,&header->key,header->maxcomps,&header->offset,header->maxcomps,&header->nvar,header->maxcomps,&header->offsetvarrel));
159566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(header->maxcomps,&cvalue->data));
1654dfd506SHong Zhang 
1754dfd506SHong Zhang   /* The size of the header is the size of struct _p_DMNetworkComponentHeader. Since the struct contains PetscInt pointers we cannot use sizeof(struct). So, we need to explicitly calculate the size.
1854dfd506SHong Zhang    If the data header struct changes then this header size calculation needs to be updated. */
1954dfd506SHong Zhang   header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5*header->maxcomps*sizeof(PetscInt);
2054dfd506SHong Zhang   header->hsize /= sizeof(DMNetworkComponentGenericDataType);
2154dfd506SHong Zhang   PetscFunctionReturn(0);
2254dfd506SHong Zhang }
2354dfd506SHong Zhang 
24*daad07d3SAidan Hamilton PetscErrorCode DMNetworkInitializeHeaderComponentData(DM dm)
25*daad07d3SAidan Hamilton {
26*daad07d3SAidan Hamilton   DM_Network *network = (DM_Network*)dm->data;
27*daad07d3SAidan Hamilton   PetscInt np,p,defaultnumcomp = DMNETWORK_MAX_COMP_AT_POINT_DEFAULT;
28*daad07d3SAidan Hamilton 
29*daad07d3SAidan Hamilton   PetscFunctionBegin;
30*daad07d3SAidan Hamilton   np = network->cloneshared->pEnd - network->cloneshared->pStart;
31*daad07d3SAidan Hamilton   if (network->header)
32*daad07d3SAidan Hamilton   for (p=0; p < np; p++) {
33*daad07d3SAidan Hamilton     network->header[p].maxcomps = defaultnumcomp;
34*daad07d3SAidan Hamilton     PetscCall(SetUpNetworkHeaderComponentValue(dm,&network->header[p],&network->cvalue[p]));
35*daad07d3SAidan Hamilton   }
36*daad07d3SAidan Hamilton 
37*daad07d3SAidan Hamilton   PetscFunctionReturn(0);
38*daad07d3SAidan Hamilton }
395f2c45f1SShri Abhyankar /*@
40556ed216SShri Abhyankar   DMNetworkGetPlex - Gets the Plex DM associated with this network DM
41556ed216SShri Abhyankar 
42556ed216SShri Abhyankar   Not collective
43556ed216SShri Abhyankar 
44556ed216SShri Abhyankar   Input Parameters:
452bf73ac6SHong Zhang . dm - the dm object
462bf73ac6SHong Zhang 
472bf73ac6SHong Zhang   Output Parameters:
482bf73ac6SHong Zhang . plexdm - the plex dm object
49556ed216SShri Abhyankar 
50556ed216SShri Abhyankar   Level: Advanced
51556ed216SShri Abhyankar 
52db781477SPatrick Sanan .seealso: `DMNetworkCreate()`
53556ed216SShri Abhyankar @*/
542bf73ac6SHong Zhang PetscErrorCode DMNetworkGetPlex(DM dm,DM *plexdm)
55556ed216SShri Abhyankar {
562bf73ac6SHong Zhang   DM_Network *network = (DM_Network*)dm->data;
57556ed216SShri Abhyankar 
58556ed216SShri Abhyankar   PetscFunctionBegin;
59556ed216SShri Abhyankar   *plexdm = network->plex;
60556ed216SShri Abhyankar   PetscFunctionReturn(0);
61556ed216SShri Abhyankar }
62556ed216SShri Abhyankar 
63556ed216SShri Abhyankar /*@
642bf73ac6SHong Zhang   DMNetworkGetNumSubNetworks - Gets the the number of subnetworks
6572c3e9f7SHong Zhang 
662bf73ac6SHong Zhang   Not collective
6772c3e9f7SHong Zhang 
68f899ff85SJose E. Roman   Input Parameter:
692bf73ac6SHong Zhang . dm - the dm object
702bf73ac6SHong Zhang 
712bf73ac6SHong Zhang   Output Parameters:
722bf73ac6SHong Zhang + nsubnet - local number of subnetworks
732bf73ac6SHong Zhang - Nsubnet - global number of subnetworks
7472c3e9f7SHong Zhang 
7597bb938eSShri Abhyankar   Level: beginner
7672c3e9f7SHong Zhang 
77db781477SPatrick Sanan .seealso: `DMNetworkCreate()`, `DMNetworkSetNumSubNetworks()`
7872c3e9f7SHong Zhang @*/
792bf73ac6SHong Zhang PetscErrorCode DMNetworkGetNumSubNetworks(DM dm,PetscInt *nsubnet,PetscInt *Nsubnet)
8072c3e9f7SHong Zhang {
812bf73ac6SHong Zhang   DM_Network *network = (DM_Network*)dm->data;
8272c3e9f7SHong Zhang 
8372c3e9f7SHong Zhang   PetscFunctionBegin;
84*daad07d3SAidan Hamilton   if (nsubnet) *nsubnet = network->cloneshared->nsubnet;
85*daad07d3SAidan Hamilton   if (Nsubnet) *Nsubnet = network->cloneshared->Nsubnet;
8672c3e9f7SHong Zhang   PetscFunctionReturn(0);
8772c3e9f7SHong Zhang }
8872c3e9f7SHong Zhang 
8972c3e9f7SHong Zhang /*@
902bf73ac6SHong Zhang   DMNetworkSetNumSubNetworks - Sets the number of subnetworks
915f2c45f1SShri Abhyankar 
92d083f849SBarry Smith   Collective on dm
935f2c45f1SShri Abhyankar 
945f2c45f1SShri Abhyankar   Input Parameters:
955f2c45f1SShri Abhyankar + dm - the dm object
962bf73ac6SHong Zhang . nsubnet - local number of subnetworks
972bf73ac6SHong Zhang - Nsubnet - global number of subnetworks
985f2c45f1SShri Abhyankar 
9997bb938eSShri Abhyankar    Level: beginner
1001b266c99SBarry Smith 
101db781477SPatrick Sanan .seealso: `DMNetworkCreate()`, `DMNetworkGetNumSubNetworks()`
1025f2c45f1SShri Abhyankar @*/
1032bf73ac6SHong Zhang PetscErrorCode DMNetworkSetNumSubNetworks(DM dm,PetscInt nsubnet,PetscInt Nsubnet)
1045f2c45f1SShri Abhyankar {
1055f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
1065f2c45f1SShri Abhyankar 
1075f2c45f1SShri Abhyankar   PetscFunctionBegin;
108*daad07d3SAidan Hamilton   PetscCheck(network->cloneshared->Nsubnet == 0,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network");
1092bf73ac6SHong Zhang 
1105f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1112bf73ac6SHong Zhang   PetscValidLogicalCollectiveInt(dm,nsubnet,2);
1122bf73ac6SHong Zhang   PetscValidLogicalCollectiveInt(dm,Nsubnet,3);
1137765340cSHong Zhang 
1142bf73ac6SHong Zhang   if (Nsubnet == PETSC_DECIDE) {
1155c6496baSHong Zhang     PetscCheck(nsubnet >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of local subnetworks %" PetscInt_FMT " cannot be less than 0",nsubnet);
1161c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&nsubnet,&Nsubnet,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm)));
117caf410d2SHong Zhang   }
1185c6496baSHong Zhang   PetscCheck(Nsubnet >= 1,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Number of global subnetworks %" PetscInt_FMT " cannot be less than 1",Nsubnet);
119caf410d2SHong Zhang 
120*daad07d3SAidan Hamilton   network->cloneshared->Nsubnet  = Nsubnet;
121*daad07d3SAidan Hamilton   network->cloneshared->nsubnet  = 0;       /* initial value; will be determind by DMNetworkAddSubnetwork() */
122*daad07d3SAidan Hamilton   PetscCall(PetscCalloc1(Nsubnet,&network->cloneshared->subnet));
123caf410d2SHong Zhang 
1242bf73ac6SHong Zhang   /* num of shared vertices */
125*daad07d3SAidan Hamilton   network->cloneshared->nsvtx = 0;
126*daad07d3SAidan Hamilton   network->cloneshared->Nsvtx = 0;
1277765340cSHong Zhang   PetscFunctionReturn(0);
1287765340cSHong Zhang }
1297765340cSHong Zhang 
1305f2c45f1SShri Abhyankar /*@
1312bf73ac6SHong Zhang   DMNetworkAddSubnetwork - Add a subnetwork
1325f2c45f1SShri Abhyankar 
1332bf73ac6SHong Zhang   Collective on dm
1345f2c45f1SShri Abhyankar 
1355f2c45f1SShri Abhyankar   Input Parameters:
136e3e68989SHong Zhang + dm - the dm object
1372bf73ac6SHong Zhang . name - name of the subnetwork
1382bf73ac6SHong Zhang . ne - number of local edges of this subnetwork
1398d1cd658SBarry Smith - edgelist - list of edges for this subnetwork, this is a one dimensional array with pairs of entries being the two vertices (in global numbering of the vertices) of each edge,
1408d1cd658SBarry Smith $            [first vertex of first edge, second vertex of first edge, first vertex of second edge, second vertex of second edge, etc]
1412bf73ac6SHong Zhang 
1422bf73ac6SHong Zhang   Output Parameters:
1432bf73ac6SHong Zhang . netnum - global index of the subnetwork
1445f2c45f1SShri Abhyankar 
1455f2c45f1SShri Abhyankar   Notes:
1465f2c45f1SShri Abhyankar   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
1472bf73ac6SHong Zhang   not be destroyed before the call to DMNetworkLayoutSetUp()
1485f2c45f1SShri Abhyankar 
1498d1cd658SBarry Smith   A network can comprise of a single subnetwork OR multiple subnetworks. For a single subnetwork, the subnetwork can be read either in serial or parallel. For a multiple subnetworks,
1508d1cd658SBarry Smith   each subnetwork topology needs to be set on a unique rank and the communicator size needs to be at least equal to the number of subnetworks.
151f11a936eSBarry Smith 
15297bb938eSShri Abhyankar   Level: beginner
1535f2c45f1SShri Abhyankar 
154e3e68989SHong Zhang   Example usage:
155f11a936eSBarry Smith   Consider the following networks:
156f11a936eSBarry Smith   1) A sigle subnetwork:
157e3e68989SHong Zhang .vb
158f11a936eSBarry Smith  network 0:
159f11a936eSBarry Smith  rank[0]:
160f11a936eSBarry Smith    v0 --> v2; v1 --> v2
161f11a936eSBarry Smith  rank[1]:
162f11a936eSBarry Smith    v3 --> v5; v4 --> v5
163e3e68989SHong Zhang .ve
164e3e68989SHong Zhang 
165e3e68989SHong Zhang  The resulting input
166f11a936eSBarry Smith  network 0:
167f11a936eSBarry Smith  rank[0]:
168f11a936eSBarry Smith    ne = 2
169f11a936eSBarry Smith    edgelist = [0 2 | 1 2]
170f11a936eSBarry Smith  rank[1]:
171f11a936eSBarry Smith    ne = 2
172f11a936eSBarry Smith    edgelist = [3 5 | 4 5]
173f11a936eSBarry Smith 
174f11a936eSBarry Smith   2) Two subnetworks:
175f11a936eSBarry Smith .vb
176f11a936eSBarry Smith  subnetwork 0:
177f11a936eSBarry Smith  rank[0]:
178f11a936eSBarry Smith    v0 --> v2; v2 --> v1; v1 --> v3;
179f11a936eSBarry Smith  subnetwork 1:
180f11a936eSBarry Smith  rank[1]:
181f11a936eSBarry Smith    v0 --> v3; v3 --> v2; v2 --> v1;
182f11a936eSBarry Smith .ve
183f11a936eSBarry Smith 
184f11a936eSBarry Smith  The resulting input
185f11a936eSBarry Smith  subnetwork 0:
186f11a936eSBarry Smith  rank[0]:
187f11a936eSBarry Smith    ne = 3
188f11a936eSBarry Smith    edgelist = [0 2 | 2 1 | 1 3]
189f11a936eSBarry Smith  rank[1]:
190f11a936eSBarry Smith    ne = 0
191f11a936eSBarry Smith    edgelist = NULL
192f11a936eSBarry Smith 
193f11a936eSBarry Smith  subnetwork 1:
194f11a936eSBarry Smith  rank[0]:
195f11a936eSBarry Smith    ne = 0
196f11a936eSBarry Smith    edgelist = NULL
197f11a936eSBarry Smith  rank[1]:
198f11a936eSBarry Smith    edgelist = [0 3 | 3 2 | 2 1]
199e3e68989SHong Zhang 
200db781477SPatrick Sanan .seealso: `DMNetworkCreate()`, `DMNetworkSetNumSubnetworks()`
2015f2c45f1SShri Abhyankar @*/
202f11a936eSBarry Smith PetscErrorCode DMNetworkAddSubnetwork(DM dm,const char* name,PetscInt ne,PetscInt edgelist[],PetscInt *netnum)
2035f2c45f1SShri Abhyankar {
2045f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
2055c6496baSHong Zhang   PetscInt       i,Nedge,j,Nvtx,nvtx,nvtx_min=-1,nvtx_max=0;
206f11a936eSBarry Smith   PetscBT        table;
2075f2c45f1SShri Abhyankar 
2085f2c45f1SShri Abhyankar   PetscFunctionBegin;
2098d1cd658SBarry Smith   for (i=0; i<ne; i++) {
2105c6496baSHong Zhang     PetscCheck(edgelist[2*i] != edgelist[2*i+1],PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Edge %" PetscInt_FMT " has the same vertex %" PetscInt_FMT " at each endpoint",i,edgelist[2*i]);
2118d1cd658SBarry Smith   }
212f11a936eSBarry Smith 
2135c6496baSHong Zhang   i = 0;
2145c6496baSHong Zhang   if (ne) nvtx_min = nvtx_max = edgelist[0];
2155c6496baSHong Zhang   for (j=0; j<ne; j++) {
2165c6496baSHong Zhang     nvtx_min = PetscMin(nvtx_min, edgelist[i]);
2175c6496baSHong Zhang     nvtx_max = PetscMax(nvtx_max, edgelist[i]);
2185c6496baSHong Zhang     i++;
2195c6496baSHong Zhang     nvtx_min = PetscMin(nvtx_min, edgelist[i]);
2205c6496baSHong Zhang     nvtx_max = PetscMax(nvtx_max, edgelist[i]);
2215c6496baSHong Zhang     i++;
2225c6496baSHong Zhang   }
2235c6496baSHong Zhang   Nvtx = nvtx_max - nvtx_min + 1; /* approximated total local nvtx for this subnet */
2245c6496baSHong Zhang 
2255c6496baSHong Zhang   /* Get exact local nvtx for this subnet: counting local values between nvtx_min and nvtx_max */
2269566063dSJacob Faibussowitsch   PetscCall(PetscBTCreate(Nvtx,&table));
2279566063dSJacob Faibussowitsch   PetscCall(PetscBTMemzero(Nvtx,table));
228f11a936eSBarry Smith   i = 0;
229f11a936eSBarry Smith   for (j=0; j<ne; j++) {
2309566063dSJacob Faibussowitsch     PetscCall(PetscBTSet(table,edgelist[i++]-nvtx_min));
2319566063dSJacob Faibussowitsch     PetscCall(PetscBTSet(table,edgelist[i++]-nvtx_min));
232f11a936eSBarry Smith   }
233f11a936eSBarry Smith   nvtx = 0;
234f11a936eSBarry Smith   for (j=0; j<Nvtx; j++) {
235f11a936eSBarry Smith     if (PetscBTLookup(table,j)) nvtx++;
236f11a936eSBarry Smith   }
2375c6496baSHong Zhang 
2385c6496baSHong Zhang   /* Get global total Nvtx = max(edgelist[])+1 for this subnet */
2391c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&nvtx_max,&Nvtx,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)dm)));
2405c6496baSHong Zhang   Nvtx++;
2419566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&table));
242f11a936eSBarry Smith 
243f11a936eSBarry Smith   /* Get global total Nedge for this subnet */
2441c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&ne,&Nedge,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm)));
245f11a936eSBarry Smith 
246*daad07d3SAidan Hamilton   i = network->cloneshared->nsubnet;
2472bf73ac6SHong Zhang   if (name) {
248*daad07d3SAidan Hamilton     PetscCall(PetscStrcpy(network->cloneshared->subnet[i].name,name));
249e3e68989SHong Zhang   }
250*daad07d3SAidan Hamilton   network->cloneshared->subnet[i].nvtx     = nvtx; /* include ghost vertices */
251*daad07d3SAidan Hamilton   network->cloneshared->subnet[i].nedge    = ne;
252*daad07d3SAidan Hamilton   network->cloneshared->subnet[i].edgelist = edgelist;
253*daad07d3SAidan Hamilton   network->cloneshared->subnet[i].Nvtx     = Nvtx;
254*daad07d3SAidan Hamilton   network->cloneshared->subnet[i].Nedge    = Nedge;
2552bf73ac6SHong Zhang 
2562bf73ac6SHong Zhang   /* ----------------------------------------------------------
2572bf73ac6SHong Zhang    p=v or e;
2582bf73ac6SHong Zhang    subnet[0].pStart   = 0
2592bf73ac6SHong Zhang    subnet[i+1].pStart = subnet[i].pEnd = subnet[i].pStart + (nE[i] or NV[i])
2602bf73ac6SHong Zhang    ----------------------------------------------------------------------- */
2612bf73ac6SHong Zhang   /* GLOBAL subnet[].vStart and vEnd, used by DMNetworkLayoutSetUp() */
262*daad07d3SAidan Hamilton   network->cloneshared->subnet[i].vStart = network->cloneshared->NVertices;
263*daad07d3SAidan Hamilton   network->cloneshared->subnet[i].vEnd   = network->cloneshared->subnet[i].vStart + network->cloneshared->subnet[i].Nvtx; /* global vEnd of subnet[i] */
2642bf73ac6SHong Zhang 
265*daad07d3SAidan Hamilton   network->cloneshared->nVertices += nvtx; /* include ghost vertices */
266*daad07d3SAidan Hamilton   network->cloneshared->NVertices += network->cloneshared->subnet[i].Nvtx;
2672bf73ac6SHong Zhang 
2682bf73ac6SHong Zhang   /* LOCAL subnet[].eStart and eEnd, used by DMNetworkLayoutSetUp() */
269*daad07d3SAidan Hamilton   network->cloneshared->subnet[i].eStart = network->cloneshared->nEdges;
270*daad07d3SAidan Hamilton   network->cloneshared->subnet[i].eEnd   = network->cloneshared->subnet[i].eStart + ne;
271*daad07d3SAidan Hamilton   network->cloneshared->nEdges += ne;
272*daad07d3SAidan Hamilton   network->cloneshared->NEdges += network->cloneshared->subnet[i].Nedge;
2732bf73ac6SHong Zhang 
274*daad07d3SAidan Hamilton   PetscCall(PetscStrcpy(network->cloneshared->subnet[i].name,name));
275*daad07d3SAidan Hamilton   if (netnum) *netnum = network->cloneshared->nsubnet;
276*daad07d3SAidan Hamilton   network->cloneshared->nsubnet++;
2772bf73ac6SHong Zhang   PetscFunctionReturn(0);
2782bf73ac6SHong Zhang }
2792bf73ac6SHong Zhang 
2805c6496baSHong Zhang /*@C
2815c6496baSHong Zhang   DMNetworkSharedVertexGetInfo - Get info of a shared vertex struct, see petsc/private/dmnetworkimpl.h
2822bf73ac6SHong Zhang 
2835c6496baSHong Zhang   Not collective
2845c6496baSHong Zhang 
2855c6496baSHong Zhang   Input Parameters:
2865c6496baSHong Zhang + dm - the DM object
2875c6496baSHong Zhang - v - vertex point
2885c6496baSHong Zhang 
2895c6496baSHong Zhang   Output Parameters:
2905c6496baSHong Zhang + gidx - global number of this shared vertex in the internal dmplex
2915c6496baSHong Zhang . n - number of subnetworks that share this vertex
2925c6496baSHong Zhang - sv - array of size n: sv[2*i,2*i+1]=(net[i], idx[i]), i=0,...,n-1
2935c6496baSHong Zhang 
2945c6496baSHong Zhang   Level: intermediate
2955c6496baSHong Zhang 
296db781477SPatrick Sanan .seealso: `DMNetworkGetSharedVertices()`
2975c6496baSHong Zhang @*/
2985c6496baSHong Zhang PetscErrorCode DMNetworkSharedVertexGetInfo(DM dm,PetscInt v,PetscInt *gidx,PetscInt *n,const PetscInt **sv)
2992bf73ac6SHong Zhang {
3005c6496baSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
301*daad07d3SAidan Hamilton   SVtx           *svtx = network->cloneshared->svtx;
3025c6496baSHong Zhang   PetscInt       i,gidx_tmp;
3035c6496baSHong Zhang 
3045c6496baSHong Zhang   PetscFunctionBegin;
3059566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetGlobalVertexIndex(dm,v,&gidx_tmp));
306*daad07d3SAidan Hamilton   PetscCall(PetscTableFind(network->cloneshared->svtable,gidx_tmp+1,&i));
3075c6496baSHong Zhang   PetscCheck(i > 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"input vertex is not a shared vertex");
3085c6496baSHong Zhang 
3095c6496baSHong Zhang   i--;
3105c6496baSHong Zhang   if (gidx) *gidx = gidx_tmp;
3115c6496baSHong Zhang   if (n)    *n    = svtx[i].n;
3125c6496baSHong Zhang   if (sv)   *sv   = svtx[i].sv;
3135c6496baSHong Zhang   PetscFunctionReturn(0);
3145c6496baSHong Zhang }
3155c6496baSHong Zhang 
3165c6496baSHong Zhang /*
3175c6496baSHong Zhang   VtxGetInfo - Get info of an input vertex=(net,idx)
3185c6496baSHong Zhang 
3195c6496baSHong Zhang   Input Parameters:
3205c6496baSHong Zhang + Nsvtx - global num of shared vertices
3215c6496baSHong Zhang . svtx - array of shared vertices (global)
3225c6496baSHong Zhang - (net,idx) - subnet number and local index for a vertex
3235c6496baSHong Zhang 
3245c6496baSHong Zhang   Output Parameters:
3255c6496baSHong Zhang + gidx - global index of (net,idx)
3265c6496baSHong Zhang . svtype - see petsc/private/dmnetworkimpl.h
3275c6496baSHong Zhang - svtx_idx - ordering in the svtx array
3285c6496baSHong Zhang */
3295c6496baSHong Zhang static inline PetscErrorCode VtxGetInfo(PetscInt Nsvtx,SVtx *svtx,PetscInt net,PetscInt idx,PetscInt *gidx,SVtxType *svtype,PetscInt *svtx_idx)
3305c6496baSHong Zhang {
3315c6496baSHong Zhang   PetscInt i,j,*svto,g_idx;
3322bf73ac6SHong Zhang   SVtxType vtype;
3332bf73ac6SHong Zhang 
3342bf73ac6SHong Zhang   PetscFunctionBegin;
3352bf73ac6SHong Zhang   if (!Nsvtx) PetscFunctionReturn(0);
3362bf73ac6SHong Zhang 
3375c6496baSHong Zhang   g_idx = -1;
3382bf73ac6SHong Zhang   vtype = SVNONE;
3395c6496baSHong Zhang 
3402bf73ac6SHong Zhang   for (i=0; i<Nsvtx; i++) {
3412bf73ac6SHong Zhang     if (net == svtx[i].sv[0] && idx == svtx[i].sv[1]) {
3425c6496baSHong Zhang       g_idx = svtx[i].gidx;
3432bf73ac6SHong Zhang       vtype = SVFROM;
3442bf73ac6SHong Zhang     } else { /* loop over svtx[i].n */
3452bf73ac6SHong Zhang       for (j=1; j<svtx[i].n; j++) {
3462bf73ac6SHong Zhang         svto = svtx[i].sv + 2*j;
3472bf73ac6SHong Zhang         if (net == svto[0] && idx == svto[1]) {
3482bf73ac6SHong Zhang           /* input vertex net.idx is a shared to_vertex, output its global index and its svtype */
3495c6496baSHong Zhang           g_idx = svtx[i].gidx; /* output gidx for to_vertex */
3502bf73ac6SHong Zhang           vtype = SVTO;
3512bf73ac6SHong Zhang         }
3522bf73ac6SHong Zhang       }
3532bf73ac6SHong Zhang     }
3542bf73ac6SHong Zhang     if (vtype != SVNONE) break;
3552bf73ac6SHong Zhang   }
3565c6496baSHong Zhang   if (gidx)     *gidx     = g_idx;
3572bf73ac6SHong Zhang   if (svtype)   *svtype   = vtype;
3582bf73ac6SHong Zhang   if (svtx_idx) *svtx_idx = i;
3592bf73ac6SHong Zhang   PetscFunctionReturn(0);
3602bf73ac6SHong Zhang }
3612bf73ac6SHong Zhang 
3622bf73ac6SHong Zhang /*
3635c6496baSHong Zhang   TableAddSVtx - Add a new shared vertice from sedgelist[k] to a ctable svta
3645c6496baSHong Zhang 
3655c6496baSHong Zhang   Input:  network, sedgelist, k, svta
3665c6496baSHong Zhang   Output: svta, tdata, ta2sv
3672bf73ac6SHong Zhang */
3685c6496baSHong Zhang static inline PetscErrorCode TableAddSVtx(DM_Network *network,PetscInt *sedgelist,PetscInt k,PetscTable svta,PetscInt* tdata,PetscInt *ta2sv)
3692bf73ac6SHong Zhang {
3705c6496baSHong Zhang   PetscInt       net,idx,gidx;
3712bf73ac6SHong Zhang 
3722bf73ac6SHong Zhang   PetscFunctionBegin;
3735c6496baSHong Zhang   net = sedgelist[k];
3745c6496baSHong Zhang   idx = sedgelist[k+1];
375*daad07d3SAidan Hamilton   gidx = network->cloneshared->subnet[net].vStart + idx;
3769566063dSJacob Faibussowitsch   PetscCall(PetscTableAdd(svta,gidx+1,*tdata+1,INSERT_VALUES));
3775c6496baSHong Zhang 
3785c6496baSHong Zhang   ta2sv[*tdata] = k; /* maps tdata to index of sedgelist */
3795c6496baSHong Zhang   (*tdata)++;
3802bf73ac6SHong Zhang   PetscFunctionReturn(0);
3812bf73ac6SHong Zhang }
3822bf73ac6SHong Zhang 
3832bf73ac6SHong Zhang /*
3845c6496baSHong Zhang   SharedVtxCreate - Create an array of global shared vertices. See SVtx and SVtxType in dmnetworkimpl.h
3852bf73ac6SHong Zhang 
3862bf73ac6SHong Zhang   Input:  dm, Nsedgelist, sedgelist
3872bf73ac6SHong Zhang 
3882bf73ac6SHong Zhang   Note: Output svtx is organized as
3892bf73ac6SHong Zhang         sv(net[0],idx[0]) --> sv(net[1],idx[1])
3902bf73ac6SHong Zhang                           --> sv(net[1],idx[1])
3912bf73ac6SHong Zhang                           ...
3922bf73ac6SHong Zhang                           --> sv(net[n-1],idx[n-1])
3932bf73ac6SHong Zhang         and net[0] < net[1] < ... < net[n-1]
3942bf73ac6SHong Zhang         where sv[0] has SVFROM type, sv[i], i>0, has SVTO type.
3952bf73ac6SHong Zhang  */
3965c6496baSHong Zhang static PetscErrorCode SharedVtxCreate(DM dm,PetscInt Nsedgelist,PetscInt *sedgelist)
3972bf73ac6SHong Zhang {
3985c6496baSHong Zhang   SVtx               *svtx = NULL;
3992bf73ac6SHong Zhang   PetscInt           *sv,k,j,nsv,*tdata,**ta2sv;
4002bf73ac6SHong Zhang   PetscTable         *svtas;
4015c6496baSHong Zhang   PetscInt           gidx,net,idx,i,nta,ita,idx_from,idx_to,n;
4022bf73ac6SHong Zhang   DM_Network         *network = (DM_Network*)dm->data;
4032bf73ac6SHong Zhang   PetscTablePosition ppos;
4042bf73ac6SHong Zhang 
4052bf73ac6SHong Zhang   PetscFunctionBegin;
4065c6496baSHong Zhang   /* (1) Crete an array of ctables svtas to map (net,idx) -> gidx; a svtas[] for a shared/merged vertex */
4079566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(Nsedgelist,&svtas,Nsedgelist,&tdata,2*Nsedgelist,&ta2sv));
4082bf73ac6SHong Zhang 
4092bf73ac6SHong Zhang   k   = 0;   /* sedgelist vertex counter j = 4*k */
4105c6496baSHong Zhang   nta = 0;   /* num of svta tables created */
4112bf73ac6SHong Zhang 
4122bf73ac6SHong Zhang   /* for j=0 */
413*daad07d3SAidan Hamilton   PetscCall(PetscTableCreate(2*Nsedgelist,network->cloneshared->NVertices+1,&svtas[nta]));
4149566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*Nsedgelist,&ta2sv[nta]));
4152bf73ac6SHong Zhang 
4169566063dSJacob Faibussowitsch   PetscCall(TableAddSVtx(network,sedgelist,k,svtas[nta],&tdata[nta],ta2sv[nta]));
4179566063dSJacob Faibussowitsch   PetscCall(TableAddSVtx(network,sedgelist,k+2,svtas[nta],&tdata[nta],ta2sv[nta]));
4182bf73ac6SHong Zhang   nta++; k += 4;
4192bf73ac6SHong Zhang 
4205c6496baSHong Zhang   for (j = 1; j < Nsedgelist; j++) { /* j: sedgelist counter */
4212bf73ac6SHong Zhang     for (ita = 0; ita < nta; ita++) {
4222bf73ac6SHong Zhang       /* vfrom */
4232bf73ac6SHong Zhang       net = sedgelist[k]; idx = sedgelist[k+1];
424*daad07d3SAidan Hamilton       gidx = network->cloneshared->subnet[net].vStart + idx; /* global index of the vertex net.idx before merging shared vertices */
4259566063dSJacob Faibussowitsch       PetscCall(PetscTableFind(svtas[ita],gidx+1,&idx_from));
4262bf73ac6SHong Zhang 
4272bf73ac6SHong Zhang       /* vto */
4282bf73ac6SHong Zhang       net = sedgelist[k+2]; idx = sedgelist[k+3];
429*daad07d3SAidan Hamilton       gidx = network->cloneshared->subnet[net].vStart + idx;
4309566063dSJacob Faibussowitsch       PetscCall(PetscTableFind(svtas[ita],gidx+1,&idx_to));
4312bf73ac6SHong Zhang 
4322bf73ac6SHong Zhang       if (idx_from || idx_to) { /* vfrom or vto is on table svtas[ita] */
4332bf73ac6SHong Zhang         idx_from--; idx_to--;
4342bf73ac6SHong Zhang         if (idx_from < 0) { /* vto is on svtas[ita] */
4359566063dSJacob Faibussowitsch           PetscCall(TableAddSVtx(network,sedgelist,k,svtas[ita],&tdata[ita],ta2sv[ita]));
4362bf73ac6SHong Zhang           break;
4372bf73ac6SHong Zhang         } else if (idx_to < 0) {
4389566063dSJacob Faibussowitsch           PetscCall(TableAddSVtx(network,sedgelist,k+2,svtas[ita],&tdata[ita],ta2sv[ita]));
4392bf73ac6SHong Zhang           break;
4402bf73ac6SHong Zhang         }
4412bf73ac6SHong Zhang       }
4422bf73ac6SHong Zhang     }
4432bf73ac6SHong Zhang 
4442bf73ac6SHong Zhang     if (ita == nta) {
445*daad07d3SAidan Hamilton       PetscCall(PetscTableCreate(2*Nsedgelist,network->cloneshared->NVertices+1,&svtas[nta]));
4469566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2*Nsedgelist, &ta2sv[nta]));
4472bf73ac6SHong Zhang 
4489566063dSJacob Faibussowitsch       PetscCall(TableAddSVtx(network,sedgelist,k,svtas[nta],&tdata[nta],ta2sv[nta]));
4499566063dSJacob Faibussowitsch       PetscCall(TableAddSVtx(network,sedgelist,k+2,svtas[nta],&tdata[nta],ta2sv[nta]));
4502bf73ac6SHong Zhang       nta++;
4512bf73ac6SHong Zhang     }
4522bf73ac6SHong Zhang     k += 4;
4532bf73ac6SHong Zhang   }
4542bf73ac6SHong Zhang 
4556aad120cSJose E. Roman   /* (2) Create svtable for query shared vertices using gidx */
456*daad07d3SAidan Hamilton   PetscCall(PetscTableCreate(nta,network->cloneshared->NVertices+1,&network->cloneshared->svtable));
4575c6496baSHong Zhang 
4585c6496baSHong Zhang   /* (3) Construct svtx from svtas
4595c6496baSHong Zhang      svtx: array of SVtx: sv[0]=(net[0],idx[0]) to vertices sv[k], k=1,...,n-1;
4604f5c2772SJose E. Roman      net[k], k=0, ...,n-1, are in ascending order */
4619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nta,&svtx));
4622bf73ac6SHong Zhang   for (nsv = 0; nsv < nta; nsv++) {
4634f5c2772SJose E. Roman     /* for a single svtx, put shared vertices in ascending order of gidx */
4649566063dSJacob Faibussowitsch     PetscCall(PetscTableGetCount(svtas[nsv],&n));
4659566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(2*n,&sv));
4665c6496baSHong Zhang     svtx[nsv].sv   = sv;
4675c6496baSHong Zhang     svtx[nsv].n    = n;
468*daad07d3SAidan Hamilton     svtx[nsv].gidx = network->cloneshared->NVertices; /* initialization */
4692bf73ac6SHong Zhang 
4709566063dSJacob Faibussowitsch     PetscCall(PetscTableGetHeadPosition(svtas[nsv],&ppos));
4714f5c2772SJose E. Roman     for (k=0; k<n; k++) { /* gidx is sorted in ascending order */
4729566063dSJacob Faibussowitsch       PetscCall(PetscTableGetNext(svtas[nsv],&ppos,&gidx,&i));
4732bf73ac6SHong Zhang       gidx--; i--;
4742bf73ac6SHong Zhang 
4755c6496baSHong Zhang       if (svtx[nsv].gidx > gidx) svtx[nsv].gidx = gidx; /*svtx[nsv].gidx = min(gidx) */
4765c6496baSHong Zhang 
4775c6496baSHong Zhang       j = ta2sv[nsv][i]; /* maps i to index of sedgelist */
4785c6496baSHong Zhang       sv[2*k]   = sedgelist[j];   /* subnet number */
4795c6496baSHong Zhang       sv[2*k+1] = sedgelist[j+1]; /* index on the subnet */
4802bf73ac6SHong Zhang     }
4815c6496baSHong Zhang 
4826aad120cSJose E. Roman     /* Setup svtable for query shared vertices */
483*daad07d3SAidan Hamilton     PetscCall(PetscTableAdd(network->cloneshared->svtable,svtx[nsv].gidx+1,nsv+1,INSERT_VALUES));
4842bf73ac6SHong Zhang   }
4852bf73ac6SHong Zhang 
4862bf73ac6SHong Zhang   for (j=0; j<nta; j++) {
4879566063dSJacob Faibussowitsch     PetscCall(PetscTableDestroy(&svtas[j]));
4889566063dSJacob Faibussowitsch     PetscCall(PetscFree(ta2sv[j]));
4892bf73ac6SHong Zhang   }
4909566063dSJacob Faibussowitsch   PetscCall(PetscFree3(svtas,tdata,ta2sv));
4912bf73ac6SHong Zhang 
492*daad07d3SAidan Hamilton   network->cloneshared->Nsvtx = nta;
493*daad07d3SAidan Hamilton   network->cloneshared->svtx  = svtx;
4942bf73ac6SHong Zhang   PetscFunctionReturn(0);
4952bf73ac6SHong Zhang }
4962bf73ac6SHong Zhang 
4975c6496baSHong Zhang /*
4985c6496baSHong Zhang   GetEdgelist_Coupling - Get an integrated edgelist for dmplex from user-provided subnet[].edgelist when subnets are coupled by shared vertices
4995c6496baSHong Zhang 
5005c6496baSHong Zhang   Input Parameters:
5015c6496baSHong Zhang . dm - the dmnetwork object
5025c6496baSHong Zhang 
5035c6496baSHong Zhang    Output Parameters:
5045c6496baSHong Zhang +  edges - the integrated edgelist for dmplex
5055c6496baSHong Zhang -  nmerged_ptr - num of vertices being merged
5065c6496baSHong Zhang */
5075c6496baSHong Zhang static PetscErrorCode GetEdgelist_Coupling(DM dm,PetscInt *edges,PetscInt *nmerged_ptr)
5082bf73ac6SHong Zhang {
5092bf73ac6SHong Zhang   MPI_Comm       comm;
5102bf73ac6SHong Zhang   PetscMPIInt    size,rank,*recvcounts=NULL,*displs=NULL;
511eac198afSGetnet   DM_Network     *network = (DM_Network*)dm->data;
512eac198afSGetnet   PetscInt       i,j,ctr,np;
513*daad07d3SAidan Hamilton   PetscInt       *vidxlTog,Nsv,Nsubnet=network->cloneshared->Nsubnet;
514*daad07d3SAidan Hamilton   PetscInt       *sedgelist=network->cloneshared->sedgelist;
5155c6496baSHong Zhang   PetscInt       net,idx,gidx,nmerged,*vrange,gidx_from,net_from,sv_idx;
5162bf73ac6SHong Zhang   SVtxType       svtype = SVNONE;
5175c6496baSHong Zhang   SVtx           *svtx;
5182bf73ac6SHong Zhang 
5192bf73ac6SHong Zhang   PetscFunctionBegin;
5209566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
5219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
5229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
5232bf73ac6SHong Zhang 
5245c6496baSHong Zhang   /* (1) Create global svtx[] from sedgelist */
5255c6496baSHong Zhang   /* --------------------------------------- */
526*daad07d3SAidan Hamilton   PetscCall(SharedVtxCreate(dm,network->cloneshared->Nsvtx,sedgelist));
527*daad07d3SAidan Hamilton   Nsv  = network->cloneshared->Nsvtx;
528*daad07d3SAidan Hamilton   svtx = network->cloneshared->svtx;
5292bf73ac6SHong Zhang 
5305c6496baSHong Zhang   /* (2) Merge shared vto vertices to their vfrom vertex with same global vetex index (gidx) */
5315c6496baSHong Zhang   /* --------------------------------------------------------------------------------------- */
5322bf73ac6SHong Zhang   /* (2.1) compute vrage[rank]: global index of 1st local vertex in proc[rank] */
533*daad07d3SAidan Hamilton   PetscCall(PetscMalloc4(size+1,&vrange,size,&displs,size,&recvcounts,network->cloneshared->nVertices,&vidxlTog));
5342bf73ac6SHong Zhang   for (i=0; i<size; i++) {displs[i] = i; recvcounts[i] = 1;}
5352bf73ac6SHong Zhang 
5362bf73ac6SHong Zhang   vrange[0] = 0;
537*daad07d3SAidan Hamilton   PetscCallMPI(MPI_Allgatherv(&network->cloneshared->nVertices,1,MPIU_INT,vrange+1,recvcounts,displs,MPIU_INT,comm));
5385c6496baSHong Zhang   for (i=2; i<size+1; i++) vrange[i] += vrange[i-1];
5392bf73ac6SHong Zhang 
5402bf73ac6SHong Zhang   /* (2.2) Create vidxlTog: maps UN-MERGED local vertex index i to global index gidx (plex, excluding ghost vertices) */
5412bf73ac6SHong Zhang   i = 0; gidx = 0;
5422bf73ac6SHong Zhang   nmerged        = 0; /* local num of merged vertices */
543*daad07d3SAidan Hamilton   network->cloneshared->nsvtx = 0; /* local num of SVtx structs, including ghosts */
5442bf73ac6SHong Zhang   for (net=0; net<Nsubnet; net++) {
545*daad07d3SAidan Hamilton     for (idx=0; idx<network->cloneshared->subnet[net].Nvtx; idx++) { /* Note: global subnet[net].Nvtx */
5469566063dSJacob Faibussowitsch       PetscCall(VtxGetInfo(Nsv,svtx,net,idx,&gidx_from,&svtype,&sv_idx));
5472bf73ac6SHong Zhang       if (svtype == SVTO) {
548*daad07d3SAidan Hamilton         if (network->cloneshared->subnet[net].nvtx) {/* this proc owns sv_to */
5495c6496baSHong Zhang           net_from = svtx[sv_idx].sv[0]; /* subnet number of its shared vertex */
550*daad07d3SAidan Hamilton           if (network->cloneshared->subnet[net_from].nvtx == 0) {
5515c6496baSHong Zhang             /* this proc does not own v_from, thus a ghost local vertex */
552*daad07d3SAidan Hamilton             network->cloneshared->nsvtx++;
5532bf73ac6SHong Zhang           }
5545c6496baSHong Zhang           vidxlTog[i++] = gidx_from; /* gidx before merging! Bug??? */
5555c6496baSHong Zhang           nmerged++; /* a shared vertex -- merged */
5562bf73ac6SHong Zhang         }
5572bf73ac6SHong Zhang       } else {
558*daad07d3SAidan Hamilton         if (svtype == SVFROM && network->cloneshared->subnet[net].nvtx) {
5595c6496baSHong Zhang           /* this proc owns this v_from, a new local shared vertex */
560*daad07d3SAidan Hamilton           network->cloneshared->nsvtx++;
5612bf73ac6SHong Zhang         }
562*daad07d3SAidan Hamilton         if (network->cloneshared->subnet[net].nvtx) vidxlTog[i++] = gidx;
5632bf73ac6SHong Zhang         gidx++;
5642bf73ac6SHong Zhang       }
5652bf73ac6SHong Zhang     }
5662bf73ac6SHong Zhang   }
5672bf73ac6SHong Zhang #if defined(PETSC_USE_DEBUG)
568*daad07d3SAidan Hamilton   PetscCheck(i == network->cloneshared->nVertices,PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"%" PetscInt_FMT " != %" PetscInt_FMT " nVertices",i,network->cloneshared->nVertices);
5692bf73ac6SHong Zhang #endif
5702bf73ac6SHong Zhang 
5715c6496baSHong Zhang   /* (2.3) Shared vertices in the subnetworks are merged, update global NVertices: np = sum(local nmerged) */
5729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&nmerged,&np,1,MPIU_INT,MPI_SUM,comm));
573*daad07d3SAidan Hamilton   network->cloneshared->NVertices -= np;
5742bf73ac6SHong Zhang 
5752bf73ac6SHong Zhang   ctr = 0;
5762bf73ac6SHong Zhang   for (net=0; net<Nsubnet; net++) {
577*daad07d3SAidan Hamilton     for (j = 0; j < network->cloneshared->subnet[net].nedge; j++) {
5782bf73ac6SHong Zhang       /* vfrom: */
579*daad07d3SAidan Hamilton       i = network->cloneshared->subnet[net].edgelist[2*j] + (network->cloneshared->subnet[net].vStart - vrange[rank]);
5802bf73ac6SHong Zhang       edges[2*ctr] = vidxlTog[i];
5812bf73ac6SHong Zhang 
5822bf73ac6SHong Zhang       /* vto */
583*daad07d3SAidan Hamilton       i = network->cloneshared->subnet[net].edgelist[2*j+1] + (network->cloneshared->subnet[net].vStart - vrange[rank]);
5842bf73ac6SHong Zhang       edges[2*ctr+1] = vidxlTog[i];
5852bf73ac6SHong Zhang       ctr++;
5862bf73ac6SHong Zhang     }
5872bf73ac6SHong Zhang   }
5889566063dSJacob Faibussowitsch   PetscCall(PetscFree4(vrange,displs,recvcounts,vidxlTog));
5899566063dSJacob Faibussowitsch   PetscCall(PetscFree(sedgelist)); /* created in DMNetworkAddSharedVertices() */
5902bf73ac6SHong Zhang 
591eac198afSGetnet   *nmerged_ptr = nmerged;
5925f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5935f2c45f1SShri Abhyankar }
5945f2c45f1SShri Abhyankar 
595*daad07d3SAidan Hamilton PetscErrorCode DMNetworkInitializeNonTopological(DM dm)
596*daad07d3SAidan Hamilton {
597*daad07d3SAidan Hamilton   DM_Network     *network = (DM_Network*)dm->data;
598*daad07d3SAidan Hamilton   PetscInt       p,pStart = network->cloneshared->pStart, pEnd = network->cloneshared->pEnd;
599*daad07d3SAidan Hamilton   MPI_Comm       comm;
600*daad07d3SAidan Hamilton 
601*daad07d3SAidan Hamilton   PetscFunctionBegin;
602*daad07d3SAidan Hamilton   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
603*daad07d3SAidan Hamilton 
604*daad07d3SAidan Hamilton   PetscCall(PetscSectionCreate(comm,&network->DataSection));
605*daad07d3SAidan Hamilton   PetscCall(PetscSectionCreate(comm,&network->DofSection));
606*daad07d3SAidan Hamilton   PetscCall(PetscSectionSetChart(network->DataSection,pStart,pEnd));
607*daad07d3SAidan Hamilton   PetscCall(PetscSectionSetChart(network->DofSection,pStart,pEnd));
608*daad07d3SAidan Hamilton 
609*daad07d3SAidan Hamilton   PetscCall(DMNetworkInitializeHeaderComponentData(dm));
610*daad07d3SAidan Hamilton 
611*daad07d3SAidan Hamilton   for (p = 0; p <pEnd-pStart; p++) {
612*daad07d3SAidan Hamilton     network->header[p].ndata           = 0;
613*daad07d3SAidan Hamilton     network->header[p].offset[0]       = 0;
614*daad07d3SAidan Hamilton     network->header[p].offsetvarrel[0] = 0;
615*daad07d3SAidan Hamilton     PetscCall(PetscSectionAddDof(network->DataSection,p,network->header[p].hsize));
616*daad07d3SAidan Hamilton   }
617*daad07d3SAidan Hamilton   PetscFunctionReturn(0);
618*daad07d3SAidan Hamilton }
619*daad07d3SAidan Hamilton 
6205f2c45f1SShri Abhyankar /*@
6215f2c45f1SShri Abhyankar   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
6225f2c45f1SShri Abhyankar 
623eac198afSGetnet   Not Collective
6245f2c45f1SShri Abhyankar 
6257a7aea1fSJed Brown   Input Parameters:
626f11a936eSBarry Smith . dm - the dmnetwork object
6275f2c45f1SShri Abhyankar 
6285f2c45f1SShri Abhyankar   Notes:
6295f2c45f1SShri Abhyankar   This routine should be called after the network sizes and edgelists have been provided. It creates
6305f2c45f1SShri Abhyankar   the bare layout of the network and sets up the network to begin insertion of components.
6315f2c45f1SShri Abhyankar 
6325f2c45f1SShri Abhyankar   All the components should be registered before calling this routine.
6335f2c45f1SShri Abhyankar 
63497bb938eSShri Abhyankar   Level: beginner
6355f2c45f1SShri Abhyankar 
636db781477SPatrick Sanan .seealso: `DMNetworkSetNumSubNetworks()`, `DMNetworkAddSubnetwork()`
6375f2c45f1SShri Abhyankar @*/
6385f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm)
6395f2c45f1SShri Abhyankar {
6405f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
641*daad07d3SAidan Hamilton   PetscInt       i,j,ctr,Nsubnet=network->cloneshared->Nsubnet,np,*edges,*subnetvtx,*subnetedge,e,v,vfrom,vto,net,globaledgeoff;
642caf410d2SHong Zhang   const PetscInt *cone;
643caf410d2SHong Zhang   MPI_Comm       comm;
6445503a911SAidan Hamilton   PetscMPIInt    size;
645eac198afSGetnet   PetscSection   sectiong;
6465c6496baSHong Zhang   PetscInt       nmerged=0;
6476500d4abSHong Zhang 
6486500d4abSHong Zhang   PetscFunctionBegin;
649e3b1a053SGetnet Betrie   PetscCall(PetscLogEventBegin(DMNetwork_LayoutSetUp,dm,0,0,0));
650*daad07d3SAidan Hamilton   PetscCheck(network->cloneshared->nsubnet == Nsubnet,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Must call DMNetworkAddSubnetwork() %" PetscInt_FMT " times",Nsubnet);
6512bf73ac6SHong Zhang 
652eac198afSGetnet   /* This implementation requires user input each subnet by a single processor when Nsubnet>1, thus subnet[net].nvtx=subnet[net].Nvtx when net>0 */
653eac198afSGetnet   for (net=1; net<Nsubnet; net++) {
654*daad07d3SAidan Hamilton     if (network->cloneshared->subnet[net].nvtx) PetscCheck(network->cloneshared->subnet[net].nvtx == network->cloneshared->subnet[net].Nvtx,PETSC_COMM_SELF,PETSC_ERR_SUP,"subnetwork %" PetscInt_FMT " local num of vertices %" PetscInt_FMT " != %" PetscInt_FMT " global num",net,network->cloneshared->subnet[net].nvtx,network->cloneshared->subnet[net].Nvtx);
655eac198afSGetnet   }
656eac198afSGetnet 
6579566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
6585503a911SAidan Hamilton   PetscCall(MPI_Comm_size(comm,&size));
6596500d4abSHong Zhang 
660f11a936eSBarry Smith   /* Create LOCAL edgelist in global vertex ordering for the network by concatenating local input edgelists of the subnetworks */
661*daad07d3SAidan Hamilton   PetscCall(PetscCalloc1(2*network->cloneshared->nEdges,&edges));
662eac198afSGetnet 
663*daad07d3SAidan Hamilton   if (network->cloneshared->Nsvtx) { /* subnetworks are coupled via shared vertices */
6649566063dSJacob Faibussowitsch     PetscCall(GetEdgelist_Coupling(dm,edges,&nmerged));
665eac198afSGetnet   } else { /* subnetworks are not coupled */
6666aad120cSJose E. Roman     /* Create a 0-size svtable for query shared vertices */
667*daad07d3SAidan Hamilton     PetscCall(PetscTableCreate(0,network->cloneshared->NVertices+1,&network->cloneshared->svtable));
668caf410d2SHong Zhang     ctr = 0;
6692bf73ac6SHong Zhang     for (i=0; i < Nsubnet; i++) {
670*daad07d3SAidan Hamilton       for (j = 0; j < network->cloneshared->subnet[i].nedge; j++) {
671*daad07d3SAidan Hamilton         edges[2*ctr]   = network->cloneshared->subnet[i].vStart + network->cloneshared->subnet[i].edgelist[2*j];
672*daad07d3SAidan Hamilton         edges[2*ctr+1] = network->cloneshared->subnet[i].vStart + network->cloneshared->subnet[i].edgelist[2*j+1];
6736500d4abSHong Zhang         ctr++;
6746500d4abSHong Zhang       }
6756500d4abSHong Zhang     }
676eac198afSGetnet   }
6777765340cSHong Zhang 
6782bf73ac6SHong Zhang   /* Create network->plex; One dimensional network, numCorners=2 */
6799566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm,&network->plex));
6809566063dSJacob Faibussowitsch   PetscCall(DMSetType(network->plex,DMPLEX));
6819566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(network->plex,1));
682eac198afSGetnet 
683*daad07d3SAidan Hamilton   if (size == 1) PetscCall(DMPlexBuildFromCellList(network->plex,network->cloneshared->nEdges,PETSC_DECIDE,2,edges));
684*daad07d3SAidan Hamilton   else PetscCall(DMPlexBuildFromCellListParallel(network->plex,network->cloneshared->nEdges,PETSC_DECIDE,PETSC_DECIDE,2,edges,NULL, NULL));
6856500d4abSHong Zhang 
686*daad07d3SAidan Hamilton   PetscCall(DMPlexGetChart(network->plex,&network->cloneshared->pStart,&network->cloneshared->pEnd));
687*daad07d3SAidan Hamilton   PetscCall(DMPlexGetHeightStratum(network->plex,0,&network->cloneshared->eStart,&network->cloneshared->eEnd));
688*daad07d3SAidan Hamilton   PetscCall(DMPlexGetHeightStratum(network->plex,1,&network->cloneshared->vStart,&network->cloneshared->vEnd));
689*daad07d3SAidan Hamilton   np = network->cloneshared->pEnd - network->cloneshared->pStart;
6909566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(np,&network->header,np,&network->cvalue));
691caf410d2SHong Zhang 
692f11a936eSBarry Smith   /* Create edge and vertex arrays for the subnetworks
693f11a936eSBarry Smith      This implementation assumes that DMNetwork reads
694f11a936eSBarry Smith      (1) a single subnetwork in parallel; or
695f11a936eSBarry Smith      (2) n subnetworks using n processors, one subnetwork/processor.
696f11a936eSBarry Smith   */
697*daad07d3SAidan Hamilton   PetscCall(PetscCalloc2(network->cloneshared->nEdges,&subnetedge,network->cloneshared->nVertices+network->cloneshared->nsvtx,&subnetvtx)); /* Maps local edge/vertex to local subnetwork's edge/vertex */
698*daad07d3SAidan Hamilton   network->cloneshared->subnetedge = subnetedge;
699*daad07d3SAidan Hamilton   network->cloneshared->subnetvtx  = subnetvtx;
7005c6496baSHong Zhang   for (j=0; j < Nsubnet; j++) {
701*daad07d3SAidan Hamilton     network->cloneshared->subnet[j].edges = subnetedge;
702*daad07d3SAidan Hamilton     subnetedge              += network->cloneshared->subnet[j].nedge;
703f11a936eSBarry Smith 
704*daad07d3SAidan Hamilton     network->cloneshared->subnet[j].vertices = subnetvtx;
705*daad07d3SAidan Hamilton     subnetvtx                  += network->cloneshared->subnet[j].nvtx;
7066500d4abSHong Zhang   }
707*daad07d3SAidan Hamilton   network->cloneshared->svertices = subnetvtx;
708caf410d2SHong Zhang 
709caf410d2SHong Zhang   /* Get edge ownership */
710*daad07d3SAidan Hamilton   np = network->cloneshared->eEnd - network->cloneshared->eStart;
7115503a911SAidan Hamilton   PetscCallMPI(MPI_Scan(&np,&globaledgeoff,1,MPIU_INT,MPI_SUM,comm));
7125503a911SAidan Hamilton   globaledgeoff -= np;
713caf410d2SHong Zhang 
714eac198afSGetnet   /* Setup local edge and vertex arrays for subnetworks */
7152bf73ac6SHong Zhang   e = 0;
7162bf73ac6SHong Zhang   for (i=0; i < Nsubnet; i++) {
717f11a936eSBarry Smith     ctr = 0;
718*daad07d3SAidan Hamilton     for (j = 0; j < network->cloneshared->subnet[i].nedge; j++) {
7192bf73ac6SHong Zhang       /* edge e */
7205503a911SAidan Hamilton       network->header[e].index    = e + globaledgeoff;   /* Global edge index */
7212bf73ac6SHong Zhang       network->header[e].subnetid = i;
722*daad07d3SAidan Hamilton       network->cloneshared->subnet[i].edges[j] = e;
7232bf73ac6SHong Zhang 
7242bf73ac6SHong Zhang       /* connected vertices */
7259566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(network->plex,e,&cone));
7262bf73ac6SHong Zhang 
7272bf73ac6SHong Zhang       /* vertex cone[0] */
728f11a936eSBarry Smith       v = cone[0];
729f11a936eSBarry Smith       network->header[v].index     = edges[2*e];  /* Global vertex index */
730f11a936eSBarry Smith       network->header[v].subnetid  = i;           /* Subnetwork id */
731f11a936eSBarry Smith       if (Nsubnet == 1) {
732*daad07d3SAidan Hamilton         network->cloneshared->subnet[i].vertices[v - network->cloneshared->vStart] = v; /* user's subnet[].idx = petsc's v */
733f11a936eSBarry Smith       } else {
734*daad07d3SAidan Hamilton         vfrom = network->cloneshared->subnet[i].edgelist[2*ctr];     /* =subnet[i].idx, Global index! */
735*daad07d3SAidan Hamilton         network->cloneshared->subnet[i].vertices[vfrom] = v; /* user's subnet[].dix = petsc's v */
736f11a936eSBarry Smith       }
7372bf73ac6SHong Zhang 
7382bf73ac6SHong Zhang       /* vertex cone[1] */
739f11a936eSBarry Smith       v = cone[1];
740f11a936eSBarry Smith       network->header[v].index    = edges[2*e+1];   /* Global vertex index */
741f11a936eSBarry Smith       network->header[v].subnetid = i;              /* Subnetwork id */
742f11a936eSBarry Smith       if (Nsubnet == 1) {
743*daad07d3SAidan Hamilton         network->cloneshared->subnet[i].vertices[v - network->cloneshared->vStart] = v; /* user's subnet[].idx = petsc's v */
744f11a936eSBarry Smith       } else {
745*daad07d3SAidan Hamilton         vto = network->cloneshared->subnet[i].edgelist[2*ctr+1];     /* =subnet[i].idx, Global index! */
746*daad07d3SAidan Hamilton         network->cloneshared->subnet[i].vertices[vto] = v; /* user's subnet[].dix = petsc's v */
747f11a936eSBarry Smith       }
7482bf73ac6SHong Zhang 
749f11a936eSBarry Smith       e++; ctr++;
7506500d4abSHong Zhang     }
7516500d4abSHong Zhang   }
7525503a911SAidan Hamilton   PetscCall(PetscFree(edges));
753caf410d2SHong Zhang 
754eac198afSGetnet   /* Set local vertex array for the subnetworks */
755eac198afSGetnet   j = 0;
756*daad07d3SAidan Hamilton   for (v = network->cloneshared->vStart; v < network->cloneshared->vEnd; v++) {
757eac198afSGetnet     /* local shared vertex */
758*daad07d3SAidan Hamilton     PetscCall(PetscTableFind(network->cloneshared->svtable,network->header[v].index+1,&i));
759*daad07d3SAidan Hamilton     if (i) network->cloneshared->svertices[j++] = v;
7606500d4abSHong Zhang   }
761eac198afSGetnet 
762eac198afSGetnet   /* Create a global section to be used by DMNetworkIsGhostVertex() which is a non-collective routine */
763eac198afSGetnet   /* see snes_tutorials_network-ex1_4 */
7649566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(network->plex,&sectiong));
765*daad07d3SAidan Hamilton   /* Initialize non-topological data structures  */
766*daad07d3SAidan Hamilton   PetscCall(DMNetworkInitializeNonTopological(dm));
767e3b1a053SGetnet Betrie   PetscCall(PetscLogEventEnd(DMNetwork_LayoutSetUp,dm,0,0,0));
7685f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7695f2c45f1SShri Abhyankar }
7705f2c45f1SShri Abhyankar 
77194ef8ddeSSatish Balay /*@C
7722bf73ac6SHong Zhang   DMNetworkGetSubnetwork - Returns the information about a requested subnetwork
7732bf73ac6SHong Zhang 
7742bf73ac6SHong Zhang   Not collective
7752727e31bSShri Abhyankar 
7767a7aea1fSJed Brown   Input Parameters:
777caf410d2SHong Zhang + dm - the DM object
7782bf73ac6SHong Zhang - netnum - the global index of the subnetwork
7792727e31bSShri Abhyankar 
7807a7aea1fSJed Brown   Output Parameters:
7812727e31bSShri Abhyankar + nv - number of vertices (local)
7822727e31bSShri Abhyankar . ne - number of edges (local)
7832bf73ac6SHong Zhang . vtx - local vertices of the subnetwork
7842bf73ac6SHong Zhang - edge - local edges of the subnetwork
7852727e31bSShri Abhyankar 
7862727e31bSShri Abhyankar   Notes:
7872727e31bSShri Abhyankar     Cannot call this routine before DMNetworkLayoutSetup()
7882727e31bSShri Abhyankar 
789eac198afSGetnet     The local vertices returned on each rank are determined by DMNetwork. The user does not have any control over what vertices are local.
790eac198afSGetnet 
79106dd6b0eSSatish Balay   Level: intermediate
79206dd6b0eSSatish Balay 
793db781477SPatrick Sanan .seealso: `DMNetworkCreate()`, `DMNetworkAddSubnetwork()`, `DMNetworkLayoutSetUp()`
7942727e31bSShri Abhyankar @*/
7952bf73ac6SHong Zhang PetscErrorCode DMNetworkGetSubnetwork(DM dm,PetscInt netnum,PetscInt *nv,PetscInt *ne,const PetscInt **vtx,const PetscInt **edge)
7962727e31bSShri Abhyankar {
797caf410d2SHong Zhang   DM_Network *network = (DM_Network*)dm->data;
7982727e31bSShri Abhyankar 
7992727e31bSShri Abhyankar   PetscFunctionBegin;
800*daad07d3SAidan Hamilton   PetscCheck(netnum < network->cloneshared->Nsubnet,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Subnet index %" PetscInt_FMT " exceeds the num of subnets %" PetscInt_FMT,netnum,network->cloneshared->Nsubnet);
801*daad07d3SAidan Hamilton   if (nv) *nv     = network->cloneshared->subnet[netnum].nvtx;
802*daad07d3SAidan Hamilton   if (ne) *ne     = network->cloneshared->subnet[netnum].nedge;
803*daad07d3SAidan Hamilton   if (vtx) *vtx   = network->cloneshared->subnet[netnum].vertices;
804*daad07d3SAidan Hamilton   if (edge) *edge = network->cloneshared->subnet[netnum].edges;
8052bf73ac6SHong Zhang   PetscFunctionReturn(0);
8062bf73ac6SHong Zhang }
8072bf73ac6SHong Zhang 
8082bf73ac6SHong Zhang /*@
8092bf73ac6SHong Zhang   DMNetworkAddSharedVertices - Add shared vertices that connect two given subnetworks
8102bf73ac6SHong Zhang 
8112bf73ac6SHong Zhang   Collective on dm
8122bf73ac6SHong Zhang 
8132bf73ac6SHong Zhang   Input Parameters:
8142bf73ac6SHong Zhang + dm - the dm object
8152bf73ac6SHong Zhang . anetnum - first subnetwork global numbering returned by DMNetworkAddSubnetwork()
8162bf73ac6SHong Zhang . bnetnum - second subnetwork global numbering returned by DMNetworkAddSubnetwork()
8172bf73ac6SHong Zhang . nsvtx - number of vertices that are shared by the two subnetworks
8182bf73ac6SHong Zhang . asvtx - vertex index in the first subnetwork
8192bf73ac6SHong Zhang - bsvtx - vertex index in the second subnetwork
8202bf73ac6SHong Zhang 
8212bf73ac6SHong Zhang   Level: beginner
8222bf73ac6SHong Zhang 
823db781477SPatrick Sanan .seealso: `DMNetworkCreate()`, `DMNetworkAddSubnetwork()`, `DMNetworkGetSharedVertices()`
8242bf73ac6SHong Zhang @*/
8252bf73ac6SHong Zhang PetscErrorCode DMNetworkAddSharedVertices(DM dm,PetscInt anetnum,PetscInt bnetnum,PetscInt nsvtx,PetscInt asvtx[],PetscInt bsvtx[])
8262bf73ac6SHong Zhang {
8272bf73ac6SHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
828*daad07d3SAidan Hamilton   PetscInt       i,nsubnet = network->cloneshared->Nsubnet,*sedgelist,Nsvtx=network->cloneshared->Nsvtx;
8292bf73ac6SHong Zhang 
8302bf73ac6SHong Zhang   PetscFunctionBegin;
8315c6496baSHong Zhang   PetscCheck(anetnum != bnetnum,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Subnetworks must have different netnum");
8325c6496baSHong Zhang   PetscCheck(anetnum >= 0 && bnetnum >= 0,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"netnum cannot be negative");
8332bf73ac6SHong Zhang   if (!Nsvtx) {
8342bf73ac6SHong Zhang     /* allocate network->sedgelist to hold at most 2*nsubnet pairs of shared vertices */
835*daad07d3SAidan Hamilton     PetscCall(PetscMalloc1(2*4*nsubnet,&network->cloneshared->sedgelist));
8362bf73ac6SHong Zhang   }
8372bf73ac6SHong Zhang 
838*daad07d3SAidan Hamilton   sedgelist = network->cloneshared->sedgelist;
8392bf73ac6SHong Zhang   for (i=0; i<nsvtx; i++) {
8402bf73ac6SHong Zhang     sedgelist[4*Nsvtx]   = anetnum; sedgelist[4*Nsvtx+1] = asvtx[i];
8412bf73ac6SHong Zhang     sedgelist[4*Nsvtx+2] = bnetnum; sedgelist[4*Nsvtx+3] = bsvtx[i];
8422bf73ac6SHong Zhang     Nsvtx++;
8432bf73ac6SHong Zhang   }
8445c6496baSHong Zhang   PetscCheck(Nsvtx <= 2*nsubnet,PETSC_COMM_SELF,PETSC_ERR_SUP,"allocate more space for coupling edgelist");
845*daad07d3SAidan Hamilton   network->cloneshared->Nsvtx = Nsvtx;
8462727e31bSShri Abhyankar   PetscFunctionReturn(0);
8472727e31bSShri Abhyankar }
8482727e31bSShri Abhyankar 
8492727e31bSShri Abhyankar /*@C
8502bf73ac6SHong Zhang   DMNetworkGetSharedVertices - Returns the info for the shared vertices
8512bf73ac6SHong Zhang 
8522bf73ac6SHong Zhang   Not collective
853caf410d2SHong Zhang 
854f899ff85SJose E. Roman   Input Parameter:
8552bf73ac6SHong Zhang . dm - the DM object
856caf410d2SHong Zhang 
8577a7aea1fSJed Brown   Output Parameters:
8582bf73ac6SHong Zhang + nsv - number of local shared vertices
8592bf73ac6SHong Zhang - svtx - local shared vertices
860caf410d2SHong Zhang 
861caf410d2SHong Zhang   Notes:
862caf410d2SHong Zhang   Cannot call this routine before DMNetworkLayoutSetup()
863caf410d2SHong Zhang 
864caf410d2SHong Zhang   Level: intermediate
865caf410d2SHong Zhang 
866db781477SPatrick Sanan .seealso: `DMNetworkGetSubnetwork()`, `DMNetworkLayoutSetUp()`, `DMNetworkAddSharedVertices()`
867caf410d2SHong Zhang @*/
8682bf73ac6SHong Zhang PetscErrorCode DMNetworkGetSharedVertices(DM dm,PetscInt *nsv,const PetscInt **svtx)
869caf410d2SHong Zhang {
870caf410d2SHong Zhang   DM_Network *net = (DM_Network*)dm->data;
871caf410d2SHong Zhang 
872caf410d2SHong Zhang   PetscFunctionBegin;
8735c6496baSHong Zhang   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
874*daad07d3SAidan Hamilton   if (nsv)  *nsv  = net->cloneshared->nsvtx;
875*daad07d3SAidan Hamilton   if (svtx) *svtx = net->cloneshared->svertices;
876caf410d2SHong Zhang   PetscFunctionReturn(0);
877caf410d2SHong Zhang }
878caf410d2SHong Zhang 
879caf410d2SHong Zhang /*@C
8805f2c45f1SShri Abhyankar   DMNetworkRegisterComponent - Registers the network component
8815f2c45f1SShri Abhyankar 
882d083f849SBarry Smith   Logically collective on dm
8835f2c45f1SShri Abhyankar 
8847a7aea1fSJed Brown   Input Parameters:
8855f2c45f1SShri Abhyankar + dm - the network object
8865f2c45f1SShri Abhyankar . name - the component name
8875f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data
8885f2c45f1SShri Abhyankar 
8897a7aea1fSJed Brown    Output Parameters:
8905f2c45f1SShri Abhyankar .  key - an integer key that defines the component
8915f2c45f1SShri Abhyankar 
8925f2c45f1SShri Abhyankar    Notes
8935f2c45f1SShri Abhyankar    This routine should be called by all processors before calling DMNetworkLayoutSetup().
8945f2c45f1SShri Abhyankar 
89597bb938eSShri Abhyankar    Level: beginner
8965f2c45f1SShri Abhyankar 
897db781477SPatrick Sanan .seealso: `DMNetworkCreate()`, `DMNetworkLayoutSetUp()`
8985f2c45f1SShri Abhyankar @*/
899caf410d2SHong Zhang PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,size_t size,PetscInt *key)
9005f2c45f1SShri Abhyankar {
9015f2c45f1SShri Abhyankar   DM_Network            *network = (DM_Network*) dm->data;
90254dfd506SHong Zhang   DMNetworkComponent    *component=NULL,*newcomponent=NULL;
9035f2c45f1SShri Abhyankar   PetscBool             flg=PETSC_FALSE;
9045f2c45f1SShri Abhyankar   PetscInt              i;
9055f2c45f1SShri Abhyankar 
9065f2c45f1SShri Abhyankar   PetscFunctionBegin;
90754dfd506SHong Zhang   if (!network->component) {
9089566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(network->max_comps_registered,&network->component));
90954dfd506SHong Zhang   }
91054dfd506SHong Zhang 
9115f2c45f1SShri Abhyankar   for (i=0; i < network->ncomponent; i++) {
9129566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(network->component[i].name,name,&flg));
9135f2c45f1SShri Abhyankar     if (flg) {
9145f2c45f1SShri Abhyankar       *key = i;
9155f2c45f1SShri Abhyankar       PetscFunctionReturn(0);
9165f2c45f1SShri Abhyankar     }
9176d64e262SShri Abhyankar   }
91854dfd506SHong Zhang 
91954dfd506SHong Zhang   if (network->ncomponent == network->max_comps_registered) {
92054dfd506SHong Zhang     /* Reached max allowed so resize component */
92154dfd506SHong Zhang     network->max_comps_registered += 2;
9229566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(network->max_comps_registered,&newcomponent));
92354dfd506SHong Zhang     /* Copy over the previous component info */
92454dfd506SHong Zhang     for (i=0; i < network->ncomponent; i++) {
9259566063dSJacob Faibussowitsch       PetscCall(PetscStrcpy(newcomponent[i].name,network->component[i].name));
92654dfd506SHong Zhang       newcomponent[i].size = network->component[i].size;
9275f2c45f1SShri Abhyankar     }
92854dfd506SHong Zhang     /* Free old one */
9299566063dSJacob Faibussowitsch     PetscCall(PetscFree(network->component));
93054dfd506SHong Zhang     /* Update pointer */
93154dfd506SHong Zhang     network->component = newcomponent;
93254dfd506SHong Zhang   }
93354dfd506SHong Zhang 
93454dfd506SHong Zhang   component = &network->component[network->ncomponent];
9355f2c45f1SShri Abhyankar 
9369566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(component->name,name));
9375f2c45f1SShri Abhyankar   component->size = size/sizeof(DMNetworkComponentGenericDataType);
9385f2c45f1SShri Abhyankar   *key = network->ncomponent;
9395f2c45f1SShri Abhyankar   network->ncomponent++;
9405f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9415f2c45f1SShri Abhyankar }
9425f2c45f1SShri Abhyankar 
9435f2c45f1SShri Abhyankar /*@
9442bf73ac6SHong Zhang   DMNetworkGetVertexRange - Get the bounds [start, end) for the local vertices
9455f2c45f1SShri Abhyankar 
9465f2c45f1SShri Abhyankar   Not Collective
9475f2c45f1SShri Abhyankar 
948f899ff85SJose E. Roman   Input Parameter:
9492bf73ac6SHong Zhang . dm - the DMNetwork object
9505f2c45f1SShri Abhyankar 
951fd292e60Sprj-   Output Parameters:
9522bf73ac6SHong Zhang + vStart - the first vertex point
9532bf73ac6SHong Zhang - vEnd - one beyond the last vertex point
9545f2c45f1SShri Abhyankar 
95597bb938eSShri Abhyankar   Level: beginner
9565f2c45f1SShri Abhyankar 
957db781477SPatrick Sanan .seealso: `DMNetworkGetEdgeRange()`
9585f2c45f1SShri Abhyankar @*/
9595f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
9605f2c45f1SShri Abhyankar {
9615f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network*)dm->data;
9625f2c45f1SShri Abhyankar 
9635f2c45f1SShri Abhyankar   PetscFunctionBegin;
964*daad07d3SAidan Hamilton   if (vStart) *vStart = network->cloneshared->vStart;
965*daad07d3SAidan Hamilton   if (vEnd) *vEnd = network->cloneshared->vEnd;
9665f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9675f2c45f1SShri Abhyankar }
9685f2c45f1SShri Abhyankar 
9695f2c45f1SShri Abhyankar /*@
9702bf73ac6SHong Zhang   DMNetworkGetEdgeRange - Get the bounds [start, end) for the local edges
9715f2c45f1SShri Abhyankar 
9725f2c45f1SShri Abhyankar   Not Collective
9735f2c45f1SShri Abhyankar 
974f899ff85SJose E. Roman   Input Parameter:
9752bf73ac6SHong Zhang . dm - the DMNetwork object
9765f2c45f1SShri Abhyankar 
977fd292e60Sprj-   Output Parameters:
9785f2c45f1SShri Abhyankar + eStart - The first edge point
9795f2c45f1SShri Abhyankar - eEnd - One beyond the last edge point
9805f2c45f1SShri Abhyankar 
98197bb938eSShri Abhyankar   Level: beginner
9825f2c45f1SShri Abhyankar 
983db781477SPatrick Sanan .seealso: `DMNetworkGetVertexRange()`
9845f2c45f1SShri Abhyankar @*/
9855f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
9865f2c45f1SShri Abhyankar {
9875f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network*)dm->data;
9885f2c45f1SShri Abhyankar 
9895f2c45f1SShri Abhyankar   PetscFunctionBegin;
9905c6496baSHong Zhang   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
991*daad07d3SAidan Hamilton   if (eStart) *eStart = network->cloneshared->eStart;
992*daad07d3SAidan Hamilton   if (eEnd)   *eEnd   = network->cloneshared->eEnd;
9935f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9945f2c45f1SShri Abhyankar }
9955f2c45f1SShri Abhyankar 
996*daad07d3SAidan Hamilton PetscErrorCode DMNetworkGetIndex(DM dm,PetscInt p,PetscInt *index)
9972bf73ac6SHong Zhang {
9982bf73ac6SHong Zhang   DM_Network *network = (DM_Network*)dm->data;
9995c6496baSHong Zhang 
10005c6496baSHong Zhang   PetscFunctionBegin;
10015c6496baSHong Zhang   if (network->header) {
10025c6496baSHong Zhang     *index = network->header[p].index;
10035c6496baSHong Zhang   } else {
10042bf73ac6SHong Zhang     PetscInt                 offsetp;
10052bf73ac6SHong Zhang     DMNetworkComponentHeader header;
10062bf73ac6SHong Zhang 
10079566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(network->DataSection,p,&offsetp));
10082bf73ac6SHong Zhang     header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
10092bf73ac6SHong Zhang     *index = header->index;
10105c6496baSHong Zhang   }
10112bf73ac6SHong Zhang   PetscFunctionReturn(0);
10122bf73ac6SHong Zhang }
10132bf73ac6SHong Zhang 
1014*daad07d3SAidan Hamilton PetscErrorCode DMNetworkGetSubnetID(DM dm,PetscInt p,PetscInt *subnetid)
1015*daad07d3SAidan Hamilton {
1016*daad07d3SAidan Hamilton   DM_Network *network = (DM_Network*)dm->data;
1017*daad07d3SAidan Hamilton 
1018*daad07d3SAidan Hamilton   PetscFunctionBegin;
1019*daad07d3SAidan Hamilton   if (network->header) {
1020*daad07d3SAidan Hamilton     *subnetid= network->header[p].subnetid;
1021*daad07d3SAidan Hamilton   } else {
1022*daad07d3SAidan Hamilton     PetscInt                 offsetp;
1023*daad07d3SAidan Hamilton     DMNetworkComponentHeader header;
1024*daad07d3SAidan Hamilton 
1025*daad07d3SAidan Hamilton     PetscCall(PetscSectionGetOffset(network->DataSection,p,&offsetp));
1026*daad07d3SAidan Hamilton     header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
1027*daad07d3SAidan Hamilton     *subnetid = header->subnetid;
1028*daad07d3SAidan Hamilton   }
1029*daad07d3SAidan Hamilton   PetscFunctionReturn(0);
1030*daad07d3SAidan Hamilton }
1031*daad07d3SAidan Hamilton 
10327b6afd5bSHong Zhang /*@
10332bf73ac6SHong Zhang   DMNetworkGetGlobalEdgeIndex - Get the global numbering for the edge on the network
10347b6afd5bSHong Zhang 
10357b6afd5bSHong Zhang   Not Collective
10367b6afd5bSHong Zhang 
10377b6afd5bSHong Zhang   Input Parameters:
10387b6afd5bSHong Zhang + dm - DMNetwork object
1039e85e6aecSHong Zhang - p - edge point
10407b6afd5bSHong Zhang 
1041fd292e60Sprj-   Output Parameters:
10422bf73ac6SHong Zhang . index - the global numbering for the edge
10437b6afd5bSHong Zhang 
10447b6afd5bSHong Zhang   Level: intermediate
10457b6afd5bSHong Zhang 
1046db781477SPatrick Sanan .seealso: `DMNetworkGetGlobalVertexIndex()`
10477b6afd5bSHong Zhang @*/
1048e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
10497b6afd5bSHong Zhang {
10507b6afd5bSHong Zhang   PetscFunctionBegin;
10519566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetIndex(dm,p,index));
10527b6afd5bSHong Zhang   PetscFunctionReturn(0);
10537b6afd5bSHong Zhang }
10547b6afd5bSHong Zhang 
10555f2c45f1SShri Abhyankar /*@
10562bf73ac6SHong Zhang   DMNetworkGetGlobalVertexIndex - Get the global numbering for the vertex on the network
1057e85e6aecSHong Zhang 
1058e85e6aecSHong Zhang   Not Collective
1059e85e6aecSHong Zhang 
1060e85e6aecSHong Zhang   Input Parameters:
1061e85e6aecSHong Zhang + dm - DMNetwork object
1062e85e6aecSHong Zhang - p  - vertex point
1063e85e6aecSHong Zhang 
1064fd292e60Sprj-   Output Parameters:
10652bf73ac6SHong Zhang . index - the global numbering for the vertex
1066e85e6aecSHong Zhang 
1067e85e6aecSHong Zhang   Level: intermediate
1068e85e6aecSHong Zhang 
1069db781477SPatrick Sanan .seealso: `DMNetworkGetGlobalEdgeIndex()`, `DMNetworkGetLocalVertexIndex()`
1070e85e6aecSHong Zhang @*/
1071e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
1072e85e6aecSHong Zhang {
1073e85e6aecSHong Zhang   PetscFunctionBegin;
10749566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetIndex(dm,p,index));
10759988b915SShri Abhyankar   PetscFunctionReturn(0);
10769988b915SShri Abhyankar }
10779988b915SShri Abhyankar 
10789988b915SShri Abhyankar /*@
10795f2c45f1SShri Abhyankar   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
10805f2c45f1SShri Abhyankar 
10815f2c45f1SShri Abhyankar   Not Collective
10825f2c45f1SShri Abhyankar 
10835f2c45f1SShri Abhyankar   Input Parameters:
10842bf73ac6SHong Zhang + dm - the DMNetwork object
1085a2b725a8SWilliam Gropp - p - vertex/edge point
10865f2c45f1SShri Abhyankar 
10875f2c45f1SShri Abhyankar   Output Parameters:
10885f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge
10895f2c45f1SShri Abhyankar 
109097bb938eSShri Abhyankar   Level: beginner
10915f2c45f1SShri Abhyankar 
1092db781477SPatrick Sanan .seealso: `DMNetworkRegisterComponent()`, `DMNetworkAddComponent()`
10935f2c45f1SShri Abhyankar @*/
10945f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
10955f2c45f1SShri Abhyankar {
10965f2c45f1SShri Abhyankar   PetscInt       offset;
10975f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
10985f2c45f1SShri Abhyankar 
10995f2c45f1SShri Abhyankar   PetscFunctionBegin;
11009566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offset));
11015f2c45f1SShri Abhyankar   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
11025f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
11035f2c45f1SShri Abhyankar }
11045f2c45f1SShri Abhyankar 
11055f2c45f1SShri Abhyankar /*@
11062bf73ac6SHong Zhang   DMNetworkGetLocalVecOffset - Get the offset for accessing the variables associated with a component at the given vertex/edge from the local vector
11075f2c45f1SShri Abhyankar 
11085f2c45f1SShri Abhyankar   Not Collective
11095f2c45f1SShri Abhyankar 
11105f2c45f1SShri Abhyankar   Input Parameters:
11112bf73ac6SHong Zhang + dm - the DMNetwork object
1112f2f58720SBarry Smith . p - the edge or vertex point
11132bf73ac6SHong Zhang - compnum - component number; use ALL_COMPONENTS if no specific component is requested
11149988b915SShri Abhyankar 
11159988b915SShri Abhyankar   Output Parameters:
11162bf73ac6SHong Zhang . offset - the local offset
11179988b915SShri Abhyankar 
1118f2f58720SBarry Smith   Notes:
1119f2f58720SBarry Smith     These offsets can be passed to MatSetValuesLocal() for matrices obtained with DMCreateMatrix().
1120f2f58720SBarry Smith 
1121f2f58720SBarry Smith     For vectors obtained with DMCreateLocalVector() the offsets can be used with VecSetValues().
1122f2f58720SBarry Smith 
1123f2f58720SBarry Smith     For vectors obtained with DMCreateLocalVector() and the array obtained with VecGetArray(vec,&array) you can access or set
1124f2f58720SBarry Smith     the vector values with array[offset].
1125f2f58720SBarry Smith 
1126f2f58720SBarry Smith     For vectors obtained with DMCreateGlobalVector() the offsets can be used with VecSetValuesLocal().
1127f2f58720SBarry Smith 
11289988b915SShri Abhyankar   Level: intermediate
11299988b915SShri Abhyankar 
1130db781477SPatrick Sanan .seealso: `DMGetLocalVector()`, `DMNetworkGetComponent()`, `DMNetworkGetGlobalVecOffset()`, `DMCreateGlobalVector()`, `VecGetArray()`, `VecSetValuesLocal()`, `MatSetValuesLocal()`
11319988b915SShri Abhyankar @*/
11322bf73ac6SHong Zhang PetscErrorCode DMNetworkGetLocalVecOffset(DM dm,PetscInt p,PetscInt compnum,PetscInt *offset)
11339988b915SShri Abhyankar {
11349988b915SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
11359988b915SShri Abhyankar   PetscInt       offsetp,offsetd;
11369988b915SShri Abhyankar   DMNetworkComponentHeader header;
11379988b915SShri Abhyankar 
11389988b915SShri Abhyankar   PetscFunctionBegin;
11399566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->plex->localSection,p,&offsetp));
11402bf73ac6SHong Zhang   if (compnum == ALL_COMPONENTS) {
11412bf73ac6SHong Zhang     *offset = offsetp;
11422bf73ac6SHong Zhang     PetscFunctionReturn(0);
11432bf73ac6SHong Zhang   }
11442bf73ac6SHong Zhang 
11459566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offsetd));
11469988b915SShri Abhyankar   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetd);
11479988b915SShri Abhyankar   *offset = offsetp + header->offsetvarrel[compnum];
11489988b915SShri Abhyankar   PetscFunctionReturn(0);
11499988b915SShri Abhyankar }
11509988b915SShri Abhyankar 
11519988b915SShri Abhyankar /*@
11522bf73ac6SHong Zhang   DMNetworkGetGlobalVecOffset - Get the global offset for accessing the variables associated with a component for the given vertex/edge from the global vector
11539988b915SShri Abhyankar 
11549988b915SShri Abhyankar   Not Collective
11559988b915SShri Abhyankar 
11569988b915SShri Abhyankar   Input Parameters:
11572bf73ac6SHong Zhang + dm - the DMNetwork object
1158f2f58720SBarry Smith . p - the edge or vertex point
11592bf73ac6SHong Zhang - compnum - component number; use ALL_COMPONENTS if no specific component is requested
11609988b915SShri Abhyankar 
11619988b915SShri Abhyankar   Output Parameters:
11629988b915SShri Abhyankar . offsetg - the global offset
11639988b915SShri Abhyankar 
1164f2f58720SBarry Smith   Notes:
1165f2f58720SBarry Smith     These offsets can be passed to MatSetValues() for matrices obtained with DMCreateMatrix().
1166f2f58720SBarry Smith 
1167f2f58720SBarry Smith     For vectors obtained with DMCreateGlobalVector() the offsets can be used with VecSetValues().
1168f2f58720SBarry Smith 
1169f2f58720SBarry Smith     For vectors obtained with DMCreateGlobalVector() and the array obtained with VecGetArray(vec,&array) you can access or set
1170f2f58720SBarry Smith     the vector values with array[offset - rstart] where restart is obtained with VecGetOwnershipRange(v,&rstart,NULL);
1171f2f58720SBarry Smith 
11729988b915SShri Abhyankar   Level: intermediate
11739988b915SShri Abhyankar 
1174db781477SPatrick Sanan .seealso: `DMNetworkGetLocalVecOffset()`, `DMGetGlobalVector()`, `DMNetworkGetComponent()`, `DMCreateGlobalVector()`, `VecGetArray()`, `VecSetValues()`, `MatSetValues()`
11759988b915SShri Abhyankar @*/
11762bf73ac6SHong Zhang PetscErrorCode DMNetworkGetGlobalVecOffset(DM dm,PetscInt p,PetscInt compnum,PetscInt *offsetg)
11779988b915SShri Abhyankar {
11789988b915SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
11799988b915SShri Abhyankar   PetscInt                 offsetp,offsetd;
11809988b915SShri Abhyankar   DMNetworkComponentHeader header;
11819988b915SShri Abhyankar 
11829988b915SShri Abhyankar   PetscFunctionBegin;
11839566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->plex->globalSection,p,&offsetp));
11842bf73ac6SHong Zhang   if (offsetp < 0) offsetp = -(offsetp + 1); /* Convert to actual global offset for ghost vertex */
11852bf73ac6SHong Zhang 
11862bf73ac6SHong Zhang   if (compnum == ALL_COMPONENTS) {
11872bf73ac6SHong Zhang     *offsetg = offsetp;
11882bf73ac6SHong Zhang     PetscFunctionReturn(0);
11892bf73ac6SHong Zhang   }
11909566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offsetd));
11919988b915SShri Abhyankar   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetd);
11929988b915SShri Abhyankar   *offsetg = offsetp + header->offsetvarrel[compnum];
11939988b915SShri Abhyankar   PetscFunctionReturn(0);
11949988b915SShri Abhyankar }
11959988b915SShri Abhyankar 
11969988b915SShri Abhyankar /*@
11972bf73ac6SHong Zhang   DMNetworkGetEdgeOffset - Get the offset for accessing the variables associated with the given edge from the local subvector
119824121865SAdrian Maldonado 
119924121865SAdrian Maldonado   Not Collective
120024121865SAdrian Maldonado 
120124121865SAdrian Maldonado   Input Parameters:
12022bf73ac6SHong Zhang + dm - the DMNetwork object
120324121865SAdrian Maldonado - p - the edge point
120424121865SAdrian Maldonado 
120524121865SAdrian Maldonado   Output Parameters:
120624121865SAdrian Maldonado . offset - the offset
120724121865SAdrian Maldonado 
120824121865SAdrian Maldonado   Level: intermediate
120924121865SAdrian Maldonado 
1210db781477SPatrick Sanan .seealso: `DMNetworkGetLocalVecOffset()`, `DMGetLocalVector()`
121124121865SAdrian Maldonado @*/
121224121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
121324121865SAdrian Maldonado {
121424121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
121524121865SAdrian Maldonado 
121624121865SAdrian Maldonado   PetscFunctionBegin;
12179566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->edge.DofSection,p,offset));
121824121865SAdrian Maldonado   PetscFunctionReturn(0);
121924121865SAdrian Maldonado }
122024121865SAdrian Maldonado 
122124121865SAdrian Maldonado /*@
12222bf73ac6SHong Zhang   DMNetworkGetVertexOffset - Get the offset for accessing the variables associated with the given vertex from the local subvector
122324121865SAdrian Maldonado 
122424121865SAdrian Maldonado   Not Collective
122524121865SAdrian Maldonado 
122624121865SAdrian Maldonado   Input Parameters:
12272bf73ac6SHong Zhang + dm - the DMNetwork object
122824121865SAdrian Maldonado - p - the vertex point
122924121865SAdrian Maldonado 
123024121865SAdrian Maldonado   Output Parameters:
123124121865SAdrian Maldonado . offset - the offset
123224121865SAdrian Maldonado 
123324121865SAdrian Maldonado   Level: intermediate
123424121865SAdrian Maldonado 
1235db781477SPatrick Sanan .seealso: `DMNetworkGetEdgeOffset()`, `DMGetLocalVector()`
123624121865SAdrian Maldonado @*/
123724121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
123824121865SAdrian Maldonado {
123924121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
124024121865SAdrian Maldonado 
124124121865SAdrian Maldonado   PetscFunctionBegin;
1242*daad07d3SAidan Hamilton   p -= network->cloneshared->vStart;
12439566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->vertex.DofSection,p,offset));
124424121865SAdrian Maldonado   PetscFunctionReturn(0);
124524121865SAdrian Maldonado }
12462bf73ac6SHong Zhang 
12475f2c45f1SShri Abhyankar /*@
12482bf73ac6SHong Zhang   DMNetworkAddComponent - Adds a network component and number of variables at the given point (vertex/edge)
12495f2c45f1SShri Abhyankar 
1250eac198afSGetnet   Collective on dm
12515f2c45f1SShri Abhyankar 
12525f2c45f1SShri Abhyankar   Input Parameters:
12534dc646bcSBarry Smith + dm - the DMNetwork
1254eac198afSGetnet . p - the vertex/edge point. These points are local indices provided by DMNetworkGetSubnetwork()
1255eac198afSGetnet . componentkey - component key returned while registering the component with DMNetworkRegisterComponent()
125624359064SBarry Smith . compvalue - pointer to the data structure for the component, or NULL if the component does not require data, this data is not copied so you cannot
125724359064SBarry Smith               free this space until after DMSetUp() is called.
1258eac198afSGetnet - nvar - number of variables for the component at the vertex/edge point, zero if the component does not introduce any degrees of freedom at the point
12595f2c45f1SShri Abhyankar 
1260eac198afSGetnet   Notes:
1261eac198afSGetnet     The owning rank and any other ranks that have this point as a ghost location must call this routine to add a component and number of variables in the same order at the given point.
1262eac198afSGetnet 
1263eac198afSGetnet     DMNetworkLayoutSetUp() must be called before this routine.
1264eac198afSGetnet 
1265eac198afSGetnet   Developer Notes:
1266eac198afSGetnet      The requirement that all the ranks with access to a vertex (as owner or as ghost) add all the components comes from a limitation of the underlying implementation based on DMPLEX.
126797bb938eSShri Abhyankar   Level: beginner
12685f2c45f1SShri Abhyankar 
1269db781477SPatrick Sanan .seealso: `DMNetworkGetComponent()`, `DMNetworkGetSubnetwork()`, `DMNetworkIsGhostVertex()`, `DMNetworkLayoutSetUp()`
12705f2c45f1SShri Abhyankar @*/
12712bf73ac6SHong Zhang PetscErrorCode DMNetworkAddComponent(DM dm,PetscInt p,PetscInt componentkey,void* compvalue,PetscInt nvar)
12725f2c45f1SShri Abhyankar {
12735f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
12742bf73ac6SHong Zhang   DMNetworkComponent       *component = &network->component[componentkey];
127554dfd506SHong Zhang   DMNetworkComponentHeader header;
127654dfd506SHong Zhang   DMNetworkComponentValue  cvalue;
127754dfd506SHong Zhang   PetscInt                 compnum;
127854dfd506SHong Zhang   PetscInt                 *compsize,*compkey,*compoffset,*compnvar,*compoffsetvarrel;
127954dfd506SHong Zhang   void*                    *compdata;
12805f2c45f1SShri Abhyankar 
12815f2c45f1SShri Abhyankar   PetscFunctionBegin;
12825c6496baSHong Zhang   PetscCheck(componentkey >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"componentkey %" PetscInt_FMT " cannot be negative. Input a component key returned while registering the component with DMNetworkRegisterComponent()",componentkey);
1283*daad07d3SAidan Hamilton   PetscCheck(network->componentsetup == PETSC_FALSE,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"The network has already finalized the components. No new components can be added.");
1284eac198afSGetnet   /* The owning rank and all ghost ranks add nvar */
12859566063dSJacob Faibussowitsch   PetscCall(PetscSectionAddDof(network->DofSection,p,nvar));
12862bf73ac6SHong Zhang 
1287eac198afSGetnet   /* The owning rank and all ghost ranks add a component, including compvalue=NULL */
128854dfd506SHong Zhang   header = &network->header[p];
128954dfd506SHong Zhang   cvalue = &network->cvalue[p];
129054dfd506SHong Zhang   if (header->ndata == header->maxcomps) {
1291f11a936eSBarry Smith     PetscInt additional_size;
1292f11a936eSBarry Smith 
129354dfd506SHong Zhang     /* Reached limit so resize header component arrays */
129454dfd506SHong Zhang     header->maxcomps += 2;
129554dfd506SHong Zhang 
129654dfd506SHong Zhang     /* Allocate arrays for component information and value */
12979566063dSJacob Faibussowitsch     PetscCall(PetscCalloc5(header->maxcomps,&compsize,header->maxcomps,&compkey,header->maxcomps,&compoffset,header->maxcomps,&compnvar,header->maxcomps,&compoffsetvarrel));
12989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(header->maxcomps,&compdata));
129954dfd506SHong Zhang 
130054dfd506SHong Zhang     /* Recalculate header size */
130154dfd506SHong Zhang     header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5*header->maxcomps*sizeof(PetscInt);
130254dfd506SHong Zhang 
130354dfd506SHong Zhang     header->hsize /= sizeof(DMNetworkComponentGenericDataType);
130454dfd506SHong Zhang 
130554dfd506SHong Zhang     /* Copy over component info */
13069566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(compsize,header->size,header->ndata*sizeof(PetscInt)));
13079566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(compkey,header->key,header->ndata*sizeof(PetscInt)));
13089566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(compoffset,header->offset,header->ndata*sizeof(PetscInt)));
13099566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(compnvar,header->nvar,header->ndata*sizeof(PetscInt)));
13109566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(compoffsetvarrel,header->offsetvarrel,header->ndata*sizeof(PetscInt)));
131154dfd506SHong Zhang 
131254dfd506SHong Zhang     /* Copy over component data pointers */
13139566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(compdata,cvalue->data,header->ndata*sizeof(void*)));
131454dfd506SHong Zhang 
131554dfd506SHong Zhang     /* Free old arrays */
13169566063dSJacob Faibussowitsch     PetscCall(PetscFree5(header->size,header->key,header->offset,header->nvar,header->offsetvarrel));
13179566063dSJacob Faibussowitsch     PetscCall(PetscFree(cvalue->data));
131854dfd506SHong Zhang 
131954dfd506SHong Zhang     /* Update pointers */
132054dfd506SHong Zhang     header->size         = compsize;
132154dfd506SHong Zhang     header->key          = compkey;
132254dfd506SHong Zhang     header->offset       = compoffset;
132354dfd506SHong Zhang     header->nvar         = compnvar;
132454dfd506SHong Zhang     header->offsetvarrel = compoffsetvarrel;
132554dfd506SHong Zhang 
132654dfd506SHong Zhang     cvalue->data = compdata;
132754dfd506SHong Zhang 
132854dfd506SHong Zhang     /* Update DataSection Dofs */
132954dfd506SHong Zhang     /* The dofs for datasection point p equals sizeof the header (i.e. header->hsize) + sizes of the components added at point p. With the resizing of the header, we need to update the dofs for point p. Hence, we add the extra size added for the header */
1330f11a936eSBarry Smith     additional_size = (5*(header->maxcomps - header->ndata)*sizeof(PetscInt))/sizeof(DMNetworkComponentGenericDataType);
13319566063dSJacob Faibussowitsch     PetscCall(PetscSectionAddDof(network->DataSection,p,additional_size));
133254dfd506SHong Zhang   }
133354dfd506SHong Zhang   header = &network->header[p];
133454dfd506SHong Zhang   cvalue = &network->cvalue[p];
133554dfd506SHong Zhang 
133654dfd506SHong Zhang   compnum = header->ndata;
13372bf73ac6SHong Zhang 
13382bf73ac6SHong Zhang   header->size[compnum] = component->size;
13399566063dSJacob Faibussowitsch   PetscCall(PetscSectionAddDof(network->DataSection,p,component->size));
13402bf73ac6SHong Zhang   header->key[compnum] = componentkey;
13412bf73ac6SHong Zhang   if (compnum != 0) header->offset[compnum] = header->offset[compnum-1] + header->size[compnum-1];
13422bf73ac6SHong Zhang   cvalue->data[compnum] = (void*)compvalue;
13432bf73ac6SHong Zhang 
13442bf73ac6SHong Zhang   /* variables */
13452bf73ac6SHong Zhang   header->nvar[compnum] += nvar;
13462bf73ac6SHong Zhang   if (compnum != 0) header->offsetvarrel[compnum] = header->offsetvarrel[compnum-1] + header->nvar[compnum-1];
13472bf73ac6SHong Zhang 
13482bf73ac6SHong Zhang   header->ndata++;
13495f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
13505f2c45f1SShri Abhyankar }
13515f2c45f1SShri Abhyankar 
135227f51fceSHong Zhang /*@
13532bf73ac6SHong Zhang   DMNetworkGetComponent - Gets the component key, the component data, and the number of variables at a given network point
135427f51fceSHong Zhang 
135527f51fceSHong Zhang   Not Collective
135627f51fceSHong Zhang 
135727f51fceSHong Zhang   Input Parameters:
13582bf73ac6SHong Zhang + dm - the DMNetwork object
13592bf73ac6SHong Zhang . p - vertex/edge point
136099c90e12SSatish Balay - compnum - component number; use ALL_COMPONENTS if sum up all the components
136127f51fceSHong Zhang 
136227f51fceSHong Zhang   Output Parameters:
13632bf73ac6SHong Zhang + compkey - the key obtained when registering the component (use NULL if not required)
13642bf73ac6SHong Zhang . component - the component data (use NULL if not required)
13652bf73ac6SHong Zhang - nvar  - number of variables (use NULL if not required)
136627f51fceSHong Zhang 
136797bb938eSShri Abhyankar   Level: beginner
136827f51fceSHong Zhang 
1369db781477SPatrick Sanan .seealso: `DMNetworkAddComponent()`, `DMNetworkGetNumComponents()`
137027f51fceSHong Zhang @*/
13712bf73ac6SHong Zhang PetscErrorCode DMNetworkGetComponent(DM dm,PetscInt p,PetscInt compnum,PetscInt *compkey,void **component,PetscInt *nvar)
137227f51fceSHong Zhang {
137327f51fceSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
13742bf73ac6SHong Zhang   PetscInt       offset = 0;
13752bf73ac6SHong Zhang   DMNetworkComponentHeader header;
137627f51fceSHong Zhang 
137727f51fceSHong Zhang   PetscFunctionBegin;
13782bf73ac6SHong Zhang   if (compnum == ALL_COMPONENTS) {
13799566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(network->DofSection,p,nvar));
138027f51fceSHong Zhang     PetscFunctionReturn(0);
138127f51fceSHong Zhang   }
138227f51fceSHong Zhang 
13839566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offset));
138442dc13f1SHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offset);
13855f2c45f1SShri Abhyankar 
13862bf73ac6SHong Zhang   if (compnum >= 0) {
13872bf73ac6SHong Zhang     if (compkey) *compkey = header->key[compnum];
13882bf73ac6SHong Zhang     if (component) {
138954dfd506SHong Zhang       offset += header->hsize+header->offset[compnum];
13902bf73ac6SHong Zhang       *component = network->componentdataarray+offset;
13912bf73ac6SHong Zhang     }
13922bf73ac6SHong Zhang   }
13935f2c45f1SShri Abhyankar 
13942bf73ac6SHong Zhang   if (nvar) *nvar = header->nvar[compnum];
139554dfd506SHong Zhang 
13965f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
13975f2c45f1SShri Abhyankar }
13985f2c45f1SShri Abhyankar 
13992bf73ac6SHong Zhang /*
14002bf73ac6SHong Zhang  Sets up the array that holds the data for all components and its associated section.
14012bf73ac6SHong Zhang  It copies the data for all components in a contiguous array called componentdataarray. The component data is stored pointwise with an additional header (metadata) stored for each point. The header has metadata information such as number of components at each point, number of variables for each component, offsets for the components data, etc.
14022bf73ac6SHong Zhang */
14035f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm)
14045f2c45f1SShri Abhyankar {
14055f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
1406f11a936eSBarry Smith   PetscInt                 arr_size,p,offset,offsetp,ncomp,i,*headerarr;
14075f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
14085f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue;
1409f11a936eSBarry Smith   DMNetworkComponentHeader headerinfo;
14105f2c45f1SShri Abhyankar   DMNetworkComponentGenericDataType *componentdataarray;
14115f2c45f1SShri Abhyankar 
14125f2c45f1SShri Abhyankar   PetscFunctionBegin;
14139566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(network->DataSection));
14149566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(network->DataSection,&arr_size));
1415f11a936eSBarry Smith   /* arr_size+1 fixes pipeline test of opensolaris-misc for src/dm/tests/ex10.c -- Do not know why */
14169566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(arr_size+1,&network->componentdataarray));
14175f2c45f1SShri Abhyankar   componentdataarray = network->componentdataarray;
1418*daad07d3SAidan Hamilton   for (p = network->cloneshared->pStart; p < network->cloneshared->pEnd; p++) {
14199566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(network->DataSection,p,&offsetp));
14205f2c45f1SShri Abhyankar     /* Copy header */
14215f2c45f1SShri Abhyankar     header = &network->header[p];
1422f11a936eSBarry Smith     headerinfo = (DMNetworkComponentHeader)(componentdataarray+offsetp);
14239566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(headerinfo,header,sizeof(struct _p_DMNetworkComponentHeader)));
1424f11a936eSBarry Smith     headerarr = (PetscInt*)(headerinfo+1);
14259566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(headerarr,header->size,header->maxcomps*sizeof(PetscInt)));
1426aeb78aa7SBarry Smith     headerinfo->size = headerarr;
142754dfd506SHong Zhang     headerarr += header->maxcomps;
14289566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(headerarr,header->key,header->maxcomps*sizeof(PetscInt)));
1429aeb78aa7SBarry Smith     headerinfo->key = headerarr;
143054dfd506SHong Zhang     headerarr += header->maxcomps;
14319566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(headerarr,header->offset,header->maxcomps*sizeof(PetscInt)));
1432aeb78aa7SBarry Smith     headerinfo->offset = headerarr;
143354dfd506SHong Zhang     headerarr += header->maxcomps;
14349566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(headerarr,header->nvar,header->maxcomps*sizeof(PetscInt)));
1435aeb78aa7SBarry Smith     headerinfo->nvar = headerarr;
143654dfd506SHong Zhang     headerarr += header->maxcomps;
14379566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(headerarr,header->offsetvarrel,header->maxcomps*sizeof(PetscInt)));
1438aeb78aa7SBarry Smith     headerinfo->offsetvarrel = headerarr;
143954dfd506SHong Zhang 
14405f2c45f1SShri Abhyankar     /* Copy data */
14415f2c45f1SShri Abhyankar     cvalue = &network->cvalue[p];
14425f2c45f1SShri Abhyankar     ncomp  = header->ndata;
14432bf73ac6SHong Zhang 
14445f2c45f1SShri Abhyankar     for (i = 0; i < ncomp; i++) {
144554dfd506SHong Zhang       offset = offsetp + header->hsize + header->offset[i];
14469566063dSJacob Faibussowitsch       PetscCall(PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType)));
14475f2c45f1SShri Abhyankar     }
14485f2c45f1SShri Abhyankar   }
1449aeb78aa7SBarry Smith 
1450*daad07d3SAidan Hamilton   for (i=network->cloneshared->pStart; i < network->cloneshared->pEnd; i++) {
1451aeb78aa7SBarry Smith     PetscCall(PetscFree5(network->header[i].size,network->header[i].key,network->header[i].offset,network->header[i].nvar,network->header[i].offsetvarrel));
1452aeb78aa7SBarry Smith     PetscCall(PetscFree(network->cvalue[i].data));
1453aeb78aa7SBarry Smith   }
1454aeb78aa7SBarry Smith   PetscCall(PetscFree2(network->header,network->cvalue));
14555f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
14565f2c45f1SShri Abhyankar }
14575f2c45f1SShri Abhyankar 
14585f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */
14592bf73ac6SHong Zhang static PetscErrorCode DMNetworkVariablesSetUp(DM dm)
14605f2c45f1SShri Abhyankar {
14615f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
14625f2c45f1SShri Abhyankar 
14635f2c45f1SShri Abhyankar   PetscFunctionBegin;
14649566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(network->DofSection));
14655f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
14665f2c45f1SShri Abhyankar }
14675f2c45f1SShri Abhyankar 
146824121865SAdrian Maldonado /* Get a subsection from a range of points */
14699dddd249SSatish Balay static PetscErrorCode DMNetworkGetSubSection_private(PetscSection main,PetscInt pstart,PetscInt pend,PetscSection *subsection)
147024121865SAdrian Maldonado {
147124121865SAdrian Maldonado   PetscInt       i, nvar;
147224121865SAdrian Maldonado 
147324121865SAdrian Maldonado   PetscFunctionBegin;
14749566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)main), subsection));
14759566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(*subsection, 0, pend - pstart));
147624121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
14779566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(main,i,&nvar));
14789566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(*subsection, i - pstart, nvar));
147924121865SAdrian Maldonado   }
148024121865SAdrian Maldonado 
14819566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(*subsection));
148224121865SAdrian Maldonado   PetscFunctionReturn(0);
148324121865SAdrian Maldonado }
148424121865SAdrian Maldonado 
148524121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */
14862bf73ac6SHong Zhang static PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
148724121865SAdrian Maldonado {
148824121865SAdrian Maldonado   PetscInt       i, *subpoints;
148924121865SAdrian Maldonado 
149024121865SAdrian Maldonado   PetscFunctionBegin;
149124121865SAdrian Maldonado   /* Create index sets to map from "points" to "subpoints" */
14929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pend - pstart, &subpoints));
149324121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
149424121865SAdrian Maldonado     subpoints[i - pstart] = i;
149524121865SAdrian Maldonado   }
14969566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map));
14979566063dSJacob Faibussowitsch   PetscCall(PetscFree(subpoints));
149824121865SAdrian Maldonado   PetscFunctionReturn(0);
149924121865SAdrian Maldonado }
150024121865SAdrian Maldonado 
150124121865SAdrian Maldonado /*@
15022bf73ac6SHong Zhang   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute
150324121865SAdrian Maldonado 
15042bf73ac6SHong Zhang   Collective on dm
150524121865SAdrian Maldonado 
150624121865SAdrian Maldonado   Input Parameters:
15072bf73ac6SHong Zhang . dm - the DMNetworkObject
150824121865SAdrian Maldonado 
150924121865SAdrian Maldonado   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
151024121865SAdrian Maldonado 
151124121865SAdrian Maldonado   points = [0 1 2 3 4 5 6]
151224121865SAdrian Maldonado 
1513caf410d2SHong Zhang   where edges = [0,1,2,3] and vertices = [4,5,6]. The new orderings will be specific to the subset (i.e vertices = [0,1,2] <- [4,5,6]).
151424121865SAdrian Maldonado 
151524121865SAdrian Maldonado   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
151624121865SAdrian Maldonado 
151724121865SAdrian Maldonado   Level: intermediate
151824121865SAdrian Maldonado 
151924121865SAdrian Maldonado @*/
152024121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
152124121865SAdrian Maldonado {
152224121865SAdrian Maldonado   MPI_Comm       comm;
1523f11a936eSBarry Smith   PetscMPIInt    size;
152424121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
152524121865SAdrian Maldonado 
1526eab1376dSHong Zhang   PetscFunctionBegin;
15279566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
15289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
152924121865SAdrian Maldonado 
153024121865SAdrian Maldonado   /* Create maps for vertices and edges */
1531*daad07d3SAidan Hamilton   PetscCall(DMNetworkSetSubMap_private(network->cloneshared->vStart,network->cloneshared->vEnd,&network->vertex.mapping));
1532*daad07d3SAidan Hamilton   PetscCall(DMNetworkSetSubMap_private(network->cloneshared->eStart,network->cloneshared->eEnd,&network->edge.mapping));
153324121865SAdrian Maldonado 
153424121865SAdrian Maldonado   /* Create local sub-sections */
1535*daad07d3SAidan Hamilton   PetscCall(DMNetworkGetSubSection_private(network->DofSection,network->cloneshared->vStart,network->cloneshared->vEnd,&network->vertex.DofSection));
1536*daad07d3SAidan Hamilton   PetscCall(DMNetworkGetSubSection_private(network->DofSection,network->cloneshared->eStart,network->cloneshared->eEnd,&network->edge.DofSection));
153724121865SAdrian Maldonado 
15389852e123SBarry Smith   if (size > 1) {
15399566063dSJacob Faibussowitsch     PetscCall(PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf));
154022bbedd7SHong Zhang 
15419566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection));
15429566063dSJacob Faibussowitsch     PetscCall(PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf));
15439566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection));
154424121865SAdrian Maldonado   } else {
154524121865SAdrian Maldonado     /* create structures for vertex */
15469566063dSJacob Faibussowitsch     PetscCall(PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection));
154724121865SAdrian Maldonado     /* create structures for edge */
15489566063dSJacob Faibussowitsch     PetscCall(PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection));
154924121865SAdrian Maldonado   }
155024121865SAdrian Maldonado 
155124121865SAdrian Maldonado   /* Add viewers */
15529566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section"));
15539566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section"));
15549566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view"));
15559566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view"));
155624121865SAdrian Maldonado   PetscFunctionReturn(0);
155724121865SAdrian Maldonado }
15587b6afd5bSHong Zhang 
1559f11a936eSBarry Smith /*
15605c6496baSHong Zhang    Setup a lookup btable for the input v's owning subnetworks
15615c6496baSHong Zhang    - add all owing subnetworks that connect to this v to the btable
1562f11a936eSBarry Smith      vertex_subnetid = supportingedge_subnetid
1563f11a936eSBarry Smith */
15649fbee547SJacob Faibussowitsch static inline PetscErrorCode SetSubnetIdLookupBT(DM dm,PetscInt v,PetscInt Nsubnet,PetscBT btable)
1565f11a936eSBarry Smith {
1566f11a936eSBarry Smith   PetscInt       e,nedges,offset;
1567f11a936eSBarry Smith   const PetscInt *edges;
1568f11a936eSBarry Smith   DM_Network     *newDMnetwork = (DM_Network*)dm->data;
1569f11a936eSBarry Smith   DMNetworkComponentHeader header;
1570f11a936eSBarry Smith 
1571f11a936eSBarry Smith   PetscFunctionBegin;
15729566063dSJacob Faibussowitsch   PetscCall(PetscBTMemzero(Nsubnet,btable));
15739566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSupportingEdges(dm,v,&nedges,&edges));
1574f11a936eSBarry Smith   for (e=0; e<nedges; e++) {
15759566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,edges[e],&offset));
1576f11a936eSBarry Smith     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
15779566063dSJacob Faibussowitsch     PetscCall(PetscBTSet(btable,header->subnetid));
1578f11a936eSBarry Smith   }
1579f11a936eSBarry Smith   PetscFunctionReturn(0);
1580f11a936eSBarry Smith }
1581f11a936eSBarry Smith 
15825f2c45f1SShri Abhyankar /*@
15832bf73ac6SHong Zhang   DMNetworkDistribute - Distributes the network and moves associated component data
15845f2c45f1SShri Abhyankar 
15855f2c45f1SShri Abhyankar   Collective
15865f2c45f1SShri Abhyankar 
1587d8d19677SJose E. Roman   Input Parameters:
1588d3464fd4SAdrian Maldonado + DM - the DMNetwork object
15892bf73ac6SHong Zhang - overlap - the overlap of partitions, 0 is the default
15905f2c45f1SShri Abhyankar 
15915b3f975aSHong Zhang   Options Database Key:
15925b3f975aSHong Zhang + -dmnetwork_view - Calls DMView() at the conclusion of DMSetUp()
15935b3f975aSHong Zhang - -dmnetwork_view_distributed - Calls DMView() at the conclusion of DMNetworkDistribute()
15945b3f975aSHong Zhang 
15955f2c45f1SShri Abhyankar   Notes:
15968b171c8eSHong Zhang   Distributes the network with <overlap>-overlapping partitioning of the edges.
15975f2c45f1SShri Abhyankar 
15985f2c45f1SShri Abhyankar   Level: intermediate
15995f2c45f1SShri Abhyankar 
1600db781477SPatrick Sanan .seealso: `DMNetworkCreate()`
16015f2c45f1SShri Abhyankar @*/
1602d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
16035f2c45f1SShri Abhyankar {
1604d3464fd4SAdrian Maldonado   MPI_Comm       comm;
1605d3464fd4SAdrian Maldonado   PetscMPIInt    size;
1606d3464fd4SAdrian Maldonado   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
1607d3464fd4SAdrian Maldonado   DM_Network     *newDMnetwork;
1608caf410d2SHong Zhang   PetscSF        pointsf=NULL;
16095f2c45f1SShri Abhyankar   DM             newDM;
1610f11a936eSBarry Smith   PetscInt       j,e,v,offset,*subnetvtx,*subnetedge,Nsubnet,gidx,svtx_idx,nv;
16115c6496baSHong Zhang   PetscInt       net,*sv;
1612f11a936eSBarry Smith   PetscBT        btable;
161351ac5effSHong Zhang   PetscPartitioner         part;
1614b9c6e19dSShri Abhyankar   DMNetworkComponentHeader header;
16155f2c45f1SShri Abhyankar 
16165f2c45f1SShri Abhyankar   PetscFunctionBegin;
16175c6496baSHong Zhang   PetscValidPointer(dm,1);
16185c6496baSHong Zhang   PetscValidHeaderSpecific(*dm,DM_CLASSID,1);
16199566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)*dm,&comm));
16209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1621d3464fd4SAdrian Maldonado   if (size == 1) PetscFunctionReturn(0);
1622d3464fd4SAdrian Maldonado 
16235c6496baSHong Zhang   PetscCheck(!overlap,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"overlap %" PetscInt_FMT " != 0 is not supported yet",overlap);
16245b3f975aSHong Zhang 
16252bf73ac6SHong Zhang   /* This routine moves the component data to the appropriate processors. It makes use of the DataSection and the componentdataarray to move the component data to appropriate processors and returns a new DataSection and new componentdataarray. */
1626e3b1a053SGetnet Betrie   PetscCall(PetscLogEventBegin(DMNetwork_Distribute,dm,0,0,0));
16279566063dSJacob Faibussowitsch   PetscCall(DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM));
16285f2c45f1SShri Abhyankar   newDMnetwork = (DM_Network*)newDM->data;
162954dfd506SHong Zhang   newDMnetwork->max_comps_registered = oldDMnetwork->max_comps_registered;
16309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(newDMnetwork->max_comps_registered,&newDMnetwork->component));
163151ac5effSHong Zhang 
163251ac5effSHong Zhang   /* Enable runtime options for petscpartitioner */
16339566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPartitioner(oldDMnetwork->plex,&part));
16349566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetFromOptions(part));
163551ac5effSHong Zhang 
16362bf73ac6SHong Zhang   /* Distribute plex dm */
16379566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex));
163851ac5effSHong Zhang 
16395f2c45f1SShri Abhyankar   /* Distribute dof section */
16409566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm,&newDMnetwork->DofSection));
16419566063dSJacob Faibussowitsch   PetscCall(PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection));
164251ac5effSHong Zhang 
16435f2c45f1SShri Abhyankar   /* Distribute data and associated section */
16449566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm,&newDMnetwork->DataSection));
16459566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray));
164624121865SAdrian Maldonado 
1647*daad07d3SAidan Hamilton   PetscCall(PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->cloneshared->pStart,&newDMnetwork->cloneshared->pEnd));
1648*daad07d3SAidan Hamilton   PetscCall(DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->cloneshared->eStart,&newDMnetwork->cloneshared->eEnd));
1649*daad07d3SAidan Hamilton   PetscCall(DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->cloneshared->vStart,&newDMnetwork->cloneshared->vEnd));
1650*daad07d3SAidan Hamilton   newDMnetwork->cloneshared->nEdges    = newDMnetwork->cloneshared->eEnd - newDMnetwork->cloneshared->eStart;
1651*daad07d3SAidan Hamilton   newDMnetwork->cloneshared->nVertices = newDMnetwork->cloneshared->vEnd - newDMnetwork->cloneshared->vStart;
1652*daad07d3SAidan Hamilton   newDMnetwork->cloneshared->NVertices = oldDMnetwork->cloneshared->NVertices;
1653*daad07d3SAidan Hamilton   newDMnetwork->cloneshared->NEdges    = oldDMnetwork->cloneshared->NEdges;
1654*daad07d3SAidan Hamilton   newDMnetwork->cloneshared->svtable   = oldDMnetwork->cloneshared->svtable; /* global table! */
1655*daad07d3SAidan Hamilton   oldDMnetwork->cloneshared->svtable   = NULL;
165624121865SAdrian Maldonado 
16571bb6d2a8SBarry Smith   /* Set Dof section as the section for dm */
16589566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(newDMnetwork->plex,newDMnetwork->DofSection));
16599566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection));
16605f2c45f1SShri Abhyankar 
1661b9c6e19dSShri Abhyankar   /* Setup subnetwork info in the newDM */
1662*daad07d3SAidan Hamilton   newDMnetwork->cloneshared->Nsubnet = oldDMnetwork->cloneshared->Nsubnet;
1663*daad07d3SAidan Hamilton   newDMnetwork->cloneshared->Nsvtx   = oldDMnetwork->cloneshared->Nsvtx;
1664*daad07d3SAidan Hamilton   oldDMnetwork->cloneshared->Nsvtx   = 0;
1665*daad07d3SAidan Hamilton   newDMnetwork->cloneshared->svtx    = oldDMnetwork->cloneshared->svtx; /* global vertices! */
1666*daad07d3SAidan Hamilton   oldDMnetwork->cloneshared->svtx    = NULL;
1667*daad07d3SAidan Hamilton   PetscCall(PetscCalloc1(newDMnetwork->cloneshared->Nsubnet,&newDMnetwork->cloneshared->subnet));
16682bf73ac6SHong Zhang 
16692bf73ac6SHong Zhang   /* Copy over the global number of vertices and edges in each subnetwork.
16702bf73ac6SHong Zhang      Note: these are calculated in DMNetworkLayoutSetUp()
1671b9c6e19dSShri Abhyankar   */
1672*daad07d3SAidan Hamilton   Nsubnet = newDMnetwork->cloneshared->Nsubnet;
1673f11a936eSBarry Smith   for (j = 0; j < Nsubnet; j++) {
1674*daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[j].Nvtx  = oldDMnetwork->cloneshared->subnet[j].Nvtx;
1675*daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[j].Nedge = oldDMnetwork->cloneshared->subnet[j].Nedge;
1676b9c6e19dSShri Abhyankar   }
1677b9c6e19dSShri Abhyankar 
1678f11a936eSBarry Smith   /* Count local nedges for subnetworks */
1679*daad07d3SAidan Hamilton   for (e = newDMnetwork->cloneshared->eStart; e < newDMnetwork->cloneshared->eEnd; e++) {
16809566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset));
168154dfd506SHong Zhang     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
16825c6496baSHong Zhang 
168354dfd506SHong Zhang     /* Update pointers */
168454dfd506SHong Zhang     header->size          = (PetscInt*)(header + 1);
168554dfd506SHong Zhang     header->key           = header->size   + header->maxcomps;
168654dfd506SHong Zhang     header->offset        = header->key    + header->maxcomps;
168754dfd506SHong Zhang     header->nvar          = header->offset + header->maxcomps;
168854dfd506SHong Zhang     header->offsetvarrel  = header->nvar   + header->maxcomps;
168954dfd506SHong Zhang 
1690*daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[header->subnetid].nedge++;
1691b9c6e19dSShri Abhyankar   }
1692b9c6e19dSShri Abhyankar 
16935c6496baSHong Zhang   /* Setup a btable to keep track subnetworks owned by this process at a shared vertex */
1694*daad07d3SAidan Hamilton   if (newDMnetwork->cloneshared->Nsvtx) {
16959566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(Nsubnet,&btable));
16965c6496baSHong Zhang   }
16975c6496baSHong Zhang 
16985c6496baSHong Zhang   /* Count local nvtx for subnetworks */
1699*daad07d3SAidan Hamilton   for (v = newDMnetwork->cloneshared->vStart; v < newDMnetwork->cloneshared->vEnd; v++) {
17009566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset));
17015f80ce2aSJacob Faibussowitsch     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
17025c6496baSHong Zhang 
170354dfd506SHong Zhang     /* Update pointers */
170454dfd506SHong Zhang     header->size          = (PetscInt*)(header + 1);
170554dfd506SHong Zhang     header->key           = header->size   + header->maxcomps;
170654dfd506SHong Zhang     header->offset        = header->key    + header->maxcomps;
170754dfd506SHong Zhang     header->nvar          = header->offset + header->maxcomps;
170854dfd506SHong Zhang     header->offsetvarrel  = header->nvar   + header->maxcomps;
170954dfd506SHong Zhang 
17102bf73ac6SHong Zhang     /* shared vertices: use gidx=header->index to check if v is a shared vertex */
17112bf73ac6SHong Zhang     gidx = header->index;
1712*daad07d3SAidan Hamilton     PetscCall(PetscTableFind(newDMnetwork->cloneshared->svtable,gidx+1,&svtx_idx));
17132bf73ac6SHong Zhang     svtx_idx--;
17142bf73ac6SHong Zhang 
17152bf73ac6SHong Zhang     if (svtx_idx < 0) { /* not a shared vertex */
1716*daad07d3SAidan Hamilton       newDMnetwork->cloneshared->subnet[header->subnetid].nvtx++;
17172bf73ac6SHong Zhang     } else { /* a shared vertex belongs to more than one subnetworks, it is being counted by multiple subnets */
17185c6496baSHong Zhang       /* Setup a lookup btable for this v's owning subnetworks */
17199566063dSJacob Faibussowitsch       PetscCall(SetSubnetIdLookupBT(newDM,v,Nsubnet,btable));
1720f11a936eSBarry Smith 
1721*daad07d3SAidan Hamilton       for (j=0; j<newDMnetwork->cloneshared->svtx[svtx_idx].n; j++) {
1722*daad07d3SAidan Hamilton         sv  = newDMnetwork->cloneshared->svtx[svtx_idx].sv + 2*j;
17235c6496baSHong Zhang         net = sv[0];
1724*daad07d3SAidan Hamilton         if (PetscBTLookup(btable,net)) newDMnetwork->cloneshared->subnet[net].nvtx++; /* sv is on net owned by this proces */
17252bf73ac6SHong Zhang       }
17262bf73ac6SHong Zhang     }
1727b9c6e19dSShri Abhyankar   }
1728b9c6e19dSShri Abhyankar 
17292bf73ac6SHong Zhang   /* Get total local nvtx for subnetworks */
17302bf73ac6SHong Zhang   nv = 0;
1731*daad07d3SAidan Hamilton   for (j=0; j<Nsubnet; j++) nv += newDMnetwork->cloneshared->subnet[j].nvtx;
1732*daad07d3SAidan Hamilton   nv += newDMnetwork->cloneshared->Nsvtx;
1733caf410d2SHong Zhang 
17342bf73ac6SHong Zhang   /* Now create the vertices and edge arrays for the subnetworks */
1735*daad07d3SAidan Hamilton   PetscCall(PetscCalloc2(newDMnetwork->cloneshared->nEdges,&subnetedge,nv,&subnetvtx)); /* Maps local vertex to local subnetwork's vertex */
1736*daad07d3SAidan Hamilton   newDMnetwork->cloneshared->subnetedge = subnetedge;
1737*daad07d3SAidan Hamilton   newDMnetwork->cloneshared->subnetvtx  = subnetvtx;
1738*daad07d3SAidan Hamilton   for (j=0; j < newDMnetwork->cloneshared->Nsubnet; j++) {
1739*daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[j].edges = subnetedge;
1740*daad07d3SAidan Hamilton     subnetedge                   += newDMnetwork->cloneshared->subnet[j].nedge;
17412bf73ac6SHong Zhang 
1742*daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[j].vertices = subnetvtx;
1743*daad07d3SAidan Hamilton     subnetvtx                       += newDMnetwork->cloneshared->subnet[j].nvtx;
1744caf410d2SHong Zhang 
17452bf73ac6SHong Zhang     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. These get updated when the vertices and edges are added. */
1746*daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[j].nvtx = newDMnetwork->cloneshared->subnet[j].nedge = 0;
1747b9c6e19dSShri Abhyankar   }
1748*daad07d3SAidan Hamilton   newDMnetwork->cloneshared->svertices = subnetvtx;
1749b9c6e19dSShri Abhyankar 
17502bf73ac6SHong Zhang   /* Set the edges and vertices in each subnetwork */
1751*daad07d3SAidan Hamilton   for (e = newDMnetwork->cloneshared->eStart; e < newDMnetwork->cloneshared->eEnd; e++) {
17529566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset));
17535f80ce2aSJacob Faibussowitsch     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1754*daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[header->subnetid].edges[newDMnetwork->cloneshared->subnet[header->subnetid].nedge++] = e;
1755b9c6e19dSShri Abhyankar   }
1756b9c6e19dSShri Abhyankar 
17572bf73ac6SHong Zhang   nv = 0;
1758*daad07d3SAidan Hamilton   for (v = newDMnetwork->cloneshared->vStart; v < newDMnetwork->cloneshared->vEnd; v++) {
17599566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset));
17605f80ce2aSJacob Faibussowitsch     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
17612bf73ac6SHong Zhang 
17622bf73ac6SHong Zhang     /* coupling vertices: use gidx = header->index to check if v is a coupling vertex */
1763*daad07d3SAidan Hamilton     PetscCall(PetscTableFind(newDMnetwork->cloneshared->svtable,header->index+1,&svtx_idx));
17642bf73ac6SHong Zhang     svtx_idx--;
17652bf73ac6SHong Zhang     if (svtx_idx < 0) {
1766*daad07d3SAidan Hamilton       newDMnetwork->cloneshared->subnet[header->subnetid].vertices[newDMnetwork->cloneshared->subnet[header->subnetid].nvtx++] = v;
17672bf73ac6SHong Zhang     } else { /* a shared vertex */
1768*daad07d3SAidan Hamilton       newDMnetwork->cloneshared->svertices[nv++] = v;
17692bf73ac6SHong Zhang 
17705c6496baSHong Zhang       /* Setup a lookup btable for this v's owning subnetworks */
17719566063dSJacob Faibussowitsch       PetscCall(SetSubnetIdLookupBT(newDM,v,Nsubnet,btable));
1772f11a936eSBarry Smith 
1773*daad07d3SAidan Hamilton       for (j=0; j<newDMnetwork->cloneshared->svtx[svtx_idx].n; j++) {
1774*daad07d3SAidan Hamilton         sv  = newDMnetwork->cloneshared->svtx[svtx_idx].sv + 2*j;
17755c6496baSHong Zhang         net = sv[0];
17765c6496baSHong Zhang         if (PetscBTLookup(btable,net))
1777*daad07d3SAidan Hamilton           newDMnetwork->cloneshared->subnet[net].vertices[newDMnetwork->cloneshared->subnet[net].nvtx++] = v;
1778b9c6e19dSShri Abhyankar       }
17792bf73ac6SHong Zhang     }
17802bf73ac6SHong Zhang   }
1781*daad07d3SAidan Hamilton   newDMnetwork->cloneshared->nsvtx = nv;   /* num of local shared vertices */
1782b9c6e19dSShri Abhyankar 
1783caf410d2SHong Zhang   newDM->setupcalled = (*dm)->setupcalled;
1784*daad07d3SAidan Hamilton   newDMnetwork->cloneshared->distributecalled = PETSC_TRUE;
1785caf410d2SHong Zhang 
17862bf73ac6SHong Zhang   /* Free spaces */
17879566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&pointsf));
17889566063dSJacob Faibussowitsch   PetscCall(DMDestroy(dm));
1789*daad07d3SAidan Hamilton   if (newDMnetwork->cloneshared->Nsvtx) {
17909566063dSJacob Faibussowitsch     PetscCall(PetscBTDestroy(&btable));
17915c6496baSHong Zhang   }
17922bf73ac6SHong Zhang 
17935b3f975aSHong Zhang   /* View distributed dmnetwork */
17949566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(newDM,NULL,"-dmnetwork_view_distributed"));
17955b3f975aSHong Zhang 
1796d3464fd4SAdrian Maldonado   *dm  = newDM;
1797e3b1a053SGetnet Betrie   PetscCall(PetscLogEventEnd(DMNetwork_Distribute,dm,0,0,0));
17985f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
17995f2c45f1SShri Abhyankar }
18005f2c45f1SShri Abhyankar 
180124121865SAdrian Maldonado /*@C
18022bf73ac6SHong Zhang   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering
18032bf73ac6SHong Zhang 
18042bf73ac6SHong Zhang  Collective
180524121865SAdrian Maldonado 
180624121865SAdrian Maldonado   Input Parameters:
18079dddd249SSatish Balay + mainSF - the original SF structure
180824121865SAdrian Maldonado - map - a ISLocalToGlobal mapping that contains the subset of points
180924121865SAdrian Maldonado 
181097bb3fdcSJose E. Roman   Output Parameter:
18119dddd249SSatish Balay . subSF - a subset of the mainSF for the desired subset.
1812ee300463SSatish Balay 
1813ee300463SSatish Balay   Level: intermediate
18147d928bffSSatish Balay @*/
18159dddd249SSatish Balay PetscErrorCode PetscSFGetSubSF(PetscSF mainsf,ISLocalToGlobalMapping map,PetscSF *subSF)
18162bf73ac6SHong Zhang {
181724121865SAdrian Maldonado   PetscInt              nroots, nleaves, *ilocal_sub;
181824121865SAdrian Maldonado   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
181924121865SAdrian Maldonado   PetscInt              *local_points, *remote_points;
182024121865SAdrian Maldonado   PetscSFNode           *iremote_sub;
182124121865SAdrian Maldonado   const PetscInt        *ilocal;
182224121865SAdrian Maldonado   const PetscSFNode     *iremote;
182324121865SAdrian Maldonado 
182424121865SAdrian Maldonado   PetscFunctionBegin;
18259566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(mainsf,&nroots,&nleaves,&ilocal,&iremote));
182624121865SAdrian Maldonado 
182724121865SAdrian Maldonado   /* Look for leaves that pertain to the subset of points. Get the local ordering */
18289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves,&ilocal_map));
18299566063dSJacob Faibussowitsch   PetscCall(ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map));
183024121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
183124121865SAdrian Maldonado     if (ilocal_map[i] != -1) nleaves_sub += 1;
183224121865SAdrian Maldonado   }
183324121865SAdrian Maldonado   /* Re-number ilocal with subset numbering. Need information from roots */
18349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nroots,&local_points,nroots,&remote_points));
183524121865SAdrian Maldonado   for (i = 0; i < nroots; i++) local_points[i] = i;
18369566063dSJacob Faibussowitsch   PetscCall(ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points));
18379566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(mainsf, MPIU_INT, local_points, remote_points,MPI_REPLACE));
18389566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(mainsf, MPIU_INT, local_points, remote_points,MPI_REPLACE));
183924121865SAdrian Maldonado   /* Fill up graph using local (that is, local to the subset) numbering. */
18409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves_sub,&ilocal_sub));
18419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves_sub,&iremote_sub));
184224121865SAdrian Maldonado   nleaves_sub = 0;
184324121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
184424121865SAdrian Maldonado     if (ilocal_map[i] != -1) {
184524121865SAdrian Maldonado       ilocal_sub[nleaves_sub] = ilocal_map[i];
18464b70a8deSAdrian Maldonado       iremote_sub[nleaves_sub].rank = iremote[i].rank;
184724121865SAdrian Maldonado       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
184824121865SAdrian Maldonado       nleaves_sub += 1;
184924121865SAdrian Maldonado     }
185024121865SAdrian Maldonado   }
18519566063dSJacob Faibussowitsch   PetscCall(PetscFree2(local_points,remote_points));
18529566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetSize(map,&nroots_sub));
185324121865SAdrian Maldonado 
185424121865SAdrian Maldonado   /* Create new subSF */
18559566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PETSC_COMM_WORLD,subSF));
18569566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(*subSF));
18579566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES));
18589566063dSJacob Faibussowitsch   PetscCall(PetscFree(ilocal_map));
18599566063dSJacob Faibussowitsch   PetscCall(PetscFree(iremote_sub));
186024121865SAdrian Maldonado   PetscFunctionReturn(0);
186124121865SAdrian Maldonado }
186224121865SAdrian Maldonado 
18635f2c45f1SShri Abhyankar /*@C
18645f2c45f1SShri Abhyankar   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
18655f2c45f1SShri Abhyankar 
18665f2c45f1SShri Abhyankar   Not Collective
18675f2c45f1SShri Abhyankar 
18685f2c45f1SShri Abhyankar   Input Parameters:
18692bf73ac6SHong Zhang + dm - the DMNetwork object
18705f2c45f1SShri Abhyankar - p  - the vertex point
18715f2c45f1SShri Abhyankar 
1872fd292e60Sprj-   Output Parameters:
18735f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point
18742bf73ac6SHong Zhang - edges  - list of edge points
18755f2c45f1SShri Abhyankar 
187697bb938eSShri Abhyankar   Level: beginner
18775f2c45f1SShri Abhyankar 
18785f2c45f1SShri Abhyankar   Fortran Notes:
18795f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
18805f2c45f1SShri Abhyankar   include petsc.h90 in your code.
18815f2c45f1SShri Abhyankar 
1882db781477SPatrick Sanan .seealso: `DMNetworkCreate()`, `DMNetworkGetConnectedVertices()`
18835f2c45f1SShri Abhyankar @*/
18845f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
18855f2c45f1SShri Abhyankar {
18865f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
18875f2c45f1SShri Abhyankar 
18885f2c45f1SShri Abhyankar   PetscFunctionBegin;
18899566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSupportSize(network->plex,vertex,nedges));
18909566063dSJacob Faibussowitsch   if (edges) PetscCall(DMPlexGetSupport(network->plex,vertex,edges));
18915f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
18925f2c45f1SShri Abhyankar }
18935f2c45f1SShri Abhyankar 
18945f2c45f1SShri Abhyankar /*@C
1895d842c372SHong Zhang   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
18965f2c45f1SShri Abhyankar 
18975f2c45f1SShri Abhyankar   Not Collective
18985f2c45f1SShri Abhyankar 
18995f2c45f1SShri Abhyankar   Input Parameters:
19002bf73ac6SHong Zhang + dm - the DMNetwork object
19015f2c45f1SShri Abhyankar - p - the edge point
19025f2c45f1SShri Abhyankar 
1903fd292e60Sprj-   Output Parameters:
19045f2c45f1SShri Abhyankar . vertices - vertices connected to this edge
19055f2c45f1SShri Abhyankar 
190697bb938eSShri Abhyankar   Level: beginner
19075f2c45f1SShri Abhyankar 
19085f2c45f1SShri Abhyankar   Fortran Notes:
19095f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
19105f2c45f1SShri Abhyankar   include petsc.h90 in your code.
19115f2c45f1SShri Abhyankar 
1912db781477SPatrick Sanan .seealso: `DMNetworkCreate()`, `DMNetworkGetSupportingEdges()`
19135f2c45f1SShri Abhyankar @*/
1914d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
19155f2c45f1SShri Abhyankar {
19165f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
19175f2c45f1SShri Abhyankar 
19185f2c45f1SShri Abhyankar   PetscFunctionBegin;
19199566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(network->plex,edge,vertices));
19205f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19215f2c45f1SShri Abhyankar }
19225f2c45f1SShri Abhyankar 
19235f2c45f1SShri Abhyankar /*@
19242bf73ac6SHong Zhang   DMNetworkIsSharedVertex - Returns TRUE if the vertex is shared by subnetworks
19252bf73ac6SHong Zhang 
19262bf73ac6SHong Zhang   Not Collective
19272bf73ac6SHong Zhang 
19282bf73ac6SHong Zhang   Input Parameters:
19292bf73ac6SHong Zhang + dm - the DMNetwork object
19302bf73ac6SHong Zhang - p - the vertex point
19312bf73ac6SHong Zhang 
19322bf73ac6SHong Zhang   Output Parameter:
19332bf73ac6SHong Zhang . flag - TRUE if the vertex is shared by subnetworks
19342bf73ac6SHong Zhang 
19352bf73ac6SHong Zhang   Level: beginner
19362bf73ac6SHong Zhang 
1937db781477SPatrick Sanan .seealso: `DMNetworkAddSharedVertices()`, `DMNetworkIsGhostVertex()`
19382bf73ac6SHong Zhang @*/
19392bf73ac6SHong Zhang PetscErrorCode DMNetworkIsSharedVertex(DM dm,PetscInt p,PetscBool *flag)
19402bf73ac6SHong Zhang {
19412bf73ac6SHong Zhang   PetscInt       i;
19422bf73ac6SHong Zhang 
19432bf73ac6SHong Zhang   PetscFunctionBegin;
19442bf73ac6SHong Zhang   *flag = PETSC_FALSE;
19452bf73ac6SHong Zhang 
19462bf73ac6SHong Zhang   if (dm->setupcalled) { /* DMNetworkGetGlobalVertexIndex() requires DMSetUp() be called */
19472bf73ac6SHong Zhang     DM_Network     *network = (DM_Network*)dm->data;
19482bf73ac6SHong Zhang     PetscInt       gidx;
19499566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVertexIndex(dm,p,&gidx));
1950*daad07d3SAidan Hamilton     PetscCall(PetscTableFind(network->cloneshared->svtable,gidx+1,&i));
19512bf73ac6SHong Zhang     if (i) *flag = PETSC_TRUE;
19522bf73ac6SHong Zhang   } else { /* would be removed? */
19532bf73ac6SHong Zhang     PetscInt       nv;
19542bf73ac6SHong Zhang     const PetscInt *vtx;
19559566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetSharedVertices(dm,&nv,&vtx));
19562bf73ac6SHong Zhang     for (i=0; i<nv; i++) {
19572bf73ac6SHong Zhang       if (p == vtx[i]) {
19582bf73ac6SHong Zhang         *flag = PETSC_TRUE;
19592bf73ac6SHong Zhang         break;
19602bf73ac6SHong Zhang       }
19612bf73ac6SHong Zhang     }
19622bf73ac6SHong Zhang   }
19632bf73ac6SHong Zhang   PetscFunctionReturn(0);
19642bf73ac6SHong Zhang }
19652bf73ac6SHong Zhang 
19662bf73ac6SHong Zhang /*@
19675f2c45f1SShri Abhyankar   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
19685f2c45f1SShri Abhyankar 
19695f2c45f1SShri Abhyankar   Not Collective
19705f2c45f1SShri Abhyankar 
19715f2c45f1SShri Abhyankar   Input Parameters:
19722bf73ac6SHong Zhang + dm - the DMNetwork object
1973a2b725a8SWilliam Gropp - p - the vertex point
19745f2c45f1SShri Abhyankar 
19755f2c45f1SShri Abhyankar   Output Parameter:
19765f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point
19775f2c45f1SShri Abhyankar 
197897bb938eSShri Abhyankar   Level: beginner
19795f2c45f1SShri Abhyankar 
1980db781477SPatrick Sanan .seealso: `DMNetworkGetConnectedVertices()`, `DMNetworkGetVertexRange()`, `DMNetworkIsSharedVertex()`
19815f2c45f1SShri Abhyankar @*/
19825f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
19835f2c45f1SShri Abhyankar {
19845f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
19855f2c45f1SShri Abhyankar   PetscInt       offsetg;
19865f2c45f1SShri Abhyankar   PetscSection   sectiong;
19875f2c45f1SShri Abhyankar 
19885f2c45f1SShri Abhyankar   PetscFunctionBegin;
19895f2c45f1SShri Abhyankar   *isghost = PETSC_FALSE;
19909566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(network->plex,&sectiong));
19919566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(sectiong,p,&offsetg));
19925f2c45f1SShri Abhyankar   if (offsetg < 0) *isghost = PETSC_TRUE;
19935f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
19945f2c45f1SShri Abhyankar }
19955f2c45f1SShri Abhyankar 
19965f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm)
19975f2c45f1SShri Abhyankar {
19985f2c45f1SShri Abhyankar   PetscFunctionBegin;
1999e3b1a053SGetnet Betrie   PetscCall(PetscLogEventBegin(DMNetwork_SetUpNetwork,dm,0,0,0));
2000*daad07d3SAidan Hamilton   PetscCall(DMNetworkFinalizeComponents(dm));
20015b3f975aSHong Zhang   /* View dmnetwork */
20029566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm,NULL,"-dmnetwork_view"));
2003e3b1a053SGetnet Betrie   PetscCall(PetscLogEventEnd(DMNetwork_SetUpNetwork,dm,0,0,0));
20045f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
20055f2c45f1SShri Abhyankar }
20065f2c45f1SShri Abhyankar 
20071ad426b7SHong Zhang /*@
200817df6e9eSHong Zhang   DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
20091ad426b7SHong Zhang       -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
20101ad426b7SHong Zhang 
20111ad426b7SHong Zhang   Collective
20121ad426b7SHong Zhang 
20131ad426b7SHong Zhang   Input Parameters:
20142bf73ac6SHong Zhang + dm - the DMNetwork object
201583b2e829SHong Zhang . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
201683b2e829SHong Zhang - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
20171ad426b7SHong Zhang 
20181ad426b7SHong Zhang  Level: intermediate
20191ad426b7SHong Zhang 
20201ad426b7SHong Zhang @*/
202183b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
20221ad426b7SHong Zhang {
20231ad426b7SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
2024*daad07d3SAidan Hamilton   PetscInt       nVertices = network->cloneshared->nVertices;
20251ad426b7SHong Zhang 
20261ad426b7SHong Zhang   PetscFunctionBegin;
202783b2e829SHong Zhang   network->userEdgeJacobian   = eflg;
202883b2e829SHong Zhang   network->userVertexJacobian = vflg;
20298675203cSHong Zhang 
20308675203cSHong Zhang   if (eflg && !network->Je) {
2031*daad07d3SAidan Hamilton     PetscCall(PetscCalloc1(3*network->cloneshared->nEdges,&network->Je));
20328675203cSHong Zhang   }
20338675203cSHong Zhang 
203466b4e0ffSHong Zhang   if (vflg && !network->Jv && nVertices) {
2035*daad07d3SAidan Hamilton     PetscInt       i,*vptr,nedges,vStart=network->cloneshared->vStart;
203666b4e0ffSHong Zhang     PetscInt       nedges_total;
20378675203cSHong Zhang     const PetscInt *edges;
20388675203cSHong Zhang 
20398675203cSHong Zhang     /* count nvertex_total */
20408675203cSHong Zhang     nedges_total = 0;
20419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nVertices+1,&vptr));
20428675203cSHong Zhang 
20438675203cSHong Zhang     vptr[0] = 0;
20448675203cSHong Zhang     for (i=0; i<nVertices; i++) {
20459566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges));
20468675203cSHong Zhang       nedges_total += nedges;
20478675203cSHong Zhang       vptr[i+1] = vptr[i] + 2*nedges + 1;
20488675203cSHong Zhang     }
20498675203cSHong Zhang 
20509566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(2*nedges_total+nVertices,&network->Jv));
20518675203cSHong Zhang     network->Jvptr = vptr;
20528675203cSHong Zhang   }
20531ad426b7SHong Zhang   PetscFunctionReturn(0);
20541ad426b7SHong Zhang }
20551ad426b7SHong Zhang 
20561ad426b7SHong Zhang /*@
205783b2e829SHong Zhang   DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
205883b2e829SHong Zhang 
205983b2e829SHong Zhang   Not Collective
206083b2e829SHong Zhang 
206183b2e829SHong Zhang   Input Parameters:
20622bf73ac6SHong Zhang + dm - the DMNetwork object
206383b2e829SHong Zhang . p - the edge point
20643e97b6e8SHong Zhang - J - array (size = 3) of Jacobian submatrices for this edge point:
20653e97b6e8SHong Zhang         J[0]: this edge
2066d842c372SHong Zhang         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()
206783b2e829SHong Zhang 
206897bb938eSShri Abhyankar   Level: advanced
206983b2e829SHong Zhang 
2070db781477SPatrick Sanan .seealso: `DMNetworkVertexSetMatrix()`
207183b2e829SHong Zhang @*/
207283b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
207383b2e829SHong Zhang {
207483b2e829SHong Zhang   DM_Network *network=(DM_Network*)dm->data;
207583b2e829SHong Zhang 
207683b2e829SHong Zhang   PetscFunctionBegin;
20775c6496baSHong Zhang   PetscCheck(network->Je,PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
20788675203cSHong Zhang 
20798675203cSHong Zhang   if (J) {
2080883e35e8SHong Zhang     network->Je[3*p]   = J[0];
2081883e35e8SHong Zhang     network->Je[3*p+1] = J[1];
2082883e35e8SHong Zhang     network->Je[3*p+2] = J[2];
20838675203cSHong Zhang   }
208483b2e829SHong Zhang   PetscFunctionReturn(0);
208583b2e829SHong Zhang }
208683b2e829SHong Zhang 
208783b2e829SHong Zhang /*@
208876ddfea5SHong Zhang   DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
20891ad426b7SHong Zhang 
20901ad426b7SHong Zhang   Not Collective
20911ad426b7SHong Zhang 
20921ad426b7SHong Zhang   Input Parameters:
20931ad426b7SHong Zhang + dm - The DMNetwork object
20941ad426b7SHong Zhang . p - the vertex point
20953e97b6e8SHong Zhang - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
20963e97b6e8SHong Zhang         J[0]:       this vertex
20973e97b6e8SHong Zhang         J[1+2*i]:   i-th supporting edge
20983e97b6e8SHong Zhang         J[1+2*i+1]: i-th connected vertex
20991ad426b7SHong Zhang 
210097bb938eSShri Abhyankar   Level: advanced
21011ad426b7SHong Zhang 
2102db781477SPatrick Sanan .seealso: `DMNetworkEdgeSetMatrix()`
21031ad426b7SHong Zhang @*/
2104883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
21055f2c45f1SShri Abhyankar {
21065f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
2107*daad07d3SAidan Hamilton   PetscInt       i,*vptr,nedges,vStart=network->cloneshared->vStart;
2108883e35e8SHong Zhang   const PetscInt *edges;
21095f2c45f1SShri Abhyankar 
21105f2c45f1SShri Abhyankar   PetscFunctionBegin;
21115c6496baSHong Zhang   PetscCheck(network->Jv,PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
2112883e35e8SHong Zhang 
21138675203cSHong Zhang   if (J) {
2114883e35e8SHong Zhang     vptr = network->Jvptr;
21153e97b6e8SHong Zhang     network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
21163e97b6e8SHong Zhang 
21173e97b6e8SHong Zhang     /* Set Jacobian for each supporting edge and connected vertex */
21189566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetSupportingEdges(dm,p,&nedges,&edges));
2119883e35e8SHong Zhang     for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
21208675203cSHong Zhang   }
2121883e35e8SHong Zhang   PetscFunctionReturn(0);
2122883e35e8SHong Zhang }
2123883e35e8SHong Zhang 
21249fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
21255cf7da58SHong Zhang {
21265cf7da58SHong Zhang   PetscInt       j;
21275cf7da58SHong Zhang   PetscScalar    val=(PetscScalar)ncols;
21285cf7da58SHong Zhang 
21295cf7da58SHong Zhang   PetscFunctionBegin;
21305cf7da58SHong Zhang   if (!ghost) {
21315cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
21329566063dSJacob Faibussowitsch       PetscCall(VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES));
21335cf7da58SHong Zhang     }
21345cf7da58SHong Zhang   } else {
21355cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
21369566063dSJacob Faibussowitsch       PetscCall(VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES));
21375cf7da58SHong Zhang     }
21385cf7da58SHong Zhang   }
21395cf7da58SHong Zhang   PetscFunctionReturn(0);
21405cf7da58SHong Zhang }
21415cf7da58SHong Zhang 
21429fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
21435cf7da58SHong Zhang {
21445cf7da58SHong Zhang   PetscInt       j,ncols_u;
21455cf7da58SHong Zhang   PetscScalar    val;
21465cf7da58SHong Zhang 
21475cf7da58SHong Zhang   PetscFunctionBegin;
21485cf7da58SHong Zhang   if (!ghost) {
21495cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
21509566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Ju,j,&ncols_u,NULL,NULL));
21515cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
21529566063dSJacob Faibussowitsch       PetscCall(VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES));
21539566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Ju,j,&ncols_u,NULL,NULL));
21545cf7da58SHong Zhang     }
21555cf7da58SHong Zhang   } else {
21565cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
21579566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Ju,j,&ncols_u,NULL,NULL));
21585cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
21599566063dSJacob Faibussowitsch       PetscCall(VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES));
21609566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Ju,j,&ncols_u,NULL,NULL));
21615cf7da58SHong Zhang     }
21625cf7da58SHong Zhang   }
21635cf7da58SHong Zhang   PetscFunctionReturn(0);
21645cf7da58SHong Zhang }
21655cf7da58SHong Zhang 
21669fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
21675cf7da58SHong Zhang {
21685cf7da58SHong Zhang   PetscFunctionBegin;
21695cf7da58SHong Zhang   if (Ju) {
21709566063dSJacob Faibussowitsch     PetscCall(MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz));
21715cf7da58SHong Zhang   } else {
21729566063dSJacob Faibussowitsch     PetscCall(MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz));
21735cf7da58SHong Zhang   }
21745cf7da58SHong Zhang   PetscFunctionReturn(0);
21755cf7da58SHong Zhang }
21765cf7da58SHong Zhang 
21779fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
2178883e35e8SHong Zhang {
2179883e35e8SHong Zhang   PetscInt       j,*cols;
2180883e35e8SHong Zhang   PetscScalar    *zeros;
2181883e35e8SHong Zhang 
2182883e35e8SHong Zhang   PetscFunctionBegin;
21839566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(ncols,&cols,nrows*ncols,&zeros));
2184883e35e8SHong Zhang   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
21859566063dSJacob Faibussowitsch   PetscCall(MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES));
21869566063dSJacob Faibussowitsch   PetscCall(PetscFree2(cols,zeros));
21871ad426b7SHong Zhang   PetscFunctionReturn(0);
21881ad426b7SHong Zhang }
2189a4e85ca8SHong Zhang 
21909fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
21913e97b6e8SHong Zhang {
21923e97b6e8SHong Zhang   PetscInt       j,M,N,row,col,ncols_u;
21933e97b6e8SHong Zhang   const PetscInt *cols;
21943e97b6e8SHong Zhang   PetscScalar    zero=0.0;
21953e97b6e8SHong Zhang 
21963e97b6e8SHong Zhang   PetscFunctionBegin;
21979566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Ju,&M,&N));
219863a3b9bcSJacob Faibussowitsch   PetscCheck(nrows == M && ncols == N,PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%" PetscInt_FMT " by %" PetscInt_FMT " must equal %" PetscInt_FMT " by %" PetscInt_FMT,nrows,ncols,M,N);
21993e97b6e8SHong Zhang 
22003e97b6e8SHong Zhang   for (row=0; row<nrows; row++) {
22019566063dSJacob Faibussowitsch     PetscCall(MatGetRow(Ju,row,&ncols_u,&cols,NULL));
22023e97b6e8SHong Zhang     for (j=0; j<ncols_u; j++) {
22033e97b6e8SHong Zhang       col = cols[j] + cstart;
22049566063dSJacob Faibussowitsch       PetscCall(MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES));
22053e97b6e8SHong Zhang     }
22069566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(Ju,row,&ncols_u,&cols,NULL));
22073e97b6e8SHong Zhang   }
22083e97b6e8SHong Zhang   PetscFunctionReturn(0);
22093e97b6e8SHong Zhang }
22101ad426b7SHong Zhang 
22119fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
2212a4e85ca8SHong Zhang {
2213a4e85ca8SHong Zhang   PetscFunctionBegin;
2214a4e85ca8SHong Zhang   if (Ju) {
22159566063dSJacob Faibussowitsch     PetscCall(MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J));
2216a4e85ca8SHong Zhang   } else {
22179566063dSJacob Faibussowitsch     PetscCall(MatSetDenseblock_private(nrows,rows,ncols,cstart,J));
2218a4e85ca8SHong Zhang   }
2219a4e85ca8SHong Zhang   PetscFunctionReturn(0);
2220a4e85ca8SHong Zhang }
2221a4e85ca8SHong Zhang 
222224121865SAdrian 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.
222324121865SAdrian Maldonado */
222424121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
222524121865SAdrian Maldonado {
222624121865SAdrian Maldonado   PetscInt       i,size,dof;
222724121865SAdrian Maldonado   PetscInt       *glob2loc;
222824121865SAdrian Maldonado 
222924121865SAdrian Maldonado   PetscFunctionBegin;
22309566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(localsec,&size));
22319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size,&glob2loc));
223224121865SAdrian Maldonado 
223324121865SAdrian Maldonado   for (i = 0; i < size; i++) {
22349566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(globalsec,i,&dof));
223524121865SAdrian Maldonado     dof = (dof >= 0) ? dof : -(dof + 1);
223624121865SAdrian Maldonado     glob2loc[i] = dof;
223724121865SAdrian Maldonado   }
223824121865SAdrian Maldonado 
22399566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog));
224024121865SAdrian Maldonado #if 0
22419566063dSJacob Faibussowitsch   PetscCall(PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD));
224224121865SAdrian Maldonado #endif
224324121865SAdrian Maldonado   PetscFunctionReturn(0);
224424121865SAdrian Maldonado }
224524121865SAdrian Maldonado 
224601ad2aeeSHong Zhang #include <petsc/private/matimpl.h>
22479e1d080bSHong Zhang 
22489e1d080bSHong Zhang PetscErrorCode DMCreateMatrix_Network_Nest(DM dm,Mat *J)
22491ad426b7SHong Zhang {
22501ad426b7SHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
225124121865SAdrian Maldonado   PetscInt       eDof,vDof;
225224121865SAdrian Maldonado   Mat            j11,j12,j21,j22,bA[2][2];
22539e1d080bSHong Zhang   MPI_Comm       comm;
225424121865SAdrian Maldonado   ISLocalToGlobalMapping eISMap,vISMap;
225524121865SAdrian Maldonado 
22569e1d080bSHong Zhang   PetscFunctionBegin;
22579566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
225824121865SAdrian Maldonado 
22599566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof));
22609566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof));
226124121865SAdrian Maldonado 
22629566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &j11));
22639566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE));
22649566063dSJacob Faibussowitsch   PetscCall(MatSetType(j11, MATMPIAIJ));
226524121865SAdrian Maldonado 
22669566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &j12));
22679566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE));
22689566063dSJacob Faibussowitsch   PetscCall(MatSetType(j12, MATMPIAIJ));
226924121865SAdrian Maldonado 
22709566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &j21));
22719566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE));
22729566063dSJacob Faibussowitsch   PetscCall(MatSetType(j21, MATMPIAIJ));
227324121865SAdrian Maldonado 
22749566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &j22));
22759566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE));
22769566063dSJacob Faibussowitsch   PetscCall(MatSetType(j22, MATMPIAIJ));
227724121865SAdrian Maldonado 
22783f6a6bdaSHong Zhang   bA[0][0] = j11;
22793f6a6bdaSHong Zhang   bA[0][1] = j12;
22803f6a6bdaSHong Zhang   bA[1][0] = j21;
22813f6a6bdaSHong Zhang   bA[1][1] = j22;
228224121865SAdrian Maldonado 
22839566063dSJacob Faibussowitsch   PetscCall(CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap));
22849566063dSJacob Faibussowitsch   PetscCall(CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap));
228524121865SAdrian Maldonado 
22869566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(j11,eISMap,eISMap));
22879566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(j12,eISMap,vISMap));
22889566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(j21,vISMap,eISMap));
22899566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(j22,vISMap,vISMap));
229024121865SAdrian Maldonado 
22919566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j11));
22929566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j12));
22939566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j21));
22949566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j22));
229524121865SAdrian Maldonado 
22969566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J));
22979566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*J));
22989566063dSJacob Faibussowitsch   PetscCall(MatNestSetVecType(*J,VECNEST));
22999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j11));
23009566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j12));
23019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j21));
23029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j22));
230324121865SAdrian Maldonado 
23049566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY));
23059566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY));
23069566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
230724121865SAdrian Maldonado 
230824121865SAdrian Maldonado   /* Free structures */
23099566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&eISMap));
23109566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&vISMap));
231124121865SAdrian Maldonado   PetscFunctionReturn(0);
23129e1d080bSHong Zhang }
23139e1d080bSHong Zhang 
23149e1d080bSHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
23159e1d080bSHong Zhang {
23169e1d080bSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
23179e1d080bSHong Zhang   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
23189e1d080bSHong Zhang   PetscInt       cstart,ncols,j,e,v;
23199e1d080bSHong Zhang   PetscBool      ghost,ghost_vc,ghost2,isNest;
23209e1d080bSHong Zhang   Mat            Juser;
23219e1d080bSHong Zhang   PetscSection   sectionGlobal;
23229e1d080bSHong Zhang   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
23239e1d080bSHong Zhang   const PetscInt *edges,*cone;
23249e1d080bSHong Zhang   MPI_Comm       comm;
23259e1d080bSHong Zhang   MatType        mtype;
23269e1d080bSHong Zhang   Vec            vd_nz,vo_nz;
23279e1d080bSHong Zhang   PetscInt       *dnnz,*onnz;
23289e1d080bSHong Zhang   PetscScalar    *vdnz,*vonz;
23299e1d080bSHong Zhang 
23309e1d080bSHong Zhang   PetscFunctionBegin;
23319e1d080bSHong Zhang   mtype = dm->mattype;
23329566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype,MATNEST,&isNest));
23339e1d080bSHong Zhang   if (isNest) {
23349566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix_Network_Nest(dm,J));
23359566063dSJacob Faibussowitsch     PetscCall(MatSetDM(*J,dm));
23369e1d080bSHong Zhang     PetscFunctionReturn(0);
23379e1d080bSHong Zhang   }
23389e1d080bSHong Zhang 
23399e1d080bSHong Zhang   if (!network->userEdgeJacobian && !network->userVertexJacobian) {
2340a4e85ca8SHong Zhang     /* user does not provide Jacobian blocks */
23419566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix_Plex(network->plex,J));
23429566063dSJacob Faibussowitsch     PetscCall(MatSetDM(*J,dm));
23431ad426b7SHong Zhang     PetscFunctionReturn(0);
23441ad426b7SHong Zhang   }
23451ad426b7SHong Zhang 
23469566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dm),J));
23479566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(network->plex,&sectionGlobal));
23489566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize));
23499566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE));
23502a945128SHong Zhang 
23519566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J,MATAIJ));
23529566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(*J));
235389898e50SHong Zhang 
235489898e50SHong Zhang   /* (1) Set matrix preallocation */
235589898e50SHong Zhang   /*------------------------------*/
23569566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
23579566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm,&vd_nz));
23589566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(vd_nz,localSize,PETSC_DECIDE));
23599566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(vd_nz));
23609566063dSJacob Faibussowitsch   PetscCall(VecSet(vd_nz,0.0));
23619566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(vd_nz,&vo_nz));
2362840c2264SHong Zhang 
236389898e50SHong Zhang   /* Set preallocation for edges */
236489898e50SHong Zhang   /*-----------------------------*/
23659566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(dm,&eStart,&eEnd));
2366840c2264SHong Zhang 
23679566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(localSize,&rows));
2368840c2264SHong Zhang   for (e=eStart; e<eEnd; e++) {
2369840c2264SHong Zhang     /* Get row indices */
23709566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVecOffset(dm,e,ALL_COMPONENTS,&rstart));
23719566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(network->DofSection,e,&nrows));
2372840c2264SHong Zhang     if (nrows) {
2373840c2264SHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
2374840c2264SHong Zhang 
2375a5b23f4aSJose E. Roman       /* Set preallocation for connected vertices */
23769566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetConnectedVertices(dm,e,&cone));
2377840c2264SHong Zhang       for (v=0; v<2; v++) {
23789566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(network->DofSection,cone[v],&ncols));
2379840c2264SHong Zhang 
23808675203cSHong Zhang         if (network->Je) {
2381840c2264SHong Zhang           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
23828675203cSHong Zhang         } else Juser = NULL;
23839566063dSJacob Faibussowitsch         PetscCall(DMNetworkIsGhostVertex(dm,cone[v],&ghost));
23849566063dSJacob Faibussowitsch         PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz));
2385840c2264SHong Zhang       }
2386840c2264SHong Zhang 
238789898e50SHong Zhang       /* Set preallocation for edge self */
2388840c2264SHong Zhang       cstart = rstart;
23898675203cSHong Zhang       if (network->Je) {
2390840c2264SHong Zhang         Juser = network->Je[3*e]; /* Jacobian(e,e) */
23918675203cSHong Zhang       } else Juser = NULL;
23929566063dSJacob Faibussowitsch       PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz));
2393840c2264SHong Zhang     }
2394840c2264SHong Zhang   }
2395840c2264SHong Zhang 
239689898e50SHong Zhang   /* Set preallocation for vertices */
239789898e50SHong Zhang   /*--------------------------------*/
23989566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(dm,&vStart,&vEnd));
23998675203cSHong Zhang   if (vEnd - vStart) vptr = network->Jvptr;
2400840c2264SHong Zhang 
2401840c2264SHong Zhang   for (v=vStart; v<vEnd; v++) {
2402840c2264SHong Zhang     /* Get row indices */
24039566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&rstart));
24049566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(network->DofSection,v,&nrows));
2405840c2264SHong Zhang     if (!nrows) continue;
2406840c2264SHong Zhang 
24079566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(dm,v,&ghost));
2408bdcb62a2SHong Zhang     if (ghost) {
24099566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nrows,&rows_v));
2410bdcb62a2SHong Zhang     } else {
2411bdcb62a2SHong Zhang       rows_v = rows;
2412bdcb62a2SHong Zhang     }
2413bdcb62a2SHong Zhang 
2414bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
2415840c2264SHong Zhang 
2416840c2264SHong Zhang     /* Get supporting edges and connected vertices */
24179566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetSupportingEdges(dm,v,&nedges,&edges));
2418840c2264SHong Zhang 
2419840c2264SHong Zhang     for (e=0; e<nedges; e++) {
2420840c2264SHong Zhang       /* Supporting edges */
24219566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetGlobalVecOffset(dm,edges[e],ALL_COMPONENTS,&cstart));
24229566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(network->DofSection,edges[e],&ncols));
2423840c2264SHong Zhang 
24248675203cSHong Zhang       if (network->Jv) {
2425840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
24268675203cSHong Zhang       } else Juser = NULL;
24279566063dSJacob Faibussowitsch       PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz));
2428840c2264SHong Zhang 
2429840c2264SHong Zhang       /* Connected vertices */
24309566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetConnectedVertices(dm,edges[e],&cone));
2431840c2264SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
24329566063dSJacob Faibussowitsch       PetscCall(DMNetworkIsGhostVertex(dm,vc,&ghost_vc));
2433840c2264SHong Zhang 
24349566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(network->DofSection,vc,&ncols));
2435840c2264SHong Zhang 
24368675203cSHong Zhang       if (network->Jv) {
2437840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
24388675203cSHong Zhang       } else Juser = NULL;
2439e102a522SHong Zhang       if (ghost_vc||ghost) {
2440e102a522SHong Zhang         ghost2 = PETSC_TRUE;
2441e102a522SHong Zhang       } else {
2442e102a522SHong Zhang         ghost2 = PETSC_FALSE;
2443e102a522SHong Zhang       }
24449566063dSJacob Faibussowitsch       PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz));
2445840c2264SHong Zhang     }
2446840c2264SHong Zhang 
244789898e50SHong Zhang     /* Set preallocation for vertex self */
24489566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(dm,v,&ghost));
2449840c2264SHong Zhang     if (!ghost) {
24509566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&cstart));
24518675203cSHong Zhang       if (network->Jv) {
2452840c2264SHong Zhang         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
24538675203cSHong Zhang       } else Juser = NULL;
24549566063dSJacob Faibussowitsch       PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz));
2455840c2264SHong Zhang     }
24561baa6e33SBarry Smith     if (ghost) PetscCall(PetscFree(rows_v));
2457840c2264SHong Zhang   }
2458840c2264SHong Zhang 
24599566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(vd_nz));
24609566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(vo_nz));
24615cf7da58SHong Zhang 
24629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(localSize,&dnnz,localSize,&onnz));
24635cf7da58SHong Zhang 
24649566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(vd_nz));
24659566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(vo_nz));
2466840c2264SHong Zhang 
24679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vd_nz,&vdnz));
24689566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vo_nz,&vonz));
2469840c2264SHong Zhang   for (j=0; j<localSize; j++) {
2470e102a522SHong Zhang     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
2471e102a522SHong Zhang     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
2472840c2264SHong Zhang   }
24739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vd_nz,&vdnz));
24749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vo_nz,&vonz));
24759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vd_nz));
24769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vo_nz));
2477840c2264SHong Zhang 
24789566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*J,0,dnnz));
24799566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz));
24809566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
24815cf7da58SHong Zhang 
24829566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dnnz,onnz));
24835cf7da58SHong Zhang 
248489898e50SHong Zhang   /* (2) Set matrix entries for edges */
248589898e50SHong Zhang   /*----------------------------------*/
24861ad426b7SHong Zhang   for (e=eStart; e<eEnd; e++) {
2487bfbc38dcSHong Zhang     /* Get row indices */
24889566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVecOffset(dm,e,ALL_COMPONENTS,&rstart));
24899566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(network->DofSection,e,&nrows));
24904b976069SHong Zhang     if (nrows) {
249117df6e9eSHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
24921ad426b7SHong Zhang 
2493a5b23f4aSJose E. Roman       /* Set matrix entries for connected vertices */
24949566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetConnectedVertices(dm,e,&cone));
2495bfbc38dcSHong Zhang       for (v=0; v<2; v++) {
24969566063dSJacob Faibussowitsch         PetscCall(DMNetworkGetGlobalVecOffset(dm,cone[v],ALL_COMPONENTS,&cstart));
24979566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(network->DofSection,cone[v],&ncols));
24983e97b6e8SHong Zhang 
24998675203cSHong Zhang         if (network->Je) {
2500a4e85ca8SHong Zhang           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
25018675203cSHong Zhang         } else Juser = NULL;
25029566063dSJacob Faibussowitsch         PetscCall(MatSetblock_private(Juser,nrows,rows,ncols,cstart,J));
2503bfbc38dcSHong Zhang       }
250417df6e9eSHong Zhang 
2505bfbc38dcSHong Zhang       /* Set matrix entries for edge self */
25063e97b6e8SHong Zhang       cstart = rstart;
25078675203cSHong Zhang       if (network->Je) {
2508a4e85ca8SHong Zhang         Juser = network->Je[3*e]; /* Jacobian(e,e) */
25098675203cSHong Zhang       } else Juser = NULL;
25109566063dSJacob Faibussowitsch       PetscCall(MatSetblock_private(Juser,nrows,rows,nrows,cstart,J));
25111ad426b7SHong Zhang     }
25124b976069SHong Zhang   }
25131ad426b7SHong Zhang 
2514bfbc38dcSHong Zhang   /* Set matrix entries for vertices */
251583b2e829SHong Zhang   /*---------------------------------*/
25161ad426b7SHong Zhang   for (v=vStart; v<vEnd; v++) {
2517bfbc38dcSHong Zhang     /* Get row indices */
25189566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&rstart));
25199566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(network->DofSection,v,&nrows));
25204b976069SHong Zhang     if (!nrows) continue;
2521596e729fSHong Zhang 
25229566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(dm,v,&ghost));
2523bdcb62a2SHong Zhang     if (ghost) {
25249566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nrows,&rows_v));
2525bdcb62a2SHong Zhang     } else {
2526bdcb62a2SHong Zhang       rows_v = rows;
2527bdcb62a2SHong Zhang     }
2528bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
2529596e729fSHong Zhang 
2530bfbc38dcSHong Zhang     /* Get supporting edges and connected vertices */
25319566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetSupportingEdges(dm,v,&nedges,&edges));
2532596e729fSHong Zhang 
2533596e729fSHong Zhang     for (e=0; e<nedges; e++) {
2534bfbc38dcSHong Zhang       /* Supporting edges */
25359566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetGlobalVecOffset(dm,edges[e],ALL_COMPONENTS,&cstart));
25369566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(network->DofSection,edges[e],&ncols));
2537596e729fSHong Zhang 
25388675203cSHong Zhang       if (network->Jv) {
2539a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
25408675203cSHong Zhang       } else Juser = NULL;
25419566063dSJacob Faibussowitsch       PetscCall(MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J));
2542596e729fSHong Zhang 
2543bfbc38dcSHong Zhang       /* Connected vertices */
25449566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetConnectedVertices(dm,edges[e],&cone));
25452a945128SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
25462a945128SHong Zhang 
25479566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetGlobalVecOffset(dm,vc,ALL_COMPONENTS,&cstart));
25489566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(network->DofSection,vc,&ncols));
2549a4e85ca8SHong Zhang 
25508675203cSHong Zhang       if (network->Jv) {
2551a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
25528675203cSHong Zhang       } else Juser = NULL;
25539566063dSJacob Faibussowitsch       PetscCall(MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J));
2554596e729fSHong Zhang     }
2555596e729fSHong Zhang 
2556bfbc38dcSHong Zhang     /* Set matrix entries for vertex self */
25571ad426b7SHong Zhang     if (!ghost) {
25589566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&cstart));
25598675203cSHong Zhang       if (network->Jv) {
2560a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
25618675203cSHong Zhang       } else Juser = NULL;
25629566063dSJacob Faibussowitsch       PetscCall(MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J));
2563bdcb62a2SHong Zhang     }
25641baa6e33SBarry Smith     if (ghost) PetscCall(PetscFree(rows_v));
25651ad426b7SHong Zhang   }
25669566063dSJacob Faibussowitsch   PetscCall(PetscFree(rows));
2567bdcb62a2SHong Zhang 
25689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY));
25699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY));
2570dd6f46cdSHong Zhang 
25719566063dSJacob Faibussowitsch   PetscCall(MatSetDM(*J,dm));
25725f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
25735f2c45f1SShri Abhyankar }
25745f2c45f1SShri Abhyankar 
2575*daad07d3SAidan Hamilton static PetscErrorCode DMNetworkDestroyComponentData(DM dm)
2576*daad07d3SAidan Hamilton {
2577*daad07d3SAidan Hamilton   DM_Network     *network = (DM_Network*)dm->data;
2578*daad07d3SAidan Hamilton   PetscInt       j,np;
2579*daad07d3SAidan Hamilton 
2580*daad07d3SAidan Hamilton   PetscFunctionBegin;
2581*daad07d3SAidan Hamilton   if (network->header) {
2582*daad07d3SAidan Hamilton     np = network->cloneshared->pEnd - network->cloneshared->pStart;
2583*daad07d3SAidan Hamilton     for (j=0; j < np; j++) {
2584*daad07d3SAidan Hamilton       PetscCall(PetscFree5(network->header[j].size,network->header[j].key,network->header[j].offset,network->header[j].nvar,network->header[j].offsetvarrel));
2585*daad07d3SAidan Hamilton       PetscCall(PetscFree(network->cvalue[j].data));
2586*daad07d3SAidan Hamilton     }
2587*daad07d3SAidan Hamilton     PetscCall(PetscFree2(network->header,network->cvalue));
2588*daad07d3SAidan Hamilton   }
2589*daad07d3SAidan Hamilton   PetscFunctionReturn(0);
2590*daad07d3SAidan Hamilton }
2591*daad07d3SAidan Hamilton 
25925f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm)
25935f2c45f1SShri Abhyankar {
25945f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
2595aeb78aa7SBarry Smith   PetscInt       j;
25965f2c45f1SShri Abhyankar 
25975f2c45f1SShri Abhyankar   PetscFunctionBegin;
2598*daad07d3SAidan Hamilton   /*
2599*daad07d3SAidan Hamilton     Developers Note: Due to the mixed topological definition of DMNetwork and data defined ON the
2600*daad07d3SAidan Hamilton     network like DofSection, DataSection, *componentdataarray, and friends, when cloning, we share
2601*daad07d3SAidan Hamilton     only the true topological data, and make our own data ON the network. Thus refct only refers
2602*daad07d3SAidan Hamilton     to the number of references to topological data, and data ON the network is always destroyed.
2603*daad07d3SAidan Hamilton     It is understood this is atypical for a DM, but is very intentional.
2604*daad07d3SAidan Hamilton   */
2605*daad07d3SAidan Hamilton 
2606*daad07d3SAidan Hamilton   /* Always destroy data ON the network */
26079566063dSJacob Faibussowitsch   PetscCall(PetscFree(network->Je));
260883b2e829SHong Zhang   if (network->Jv) {
26099566063dSJacob Faibussowitsch     PetscCall(PetscFree(network->Jvptr));
26109566063dSJacob Faibussowitsch     PetscCall(PetscFree(network->Jv));
26111ad426b7SHong Zhang   }
2612*daad07d3SAidan Hamilton   PetscCall(PetscSectionDestroy(&network->DataSection));
2613*daad07d3SAidan Hamilton   PetscCall(PetscSectionDestroy(&network->DofSection));
2614*daad07d3SAidan Hamilton   PetscCall(PetscFree(network->component));
2615*daad07d3SAidan Hamilton   PetscCall(PetscFree(network->componentdataarray));
2616*daad07d3SAidan Hamilton   PetscCall(DMNetworkDestroyComponentData(dm));
261713c2a604SAdrian Maldonado 
2618*daad07d3SAidan Hamilton   PetscCall(DMDestroy(&network->plex)); /* this is cloned in DMClone_Network, so safe to destroy */
2619*daad07d3SAidan Hamilton 
2620*daad07d3SAidan Hamilton   /*
2621*daad07d3SAidan Hamilton   Developers Note: The DMNetworkVertexInfo and DMNetworkEdgeInfo data structures are completely
2622*daad07d3SAidan Hamilton   destroyed as they are again a mix of topological data:
2623*daad07d3SAidan Hamilton     ISLocalToGlobalMapping            mapping;
2624*daad07d3SAidan Hamilton     PetscSF                           sf;
2625*daad07d3SAidan Hamilton   and data ON the network
2626*daad07d3SAidan Hamilton     PetscSection                      DofSection;
2627*daad07d3SAidan Hamilton     PetscSection                      GlobalDofSection;
2628*daad07d3SAidan Hamilton   And the only way to access them currently is through DMNetworkAssembleGraphStructures which assembles
2629*daad07d3SAidan Hamilton   everything. So we must destroy everything and require DMNetworkAssembleGraphStructures is called again
2630*daad07d3SAidan Hamilton   for any clone.
2631*daad07d3SAidan Hamilton   */
26329566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&network->vertex.mapping));
26339566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&network->vertex.DofSection));
26349566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&network->vertex.GlobalDofSection));
26359566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&network->vertex.sf));
263613c2a604SAdrian Maldonado   /* edge */
26379566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&network->edge.mapping));
26389566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&network->edge.DofSection));
26399566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&network->edge.GlobalDofSection));
26409566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&network->edge.sf));
2641*daad07d3SAidan Hamilton   /* Destroy the potentially cloneshared data */
2642*daad07d3SAidan Hamilton   if (--network->cloneshared->refct <= 0) {
2643*daad07d3SAidan Hamilton   /* Developer Note: I'm not sure if vltog can be reused or not, as I'm not sure what it's purpose is. I
2644*daad07d3SAidan Hamilton      naively think it can be reused. */
2645*daad07d3SAidan Hamilton     PetscCall(PetscFree(network->cloneshared->vltog));
2646*daad07d3SAidan Hamilton     for (j=0; j<network->cloneshared->Nsvtx; j++) PetscCall(PetscFree(network->cloneshared->svtx[j].sv));
2647*daad07d3SAidan Hamilton     PetscCall(PetscFree(network->cloneshared->svtx));
2648*daad07d3SAidan Hamilton     PetscCall(PetscFree2(network->cloneshared->subnetedge,network->cloneshared->subnetvtx));
2649*daad07d3SAidan Hamilton     PetscCall(PetscTableDestroy(&network->cloneshared->svtable));
2650*daad07d3SAidan Hamilton     PetscCall(PetscFree(network->cloneshared->subnet));
2651*daad07d3SAidan Hamilton     PetscCall(PetscFree(network->cloneshared));
2652*daad07d3SAidan Hamilton   }
2653*daad07d3SAidan Hamilton   PetscCall(PetscFree(network)); /* Always freed as this structure is copied in a clone, not cloneshared */
26545f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
26555f2c45f1SShri Abhyankar }
26565f2c45f1SShri Abhyankar 
26575f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm,PetscViewer viewer)
26585f2c45f1SShri Abhyankar {
2659caf410d2SHong Zhang   PetscBool      iascii;
2660caf410d2SHong Zhang   PetscMPIInt    rank;
26615f2c45f1SShri Abhyankar 
26625f2c45f1SShri Abhyankar   PetscFunctionBegin;
2663caf410d2SHong Zhang   PetscValidHeaderSpecific(dm,DM_CLASSID, 1);
2664caf410d2SHong Zhang   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
26659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
26669566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
2667caf410d2SHong Zhang   if (iascii) {
2668caf410d2SHong Zhang     const PetscInt *cone,*vtx,*edges;
26695c6496baSHong Zhang     PetscInt       vfrom,vto,i,j,nv,ne,nsv,p,nsubnet;
26702bf73ac6SHong Zhang     DM_Network     *network = (DM_Network*)dm->data;
2671caf410d2SHong Zhang 
2672*daad07d3SAidan Hamilton     nsubnet = network->cloneshared->Nsubnet; /* num of subnetworks */
2673dd400576SPatrick Sanan     if (rank == 0) {
2674*daad07d3SAidan Hamilton       PetscCall(PetscPrintf(PETSC_COMM_SELF,"  NSubnets: %" PetscInt_FMT "; NEdges: %" PetscInt_FMT "; NVertices: %" PetscInt_FMT "; NSharedVertices: %" PetscInt_FMT ".\n",nsubnet,network->cloneshared->NEdges,network->cloneshared->NVertices,network->cloneshared->Nsvtx));
26752bf73ac6SHong Zhang     }
26762bf73ac6SHong Zhang 
26779566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetSharedVertices(dm,&nsv,NULL));
26789566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
2679*daad07d3SAidan Hamilton     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] nEdges: %" PetscInt_FMT "; nVertices: %" PetscInt_FMT "; nSharedVertices: %" PetscInt_FMT "\n",rank,network->cloneshared->nEdges,network->cloneshared->nVertices,nsv));
2680caf410d2SHong Zhang 
2681caf410d2SHong Zhang     for (i=0; i<nsubnet; i++) {
26829566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetSubnetwork(dm,i,&nv,&ne,&vtx,&edges));
2683caf410d2SHong Zhang       if (ne) {
26849566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "     Subnet %" PetscInt_FMT ": nEdges %" PetscInt_FMT ", nVertices(include shared vertices) %" PetscInt_FMT "\n",i,ne,nv));
2685caf410d2SHong Zhang         for (j=0; j<ne; j++) {
2686caf410d2SHong Zhang           p = edges[j];
26879566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetConnectedVertices(dm,p,&cone));
26889566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetGlobalVertexIndex(dm,cone[0],&vfrom));
26899566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetGlobalVertexIndex(dm,cone[1],&vto));
26909566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetGlobalEdgeIndex(dm,edges[j],&p));
26919566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "       edge %" PetscInt_FMT ": %" PetscInt_FMT " ----> %" PetscInt_FMT "\n",p,vfrom,vto));
2692caf410d2SHong Zhang         }
2693caf410d2SHong Zhang       }
2694caf410d2SHong Zhang     }
26952bf73ac6SHong Zhang 
26962bf73ac6SHong Zhang     /* Shared vertices */
26979566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetSharedVertices(dm,NULL,&vtx));
26985c6496baSHong Zhang     if (nsv) {
26995c6496baSHong Zhang       PetscInt       gidx;
27002bf73ac6SHong Zhang       PetscBool      ghost;
27015c6496baSHong Zhang       const PetscInt *sv=NULL;
27025c6496baSHong Zhang 
27039566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "     SharedVertices:\n"));
27045c6496baSHong Zhang       for (i=0; i<nsv; i++) {
27059566063dSJacob Faibussowitsch         PetscCall(DMNetworkIsGhostVertex(dm,vtx[i],&ghost));
27062bf73ac6SHong Zhang         if (ghost) continue;
27072bf73ac6SHong Zhang 
27089566063dSJacob Faibussowitsch         PetscCall(DMNetworkSharedVertexGetInfo(dm,vtx[i],&gidx,&nv,&sv));
27099566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "       svtx %" PetscInt_FMT ": global index %" PetscInt_FMT ", subnet[%" PetscInt_FMT "].%" PetscInt_FMT " ---->\n",i,gidx,sv[0],sv[1]));
27105c6496baSHong Zhang         for (j=1; j<nv; j++) {
27119566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "                                           ----> subnet[%" PetscInt_FMT "].%" PetscInt_FMT "\n",sv[2*j],sv[2*j+1]));
2712caf410d2SHong Zhang         }
2713caf410d2SHong Zhang       }
2714caf410d2SHong Zhang     }
27159566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
27169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
27175c6496baSHong Zhang   } else PetscCheck(iascii,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Viewer type %s not yet supported for DMNetwork writing",((PetscObject)viewer)->type_name);
27185f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
27195f2c45f1SShri Abhyankar }
27205f2c45f1SShri Abhyankar 
27215f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
27225f2c45f1SShri Abhyankar {
27235f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
27245f2c45f1SShri Abhyankar 
27255f2c45f1SShri Abhyankar   PetscFunctionBegin;
27269566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(network->plex,g,mode,l));
27275f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
27285f2c45f1SShri Abhyankar }
27295f2c45f1SShri Abhyankar 
27305f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
27315f2c45f1SShri Abhyankar {
27325f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
27335f2c45f1SShri Abhyankar 
27345f2c45f1SShri Abhyankar   PetscFunctionBegin;
27359566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(network->plex,g,mode,l));
27365f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
27375f2c45f1SShri Abhyankar }
27385f2c45f1SShri Abhyankar 
27395f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
27405f2c45f1SShri Abhyankar {
27415f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
27425f2c45f1SShri Abhyankar 
27435f2c45f1SShri Abhyankar   PetscFunctionBegin;
27449566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(network->plex,l,mode,g));
27455f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
27465f2c45f1SShri Abhyankar }
27475f2c45f1SShri Abhyankar 
27485f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
27495f2c45f1SShri Abhyankar {
27505f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
27515f2c45f1SShri Abhyankar 
27525f2c45f1SShri Abhyankar   PetscFunctionBegin;
27539566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(network->plex,l,mode,g));
27545f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
27555f2c45f1SShri Abhyankar }
275622bbedd7SHong Zhang 
275722bbedd7SHong Zhang /*@
275864238763SRodolfo Oliveira   DMNetworkGetVertexLocalToGlobalOrdering - Get vertex global index
275922bbedd7SHong Zhang 
276022bbedd7SHong Zhang   Not collective
276122bbedd7SHong Zhang 
27627a7aea1fSJed Brown   Input Parameters:
276322bbedd7SHong Zhang + dm - the dm object
276422bbedd7SHong Zhang - vloc - local vertex ordering, start from 0
276522bbedd7SHong Zhang 
27667a7aea1fSJed Brown   Output Parameters:
2767f0fc11ceSJed Brown .  vg  - global vertex ordering, start from 0
276822bbedd7SHong Zhang 
276997bb938eSShri Abhyankar   Level: advanced
277022bbedd7SHong Zhang 
2771db781477SPatrick Sanan .seealso: `DMNetworkSetVertexLocalToGlobalOrdering()`
277222bbedd7SHong Zhang @*/
277322bbedd7SHong Zhang PetscErrorCode DMNetworkGetVertexLocalToGlobalOrdering(DM dm,PetscInt vloc,PetscInt *vg)
277422bbedd7SHong Zhang {
277522bbedd7SHong Zhang   DM_Network  *network = (DM_Network*)dm->data;
2776*daad07d3SAidan Hamilton   PetscInt    *vltog = network->cloneshared->vltog;
277722bbedd7SHong Zhang 
277822bbedd7SHong Zhang   PetscFunctionBegin;
27795c6496baSHong Zhang   PetscCheck(vltog,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Must call DMNetworkSetVertexLocalToGlobalOrdering() first");
278022bbedd7SHong Zhang   *vg = vltog[vloc];
278122bbedd7SHong Zhang   PetscFunctionReturn(0);
278222bbedd7SHong Zhang }
278322bbedd7SHong Zhang 
278422bbedd7SHong Zhang /*@
278564238763SRodolfo Oliveira   DMNetworkSetVertexLocalToGlobalOrdering - Create and setup vertex local to global map
278622bbedd7SHong Zhang 
278722bbedd7SHong Zhang   Collective
278822bbedd7SHong Zhang 
278922bbedd7SHong Zhang   Input Parameters:
2790f0fc11ceSJed Brown . dm - the dm object
279122bbedd7SHong Zhang 
279297bb938eSShri Abhyankar   Level: advanced
279322bbedd7SHong Zhang 
2794db781477SPatrick Sanan .seealso: `DMNetworkGetGlobalVertexIndex()`
279522bbedd7SHong Zhang @*/
279622bbedd7SHong Zhang PetscErrorCode DMNetworkSetVertexLocalToGlobalOrdering(DM dm)
279722bbedd7SHong Zhang {
279822bbedd7SHong Zhang   DM_Network        *network=(DM_Network*)dm->data;
279922bbedd7SHong Zhang   MPI_Comm          comm;
28002bf73ac6SHong Zhang   PetscMPIInt       rank,size,*displs=NULL,*recvcounts=NULL,remoterank;
280122bbedd7SHong Zhang   PetscBool         ghost;
280263029d29SHong Zhang   PetscInt          *vltog,nroots,nleaves,i,*vrange,k,N,lidx;
280322bbedd7SHong Zhang   const PetscSFNode *iremote;
280422bbedd7SHong Zhang   PetscSF           vsf;
280563029d29SHong Zhang   Vec               Vleaves,Vleaves_seq;
280663029d29SHong Zhang   VecScatter        ctx;
280763029d29SHong Zhang   PetscScalar       *varr,val;
280863029d29SHong Zhang   const PetscScalar *varr_read;
280922bbedd7SHong Zhang 
281022bbedd7SHong Zhang   PetscFunctionBegin;
28119566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
28129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
28139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
281422bbedd7SHong Zhang 
281522bbedd7SHong Zhang   if (size == 1) {
2816*daad07d3SAidan Hamilton     nroots = network->cloneshared->vEnd - network->cloneshared->vStart;
28179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nroots, &vltog));
281822bbedd7SHong Zhang     for (i=0; i<nroots; i++) vltog[i] = i;
2819*daad07d3SAidan Hamilton     network->cloneshared->vltog = vltog;
282022bbedd7SHong Zhang     PetscFunctionReturn(0);
282122bbedd7SHong Zhang   }
282222bbedd7SHong Zhang 
2823*daad07d3SAidan Hamilton   PetscCheck(network->cloneshared->distributecalled,comm, PETSC_ERR_ARG_WRONGSTATE,"Must call DMNetworkDistribute() first");
2824*daad07d3SAidan Hamilton   if (network->cloneshared->vltog) PetscCall(PetscFree(network->cloneshared->vltog));
282522bbedd7SHong Zhang 
2826*daad07d3SAidan Hamilton   PetscCall(DMNetworkSetSubMap_private(network->cloneshared->vStart,network->cloneshared->vEnd,&network->vertex.mapping));
28279566063dSJacob Faibussowitsch   PetscCall(PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf));
282822bbedd7SHong Zhang   vsf = network->vertex.sf;
282922bbedd7SHong Zhang 
28309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(size+1,&vrange,size,&displs,size,&recvcounts));
28319566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(vsf,&nroots,&nleaves,NULL,&iremote));
283222bbedd7SHong Zhang 
283322bbedd7SHong Zhang   for (i=0; i<size; i++) { displs[i] = i; recvcounts[i] = 1;}
283422bbedd7SHong Zhang 
283522bbedd7SHong Zhang   i         = nroots - nleaves; /* local number of vertices, excluding ghosts */
283622bbedd7SHong Zhang   vrange[0] = 0;
28379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgatherv(&i,1,MPIU_INT,vrange+1,recvcounts,displs,MPIU_INT,comm));
283867dd800bSHong Zhang   for (i=2; i<size+1; i++) {vrange[i] += vrange[i-1];}
283922bbedd7SHong Zhang 
28409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nroots, &vltog));
2841*daad07d3SAidan Hamilton   network->cloneshared->vltog = vltog;
284222bbedd7SHong Zhang 
284363029d29SHong Zhang   /* Set vltog for non-ghost vertices */
284463029d29SHong Zhang   k = 0;
284522bbedd7SHong Zhang   for (i=0; i<nroots; i++) {
2846*daad07d3SAidan Hamilton     PetscCall(DMNetworkIsGhostVertex(dm,i+network->cloneshared->vStart,&ghost));
284763029d29SHong Zhang     if (ghost) continue;
284863029d29SHong Zhang     vltog[i] = vrange[rank] + k++;
284922bbedd7SHong Zhang   }
28509566063dSJacob Faibussowitsch   PetscCall(PetscFree3(vrange,displs,recvcounts));
285163029d29SHong Zhang 
285263029d29SHong Zhang   /* Set vltog for ghost vertices */
285363029d29SHong Zhang   /* (a) create parallel Vleaves and sequential Vleaves_seq to convert local iremote[*].index to global index */
28549566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm,&Vleaves));
28559566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(Vleaves,2*nleaves,PETSC_DETERMINE));
28569566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(Vleaves));
28579566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Vleaves,&varr));
285863029d29SHong Zhang   for (i=0; i<nleaves; i++) {
285963029d29SHong Zhang     varr[2*i]   = (PetscScalar)(iremote[i].rank);  /* rank of remote process */
286063029d29SHong Zhang     varr[2*i+1] = (PetscScalar)(iremote[i].index); /* local index in remote process */
286163029d29SHong Zhang   }
28629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Vleaves,&varr));
286363029d29SHong Zhang 
286463029d29SHong Zhang   /* (b) scatter local info to remote processes via VecScatter() */
28659566063dSJacob Faibussowitsch   PetscCall(VecScatterCreateToAll(Vleaves,&ctx,&Vleaves_seq));
28669566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD));
28679566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD));
286863029d29SHong Zhang 
286963029d29SHong Zhang   /* (c) convert local indices to global indices in parallel vector Vleaves */
28709566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Vleaves_seq,&N));
28719566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vleaves_seq,&varr_read));
287263029d29SHong Zhang   for (i=0; i<N; i+=2) {
28732e4cff2eSHong Zhang     remoterank = (PetscMPIInt)PetscRealPart(varr_read[i]);
287463029d29SHong Zhang     if (remoterank == rank) {
287563029d29SHong Zhang       k = i+1; /* row number */
28762e4cff2eSHong Zhang       lidx = (PetscInt)PetscRealPart(varr_read[i+1]);
287763029d29SHong Zhang       val  = (PetscScalar)vltog[lidx]; /* global index for non-ghost vertex computed above */
28789566063dSJacob Faibussowitsch       PetscCall(VecSetValues(Vleaves,1,&k,&val,INSERT_VALUES));
287963029d29SHong Zhang     }
288063029d29SHong Zhang   }
28819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vleaves_seq,&varr_read));
28829566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Vleaves));
28839566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Vleaves));
288463029d29SHong Zhang 
288563029d29SHong Zhang   /* (d) Set vltog for ghost vertices by copying local values of Vleaves */
28869566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vleaves,&varr_read));
288763029d29SHong Zhang   k = 0;
288863029d29SHong Zhang   for (i=0; i<nroots; i++) {
2889*daad07d3SAidan Hamilton     PetscCall(DMNetworkIsGhostVertex(dm,i+network->cloneshared->vStart,&ghost));
289063029d29SHong Zhang     if (!ghost) continue;
28912e4cff2eSHong Zhang     vltog[i] = (PetscInt)PetscRealPart(varr_read[2*k+1]); k++;
289263029d29SHong Zhang   }
28939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vleaves,&varr_read));
289463029d29SHong Zhang 
28959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Vleaves));
28969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Vleaves_seq));
28979566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&ctx));
289822bbedd7SHong Zhang   PetscFunctionReturn(0);
289922bbedd7SHong Zhang }
290042dc13f1SHong Zhang 
29019fbee547SJacob Faibussowitsch static inline PetscErrorCode DMISAddSize_private(DM_Network *network,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *nidx)
290242dc13f1SHong Zhang {
290342dc13f1SHong Zhang   PetscInt                 i,j,ncomps,nvar,key,offset=0;
290442dc13f1SHong Zhang   DMNetworkComponentHeader header;
290542dc13f1SHong Zhang 
290642dc13f1SHong Zhang   PetscFunctionBegin;
29079566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offset));
290842dc13f1SHong Zhang   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
290942dc13f1SHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offset);
291042dc13f1SHong Zhang 
291142dc13f1SHong Zhang   for (i=0; i<ncomps; i++) {
291242dc13f1SHong Zhang     key  = header->key[i];
291342dc13f1SHong Zhang     nvar = header->nvar[i];
291442dc13f1SHong Zhang     for (j=0; j<numkeys; j++) {
291542dc13f1SHong Zhang       if (key == keys[j]) {
291642dc13f1SHong Zhang         if (!blocksize || blocksize[j] == -1) {
291742dc13f1SHong Zhang           *nidx += nvar;
291842dc13f1SHong Zhang         } else {
291942dc13f1SHong Zhang           *nidx += nselectedvar[j]*nvar/blocksize[j];
292042dc13f1SHong Zhang         }
292142dc13f1SHong Zhang       }
292242dc13f1SHong Zhang     }
292342dc13f1SHong Zhang   }
292442dc13f1SHong Zhang   PetscFunctionReturn(0);
292542dc13f1SHong Zhang }
292642dc13f1SHong Zhang 
29279fbee547SJacob Faibussowitsch static inline PetscErrorCode DMISComputeIdx_private(DM dm,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],PetscInt *ii,PetscInt *idx)
292842dc13f1SHong Zhang {
292942dc13f1SHong Zhang   PetscInt                 i,j,ncomps,nvar,key,offsetg,k,k1,offset=0;
293042dc13f1SHong Zhang   DM_Network               *network = (DM_Network*)dm->data;
293142dc13f1SHong Zhang   DMNetworkComponentHeader header;
293242dc13f1SHong Zhang 
293342dc13f1SHong Zhang   PetscFunctionBegin;
29349566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offset));
293542dc13f1SHong Zhang   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
293642dc13f1SHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offset);
293742dc13f1SHong Zhang 
293842dc13f1SHong Zhang   for (i=0; i<ncomps; i++) {
293942dc13f1SHong Zhang     key  = header->key[i];
294042dc13f1SHong Zhang     nvar = header->nvar[i];
294142dc13f1SHong Zhang     for (j=0; j<numkeys; j++) {
294242dc13f1SHong Zhang       if (key != keys[j]) continue;
294342dc13f1SHong Zhang 
29449566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetGlobalVecOffset(dm,p,i,&offsetg));
294542dc13f1SHong Zhang       if (!blocksize || blocksize[j] == -1) {
294642dc13f1SHong Zhang         for (k=0; k<nvar; k++) idx[(*ii)++] = offsetg + k;
294742dc13f1SHong Zhang       } else {
294842dc13f1SHong Zhang         for (k=0; k<nvar; k+=blocksize[j]) {
294942dc13f1SHong Zhang           for (k1=0; k1<nselectedvar[j]; k1++) idx[(*ii)++] = offsetg + k + selectedvar[j][k1];
295042dc13f1SHong Zhang         }
295142dc13f1SHong Zhang       }
295242dc13f1SHong Zhang     }
295342dc13f1SHong Zhang   }
295442dc13f1SHong Zhang   PetscFunctionReturn(0);
295542dc13f1SHong Zhang }
295642dc13f1SHong Zhang 
295742dc13f1SHong Zhang /*@
295842dc13f1SHong Zhang   DMNetworkCreateIS - Create an index set object from the global vector of the network
295942dc13f1SHong Zhang 
296042dc13f1SHong Zhang   Collective
296142dc13f1SHong Zhang 
296242dc13f1SHong Zhang   Input Parameters:
296342dc13f1SHong Zhang + dm - DMNetwork object
296442dc13f1SHong Zhang . numkeys - number of keys
296542dc13f1SHong Zhang . keys - array of keys that define the components of the variables you wish to extract
296642dc13f1SHong Zhang . blocksize - block size of the variables associated to the component
296742dc13f1SHong Zhang . nselectedvar - number of variables in each block to select
296842dc13f1SHong Zhang - selectedvar - the offset into the block of each variable in each block to select
296942dc13f1SHong Zhang 
297042dc13f1SHong Zhang   Output Parameters:
297142dc13f1SHong Zhang . is - the index set
297242dc13f1SHong Zhang 
297342dc13f1SHong Zhang   Level: Advanced
297442dc13f1SHong Zhang 
297542dc13f1SHong Zhang   Notes:
297642dc13f1SHong Zhang     Use blocksize[i] of -1 to indicate select all the variables of the i-th component, for which nselectvar[i] and selectedvar[i] are ignored. Use NULL, NULL, NULL to indicate for all selected components one wishes to obtain all the values of that component. For example, DMNetworkCreateIS(dm,1,&key,NULL,NULL,NULL,&is) will return an is that extracts all the variables for the 0-th component.
297742dc13f1SHong Zhang 
2978db781477SPatrick Sanan .seealso: `DMNetworkCreate()`, `ISCreateGeneral()`, `DMNetworkCreateLocalIS()`
297942dc13f1SHong Zhang @*/
298042dc13f1SHong Zhang PetscErrorCode DMNetworkCreateIS(DM dm,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],IS *is)
298142dc13f1SHong Zhang {
298242dc13f1SHong Zhang   MPI_Comm       comm;
298342dc13f1SHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
298442dc13f1SHong Zhang   PetscInt       i,p,estart,eend,vstart,vend,nidx,*idx;
298542dc13f1SHong Zhang   PetscBool      ghost;
298642dc13f1SHong Zhang 
298742dc13f1SHong Zhang   PetscFunctionBegin;
29889566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
298942dc13f1SHong Zhang 
299042dc13f1SHong Zhang   /* Check input parameters */
299142dc13f1SHong Zhang   for (i=0; i<numkeys; i++) {
299242dc13f1SHong Zhang     if (!blocksize || blocksize[i] == -1) continue;
299363a3b9bcSJacob Faibussowitsch     PetscCheck(nselectedvar[i] <= blocksize[i],PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"number of selectedvariables %" PetscInt_FMT " cannot be larger than blocksize %" PetscInt_FMT,nselectedvar[i],blocksize[i]);
299442dc13f1SHong Zhang   }
299542dc13f1SHong Zhang 
29969566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(dm,&estart,&eend));
29979566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(dm,&vstart,&vend));
299842dc13f1SHong Zhang 
299942dc13f1SHong Zhang   /* Get local number of idx */
300042dc13f1SHong Zhang   nidx = 0;
300142dc13f1SHong Zhang   for (p=estart; p<eend; p++) {
30029566063dSJacob Faibussowitsch     PetscCall(DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx));
300342dc13f1SHong Zhang   }
300442dc13f1SHong Zhang   for (p=vstart; p<vend; p++) {
30059566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(dm,p,&ghost));
300642dc13f1SHong Zhang     if (ghost) continue;
30079566063dSJacob Faibussowitsch     PetscCall(DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx));
300842dc13f1SHong Zhang   }
300942dc13f1SHong Zhang 
301042dc13f1SHong Zhang   /* Compute idx */
30119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nidx,&idx));
301242dc13f1SHong Zhang   i = 0;
301342dc13f1SHong Zhang   for (p=estart; p<eend; p++) {
30149566063dSJacob Faibussowitsch     PetscCall(DMISComputeIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx));
301542dc13f1SHong Zhang   }
301642dc13f1SHong Zhang   for (p=vstart; p<vend; p++) {
30179566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(dm,p,&ghost));
301842dc13f1SHong Zhang     if (ghost) continue;
30199566063dSJacob Faibussowitsch     PetscCall(DMISComputeIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx));
302042dc13f1SHong Zhang   }
302142dc13f1SHong Zhang 
302242dc13f1SHong Zhang   /* Create is */
30239566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm,nidx,idx,PETSC_COPY_VALUES,is));
30249566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
302542dc13f1SHong Zhang   PetscFunctionReturn(0);
302642dc13f1SHong Zhang }
302742dc13f1SHong Zhang 
30289fbee547SJacob Faibussowitsch static inline PetscErrorCode DMISComputeLocalIdx_private(DM dm,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],PetscInt *ii,PetscInt *idx)
302942dc13f1SHong Zhang {
303042dc13f1SHong Zhang   PetscInt                 i,j,ncomps,nvar,key,offsetl,k,k1,offset=0;
303142dc13f1SHong Zhang   DM_Network               *network = (DM_Network*)dm->data;
303242dc13f1SHong Zhang   DMNetworkComponentHeader header;
303342dc13f1SHong Zhang 
303442dc13f1SHong Zhang   PetscFunctionBegin;
30359566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offset));
303642dc13f1SHong Zhang   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
303742dc13f1SHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offset);
303842dc13f1SHong Zhang 
303942dc13f1SHong Zhang   for (i=0; i<ncomps; i++) {
304042dc13f1SHong Zhang     key  = header->key[i];
304142dc13f1SHong Zhang     nvar = header->nvar[i];
304242dc13f1SHong Zhang     for (j=0; j<numkeys; j++) {
304342dc13f1SHong Zhang       if (key != keys[j]) continue;
304442dc13f1SHong Zhang 
30459566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetLocalVecOffset(dm,p,i,&offsetl));
304642dc13f1SHong Zhang       if (!blocksize || blocksize[j] == -1) {
304742dc13f1SHong Zhang         for (k=0; k<nvar; k++) idx[(*ii)++] = offsetl + k;
304842dc13f1SHong Zhang       } else {
304942dc13f1SHong Zhang         for (k=0; k<nvar; k+=blocksize[j]) {
305042dc13f1SHong Zhang           for (k1=0; k1<nselectedvar[j]; k1++) idx[(*ii)++] = offsetl + k + selectedvar[j][k1];
305142dc13f1SHong Zhang         }
305242dc13f1SHong Zhang       }
305342dc13f1SHong Zhang     }
305442dc13f1SHong Zhang   }
305542dc13f1SHong Zhang   PetscFunctionReturn(0);
305642dc13f1SHong Zhang }
305742dc13f1SHong Zhang 
305842dc13f1SHong Zhang /*@
305942dc13f1SHong Zhang   DMNetworkCreateLocalIS - Create an index set object from the local vector of the network
306042dc13f1SHong Zhang 
306142dc13f1SHong Zhang   Not collective
306242dc13f1SHong Zhang 
306342dc13f1SHong Zhang   Input Parameters:
306442dc13f1SHong Zhang + dm - DMNetwork object
306542dc13f1SHong Zhang . numkeys - number of keys
306642dc13f1SHong Zhang . keys - array of keys that define the components of the variables you wish to extract
306742dc13f1SHong Zhang . blocksize - block size of the variables associated to the component
306842dc13f1SHong Zhang . nselectedvar - number of variables in each block to select
306942dc13f1SHong Zhang - selectedvar - the offset into the block of each variable in each block to select
307042dc13f1SHong Zhang 
307142dc13f1SHong Zhang   Output Parameters:
307242dc13f1SHong Zhang . is - the index set
307342dc13f1SHong Zhang 
307442dc13f1SHong Zhang   Level: Advanced
307542dc13f1SHong Zhang 
307642dc13f1SHong Zhang   Notes:
307742dc13f1SHong Zhang     Use blocksize[i] of -1 to indicate select all the variables of the i-th component, for which nselectvar[i] and selectedvar[i] are ignored. Use NULL, NULL, NULL to indicate for all selected components one wishes to obtain all the values of that component. For example, DMNetworkCreateLocalIS(dm,1,&key,NULL,NULL,NULL,&is) will return an is that extracts all the variables for the 0-th component.
307842dc13f1SHong Zhang 
3079db781477SPatrick Sanan .seealso: `DMNetworkCreate()`, `DMNetworkCreateIS`, `ISCreateGeneral()`
308042dc13f1SHong Zhang @*/
308142dc13f1SHong Zhang PetscErrorCode DMNetworkCreateLocalIS(DM dm,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],IS *is)
308242dc13f1SHong Zhang {
308342dc13f1SHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
308442dc13f1SHong Zhang   PetscInt       i,p,pstart,pend,nidx,*idx;
308542dc13f1SHong Zhang 
308642dc13f1SHong Zhang   PetscFunctionBegin;
308742dc13f1SHong Zhang   /* Check input parameters */
308842dc13f1SHong Zhang   for (i=0; i<numkeys; i++) {
308942dc13f1SHong Zhang     if (!blocksize || blocksize[i] == -1) continue;
309063a3b9bcSJacob Faibussowitsch     PetscCheck(nselectedvar[i] <= blocksize[i],PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"number of selectedvariables %" PetscInt_FMT " cannot be larger than blocksize %" PetscInt_FMT,nselectedvar[i],blocksize[i]);
309142dc13f1SHong Zhang   }
309242dc13f1SHong Zhang 
3093*daad07d3SAidan Hamilton   pstart = network->cloneshared->pStart;
3094*daad07d3SAidan Hamilton   pend   = network->cloneshared->pEnd;
309542dc13f1SHong Zhang 
309642dc13f1SHong Zhang   /* Get local number of idx */
309742dc13f1SHong Zhang   nidx = 0;
309842dc13f1SHong Zhang   for (p=pstart; p<pend; p++) {
30999566063dSJacob Faibussowitsch     PetscCall(DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx));
310042dc13f1SHong Zhang   }
310142dc13f1SHong Zhang 
310242dc13f1SHong Zhang   /* Compute local idx */
31039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nidx,&idx));
310442dc13f1SHong Zhang   i = 0;
310542dc13f1SHong Zhang   for (p=pstart; p<pend; p++) {
31069566063dSJacob Faibussowitsch     PetscCall(DMISComputeLocalIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx));
310742dc13f1SHong Zhang   }
310842dc13f1SHong Zhang 
310942dc13f1SHong Zhang   /* Create is */
31109566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,is));
31119566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
311242dc13f1SHong Zhang   PetscFunctionReturn(0);
311342dc13f1SHong Zhang }
3114*daad07d3SAidan Hamilton /*@
3115*daad07d3SAidan Hamilton   DMNetworkFinalizeComponents - Sets up interal data structures for the sections and components. It is called after registering new components and adding all components
3116*daad07d3SAidan Hamilton   to the cloned network. After calling this subroutine, no new components can be added to the network.
3117*daad07d3SAidan Hamilton 
3118*daad07d3SAidan Hamilton   Collective
3119*daad07d3SAidan Hamilton 
3120*daad07d3SAidan Hamilton   Input Parameters:
3121*daad07d3SAidan Hamilton . dm - the dmnetwork object
3122*daad07d3SAidan Hamilton 
3123*daad07d3SAidan Hamilton   Level: beginner
3124*daad07d3SAidan Hamilton 
3125*daad07d3SAidan Hamilton .seealso: `DMNetworkAddComponent()`, `DMNetworkRegisterComponent()`, `DMSetUp()`
3126*daad07d3SAidan Hamilton @*/
3127*daad07d3SAidan Hamilton PetscErrorCode DMNetworkFinalizeComponents(DM dm)
3128*daad07d3SAidan Hamilton {
3129*daad07d3SAidan Hamilton   DM_Network     *network = (DM_Network*)dm->data;
3130*daad07d3SAidan Hamilton 
3131*daad07d3SAidan Hamilton   PetscFunctionBegin;
3132*daad07d3SAidan Hamilton   if (network->componentsetup) PetscFunctionReturn(0);
3133*daad07d3SAidan Hamilton   PetscCall(DMNetworkComponentSetUp(dm));
3134*daad07d3SAidan Hamilton   PetscCall(DMNetworkVariablesSetUp(dm));
3135*daad07d3SAidan Hamilton   PetscCall(DMSetLocalSection(network->plex,network->DofSection));
3136*daad07d3SAidan Hamilton   PetscCall(DMGetGlobalSection(network->plex,&network->GlobalDofSection));
3137*daad07d3SAidan Hamilton   network->componentsetup = PETSC_TRUE;
3138*daad07d3SAidan Hamilton   PetscFunctionReturn(0);
3139*daad07d3SAidan Hamilton }
3140