xref: /petsc/src/dm/impls/network/network.c (revision 60225df5d8469840be2bf9c1f64795a92b19f3c2)
1af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h> /*I  "petscdmnetwork.h"  I*/
25f25b224SDuncan Campbell #include "petscis.h"
35f2c45f1SShri Abhyankar 
4e3b1a053SGetnet Betrie PetscLogEvent DMNetwork_LayoutSetUp;
5e3b1a053SGetnet Betrie PetscLogEvent DMNetwork_SetUpNetwork;
6e3b1a053SGetnet Betrie PetscLogEvent DMNetwork_Distribute;
7e3b1a053SGetnet Betrie 
854dfd506SHong Zhang /*
954dfd506SHong Zhang   Creates the component header and value objects for a network point
1054dfd506SHong Zhang */
11d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetUpNetworkHeaderComponentValue(DM dm, DMNetworkComponentHeader header, DMNetworkComponentValue cvalue)
12d71ae5a4SJacob Faibussowitsch {
1354dfd506SHong Zhang   PetscFunctionBegin;
1454dfd506SHong Zhang   /* Allocate arrays for component information */
159566063dSJacob Faibussowitsch   PetscCall(PetscCalloc5(header->maxcomps, &header->size, header->maxcomps, &header->key, header->maxcomps, &header->offset, header->maxcomps, &header->nvar, header->maxcomps, &header->offsetvarrel));
169566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(header->maxcomps, &cvalue->data));
1754dfd506SHong Zhang 
1854dfd506SHong 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.
1954dfd506SHong Zhang    If the data header struct changes then this header size calculation needs to be updated. */
2054dfd506SHong Zhang   header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5 * header->maxcomps * sizeof(PetscInt);
2154dfd506SHong Zhang   header->hsize /= sizeof(DMNetworkComponentGenericDataType);
22dd5bbc93SStefano Zampini #if defined(__NEC__)
23dd5bbc93SStefano Zampini   /* NEC/LG: quick hack to keep data aligned on 8 bytes. */
24dd5bbc93SStefano Zampini   header->hsize = (header->hsize + (8 - 1)) & ~(8 - 1);
25dd5bbc93SStefano Zampini #endif
263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2754dfd506SHong Zhang }
2854dfd506SHong Zhang 
29d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkInitializeHeaderComponentData(DM dm)
30d71ae5a4SJacob Faibussowitsch {
31daad07d3SAidan Hamilton   DM_Network *network = (DM_Network *)dm->data;
32daad07d3SAidan Hamilton   PetscInt    np, p, defaultnumcomp = DMNETWORK_MAX_COMP_AT_POINT_DEFAULT;
33daad07d3SAidan Hamilton 
34daad07d3SAidan Hamilton   PetscFunctionBegin;
35daad07d3SAidan Hamilton   np = network->cloneshared->pEnd - network->cloneshared->pStart;
36daad07d3SAidan Hamilton   if (network->header)
37daad07d3SAidan Hamilton     for (p = 0; p < np; p++) {
38daad07d3SAidan Hamilton       network->header[p].maxcomps = defaultnumcomp;
39daad07d3SAidan Hamilton       PetscCall(SetUpNetworkHeaderComponentValue(dm, &network->header[p], &network->cvalue[p]));
40daad07d3SAidan Hamilton     }
41daad07d3SAidan Hamilton 
423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43daad07d3SAidan Hamilton }
445f25b224SDuncan Campbell 
455f2c45f1SShri Abhyankar /*@
4620f4b53cSBarry Smith   DMNetworkGetPlex - Gets the `DMPLEX` associated with this `DMNETWORK`
47556ed216SShri Abhyankar 
4820f4b53cSBarry Smith   Not Collective
49556ed216SShri Abhyankar 
502fe279fdSBarry Smith   Input Parameter:
5120f4b53cSBarry Smith . dm - the `DMNETWORK` object
522bf73ac6SHong Zhang 
532fe279fdSBarry Smith   Output Parameter:
5420f4b53cSBarry Smith . plexdm - the `DMPLEX` object
55556ed216SShri Abhyankar 
5620f4b53cSBarry Smith   Level: advanced
57556ed216SShri Abhyankar 
5820f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMPLEX`, `DMNetworkCreate()`
59556ed216SShri Abhyankar @*/
60d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetPlex(DM dm, DM *plexdm)
61d71ae5a4SJacob Faibussowitsch {
622bf73ac6SHong Zhang   DM_Network *network = (DM_Network *)dm->data;
63556ed216SShri Abhyankar 
64556ed216SShri Abhyankar   PetscFunctionBegin;
65556ed216SShri Abhyankar   *plexdm = network->plex;
663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
67556ed216SShri Abhyankar }
68556ed216SShri Abhyankar 
69556ed216SShri Abhyankar /*@
702bf73ac6SHong Zhang   DMNetworkGetNumSubNetworks - Gets the the number of subnetworks
7172c3e9f7SHong Zhang 
7220f4b53cSBarry Smith   Not Collective
7372c3e9f7SHong Zhang 
74f899ff85SJose E. Roman   Input Parameter:
7520f4b53cSBarry Smith . dm - the `DMNETWORK` object
762bf73ac6SHong Zhang 
772bf73ac6SHong Zhang   Output Parameters:
782bf73ac6SHong Zhang + nsubnet - local number of subnetworks
792bf73ac6SHong Zhang - Nsubnet - global number of subnetworks
8072c3e9f7SHong Zhang 
8197bb938eSShri Abhyankar   Level: beginner
8272c3e9f7SHong Zhang 
8320f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkCreate()`, `DMNetworkSetNumSubNetworks()`
8472c3e9f7SHong Zhang @*/
85d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetNumSubNetworks(DM dm, PetscInt *nsubnet, PetscInt *Nsubnet)
86d71ae5a4SJacob Faibussowitsch {
872bf73ac6SHong Zhang   DM_Network *network = (DM_Network *)dm->data;
8872c3e9f7SHong Zhang 
8972c3e9f7SHong Zhang   PetscFunctionBegin;
90daad07d3SAidan Hamilton   if (nsubnet) *nsubnet = network->cloneshared->nsubnet;
91daad07d3SAidan Hamilton   if (Nsubnet) *Nsubnet = network->cloneshared->Nsubnet;
923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9372c3e9f7SHong Zhang }
9472c3e9f7SHong Zhang 
9572c3e9f7SHong Zhang /*@
962bf73ac6SHong Zhang   DMNetworkSetNumSubNetworks - Sets the number of subnetworks
975f2c45f1SShri Abhyankar 
9820f4b53cSBarry Smith   Collective
995f2c45f1SShri Abhyankar 
1005f2c45f1SShri Abhyankar   Input Parameters:
10120f4b53cSBarry Smith + dm      - the `DMNETWORK` object
1022bf73ac6SHong Zhang . nsubnet - local number of subnetworks
1032bf73ac6SHong Zhang - Nsubnet - global number of subnetworks
1045f2c45f1SShri Abhyankar 
10597bb938eSShri Abhyankar   Level: beginner
1061b266c99SBarry Smith 
10720f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkCreate()`, `DMNetworkGetNumSubNetworks()`
1085f2c45f1SShri Abhyankar @*/
109d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkSetNumSubNetworks(DM dm, PetscInt nsubnet, PetscInt Nsubnet)
110d71ae5a4SJacob Faibussowitsch {
1115f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
1125f2c45f1SShri Abhyankar 
1135f2c45f1SShri Abhyankar   PetscFunctionBegin;
114d5b43468SJose E. Roman   PetscCheck(network->cloneshared->Nsubnet == 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "Network sizes already set, cannot resize the network");
1152bf73ac6SHong Zhang 
1165f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1172bf73ac6SHong Zhang   PetscValidLogicalCollectiveInt(dm, nsubnet, 2);
1182bf73ac6SHong Zhang   PetscValidLogicalCollectiveInt(dm, Nsubnet, 3);
1197765340cSHong Zhang 
1202bf73ac6SHong Zhang   if (Nsubnet == PETSC_DECIDE) {
1215c6496baSHong Zhang     PetscCheck(nsubnet >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of local subnetworks %" PetscInt_FMT " cannot be less than 0", nsubnet);
1221c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&nsubnet, &Nsubnet, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
123caf410d2SHong Zhang   }
1245c6496baSHong Zhang   PetscCheck(Nsubnet >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "Number of global subnetworks %" PetscInt_FMT " cannot be less than 1", Nsubnet);
125caf410d2SHong Zhang 
126daad07d3SAidan Hamilton   network->cloneshared->Nsubnet = Nsubnet;
127da81f932SPierre Jolivet   network->cloneshared->nsubnet = 0; /* initial value; will be determined by DMNetworkAddSubnetwork() */
128daad07d3SAidan Hamilton   PetscCall(PetscCalloc1(Nsubnet, &network->cloneshared->subnet));
129caf410d2SHong Zhang 
1302bf73ac6SHong Zhang   /* num of shared vertices */
131daad07d3SAidan Hamilton   network->cloneshared->nsvtx = 0;
132daad07d3SAidan Hamilton   network->cloneshared->Nsvtx = 0;
1333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1347765340cSHong Zhang }
1357765340cSHong Zhang 
136*60225df5SJacob Faibussowitsch /*@C
1372bf73ac6SHong Zhang   DMNetworkAddSubnetwork - Add a subnetwork
1385f2c45f1SShri Abhyankar 
13920f4b53cSBarry Smith   Collective
1405f2c45f1SShri Abhyankar 
1415f2c45f1SShri Abhyankar   Input Parameters:
14220f4b53cSBarry Smith + dm       - the `DMNETWORK`  object
1432bf73ac6SHong Zhang . name     - name of the subnetwork
1442bf73ac6SHong Zhang . ne       - number of local edges of this subnetwork
14520f4b53cSBarry 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
14620f4b53cSBarry Smith               of the vertices) of each edge,
1478d1cd658SBarry Smith $            [first vertex of first edge, second vertex of first edge, first vertex of second edge, second vertex of second edge, etc]
1482bf73ac6SHong Zhang 
1492fe279fdSBarry Smith   Output Parameter:
1502bf73ac6SHong Zhang . netnum - global index of the subnetwork
1515f2c45f1SShri Abhyankar 
15220f4b53cSBarry Smith   Level: beginner
15320f4b53cSBarry Smith 
1545f2c45f1SShri Abhyankar   Notes:
1555f2c45f1SShri Abhyankar   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
15620f4b53cSBarry Smith   not be destroyed before the call to `DMNetworkLayoutSetUp()`
1575f2c45f1SShri Abhyankar 
15820f4b53cSBarry 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.
15920f4b53cSBarry Smith   For a multiple subnetworks,
1608d1cd658SBarry 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.
161f11a936eSBarry Smith 
162e3e68989SHong Zhang   Example usage:
163f11a936eSBarry Smith   Consider the following networks:
164d5b43468SJose E. Roman   1) A single subnetwork:
165e3e68989SHong Zhang .vb
166f11a936eSBarry Smith  network 0:
167f11a936eSBarry Smith  rank[0]:
168f11a936eSBarry Smith    v0 --> v2; v1 --> v2
169f11a936eSBarry Smith  rank[1]:
170f11a936eSBarry Smith    v3 --> v5; v4 --> v5
171e3e68989SHong Zhang .ve
172e3e68989SHong Zhang 
173*60225df5SJacob Faibussowitsch   The resulting input network 0\:
174*60225df5SJacob Faibussowitsch .vb
175f11a936eSBarry Smith   rank[0]:
176f11a936eSBarry Smith   ne = 2
177f11a936eSBarry Smith   edgelist = [0 2 | 1 2]
178*60225df5SJacob Faibussowitsch 
179f11a936eSBarry Smith   rank[1]:
180f11a936eSBarry Smith   ne = 2
181f11a936eSBarry Smith   edgelist = [3 5 | 4 5]
182*60225df5SJacob Faibussowitsch .ve
183*60225df5SJacob Faibussowitsch   2) Two subnetworks\:
184f11a936eSBarry Smith .vb
185f11a936eSBarry Smith  subnetwork 0:
186f11a936eSBarry Smith  rank[0]:
187f11a936eSBarry Smith    v0 --> v2; v2 --> v1; v1 --> v3;
188f11a936eSBarry Smith  subnetwork 1:
189f11a936eSBarry Smith  rank[1]:
190f11a936eSBarry Smith    v0 --> v3; v3 --> v2; v2 --> v1;
191f11a936eSBarry Smith .ve
192f11a936eSBarry Smith 
193*60225df5SJacob Faibussowitsch   The resulting input subnetwork 0\:
194*60225df5SJacob Faibussowitsch .vb
195f11a936eSBarry Smith   rank[0]:
196f11a936eSBarry Smith   ne = 3
197f11a936eSBarry Smith   edgelist = [0 2 | 2 1 | 1 3]
198*60225df5SJacob Faibussowitsch 
199f11a936eSBarry Smith   rank[1]:
200f11a936eSBarry Smith   ne = 0
201f11a936eSBarry Smith   edgelist = NULL
202*60225df5SJacob Faibussowitsch .ve
203*60225df5SJacob Faibussowitsch   subnetwork 1\:
204*60225df5SJacob Faibussowitsch .vb
205f11a936eSBarry Smith   rank[0]:
206f11a936eSBarry Smith   ne = 0
207f11a936eSBarry Smith   edgelist = NULL
208*60225df5SJacob Faibussowitsch 
209f11a936eSBarry Smith   rank[1]:
210f11a936eSBarry Smith   edgelist = [0 3 | 3 2 | 2 1]
211*60225df5SJacob Faibussowitsch .ve
212e3e68989SHong Zhang 
21320f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkCreate()`, `DMNetworkSetNumSubnetworks()`
2145f2c45f1SShri Abhyankar @*/
215d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkAddSubnetwork(DM dm, const char *name, PetscInt ne, PetscInt edgelist[], PetscInt *netnum)
216d71ae5a4SJacob Faibussowitsch {
2175f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
2185c6496baSHong Zhang   PetscInt    i, Nedge, j, Nvtx, nvtx, nvtx_min = -1, nvtx_max = 0;
219f11a936eSBarry Smith   PetscBT     table;
2205f2c45f1SShri Abhyankar 
2215f2c45f1SShri Abhyankar   PetscFunctionBegin;
222ad540459SPierre Jolivet   for (i = 0; i < ne; i++) 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]);
223f11a936eSBarry Smith 
2245c6496baSHong Zhang   i = 0;
2255c6496baSHong Zhang   if (ne) nvtx_min = nvtx_max = edgelist[0];
2265c6496baSHong Zhang   for (j = 0; j < ne; j++) {
2275c6496baSHong Zhang     nvtx_min = PetscMin(nvtx_min, edgelist[i]);
2285c6496baSHong Zhang     nvtx_max = PetscMax(nvtx_max, edgelist[i]);
2295c6496baSHong Zhang     i++;
2305c6496baSHong Zhang     nvtx_min = PetscMin(nvtx_min, edgelist[i]);
2315c6496baSHong Zhang     nvtx_max = PetscMax(nvtx_max, edgelist[i]);
2325c6496baSHong Zhang     i++;
2335c6496baSHong Zhang   }
2345c6496baSHong Zhang   Nvtx = nvtx_max - nvtx_min + 1; /* approximated total local nvtx for this subnet */
2355c6496baSHong Zhang 
2365c6496baSHong Zhang   /* Get exact local nvtx for this subnet: counting local values between nvtx_min and nvtx_max */
2379566063dSJacob Faibussowitsch   PetscCall(PetscBTCreate(Nvtx, &table));
2389566063dSJacob Faibussowitsch   PetscCall(PetscBTMemzero(Nvtx, table));
239f11a936eSBarry Smith   i = 0;
240f11a936eSBarry Smith   for (j = 0; j < ne; j++) {
2419566063dSJacob Faibussowitsch     PetscCall(PetscBTSet(table, edgelist[i++] - nvtx_min));
2429566063dSJacob Faibussowitsch     PetscCall(PetscBTSet(table, edgelist[i++] - nvtx_min));
243f11a936eSBarry Smith   }
244f11a936eSBarry Smith   nvtx = 0;
245f11a936eSBarry Smith   for (j = 0; j < Nvtx; j++) {
246f11a936eSBarry Smith     if (PetscBTLookup(table, j)) nvtx++;
247f11a936eSBarry Smith   }
2485c6496baSHong Zhang 
2495c6496baSHong Zhang   /* Get global total Nvtx = max(edgelist[])+1 for this subnet */
2501c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&nvtx_max, &Nvtx, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
2515c6496baSHong Zhang   Nvtx++;
2529566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&table));
253f11a936eSBarry Smith 
254f11a936eSBarry Smith   /* Get global total Nedge for this subnet */
2551c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&ne, &Nedge, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
256f11a936eSBarry Smith 
257daad07d3SAidan Hamilton   i = network->cloneshared->nsubnet;
258c6a7a370SJeremy L Thompson   if (name) PetscCall(PetscStrncpy(network->cloneshared->subnet[i].name, name, sizeof(network->cloneshared->subnet[i].name)));
259daad07d3SAidan Hamilton   network->cloneshared->subnet[i].nvtx     = nvtx; /* include ghost vertices */
260daad07d3SAidan Hamilton   network->cloneshared->subnet[i].nedge    = ne;
261daad07d3SAidan Hamilton   network->cloneshared->subnet[i].edgelist = edgelist;
262daad07d3SAidan Hamilton   network->cloneshared->subnet[i].Nvtx     = Nvtx;
263daad07d3SAidan Hamilton   network->cloneshared->subnet[i].Nedge    = Nedge;
2642bf73ac6SHong Zhang 
2652bf73ac6SHong Zhang   /* ----------------------------------------------------------
2662bf73ac6SHong Zhang    p=v or e;
2672bf73ac6SHong Zhang    subnet[0].pStart   = 0
2682bf73ac6SHong Zhang    subnet[i+1].pStart = subnet[i].pEnd = subnet[i].pStart + (nE[i] or NV[i])
2692bf73ac6SHong Zhang    ----------------------------------------------------------------------- */
2702bf73ac6SHong Zhang   /* GLOBAL subnet[].vStart and vEnd, used by DMNetworkLayoutSetUp() */
271daad07d3SAidan Hamilton   network->cloneshared->subnet[i].vStart = network->cloneshared->NVertices;
272daad07d3SAidan Hamilton   network->cloneshared->subnet[i].vEnd   = network->cloneshared->subnet[i].vStart + network->cloneshared->subnet[i].Nvtx; /* global vEnd of subnet[i] */
2732bf73ac6SHong Zhang 
274daad07d3SAidan Hamilton   network->cloneshared->nVertices += nvtx; /* include ghost vertices */
275daad07d3SAidan Hamilton   network->cloneshared->NVertices += network->cloneshared->subnet[i].Nvtx;
2762bf73ac6SHong Zhang 
2772bf73ac6SHong Zhang   /* LOCAL subnet[].eStart and eEnd, used by DMNetworkLayoutSetUp() */
278daad07d3SAidan Hamilton   network->cloneshared->subnet[i].eStart = network->cloneshared->nEdges;
279daad07d3SAidan Hamilton   network->cloneshared->subnet[i].eEnd   = network->cloneshared->subnet[i].eStart + ne;
280daad07d3SAidan Hamilton   network->cloneshared->nEdges += ne;
281daad07d3SAidan Hamilton   network->cloneshared->NEdges += network->cloneshared->subnet[i].Nedge;
2822bf73ac6SHong Zhang 
283c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(network->cloneshared->subnet[i].name, name, sizeof(network->cloneshared->subnet[i].name)));
284daad07d3SAidan Hamilton   if (netnum) *netnum = network->cloneshared->nsubnet;
285daad07d3SAidan Hamilton   network->cloneshared->nsubnet++;
2863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2872bf73ac6SHong Zhang }
2882bf73ac6SHong Zhang 
2895c6496baSHong Zhang /*@C
2905c6496baSHong Zhang   DMNetworkSharedVertexGetInfo - Get info of a shared vertex struct, see petsc/private/dmnetworkimpl.h
2912bf73ac6SHong Zhang 
29220f4b53cSBarry Smith   Not Collective
2935c6496baSHong Zhang 
2945c6496baSHong Zhang   Input Parameters:
29520f4b53cSBarry Smith + dm - the `DM` object
2965c6496baSHong Zhang - v  - vertex point
2975c6496baSHong Zhang 
2985c6496baSHong Zhang   Output Parameters:
2995c6496baSHong Zhang + gidx - global number of this shared vertex in the internal dmplex
3005c6496baSHong Zhang . n    - number of subnetworks that share this vertex
3015c6496baSHong Zhang - sv   - array of size n: sv[2*i,2*i+1]=(net[i], idx[i]), i=0,...,n-1
3025c6496baSHong Zhang 
3035c6496baSHong Zhang   Level: intermediate
3045c6496baSHong Zhang 
30520f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetSharedVertices()`
3065c6496baSHong Zhang @*/
307d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkSharedVertexGetInfo(DM dm, PetscInt v, PetscInt *gidx, PetscInt *n, const PetscInt **sv)
308d71ae5a4SJacob Faibussowitsch {
3095c6496baSHong Zhang   DM_Network *network = (DM_Network *)dm->data;
310daad07d3SAidan Hamilton   SVtx       *svtx    = network->cloneshared->svtx;
3115c6496baSHong Zhang   PetscInt    i, gidx_tmp;
3125c6496baSHong Zhang 
3135c6496baSHong Zhang   PetscFunctionBegin;
3149566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetGlobalVertexIndex(dm, v, &gidx_tmp));
315eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapIGetWithDefault(network->cloneshared->svtable, gidx_tmp + 1, 0, &i));
3165c6496baSHong Zhang   PetscCheck(i > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "input vertex is not a shared vertex");
3175c6496baSHong Zhang 
3185c6496baSHong Zhang   i--;
3195c6496baSHong Zhang   if (gidx) *gidx = gidx_tmp;
3205c6496baSHong Zhang   if (n) *n = svtx[i].n;
3215c6496baSHong Zhang   if (sv) *sv = svtx[i].sv;
3223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3235c6496baSHong Zhang }
3245c6496baSHong Zhang 
3255c6496baSHong Zhang /*
3265c6496baSHong Zhang   VtxGetInfo - Get info of an input vertex=(net,idx)
3275c6496baSHong Zhang 
3285c6496baSHong Zhang   Input Parameters:
3295c6496baSHong Zhang + Nsvtx - global num of shared vertices
3305c6496baSHong Zhang . svtx - array of shared vertices (global)
3315c6496baSHong Zhang - (net,idx) - subnet number and local index for a vertex
3325c6496baSHong Zhang 
3335c6496baSHong Zhang   Output Parameters:
3345c6496baSHong Zhang + gidx - global index of (net,idx)
3355c6496baSHong Zhang . svtype - see petsc/private/dmnetworkimpl.h
3365c6496baSHong Zhang - svtx_idx - ordering in the svtx array
3375c6496baSHong Zhang */
338d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode VtxGetInfo(PetscInt Nsvtx, SVtx *svtx, PetscInt net, PetscInt idx, PetscInt *gidx, SVtxType *svtype, PetscInt *svtx_idx)
339d71ae5a4SJacob Faibussowitsch {
3405c6496baSHong Zhang   PetscInt i, j, *svto, g_idx;
3412bf73ac6SHong Zhang   SVtxType vtype;
3422bf73ac6SHong Zhang 
3432bf73ac6SHong Zhang   PetscFunctionBegin;
3443ba16761SJacob Faibussowitsch   if (!Nsvtx) PetscFunctionReturn(PETSC_SUCCESS);
3452bf73ac6SHong Zhang 
3465c6496baSHong Zhang   g_idx = -1;
3472bf73ac6SHong Zhang   vtype = SVNONE;
3485c6496baSHong Zhang 
3492bf73ac6SHong Zhang   for (i = 0; i < Nsvtx; i++) {
3502bf73ac6SHong Zhang     if (net == svtx[i].sv[0] && idx == svtx[i].sv[1]) {
3515c6496baSHong Zhang       g_idx = svtx[i].gidx;
3522bf73ac6SHong Zhang       vtype = SVFROM;
3532bf73ac6SHong Zhang     } else { /* loop over svtx[i].n */
3542bf73ac6SHong Zhang       for (j = 1; j < svtx[i].n; j++) {
3552bf73ac6SHong Zhang         svto = svtx[i].sv + 2 * j;
3562bf73ac6SHong Zhang         if (net == svto[0] && idx == svto[1]) {
3572bf73ac6SHong Zhang           /* input vertex net.idx is a shared to_vertex, output its global index and its svtype */
3585c6496baSHong Zhang           g_idx = svtx[i].gidx; /* output gidx for to_vertex */
3592bf73ac6SHong Zhang           vtype = SVTO;
3602bf73ac6SHong Zhang         }
3612bf73ac6SHong Zhang       }
3622bf73ac6SHong Zhang     }
3632bf73ac6SHong Zhang     if (vtype != SVNONE) break;
3642bf73ac6SHong Zhang   }
3655c6496baSHong Zhang   if (gidx) *gidx = g_idx;
3662bf73ac6SHong Zhang   if (svtype) *svtype = vtype;
3672bf73ac6SHong Zhang   if (svtx_idx) *svtx_idx = i;
3683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3692bf73ac6SHong Zhang }
3702bf73ac6SHong Zhang 
3712bf73ac6SHong Zhang /*
3725c6496baSHong Zhang   TableAddSVtx - Add a new shared vertice from sedgelist[k] to a ctable svta
3735c6496baSHong Zhang 
3745c6496baSHong Zhang   Input:  network, sedgelist, k, svta
3755c6496baSHong Zhang   Output: svta, tdata, ta2sv
3762bf73ac6SHong Zhang */
377eec179cfSJacob Faibussowitsch static inline PetscErrorCode TableAddSVtx(DM_Network *network, PetscInt *sedgelist, PetscInt k, PetscHMapI svta, PetscInt *tdata, PetscInt *ta2sv)
378d71ae5a4SJacob Faibussowitsch {
3795c6496baSHong Zhang   PetscInt net, idx, gidx;
3802bf73ac6SHong Zhang 
3812bf73ac6SHong Zhang   PetscFunctionBegin;
3825c6496baSHong Zhang   net  = sedgelist[k];
3835c6496baSHong Zhang   idx  = sedgelist[k + 1];
384daad07d3SAidan Hamilton   gidx = network->cloneshared->subnet[net].vStart + idx;
385c76ffc5fSJacob Faibussowitsch   PetscCall(PetscHMapISet(svta, gidx + 1, *tdata + 1));
3865c6496baSHong Zhang 
3875c6496baSHong Zhang   ta2sv[*tdata] = k; /* maps tdata to index of sedgelist */
3885c6496baSHong Zhang   (*tdata)++;
3893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3902bf73ac6SHong Zhang }
3912bf73ac6SHong Zhang 
3922bf73ac6SHong Zhang /*
3935c6496baSHong Zhang   SharedVtxCreate - Create an array of global shared vertices. See SVtx and SVtxType in dmnetworkimpl.h
3942bf73ac6SHong Zhang 
3952bf73ac6SHong Zhang   Input:  dm, Nsedgelist, sedgelist
3962bf73ac6SHong Zhang 
3972bf73ac6SHong Zhang   Note: Output svtx is organized as
3982bf73ac6SHong Zhang         sv(net[0],idx[0]) --> sv(net[1],idx[1])
3992bf73ac6SHong Zhang                           --> sv(net[1],idx[1])
4002bf73ac6SHong Zhang                           ...
4012bf73ac6SHong Zhang                           --> sv(net[n-1],idx[n-1])
4022bf73ac6SHong Zhang         and net[0] < net[1] < ... < net[n-1]
4032bf73ac6SHong Zhang         where sv[0] has SVFROM type, sv[i], i>0, has SVTO type.
4042bf73ac6SHong Zhang  */
405d71ae5a4SJacob Faibussowitsch static PetscErrorCode SharedVtxCreate(DM dm, PetscInt Nsedgelist, PetscInt *sedgelist)
406d71ae5a4SJacob Faibussowitsch {
4075c6496baSHong Zhang   SVtx         *svtx = NULL;
4082bf73ac6SHong Zhang   PetscInt     *sv, k, j, nsv, *tdata, **ta2sv;
409eec179cfSJacob Faibussowitsch   PetscHMapI   *svtas;
4107340ca2cSHong Zhang   PetscInt      gidx, net, idx, i, nta, ita, idx_from, idx_to, n, *net_tmp, *idx_tmp, *gidx_tmp;
4112bf73ac6SHong Zhang   DM_Network   *network = (DM_Network *)dm->data;
412eec179cfSJacob Faibussowitsch   PetscHashIter ppos;
4132bf73ac6SHong Zhang 
4142bf73ac6SHong Zhang   PetscFunctionBegin;
4155c6496baSHong Zhang   /* (1) Crete an array of ctables svtas to map (net,idx) -> gidx; a svtas[] for a shared/merged vertex */
4169566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(Nsedgelist, &svtas, Nsedgelist, &tdata, 2 * Nsedgelist, &ta2sv));
4172bf73ac6SHong Zhang 
4182bf73ac6SHong Zhang   k   = 0; /* sedgelist vertex counter j = 4*k */
4197340ca2cSHong Zhang   nta = 0; /* num of svta tables created = num of groups of shared vertices */
4202bf73ac6SHong Zhang 
4212bf73ac6SHong Zhang   /* for j=0 */
422eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapICreateWithSize(2 * Nsedgelist, svtas + nta));
4239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * Nsedgelist, &ta2sv[nta]));
4242bf73ac6SHong Zhang 
4259566063dSJacob Faibussowitsch   PetscCall(TableAddSVtx(network, sedgelist, k, svtas[nta], &tdata[nta], ta2sv[nta]));
4269566063dSJacob Faibussowitsch   PetscCall(TableAddSVtx(network, sedgelist, k + 2, svtas[nta], &tdata[nta], ta2sv[nta]));
4279371c9d4SSatish Balay   nta++;
4289371c9d4SSatish Balay   k += 4;
4292bf73ac6SHong Zhang 
4305c6496baSHong Zhang   for (j = 1; j < Nsedgelist; j++) { /* j: sedgelist counter */
4312bf73ac6SHong Zhang     for (ita = 0; ita < nta; ita++) {
4322bf73ac6SHong Zhang       /* vfrom */
4339371c9d4SSatish Balay       net  = sedgelist[k];
4349371c9d4SSatish Balay       idx  = sedgelist[k + 1];
435daad07d3SAidan Hamilton       gidx = network->cloneshared->subnet[net].vStart + idx; /* global index of the vertex net.idx before merging shared vertices */
436eec179cfSJacob Faibussowitsch       PetscCall(PetscHMapIGetWithDefault(svtas[ita], gidx + 1, 0, &idx_from));
4372bf73ac6SHong Zhang 
4382bf73ac6SHong Zhang       /* vto */
4399371c9d4SSatish Balay       net  = sedgelist[k + 2];
4409371c9d4SSatish Balay       idx  = sedgelist[k + 3];
441daad07d3SAidan Hamilton       gidx = network->cloneshared->subnet[net].vStart + idx;
442eec179cfSJacob Faibussowitsch       PetscCall(PetscHMapIGetWithDefault(svtas[ita], gidx + 1, 0, &idx_to));
4432bf73ac6SHong Zhang 
4442bf73ac6SHong Zhang       if (idx_from || idx_to) { /* vfrom or vto is on table svtas[ita] */
4459371c9d4SSatish Balay         idx_from--;
4469371c9d4SSatish Balay         idx_to--;
4472bf73ac6SHong Zhang         if (idx_from < 0) { /* vto is on svtas[ita] */
4489566063dSJacob Faibussowitsch           PetscCall(TableAddSVtx(network, sedgelist, k, svtas[ita], &tdata[ita], ta2sv[ita]));
4492bf73ac6SHong Zhang           break;
4502bf73ac6SHong Zhang         } else if (idx_to < 0) {
4519566063dSJacob Faibussowitsch           PetscCall(TableAddSVtx(network, sedgelist, k + 2, svtas[ita], &tdata[ita], ta2sv[ita]));
4522bf73ac6SHong Zhang           break;
4532bf73ac6SHong Zhang         }
4542bf73ac6SHong Zhang       }
4552bf73ac6SHong Zhang     }
4562bf73ac6SHong Zhang 
4572bf73ac6SHong Zhang     if (ita == nta) {
458eec179cfSJacob Faibussowitsch       PetscCall(PetscHMapICreateWithSize(2 * Nsedgelist, svtas + nta));
4599566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2 * Nsedgelist, &ta2sv[nta]));
4602bf73ac6SHong Zhang 
4619566063dSJacob Faibussowitsch       PetscCall(TableAddSVtx(network, sedgelist, k, svtas[nta], &tdata[nta], ta2sv[nta]));
4629566063dSJacob Faibussowitsch       PetscCall(TableAddSVtx(network, sedgelist, k + 2, svtas[nta], &tdata[nta], ta2sv[nta]));
4632bf73ac6SHong Zhang       nta++;
4642bf73ac6SHong Zhang     }
4652bf73ac6SHong Zhang     k += 4;
4662bf73ac6SHong Zhang   }
4672bf73ac6SHong Zhang 
4686aad120cSJose E. Roman   /* (2) Create svtable for query shared vertices using gidx */
469eec179cfSJacob Faibussowitsch   PetscCall(PetscHMapICreateWithSize(nta, &network->cloneshared->svtable));
4705c6496baSHong Zhang 
4715c6496baSHong Zhang   /* (3) Construct svtx from svtas
4727340ca2cSHong Zhang    svtx: array of SVtx: sv[0]=(net[0],idx[0]) to vertices sv[k], k=1,...,n-1. */
4739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nta, &svtx));
4742bf73ac6SHong Zhang   for (nsv = 0; nsv < nta; nsv++) {
4754f5c2772SJose E. Roman     /* for a single svtx, put shared vertices in ascending order of gidx */
476eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIGetSize(svtas[nsv], &n));
4779566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(2 * n, &sv));
4787340ca2cSHong Zhang     PetscCall(PetscMalloc3(n, &gidx_tmp, n, &net_tmp, n, &idx_tmp));
4795c6496baSHong Zhang     svtx[nsv].sv   = sv;
4805c6496baSHong Zhang     svtx[nsv].n    = n;
481daad07d3SAidan Hamilton     svtx[nsv].gidx = network->cloneshared->NVertices; /* initialization */
4822bf73ac6SHong Zhang 
483eec179cfSJacob Faibussowitsch     PetscHashIterBegin(svtas[nsv], ppos);
4844f5c2772SJose E. Roman     for (k = 0; k < n; k++) { /* gidx is sorted in ascending order */
485eec179cfSJacob Faibussowitsch       PetscHashIterGetKey(svtas[nsv], ppos, gidx);
486eec179cfSJacob Faibussowitsch       PetscHashIterGetVal(svtas[nsv], ppos, i);
487eec179cfSJacob Faibussowitsch       PetscHashIterNext(svtas[nsv], ppos);
4889371c9d4SSatish Balay       gidx--;
4899371c9d4SSatish Balay       i--;
4905c6496baSHong Zhang       j           = ta2sv[nsv][i];    /* maps i to index of sedgelist */
4917340ca2cSHong Zhang       net_tmp[k]  = sedgelist[j];     /* subnet number */
4927340ca2cSHong Zhang       idx_tmp[k]  = sedgelist[j + 1]; /* index on the subnet */
4937340ca2cSHong Zhang       gidx_tmp[k] = gidx;             /* gidx in un-merged dmnetwork */
4942bf73ac6SHong Zhang     }
4955c6496baSHong Zhang 
4967340ca2cSHong Zhang     /* current implementation requires sv[]=[net,idx] in ascending order of its gidx in un-merged dmnetwork */
4977340ca2cSHong Zhang     PetscCall(PetscSortIntWithArrayPair(n, gidx_tmp, net_tmp, idx_tmp));
4987340ca2cSHong Zhang     svtx[nsv].gidx = gidx_tmp[0]; /* = min(gidx) */
4997340ca2cSHong Zhang     for (k = 0; k < n; k++) {
5007340ca2cSHong Zhang       sv[2 * k]     = net_tmp[k];
5017340ca2cSHong Zhang       sv[2 * k + 1] = idx_tmp[k];
5027340ca2cSHong Zhang     }
5037340ca2cSHong Zhang     PetscCall(PetscFree3(gidx_tmp, net_tmp, idx_tmp));
5047340ca2cSHong Zhang 
5056aad120cSJose E. Roman     /* Setup svtable for query shared vertices */
506c76ffc5fSJacob Faibussowitsch     PetscCall(PetscHMapISet(network->cloneshared->svtable, svtx[nsv].gidx + 1, nsv + 1));
5072bf73ac6SHong Zhang   }
5082bf73ac6SHong Zhang 
5092bf73ac6SHong Zhang   for (j = 0; j < nta; j++) {
510eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIDestroy(svtas + j));
5119566063dSJacob Faibussowitsch     PetscCall(PetscFree(ta2sv[j]));
5122bf73ac6SHong Zhang   }
5139566063dSJacob Faibussowitsch   PetscCall(PetscFree3(svtas, tdata, ta2sv));
5142bf73ac6SHong Zhang 
515daad07d3SAidan Hamilton   network->cloneshared->Nsvtx = nta;
516daad07d3SAidan Hamilton   network->cloneshared->svtx  = svtx;
5173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5182bf73ac6SHong Zhang }
5192bf73ac6SHong Zhang 
5205c6496baSHong Zhang /*
5215c6496baSHong Zhang   GetEdgelist_Coupling - Get an integrated edgelist for dmplex from user-provided subnet[].edgelist when subnets are coupled by shared vertices
5225c6496baSHong Zhang 
5235c6496baSHong Zhang   Input Parameters:
5245c6496baSHong Zhang . dm - the dmnetwork object
5255c6496baSHong Zhang 
5265c6496baSHong Zhang    Output Parameters:
5275c6496baSHong Zhang +  edges - the integrated edgelist for dmplex
5285c6496baSHong Zhang -  nmerged_ptr - num of vertices being merged
5295c6496baSHong Zhang */
530d71ae5a4SJacob Faibussowitsch static PetscErrorCode GetEdgelist_Coupling(DM dm, PetscInt *edges, PetscInt *nmerged_ptr)
531d71ae5a4SJacob Faibussowitsch {
5322bf73ac6SHong Zhang   MPI_Comm    comm;
5335c795782SAidan Hamilton   PetscMPIInt size, rank;
534eac198afSGetnet   DM_Network *network = (DM_Network *)dm->data;
535eac198afSGetnet   PetscInt    i, j, ctr, np;
536daad07d3SAidan Hamilton   PetscInt   *vidxlTog, Nsv, Nsubnet = network->cloneshared->Nsubnet;
5375c795782SAidan Hamilton   PetscInt   *sedgelist = network->cloneshared->sedgelist, vrange;
5385c795782SAidan Hamilton   PetscInt    net, idx, gidx, nmerged, gidx_from, net_from, sv_idx;
5392bf73ac6SHong Zhang   SVtxType    svtype = SVNONE;
5405c6496baSHong Zhang   SVtx       *svtx;
5412bf73ac6SHong Zhang 
5422bf73ac6SHong Zhang   PetscFunctionBegin;
5439566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
5462bf73ac6SHong Zhang 
5475c6496baSHong Zhang   /* (1) Create global svtx[] from sedgelist */
5485c6496baSHong Zhang   /* --------------------------------------- */
549daad07d3SAidan Hamilton   PetscCall(SharedVtxCreate(dm, network->cloneshared->Nsvtx, sedgelist));
550daad07d3SAidan Hamilton   Nsv  = network->cloneshared->Nsvtx;
551daad07d3SAidan Hamilton   svtx = network->cloneshared->svtx;
5522bf73ac6SHong Zhang 
553d5b43468SJose E. Roman   /* (2) Merge shared vto vertices to their vfrom vertex with same global vertex index (gidx) */
5545c6496baSHong Zhang   /* --------------------------------------------------------------------------------------- */
5552bf73ac6SHong Zhang   /* (2.1) compute vrage[rank]: global index of 1st local vertex in proc[rank] */
5565c795782SAidan Hamilton   PetscCall(PetscMalloc1(network->cloneshared->nVertices, &vidxlTog));
5572bf73ac6SHong Zhang 
5585c795782SAidan Hamilton   PetscCallMPI(MPI_Scan(&network->cloneshared->nVertices, &vrange, 1, MPIU_INT, MPI_SUM, comm));
5595c795782SAidan Hamilton   vrange -= network->cloneshared->nVertices;
5602bf73ac6SHong Zhang 
5612bf73ac6SHong Zhang   /* (2.2) Create vidxlTog: maps UN-MERGED local vertex index i to global index gidx (plex, excluding ghost vertices) */
5629371c9d4SSatish Balay   i                           = 0;
5639371c9d4SSatish Balay   gidx                        = 0;
5642bf73ac6SHong Zhang   nmerged                     = 0; /* local num of merged vertices */
565daad07d3SAidan Hamilton   network->cloneshared->nsvtx = 0; /* local num of SVtx structs, including ghosts */
5662bf73ac6SHong Zhang   for (net = 0; net < Nsubnet; net++) {
567daad07d3SAidan Hamilton     for (idx = 0; idx < network->cloneshared->subnet[net].Nvtx; idx++) { /* Note: global subnet[net].Nvtx */
5689566063dSJacob Faibussowitsch       PetscCall(VtxGetInfo(Nsv, svtx, net, idx, &gidx_from, &svtype, &sv_idx));
5692bf73ac6SHong Zhang       if (svtype == SVTO) {
570daad07d3SAidan Hamilton         if (network->cloneshared->subnet[net].nvtx) { /* this proc owns sv_to */
5715c6496baSHong Zhang           net_from = svtx[sv_idx].sv[0];              /* subnet number of its shared vertex */
572daad07d3SAidan Hamilton           if (network->cloneshared->subnet[net_from].nvtx == 0) {
5735c6496baSHong Zhang             /* this proc does not own v_from, thus a ghost local vertex */
574daad07d3SAidan Hamilton             network->cloneshared->nsvtx++;
5752bf73ac6SHong Zhang           }
5765c6496baSHong Zhang           vidxlTog[i++] = gidx_from; /* gidx before merging! Bug??? */
5775c6496baSHong Zhang           nmerged++;                 /* a shared vertex -- merged */
5782bf73ac6SHong Zhang         }
5792bf73ac6SHong Zhang       } else {
580daad07d3SAidan Hamilton         if (svtype == SVFROM && network->cloneshared->subnet[net].nvtx) {
5815c6496baSHong Zhang           /* this proc owns this v_from, a new local shared vertex */
582daad07d3SAidan Hamilton           network->cloneshared->nsvtx++;
5832bf73ac6SHong Zhang         }
584daad07d3SAidan Hamilton         if (network->cloneshared->subnet[net].nvtx) vidxlTog[i++] = gidx;
5852bf73ac6SHong Zhang         gidx++;
5862bf73ac6SHong Zhang       }
5872bf73ac6SHong Zhang     }
5882bf73ac6SHong Zhang   }
5892bf73ac6SHong Zhang #if defined(PETSC_USE_DEBUG)
590daad07d3SAidan Hamilton   PetscCheck(i == network->cloneshared->nVertices, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "%" PetscInt_FMT " != %" PetscInt_FMT " nVertices", i, network->cloneshared->nVertices);
5912bf73ac6SHong Zhang #endif
5922bf73ac6SHong Zhang 
5935c6496baSHong Zhang   /* (2.3) Shared vertices in the subnetworks are merged, update global NVertices: np = sum(local nmerged) */
594712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(&nmerged, &np, 1, MPIU_INT, MPI_SUM, comm));
595daad07d3SAidan Hamilton   network->cloneshared->NVertices -= np;
5962bf73ac6SHong Zhang 
5972bf73ac6SHong Zhang   ctr = 0;
5982bf73ac6SHong Zhang   for (net = 0; net < Nsubnet; net++) {
599daad07d3SAidan Hamilton     for (j = 0; j < network->cloneshared->subnet[net].nedge; j++) {
6002bf73ac6SHong Zhang       /* vfrom: */
6015c795782SAidan Hamilton       i              = network->cloneshared->subnet[net].edgelist[2 * j] + (network->cloneshared->subnet[net].vStart - vrange);
6022bf73ac6SHong Zhang       edges[2 * ctr] = vidxlTog[i];
6032bf73ac6SHong Zhang 
6042bf73ac6SHong Zhang       /* vto */
6055c795782SAidan Hamilton       i                  = network->cloneshared->subnet[net].edgelist[2 * j + 1] + (network->cloneshared->subnet[net].vStart - vrange);
6062bf73ac6SHong Zhang       edges[2 * ctr + 1] = vidxlTog[i];
6072bf73ac6SHong Zhang       ctr++;
6082bf73ac6SHong Zhang     }
6092bf73ac6SHong Zhang   }
6105c795782SAidan Hamilton   PetscCall(PetscFree(vidxlTog));
6119566063dSJacob Faibussowitsch   PetscCall(PetscFree(sedgelist)); /* created in DMNetworkAddSharedVertices() */
6122bf73ac6SHong Zhang 
613eac198afSGetnet   *nmerged_ptr = nmerged;
6143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6155f2c45f1SShri Abhyankar }
6165f2c45f1SShri Abhyankar 
617d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkInitializeNonTopological(DM dm)
618d71ae5a4SJacob Faibussowitsch {
619daad07d3SAidan Hamilton   DM_Network *network = (DM_Network *)dm->data;
620daad07d3SAidan Hamilton   PetscInt    p, pStart = network->cloneshared->pStart, pEnd = network->cloneshared->pEnd;
621daad07d3SAidan Hamilton   MPI_Comm    comm;
622daad07d3SAidan Hamilton 
623daad07d3SAidan Hamilton   PetscFunctionBegin;
624daad07d3SAidan Hamilton   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
625daad07d3SAidan Hamilton 
626daad07d3SAidan Hamilton   PetscCall(PetscSectionCreate(comm, &network->DataSection));
627daad07d3SAidan Hamilton   PetscCall(PetscSectionCreate(comm, &network->DofSection));
628daad07d3SAidan Hamilton   PetscCall(PetscSectionSetChart(network->DataSection, pStart, pEnd));
629daad07d3SAidan Hamilton   PetscCall(PetscSectionSetChart(network->DofSection, pStart, pEnd));
630daad07d3SAidan Hamilton 
631daad07d3SAidan Hamilton   PetscCall(DMNetworkInitializeHeaderComponentData(dm));
632daad07d3SAidan Hamilton 
633daad07d3SAidan Hamilton   for (p = 0; p < pEnd - pStart; p++) {
634daad07d3SAidan Hamilton     network->header[p].ndata           = 0;
635daad07d3SAidan Hamilton     network->header[p].offset[0]       = 0;
636daad07d3SAidan Hamilton     network->header[p].offsetvarrel[0] = 0;
637daad07d3SAidan Hamilton     PetscCall(PetscSectionAddDof(network->DataSection, p, network->header[p].hsize));
638daad07d3SAidan Hamilton   }
6393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
640daad07d3SAidan Hamilton }
641daad07d3SAidan Hamilton 
6425f2c45f1SShri Abhyankar /*@
6435f2c45f1SShri Abhyankar   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
6445f2c45f1SShri Abhyankar 
645eac198afSGetnet   Not Collective
6465f2c45f1SShri Abhyankar 
6472fe279fdSBarry Smith   Input Parameter:
64820f4b53cSBarry Smith . dm - the `DMNETWORK` object
64920f4b53cSBarry Smith 
65020f4b53cSBarry Smith   Level: beginner
6515f2c45f1SShri Abhyankar 
6525f2c45f1SShri Abhyankar   Notes:
6535f2c45f1SShri Abhyankar   This routine should be called after the network sizes and edgelists have been provided. It creates
6545f2c45f1SShri Abhyankar   the bare layout of the network and sets up the network to begin insertion of components.
6555f2c45f1SShri Abhyankar 
6565f2c45f1SShri Abhyankar   All the components should be registered before calling this routine.
6575f2c45f1SShri Abhyankar 
65820f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkSetNumSubNetworks()`, `DMNetworkAddSubnetwork()`
6595f2c45f1SShri Abhyankar @*/
660d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkLayoutSetUp(DM dm)
661d71ae5a4SJacob Faibussowitsch {
6625f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network *)dm->data;
663daad07d3SAidan Hamilton   PetscInt        i, j, ctr, Nsubnet = network->cloneshared->Nsubnet, np, *edges, *subnetvtx, *subnetedge, e, v, vfrom, vto, net, globaledgeoff;
664caf410d2SHong Zhang   const PetscInt *cone;
665caf410d2SHong Zhang   MPI_Comm        comm;
6665503a911SAidan Hamilton   PetscMPIInt     size;
667eac198afSGetnet   PetscSection    sectiong;
6685c6496baSHong Zhang   PetscInt        nmerged = 0;
6696500d4abSHong Zhang 
6706500d4abSHong Zhang   PetscFunctionBegin;
671e3b1a053SGetnet Betrie   PetscCall(PetscLogEventBegin(DMNetwork_LayoutSetUp, dm, 0, 0, 0));
672daad07d3SAidan Hamilton   PetscCheck(network->cloneshared->nsubnet == Nsubnet, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Must call DMNetworkAddSubnetwork() %" PetscInt_FMT " times", Nsubnet);
6732bf73ac6SHong Zhang 
674eac198afSGetnet   /* This implementation requires user input each subnet by a single processor when Nsubnet>1, thus subnet[net].nvtx=subnet[net].Nvtx when net>0 */
675eac198afSGetnet   for (net = 1; net < Nsubnet; net++) {
6769371c9d4SSatish Balay     if (network->cloneshared->subnet[net].nvtx)
6779371c9d4SSatish Balay       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,
6789371c9d4SSatish Balay                  network->cloneshared->subnet[net].nvtx, network->cloneshared->subnet[net].Nvtx);
679eac198afSGetnet   }
680eac198afSGetnet 
6819566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
6823ba16761SJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
6836500d4abSHong Zhang 
684f11a936eSBarry Smith   /* Create LOCAL edgelist in global vertex ordering for the network by concatenating local input edgelists of the subnetworks */
685daad07d3SAidan Hamilton   PetscCall(PetscCalloc1(2 * network->cloneshared->nEdges, &edges));
686eac198afSGetnet 
687daad07d3SAidan Hamilton   if (network->cloneshared->Nsvtx) { /* subnetworks are coupled via shared vertices */
6889566063dSJacob Faibussowitsch     PetscCall(GetEdgelist_Coupling(dm, edges, &nmerged));
689eac198afSGetnet   } else { /* subnetworks are not coupled */
6906aad120cSJose E. Roman     /* Create a 0-size svtable for query shared vertices */
691eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&network->cloneshared->svtable));
692caf410d2SHong Zhang     ctr = 0;
6932bf73ac6SHong Zhang     for (i = 0; i < Nsubnet; i++) {
694daad07d3SAidan Hamilton       for (j = 0; j < network->cloneshared->subnet[i].nedge; j++) {
695daad07d3SAidan Hamilton         edges[2 * ctr]     = network->cloneshared->subnet[i].vStart + network->cloneshared->subnet[i].edgelist[2 * j];
696daad07d3SAidan Hamilton         edges[2 * ctr + 1] = network->cloneshared->subnet[i].vStart + network->cloneshared->subnet[i].edgelist[2 * j + 1];
6976500d4abSHong Zhang         ctr++;
6986500d4abSHong Zhang       }
6996500d4abSHong Zhang     }
700eac198afSGetnet   }
7017765340cSHong Zhang 
7022bf73ac6SHong Zhang   /* Create network->plex; One dimensional network, numCorners=2 */
7039566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, &network->plex));
7049566063dSJacob Faibussowitsch   PetscCall(DMSetType(network->plex, DMPLEX));
7059566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(network->plex, 1));
706eac198afSGetnet 
707daad07d3SAidan Hamilton   if (size == 1) PetscCall(DMPlexBuildFromCellList(network->plex, network->cloneshared->nEdges, PETSC_DECIDE, 2, edges));
708daad07d3SAidan Hamilton   else PetscCall(DMPlexBuildFromCellListParallel(network->plex, network->cloneshared->nEdges, PETSC_DECIDE, PETSC_DECIDE, 2, edges, NULL, NULL));
7096500d4abSHong Zhang 
710daad07d3SAidan Hamilton   PetscCall(DMPlexGetChart(network->plex, &network->cloneshared->pStart, &network->cloneshared->pEnd));
711daad07d3SAidan Hamilton   PetscCall(DMPlexGetHeightStratum(network->plex, 0, &network->cloneshared->eStart, &network->cloneshared->eEnd));
712daad07d3SAidan Hamilton   PetscCall(DMPlexGetHeightStratum(network->plex, 1, &network->cloneshared->vStart, &network->cloneshared->vEnd));
713daad07d3SAidan Hamilton   np = network->cloneshared->pEnd - network->cloneshared->pStart;
7149566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(np, &network->header, np, &network->cvalue));
715caf410d2SHong Zhang 
716f11a936eSBarry Smith   /* Create edge and vertex arrays for the subnetworks
717f11a936eSBarry Smith      This implementation assumes that DMNetwork reads
718f11a936eSBarry Smith      (1) a single subnetwork in parallel; or
719f11a936eSBarry Smith      (2) n subnetworks using n processors, one subnetwork/processor.
720f11a936eSBarry Smith   */
721daad07d3SAidan Hamilton   PetscCall(PetscCalloc2(network->cloneshared->nEdges, &subnetedge, network->cloneshared->nVertices + network->cloneshared->nsvtx, &subnetvtx)); /* Maps local edge/vertex to local subnetwork's edge/vertex */
722daad07d3SAidan Hamilton   network->cloneshared->subnetedge = subnetedge;
723daad07d3SAidan Hamilton   network->cloneshared->subnetvtx  = subnetvtx;
7245c6496baSHong Zhang   for (j = 0; j < Nsubnet; j++) {
725daad07d3SAidan Hamilton     network->cloneshared->subnet[j].edges = subnetedge;
726daad07d3SAidan Hamilton     subnetedge += network->cloneshared->subnet[j].nedge;
727f11a936eSBarry Smith 
728daad07d3SAidan Hamilton     network->cloneshared->subnet[j].vertices = subnetvtx;
729daad07d3SAidan Hamilton     subnetvtx += network->cloneshared->subnet[j].nvtx;
7306500d4abSHong Zhang   }
731daad07d3SAidan Hamilton   network->cloneshared->svertices = subnetvtx;
732caf410d2SHong Zhang 
733caf410d2SHong Zhang   /* Get edge ownership */
734daad07d3SAidan Hamilton   np = network->cloneshared->eEnd - network->cloneshared->eStart;
7355503a911SAidan Hamilton   PetscCallMPI(MPI_Scan(&np, &globaledgeoff, 1, MPIU_INT, MPI_SUM, comm));
7365503a911SAidan Hamilton   globaledgeoff -= np;
737caf410d2SHong Zhang 
738eac198afSGetnet   /* Setup local edge and vertex arrays for subnetworks */
7392bf73ac6SHong Zhang   e = 0;
7402bf73ac6SHong Zhang   for (i = 0; i < Nsubnet; i++) {
741f11a936eSBarry Smith     ctr = 0;
742daad07d3SAidan Hamilton     for (j = 0; j < network->cloneshared->subnet[i].nedge; j++) {
7432bf73ac6SHong Zhang       /* edge e */
7445503a911SAidan Hamilton       network->header[e].index                 = e + globaledgeoff; /* Global edge index */
7452bf73ac6SHong Zhang       network->header[e].subnetid              = i;
746daad07d3SAidan Hamilton       network->cloneshared->subnet[i].edges[j] = e;
7472bf73ac6SHong Zhang 
7482bf73ac6SHong Zhang       /* connected vertices */
7499566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(network->plex, e, &cone));
7502bf73ac6SHong Zhang 
7512bf73ac6SHong Zhang       /* vertex cone[0] */
752f11a936eSBarry Smith       v                           = cone[0];
753f11a936eSBarry Smith       network->header[v].index    = edges[2 * e]; /* Global vertex index */
754f11a936eSBarry Smith       network->header[v].subnetid = i;            /* Subnetwork id */
755f11a936eSBarry Smith       if (Nsubnet == 1) {
756daad07d3SAidan Hamilton         network->cloneshared->subnet[i].vertices[v - network->cloneshared->vStart] = v; /* user's subnet[].idx = petsc's v */
757f11a936eSBarry Smith       } else {
758daad07d3SAidan Hamilton         vfrom                                           = network->cloneshared->subnet[i].edgelist[2 * ctr]; /* =subnet[i].idx, Global index! */
759daad07d3SAidan Hamilton         network->cloneshared->subnet[i].vertices[vfrom] = v;                                                 /* user's subnet[].dix = petsc's v */
760f11a936eSBarry Smith       }
7612bf73ac6SHong Zhang 
7622bf73ac6SHong Zhang       /* vertex cone[1] */
763f11a936eSBarry Smith       v                           = cone[1];
764f11a936eSBarry Smith       network->header[v].index    = edges[2 * e + 1]; /* Global vertex index */
765f11a936eSBarry Smith       network->header[v].subnetid = i;                /* Subnetwork id */
766f11a936eSBarry Smith       if (Nsubnet == 1) {
767daad07d3SAidan Hamilton         network->cloneshared->subnet[i].vertices[v - network->cloneshared->vStart] = v; /* user's subnet[].idx = petsc's v */
768f11a936eSBarry Smith       } else {
769daad07d3SAidan Hamilton         vto                                           = network->cloneshared->subnet[i].edgelist[2 * ctr + 1]; /* =subnet[i].idx, Global index! */
770daad07d3SAidan Hamilton         network->cloneshared->subnet[i].vertices[vto] = v;                                                     /* user's subnet[].dix = petsc's v */
771f11a936eSBarry Smith       }
7722bf73ac6SHong Zhang 
7739371c9d4SSatish Balay       e++;
7749371c9d4SSatish Balay       ctr++;
7756500d4abSHong Zhang     }
7766500d4abSHong Zhang   }
7775503a911SAidan Hamilton   PetscCall(PetscFree(edges));
778caf410d2SHong Zhang 
779eac198afSGetnet   /* Set local vertex array for the subnetworks */
780eac198afSGetnet   j = 0;
781daad07d3SAidan Hamilton   for (v = network->cloneshared->vStart; v < network->cloneshared->vEnd; v++) {
782eac198afSGetnet     /* local shared vertex */
783eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIGetWithDefault(network->cloneshared->svtable, network->header[v].index + 1, 0, &i));
784daad07d3SAidan Hamilton     if (i) network->cloneshared->svertices[j++] = v;
7856500d4abSHong Zhang   }
786eac198afSGetnet 
787eac198afSGetnet   /* Create a global section to be used by DMNetworkIsGhostVertex() which is a non-collective routine */
788eac198afSGetnet   /* see snes_tutorials_network-ex1_4 */
7899566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(network->plex, &sectiong));
790daad07d3SAidan Hamilton   /* Initialize non-topological data structures  */
791daad07d3SAidan Hamilton   PetscCall(DMNetworkInitializeNonTopological(dm));
792e3b1a053SGetnet Betrie   PetscCall(PetscLogEventEnd(DMNetwork_LayoutSetUp, dm, 0, 0, 0));
7933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7945f2c45f1SShri Abhyankar }
7955f2c45f1SShri Abhyankar 
79694ef8ddeSSatish Balay /*@C
7972bf73ac6SHong Zhang   DMNetworkGetSubnetwork - Returns the information about a requested subnetwork
7982bf73ac6SHong Zhang 
79920f4b53cSBarry Smith   Not Collective
8002727e31bSShri Abhyankar 
8017a7aea1fSJed Brown   Input Parameters:
80220f4b53cSBarry Smith + dm     - the `DMNETWORK` object
8032bf73ac6SHong Zhang - netnum - the global index of the subnetwork
8042727e31bSShri Abhyankar 
8057a7aea1fSJed Brown   Output Parameters:
8062727e31bSShri Abhyankar + nv   - number of vertices (local)
8072727e31bSShri Abhyankar . ne   - number of edges (local)
8082bf73ac6SHong Zhang . vtx  - local vertices of the subnetwork
8092bf73ac6SHong Zhang - edge - local edges of the subnetwork
8102727e31bSShri Abhyankar 
81106dd6b0eSSatish Balay   Level: intermediate
81206dd6b0eSSatish Balay 
81320f4b53cSBarry Smith   Notes:
81420f4b53cSBarry Smith   Cannot call this routine before `DMNetworkLayoutSetup()`
81520f4b53cSBarry Smith 
81620f4b53cSBarry Smith   The local vertices returned on each rank are determined by `DMNETWORK`. The user does not have any control over what vertices are local.
81720f4b53cSBarry Smith 
81820f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkCreate()`, `DMNetworkAddSubnetwork()`, `DMNetworkLayoutSetUp()`
8192727e31bSShri Abhyankar @*/
820d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetSubnetwork(DM dm, PetscInt netnum, PetscInt *nv, PetscInt *ne, const PetscInt **vtx, const PetscInt **edge)
821d71ae5a4SJacob Faibussowitsch {
822caf410d2SHong Zhang   DM_Network *network = (DM_Network *)dm->data;
8232727e31bSShri Abhyankar 
8242727e31bSShri Abhyankar   PetscFunctionBegin;
825daad07d3SAidan 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);
826daad07d3SAidan Hamilton   if (nv) *nv = network->cloneshared->subnet[netnum].nvtx;
827daad07d3SAidan Hamilton   if (ne) *ne = network->cloneshared->subnet[netnum].nedge;
828daad07d3SAidan Hamilton   if (vtx) *vtx = network->cloneshared->subnet[netnum].vertices;
829daad07d3SAidan Hamilton   if (edge) *edge = network->cloneshared->subnet[netnum].edges;
8303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8312bf73ac6SHong Zhang }
8322bf73ac6SHong Zhang 
8332bf73ac6SHong Zhang /*@
8342bf73ac6SHong Zhang   DMNetworkAddSharedVertices - Add shared vertices that connect two given subnetworks
8352bf73ac6SHong Zhang 
83620f4b53cSBarry Smith   Collective
8372bf73ac6SHong Zhang 
8382bf73ac6SHong Zhang   Input Parameters:
83920f4b53cSBarry Smith + dm      - the `DMNETWORK` object
84020f4b53cSBarry Smith . anetnum - first subnetwork global numbering returned by `DMNetworkAddSubnetwork()`
84120f4b53cSBarry Smith . bnetnum - second subnetwork global numbering returned by `DMNetworkAddSubnetwork()`
8422bf73ac6SHong Zhang . nsvtx   - number of vertices that are shared by the two subnetworks
8432bf73ac6SHong Zhang . asvtx   - vertex index in the first subnetwork
8442bf73ac6SHong Zhang - bsvtx   - vertex index in the second subnetwork
8452bf73ac6SHong Zhang 
8462bf73ac6SHong Zhang   Level: beginner
8472bf73ac6SHong Zhang 
84820f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkCreate()`, `DMNetworkAddSubnetwork()`, `DMNetworkGetSharedVertices()`
8492bf73ac6SHong Zhang @*/
850d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkAddSharedVertices(DM dm, PetscInt anetnum, PetscInt bnetnum, PetscInt nsvtx, PetscInt asvtx[], PetscInt bsvtx[])
851d71ae5a4SJacob Faibussowitsch {
8522bf73ac6SHong Zhang   DM_Network *network = (DM_Network *)dm->data;
853daad07d3SAidan Hamilton   PetscInt    i, nsubnet = network->cloneshared->Nsubnet, *sedgelist, Nsvtx = network->cloneshared->Nsvtx;
8542bf73ac6SHong Zhang 
8552bf73ac6SHong Zhang   PetscFunctionBegin;
8565c6496baSHong Zhang   PetscCheck(anetnum != bnetnum, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Subnetworks must have different netnum");
8575c6496baSHong Zhang   PetscCheck(anetnum >= 0 && bnetnum >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "netnum cannot be negative");
8582bf73ac6SHong Zhang   if (!Nsvtx) {
8592bf73ac6SHong Zhang     /* allocate network->sedgelist to hold at most 2*nsubnet pairs of shared vertices */
860daad07d3SAidan Hamilton     PetscCall(PetscMalloc1(2 * 4 * nsubnet, &network->cloneshared->sedgelist));
8612bf73ac6SHong Zhang   }
8622bf73ac6SHong Zhang 
863daad07d3SAidan Hamilton   sedgelist = network->cloneshared->sedgelist;
8642bf73ac6SHong Zhang   for (i = 0; i < nsvtx; i++) {
8659371c9d4SSatish Balay     sedgelist[4 * Nsvtx]     = anetnum;
8669371c9d4SSatish Balay     sedgelist[4 * Nsvtx + 1] = asvtx[i];
8679371c9d4SSatish Balay     sedgelist[4 * Nsvtx + 2] = bnetnum;
8689371c9d4SSatish Balay     sedgelist[4 * Nsvtx + 3] = bsvtx[i];
8692bf73ac6SHong Zhang     Nsvtx++;
8702bf73ac6SHong Zhang   }
8715c6496baSHong Zhang   PetscCheck(Nsvtx <= 2 * nsubnet, PETSC_COMM_SELF, PETSC_ERR_SUP, "allocate more space for coupling edgelist");
872daad07d3SAidan Hamilton   network->cloneshared->Nsvtx = Nsvtx;
8733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8742727e31bSShri Abhyankar }
8752727e31bSShri Abhyankar 
8762727e31bSShri Abhyankar /*@C
8772bf73ac6SHong Zhang   DMNetworkGetSharedVertices - Returns the info for the shared vertices
8782bf73ac6SHong Zhang 
87920f4b53cSBarry Smith   Not Collective
880caf410d2SHong Zhang 
881f899ff85SJose E. Roman   Input Parameter:
88220f4b53cSBarry Smith . dm - the `DMNETWORK` object
883caf410d2SHong Zhang 
8847a7aea1fSJed Brown   Output Parameters:
8852bf73ac6SHong Zhang + nsv  - number of local shared vertices
8862bf73ac6SHong Zhang - svtx - local shared vertices
887caf410d2SHong Zhang 
888caf410d2SHong Zhang   Level: intermediate
889caf410d2SHong Zhang 
89020f4b53cSBarry Smith   Notes:
89120f4b53cSBarry Smith   Cannot call this routine before `DMNetworkLayoutSetup()`
89220f4b53cSBarry Smith 
89320f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetSubnetwork()`, `DMNetworkLayoutSetUp()`, `DMNetworkAddSharedVertices()`
894caf410d2SHong Zhang @*/
895d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetSharedVertices(DM dm, PetscInt *nsv, const PetscInt **svtx)
896d71ae5a4SJacob Faibussowitsch {
897caf410d2SHong Zhang   DM_Network *net = (DM_Network *)dm->data;
898caf410d2SHong Zhang 
899caf410d2SHong Zhang   PetscFunctionBegin;
9005c6496baSHong Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
901daad07d3SAidan Hamilton   if (nsv) *nsv = net->cloneshared->nsvtx;
902daad07d3SAidan Hamilton   if (svtx) *svtx = net->cloneshared->svertices;
9033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
904caf410d2SHong Zhang }
905caf410d2SHong Zhang 
906caf410d2SHong Zhang /*@C
9075f2c45f1SShri Abhyankar   DMNetworkRegisterComponent - Registers the network component
9085f2c45f1SShri Abhyankar 
90920f4b53cSBarry Smith   Logically Collective
9105f2c45f1SShri Abhyankar 
9117a7aea1fSJed Brown   Input Parameters:
91220f4b53cSBarry Smith + dm   - the `DMNETWORK` object
9135f2c45f1SShri Abhyankar . name - the component name
9145f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data
9155f2c45f1SShri Abhyankar 
9162fe279fdSBarry Smith   Output Parameter:
9175f2c45f1SShri Abhyankar . key - an integer key that defines the component
9185f2c45f1SShri Abhyankar 
91920f4b53cSBarry Smith   Level: beginner
92020f4b53cSBarry Smith 
92187497f52SBarry Smith   Note:
92287497f52SBarry Smith   This routine should be called by all processors before calling `DMNetworkLayoutSetup()`.
9235f2c45f1SShri Abhyankar 
92420f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkCreate()`, `DMNetworkLayoutSetUp()`
9255f2c45f1SShri Abhyankar @*/
926d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkRegisterComponent(DM dm, const char *name, size_t size, PetscInt *key)
927d71ae5a4SJacob Faibussowitsch {
9285f2c45f1SShri Abhyankar   DM_Network         *network   = (DM_Network *)dm->data;
92954dfd506SHong Zhang   DMNetworkComponent *component = NULL, *newcomponent = NULL;
9305f2c45f1SShri Abhyankar   PetscBool           flg = PETSC_FALSE;
9315f2c45f1SShri Abhyankar   PetscInt            i;
9325f2c45f1SShri Abhyankar 
9335f2c45f1SShri Abhyankar   PetscFunctionBegin;
93448a46eb9SPierre Jolivet   if (!network->component) PetscCall(PetscCalloc1(network->max_comps_registered, &network->component));
93554dfd506SHong Zhang 
9365f2c45f1SShri Abhyankar   for (i = 0; i < network->ncomponent; i++) {
9379566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(network->component[i].name, name, &flg));
9385f2c45f1SShri Abhyankar     if (flg) {
9395f2c45f1SShri Abhyankar       *key = i;
9403ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
9415f2c45f1SShri Abhyankar     }
9426d64e262SShri Abhyankar   }
94354dfd506SHong Zhang 
94454dfd506SHong Zhang   if (network->ncomponent == network->max_comps_registered) {
94554dfd506SHong Zhang     /* Reached max allowed so resize component */
94654dfd506SHong Zhang     network->max_comps_registered += 2;
9479566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(network->max_comps_registered, &newcomponent));
94854dfd506SHong Zhang     /* Copy over the previous component info */
94954dfd506SHong Zhang     for (i = 0; i < network->ncomponent; i++) {
950c6a7a370SJeremy L Thompson       PetscCall(PetscStrncpy(newcomponent[i].name, network->component[i].name, sizeof(newcomponent[i].name)));
95154dfd506SHong Zhang       newcomponent[i].size = network->component[i].size;
9525f2c45f1SShri Abhyankar     }
95354dfd506SHong Zhang     /* Free old one */
9549566063dSJacob Faibussowitsch     PetscCall(PetscFree(network->component));
95554dfd506SHong Zhang     /* Update pointer */
95654dfd506SHong Zhang     network->component = newcomponent;
95754dfd506SHong Zhang   }
95854dfd506SHong Zhang 
95954dfd506SHong Zhang   component = &network->component[network->ncomponent];
9605f2c45f1SShri Abhyankar 
961c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(component->name, name, sizeof(component->name)));
9625f2c45f1SShri Abhyankar   component->size = size / sizeof(DMNetworkComponentGenericDataType);
9635f2c45f1SShri Abhyankar   *key            = network->ncomponent;
9645f2c45f1SShri Abhyankar   network->ncomponent++;
9653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9665f2c45f1SShri Abhyankar }
9675f2c45f1SShri Abhyankar 
9685f2c45f1SShri Abhyankar /*@
9698afb7921SAidan Hamilton   DMNetworkGetNumVertices - Get the local and global number of vertices for the entire network.
9708afb7921SAidan Hamilton 
9718afb7921SAidan Hamilton   Not Collective
9728afb7921SAidan Hamilton 
9738afb7921SAidan Hamilton   Input Parameter:
97420f4b53cSBarry Smith . dm - the `DMNETWORK` object
9758afb7921SAidan Hamilton 
9768afb7921SAidan Hamilton   Output Parameters:
9778afb7921SAidan Hamilton + nVertices - the local number of vertices
9788afb7921SAidan Hamilton - NVertices - the global number of vertices
9798afb7921SAidan Hamilton 
9808afb7921SAidan Hamilton   Level: beginner
9818afb7921SAidan Hamilton 
98220f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetNumEdges()`
9838afb7921SAidan Hamilton @*/
9848afb7921SAidan Hamilton PetscErrorCode DMNetworkGetNumVertices(DM dm, PetscInt *nVertices, PetscInt *NVertices)
9858afb7921SAidan Hamilton {
9868afb7921SAidan Hamilton   DM_Network *network = (DM_Network *)dm->data;
9878afb7921SAidan Hamilton 
9888afb7921SAidan Hamilton   PetscFunctionBegin;
9898afb7921SAidan Hamilton   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK);
9908afb7921SAidan Hamilton   if (nVertices) {
9918afb7921SAidan Hamilton     PetscValidIntPointer(nVertices, 2);
9928afb7921SAidan Hamilton     *nVertices = network->cloneshared->nVertices;
9938afb7921SAidan Hamilton   }
9948afb7921SAidan Hamilton   if (NVertices) {
9958afb7921SAidan Hamilton     PetscValidIntPointer(NVertices, 3);
9968afb7921SAidan Hamilton     *NVertices = network->cloneshared->NVertices;
9978afb7921SAidan Hamilton   }
9983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9998afb7921SAidan Hamilton }
10008afb7921SAidan Hamilton 
10018afb7921SAidan Hamilton /*@
10028afb7921SAidan Hamilton   DMNetworkGetNumEdges - Get the local and global number of edges for the entire network.
10038afb7921SAidan Hamilton 
10048afb7921SAidan Hamilton   Not Collective
10058afb7921SAidan Hamilton 
10068afb7921SAidan Hamilton   Input Parameter:
100720f4b53cSBarry Smith . dm - the `DMNETWORK` object
10088afb7921SAidan Hamilton 
10098afb7921SAidan Hamilton   Output Parameters:
10108afb7921SAidan Hamilton + nEdges - the local number of edges
10118afb7921SAidan Hamilton - NEdges - the global number of edges
10128afb7921SAidan Hamilton 
10138afb7921SAidan Hamilton   Level: beginner
10148afb7921SAidan Hamilton 
101520f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetNumVertices()`
10168afb7921SAidan Hamilton @*/
10178afb7921SAidan Hamilton PetscErrorCode DMNetworkGetNumEdges(DM dm, PetscInt *nEdges, PetscInt *NEdges)
10188afb7921SAidan Hamilton {
10198afb7921SAidan Hamilton   DM_Network *network = (DM_Network *)dm->data;
10208afb7921SAidan Hamilton 
10218afb7921SAidan Hamilton   PetscFunctionBegin;
10228afb7921SAidan Hamilton   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMNETWORK);
10238afb7921SAidan Hamilton   if (nEdges) {
10248afb7921SAidan Hamilton     PetscValidIntPointer(nEdges, 2);
10258afb7921SAidan Hamilton     *nEdges = network->cloneshared->nEdges;
10268afb7921SAidan Hamilton   }
10278afb7921SAidan Hamilton   if (NEdges) {
10288afb7921SAidan Hamilton     PetscValidIntPointer(NEdges, 3);
10298afb7921SAidan Hamilton     *NEdges = network->cloneshared->NEdges;
10308afb7921SAidan Hamilton   }
10313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10328afb7921SAidan Hamilton }
10338afb7921SAidan Hamilton 
10348afb7921SAidan Hamilton /*@
10352bf73ac6SHong Zhang   DMNetworkGetVertexRange - Get the bounds [start, end) for the local vertices
10365f2c45f1SShri Abhyankar 
10375f2c45f1SShri Abhyankar   Not Collective
10385f2c45f1SShri Abhyankar 
1039f899ff85SJose E. Roman   Input Parameter:
104020f4b53cSBarry Smith . dm - the `DMNETWORK` object
10415f2c45f1SShri Abhyankar 
1042fd292e60Sprj-   Output Parameters:
10432bf73ac6SHong Zhang + vStart - the first vertex point
10442bf73ac6SHong Zhang - vEnd   - one beyond the last vertex point
10455f2c45f1SShri Abhyankar 
104697bb938eSShri Abhyankar   Level: beginner
10475f2c45f1SShri Abhyankar 
104820f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetEdgeRange()`
10495f2c45f1SShri Abhyankar @*/
1050d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetVertexRange(DM dm, PetscInt *vStart, PetscInt *vEnd)
1051d71ae5a4SJacob Faibussowitsch {
10525f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
10535f2c45f1SShri Abhyankar 
10545f2c45f1SShri Abhyankar   PetscFunctionBegin;
1055daad07d3SAidan Hamilton   if (vStart) *vStart = network->cloneshared->vStart;
1056daad07d3SAidan Hamilton   if (vEnd) *vEnd = network->cloneshared->vEnd;
10573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10585f2c45f1SShri Abhyankar }
10595f2c45f1SShri Abhyankar 
10605f2c45f1SShri Abhyankar /*@
10612bf73ac6SHong Zhang   DMNetworkGetEdgeRange - Get the bounds [start, end) for the local edges
10625f2c45f1SShri Abhyankar 
10635f2c45f1SShri Abhyankar   Not Collective
10645f2c45f1SShri Abhyankar 
1065f899ff85SJose E. Roman   Input Parameter:
106620f4b53cSBarry Smith . dm - the `DMNETWORK` object
10675f2c45f1SShri Abhyankar 
1068fd292e60Sprj-   Output Parameters:
10695f2c45f1SShri Abhyankar + eStart - The first edge point
10705f2c45f1SShri Abhyankar - eEnd   - One beyond the last edge point
10715f2c45f1SShri Abhyankar 
107297bb938eSShri Abhyankar   Level: beginner
10735f2c45f1SShri Abhyankar 
107420f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetVertexRange()`
10755f2c45f1SShri Abhyankar @*/
1076d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetEdgeRange(DM dm, PetscInt *eStart, PetscInt *eEnd)
1077d71ae5a4SJacob Faibussowitsch {
10785f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
10795f2c45f1SShri Abhyankar 
10805f2c45f1SShri Abhyankar   PetscFunctionBegin;
10815c6496baSHong Zhang   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1082daad07d3SAidan Hamilton   if (eStart) *eStart = network->cloneshared->eStart;
1083daad07d3SAidan Hamilton   if (eEnd) *eEnd = network->cloneshared->eEnd;
10843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10855f2c45f1SShri Abhyankar }
10865f2c45f1SShri Abhyankar 
1087d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetIndex(DM dm, PetscInt p, PetscInt *index)
1088d71ae5a4SJacob Faibussowitsch {
10892bf73ac6SHong Zhang   DM_Network *network = (DM_Network *)dm->data;
10905c6496baSHong Zhang 
10915c6496baSHong Zhang   PetscFunctionBegin;
10925c6496baSHong Zhang   if (network->header) {
10935c6496baSHong Zhang     *index = network->header[p].index;
10945c6496baSHong Zhang   } else {
10952bf73ac6SHong Zhang     PetscInt                 offsetp;
10962bf73ac6SHong Zhang     DMNetworkComponentHeader header;
10972bf73ac6SHong Zhang 
10989566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(network->DataSection, p, &offsetp));
10992bf73ac6SHong Zhang     header = (DMNetworkComponentHeader)(network->componentdataarray + offsetp);
11002bf73ac6SHong Zhang     *index = header->index;
11015c6496baSHong Zhang   }
11023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11032bf73ac6SHong Zhang }
11042bf73ac6SHong Zhang 
1105d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetSubnetID(DM dm, PetscInt p, PetscInt *subnetid)
1106d71ae5a4SJacob Faibussowitsch {
1107daad07d3SAidan Hamilton   DM_Network *network = (DM_Network *)dm->data;
1108daad07d3SAidan Hamilton 
1109daad07d3SAidan Hamilton   PetscFunctionBegin;
1110daad07d3SAidan Hamilton   if (network->header) {
1111daad07d3SAidan Hamilton     *subnetid = network->header[p].subnetid;
1112daad07d3SAidan Hamilton   } else {
1113daad07d3SAidan Hamilton     PetscInt                 offsetp;
1114daad07d3SAidan Hamilton     DMNetworkComponentHeader header;
1115daad07d3SAidan Hamilton 
1116daad07d3SAidan Hamilton     PetscCall(PetscSectionGetOffset(network->DataSection, p, &offsetp));
1117daad07d3SAidan Hamilton     header    = (DMNetworkComponentHeader)(network->componentdataarray + offsetp);
1118daad07d3SAidan Hamilton     *subnetid = header->subnetid;
1119daad07d3SAidan Hamilton   }
11203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1121daad07d3SAidan Hamilton }
1122daad07d3SAidan Hamilton 
11237b6afd5bSHong Zhang /*@
11242bf73ac6SHong Zhang   DMNetworkGetGlobalEdgeIndex - Get the global numbering for the edge on the network
11257b6afd5bSHong Zhang 
11267b6afd5bSHong Zhang   Not Collective
11277b6afd5bSHong Zhang 
11287b6afd5bSHong Zhang   Input Parameters:
112920f4b53cSBarry Smith + dm - `DMNETWORK` object
1130e85e6aecSHong Zhang - p  - edge point
11317b6afd5bSHong Zhang 
11322fe279fdSBarry Smith   Output Parameter:
11332bf73ac6SHong Zhang . index - the global numbering for the edge
11347b6afd5bSHong Zhang 
11357b6afd5bSHong Zhang   Level: intermediate
11367b6afd5bSHong Zhang 
113720f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetGlobalVertexIndex()`
11387b6afd5bSHong Zhang @*/
1139d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm, PetscInt p, PetscInt *index)
1140d71ae5a4SJacob Faibussowitsch {
11417b6afd5bSHong Zhang   PetscFunctionBegin;
11429566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetIndex(dm, p, index));
11433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11447b6afd5bSHong Zhang }
11457b6afd5bSHong Zhang 
11465f2c45f1SShri Abhyankar /*@
11472bf73ac6SHong Zhang   DMNetworkGetGlobalVertexIndex - Get the global numbering for the vertex on the network
1148e85e6aecSHong Zhang 
1149e85e6aecSHong Zhang   Not Collective
1150e85e6aecSHong Zhang 
1151e85e6aecSHong Zhang   Input Parameters:
115220f4b53cSBarry Smith + dm - `DMNETWORK` object
1153e85e6aecSHong Zhang - p  - vertex point
1154e85e6aecSHong Zhang 
11552fe279fdSBarry Smith   Output Parameter:
11562bf73ac6SHong Zhang . index - the global numbering for the vertex
1157e85e6aecSHong Zhang 
1158e85e6aecSHong Zhang   Level: intermediate
1159e85e6aecSHong Zhang 
116020f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetGlobalEdgeIndex()`, `DMNetworkGetLocalVertexIndex()`
1161e85e6aecSHong Zhang @*/
1162d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm, PetscInt p, PetscInt *index)
1163d71ae5a4SJacob Faibussowitsch {
1164e85e6aecSHong Zhang   PetscFunctionBegin;
11659566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetIndex(dm, p, index));
11663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11679988b915SShri Abhyankar }
11689988b915SShri Abhyankar 
11699988b915SShri Abhyankar /*@
11705f2c45f1SShri Abhyankar   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
11715f2c45f1SShri Abhyankar 
11725f2c45f1SShri Abhyankar   Not Collective
11735f2c45f1SShri Abhyankar 
11745f2c45f1SShri Abhyankar   Input Parameters:
117520f4b53cSBarry Smith + dm - the `DMNETWORK` object
1176a2b725a8SWilliam Gropp - p  - vertex/edge point
11775f2c45f1SShri Abhyankar 
11782fe279fdSBarry Smith   Output Parameter:
11795f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge
11805f2c45f1SShri Abhyankar 
118197bb938eSShri Abhyankar   Level: beginner
11825f2c45f1SShri Abhyankar 
118320f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkRegisterComponent()`, `DMNetworkAddComponent()`
11845f2c45f1SShri Abhyankar @*/
1185d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetNumComponents(DM dm, PetscInt p, PetscInt *numcomponents)
1186d71ae5a4SJacob Faibussowitsch {
11875f2c45f1SShri Abhyankar   PetscInt    offset;
11885f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
11895f2c45f1SShri Abhyankar 
11905f2c45f1SShri Abhyankar   PetscFunctionBegin;
11919566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection, p, &offset));
11925f2c45f1SShri Abhyankar   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray + offset))->ndata;
11933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11945f2c45f1SShri Abhyankar }
11955f2c45f1SShri Abhyankar 
11965f2c45f1SShri Abhyankar /*@
11972bf73ac6SHong Zhang   DMNetworkGetLocalVecOffset - Get the offset for accessing the variables associated with a component at the given vertex/edge from the local vector
11985f2c45f1SShri Abhyankar 
11995f2c45f1SShri Abhyankar   Not Collective
12005f2c45f1SShri Abhyankar 
12015f2c45f1SShri Abhyankar   Input Parameters:
120220f4b53cSBarry Smith + dm      - the `DMNETWORK` object
1203f2f58720SBarry Smith . p       - the edge or vertex point
12042bf73ac6SHong Zhang - compnum - component number; use ALL_COMPONENTS if no specific component is requested
12059988b915SShri Abhyankar 
12062fe279fdSBarry Smith   Output Parameter:
12072bf73ac6SHong Zhang . offset - the local offset
12089988b915SShri Abhyankar 
12099988b915SShri Abhyankar   Level: intermediate
12109988b915SShri Abhyankar 
121120f4b53cSBarry Smith   Notes:
121220f4b53cSBarry Smith   These offsets can be passed to `MatSetValuesLocal()` for matrices obtained with `DMCreateMatrix()`.
121320f4b53cSBarry Smith 
121420f4b53cSBarry Smith   For vectors obtained with `DMCreateLocalVector()` the offsets can be used with `VecSetValues()`.
121520f4b53cSBarry Smith 
121620f4b53cSBarry Smith   For vectors obtained with `DMCreateLocalVector()` and the array obtained with `VecGetArray`(vec,&array) you can access or set
121720f4b53cSBarry Smith   the vector values with array[offset].
121820f4b53cSBarry Smith 
121920f4b53cSBarry Smith   For vectors obtained with `DMCreateGlobalVector()` the offsets can be used with `VecSetValuesLocal()`.
122020f4b53cSBarry Smith 
122120f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMGetLocalVector()`, `DMNetworkGetComponent()`, `DMNetworkGetGlobalVecOffset()`, `DMCreateGlobalVector()`, `VecGetArray()`, `VecSetValuesLocal()`, `MatSetValuesLocal()`
12229988b915SShri Abhyankar @*/
1223d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetLocalVecOffset(DM dm, PetscInt p, PetscInt compnum, PetscInt *offset)
1224d71ae5a4SJacob Faibussowitsch {
12259988b915SShri Abhyankar   DM_Network              *network = (DM_Network *)dm->data;
12269988b915SShri Abhyankar   PetscInt                 offsetp, offsetd;
12279988b915SShri Abhyankar   DMNetworkComponentHeader header;
12289988b915SShri Abhyankar 
12299988b915SShri Abhyankar   PetscFunctionBegin;
12309566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->plex->localSection, p, &offsetp));
12312bf73ac6SHong Zhang   if (compnum == ALL_COMPONENTS) {
12322bf73ac6SHong Zhang     *offset = offsetp;
12333ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
12342bf73ac6SHong Zhang   }
12352bf73ac6SHong Zhang 
12369566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection, p, &offsetd));
12379988b915SShri Abhyankar   header  = (DMNetworkComponentHeader)(network->componentdataarray + offsetd);
12389988b915SShri Abhyankar   *offset = offsetp + header->offsetvarrel[compnum];
12393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12409988b915SShri Abhyankar }
12419988b915SShri Abhyankar 
12429988b915SShri Abhyankar /*@
12432bf73ac6SHong Zhang   DMNetworkGetGlobalVecOffset - Get the global offset for accessing the variables associated with a component for the given vertex/edge from the global vector
12449988b915SShri Abhyankar 
12459988b915SShri Abhyankar   Not Collective
12469988b915SShri Abhyankar 
12479988b915SShri Abhyankar   Input Parameters:
124820f4b53cSBarry Smith + dm      - the `DMNETWORK` object
1249f2f58720SBarry Smith . p       - the edge or vertex point
12502bf73ac6SHong Zhang - compnum - component number; use ALL_COMPONENTS if no specific component is requested
12519988b915SShri Abhyankar 
12522fe279fdSBarry Smith   Output Parameter:
12539988b915SShri Abhyankar . offsetg - the global offset
12549988b915SShri Abhyankar 
12559988b915SShri Abhyankar   Level: intermediate
12569988b915SShri Abhyankar 
125720f4b53cSBarry Smith   Notes:
125820f4b53cSBarry Smith   These offsets can be passed to `MatSetValues()` for matrices obtained with `DMCreateMatrix()`.
125920f4b53cSBarry Smith 
126020f4b53cSBarry Smith   For vectors obtained with `DMCreateGlobalVector()` the offsets can be used with `VecSetValues()`.
126120f4b53cSBarry Smith 
126220f4b53cSBarry Smith   For vectors obtained with `DMCreateGlobalVector()` and the array obtained with `VecGetArray`(vec,&array) you can access or set
126320f4b53cSBarry Smith   the vector values with array[offset - rstart] where restart is obtained with `VecGetOwnershipRange`(v,&rstart,`NULL`);
126420f4b53cSBarry Smith 
126520f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetLocalVecOffset()`, `DMGetGlobalVector()`, `DMNetworkGetComponent()`, `DMCreateGlobalVector()`, `VecGetArray()`, `VecSetValues()`, `MatSetValues()`
12669988b915SShri Abhyankar @*/
1267d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetGlobalVecOffset(DM dm, PetscInt p, PetscInt compnum, PetscInt *offsetg)
1268d71ae5a4SJacob Faibussowitsch {
12699988b915SShri Abhyankar   DM_Network              *network = (DM_Network *)dm->data;
12709988b915SShri Abhyankar   PetscInt                 offsetp, offsetd;
12719988b915SShri Abhyankar   DMNetworkComponentHeader header;
12729988b915SShri Abhyankar 
12739988b915SShri Abhyankar   PetscFunctionBegin;
12749566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->plex->globalSection, p, &offsetp));
12752bf73ac6SHong Zhang   if (offsetp < 0) offsetp = -(offsetp + 1); /* Convert to actual global offset for ghost vertex */
12762bf73ac6SHong Zhang 
12772bf73ac6SHong Zhang   if (compnum == ALL_COMPONENTS) {
12782bf73ac6SHong Zhang     *offsetg = offsetp;
12793ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
12802bf73ac6SHong Zhang   }
12819566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection, p, &offsetd));
12829988b915SShri Abhyankar   header   = (DMNetworkComponentHeader)(network->componentdataarray + offsetd);
12839988b915SShri Abhyankar   *offsetg = offsetp + header->offsetvarrel[compnum];
12843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12859988b915SShri Abhyankar }
12869988b915SShri Abhyankar 
12879988b915SShri Abhyankar /*@
12882bf73ac6SHong Zhang   DMNetworkGetEdgeOffset - Get the offset for accessing the variables associated with the given edge from the local subvector
128924121865SAdrian Maldonado 
129024121865SAdrian Maldonado   Not Collective
129124121865SAdrian Maldonado 
129224121865SAdrian Maldonado   Input Parameters:
129320f4b53cSBarry Smith + dm - the `DMNETWORK` object
129424121865SAdrian Maldonado - p  - the edge point
129524121865SAdrian Maldonado 
12962fe279fdSBarry Smith   Output Parameter:
129724121865SAdrian Maldonado . offset - the offset
129824121865SAdrian Maldonado 
129924121865SAdrian Maldonado   Level: intermediate
130024121865SAdrian Maldonado 
130120f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetLocalVecOffset()`, `DMGetLocalVector()`
130224121865SAdrian Maldonado @*/
1303d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetEdgeOffset(DM dm, PetscInt p, PetscInt *offset)
1304d71ae5a4SJacob Faibussowitsch {
130524121865SAdrian Maldonado   DM_Network *network = (DM_Network *)dm->data;
130624121865SAdrian Maldonado 
130724121865SAdrian Maldonado   PetscFunctionBegin;
13089566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->edge.DofSection, p, offset));
13093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
131024121865SAdrian Maldonado }
131124121865SAdrian Maldonado 
131224121865SAdrian Maldonado /*@
13132bf73ac6SHong Zhang   DMNetworkGetVertexOffset - Get the offset for accessing the variables associated with the given vertex from the local subvector
131424121865SAdrian Maldonado 
131524121865SAdrian Maldonado   Not Collective
131624121865SAdrian Maldonado 
131724121865SAdrian Maldonado   Input Parameters:
131820f4b53cSBarry Smith + dm - the `DMNETWORK` object
131924121865SAdrian Maldonado - p  - the vertex point
132024121865SAdrian Maldonado 
13212fe279fdSBarry Smith   Output Parameter:
132224121865SAdrian Maldonado . offset - the offset
132324121865SAdrian Maldonado 
132424121865SAdrian Maldonado   Level: intermediate
132524121865SAdrian Maldonado 
132620f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetEdgeOffset()`, `DMGetLocalVector()`
132724121865SAdrian Maldonado @*/
1328d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetVertexOffset(DM dm, PetscInt p, PetscInt *offset)
1329d71ae5a4SJacob Faibussowitsch {
133024121865SAdrian Maldonado   DM_Network *network = (DM_Network *)dm->data;
133124121865SAdrian Maldonado 
133224121865SAdrian Maldonado   PetscFunctionBegin;
1333daad07d3SAidan Hamilton   p -= network->cloneshared->vStart;
13349566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->vertex.DofSection, p, offset));
13353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
133624121865SAdrian Maldonado }
13372bf73ac6SHong Zhang 
13385f2c45f1SShri Abhyankar /*@
13392bf73ac6SHong Zhang   DMNetworkAddComponent - Adds a network component and number of variables at the given point (vertex/edge)
13405f2c45f1SShri Abhyankar 
134120f4b53cSBarry Smith   Collective
13425f2c45f1SShri Abhyankar 
13435f2c45f1SShri Abhyankar   Input Parameters:
13444dc646bcSBarry Smith + dm           - the DMNetwork
134520f4b53cSBarry Smith . p            - the vertex/edge point. These points are local indices provided by `DMNetworkGetSubnetwork()`
134620f4b53cSBarry Smith . componentkey - component key returned while registering the component with `DMNetworkRegisterComponent()`
134720f4b53cSBarry 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
134820f4b53cSBarry Smith               free this space until after `DMSetUp()` is called.
1349eac198afSGetnet - 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
13505f2c45f1SShri Abhyankar 
135120f4b53cSBarry Smith   Level: beginner
135220f4b53cSBarry Smith 
1353eac198afSGetnet   Notes:
1354eac198afSGetnet   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.
1355eac198afSGetnet 
135620f4b53cSBarry Smith   `DMNetworkLayoutSetUp()` must be called before this routine.
1357eac198afSGetnet 
1358*60225df5SJacob Faibussowitsch   Developer Notes:
135920f4b53cSBarry Smith   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`.
13605f2c45f1SShri Abhyankar 
136120f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetComponent()`, `DMNetworkGetSubnetwork()`, `DMNetworkIsGhostVertex()`, `DMNetworkLayoutSetUp()`
13625f2c45f1SShri Abhyankar @*/
1363d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p, PetscInt componentkey, void *compvalue, PetscInt nvar)
1364d71ae5a4SJacob Faibussowitsch {
13655f2c45f1SShri Abhyankar   DM_Network              *network   = (DM_Network *)dm->data;
13662bf73ac6SHong Zhang   DMNetworkComponent      *component = &network->component[componentkey];
136754dfd506SHong Zhang   DMNetworkComponentHeader header;
136854dfd506SHong Zhang   DMNetworkComponentValue  cvalue;
136954dfd506SHong Zhang   PetscInt                 compnum;
137054dfd506SHong Zhang   PetscInt                *compsize, *compkey, *compoffset, *compnvar, *compoffsetvarrel;
137154dfd506SHong Zhang   void                   **compdata;
13725f2c45f1SShri Abhyankar 
13735f2c45f1SShri Abhyankar   PetscFunctionBegin;
13745c6496baSHong 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);
1375daad07d3SAidan 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.");
1376eac198afSGetnet   /* The owning rank and all ghost ranks add nvar */
13779566063dSJacob Faibussowitsch   PetscCall(PetscSectionAddDof(network->DofSection, p, nvar));
13782bf73ac6SHong Zhang 
1379eac198afSGetnet   /* The owning rank and all ghost ranks add a component, including compvalue=NULL */
138054dfd506SHong Zhang   header = &network->header[p];
138154dfd506SHong Zhang   cvalue = &network->cvalue[p];
138254dfd506SHong Zhang   if (header->ndata == header->maxcomps) {
1383f11a936eSBarry Smith     PetscInt additional_size;
1384f11a936eSBarry Smith 
138554dfd506SHong Zhang     /* Reached limit so resize header component arrays */
138654dfd506SHong Zhang     header->maxcomps += 2;
138754dfd506SHong Zhang 
138854dfd506SHong Zhang     /* Allocate arrays for component information and value */
13899566063dSJacob Faibussowitsch     PetscCall(PetscCalloc5(header->maxcomps, &compsize, header->maxcomps, &compkey, header->maxcomps, &compoffset, header->maxcomps, &compnvar, header->maxcomps, &compoffsetvarrel));
13909566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(header->maxcomps, &compdata));
139154dfd506SHong Zhang 
139254dfd506SHong Zhang     /* Recalculate header size */
139354dfd506SHong Zhang     header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5 * header->maxcomps * sizeof(PetscInt);
139454dfd506SHong Zhang     header->hsize /= sizeof(DMNetworkComponentGenericDataType);
1395dd5bbc93SStefano Zampini #if defined(__NEC__)
1396dd5bbc93SStefano Zampini     /* NEC/LG: quick hack to keep data aligned on 8 bytes. */
1397dd5bbc93SStefano Zampini     header->hsize = (header->hsize + (8 - 1)) & ~(8 - 1);
1398dd5bbc93SStefano Zampini #endif
139954dfd506SHong Zhang 
140054dfd506SHong Zhang     /* Copy over component info */
14019566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(compsize, header->size, header->ndata * sizeof(PetscInt)));
14029566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(compkey, header->key, header->ndata * sizeof(PetscInt)));
14039566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(compoffset, header->offset, header->ndata * sizeof(PetscInt)));
14049566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(compnvar, header->nvar, header->ndata * sizeof(PetscInt)));
14059566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(compoffsetvarrel, header->offsetvarrel, header->ndata * sizeof(PetscInt)));
140654dfd506SHong Zhang 
140754dfd506SHong Zhang     /* Copy over component data pointers */
14089566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(compdata, cvalue->data, header->ndata * sizeof(void *)));
140954dfd506SHong Zhang 
141054dfd506SHong Zhang     /* Free old arrays */
14119566063dSJacob Faibussowitsch     PetscCall(PetscFree5(header->size, header->key, header->offset, header->nvar, header->offsetvarrel));
14129566063dSJacob Faibussowitsch     PetscCall(PetscFree(cvalue->data));
141354dfd506SHong Zhang 
141454dfd506SHong Zhang     /* Update pointers */
141554dfd506SHong Zhang     header->size         = compsize;
141654dfd506SHong Zhang     header->key          = compkey;
141754dfd506SHong Zhang     header->offset       = compoffset;
141854dfd506SHong Zhang     header->nvar         = compnvar;
141954dfd506SHong Zhang     header->offsetvarrel = compoffsetvarrel;
142054dfd506SHong Zhang 
142154dfd506SHong Zhang     cvalue->data = compdata;
142254dfd506SHong Zhang 
142354dfd506SHong Zhang     /* Update DataSection Dofs */
142454dfd506SHong 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 */
1425f11a936eSBarry Smith     additional_size = (5 * (header->maxcomps - header->ndata) * sizeof(PetscInt)) / sizeof(DMNetworkComponentGenericDataType);
14269566063dSJacob Faibussowitsch     PetscCall(PetscSectionAddDof(network->DataSection, p, additional_size));
142754dfd506SHong Zhang   }
142854dfd506SHong Zhang   header = &network->header[p];
142954dfd506SHong Zhang   cvalue = &network->cvalue[p];
143054dfd506SHong Zhang 
143154dfd506SHong Zhang   compnum = header->ndata;
14322bf73ac6SHong Zhang 
14332bf73ac6SHong Zhang   header->size[compnum] = component->size;
14349566063dSJacob Faibussowitsch   PetscCall(PetscSectionAddDof(network->DataSection, p, component->size));
14352bf73ac6SHong Zhang   header->key[compnum] = componentkey;
14362bf73ac6SHong Zhang   if (compnum != 0) header->offset[compnum] = header->offset[compnum - 1] + header->size[compnum - 1];
14372bf73ac6SHong Zhang   cvalue->data[compnum] = (void *)compvalue;
14382bf73ac6SHong Zhang 
14392bf73ac6SHong Zhang   /* variables */
14402bf73ac6SHong Zhang   header->nvar[compnum] += nvar;
14412bf73ac6SHong Zhang   if (compnum != 0) header->offsetvarrel[compnum] = header->offsetvarrel[compnum - 1] + header->nvar[compnum - 1];
14422bf73ac6SHong Zhang 
14432bf73ac6SHong Zhang   header->ndata++;
14443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14455f2c45f1SShri Abhyankar }
14465f2c45f1SShri Abhyankar 
144727f51fceSHong Zhang /*@
14482bf73ac6SHong Zhang   DMNetworkGetComponent - Gets the component key, the component data, and the number of variables at a given network point
144927f51fceSHong Zhang 
145027f51fceSHong Zhang   Not Collective
145127f51fceSHong Zhang 
145227f51fceSHong Zhang   Input Parameters:
145320f4b53cSBarry Smith + dm      - the `DMNETWORK` object
14542bf73ac6SHong Zhang . p       - vertex/edge point
145599c90e12SSatish Balay - compnum - component number; use ALL_COMPONENTS if sum up all the components
145627f51fceSHong Zhang 
145727f51fceSHong Zhang   Output Parameters:
145820f4b53cSBarry Smith + compkey   - the key obtained when registering the component (use `NULL` if not required)
145920f4b53cSBarry Smith . component - the component data (use `NULL` if not required)
146020f4b53cSBarry Smith - nvar      - number of variables (use `NULL` if not required)
146127f51fceSHong Zhang 
146297bb938eSShri Abhyankar   Level: beginner
146327f51fceSHong Zhang 
146420f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkAddComponent()`, `DMNetworkGetNumComponents()`
146527f51fceSHong Zhang @*/
1466d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *compkey, void **component, PetscInt *nvar)
1467d71ae5a4SJacob Faibussowitsch {
146827f51fceSHong Zhang   DM_Network              *network = (DM_Network *)dm->data;
14692bf73ac6SHong Zhang   PetscInt                 offset  = 0;
14702bf73ac6SHong Zhang   DMNetworkComponentHeader header;
147127f51fceSHong Zhang 
147227f51fceSHong Zhang   PetscFunctionBegin;
14732bf73ac6SHong Zhang   if (compnum == ALL_COMPONENTS) {
14749566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(network->DofSection, p, nvar));
14753ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
147627f51fceSHong Zhang   }
147727f51fceSHong Zhang 
14789566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection, p, &offset));
147942dc13f1SHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray + offset);
14805f2c45f1SShri Abhyankar 
14812bf73ac6SHong Zhang   if (compnum >= 0) {
14822bf73ac6SHong Zhang     if (compkey) *compkey = header->key[compnum];
14832bf73ac6SHong Zhang     if (component) {
148454dfd506SHong Zhang       offset += header->hsize + header->offset[compnum];
14852bf73ac6SHong Zhang       *component = network->componentdataarray + offset;
14862bf73ac6SHong Zhang     }
14872bf73ac6SHong Zhang   }
14885f2c45f1SShri Abhyankar 
14892bf73ac6SHong Zhang   if (nvar) *nvar = header->nvar[compnum];
149054dfd506SHong Zhang 
14913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14925f2c45f1SShri Abhyankar }
14935f2c45f1SShri Abhyankar 
14942bf73ac6SHong Zhang /*
14952bf73ac6SHong Zhang  Sets up the array that holds the data for all components and its associated section.
14962bf73ac6SHong 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.
14972bf73ac6SHong Zhang */
1498d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkComponentSetUp(DM dm)
1499d71ae5a4SJacob Faibussowitsch {
15005f2c45f1SShri Abhyankar   DM_Network                        *network = (DM_Network *)dm->data;
1501f11a936eSBarry Smith   PetscInt                           arr_size, p, offset, offsetp, ncomp, i, *headerarr;
15025f2c45f1SShri Abhyankar   DMNetworkComponentHeader           header;
15035f2c45f1SShri Abhyankar   DMNetworkComponentValue            cvalue;
1504f11a936eSBarry Smith   DMNetworkComponentHeader           headerinfo;
15055f2c45f1SShri Abhyankar   DMNetworkComponentGenericDataType *componentdataarray;
15065f2c45f1SShri Abhyankar 
15075f2c45f1SShri Abhyankar   PetscFunctionBegin;
15089566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(network->DataSection));
15099566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(network->DataSection, &arr_size));
1510f11a936eSBarry Smith   /* arr_size+1 fixes pipeline test of opensolaris-misc for src/dm/tests/ex10.c -- Do not know why */
15119566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(arr_size + 1, &network->componentdataarray));
15125f2c45f1SShri Abhyankar   componentdataarray = network->componentdataarray;
1513daad07d3SAidan Hamilton   for (p = network->cloneshared->pStart; p < network->cloneshared->pEnd; p++) {
15149566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(network->DataSection, p, &offsetp));
15155f2c45f1SShri Abhyankar     /* Copy header */
15165f2c45f1SShri Abhyankar     header     = &network->header[p];
1517f11a936eSBarry Smith     headerinfo = (DMNetworkComponentHeader)(componentdataarray + offsetp);
15189566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(headerinfo, header, sizeof(struct _p_DMNetworkComponentHeader)));
1519f11a936eSBarry Smith     headerarr = (PetscInt *)(headerinfo + 1);
15209566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(headerarr, header->size, header->maxcomps * sizeof(PetscInt)));
1521aeb78aa7SBarry Smith     headerinfo->size = headerarr;
152254dfd506SHong Zhang     headerarr += header->maxcomps;
15239566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(headerarr, header->key, header->maxcomps * sizeof(PetscInt)));
1524aeb78aa7SBarry Smith     headerinfo->key = headerarr;
152554dfd506SHong Zhang     headerarr += header->maxcomps;
15269566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(headerarr, header->offset, header->maxcomps * sizeof(PetscInt)));
1527aeb78aa7SBarry Smith     headerinfo->offset = headerarr;
152854dfd506SHong Zhang     headerarr += header->maxcomps;
15299566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(headerarr, header->nvar, header->maxcomps * sizeof(PetscInt)));
1530aeb78aa7SBarry Smith     headerinfo->nvar = headerarr;
153154dfd506SHong Zhang     headerarr += header->maxcomps;
15329566063dSJacob Faibussowitsch     PetscCall(PetscMemcpy(headerarr, header->offsetvarrel, header->maxcomps * sizeof(PetscInt)));
1533aeb78aa7SBarry Smith     headerinfo->offsetvarrel = headerarr;
153454dfd506SHong Zhang 
15355f2c45f1SShri Abhyankar     /* Copy data */
15365f2c45f1SShri Abhyankar     cvalue = &network->cvalue[p];
15375f2c45f1SShri Abhyankar     ncomp  = header->ndata;
15382bf73ac6SHong Zhang 
15395f2c45f1SShri Abhyankar     for (i = 0; i < ncomp; i++) {
154054dfd506SHong Zhang       offset = offsetp + header->hsize + header->offset[i];
15419566063dSJacob Faibussowitsch       PetscCall(PetscMemcpy(componentdataarray + offset, cvalue->data[i], header->size[i] * sizeof(DMNetworkComponentGenericDataType)));
15425f2c45f1SShri Abhyankar     }
15435f2c45f1SShri Abhyankar   }
1544aeb78aa7SBarry Smith 
1545daad07d3SAidan Hamilton   for (i = network->cloneshared->pStart; i < network->cloneshared->pEnd; i++) {
1546aeb78aa7SBarry Smith     PetscCall(PetscFree5(network->header[i].size, network->header[i].key, network->header[i].offset, network->header[i].nvar, network->header[i].offsetvarrel));
1547aeb78aa7SBarry Smith     PetscCall(PetscFree(network->cvalue[i].data));
1548aeb78aa7SBarry Smith   }
1549aeb78aa7SBarry Smith   PetscCall(PetscFree2(network->header, network->cvalue));
15503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15515f2c45f1SShri Abhyankar }
15525f2c45f1SShri Abhyankar 
15535f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */
1554d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMNetworkVariablesSetUp(DM dm)
1555d71ae5a4SJacob Faibussowitsch {
15565f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
15575f2c45f1SShri Abhyankar 
15585f2c45f1SShri Abhyankar   PetscFunctionBegin;
15599566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(network->DofSection));
15603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15615f2c45f1SShri Abhyankar }
15625f2c45f1SShri Abhyankar 
156324121865SAdrian Maldonado /* Get a subsection from a range of points */
1564d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMNetworkGetSubSection_private(PetscSection main, PetscInt pstart, PetscInt pend, PetscSection *subsection)
1565d71ae5a4SJacob Faibussowitsch {
156624121865SAdrian Maldonado   PetscInt i, nvar;
156724121865SAdrian Maldonado 
156824121865SAdrian Maldonado   PetscFunctionBegin;
15699566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)main), subsection));
15709566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(*subsection, 0, pend - pstart));
157124121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
15729566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(main, i, &nvar));
15739566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(*subsection, i - pstart, nvar));
157424121865SAdrian Maldonado   }
157524121865SAdrian Maldonado 
15769566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(*subsection));
15773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
157824121865SAdrian Maldonado }
157924121865SAdrian Maldonado 
158024121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */
1581e66d8692SHong Zhang static PetscErrorCode DMNetworkSetSubMap_private(DM dm, PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
1582d71ae5a4SJacob Faibussowitsch {
158324121865SAdrian Maldonado   PetscInt i, *subpoints;
158424121865SAdrian Maldonado 
158524121865SAdrian Maldonado   PetscFunctionBegin;
158624121865SAdrian Maldonado   /* Create index sets to map from "points" to "subpoints" */
15879566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pend - pstart, &subpoints));
1588ad540459SPierre Jolivet   for (i = pstart; i < pend; i++) subpoints[i - pstart] = i;
1589e66d8692SHong Zhang   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), 1, pend - pstart, subpoints, PETSC_COPY_VALUES, map));
15909566063dSJacob Faibussowitsch   PetscCall(PetscFree(subpoints));
15913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
159224121865SAdrian Maldonado }
159324121865SAdrian Maldonado 
159424121865SAdrian Maldonado /*@
159520f4b53cSBarry Smith   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after `DMNetworkDistribute()`
159624121865SAdrian Maldonado 
159720f4b53cSBarry Smith   Collective
159824121865SAdrian Maldonado 
15992fe279fdSBarry Smith   Input Parameter:
160020f4b53cSBarry Smith . dm - the `DMNETWORK` Object
160124121865SAdrian Maldonado 
160220f4b53cSBarry Smith   Level: intermediate
160320f4b53cSBarry Smith 
160420f4b53cSBarry Smith   Note:
1605*60225df5SJacob Faibussowitsch 
160620f4b53cSBarry Smith   the routine will create alternative orderings for the vertices and edges. Assume global network points are:
160724121865SAdrian Maldonado 
160824121865SAdrian Maldonado   points = [0 1 2 3 4 5 6]
160924121865SAdrian Maldonado 
1610caf410d2SHong 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]).
161124121865SAdrian Maldonado 
161220f4b53cSBarry Smith   With this new ordering a local `PetscSection`, global `PetscSection` and` PetscSF` will be created specific to the subset.
161324121865SAdrian Maldonado 
161424121865SAdrian Maldonado @*/
1615d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
1616d71ae5a4SJacob Faibussowitsch {
161724121865SAdrian Maldonado   MPI_Comm    comm;
1618f11a936eSBarry Smith   PetscMPIInt size;
161924121865SAdrian Maldonado   DM_Network *network = (DM_Network *)dm->data;
162024121865SAdrian Maldonado 
1621eab1376dSHong Zhang   PetscFunctionBegin;
16229566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
16239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
162424121865SAdrian Maldonado 
162524121865SAdrian Maldonado   /* Create maps for vertices and edges */
1626e66d8692SHong Zhang   PetscCall(DMNetworkSetSubMap_private(dm, network->cloneshared->vStart, network->cloneshared->vEnd, &network->vertex.mapping));
1627e66d8692SHong Zhang   PetscCall(DMNetworkSetSubMap_private(dm, network->cloneshared->eStart, network->cloneshared->eEnd, &network->edge.mapping));
162824121865SAdrian Maldonado 
162924121865SAdrian Maldonado   /* Create local sub-sections */
1630daad07d3SAidan Hamilton   PetscCall(DMNetworkGetSubSection_private(network->DofSection, network->cloneshared->vStart, network->cloneshared->vEnd, &network->vertex.DofSection));
1631daad07d3SAidan Hamilton   PetscCall(DMNetworkGetSubSection_private(network->DofSection, network->cloneshared->eStart, network->cloneshared->eEnd, &network->edge.DofSection));
163224121865SAdrian Maldonado 
16339852e123SBarry Smith   if (size > 1) {
16349566063dSJacob Faibussowitsch     PetscCall(PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf));
163522bbedd7SHong Zhang 
16369566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection));
16379566063dSJacob Faibussowitsch     PetscCall(PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf));
16389566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection));
163924121865SAdrian Maldonado   } else {
164024121865SAdrian Maldonado     /* create structures for vertex */
16419566063dSJacob Faibussowitsch     PetscCall(PetscSectionClone(network->vertex.DofSection, &network->vertex.GlobalDofSection));
164224121865SAdrian Maldonado     /* create structures for edge */
16439566063dSJacob Faibussowitsch     PetscCall(PetscSectionClone(network->edge.DofSection, &network->edge.GlobalDofSection));
164424121865SAdrian Maldonado   }
164524121865SAdrian Maldonado 
164624121865SAdrian Maldonado   /* Add viewers */
16479566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)network->edge.GlobalDofSection, "Global edge dof section"));
16489566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection, "Global vertex dof section"));
16499566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view"));
16509566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view"));
16513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
165224121865SAdrian Maldonado }
16537b6afd5bSHong Zhang 
1654f11a936eSBarry Smith /*
16555c6496baSHong Zhang    Setup a lookup btable for the input v's owning subnetworks
16565c6496baSHong Zhang    - add all owing subnetworks that connect to this v to the btable
1657f11a936eSBarry Smith      vertex_subnetid = supportingedge_subnetid
1658f11a936eSBarry Smith */
1659d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode SetSubnetIdLookupBT(DM dm, PetscInt v, PetscInt Nsubnet, PetscBT btable)
1660d71ae5a4SJacob Faibussowitsch {
1661f11a936eSBarry Smith   PetscInt                 e, nedges, offset;
1662f11a936eSBarry Smith   const PetscInt          *edges;
1663f11a936eSBarry Smith   DM_Network              *newDMnetwork = (DM_Network *)dm->data;
1664f11a936eSBarry Smith   DMNetworkComponentHeader header;
1665f11a936eSBarry Smith 
1666f11a936eSBarry Smith   PetscFunctionBegin;
16679566063dSJacob Faibussowitsch   PetscCall(PetscBTMemzero(Nsubnet, btable));
16689566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSupportingEdges(dm, v, &nedges, &edges));
1669f11a936eSBarry Smith   for (e = 0; e < nedges; e++) {
16709566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection, edges[e], &offset));
1671f11a936eSBarry Smith     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);
16729566063dSJacob Faibussowitsch     PetscCall(PetscBTSet(btable, header->subnetid));
1673f11a936eSBarry Smith   }
16743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1675f11a936eSBarry Smith }
1676f11a936eSBarry Smith 
16778afb7921SAidan Hamilton /*
16788afb7921SAidan Hamilton   DMNetworkDistributeCoordinates - Internal function to distribute the coordinate network and coordinates.
16798afb7921SAidan Hamilton 
16808afb7921SAidan Hamilton   Collective
16818afb7921SAidan Hamilton 
16828afb7921SAidan Hamilton   Input Parameters:
168320f4b53cSBarry Smith   + dm - The original `DMNETWORK` object
1684aaa8cc7dSPierre Jolivet   - migrationSF - The `PetscSF` describing the migration from dm to dmnew
16858afb7921SAidan Hamilton   - newDM - The new distributed dmnetwork object.
16868afb7921SAidan Hamilton */
16878afb7921SAidan Hamilton 
16888afb7921SAidan Hamilton static PetscErrorCode DMNetworkDistributeCoordinates(DM dm, PetscSF migrationSF, DM newDM)
16898afb7921SAidan Hamilton {
16908afb7921SAidan Hamilton   DM_Network              *newDMnetwork = (DM_Network *)((newDM)->data), *newCoordnetwork, *oldCoordnetwork;
16918afb7921SAidan Hamilton   DM                       cdm, newcdm;
16928afb7921SAidan Hamilton   PetscInt                 cdim, bs, p, pStart, pEnd, offset;
16938afb7921SAidan Hamilton   Vec                      oldCoord, newCoord;
16948afb7921SAidan Hamilton   DMNetworkComponentHeader header;
16958afb7921SAidan Hamilton   const char              *name;
16968afb7921SAidan Hamilton 
16978afb7921SAidan Hamilton   PetscFunctionBegin;
16988afb7921SAidan Hamilton   /* Distribute the coordinate network and coordinates */
16998afb7921SAidan Hamilton   PetscCall(DMGetCoordinateDim(dm, &cdim));
17008afb7921SAidan Hamilton   PetscCall(DMSetCoordinateDim(newDM, cdim));
17018afb7921SAidan Hamilton 
17028afb7921SAidan Hamilton   /* Migrate only if original network had coordinates */
17038afb7921SAidan Hamilton   PetscCall(DMGetCoordinatesLocal(dm, &oldCoord));
17048afb7921SAidan Hamilton   if (oldCoord) {
17058afb7921SAidan Hamilton     PetscCall(DMGetCoordinateDM(dm, &cdm));
17068afb7921SAidan Hamilton     PetscCall(DMGetCoordinateDM(newDM, &newcdm));
17078afb7921SAidan Hamilton     newCoordnetwork = (DM_Network *)newcdm->data;
17088afb7921SAidan Hamilton     oldCoordnetwork = (DM_Network *)cdm->data;
17098afb7921SAidan Hamilton 
17108afb7921SAidan Hamilton     PetscCall(VecCreate(PETSC_COMM_SELF, &newCoord));
17118afb7921SAidan Hamilton     PetscCall(PetscObjectGetName((PetscObject)oldCoord, &name));
17128afb7921SAidan Hamilton     PetscCall(PetscObjectSetName((PetscObject)newCoord, name));
17138afb7921SAidan Hamilton     PetscCall(VecGetBlockSize(oldCoord, &bs));
17148afb7921SAidan Hamilton     PetscCall(VecSetBlockSize(newCoord, bs));
17158afb7921SAidan Hamilton 
17168afb7921SAidan Hamilton     PetscCall(DMPlexDistributeField(newDMnetwork->plex, migrationSF, oldCoordnetwork->DofSection, oldCoord, newCoordnetwork->DofSection, newCoord));
17178afb7921SAidan Hamilton     PetscCall(DMSetCoordinatesLocal(newDM, newCoord));
17188afb7921SAidan Hamilton 
17198afb7921SAidan Hamilton     PetscCall(VecDestroy(&newCoord));
1720aaa8cc7dSPierre Jolivet     /* Migrate the components from the original coordinate network to the new coordinate network */
17218afb7921SAidan Hamilton     PetscCall(DMPlexDistributeData(newDMnetwork->plex, migrationSF, oldCoordnetwork->DataSection, MPIU_INT, (void *)oldCoordnetwork->componentdataarray, newCoordnetwork->DataSection, (void **)&newCoordnetwork->componentdataarray));
17228afb7921SAidan Hamilton     /* update the header pointers in the new coordinate network components */
17238afb7921SAidan Hamilton     PetscCall(PetscSectionGetChart(newCoordnetwork->DataSection, &pStart, &pEnd));
17248afb7921SAidan Hamilton     for (p = pStart; p < pEnd; p++) {
17258afb7921SAidan Hamilton       PetscCall(PetscSectionGetOffset(newCoordnetwork->DataSection, p, &offset));
17268afb7921SAidan Hamilton       header = (DMNetworkComponentHeader)(newCoordnetwork->componentdataarray + offset);
17278afb7921SAidan Hamilton       /* Update pointers */
17288afb7921SAidan Hamilton       header->size         = (PetscInt *)(header + 1);
17298afb7921SAidan Hamilton       header->key          = header->size + header->maxcomps;
17308afb7921SAidan Hamilton       header->offset       = header->key + header->maxcomps;
17318afb7921SAidan Hamilton       header->nvar         = header->offset + header->maxcomps;
17328afb7921SAidan Hamilton       header->offsetvarrel = header->nvar + header->maxcomps;
17338afb7921SAidan Hamilton     }
17348afb7921SAidan Hamilton 
17358afb7921SAidan Hamilton     PetscCall(DMSetLocalSection(newCoordnetwork->plex, newCoordnetwork->DofSection));
17368afb7921SAidan Hamilton     PetscCall(DMGetGlobalSection(newCoordnetwork->plex, &newCoordnetwork->GlobalDofSection));
17378afb7921SAidan Hamilton     newCoordnetwork->componentsetup = PETSC_TRUE;
17388afb7921SAidan Hamilton   }
17393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17408afb7921SAidan Hamilton }
17418afb7921SAidan Hamilton 
17425f2c45f1SShri Abhyankar /*@
17432bf73ac6SHong Zhang   DMNetworkDistribute - Distributes the network and moves associated component data
17445f2c45f1SShri Abhyankar 
17455f2c45f1SShri Abhyankar   Collective
17465f2c45f1SShri Abhyankar 
1747d8d19677SJose E. Roman   Input Parameters:
1748*60225df5SJacob Faibussowitsch + dm      - the `DMNETWORK` object
17492bf73ac6SHong Zhang - overlap - the overlap of partitions, 0 is the default
17505f2c45f1SShri Abhyankar 
175120f4b53cSBarry Smith   Options Database Keys:
17526e4289a0SDuncan Campbell + -dmnetwork_view              - Calls `DMView()` at the conclusion of `DMSetUp()`
17536e4289a0SDuncan Campbell . -dmnetwork_view_distributed  - Calls `DMView()` at the conclusion of `DMNetworkDistribute()`
17545f25b224SDuncan Campbell . -dmnetwork_view_tmpdir       - Sets the temporary directory to use when viewing with the `draw` option
17555f25b224SDuncan Campbell . -dmnetwork_view_all_ranks    - Displays all of the subnetworks for each MPI rank
17565f25b224SDuncan Campbell . -dmnetwork_view_rank_range   - Displays the subnetworks for the ranks in a comma-separated list
17575f25b224SDuncan Campbell . -dmnetwork_view_no_vertices  - Disables displaying the vertices in the network visualization
17585f25b224SDuncan Campbell - -dmnetwork_view_no_numbering - Disables displaying the numbering of edges and vertices in the network visualization
17595b3f975aSHong Zhang 
17605f2c45f1SShri Abhyankar   Level: intermediate
17615f2c45f1SShri Abhyankar 
176220f4b53cSBarry Smith   Note:
176320f4b53cSBarry Smith   Distributes the network with <overlap>-overlapping partitioning of the edges.
176420f4b53cSBarry Smith 
176520f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkCreate()`
17665f2c45f1SShri Abhyankar @*/
1767d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkDistribute(DM *dm, PetscInt overlap)
1768d71ae5a4SJacob Faibussowitsch {
1769d3464fd4SAdrian Maldonado   MPI_Comm                 comm;
1770d3464fd4SAdrian Maldonado   PetscMPIInt              size;
17718afb7921SAidan Hamilton   DM_Network              *oldDMnetwork = (DM_Network *)((*dm)->data), *newDMnetwork;
1772caf410d2SHong Zhang   PetscSF                  pointsf      = NULL;
17735f2c45f1SShri Abhyankar   DM                       newDM;
17748afb7921SAidan Hamilton   PetscInt                 j, e, v, offset, *subnetvtx, *subnetedge, Nsubnet, gidx, svtx_idx, nv, net;
17758afb7921SAidan Hamilton   PetscInt                *sv;
1776f11a936eSBarry Smith   PetscBT                  btable;
177751ac5effSHong Zhang   PetscPartitioner         part;
1778b9c6e19dSShri Abhyankar   DMNetworkComponentHeader header;
17795f2c45f1SShri Abhyankar 
17805f2c45f1SShri Abhyankar   PetscFunctionBegin;
17815c6496baSHong Zhang   PetscValidPointer(dm, 1);
17825c6496baSHong Zhang   PetscValidHeaderSpecific(*dm, DM_CLASSID, 1);
17839566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)*dm, &comm));
17849566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1785a9750deaSHong Zhang   if (size == 1) {
1786a9750deaSHong Zhang     oldDMnetwork->cloneshared->distributecalled = PETSC_TRUE;
1787a9750deaSHong Zhang     PetscFunctionReturn(PETSC_SUCCESS);
1788a9750deaSHong Zhang   }
1789d3464fd4SAdrian Maldonado 
17905c6496baSHong Zhang   PetscCheck(!overlap, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "overlap %" PetscInt_FMT " != 0 is not supported yet", overlap);
17915b3f975aSHong Zhang 
17922bf73ac6SHong 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. */
1793e3b1a053SGetnet Betrie   PetscCall(PetscLogEventBegin(DMNetwork_Distribute, dm, 0, 0, 0));
17949566063dSJacob Faibussowitsch   PetscCall(DMNetworkCreate(PetscObjectComm((PetscObject)*dm), &newDM));
17955f2c45f1SShri Abhyankar   newDMnetwork                       = (DM_Network *)newDM->data;
179654dfd506SHong Zhang   newDMnetwork->max_comps_registered = oldDMnetwork->max_comps_registered;
17979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(newDMnetwork->max_comps_registered, &newDMnetwork->component));
179851ac5effSHong Zhang 
179951ac5effSHong Zhang   /* Enable runtime options for petscpartitioner */
18009566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPartitioner(oldDMnetwork->plex, &part));
18019566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetFromOptions(part));
180251ac5effSHong Zhang 
18032bf73ac6SHong Zhang   /* Distribute plex dm */
18049566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(oldDMnetwork->plex, overlap, &pointsf, &newDMnetwork->plex));
180551ac5effSHong Zhang 
18065f2c45f1SShri Abhyankar   /* Distribute dof section */
18079566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &newDMnetwork->DofSection));
18089566063dSJacob Faibussowitsch   PetscCall(PetscSFDistributeSection(pointsf, oldDMnetwork->DofSection, NULL, newDMnetwork->DofSection));
180951ac5effSHong Zhang 
18105f2c45f1SShri Abhyankar   /* Distribute data and associated section */
18119566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &newDMnetwork->DataSection));
18129566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeData(newDMnetwork->plex, pointsf, oldDMnetwork->DataSection, MPIU_INT, (void *)oldDMnetwork->componentdataarray, newDMnetwork->DataSection, (void **)&newDMnetwork->componentdataarray));
181324121865SAdrian Maldonado 
1814daad07d3SAidan Hamilton   PetscCall(PetscSectionGetChart(newDMnetwork->DataSection, &newDMnetwork->cloneshared->pStart, &newDMnetwork->cloneshared->pEnd));
1815daad07d3SAidan Hamilton   PetscCall(DMPlexGetHeightStratum(newDMnetwork->plex, 0, &newDMnetwork->cloneshared->eStart, &newDMnetwork->cloneshared->eEnd));
1816daad07d3SAidan Hamilton   PetscCall(DMPlexGetHeightStratum(newDMnetwork->plex, 1, &newDMnetwork->cloneshared->vStart, &newDMnetwork->cloneshared->vEnd));
1817daad07d3SAidan Hamilton   newDMnetwork->cloneshared->nEdges    = newDMnetwork->cloneshared->eEnd - newDMnetwork->cloneshared->eStart;
1818daad07d3SAidan Hamilton   newDMnetwork->cloneshared->nVertices = newDMnetwork->cloneshared->vEnd - newDMnetwork->cloneshared->vStart;
1819daad07d3SAidan Hamilton   newDMnetwork->cloneshared->NVertices = oldDMnetwork->cloneshared->NVertices;
1820daad07d3SAidan Hamilton   newDMnetwork->cloneshared->NEdges    = oldDMnetwork->cloneshared->NEdges;
1821daad07d3SAidan Hamilton   newDMnetwork->cloneshared->svtable   = oldDMnetwork->cloneshared->svtable; /* global table! */
1822daad07d3SAidan Hamilton   oldDMnetwork->cloneshared->svtable   = NULL;
182324121865SAdrian Maldonado 
18241bb6d2a8SBarry Smith   /* Set Dof section as the section for dm */
18259566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(newDMnetwork->plex, newDMnetwork->DofSection));
18269566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(newDMnetwork->plex, &newDMnetwork->GlobalDofSection));
18275f2c45f1SShri Abhyankar 
1828b9c6e19dSShri Abhyankar   /* Setup subnetwork info in the newDM */
1829daad07d3SAidan Hamilton   newDMnetwork->cloneshared->Nsubnet = oldDMnetwork->cloneshared->Nsubnet;
1830daad07d3SAidan Hamilton   newDMnetwork->cloneshared->Nsvtx   = oldDMnetwork->cloneshared->Nsvtx;
1831daad07d3SAidan Hamilton   oldDMnetwork->cloneshared->Nsvtx   = 0;
1832daad07d3SAidan Hamilton   newDMnetwork->cloneshared->svtx    = oldDMnetwork->cloneshared->svtx; /* global vertices! */
1833daad07d3SAidan Hamilton   oldDMnetwork->cloneshared->svtx    = NULL;
1834daad07d3SAidan Hamilton   PetscCall(PetscCalloc1(newDMnetwork->cloneshared->Nsubnet, &newDMnetwork->cloneshared->subnet));
18352bf73ac6SHong Zhang 
18362bf73ac6SHong Zhang   /* Copy over the global number of vertices and edges in each subnetwork.
18372bf73ac6SHong Zhang      Note: these are calculated in DMNetworkLayoutSetUp()
1838b9c6e19dSShri Abhyankar   */
1839daad07d3SAidan Hamilton   Nsubnet = newDMnetwork->cloneshared->Nsubnet;
1840f11a936eSBarry Smith   for (j = 0; j < Nsubnet; j++) {
1841daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[j].Nvtx  = oldDMnetwork->cloneshared->subnet[j].Nvtx;
1842daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[j].Nedge = oldDMnetwork->cloneshared->subnet[j].Nedge;
1843b9c6e19dSShri Abhyankar   }
1844b9c6e19dSShri Abhyankar 
1845f11a936eSBarry Smith   /* Count local nedges for subnetworks */
1846daad07d3SAidan Hamilton   for (e = newDMnetwork->cloneshared->eStart; e < newDMnetwork->cloneshared->eEnd; e++) {
18479566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection, e, &offset));
184854dfd506SHong Zhang     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);
18495c6496baSHong Zhang 
185054dfd506SHong Zhang     /* Update pointers */
185154dfd506SHong Zhang     header->size         = (PetscInt *)(header + 1);
185254dfd506SHong Zhang     header->key          = header->size + header->maxcomps;
185354dfd506SHong Zhang     header->offset       = header->key + header->maxcomps;
185454dfd506SHong Zhang     header->nvar         = header->offset + header->maxcomps;
185554dfd506SHong Zhang     header->offsetvarrel = header->nvar + header->maxcomps;
185654dfd506SHong Zhang 
1857daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[header->subnetid].nedge++;
1858b9c6e19dSShri Abhyankar   }
1859b9c6e19dSShri Abhyankar 
18605c6496baSHong Zhang   /* Setup a btable to keep track subnetworks owned by this process at a shared vertex */
186148a46eb9SPierre Jolivet   if (newDMnetwork->cloneshared->Nsvtx) PetscCall(PetscBTCreate(Nsubnet, &btable));
18625c6496baSHong Zhang 
18635c6496baSHong Zhang   /* Count local nvtx for subnetworks */
1864daad07d3SAidan Hamilton   for (v = newDMnetwork->cloneshared->vStart; v < newDMnetwork->cloneshared->vEnd; v++) {
18659566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection, v, &offset));
18665f80ce2aSJacob Faibussowitsch     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);
18675c6496baSHong Zhang 
186854dfd506SHong Zhang     /* Update pointers */
186954dfd506SHong Zhang     header->size         = (PetscInt *)(header + 1);
187054dfd506SHong Zhang     header->key          = header->size + header->maxcomps;
187154dfd506SHong Zhang     header->offset       = header->key + header->maxcomps;
187254dfd506SHong Zhang     header->nvar         = header->offset + header->maxcomps;
187354dfd506SHong Zhang     header->offsetvarrel = header->nvar + header->maxcomps;
187454dfd506SHong Zhang 
18752bf73ac6SHong Zhang     /* shared vertices: use gidx=header->index to check if v is a shared vertex */
18762bf73ac6SHong Zhang     gidx = header->index;
1877eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIGetWithDefault(newDMnetwork->cloneshared->svtable, gidx + 1, 0, &svtx_idx));
18782bf73ac6SHong Zhang     svtx_idx--;
18792bf73ac6SHong Zhang 
18802bf73ac6SHong Zhang     if (svtx_idx < 0) { /* not a shared vertex */
1881daad07d3SAidan Hamilton       newDMnetwork->cloneshared->subnet[header->subnetid].nvtx++;
18822bf73ac6SHong Zhang     } else { /* a shared vertex belongs to more than one subnetworks, it is being counted by multiple subnets */
18835c6496baSHong Zhang       /* Setup a lookup btable for this v's owning subnetworks */
18849566063dSJacob Faibussowitsch       PetscCall(SetSubnetIdLookupBT(newDM, v, Nsubnet, btable));
1885f11a936eSBarry Smith 
1886daad07d3SAidan Hamilton       for (j = 0; j < newDMnetwork->cloneshared->svtx[svtx_idx].n; j++) {
1887daad07d3SAidan Hamilton         sv  = newDMnetwork->cloneshared->svtx[svtx_idx].sv + 2 * j;
18885c6496baSHong Zhang         net = sv[0];
188935cb6cd3SPierre Jolivet         if (PetscBTLookup(btable, net)) newDMnetwork->cloneshared->subnet[net].nvtx++; /* sv is on net owned by this process */
18902bf73ac6SHong Zhang       }
18912bf73ac6SHong Zhang     }
1892b9c6e19dSShri Abhyankar   }
1893b9c6e19dSShri Abhyankar 
18942bf73ac6SHong Zhang   /* Get total local nvtx for subnetworks */
18952bf73ac6SHong Zhang   nv = 0;
1896daad07d3SAidan Hamilton   for (j = 0; j < Nsubnet; j++) nv += newDMnetwork->cloneshared->subnet[j].nvtx;
1897daad07d3SAidan Hamilton   nv += newDMnetwork->cloneshared->Nsvtx;
1898caf410d2SHong Zhang 
18992bf73ac6SHong Zhang   /* Now create the vertices and edge arrays for the subnetworks */
1900daad07d3SAidan Hamilton   PetscCall(PetscCalloc2(newDMnetwork->cloneshared->nEdges, &subnetedge, nv, &subnetvtx)); /* Maps local vertex to local subnetwork's vertex */
1901daad07d3SAidan Hamilton   newDMnetwork->cloneshared->subnetedge = subnetedge;
1902daad07d3SAidan Hamilton   newDMnetwork->cloneshared->subnetvtx  = subnetvtx;
1903daad07d3SAidan Hamilton   for (j = 0; j < newDMnetwork->cloneshared->Nsubnet; j++) {
1904daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[j].edges = subnetedge;
1905daad07d3SAidan Hamilton     subnetedge += newDMnetwork->cloneshared->subnet[j].nedge;
19062bf73ac6SHong Zhang 
1907daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[j].vertices = subnetvtx;
1908daad07d3SAidan Hamilton     subnetvtx += newDMnetwork->cloneshared->subnet[j].nvtx;
1909caf410d2SHong Zhang 
19102bf73ac6SHong 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. */
1911daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[j].nvtx = newDMnetwork->cloneshared->subnet[j].nedge = 0;
1912b9c6e19dSShri Abhyankar   }
1913daad07d3SAidan Hamilton   newDMnetwork->cloneshared->svertices = subnetvtx;
1914b9c6e19dSShri Abhyankar 
19152bf73ac6SHong Zhang   /* Set the edges and vertices in each subnetwork */
1916daad07d3SAidan Hamilton   for (e = newDMnetwork->cloneshared->eStart; e < newDMnetwork->cloneshared->eEnd; e++) {
19179566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection, e, &offset));
19185f80ce2aSJacob Faibussowitsch     header                                                                                                                 = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);
1919daad07d3SAidan Hamilton     newDMnetwork->cloneshared->subnet[header->subnetid].edges[newDMnetwork->cloneshared->subnet[header->subnetid].nedge++] = e;
1920b9c6e19dSShri Abhyankar   }
1921b9c6e19dSShri Abhyankar 
19222bf73ac6SHong Zhang   nv = 0;
1923daad07d3SAidan Hamilton   for (v = newDMnetwork->cloneshared->vStart; v < newDMnetwork->cloneshared->vEnd; v++) {
19249566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection, v, &offset));
19255f80ce2aSJacob Faibussowitsch     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);
19262bf73ac6SHong Zhang 
19272bf73ac6SHong Zhang     /* coupling vertices: use gidx = header->index to check if v is a coupling vertex */
1928eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIGetWithDefault(newDMnetwork->cloneshared->svtable, header->index + 1, 0, &svtx_idx));
19292bf73ac6SHong Zhang     svtx_idx--;
19302bf73ac6SHong Zhang     if (svtx_idx < 0) {
1931daad07d3SAidan Hamilton       newDMnetwork->cloneshared->subnet[header->subnetid].vertices[newDMnetwork->cloneshared->subnet[header->subnetid].nvtx++] = v;
19322bf73ac6SHong Zhang     } else { /* a shared vertex */
1933daad07d3SAidan Hamilton       newDMnetwork->cloneshared->svertices[nv++] = v;
19342bf73ac6SHong Zhang 
19355c6496baSHong Zhang       /* Setup a lookup btable for this v's owning subnetworks */
19369566063dSJacob Faibussowitsch       PetscCall(SetSubnetIdLookupBT(newDM, v, Nsubnet, btable));
1937f11a936eSBarry Smith 
1938daad07d3SAidan Hamilton       for (j = 0; j < newDMnetwork->cloneshared->svtx[svtx_idx].n; j++) {
1939daad07d3SAidan Hamilton         sv  = newDMnetwork->cloneshared->svtx[svtx_idx].sv + 2 * j;
19405c6496baSHong Zhang         net = sv[0];
19419371c9d4SSatish Balay         if (PetscBTLookup(btable, net)) newDMnetwork->cloneshared->subnet[net].vertices[newDMnetwork->cloneshared->subnet[net].nvtx++] = v;
1942b9c6e19dSShri Abhyankar       }
19432bf73ac6SHong Zhang     }
19442bf73ac6SHong Zhang   }
1945daad07d3SAidan Hamilton   newDMnetwork->cloneshared->nsvtx = nv; /* num of local shared vertices */
1946b9c6e19dSShri Abhyankar 
19478afb7921SAidan Hamilton   PetscCall(DMNetworkDistributeCoordinates(*dm, pointsf, newDM));
1948caf410d2SHong Zhang   newDM->setupcalled                          = (*dm)->setupcalled;
1949daad07d3SAidan Hamilton   newDMnetwork->cloneshared->distributecalled = PETSC_TRUE;
1950caf410d2SHong Zhang 
19512bf73ac6SHong Zhang   /* Free spaces */
19529566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&pointsf));
19539566063dSJacob Faibussowitsch   PetscCall(DMDestroy(dm));
195448a46eb9SPierre Jolivet   if (newDMnetwork->cloneshared->Nsvtx) PetscCall(PetscBTDestroy(&btable));
19552bf73ac6SHong Zhang 
19565b3f975aSHong Zhang   /* View distributed dmnetwork */
19579566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(newDM, NULL, "-dmnetwork_view_distributed"));
19585b3f975aSHong Zhang 
1959d3464fd4SAdrian Maldonado   *dm = newDM;
1960e3b1a053SGetnet Betrie   PetscCall(PetscLogEventEnd(DMNetwork_Distribute, dm, 0, 0, 0));
19613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19625f2c45f1SShri Abhyankar }
19635f2c45f1SShri Abhyankar 
196424121865SAdrian Maldonado /*@C
196520f4b53cSBarry Smith   PetscSFGetSubSF - Returns an `PetscSF` for a specific subset of points. Leaves are re-numbered to reflect the new ordering
19662bf73ac6SHong Zhang 
19672bf73ac6SHong Zhang   Collective
196824121865SAdrian Maldonado 
196924121865SAdrian Maldonado   Input Parameters:
1970*60225df5SJacob Faibussowitsch + mainsf - `PetscSF` structure
197120f4b53cSBarry Smith - map    - a `ISLocalToGlobalMapping` that contains the subset of points
197224121865SAdrian Maldonado 
197397bb3fdcSJose E. Roman   Output Parameter:
197420f4b53cSBarry Smith . subSF - a subset of the `mainSF` for the desired subset.
1975ee300463SSatish Balay 
1976ee300463SSatish Balay   Level: intermediate
197720f4b53cSBarry Smith 
197820f4b53cSBarry Smith .seealso: `PetscSF`
19797d928bffSSatish Balay @*/
1980d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSFGetSubSF(PetscSF mainsf, ISLocalToGlobalMapping map, PetscSF *subSF)
1981d71ae5a4SJacob Faibussowitsch {
198224121865SAdrian Maldonado   PetscInt           nroots, nleaves, *ilocal_sub;
198324121865SAdrian Maldonado   PetscInt           i, *ilocal_map, nroots_sub, nleaves_sub = 0;
198424121865SAdrian Maldonado   PetscInt          *local_points, *remote_points;
198524121865SAdrian Maldonado   PetscSFNode       *iremote_sub;
198624121865SAdrian Maldonado   const PetscInt    *ilocal;
198724121865SAdrian Maldonado   const PetscSFNode *iremote;
198824121865SAdrian Maldonado 
198924121865SAdrian Maldonado   PetscFunctionBegin;
19909566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(mainsf, &nroots, &nleaves, &ilocal, &iremote));
199124121865SAdrian Maldonado 
199224121865SAdrian Maldonado   /* Look for leaves that pertain to the subset of points. Get the local ordering */
19939566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &ilocal_map));
19949566063dSJacob Faibussowitsch   PetscCall(ISGlobalToLocalMappingApply(map, IS_GTOLM_MASK, nleaves, ilocal, NULL, ilocal_map));
199524121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
199624121865SAdrian Maldonado     if (ilocal_map[i] != -1) nleaves_sub += 1;
199724121865SAdrian Maldonado   }
199824121865SAdrian Maldonado   /* Re-number ilocal with subset numbering. Need information from roots */
19999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nroots, &local_points, nroots, &remote_points));
200024121865SAdrian Maldonado   for (i = 0; i < nroots; i++) local_points[i] = i;
20019566063dSJacob Faibussowitsch   PetscCall(ISGlobalToLocalMappingApply(map, IS_GTOLM_MASK, nroots, local_points, NULL, local_points));
20029566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(mainsf, MPIU_INT, local_points, remote_points, MPI_REPLACE));
20039566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(mainsf, MPIU_INT, local_points, remote_points, MPI_REPLACE));
200424121865SAdrian Maldonado   /* Fill up graph using local (that is, local to the subset) numbering. */
20059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves_sub, &ilocal_sub));
20069566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves_sub, &iremote_sub));
200724121865SAdrian Maldonado   nleaves_sub = 0;
200824121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
200924121865SAdrian Maldonado     if (ilocal_map[i] != -1) {
201024121865SAdrian Maldonado       ilocal_sub[nleaves_sub]        = ilocal_map[i];
20114b70a8deSAdrian Maldonado       iremote_sub[nleaves_sub].rank  = iremote[i].rank;
201224121865SAdrian Maldonado       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
201324121865SAdrian Maldonado       nleaves_sub += 1;
201424121865SAdrian Maldonado     }
201524121865SAdrian Maldonado   }
20169566063dSJacob Faibussowitsch   PetscCall(PetscFree2(local_points, remote_points));
20179566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetSize(map, &nroots_sub));
201824121865SAdrian Maldonado 
201924121865SAdrian Maldonado   /* Create new subSF */
2020e66d8692SHong Zhang   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)mainsf), subSF));
20219566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(*subSF));
20229566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*subSF, nroots_sub, nleaves_sub, ilocal_sub, PETSC_OWN_POINTER, iremote_sub, PETSC_COPY_VALUES));
20239566063dSJacob Faibussowitsch   PetscCall(PetscFree(ilocal_map));
20249566063dSJacob Faibussowitsch   PetscCall(PetscFree(iremote_sub));
20253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
202624121865SAdrian Maldonado }
202724121865SAdrian Maldonado 
20285f2c45f1SShri Abhyankar /*@C
20295f2c45f1SShri Abhyankar   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
20305f2c45f1SShri Abhyankar 
20315f2c45f1SShri Abhyankar   Not Collective
20325f2c45f1SShri Abhyankar 
20335f2c45f1SShri Abhyankar   Input Parameters:
203420f4b53cSBarry Smith + dm     - the `DMNETWORK` object
2035*60225df5SJacob Faibussowitsch - vertex - the vertex point
20365f2c45f1SShri Abhyankar 
2037fd292e60Sprj-   Output Parameters:
20385f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point
20392bf73ac6SHong Zhang - edges  - list of edge points
20405f2c45f1SShri Abhyankar 
204197bb938eSShri Abhyankar   Level: beginner
20425f2c45f1SShri Abhyankar 
204320f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkCreate()`, `DMNetworkGetConnectedVertices()`
20445f2c45f1SShri Abhyankar @*/
2045d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetSupportingEdges(DM dm, PetscInt vertex, PetscInt *nedges, const PetscInt *edges[])
2046d71ae5a4SJacob Faibussowitsch {
20475f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
20485f2c45f1SShri Abhyankar 
20495f2c45f1SShri Abhyankar   PetscFunctionBegin;
20509566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSupportSize(network->plex, vertex, nedges));
20519566063dSJacob Faibussowitsch   if (edges) PetscCall(DMPlexGetSupport(network->plex, vertex, edges));
20523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20535f2c45f1SShri Abhyankar }
20545f2c45f1SShri Abhyankar 
20555f2c45f1SShri Abhyankar /*@C
2056d842c372SHong Zhang   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
20575f2c45f1SShri Abhyankar 
20585f2c45f1SShri Abhyankar   Not Collective
20595f2c45f1SShri Abhyankar 
20605f2c45f1SShri Abhyankar   Input Parameters:
206120f4b53cSBarry Smith + dm   - the `DMNETWORK` object
2062*60225df5SJacob Faibussowitsch - edge - the edge point
20635f2c45f1SShri Abhyankar 
20642fe279fdSBarry Smith   Output Parameter:
20655f2c45f1SShri Abhyankar . vertices - vertices connected to this edge
20665f2c45f1SShri Abhyankar 
206797bb938eSShri Abhyankar   Level: beginner
20685f2c45f1SShri Abhyankar 
206920f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkCreate()`, `DMNetworkGetSupportingEdges()`
20705f2c45f1SShri Abhyankar @*/
2071d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetConnectedVertices(DM dm, PetscInt edge, const PetscInt *vertices[])
2072d71ae5a4SJacob Faibussowitsch {
20735f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
20745f2c45f1SShri Abhyankar 
20755f2c45f1SShri Abhyankar   PetscFunctionBegin;
20769566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(network->plex, edge, vertices));
20773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20785f2c45f1SShri Abhyankar }
20795f2c45f1SShri Abhyankar 
20805f2c45f1SShri Abhyankar /*@
208120f4b53cSBarry Smith   DMNetworkIsSharedVertex - Returns `PETSC_TRUE` if the vertex is shared by subnetworks
20822bf73ac6SHong Zhang 
20832bf73ac6SHong Zhang   Not Collective
20842bf73ac6SHong Zhang 
20852bf73ac6SHong Zhang   Input Parameters:
208620f4b53cSBarry Smith + dm - the `DMNETWORK` object
20872bf73ac6SHong Zhang - p  - the vertex point
20882bf73ac6SHong Zhang 
20892bf73ac6SHong Zhang   Output Parameter:
209020f4b53cSBarry Smith . flag - `PETSC_TRUE` if the vertex is shared by subnetworks
20912bf73ac6SHong Zhang 
20922bf73ac6SHong Zhang   Level: beginner
20932bf73ac6SHong Zhang 
209420f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkAddSharedVertices()`, `DMNetworkIsGhostVertex()`
20952bf73ac6SHong Zhang @*/
2096d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkIsSharedVertex(DM dm, PetscInt p, PetscBool *flag)
2097d71ae5a4SJacob Faibussowitsch {
20982bf73ac6SHong Zhang   PetscFunctionBegin;
2099eec179cfSJacob Faibussowitsch   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2100eec179cfSJacob Faibussowitsch   PetscValidBoolPointer(flag, 3);
21012bf73ac6SHong Zhang   if (dm->setupcalled) { /* DMNetworkGetGlobalVertexIndex() requires DMSetUp() be called */
21022bf73ac6SHong Zhang     DM_Network *network = (DM_Network *)dm->data;
21032bf73ac6SHong Zhang     PetscInt    gidx;
2104eec179cfSJacob Faibussowitsch 
21059566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVertexIndex(dm, p, &gidx));
2106eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIHas(network->cloneshared->svtable, gidx + 1, flag));
21072bf73ac6SHong Zhang   } else { /* would be removed? */
21082bf73ac6SHong Zhang     PetscInt        nv;
21092bf73ac6SHong Zhang     const PetscInt *vtx;
2110eec179cfSJacob Faibussowitsch 
21119566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetSharedVertices(dm, &nv, &vtx));
2112eec179cfSJacob Faibussowitsch     for (PetscInt i = 0; i < nv; i++) {
21132bf73ac6SHong Zhang       if (p == vtx[i]) {
21142bf73ac6SHong Zhang         *flag = PETSC_TRUE;
21153ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
21162bf73ac6SHong Zhang       }
21172bf73ac6SHong Zhang     }
2118eec179cfSJacob Faibussowitsch     *flag = PETSC_FALSE;
21192bf73ac6SHong Zhang   }
21203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21212bf73ac6SHong Zhang }
21222bf73ac6SHong Zhang 
21232bf73ac6SHong Zhang /*@
212420f4b53cSBarry Smith   DMNetworkIsGhostVertex - Returns `PETSC_TRUE` if the vertex is a ghost vertex
21255f2c45f1SShri Abhyankar 
21265f2c45f1SShri Abhyankar   Not Collective
21275f2c45f1SShri Abhyankar 
21285f2c45f1SShri Abhyankar   Input Parameters:
212920f4b53cSBarry Smith + dm - the `DMNETWORK` object
2130a2b725a8SWilliam Gropp - p  - the vertex point
21315f2c45f1SShri Abhyankar 
21325f2c45f1SShri Abhyankar   Output Parameter:
213320f4b53cSBarry Smith . isghost - `PETSC_TRUE` if the vertex is a ghost point
21345f2c45f1SShri Abhyankar 
213597bb938eSShri Abhyankar   Level: beginner
21365f2c45f1SShri Abhyankar 
213720f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetConnectedVertices()`, `DMNetworkGetVertexRange()`, `DMNetworkIsSharedVertex()`
21385f2c45f1SShri Abhyankar @*/
2139d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkIsGhostVertex(DM dm, PetscInt p, PetscBool *isghost)
2140d71ae5a4SJacob Faibussowitsch {
21415f2c45f1SShri Abhyankar   DM_Network  *network = (DM_Network *)dm->data;
21425f2c45f1SShri Abhyankar   PetscInt     offsetg;
21435f2c45f1SShri Abhyankar   PetscSection sectiong;
21445f2c45f1SShri Abhyankar 
21455f2c45f1SShri Abhyankar   PetscFunctionBegin;
21465f2c45f1SShri Abhyankar   *isghost = PETSC_FALSE;
21479566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(network->plex, &sectiong));
21489566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(sectiong, p, &offsetg));
21495f2c45f1SShri Abhyankar   if (offsetg < 0) *isghost = PETSC_TRUE;
21503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21515f2c45f1SShri Abhyankar }
21525f2c45f1SShri Abhyankar 
2153d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetUp_Network(DM dm)
2154d71ae5a4SJacob Faibussowitsch {
21555f2c45f1SShri Abhyankar   PetscFunctionBegin;
2156e3b1a053SGetnet Betrie   PetscCall(PetscLogEventBegin(DMNetwork_SetUpNetwork, dm, 0, 0, 0));
2157daad07d3SAidan Hamilton   PetscCall(DMNetworkFinalizeComponents(dm));
21585b3f975aSHong Zhang   /* View dmnetwork */
21599566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dmnetwork_view"));
2160e3b1a053SGetnet Betrie   PetscCall(PetscLogEventEnd(DMNetwork_SetUpNetwork, dm, 0, 0, 0));
21613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21625f2c45f1SShri Abhyankar }
21635f2c45f1SShri Abhyankar 
21641ad426b7SHong Zhang /*@
216517df6e9eSHong Zhang   DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
216620f4b53cSBarry Smith   -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TRUE)?
21671ad426b7SHong Zhang 
21681ad426b7SHong Zhang   Collective
21691ad426b7SHong Zhang 
21701ad426b7SHong Zhang   Input Parameters:
217120f4b53cSBarry Smith + dm   - the `DMNETWORK` object
217220f4b53cSBarry Smith . eflg - turn the option on (`PETSC_TRUE`) or off (`PETSC_FALSE`) if user provides Jacobian for edges
217320f4b53cSBarry Smith - vflg - turn the option on (`PETSC_TRUE`) or off (`PETSC_FALSE`) if user provides Jacobian for vertices
21741ad426b7SHong Zhang 
21751ad426b7SHong Zhang   Level: intermediate
21761ad426b7SHong Zhang 
21771ad426b7SHong Zhang @*/
2178d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkHasJacobian(DM dm, PetscBool eflg, PetscBool vflg)
2179d71ae5a4SJacob Faibussowitsch {
21801ad426b7SHong Zhang   DM_Network *network   = (DM_Network *)dm->data;
2181daad07d3SAidan Hamilton   PetscInt    nVertices = network->cloneshared->nVertices;
21821ad426b7SHong Zhang 
21831ad426b7SHong Zhang   PetscFunctionBegin;
218483b2e829SHong Zhang   network->userEdgeJacobian   = eflg;
218583b2e829SHong Zhang   network->userVertexJacobian = vflg;
21868675203cSHong Zhang 
218748a46eb9SPierre Jolivet   if (eflg && !network->Je) PetscCall(PetscCalloc1(3 * network->cloneshared->nEdges, &network->Je));
21888675203cSHong Zhang 
218966b4e0ffSHong Zhang   if (vflg && !network->Jv && nVertices) {
2190daad07d3SAidan Hamilton     PetscInt        i, *vptr, nedges, vStart = network->cloneshared->vStart;
219166b4e0ffSHong Zhang     PetscInt        nedges_total;
21928675203cSHong Zhang     const PetscInt *edges;
21938675203cSHong Zhang 
21948675203cSHong Zhang     /* count nvertex_total */
21958675203cSHong Zhang     nedges_total = 0;
21969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nVertices + 1, &vptr));
21978675203cSHong Zhang 
21988675203cSHong Zhang     vptr[0] = 0;
21998675203cSHong Zhang     for (i = 0; i < nVertices; i++) {
22009566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetSupportingEdges(dm, i + vStart, &nedges, &edges));
22018675203cSHong Zhang       nedges_total += nedges;
22028675203cSHong Zhang       vptr[i + 1] = vptr[i] + 2 * nedges + 1;
22038675203cSHong Zhang     }
22048675203cSHong Zhang 
22059566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(2 * nedges_total + nVertices, &network->Jv));
22068675203cSHong Zhang     network->Jvptr = vptr;
22078675203cSHong Zhang   }
22083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22091ad426b7SHong Zhang }
22101ad426b7SHong Zhang 
22111ad426b7SHong Zhang /*@
221283b2e829SHong Zhang   DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
221383b2e829SHong Zhang 
221483b2e829SHong Zhang   Not Collective
221583b2e829SHong Zhang 
221683b2e829SHong Zhang   Input Parameters:
221720f4b53cSBarry Smith + dm - the `DMNETWORK` object
221883b2e829SHong Zhang . p  - the edge point
22193e97b6e8SHong Zhang - J  - array (size = 3) of Jacobian submatrices for this edge point:
22203e97b6e8SHong Zhang         J[0]: this edge
222120f4b53cSBarry Smith         J[1] and J[2]: connected vertices, obtained by calling `DMNetworkGetConnectedVertices()`
222283b2e829SHong Zhang 
222397bb938eSShri Abhyankar   Level: advanced
222483b2e829SHong Zhang 
222520f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkVertexSetMatrix()`
222683b2e829SHong Zhang @*/
2227d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkEdgeSetMatrix(DM dm, PetscInt p, Mat J[])
2228d71ae5a4SJacob Faibussowitsch {
222983b2e829SHong Zhang   DM_Network *network = (DM_Network *)dm->data;
223083b2e829SHong Zhang 
223183b2e829SHong Zhang   PetscFunctionBegin;
22325c6496baSHong Zhang   PetscCheck(network->Je, PetscObjectComm((PetscObject)dm), PETSC_ERR_ORDER, "Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
22338675203cSHong Zhang 
22348675203cSHong Zhang   if (J) {
2235883e35e8SHong Zhang     network->Je[3 * p]     = J[0];
2236883e35e8SHong Zhang     network->Je[3 * p + 1] = J[1];
2237883e35e8SHong Zhang     network->Je[3 * p + 2] = J[2];
22388675203cSHong Zhang   }
22393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
224083b2e829SHong Zhang }
224183b2e829SHong Zhang 
224283b2e829SHong Zhang /*@
224376ddfea5SHong Zhang   DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
22441ad426b7SHong Zhang 
22451ad426b7SHong Zhang   Not Collective
22461ad426b7SHong Zhang 
22471ad426b7SHong Zhang   Input Parameters:
224820f4b53cSBarry Smith + dm - The `DMNETWORK` object
22491ad426b7SHong Zhang . p  - the vertex point
22503e97b6e8SHong Zhang - J  - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
22513e97b6e8SHong Zhang         J[0]:       this vertex
22523e97b6e8SHong Zhang         J[1+2*i]:   i-th supporting edge
22533e97b6e8SHong Zhang         J[1+2*i+1]: i-th connected vertex
22541ad426b7SHong Zhang 
225597bb938eSShri Abhyankar   Level: advanced
22561ad426b7SHong Zhang 
225720f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkEdgeSetMatrix()`
22581ad426b7SHong Zhang @*/
2259d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkVertexSetMatrix(DM dm, PetscInt p, Mat J[])
2260d71ae5a4SJacob Faibussowitsch {
22615f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network *)dm->data;
2262daad07d3SAidan Hamilton   PetscInt        i, *vptr, nedges, vStart = network->cloneshared->vStart;
2263883e35e8SHong Zhang   const PetscInt *edges;
22645f2c45f1SShri Abhyankar 
22655f2c45f1SShri Abhyankar   PetscFunctionBegin;
22665c6496baSHong Zhang   PetscCheck(network->Jv, PetscObjectComm((PetscObject)dm), PETSC_ERR_ORDER, "Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
2267883e35e8SHong Zhang 
22688675203cSHong Zhang   if (J) {
2269883e35e8SHong Zhang     vptr                          = network->Jvptr;
22703e97b6e8SHong Zhang     network->Jv[vptr[p - vStart]] = J[0]; /* Set Jacobian for this vertex */
22713e97b6e8SHong Zhang 
22723e97b6e8SHong Zhang     /* Set Jacobian for each supporting edge and connected vertex */
22739566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetSupportingEdges(dm, p, &nedges, &edges));
2274883e35e8SHong Zhang     for (i = 1; i <= 2 * nedges; i++) network->Jv[vptr[p - vStart] + i] = J[i];
22758675203cSHong Zhang   }
22763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2277883e35e8SHong Zhang }
2278883e35e8SHong Zhang 
2279d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscBool ghost, Vec vdnz, Vec vonz)
2280d71ae5a4SJacob Faibussowitsch {
22815cf7da58SHong Zhang   PetscInt    j;
22825cf7da58SHong Zhang   PetscScalar val = (PetscScalar)ncols;
22835cf7da58SHong Zhang 
22845cf7da58SHong Zhang   PetscFunctionBegin;
22855cf7da58SHong Zhang   if (!ghost) {
228648a46eb9SPierre Jolivet     for (j = 0; j < nrows; j++) PetscCall(VecSetValues(vdnz, 1, &rows[j], &val, ADD_VALUES));
22875cf7da58SHong Zhang   } else {
228848a46eb9SPierre Jolivet     for (j = 0; j < nrows; j++) PetscCall(VecSetValues(vonz, 1, &rows[j], &val, ADD_VALUES));
22895cf7da58SHong Zhang   }
22903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22915cf7da58SHong Zhang }
22925cf7da58SHong Zhang 
2293d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju, PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscBool ghost, Vec vdnz, Vec vonz)
2294d71ae5a4SJacob Faibussowitsch {
22955cf7da58SHong Zhang   PetscInt    j, ncols_u;
22965cf7da58SHong Zhang   PetscScalar val;
22975cf7da58SHong Zhang 
22985cf7da58SHong Zhang   PetscFunctionBegin;
22995cf7da58SHong Zhang   if (!ghost) {
23005cf7da58SHong Zhang     for (j = 0; j < nrows; j++) {
23019566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Ju, j, &ncols_u, NULL, NULL));
23025cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
23039566063dSJacob Faibussowitsch       PetscCall(VecSetValues(vdnz, 1, &rows[j], &val, ADD_VALUES));
23049566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Ju, j, &ncols_u, NULL, NULL));
23055cf7da58SHong Zhang     }
23065cf7da58SHong Zhang   } else {
23075cf7da58SHong Zhang     for (j = 0; j < nrows; j++) {
23089566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Ju, j, &ncols_u, NULL, NULL));
23095cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
23109566063dSJacob Faibussowitsch       PetscCall(VecSetValues(vonz, 1, &rows[j], &val, ADD_VALUES));
23119566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Ju, j, &ncols_u, NULL, NULL));
23125cf7da58SHong Zhang     }
23135cf7da58SHong Zhang   }
23143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23155cf7da58SHong Zhang }
23165cf7da58SHong Zhang 
2317d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode MatSetPreallocationblock_private(Mat Ju, PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscBool ghost, Vec vdnz, Vec vonz)
2318d71ae5a4SJacob Faibussowitsch {
23195cf7da58SHong Zhang   PetscFunctionBegin;
23205cf7da58SHong Zhang   if (Ju) {
23219566063dSJacob Faibussowitsch     PetscCall(MatSetPreallocationUserblock_private(Ju, nrows, rows, ncols, ghost, vdnz, vonz));
23225cf7da58SHong Zhang   } else {
23239566063dSJacob Faibussowitsch     PetscCall(MatSetPreallocationDenseblock_private(nrows, rows, ncols, ghost, vdnz, vonz));
23245cf7da58SHong Zhang   }
23253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23265cf7da58SHong Zhang }
23275cf7da58SHong Zhang 
2328d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode MatSetDenseblock_private(PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscInt cstart, Mat *J)
2329d71ae5a4SJacob Faibussowitsch {
2330883e35e8SHong Zhang   PetscInt     j, *cols;
2331883e35e8SHong Zhang   PetscScalar *zeros;
2332883e35e8SHong Zhang 
2333883e35e8SHong Zhang   PetscFunctionBegin;
23349566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(ncols, &cols, nrows * ncols, &zeros));
2335883e35e8SHong Zhang   for (j = 0; j < ncols; j++) cols[j] = j + cstart;
23369566063dSJacob Faibussowitsch   PetscCall(MatSetValues(*J, nrows, rows, ncols, cols, zeros, INSERT_VALUES));
23379566063dSJacob Faibussowitsch   PetscCall(PetscFree2(cols, zeros));
23383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23391ad426b7SHong Zhang }
2340a4e85ca8SHong Zhang 
2341d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode MatSetUserblock_private(Mat Ju, PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscInt cstart, Mat *J)
2342d71ae5a4SJacob Faibussowitsch {
23433e97b6e8SHong Zhang   PetscInt        j, M, N, row, col, ncols_u;
23443e97b6e8SHong Zhang   const PetscInt *cols;
23453e97b6e8SHong Zhang   PetscScalar     zero = 0.0;
23463e97b6e8SHong Zhang 
23473e97b6e8SHong Zhang   PetscFunctionBegin;
23489566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Ju, &M, &N));
234963a3b9bcSJacob 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);
23503e97b6e8SHong Zhang 
23513e97b6e8SHong Zhang   for (row = 0; row < nrows; row++) {
23529566063dSJacob Faibussowitsch     PetscCall(MatGetRow(Ju, row, &ncols_u, &cols, NULL));
23533e97b6e8SHong Zhang     for (j = 0; j < ncols_u; j++) {
23543e97b6e8SHong Zhang       col = cols[j] + cstart;
23559566063dSJacob Faibussowitsch       PetscCall(MatSetValues(*J, 1, &rows[row], 1, &col, &zero, INSERT_VALUES));
23563e97b6e8SHong Zhang     }
23579566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(Ju, row, &ncols_u, &cols, NULL));
23583e97b6e8SHong Zhang   }
23593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23603e97b6e8SHong Zhang }
23611ad426b7SHong Zhang 
2362d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode MatSetblock_private(Mat Ju, PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscInt cstart, Mat *J)
2363d71ae5a4SJacob Faibussowitsch {
2364a4e85ca8SHong Zhang   PetscFunctionBegin;
2365a4e85ca8SHong Zhang   if (Ju) {
23669566063dSJacob Faibussowitsch     PetscCall(MatSetUserblock_private(Ju, nrows, rows, ncols, cstart, J));
2367a4e85ca8SHong Zhang   } else {
23689566063dSJacob Faibussowitsch     PetscCall(MatSetDenseblock_private(nrows, rows, ncols, cstart, J));
2369a4e85ca8SHong Zhang   }
23703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2371a4e85ca8SHong Zhang }
2372a4e85ca8SHong Zhang 
237324121865SAdrian 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.
237424121865SAdrian Maldonado */
2375d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
2376d71ae5a4SJacob Faibussowitsch {
237724121865SAdrian Maldonado   PetscInt  i, size, dof;
237824121865SAdrian Maldonado   PetscInt *glob2loc;
237924121865SAdrian Maldonado 
238024121865SAdrian Maldonado   PetscFunctionBegin;
23819566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(localsec, &size));
23829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &glob2loc));
238324121865SAdrian Maldonado 
238424121865SAdrian Maldonado   for (i = 0; i < size; i++) {
23859566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(globalsec, i, &dof));
238624121865SAdrian Maldonado     dof         = (dof >= 0) ? dof : -(dof + 1);
238724121865SAdrian Maldonado     glob2loc[i] = dof;
238824121865SAdrian Maldonado   }
238924121865SAdrian Maldonado 
2390e66d8692SHong Zhang   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)globalsec), 1, size, glob2loc, PETSC_OWN_POINTER, ltog));
239124121865SAdrian Maldonado #if 0
23929566063dSJacob Faibussowitsch   PetscCall(PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD));
239324121865SAdrian Maldonado #endif
23943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
239524121865SAdrian Maldonado }
239624121865SAdrian Maldonado 
239701ad2aeeSHong Zhang #include <petsc/private/matimpl.h>
23989e1d080bSHong Zhang 
2399d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateMatrix_Network_Nest(DM dm, Mat *J)
2400d71ae5a4SJacob Faibussowitsch {
24011ad426b7SHong Zhang   DM_Network            *network = (DM_Network *)dm->data;
240224121865SAdrian Maldonado   PetscInt               eDof, vDof;
240324121865SAdrian Maldonado   Mat                    j11, j12, j21, j22, bA[2][2];
24049e1d080bSHong Zhang   MPI_Comm               comm;
240524121865SAdrian Maldonado   ISLocalToGlobalMapping eISMap, vISMap;
240624121865SAdrian Maldonado 
24079e1d080bSHong Zhang   PetscFunctionBegin;
24089566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
240924121865SAdrian Maldonado 
24109566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection, &eDof));
24119566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection, &vDof));
241224121865SAdrian Maldonado 
24139566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &j11));
24149566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE));
24159566063dSJacob Faibussowitsch   PetscCall(MatSetType(j11, MATMPIAIJ));
241624121865SAdrian Maldonado 
24179566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &j12));
24189566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE));
24199566063dSJacob Faibussowitsch   PetscCall(MatSetType(j12, MATMPIAIJ));
242024121865SAdrian Maldonado 
24219566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &j21));
24229566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE));
24239566063dSJacob Faibussowitsch   PetscCall(MatSetType(j21, MATMPIAIJ));
242424121865SAdrian Maldonado 
24259566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &j22));
24269566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE));
24279566063dSJacob Faibussowitsch   PetscCall(MatSetType(j22, MATMPIAIJ));
242824121865SAdrian Maldonado 
24293f6a6bdaSHong Zhang   bA[0][0] = j11;
24303f6a6bdaSHong Zhang   bA[0][1] = j12;
24313f6a6bdaSHong Zhang   bA[1][0] = j21;
24323f6a6bdaSHong Zhang   bA[1][1] = j22;
243324121865SAdrian Maldonado 
24349566063dSJacob Faibussowitsch   PetscCall(CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection, network->edge.DofSection, &eISMap));
24359566063dSJacob Faibussowitsch   PetscCall(CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection, network->vertex.DofSection, &vISMap));
243624121865SAdrian Maldonado 
24379566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(j11, eISMap, eISMap));
24389566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(j12, eISMap, vISMap));
24399566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(j21, vISMap, eISMap));
24409566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(j22, vISMap, vISMap));
244124121865SAdrian Maldonado 
24429566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j11));
24439566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j12));
24449566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j21));
24459566063dSJacob Faibussowitsch   PetscCall(MatSetUp(j22));
244624121865SAdrian Maldonado 
24479566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(comm, 2, NULL, 2, NULL, &bA[0][0], J));
24489566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*J));
24499566063dSJacob Faibussowitsch   PetscCall(MatNestSetVecType(*J, VECNEST));
24509566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j11));
24519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j12));
24529566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j21));
24539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&j22));
245424121865SAdrian Maldonado 
24559566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY));
24569566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY));
24579566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*J, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
245824121865SAdrian Maldonado 
245924121865SAdrian Maldonado   /* Free structures */
24609566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&eISMap));
24619566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&vISMap));
24623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24639e1d080bSHong Zhang }
24649e1d080bSHong Zhang 
2465d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateMatrix_Network(DM dm, Mat *J)
2466d71ae5a4SJacob Faibussowitsch {
24679e1d080bSHong Zhang   DM_Network     *network = (DM_Network *)dm->data;
24689e1d080bSHong Zhang   PetscInt        eStart, eEnd, vStart, vEnd, rstart, nrows, *rows, localSize;
24699e1d080bSHong Zhang   PetscInt        cstart, ncols, j, e, v;
24709e1d080bSHong Zhang   PetscBool       ghost, ghost_vc, ghost2, isNest;
24719e1d080bSHong Zhang   Mat             Juser;
24729e1d080bSHong Zhang   PetscSection    sectionGlobal;
24739e1d080bSHong Zhang   PetscInt        nedges, *vptr = NULL, vc, *rows_v; /* suppress maybe-uninitialized warning */
24749e1d080bSHong Zhang   const PetscInt *edges, *cone;
24759e1d080bSHong Zhang   MPI_Comm        comm;
24769e1d080bSHong Zhang   MatType         mtype;
24779e1d080bSHong Zhang   Vec             vd_nz, vo_nz;
24789e1d080bSHong Zhang   PetscInt       *dnnz, *onnz;
24799e1d080bSHong Zhang   PetscScalar    *vdnz, *vonz;
24809e1d080bSHong Zhang 
24819e1d080bSHong Zhang   PetscFunctionBegin;
24829e1d080bSHong Zhang   mtype = dm->mattype;
24839566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATNEST, &isNest));
24849e1d080bSHong Zhang   if (isNest) {
24859566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix_Network_Nest(dm, J));
24869566063dSJacob Faibussowitsch     PetscCall(MatSetDM(*J, dm));
24873ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
24889e1d080bSHong Zhang   }
24899e1d080bSHong Zhang 
24909e1d080bSHong Zhang   if (!network->userEdgeJacobian && !network->userVertexJacobian) {
2491a4e85ca8SHong Zhang     /* user does not provide Jacobian blocks */
24929566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix_Plex(network->plex, J));
24939566063dSJacob Faibussowitsch     PetscCall(MatSetDM(*J, dm));
24943ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
24951ad426b7SHong Zhang   }
24961ad426b7SHong Zhang 
24979566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J));
24989566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(network->plex, &sectionGlobal));
24999566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize));
25009566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, localSize, localSize, PETSC_DETERMINE, PETSC_DETERMINE));
25012a945128SHong Zhang 
25029566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, MATAIJ));
25039566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(*J));
250489898e50SHong Zhang 
250589898e50SHong Zhang   /* (1) Set matrix preallocation */
250689898e50SHong Zhang   /*------------------------------*/
25079566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
25089566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &vd_nz));
25099566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(vd_nz, localSize, PETSC_DECIDE));
25109566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(vd_nz));
25119566063dSJacob Faibussowitsch   PetscCall(VecSet(vd_nz, 0.0));
25129566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(vd_nz, &vo_nz));
2513840c2264SHong Zhang 
251489898e50SHong Zhang   /* Set preallocation for edges */
251589898e50SHong Zhang   /*-----------------------------*/
25169566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(dm, &eStart, &eEnd));
2517840c2264SHong Zhang 
25189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(localSize, &rows));
2519840c2264SHong Zhang   for (e = eStart; e < eEnd; e++) {
2520840c2264SHong Zhang     /* Get row indices */
25219566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVecOffset(dm, e, ALL_COMPONENTS, &rstart));
25229566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(network->DofSection, e, &nrows));
2523840c2264SHong Zhang     if (nrows) {
2524840c2264SHong Zhang       for (j = 0; j < nrows; j++) rows[j] = j + rstart;
2525840c2264SHong Zhang 
2526a5b23f4aSJose E. Roman       /* Set preallocation for connected vertices */
25279566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetConnectedVertices(dm, e, &cone));
2528840c2264SHong Zhang       for (v = 0; v < 2; v++) {
25299566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(network->DofSection, cone[v], &ncols));
2530840c2264SHong Zhang 
25318675203cSHong Zhang         if (network->Je) {
2532840c2264SHong Zhang           Juser = network->Je[3 * e + 1 + v]; /* Jacobian(e,v) */
25338675203cSHong Zhang         } else Juser = NULL;
25349566063dSJacob Faibussowitsch         PetscCall(DMNetworkIsGhostVertex(dm, cone[v], &ghost));
25359566063dSJacob Faibussowitsch         PetscCall(MatSetPreallocationblock_private(Juser, nrows, rows, ncols, ghost, vd_nz, vo_nz));
2536840c2264SHong Zhang       }
2537840c2264SHong Zhang 
253889898e50SHong Zhang       /* Set preallocation for edge self */
2539840c2264SHong Zhang       cstart = rstart;
25408675203cSHong Zhang       if (network->Je) {
2541840c2264SHong Zhang         Juser = network->Je[3 * e]; /* Jacobian(e,e) */
25428675203cSHong Zhang       } else Juser = NULL;
25439566063dSJacob Faibussowitsch       PetscCall(MatSetPreallocationblock_private(Juser, nrows, rows, nrows, PETSC_FALSE, vd_nz, vo_nz));
2544840c2264SHong Zhang     }
2545840c2264SHong Zhang   }
2546840c2264SHong Zhang 
254789898e50SHong Zhang   /* Set preallocation for vertices */
254889898e50SHong Zhang   /*--------------------------------*/
25499566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(dm, &vStart, &vEnd));
25508675203cSHong Zhang   if (vEnd - vStart) vptr = network->Jvptr;
2551840c2264SHong Zhang 
2552840c2264SHong Zhang   for (v = vStart; v < vEnd; v++) {
2553840c2264SHong Zhang     /* Get row indices */
25549566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVecOffset(dm, v, ALL_COMPONENTS, &rstart));
25559566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(network->DofSection, v, &nrows));
2556840c2264SHong Zhang     if (!nrows) continue;
2557840c2264SHong Zhang 
25589566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(dm, v, &ghost));
2559bdcb62a2SHong Zhang     if (ghost) {
25609566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nrows, &rows_v));
2561bdcb62a2SHong Zhang     } else {
2562bdcb62a2SHong Zhang       rows_v = rows;
2563bdcb62a2SHong Zhang     }
2564bdcb62a2SHong Zhang 
2565bdcb62a2SHong Zhang     for (j = 0; j < nrows; j++) rows_v[j] = j + rstart;
2566840c2264SHong Zhang 
2567840c2264SHong Zhang     /* Get supporting edges and connected vertices */
25689566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetSupportingEdges(dm, v, &nedges, &edges));
2569840c2264SHong Zhang 
2570840c2264SHong Zhang     for (e = 0; e < nedges; e++) {
2571840c2264SHong Zhang       /* Supporting edges */
25729566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetGlobalVecOffset(dm, edges[e], ALL_COMPONENTS, &cstart));
25739566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(network->DofSection, edges[e], &ncols));
2574840c2264SHong Zhang 
25758675203cSHong Zhang       if (network->Jv) {
2576840c2264SHong Zhang         Juser = network->Jv[vptr[v - vStart] + 2 * e + 1]; /* Jacobian(v,e) */
25778675203cSHong Zhang       } else Juser = NULL;
25789566063dSJacob Faibussowitsch       PetscCall(MatSetPreallocationblock_private(Juser, nrows, rows_v, ncols, ghost, vd_nz, vo_nz));
2579840c2264SHong Zhang 
2580840c2264SHong Zhang       /* Connected vertices */
25819566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetConnectedVertices(dm, edges[e], &cone));
2582840c2264SHong Zhang       vc = (v == cone[0]) ? cone[1] : cone[0];
25839566063dSJacob Faibussowitsch       PetscCall(DMNetworkIsGhostVertex(dm, vc, &ghost_vc));
2584840c2264SHong Zhang 
25859566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(network->DofSection, vc, &ncols));
2586840c2264SHong Zhang 
25878675203cSHong Zhang       if (network->Jv) {
2588840c2264SHong Zhang         Juser = network->Jv[vptr[v - vStart] + 2 * e + 2]; /* Jacobian(v,vc) */
25898675203cSHong Zhang       } else Juser = NULL;
2590e102a522SHong Zhang       if (ghost_vc || ghost) {
2591e102a522SHong Zhang         ghost2 = PETSC_TRUE;
2592e102a522SHong Zhang       } else {
2593e102a522SHong Zhang         ghost2 = PETSC_FALSE;
2594e102a522SHong Zhang       }
25959566063dSJacob Faibussowitsch       PetscCall(MatSetPreallocationblock_private(Juser, nrows, rows_v, ncols, ghost2, vd_nz, vo_nz));
2596840c2264SHong Zhang     }
2597840c2264SHong Zhang 
259889898e50SHong Zhang     /* Set preallocation for vertex self */
25999566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(dm, v, &ghost));
2600840c2264SHong Zhang     if (!ghost) {
26019566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetGlobalVecOffset(dm, v, ALL_COMPONENTS, &cstart));
26028675203cSHong Zhang       if (network->Jv) {
2603840c2264SHong Zhang         Juser = network->Jv[vptr[v - vStart]]; /* Jacobian(v,v) */
26048675203cSHong Zhang       } else Juser = NULL;
26059566063dSJacob Faibussowitsch       PetscCall(MatSetPreallocationblock_private(Juser, nrows, rows_v, nrows, PETSC_FALSE, vd_nz, vo_nz));
2606840c2264SHong Zhang     }
26071baa6e33SBarry Smith     if (ghost) PetscCall(PetscFree(rows_v));
2608840c2264SHong Zhang   }
2609840c2264SHong Zhang 
26109566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(vd_nz));
26119566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(vo_nz));
26125cf7da58SHong Zhang 
26139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(localSize, &dnnz, localSize, &onnz));
26145cf7da58SHong Zhang 
26159566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(vd_nz));
26169566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(vo_nz));
2617840c2264SHong Zhang 
26189566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vd_nz, &vdnz));
26199566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vo_nz, &vonz));
2620840c2264SHong Zhang   for (j = 0; j < localSize; j++) {
2621e102a522SHong Zhang     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
2622e102a522SHong Zhang     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
2623840c2264SHong Zhang   }
26249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vd_nz, &vdnz));
26259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vo_nz, &vonz));
26269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vd_nz));
26279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vo_nz));
2628840c2264SHong Zhang 
26299566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*J, 0, dnnz));
26309566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*J, 0, dnnz, 0, onnz));
26319566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*J, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
26325cf7da58SHong Zhang 
26339566063dSJacob Faibussowitsch   PetscCall(PetscFree2(dnnz, onnz));
26345cf7da58SHong Zhang 
263589898e50SHong Zhang   /* (2) Set matrix entries for edges */
263689898e50SHong Zhang   /*----------------------------------*/
26371ad426b7SHong Zhang   for (e = eStart; e < eEnd; e++) {
2638bfbc38dcSHong Zhang     /* Get row indices */
26399566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVecOffset(dm, e, ALL_COMPONENTS, &rstart));
26409566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(network->DofSection, e, &nrows));
26414b976069SHong Zhang     if (nrows) {
264217df6e9eSHong Zhang       for (j = 0; j < nrows; j++) rows[j] = j + rstart;
26431ad426b7SHong Zhang 
2644a5b23f4aSJose E. Roman       /* Set matrix entries for connected vertices */
26459566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetConnectedVertices(dm, e, &cone));
2646bfbc38dcSHong Zhang       for (v = 0; v < 2; v++) {
26479566063dSJacob Faibussowitsch         PetscCall(DMNetworkGetGlobalVecOffset(dm, cone[v], ALL_COMPONENTS, &cstart));
26489566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(network->DofSection, cone[v], &ncols));
26493e97b6e8SHong Zhang 
26508675203cSHong Zhang         if (network->Je) {
2651a4e85ca8SHong Zhang           Juser = network->Je[3 * e + 1 + v]; /* Jacobian(e,v) */
26528675203cSHong Zhang         } else Juser = NULL;
26539566063dSJacob Faibussowitsch         PetscCall(MatSetblock_private(Juser, nrows, rows, ncols, cstart, J));
2654bfbc38dcSHong Zhang       }
265517df6e9eSHong Zhang 
2656bfbc38dcSHong Zhang       /* Set matrix entries for edge self */
26573e97b6e8SHong Zhang       cstart = rstart;
26588675203cSHong Zhang       if (network->Je) {
2659a4e85ca8SHong Zhang         Juser = network->Je[3 * e]; /* Jacobian(e,e) */
26608675203cSHong Zhang       } else Juser = NULL;
26619566063dSJacob Faibussowitsch       PetscCall(MatSetblock_private(Juser, nrows, rows, nrows, cstart, J));
26621ad426b7SHong Zhang     }
26634b976069SHong Zhang   }
26641ad426b7SHong Zhang 
2665bfbc38dcSHong Zhang   /* Set matrix entries for vertices */
266683b2e829SHong Zhang   /*---------------------------------*/
26671ad426b7SHong Zhang   for (v = vStart; v < vEnd; v++) {
2668bfbc38dcSHong Zhang     /* Get row indices */
26699566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetGlobalVecOffset(dm, v, ALL_COMPONENTS, &rstart));
26709566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(network->DofSection, v, &nrows));
26714b976069SHong Zhang     if (!nrows) continue;
2672596e729fSHong Zhang 
26739566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(dm, v, &ghost));
2674bdcb62a2SHong Zhang     if (ghost) {
26759566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nrows, &rows_v));
2676bdcb62a2SHong Zhang     } else {
2677bdcb62a2SHong Zhang       rows_v = rows;
2678bdcb62a2SHong Zhang     }
2679bdcb62a2SHong Zhang     for (j = 0; j < nrows; j++) rows_v[j] = j + rstart;
2680596e729fSHong Zhang 
2681bfbc38dcSHong Zhang     /* Get supporting edges and connected vertices */
26829566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetSupportingEdges(dm, v, &nedges, &edges));
2683596e729fSHong Zhang 
2684596e729fSHong Zhang     for (e = 0; e < nedges; e++) {
2685bfbc38dcSHong Zhang       /* Supporting edges */
26869566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetGlobalVecOffset(dm, edges[e], ALL_COMPONENTS, &cstart));
26879566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(network->DofSection, edges[e], &ncols));
2688596e729fSHong Zhang 
26898675203cSHong Zhang       if (network->Jv) {
2690a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v - vStart] + 2 * e + 1]; /* Jacobian(v,e) */
26918675203cSHong Zhang       } else Juser = NULL;
26929566063dSJacob Faibussowitsch       PetscCall(MatSetblock_private(Juser, nrows, rows_v, ncols, cstart, J));
2693596e729fSHong Zhang 
2694bfbc38dcSHong Zhang       /* Connected vertices */
26959566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetConnectedVertices(dm, edges[e], &cone));
26962a945128SHong Zhang       vc = (v == cone[0]) ? cone[1] : cone[0];
26972a945128SHong Zhang 
26989566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetGlobalVecOffset(dm, vc, ALL_COMPONENTS, &cstart));
26999566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(network->DofSection, vc, &ncols));
2700a4e85ca8SHong Zhang 
27018675203cSHong Zhang       if (network->Jv) {
2702a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v - vStart] + 2 * e + 2]; /* Jacobian(v,vc) */
27038675203cSHong Zhang       } else Juser = NULL;
27049566063dSJacob Faibussowitsch       PetscCall(MatSetblock_private(Juser, nrows, rows_v, ncols, cstart, J));
2705596e729fSHong Zhang     }
2706596e729fSHong Zhang 
2707bfbc38dcSHong Zhang     /* Set matrix entries for vertex self */
27081ad426b7SHong Zhang     if (!ghost) {
27099566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetGlobalVecOffset(dm, v, ALL_COMPONENTS, &cstart));
27108675203cSHong Zhang       if (network->Jv) {
2711a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v - vStart]]; /* Jacobian(v,v) */
27128675203cSHong Zhang       } else Juser = NULL;
27139566063dSJacob Faibussowitsch       PetscCall(MatSetblock_private(Juser, nrows, rows_v, nrows, cstart, J));
2714bdcb62a2SHong Zhang     }
27151baa6e33SBarry Smith     if (ghost) PetscCall(PetscFree(rows_v));
27161ad426b7SHong Zhang   }
27179566063dSJacob Faibussowitsch   PetscCall(PetscFree(rows));
2718bdcb62a2SHong Zhang 
27199566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY));
27209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY));
2721dd6f46cdSHong Zhang 
27229566063dSJacob Faibussowitsch   PetscCall(MatSetDM(*J, dm));
27233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27245f2c45f1SShri Abhyankar }
27255f2c45f1SShri Abhyankar 
2726d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMNetworkDestroyComponentData(DM dm)
2727d71ae5a4SJacob Faibussowitsch {
2728daad07d3SAidan Hamilton   DM_Network *network = (DM_Network *)dm->data;
2729daad07d3SAidan Hamilton   PetscInt    j, np;
2730daad07d3SAidan Hamilton 
2731daad07d3SAidan Hamilton   PetscFunctionBegin;
2732daad07d3SAidan Hamilton   if (network->header) {
2733daad07d3SAidan Hamilton     np = network->cloneshared->pEnd - network->cloneshared->pStart;
2734daad07d3SAidan Hamilton     for (j = 0; j < np; j++) {
2735daad07d3SAidan Hamilton       PetscCall(PetscFree5(network->header[j].size, network->header[j].key, network->header[j].offset, network->header[j].nvar, network->header[j].offsetvarrel));
2736daad07d3SAidan Hamilton       PetscCall(PetscFree(network->cvalue[j].data));
2737daad07d3SAidan Hamilton     }
2738daad07d3SAidan Hamilton     PetscCall(PetscFree2(network->header, network->cvalue));
2739daad07d3SAidan Hamilton   }
27403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2741daad07d3SAidan Hamilton }
2742daad07d3SAidan Hamilton 
2743d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDestroy_Network(DM dm)
2744d71ae5a4SJacob Faibussowitsch {
27455f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
2746aeb78aa7SBarry Smith   PetscInt    j;
27475f2c45f1SShri Abhyankar 
27485f2c45f1SShri Abhyankar   PetscFunctionBegin;
2749daad07d3SAidan Hamilton   /*
2750daad07d3SAidan Hamilton     Developers Note: Due to the mixed topological definition of DMNetwork and data defined ON the
2751daad07d3SAidan Hamilton     network like DofSection, DataSection, *componentdataarray, and friends, when cloning, we share
2752daad07d3SAidan Hamilton     only the true topological data, and make our own data ON the network. Thus refct only refers
2753daad07d3SAidan Hamilton     to the number of references to topological data, and data ON the network is always destroyed.
2754daad07d3SAidan Hamilton     It is understood this is atypical for a DM, but is very intentional.
2755daad07d3SAidan Hamilton   */
2756daad07d3SAidan Hamilton 
2757daad07d3SAidan Hamilton   /* Always destroy data ON the network */
27589566063dSJacob Faibussowitsch   PetscCall(PetscFree(network->Je));
275983b2e829SHong Zhang   if (network->Jv) {
27609566063dSJacob Faibussowitsch     PetscCall(PetscFree(network->Jvptr));
27619566063dSJacob Faibussowitsch     PetscCall(PetscFree(network->Jv));
27621ad426b7SHong Zhang   }
2763daad07d3SAidan Hamilton   PetscCall(PetscSectionDestroy(&network->DataSection));
2764daad07d3SAidan Hamilton   PetscCall(PetscSectionDestroy(&network->DofSection));
2765daad07d3SAidan Hamilton   PetscCall(PetscFree(network->component));
2766daad07d3SAidan Hamilton   PetscCall(PetscFree(network->componentdataarray));
2767daad07d3SAidan Hamilton   PetscCall(DMNetworkDestroyComponentData(dm));
276813c2a604SAdrian Maldonado 
2769daad07d3SAidan Hamilton   PetscCall(DMDestroy(&network->plex)); /* this is cloned in DMClone_Network, so safe to destroy */
2770daad07d3SAidan Hamilton 
2771daad07d3SAidan Hamilton   /*
2772daad07d3SAidan Hamilton   Developers Note: The DMNetworkVertexInfo and DMNetworkEdgeInfo data structures are completely
2773daad07d3SAidan Hamilton   destroyed as they are again a mix of topological data:
2774daad07d3SAidan Hamilton     ISLocalToGlobalMapping            mapping;
2775daad07d3SAidan Hamilton     PetscSF                           sf;
2776daad07d3SAidan Hamilton   and data ON the network
2777daad07d3SAidan Hamilton     PetscSection                      DofSection;
2778daad07d3SAidan Hamilton     PetscSection                      GlobalDofSection;
2779daad07d3SAidan Hamilton   And the only way to access them currently is through DMNetworkAssembleGraphStructures which assembles
2780daad07d3SAidan Hamilton   everything. So we must destroy everything and require DMNetworkAssembleGraphStructures is called again
2781daad07d3SAidan Hamilton   for any clone.
2782daad07d3SAidan Hamilton   */
27839566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&network->vertex.mapping));
27849566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&network->vertex.DofSection));
27859566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&network->vertex.GlobalDofSection));
27869566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&network->vertex.sf));
278713c2a604SAdrian Maldonado   /* edge */
27889566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&network->edge.mapping));
27899566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&network->edge.DofSection));
27909566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&network->edge.GlobalDofSection));
27919566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&network->edge.sf));
27925f25b224SDuncan Campbell   /* viewer options */
27935f25b224SDuncan Campbell   PetscCall(ISDestroy(&network->vieweroptions.viewranks));
2794daad07d3SAidan Hamilton   /* Destroy the potentially cloneshared data */
2795daad07d3SAidan Hamilton   if (--network->cloneshared->refct <= 0) {
2796daad07d3SAidan 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
2797daad07d3SAidan Hamilton      naively think it can be reused. */
2798daad07d3SAidan Hamilton     PetscCall(PetscFree(network->cloneshared->vltog));
2799daad07d3SAidan Hamilton     for (j = 0; j < network->cloneshared->Nsvtx; j++) PetscCall(PetscFree(network->cloneshared->svtx[j].sv));
2800daad07d3SAidan Hamilton     PetscCall(PetscFree(network->cloneshared->svtx));
2801daad07d3SAidan Hamilton     PetscCall(PetscFree2(network->cloneshared->subnetedge, network->cloneshared->subnetvtx));
2802eec179cfSJacob Faibussowitsch     PetscCall(PetscHMapIDestroy(&network->cloneshared->svtable));
2803daad07d3SAidan Hamilton     PetscCall(PetscFree(network->cloneshared->subnet));
2804daad07d3SAidan Hamilton     PetscCall(PetscFree(network->cloneshared));
2805daad07d3SAidan Hamilton   }
2806daad07d3SAidan Hamilton   PetscCall(PetscFree(network)); /* Always freed as this structure is copied in a clone, not cloneshared */
28073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28085f2c45f1SShri Abhyankar }
28095f2c45f1SShri Abhyankar 
2810d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
2811d71ae5a4SJacob Faibussowitsch {
28125f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
28135f2c45f1SShri Abhyankar 
28145f2c45f1SShri Abhyankar   PetscFunctionBegin;
28159566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(network->plex, g, mode, l));
28163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28175f2c45f1SShri Abhyankar }
28185f2c45f1SShri Abhyankar 
2819d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
2820d71ae5a4SJacob Faibussowitsch {
28215f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
28225f2c45f1SShri Abhyankar 
28235f2c45f1SShri Abhyankar   PetscFunctionBegin;
28249566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(network->plex, g, mode, l));
28253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28265f2c45f1SShri Abhyankar }
28275f2c45f1SShri Abhyankar 
2828d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
2829d71ae5a4SJacob Faibussowitsch {
28305f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
28315f2c45f1SShri Abhyankar 
28325f2c45f1SShri Abhyankar   PetscFunctionBegin;
28339566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(network->plex, l, mode, g));
28343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28355f2c45f1SShri Abhyankar }
28365f2c45f1SShri Abhyankar 
2837d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
2838d71ae5a4SJacob Faibussowitsch {
28395f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network *)dm->data;
28405f2c45f1SShri Abhyankar 
28415f2c45f1SShri Abhyankar   PetscFunctionBegin;
28429566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(network->plex, l, mode, g));
28433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28445f2c45f1SShri Abhyankar }
284522bbedd7SHong Zhang 
284622bbedd7SHong Zhang /*@
284764238763SRodolfo Oliveira   DMNetworkGetVertexLocalToGlobalOrdering - Get vertex global index
284822bbedd7SHong Zhang 
284920f4b53cSBarry Smith   Not Collective
285022bbedd7SHong Zhang 
28517a7aea1fSJed Brown   Input Parameters:
285220f4b53cSBarry Smith + dm   - the `DMNETWORK` object
285322bbedd7SHong Zhang - vloc - local vertex ordering, start from 0
285422bbedd7SHong Zhang 
28552fe279fdSBarry Smith   Output Parameter:
2856f0fc11ceSJed Brown . vg - global vertex ordering, start from 0
285722bbedd7SHong Zhang 
285897bb938eSShri Abhyankar   Level: advanced
285922bbedd7SHong Zhang 
286020f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkSetVertexLocalToGlobalOrdering()`
286122bbedd7SHong Zhang @*/
2862d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkGetVertexLocalToGlobalOrdering(DM dm, PetscInt vloc, PetscInt *vg)
2863d71ae5a4SJacob Faibussowitsch {
286422bbedd7SHong Zhang   DM_Network *network = (DM_Network *)dm->data;
2865daad07d3SAidan Hamilton   PetscInt   *vltog   = network->cloneshared->vltog;
286622bbedd7SHong Zhang 
286722bbedd7SHong Zhang   PetscFunctionBegin;
28685c6496baSHong Zhang   PetscCheck(vltog, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Must call DMNetworkSetVertexLocalToGlobalOrdering() first");
286922bbedd7SHong Zhang   *vg = vltog[vloc];
28703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
287122bbedd7SHong Zhang }
287222bbedd7SHong Zhang 
287322bbedd7SHong Zhang /*@
287464238763SRodolfo Oliveira   DMNetworkSetVertexLocalToGlobalOrdering - Create and setup vertex local to global map
287522bbedd7SHong Zhang 
287622bbedd7SHong Zhang   Collective
287722bbedd7SHong Zhang 
287822bbedd7SHong Zhang   Input Parameters:
287920f4b53cSBarry Smith . dm - the `DMNETWORK` object
288022bbedd7SHong Zhang 
288197bb938eSShri Abhyankar   Level: advanced
288222bbedd7SHong Zhang 
288320f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkGetGlobalVertexIndex()`
288422bbedd7SHong Zhang @*/
2885d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkSetVertexLocalToGlobalOrdering(DM dm)
2886d71ae5a4SJacob Faibussowitsch {
288722bbedd7SHong Zhang   DM_Network        *network = (DM_Network *)dm->data;
288822bbedd7SHong Zhang   MPI_Comm           comm;
28892bf73ac6SHong Zhang   PetscMPIInt        rank, size, *displs = NULL, *recvcounts = NULL, remoterank;
289022bbedd7SHong Zhang   PetscBool          ghost;
289163029d29SHong Zhang   PetscInt          *vltog, nroots, nleaves, i, *vrange, k, N, lidx;
289222bbedd7SHong Zhang   const PetscSFNode *iremote;
289322bbedd7SHong Zhang   PetscSF            vsf;
289463029d29SHong Zhang   Vec                Vleaves, Vleaves_seq;
289563029d29SHong Zhang   VecScatter         ctx;
289663029d29SHong Zhang   PetscScalar       *varr, val;
289763029d29SHong Zhang   const PetscScalar *varr_read;
289822bbedd7SHong Zhang 
289922bbedd7SHong Zhang   PetscFunctionBegin;
29009566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
29019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
29029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
290322bbedd7SHong Zhang 
290422bbedd7SHong Zhang   if (size == 1) {
2905daad07d3SAidan Hamilton     nroots = network->cloneshared->vEnd - network->cloneshared->vStart;
29069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nroots, &vltog));
290722bbedd7SHong Zhang     for (i = 0; i < nroots; i++) vltog[i] = i;
2908daad07d3SAidan Hamilton     network->cloneshared->vltog = vltog;
29093ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
291022bbedd7SHong Zhang   }
291122bbedd7SHong Zhang 
2912daad07d3SAidan Hamilton   PetscCheck(network->cloneshared->distributecalled, comm, PETSC_ERR_ARG_WRONGSTATE, "Must call DMNetworkDistribute() first");
2913daad07d3SAidan Hamilton   if (network->cloneshared->vltog) PetscCall(PetscFree(network->cloneshared->vltog));
291422bbedd7SHong Zhang 
2915e66d8692SHong Zhang   PetscCall(DMNetworkSetSubMap_private(dm, network->cloneshared->vStart, network->cloneshared->vEnd, &network->vertex.mapping));
29169566063dSJacob Faibussowitsch   PetscCall(PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf));
291722bbedd7SHong Zhang   vsf = network->vertex.sf;
291822bbedd7SHong Zhang 
29199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(size + 1, &vrange, size, &displs, size, &recvcounts));
29209566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(vsf, &nroots, &nleaves, NULL, &iremote));
292122bbedd7SHong Zhang 
29229371c9d4SSatish Balay   for (i = 0; i < size; i++) {
29239371c9d4SSatish Balay     displs[i]     = i;
29249371c9d4SSatish Balay     recvcounts[i] = 1;
29259371c9d4SSatish Balay   }
292622bbedd7SHong Zhang 
292722bbedd7SHong Zhang   i         = nroots - nleaves; /* local number of vertices, excluding ghosts */
292822bbedd7SHong Zhang   vrange[0] = 0;
29299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgatherv(&i, 1, MPIU_INT, vrange + 1, recvcounts, displs, MPIU_INT, comm));
2930ad540459SPierre Jolivet   for (i = 2; i < size + 1; i++) vrange[i] += vrange[i - 1];
293122bbedd7SHong Zhang 
29329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nroots, &vltog));
2933daad07d3SAidan Hamilton   network->cloneshared->vltog = vltog;
293422bbedd7SHong Zhang 
293563029d29SHong Zhang   /* Set vltog for non-ghost vertices */
293663029d29SHong Zhang   k = 0;
293722bbedd7SHong Zhang   for (i = 0; i < nroots; i++) {
2938daad07d3SAidan Hamilton     PetscCall(DMNetworkIsGhostVertex(dm, i + network->cloneshared->vStart, &ghost));
293963029d29SHong Zhang     if (ghost) continue;
294063029d29SHong Zhang     vltog[i] = vrange[rank] + k++;
294122bbedd7SHong Zhang   }
29429566063dSJacob Faibussowitsch   PetscCall(PetscFree3(vrange, displs, recvcounts));
294363029d29SHong Zhang 
294463029d29SHong Zhang   /* Set vltog for ghost vertices */
294563029d29SHong Zhang   /* (a) create parallel Vleaves and sequential Vleaves_seq to convert local iremote[*].index to global index */
29469566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &Vleaves));
29479566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(Vleaves, 2 * nleaves, PETSC_DETERMINE));
29489566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(Vleaves));
29499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Vleaves, &varr));
295063029d29SHong Zhang   for (i = 0; i < nleaves; i++) {
295163029d29SHong Zhang     varr[2 * i]     = (PetscScalar)(iremote[i].rank);  /* rank of remote process */
295263029d29SHong Zhang     varr[2 * i + 1] = (PetscScalar)(iremote[i].index); /* local index in remote process */
295363029d29SHong Zhang   }
29549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Vleaves, &varr));
295563029d29SHong Zhang 
295663029d29SHong Zhang   /* (b) scatter local info to remote processes via VecScatter() */
29579566063dSJacob Faibussowitsch   PetscCall(VecScatterCreateToAll(Vleaves, &ctx, &Vleaves_seq));
29589566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(ctx, Vleaves, Vleaves_seq, INSERT_VALUES, SCATTER_FORWARD));
29599566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(ctx, Vleaves, Vleaves_seq, INSERT_VALUES, SCATTER_FORWARD));
296063029d29SHong Zhang 
296163029d29SHong Zhang   /* (c) convert local indices to global indices in parallel vector Vleaves */
29629566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Vleaves_seq, &N));
29639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vleaves_seq, &varr_read));
296463029d29SHong Zhang   for (i = 0; i < N; i += 2) {
29652e4cff2eSHong Zhang     remoterank = (PetscMPIInt)PetscRealPart(varr_read[i]);
296663029d29SHong Zhang     if (remoterank == rank) {
296763029d29SHong Zhang       k    = i + 1; /* row number */
29682e4cff2eSHong Zhang       lidx = (PetscInt)PetscRealPart(varr_read[i + 1]);
296963029d29SHong Zhang       val  = (PetscScalar)vltog[lidx]; /* global index for non-ghost vertex computed above */
29709566063dSJacob Faibussowitsch       PetscCall(VecSetValues(Vleaves, 1, &k, &val, INSERT_VALUES));
297163029d29SHong Zhang     }
297263029d29SHong Zhang   }
29739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vleaves_seq, &varr_read));
29749566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(Vleaves));
29759566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(Vleaves));
297663029d29SHong Zhang 
297763029d29SHong Zhang   /* (d) Set vltog for ghost vertices by copying local values of Vleaves */
29789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vleaves, &varr_read));
297963029d29SHong Zhang   k = 0;
298063029d29SHong Zhang   for (i = 0; i < nroots; i++) {
2981daad07d3SAidan Hamilton     PetscCall(DMNetworkIsGhostVertex(dm, i + network->cloneshared->vStart, &ghost));
298263029d29SHong Zhang     if (!ghost) continue;
29839371c9d4SSatish Balay     vltog[i] = (PetscInt)PetscRealPart(varr_read[2 * k + 1]);
29849371c9d4SSatish Balay     k++;
298563029d29SHong Zhang   }
29869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vleaves, &varr_read));
298763029d29SHong Zhang 
29889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Vleaves));
29899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Vleaves_seq));
29909566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&ctx));
29913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
299222bbedd7SHong Zhang }
299342dc13f1SHong Zhang 
2994d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode DMISAddSize_private(DM_Network *network, PetscInt p, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *nidx)
2995d71ae5a4SJacob Faibussowitsch {
299642dc13f1SHong Zhang   PetscInt                 i, j, ncomps, nvar, key, offset = 0;
299742dc13f1SHong Zhang   DMNetworkComponentHeader header;
299842dc13f1SHong Zhang 
299942dc13f1SHong Zhang   PetscFunctionBegin;
30009566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection, p, &offset));
300142dc13f1SHong Zhang   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray + offset))->ndata;
300242dc13f1SHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray + offset);
300342dc13f1SHong Zhang 
300442dc13f1SHong Zhang   for (i = 0; i < ncomps; i++) {
300542dc13f1SHong Zhang     key  = header->key[i];
300642dc13f1SHong Zhang     nvar = header->nvar[i];
300742dc13f1SHong Zhang     for (j = 0; j < numkeys; j++) {
300842dc13f1SHong Zhang       if (key == keys[j]) {
300942dc13f1SHong Zhang         if (!blocksize || blocksize[j] == -1) {
301042dc13f1SHong Zhang           *nidx += nvar;
301142dc13f1SHong Zhang         } else {
301242dc13f1SHong Zhang           *nidx += nselectedvar[j] * nvar / blocksize[j];
301342dc13f1SHong Zhang         }
301442dc13f1SHong Zhang       }
301542dc13f1SHong Zhang     }
301642dc13f1SHong Zhang   }
30173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
301842dc13f1SHong Zhang }
301942dc13f1SHong Zhang 
3020d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode DMISComputeIdx_private(DM dm, PetscInt p, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *selectedvar[], PetscInt *ii, PetscInt *idx)
3021d71ae5a4SJacob Faibussowitsch {
302242dc13f1SHong Zhang   PetscInt                 i, j, ncomps, nvar, key, offsetg, k, k1, offset = 0;
302342dc13f1SHong Zhang   DM_Network              *network = (DM_Network *)dm->data;
302442dc13f1SHong Zhang   DMNetworkComponentHeader header;
302542dc13f1SHong Zhang 
302642dc13f1SHong Zhang   PetscFunctionBegin;
30279566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection, p, &offset));
302842dc13f1SHong Zhang   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray + offset))->ndata;
302942dc13f1SHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray + offset);
303042dc13f1SHong Zhang 
303142dc13f1SHong Zhang   for (i = 0; i < ncomps; i++) {
303242dc13f1SHong Zhang     key  = header->key[i];
303342dc13f1SHong Zhang     nvar = header->nvar[i];
303442dc13f1SHong Zhang     for (j = 0; j < numkeys; j++) {
303542dc13f1SHong Zhang       if (key != keys[j]) continue;
303642dc13f1SHong Zhang 
30379566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetGlobalVecOffset(dm, p, i, &offsetg));
303842dc13f1SHong Zhang       if (!blocksize || blocksize[j] == -1) {
303942dc13f1SHong Zhang         for (k = 0; k < nvar; k++) idx[(*ii)++] = offsetg + k;
304042dc13f1SHong Zhang       } else {
304142dc13f1SHong Zhang         for (k = 0; k < nvar; k += blocksize[j]) {
304242dc13f1SHong Zhang           for (k1 = 0; k1 < nselectedvar[j]; k1++) idx[(*ii)++] = offsetg + k + selectedvar[j][k1];
304342dc13f1SHong Zhang         }
304442dc13f1SHong Zhang       }
304542dc13f1SHong Zhang     }
304642dc13f1SHong Zhang   }
30473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
304842dc13f1SHong Zhang }
304942dc13f1SHong Zhang 
305042dc13f1SHong Zhang /*@
305142dc13f1SHong Zhang   DMNetworkCreateIS - Create an index set object from the global vector of the network
305242dc13f1SHong Zhang 
305342dc13f1SHong Zhang   Collective
305442dc13f1SHong Zhang 
305542dc13f1SHong Zhang   Input Parameters:
305620f4b53cSBarry Smith + dm           - `DMNETWORK` object
305742dc13f1SHong Zhang . numkeys      - number of keys
305842dc13f1SHong Zhang . keys         - array of keys that define the components of the variables you wish to extract
305942dc13f1SHong Zhang . blocksize    - block size of the variables associated to the component
306042dc13f1SHong Zhang . nselectedvar - number of variables in each block to select
306142dc13f1SHong Zhang - selectedvar  - the offset into the block of each variable in each block to select
306242dc13f1SHong Zhang 
30632fe279fdSBarry Smith   Output Parameter:
306442dc13f1SHong Zhang . is - the index set
306542dc13f1SHong Zhang 
306620f4b53cSBarry Smith   Level: advanced
306742dc13f1SHong Zhang 
306842dc13f1SHong Zhang   Notes:
306920f4b53cSBarry Smith   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.
307042dc13f1SHong Zhang 
307120f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkCreate()`, `ISCreateGeneral()`, `DMNetworkCreateLocalIS()`
307242dc13f1SHong Zhang @*/
3073d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkCreateIS(DM dm, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *selectedvar[], IS *is)
3074d71ae5a4SJacob Faibussowitsch {
307542dc13f1SHong Zhang   MPI_Comm    comm;
307642dc13f1SHong Zhang   DM_Network *network = (DM_Network *)dm->data;
307742dc13f1SHong Zhang   PetscInt    i, p, estart, eend, vstart, vend, nidx, *idx;
307842dc13f1SHong Zhang   PetscBool   ghost;
307942dc13f1SHong Zhang 
308042dc13f1SHong Zhang   PetscFunctionBegin;
30819566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
308242dc13f1SHong Zhang 
308342dc13f1SHong Zhang   /* Check input parameters */
308442dc13f1SHong Zhang   for (i = 0; i < numkeys; i++) {
308542dc13f1SHong Zhang     if (!blocksize || blocksize[i] == -1) continue;
308663a3b9bcSJacob 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]);
308742dc13f1SHong Zhang   }
308842dc13f1SHong Zhang 
30899566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetEdgeRange(dm, &estart, &eend));
30909566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(dm, &vstart, &vend));
309142dc13f1SHong Zhang 
309242dc13f1SHong Zhang   /* Get local number of idx */
309342dc13f1SHong Zhang   nidx = 0;
309448a46eb9SPierre Jolivet   for (p = estart; p < eend; p++) PetscCall(DMISAddSize_private(network, p, numkeys, keys, blocksize, nselectedvar, &nidx));
309542dc13f1SHong Zhang   for (p = vstart; p < vend; p++) {
30969566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(dm, p, &ghost));
309742dc13f1SHong Zhang     if (ghost) continue;
30989566063dSJacob Faibussowitsch     PetscCall(DMISAddSize_private(network, p, numkeys, keys, blocksize, nselectedvar, &nidx));
309942dc13f1SHong Zhang   }
310042dc13f1SHong Zhang 
310142dc13f1SHong Zhang   /* Compute idx */
31029566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nidx, &idx));
310342dc13f1SHong Zhang   i = 0;
310448a46eb9SPierre Jolivet   for (p = estart; p < eend; p++) PetscCall(DMISComputeIdx_private(dm, p, numkeys, keys, blocksize, nselectedvar, selectedvar, &i, idx));
310542dc13f1SHong Zhang   for (p = vstart; p < vend; p++) {
31069566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(dm, p, &ghost));
310742dc13f1SHong Zhang     if (ghost) continue;
31089566063dSJacob Faibussowitsch     PetscCall(DMISComputeIdx_private(dm, p, numkeys, keys, blocksize, nselectedvar, selectedvar, &i, idx));
310942dc13f1SHong Zhang   }
311042dc13f1SHong Zhang 
311142dc13f1SHong Zhang   /* Create is */
31129566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, nidx, idx, PETSC_COPY_VALUES, is));
31139566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
31143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
311542dc13f1SHong Zhang }
311642dc13f1SHong Zhang 
3117d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode DMISComputeLocalIdx_private(DM dm, PetscInt p, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *selectedvar[], PetscInt *ii, PetscInt *idx)
3118d71ae5a4SJacob Faibussowitsch {
311942dc13f1SHong Zhang   PetscInt                 i, j, ncomps, nvar, key, offsetl, k, k1, offset = 0;
312042dc13f1SHong Zhang   DM_Network              *network = (DM_Network *)dm->data;
312142dc13f1SHong Zhang   DMNetworkComponentHeader header;
312242dc13f1SHong Zhang 
312342dc13f1SHong Zhang   PetscFunctionBegin;
31249566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffset(network->DataSection, p, &offset));
312542dc13f1SHong Zhang   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray + offset))->ndata;
312642dc13f1SHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray + offset);
312742dc13f1SHong Zhang 
312842dc13f1SHong Zhang   for (i = 0; i < ncomps; i++) {
312942dc13f1SHong Zhang     key  = header->key[i];
313042dc13f1SHong Zhang     nvar = header->nvar[i];
313142dc13f1SHong Zhang     for (j = 0; j < numkeys; j++) {
313242dc13f1SHong Zhang       if (key != keys[j]) continue;
313342dc13f1SHong Zhang 
31349566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetLocalVecOffset(dm, p, i, &offsetl));
313542dc13f1SHong Zhang       if (!blocksize || blocksize[j] == -1) {
313642dc13f1SHong Zhang         for (k = 0; k < nvar; k++) idx[(*ii)++] = offsetl + k;
313742dc13f1SHong Zhang       } else {
313842dc13f1SHong Zhang         for (k = 0; k < nvar; k += blocksize[j]) {
313942dc13f1SHong Zhang           for (k1 = 0; k1 < nselectedvar[j]; k1++) idx[(*ii)++] = offsetl + k + selectedvar[j][k1];
314042dc13f1SHong Zhang         }
314142dc13f1SHong Zhang       }
314242dc13f1SHong Zhang     }
314342dc13f1SHong Zhang   }
31443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
314542dc13f1SHong Zhang }
314642dc13f1SHong Zhang 
314742dc13f1SHong Zhang /*@
314842dc13f1SHong Zhang   DMNetworkCreateLocalIS - Create an index set object from the local vector of the network
314942dc13f1SHong Zhang 
315020f4b53cSBarry Smith   Not Collective
315142dc13f1SHong Zhang 
315242dc13f1SHong Zhang   Input Parameters:
315320f4b53cSBarry Smith + dm           - `DMNETWORK` object
315442dc13f1SHong Zhang . numkeys      - number of keys
315542dc13f1SHong Zhang . keys         - array of keys that define the components of the variables you wish to extract
315642dc13f1SHong Zhang . blocksize    - block size of the variables associated to the component
315742dc13f1SHong Zhang . nselectedvar - number of variables in each block to select
315842dc13f1SHong Zhang - selectedvar  - the offset into the block of each variable in each block to select
315942dc13f1SHong Zhang 
31602fe279fdSBarry Smith   Output Parameter:
316142dc13f1SHong Zhang . is - the index set
316242dc13f1SHong Zhang 
316320f4b53cSBarry Smith   Level: advanced
316442dc13f1SHong Zhang 
316542dc13f1SHong Zhang   Notes:
316620f4b53cSBarry Smith   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.
316742dc13f1SHong Zhang 
31685f25b224SDuncan Campbell .seealso: `DM`, `DMNETWORK`, `DMNetworkCreate()`, `DMNetworkCreateIS()`, `ISCreateGeneral()`
316942dc13f1SHong Zhang @*/
3170d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkCreateLocalIS(DM dm, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *selectedvar[], IS *is)
3171d71ae5a4SJacob Faibussowitsch {
317242dc13f1SHong Zhang   DM_Network *network = (DM_Network *)dm->data;
317342dc13f1SHong Zhang   PetscInt    i, p, pstart, pend, nidx, *idx;
317442dc13f1SHong Zhang 
317542dc13f1SHong Zhang   PetscFunctionBegin;
317642dc13f1SHong Zhang   /* Check input parameters */
317742dc13f1SHong Zhang   for (i = 0; i < numkeys; i++) {
317842dc13f1SHong Zhang     if (!blocksize || blocksize[i] == -1) continue;
317963a3b9bcSJacob 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]);
318042dc13f1SHong Zhang   }
318142dc13f1SHong Zhang 
3182daad07d3SAidan Hamilton   pstart = network->cloneshared->pStart;
3183daad07d3SAidan Hamilton   pend   = network->cloneshared->pEnd;
318442dc13f1SHong Zhang 
318542dc13f1SHong Zhang   /* Get local number of idx */
318642dc13f1SHong Zhang   nidx = 0;
318748a46eb9SPierre Jolivet   for (p = pstart; p < pend; p++) PetscCall(DMISAddSize_private(network, p, numkeys, keys, blocksize, nselectedvar, &nidx));
318842dc13f1SHong Zhang 
318942dc13f1SHong Zhang   /* Compute local idx */
31909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nidx, &idx));
319142dc13f1SHong Zhang   i = 0;
319248a46eb9SPierre Jolivet   for (p = pstart; p < pend; p++) PetscCall(DMISComputeLocalIdx_private(dm, p, numkeys, keys, blocksize, nselectedvar, selectedvar, &i, idx));
319342dc13f1SHong Zhang 
319442dc13f1SHong Zhang   /* Create is */
31959566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nidx, idx, PETSC_COPY_VALUES, is));
31969566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
31973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
319842dc13f1SHong Zhang }
3199daad07d3SAidan Hamilton /*@
3200d5b43468SJose E. Roman   DMNetworkFinalizeComponents - Sets up internal data structures for the sections and components. It is called after registering new components and adding all components
3201daad07d3SAidan Hamilton   to the cloned network. After calling this subroutine, no new components can be added to the network.
3202daad07d3SAidan Hamilton 
3203daad07d3SAidan Hamilton   Collective
3204daad07d3SAidan Hamilton 
32052fe279fdSBarry Smith   Input Parameter:
320620f4b53cSBarry Smith . dm - the `DMNETWORK` object
3207daad07d3SAidan Hamilton 
3208daad07d3SAidan Hamilton   Level: beginner
3209daad07d3SAidan Hamilton 
321020f4b53cSBarry Smith .seealso: `DM`, `DMNETWORK`, `DMNetworkAddComponent()`, `DMNetworkRegisterComponent()`, `DMSetUp()`
3211daad07d3SAidan Hamilton @*/
3212d71ae5a4SJacob Faibussowitsch PetscErrorCode DMNetworkFinalizeComponents(DM dm)
3213d71ae5a4SJacob Faibussowitsch {
3214daad07d3SAidan Hamilton   DM_Network *network = (DM_Network *)dm->data;
3215daad07d3SAidan Hamilton 
3216daad07d3SAidan Hamilton   PetscFunctionBegin;
32173ba16761SJacob Faibussowitsch   if (network->componentsetup) PetscFunctionReturn(PETSC_SUCCESS);
3218daad07d3SAidan Hamilton   PetscCall(DMNetworkComponentSetUp(dm));
3219daad07d3SAidan Hamilton   PetscCall(DMNetworkVariablesSetUp(dm));
3220daad07d3SAidan Hamilton   PetscCall(DMSetLocalSection(network->plex, network->DofSection));
3221daad07d3SAidan Hamilton   PetscCall(DMGetGlobalSection(network->plex, &network->GlobalDofSection));
3222daad07d3SAidan Hamilton   network->componentsetup = PETSC_TRUE;
32233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3224daad07d3SAidan Hamilton }
3225