Actual source code: network.c
1: #include <petsc/private/dmnetworkimpl.h>
3: PetscLogEvent DMNetwork_LayoutSetUp;
4: PetscLogEvent DMNetwork_SetUpNetwork;
5: PetscLogEvent DMNetwork_Distribute;
7: /*
8: Creates the component header and value objects for a network point
9: */
10: static PetscErrorCode SetUpNetworkHeaderComponentValue(DM dm, DMNetworkComponentHeader header, DMNetworkComponentValue cvalue)
11: {
12: PetscFunctionBegin;
13: /* Allocate arrays for component information */
14: PetscCall(PetscCalloc5(header->maxcomps, &header->size, header->maxcomps, &header->key, header->maxcomps, &header->offset, header->maxcomps, &header->nvar, header->maxcomps, &header->offsetvarrel));
15: PetscCall(PetscCalloc1(header->maxcomps, &cvalue->data));
17: /* 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.
18: If the data header struct changes then this header size calculation needs to be updated. */
19: header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5 * header->maxcomps * sizeof(PetscInt);
20: header->hsize /= sizeof(DMNetworkComponentGenericDataType);
21: PetscFunctionReturn(PETSC_SUCCESS);
22: }
24: PetscErrorCode DMNetworkInitializeHeaderComponentData(DM dm)
25: {
26: DM_Network *network = (DM_Network *)dm->data;
27: PetscInt np, p, defaultnumcomp = DMNETWORK_MAX_COMP_AT_POINT_DEFAULT;
29: PetscFunctionBegin;
30: np = network->cloneshared->pEnd - network->cloneshared->pStart;
31: if (network->header)
32: for (p = 0; p < np; p++) {
33: network->header[p].maxcomps = defaultnumcomp;
34: PetscCall(SetUpNetworkHeaderComponentValue(dm, &network->header[p], &network->cvalue[p]));
35: }
37: PetscFunctionReturn(PETSC_SUCCESS);
38: }
39: /*@
40: DMNetworkGetPlex - Gets the Plex DM associated with this network DM
42: Not collective
44: Input Parameters:
45: . dm - the dm object
47: Output Parameters:
48: . plexdm - the plex dm object
50: Level: Advanced
52: .seealso: `DMNetworkCreate()`
53: @*/
54: PetscErrorCode DMNetworkGetPlex(DM dm, DM *plexdm)
55: {
56: DM_Network *network = (DM_Network *)dm->data;
58: PetscFunctionBegin;
59: *plexdm = network->plex;
60: PetscFunctionReturn(PETSC_SUCCESS);
61: }
63: /*@
64: DMNetworkGetNumSubNetworks - Gets the the number of subnetworks
66: Not collective
68: Input Parameter:
69: . dm - the dm object
71: Output Parameters:
72: + nsubnet - local number of subnetworks
73: - Nsubnet - global number of subnetworks
75: Level: beginner
77: .seealso: `DMNetworkCreate()`, `DMNetworkSetNumSubNetworks()`
78: @*/
79: PetscErrorCode DMNetworkGetNumSubNetworks(DM dm, PetscInt *nsubnet, PetscInt *Nsubnet)
80: {
81: DM_Network *network = (DM_Network *)dm->data;
83: PetscFunctionBegin;
84: if (nsubnet) *nsubnet = network->cloneshared->nsubnet;
85: if (Nsubnet) *Nsubnet = network->cloneshared->Nsubnet;
86: PetscFunctionReturn(PETSC_SUCCESS);
87: }
89: /*@
90: DMNetworkSetNumSubNetworks - Sets the number of subnetworks
92: Collective on dm
94: Input Parameters:
95: + dm - the dm object
96: . nsubnet - local number of subnetworks
97: - Nsubnet - global number of subnetworks
99: Level: beginner
101: .seealso: `DMNetworkCreate()`, `DMNetworkGetNumSubNetworks()`
102: @*/
103: PetscErrorCode DMNetworkSetNumSubNetworks(DM dm, PetscInt nsubnet, PetscInt Nsubnet)
104: {
105: DM_Network *network = (DM_Network *)dm->data;
107: PetscFunctionBegin;
108: PetscCheck(network->cloneshared->Nsubnet == 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "Network sizes already set, cannot resize the network");
114: if (Nsubnet == PETSC_DECIDE) {
115: PetscCheck(nsubnet >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of local subnetworks %" PetscInt_FMT " cannot be less than 0", nsubnet);
116: PetscCall(MPIU_Allreduce(&nsubnet, &Nsubnet, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
117: }
118: PetscCheck(Nsubnet >= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "Number of global subnetworks %" PetscInt_FMT " cannot be less than 1", Nsubnet);
120: network->cloneshared->Nsubnet = Nsubnet;
121: network->cloneshared->nsubnet = 0; /* initial value; will be determined by DMNetworkAddSubnetwork() */
122: PetscCall(PetscCalloc1(Nsubnet, &network->cloneshared->subnet));
124: /* num of shared vertices */
125: network->cloneshared->nsvtx = 0;
126: network->cloneshared->Nsvtx = 0;
127: PetscFunctionReturn(PETSC_SUCCESS);
128: }
130: /*@
131: DMNetworkAddSubnetwork - Add a subnetwork
133: Collective on dm
135: Input Parameters:
136: + dm - the dm object
137: . name - name of the subnetwork
138: . ne - number of local edges of this subnetwork
139: - edgelist - list of edges for this subnetwork, this is a one dimensional array with pairs of entries being the two vertices (in global numbering of the vertices) of each edge,
140: $ [first vertex of first edge, second vertex of first edge, first vertex of second edge, second vertex of second edge, etc]
142: Output Parameters:
143: . netnum - global index of the subnetwork
145: Notes:
146: There is no copy involved in this operation, only the pointer is referenced. The edgelist should
147: not be destroyed before the call to DMNetworkLayoutSetUp()
149: A network can comprise of a single subnetwork OR multiple subnetworks. For a single subnetwork, the subnetwork can be read either in serial or parallel. For a multiple subnetworks,
150: 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.
152: Level: beginner
154: Example usage:
155: Consider the following networks:
156: 1) A single subnetwork:
157: .vb
158: network 0:
159: rank[0]:
160: v0 --> v2; v1 --> v2
161: rank[1]:
162: v3 --> v5; v4 --> v5
163: .ve
165: The resulting input
166: network 0:
167: rank[0]:
168: ne = 2
169: edgelist = [0 2 | 1 2]
170: rank[1]:
171: ne = 2
172: edgelist = [3 5 | 4 5]
174: 2) Two subnetworks:
175: .vb
176: subnetwork 0:
177: rank[0]:
178: v0 --> v2; v2 --> v1; v1 --> v3;
179: subnetwork 1:
180: rank[1]:
181: v0 --> v3; v3 --> v2; v2 --> v1;
182: .ve
184: The resulting input
185: subnetwork 0:
186: rank[0]:
187: ne = 3
188: edgelist = [0 2 | 2 1 | 1 3]
189: rank[1]:
190: ne = 0
191: edgelist = NULL
193: subnetwork 1:
194: rank[0]:
195: ne = 0
196: edgelist = NULL
197: rank[1]:
198: edgelist = [0 3 | 3 2 | 2 1]
200: .seealso: `DMNetworkCreate()`, `DMNetworkSetNumSubnetworks()`
201: @*/
202: PetscErrorCode DMNetworkAddSubnetwork(DM dm, const char *name, PetscInt ne, PetscInt edgelist[], PetscInt *netnum)
203: {
204: DM_Network *network = (DM_Network *)dm->data;
205: PetscInt i, Nedge, j, Nvtx, nvtx, nvtx_min = -1, nvtx_max = 0;
206: PetscBT table;
208: PetscFunctionBegin;
209: 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]);
211: i = 0;
212: if (ne) nvtx_min = nvtx_max = edgelist[0];
213: for (j = 0; j < ne; j++) {
214: nvtx_min = PetscMin(nvtx_min, edgelist[i]);
215: nvtx_max = PetscMax(nvtx_max, edgelist[i]);
216: i++;
217: nvtx_min = PetscMin(nvtx_min, edgelist[i]);
218: nvtx_max = PetscMax(nvtx_max, edgelist[i]);
219: i++;
220: }
221: Nvtx = nvtx_max - nvtx_min + 1; /* approximated total local nvtx for this subnet */
223: /* Get exact local nvtx for this subnet: counting local values between nvtx_min and nvtx_max */
224: PetscCall(PetscBTCreate(Nvtx, &table));
225: PetscCall(PetscBTMemzero(Nvtx, table));
226: i = 0;
227: for (j = 0; j < ne; j++) {
228: PetscCall(PetscBTSet(table, edgelist[i++] - nvtx_min));
229: PetscCall(PetscBTSet(table, edgelist[i++] - nvtx_min));
230: }
231: nvtx = 0;
232: for (j = 0; j < Nvtx; j++) {
233: if (PetscBTLookup(table, j)) nvtx++;
234: }
236: /* Get global total Nvtx = max(edgelist[])+1 for this subnet */
237: PetscCall(MPIU_Allreduce(&nvtx_max, &Nvtx, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
238: Nvtx++;
239: PetscCall(PetscBTDestroy(&table));
241: /* Get global total Nedge for this subnet */
242: PetscCall(MPIU_Allreduce(&ne, &Nedge, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
244: i = network->cloneshared->nsubnet;
245: if (name) PetscCall(PetscStrcpy(network->cloneshared->subnet[i].name, name));
246: network->cloneshared->subnet[i].nvtx = nvtx; /* include ghost vertices */
247: network->cloneshared->subnet[i].nedge = ne;
248: network->cloneshared->subnet[i].edgelist = edgelist;
249: network->cloneshared->subnet[i].Nvtx = Nvtx;
250: network->cloneshared->subnet[i].Nedge = Nedge;
252: /* ----------------------------------------------------------
253: p=v or e;
254: subnet[0].pStart = 0
255: subnet[i+1].pStart = subnet[i].pEnd = subnet[i].pStart + (nE[i] or NV[i])
256: ----------------------------------------------------------------------- */
257: /* GLOBAL subnet[].vStart and vEnd, used by DMNetworkLayoutSetUp() */
258: network->cloneshared->subnet[i].vStart = network->cloneshared->NVertices;
259: network->cloneshared->subnet[i].vEnd = network->cloneshared->subnet[i].vStart + network->cloneshared->subnet[i].Nvtx; /* global vEnd of subnet[i] */
261: network->cloneshared->nVertices += nvtx; /* include ghost vertices */
262: network->cloneshared->NVertices += network->cloneshared->subnet[i].Nvtx;
264: /* LOCAL subnet[].eStart and eEnd, used by DMNetworkLayoutSetUp() */
265: network->cloneshared->subnet[i].eStart = network->cloneshared->nEdges;
266: network->cloneshared->subnet[i].eEnd = network->cloneshared->subnet[i].eStart + ne;
267: network->cloneshared->nEdges += ne;
268: network->cloneshared->NEdges += network->cloneshared->subnet[i].Nedge;
270: PetscCall(PetscStrcpy(network->cloneshared->subnet[i].name, name));
271: if (netnum) *netnum = network->cloneshared->nsubnet;
272: network->cloneshared->nsubnet++;
273: PetscFunctionReturn(PETSC_SUCCESS);
274: }
276: /*@C
277: DMNetworkSharedVertexGetInfo - Get info of a shared vertex struct, see petsc/private/dmnetworkimpl.h
279: Not collective
281: Input Parameters:
282: + dm - the DM object
283: - v - vertex point
285: Output Parameters:
286: + gidx - global number of this shared vertex in the internal dmplex
287: . n - number of subnetworks that share this vertex
288: - sv - array of size n: sv[2*i,2*i+1]=(net[i], idx[i]), i=0,...,n-1
290: Level: intermediate
292: .seealso: `DMNetworkGetSharedVertices()`
293: @*/
294: PetscErrorCode DMNetworkSharedVertexGetInfo(DM dm, PetscInt v, PetscInt *gidx, PetscInt *n, const PetscInt **sv)
295: {
296: DM_Network *network = (DM_Network *)dm->data;
297: SVtx *svtx = network->cloneshared->svtx;
298: PetscInt i, gidx_tmp;
300: PetscFunctionBegin;
301: PetscCall(DMNetworkGetGlobalVertexIndex(dm, v, &gidx_tmp));
302: PetscCall(PetscHMapIGetWithDefault(network->cloneshared->svtable, gidx_tmp + 1, 0, &i));
303: PetscCheck(i > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "input vertex is not a shared vertex");
305: i--;
306: if (gidx) *gidx = gidx_tmp;
307: if (n) *n = svtx[i].n;
308: if (sv) *sv = svtx[i].sv;
309: PetscFunctionReturn(PETSC_SUCCESS);
310: }
312: /*
313: VtxGetInfo - Get info of an input vertex=(net,idx)
315: Input Parameters:
316: + Nsvtx - global num of shared vertices
317: . svtx - array of shared vertices (global)
318: - (net,idx) - subnet number and local index for a vertex
320: Output Parameters:
321: + gidx - global index of (net,idx)
322: . svtype - see petsc/private/dmnetworkimpl.h
323: - svtx_idx - ordering in the svtx array
324: */
325: static inline PetscErrorCode VtxGetInfo(PetscInt Nsvtx, SVtx *svtx, PetscInt net, PetscInt idx, PetscInt *gidx, SVtxType *svtype, PetscInt *svtx_idx)
326: {
327: PetscInt i, j, *svto, g_idx;
328: SVtxType vtype;
330: PetscFunctionBegin;
331: if (!Nsvtx) PetscFunctionReturn(PETSC_SUCCESS);
333: g_idx = -1;
334: vtype = SVNONE;
336: for (i = 0; i < Nsvtx; i++) {
337: if (net == svtx[i].sv[0] && idx == svtx[i].sv[1]) {
338: g_idx = svtx[i].gidx;
339: vtype = SVFROM;
340: } else { /* loop over svtx[i].n */
341: for (j = 1; j < svtx[i].n; j++) {
342: svto = svtx[i].sv + 2 * j;
343: if (net == svto[0] && idx == svto[1]) {
344: /* input vertex net.idx is a shared to_vertex, output its global index and its svtype */
345: g_idx = svtx[i].gidx; /* output gidx for to_vertex */
346: vtype = SVTO;
347: }
348: }
349: }
350: if (vtype != SVNONE) break;
351: }
352: if (gidx) *gidx = g_idx;
353: if (svtype) *svtype = vtype;
354: if (svtx_idx) *svtx_idx = i;
355: PetscFunctionReturn(PETSC_SUCCESS);
356: }
358: /*
359: TableAddSVtx - Add a new shared vertice from sedgelist[k] to a ctable svta
361: Input: network, sedgelist, k, svta
362: Output: svta, tdata, ta2sv
363: */
364: static inline PetscErrorCode TableAddSVtx(DM_Network *network, PetscInt *sedgelist, PetscInt k, PetscHMapI svta, PetscInt *tdata, PetscInt *ta2sv)
365: {
366: PetscInt net, idx, gidx;
368: PetscFunctionBegin;
369: net = sedgelist[k];
370: idx = sedgelist[k + 1];
371: gidx = network->cloneshared->subnet[net].vStart + idx;
372: PetscCall(PetscHMapISet(svta, gidx + 1, *tdata + 1));
374: ta2sv[*tdata] = k; /* maps tdata to index of sedgelist */
375: (*tdata)++;
376: PetscFunctionReturn(PETSC_SUCCESS);
377: }
379: /*
380: SharedVtxCreate - Create an array of global shared vertices. See SVtx and SVtxType in dmnetworkimpl.h
382: Input: dm, Nsedgelist, sedgelist
384: Note: Output svtx is organized as
385: sv(net[0],idx[0]) --> sv(net[1],idx[1])
386: --> sv(net[1],idx[1])
387: ...
388: --> sv(net[n-1],idx[n-1])
389: and net[0] < net[1] < ... < net[n-1]
390: where sv[0] has SVFROM type, sv[i], i>0, has SVTO type.
391: */
392: static PetscErrorCode SharedVtxCreate(DM dm, PetscInt Nsedgelist, PetscInt *sedgelist)
393: {
394: SVtx *svtx = NULL;
395: PetscInt *sv, k, j, nsv, *tdata, **ta2sv;
396: PetscHMapI *svtas;
397: PetscInt gidx, net, idx, i, nta, ita, idx_from, idx_to, n, *net_tmp, *idx_tmp, *gidx_tmp;
398: DM_Network *network = (DM_Network *)dm->data;
399: PetscHashIter ppos;
401: PetscFunctionBegin;
402: /* (1) Crete an array of ctables svtas to map (net,idx) -> gidx; a svtas[] for a shared/merged vertex */
403: PetscCall(PetscCalloc3(Nsedgelist, &svtas, Nsedgelist, &tdata, 2 * Nsedgelist, &ta2sv));
405: k = 0; /* sedgelist vertex counter j = 4*k */
406: nta = 0; /* num of svta tables created = num of groups of shared vertices */
408: /* for j=0 */
409: PetscCall(PetscHMapICreateWithSize(2 * Nsedgelist, svtas + nta));
410: PetscCall(PetscMalloc1(2 * Nsedgelist, &ta2sv[nta]));
412: PetscCall(TableAddSVtx(network, sedgelist, k, svtas[nta], &tdata[nta], ta2sv[nta]));
413: PetscCall(TableAddSVtx(network, sedgelist, k + 2, svtas[nta], &tdata[nta], ta2sv[nta]));
414: nta++;
415: k += 4;
417: for (j = 1; j < Nsedgelist; j++) { /* j: sedgelist counter */
418: for (ita = 0; ita < nta; ita++) {
419: /* vfrom */
420: net = sedgelist[k];
421: idx = sedgelist[k + 1];
422: gidx = network->cloneshared->subnet[net].vStart + idx; /* global index of the vertex net.idx before merging shared vertices */
423: PetscCall(PetscHMapIGetWithDefault(svtas[ita], gidx + 1, 0, &idx_from));
425: /* vto */
426: net = sedgelist[k + 2];
427: idx = sedgelist[k + 3];
428: gidx = network->cloneshared->subnet[net].vStart + idx;
429: PetscCall(PetscHMapIGetWithDefault(svtas[ita], gidx + 1, 0, &idx_to));
431: if (idx_from || idx_to) { /* vfrom or vto is on table svtas[ita] */
432: idx_from--;
433: idx_to--;
434: if (idx_from < 0) { /* vto is on svtas[ita] */
435: PetscCall(TableAddSVtx(network, sedgelist, k, svtas[ita], &tdata[ita], ta2sv[ita]));
436: break;
437: } else if (idx_to < 0) {
438: PetscCall(TableAddSVtx(network, sedgelist, k + 2, svtas[ita], &tdata[ita], ta2sv[ita]));
439: break;
440: }
441: }
442: }
444: if (ita == nta) {
445: PetscCall(PetscHMapICreateWithSize(2 * Nsedgelist, svtas + nta));
446: PetscCall(PetscMalloc1(2 * Nsedgelist, &ta2sv[nta]));
448: PetscCall(TableAddSVtx(network, sedgelist, k, svtas[nta], &tdata[nta], ta2sv[nta]));
449: PetscCall(TableAddSVtx(network, sedgelist, k + 2, svtas[nta], &tdata[nta], ta2sv[nta]));
450: nta++;
451: }
452: k += 4;
453: }
455: /* (2) Create svtable for query shared vertices using gidx */
456: PetscCall(PetscHMapICreateWithSize(nta, &network->cloneshared->svtable));
458: /* (3) Construct svtx from svtas
459: svtx: array of SVtx: sv[0]=(net[0],idx[0]) to vertices sv[k], k=1,...,n-1. */
460: PetscCall(PetscMalloc1(nta, &svtx));
461: for (nsv = 0; nsv < nta; nsv++) {
462: /* for a single svtx, put shared vertices in ascending order of gidx */
463: PetscCall(PetscHMapIGetSize(svtas[nsv], &n));
464: PetscCall(PetscCalloc1(2 * n, &sv));
465: PetscCall(PetscMalloc3(n, &gidx_tmp, n, &net_tmp, n, &idx_tmp));
466: svtx[nsv].sv = sv;
467: svtx[nsv].n = n;
468: svtx[nsv].gidx = network->cloneshared->NVertices; /* initialization */
470: PetscHashIterBegin(svtas[nsv], ppos);
471: for (k = 0; k < n; k++) { /* gidx is sorted in ascending order */
472: PetscHashIterGetKey(svtas[nsv], ppos, gidx);
473: PetscHashIterGetVal(svtas[nsv], ppos, i);
474: PetscHashIterNext(svtas[nsv], ppos);
475: gidx--;
476: i--;
477: j = ta2sv[nsv][i]; /* maps i to index of sedgelist */
478: net_tmp[k] = sedgelist[j]; /* subnet number */
479: idx_tmp[k] = sedgelist[j + 1]; /* index on the subnet */
480: gidx_tmp[k] = gidx; /* gidx in un-merged dmnetwork */
481: }
483: /* current implementation requires sv[]=[net,idx] in ascending order of its gidx in un-merged dmnetwork */
484: PetscCall(PetscSortIntWithArrayPair(n, gidx_tmp, net_tmp, idx_tmp));
485: svtx[nsv].gidx = gidx_tmp[0]; /* = min(gidx) */
486: for (k = 0; k < n; k++) {
487: sv[2 * k] = net_tmp[k];
488: sv[2 * k + 1] = idx_tmp[k];
489: }
490: PetscCall(PetscFree3(gidx_tmp, net_tmp, idx_tmp));
492: /* Setup svtable for query shared vertices */
493: PetscCall(PetscHMapISet(network->cloneshared->svtable, svtx[nsv].gidx + 1, nsv + 1));
494: }
496: for (j = 0; j < nta; j++) {
497: PetscCall(PetscHMapIDestroy(svtas + j));
498: PetscCall(PetscFree(ta2sv[j]));
499: }
500: PetscCall(PetscFree3(svtas, tdata, ta2sv));
502: network->cloneshared->Nsvtx = nta;
503: network->cloneshared->svtx = svtx;
504: PetscFunctionReturn(PETSC_SUCCESS);
505: }
507: /*
508: GetEdgelist_Coupling - Get an integrated edgelist for dmplex from user-provided subnet[].edgelist when subnets are coupled by shared vertices
510: Input Parameters:
511: . dm - the dmnetwork object
513: Output Parameters:
514: + edges - the integrated edgelist for dmplex
515: - nmerged_ptr - num of vertices being merged
516: */
517: static PetscErrorCode GetEdgelist_Coupling(DM dm, PetscInt *edges, PetscInt *nmerged_ptr)
518: {
519: MPI_Comm comm;
520: PetscMPIInt size, rank, *recvcounts = NULL, *displs = NULL;
521: DM_Network *network = (DM_Network *)dm->data;
522: PetscInt i, j, ctr, np;
523: PetscInt *vidxlTog, Nsv, Nsubnet = network->cloneshared->Nsubnet;
524: PetscInt *sedgelist = network->cloneshared->sedgelist;
525: PetscInt net, idx, gidx, nmerged, *vrange, gidx_from, net_from, sv_idx;
526: SVtxType svtype = SVNONE;
527: SVtx *svtx;
529: PetscFunctionBegin;
530: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
531: PetscCallMPI(MPI_Comm_rank(comm, &rank));
532: PetscCallMPI(MPI_Comm_size(comm, &size));
534: /* (1) Create global svtx[] from sedgelist */
535: /* --------------------------------------- */
536: PetscCall(SharedVtxCreate(dm, network->cloneshared->Nsvtx, sedgelist));
537: Nsv = network->cloneshared->Nsvtx;
538: svtx = network->cloneshared->svtx;
540: /* (2) Merge shared vto vertices to their vfrom vertex with same global vertex index (gidx) */
541: /* --------------------------------------------------------------------------------------- */
542: /* (2.1) compute vrage[rank]: global index of 1st local vertex in proc[rank] */
543: PetscCall(PetscMalloc4(size + 1, &vrange, size, &displs, size, &recvcounts, network->cloneshared->nVertices, &vidxlTog));
544: for (i = 0; i < size; i++) {
545: displs[i] = i;
546: recvcounts[i] = 1;
547: }
549: vrange[0] = 0;
550: PetscCallMPI(MPI_Allgatherv(&network->cloneshared->nVertices, 1, MPIU_INT, vrange + 1, recvcounts, displs, MPIU_INT, comm));
551: for (i = 2; i < size + 1; i++) vrange[i] += vrange[i - 1];
553: /* (2.2) Create vidxlTog: maps UN-MERGED local vertex index i to global index gidx (plex, excluding ghost vertices) */
554: i = 0;
555: gidx = 0;
556: nmerged = 0; /* local num of merged vertices */
557: network->cloneshared->nsvtx = 0; /* local num of SVtx structs, including ghosts */
558: for (net = 0; net < Nsubnet; net++) {
559: for (idx = 0; idx < network->cloneshared->subnet[net].Nvtx; idx++) { /* Note: global subnet[net].Nvtx */
560: PetscCall(VtxGetInfo(Nsv, svtx, net, idx, &gidx_from, &svtype, &sv_idx));
561: if (svtype == SVTO) {
562: if (network->cloneshared->subnet[net].nvtx) { /* this proc owns sv_to */
563: net_from = svtx[sv_idx].sv[0]; /* subnet number of its shared vertex */
564: if (network->cloneshared->subnet[net_from].nvtx == 0) {
565: /* this proc does not own v_from, thus a ghost local vertex */
566: network->cloneshared->nsvtx++;
567: }
568: vidxlTog[i++] = gidx_from; /* gidx before merging! Bug??? */
569: nmerged++; /* a shared vertex -- merged */
570: }
571: } else {
572: if (svtype == SVFROM && network->cloneshared->subnet[net].nvtx) {
573: /* this proc owns this v_from, a new local shared vertex */
574: network->cloneshared->nsvtx++;
575: }
576: if (network->cloneshared->subnet[net].nvtx) vidxlTog[i++] = gidx;
577: gidx++;
578: }
579: }
580: }
581: #if defined(PETSC_USE_DEBUG)
582: PetscCheck(i == network->cloneshared->nVertices, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "%" PetscInt_FMT " != %" PetscInt_FMT " nVertices", i, network->cloneshared->nVertices);
583: #endif
585: /* (2.3) Shared vertices in the subnetworks are merged, update global NVertices: np = sum(local nmerged) */
586: PetscCallMPI(MPI_Allreduce(&nmerged, &np, 1, MPIU_INT, MPI_SUM, comm));
587: network->cloneshared->NVertices -= np;
589: ctr = 0;
590: for (net = 0; net < Nsubnet; net++) {
591: for (j = 0; j < network->cloneshared->subnet[net].nedge; j++) {
592: /* vfrom: */
593: i = network->cloneshared->subnet[net].edgelist[2 * j] + (network->cloneshared->subnet[net].vStart - vrange[rank]);
594: edges[2 * ctr] = vidxlTog[i];
596: /* vto */
597: i = network->cloneshared->subnet[net].edgelist[2 * j + 1] + (network->cloneshared->subnet[net].vStart - vrange[rank]);
598: edges[2 * ctr + 1] = vidxlTog[i];
599: ctr++;
600: }
601: }
602: PetscCall(PetscFree4(vrange, displs, recvcounts, vidxlTog));
603: PetscCall(PetscFree(sedgelist)); /* created in DMNetworkAddSharedVertices() */
605: *nmerged_ptr = nmerged;
606: PetscFunctionReturn(PETSC_SUCCESS);
607: }
609: PetscErrorCode DMNetworkInitializeNonTopological(DM dm)
610: {
611: DM_Network *network = (DM_Network *)dm->data;
612: PetscInt p, pStart = network->cloneshared->pStart, pEnd = network->cloneshared->pEnd;
613: MPI_Comm comm;
615: PetscFunctionBegin;
616: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
618: PetscCall(PetscSectionCreate(comm, &network->DataSection));
619: PetscCall(PetscSectionCreate(comm, &network->DofSection));
620: PetscCall(PetscSectionSetChart(network->DataSection, pStart, pEnd));
621: PetscCall(PetscSectionSetChart(network->DofSection, pStart, pEnd));
623: PetscCall(DMNetworkInitializeHeaderComponentData(dm));
625: for (p = 0; p < pEnd - pStart; p++) {
626: network->header[p].ndata = 0;
627: network->header[p].offset[0] = 0;
628: network->header[p].offsetvarrel[0] = 0;
629: PetscCall(PetscSectionAddDof(network->DataSection, p, network->header[p].hsize));
630: }
631: PetscFunctionReturn(PETSC_SUCCESS);
632: }
634: /*@
635: DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
637: Not Collective
639: Input Parameters:
640: . dm - the dmnetwork object
642: Notes:
643: This routine should be called after the network sizes and edgelists have been provided. It creates
644: the bare layout of the network and sets up the network to begin insertion of components.
646: All the components should be registered before calling this routine.
648: Level: beginner
650: .seealso: `DMNetworkSetNumSubNetworks()`, `DMNetworkAddSubnetwork()`
651: @*/
652: PetscErrorCode DMNetworkLayoutSetUp(DM dm)
653: {
654: DM_Network *network = (DM_Network *)dm->data;
655: PetscInt i, j, ctr, Nsubnet = network->cloneshared->Nsubnet, np, *edges, *subnetvtx, *subnetedge, e, v, vfrom, vto, net, globaledgeoff;
656: const PetscInt *cone;
657: MPI_Comm comm;
658: PetscMPIInt size;
659: PetscSection sectiong;
660: PetscInt nmerged = 0;
662: PetscFunctionBegin;
663: PetscCall(PetscLogEventBegin(DMNetwork_LayoutSetUp, dm, 0, 0, 0));
664: PetscCheck(network->cloneshared->nsubnet == Nsubnet, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Must call DMNetworkAddSubnetwork() %" PetscInt_FMT " times", Nsubnet);
666: /* This implementation requires user input each subnet by a single processor when Nsubnet>1, thus subnet[net].nvtx=subnet[net].Nvtx when net>0 */
667: for (net = 1; net < Nsubnet; net++) {
668: if (network->cloneshared->subnet[net].nvtx)
669: 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,
670: network->cloneshared->subnet[net].nvtx, network->cloneshared->subnet[net].Nvtx);
671: }
673: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
674: PetscCallMPI(MPI_Comm_size(comm, &size));
676: /* Create LOCAL edgelist in global vertex ordering for the network by concatenating local input edgelists of the subnetworks */
677: PetscCall(PetscCalloc1(2 * network->cloneshared->nEdges, &edges));
679: if (network->cloneshared->Nsvtx) { /* subnetworks are coupled via shared vertices */
680: PetscCall(GetEdgelist_Coupling(dm, edges, &nmerged));
681: } else { /* subnetworks are not coupled */
682: /* Create a 0-size svtable for query shared vertices */
683: PetscCall(PetscHMapICreate(&network->cloneshared->svtable));
684: ctr = 0;
685: for (i = 0; i < Nsubnet; i++) {
686: for (j = 0; j < network->cloneshared->subnet[i].nedge; j++) {
687: edges[2 * ctr] = network->cloneshared->subnet[i].vStart + network->cloneshared->subnet[i].edgelist[2 * j];
688: edges[2 * ctr + 1] = network->cloneshared->subnet[i].vStart + network->cloneshared->subnet[i].edgelist[2 * j + 1];
689: ctr++;
690: }
691: }
692: }
694: /* Create network->plex; One dimensional network, numCorners=2 */
695: PetscCall(DMCreate(comm, &network->plex));
696: PetscCall(DMSetType(network->plex, DMPLEX));
697: PetscCall(DMSetDimension(network->plex, 1));
699: if (size == 1) PetscCall(DMPlexBuildFromCellList(network->plex, network->cloneshared->nEdges, PETSC_DECIDE, 2, edges));
700: else PetscCall(DMPlexBuildFromCellListParallel(network->plex, network->cloneshared->nEdges, PETSC_DECIDE, PETSC_DECIDE, 2, edges, NULL, NULL));
702: PetscCall(DMPlexGetChart(network->plex, &network->cloneshared->pStart, &network->cloneshared->pEnd));
703: PetscCall(DMPlexGetHeightStratum(network->plex, 0, &network->cloneshared->eStart, &network->cloneshared->eEnd));
704: PetscCall(DMPlexGetHeightStratum(network->plex, 1, &network->cloneshared->vStart, &network->cloneshared->vEnd));
705: np = network->cloneshared->pEnd - network->cloneshared->pStart;
706: PetscCall(PetscCalloc2(np, &network->header, np, &network->cvalue));
708: /* Create edge and vertex arrays for the subnetworks
709: This implementation assumes that DMNetwork reads
710: (1) a single subnetwork in parallel; or
711: (2) n subnetworks using n processors, one subnetwork/processor.
712: */
713: PetscCall(PetscCalloc2(network->cloneshared->nEdges, &subnetedge, network->cloneshared->nVertices + network->cloneshared->nsvtx, &subnetvtx)); /* Maps local edge/vertex to local subnetwork's edge/vertex */
714: network->cloneshared->subnetedge = subnetedge;
715: network->cloneshared->subnetvtx = subnetvtx;
716: for (j = 0; j < Nsubnet; j++) {
717: network->cloneshared->subnet[j].edges = subnetedge;
718: subnetedge += network->cloneshared->subnet[j].nedge;
720: network->cloneshared->subnet[j].vertices = subnetvtx;
721: subnetvtx += network->cloneshared->subnet[j].nvtx;
722: }
723: network->cloneshared->svertices = subnetvtx;
725: /* Get edge ownership */
726: np = network->cloneshared->eEnd - network->cloneshared->eStart;
727: PetscCallMPI(MPI_Scan(&np, &globaledgeoff, 1, MPIU_INT, MPI_SUM, comm));
728: globaledgeoff -= np;
730: /* Setup local edge and vertex arrays for subnetworks */
731: e = 0;
732: for (i = 0; i < Nsubnet; i++) {
733: ctr = 0;
734: for (j = 0; j < network->cloneshared->subnet[i].nedge; j++) {
735: /* edge e */
736: network->header[e].index = e + globaledgeoff; /* Global edge index */
737: network->header[e].subnetid = i;
738: network->cloneshared->subnet[i].edges[j] = e;
740: /* connected vertices */
741: PetscCall(DMPlexGetCone(network->plex, e, &cone));
743: /* vertex cone[0] */
744: v = cone[0];
745: network->header[v].index = edges[2 * e]; /* Global vertex index */
746: network->header[v].subnetid = i; /* Subnetwork id */
747: if (Nsubnet == 1) {
748: network->cloneshared->subnet[i].vertices[v - network->cloneshared->vStart] = v; /* user's subnet[].idx = petsc's v */
749: } else {
750: vfrom = network->cloneshared->subnet[i].edgelist[2 * ctr]; /* =subnet[i].idx, Global index! */
751: network->cloneshared->subnet[i].vertices[vfrom] = v; /* user's subnet[].dix = petsc's v */
752: }
754: /* vertex cone[1] */
755: v = cone[1];
756: network->header[v].index = edges[2 * e + 1]; /* Global vertex index */
757: network->header[v].subnetid = i; /* Subnetwork id */
758: if (Nsubnet == 1) {
759: network->cloneshared->subnet[i].vertices[v - network->cloneshared->vStart] = v; /* user's subnet[].idx = petsc's v */
760: } else {
761: vto = network->cloneshared->subnet[i].edgelist[2 * ctr + 1]; /* =subnet[i].idx, Global index! */
762: network->cloneshared->subnet[i].vertices[vto] = v; /* user's subnet[].dix = petsc's v */
763: }
765: e++;
766: ctr++;
767: }
768: }
769: PetscCall(PetscFree(edges));
771: /* Set local vertex array for the subnetworks */
772: j = 0;
773: for (v = network->cloneshared->vStart; v < network->cloneshared->vEnd; v++) {
774: /* local shared vertex */
775: PetscCall(PetscHMapIGetWithDefault(network->cloneshared->svtable, network->header[v].index + 1, 0, &i));
776: if (i) network->cloneshared->svertices[j++] = v;
777: }
779: /* Create a global section to be used by DMNetworkIsGhostVertex() which is a non-collective routine */
780: /* see snes_tutorials_network-ex1_4 */
781: PetscCall(DMGetGlobalSection(network->plex, §iong));
782: /* Initialize non-topological data structures */
783: PetscCall(DMNetworkInitializeNonTopological(dm));
784: PetscCall(PetscLogEventEnd(DMNetwork_LayoutSetUp, dm, 0, 0, 0));
785: PetscFunctionReturn(PETSC_SUCCESS);
786: }
788: /*@C
789: DMNetworkGetSubnetwork - Returns the information about a requested subnetwork
791: Not collective
793: Input Parameters:
794: + dm - the DM object
795: - netnum - the global index of the subnetwork
797: Output Parameters:
798: + nv - number of vertices (local)
799: . ne - number of edges (local)
800: . vtx - local vertices of the subnetwork
801: - edge - local edges of the subnetwork
803: Notes:
804: Cannot call this routine before DMNetworkLayoutSetup()
806: The local vertices returned on each rank are determined by DMNetwork. The user does not have any control over what vertices are local.
808: Level: intermediate
810: .seealso: `DMNetworkCreate()`, `DMNetworkAddSubnetwork()`, `DMNetworkLayoutSetUp()`
811: @*/
812: PetscErrorCode DMNetworkGetSubnetwork(DM dm, PetscInt netnum, PetscInt *nv, PetscInt *ne, const PetscInt **vtx, const PetscInt **edge)
813: {
814: DM_Network *network = (DM_Network *)dm->data;
816: PetscFunctionBegin;
817: 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);
818: if (nv) *nv = network->cloneshared->subnet[netnum].nvtx;
819: if (ne) *ne = network->cloneshared->subnet[netnum].nedge;
820: if (vtx) *vtx = network->cloneshared->subnet[netnum].vertices;
821: if (edge) *edge = network->cloneshared->subnet[netnum].edges;
822: PetscFunctionReturn(PETSC_SUCCESS);
823: }
825: /*@
826: DMNetworkAddSharedVertices - Add shared vertices that connect two given subnetworks
828: Collective on dm
830: Input Parameters:
831: + dm - the dm object
832: . anetnum - first subnetwork global numbering returned by DMNetworkAddSubnetwork()
833: . bnetnum - second subnetwork global numbering returned by DMNetworkAddSubnetwork()
834: . nsvtx - number of vertices that are shared by the two subnetworks
835: . asvtx - vertex index in the first subnetwork
836: - bsvtx - vertex index in the second subnetwork
838: Level: beginner
840: .seealso: `DMNetworkCreate()`, `DMNetworkAddSubnetwork()`, `DMNetworkGetSharedVertices()`
841: @*/
842: PetscErrorCode DMNetworkAddSharedVertices(DM dm, PetscInt anetnum, PetscInt bnetnum, PetscInt nsvtx, PetscInt asvtx[], PetscInt bsvtx[])
843: {
844: DM_Network *network = (DM_Network *)dm->data;
845: PetscInt i, nsubnet = network->cloneshared->Nsubnet, *sedgelist, Nsvtx = network->cloneshared->Nsvtx;
847: PetscFunctionBegin;
848: PetscCheck(anetnum != bnetnum, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Subnetworks must have different netnum");
849: PetscCheck(anetnum >= 0 && bnetnum >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "netnum cannot be negative");
850: if (!Nsvtx) {
851: /* allocate network->sedgelist to hold at most 2*nsubnet pairs of shared vertices */
852: PetscCall(PetscMalloc1(2 * 4 * nsubnet, &network->cloneshared->sedgelist));
853: }
855: sedgelist = network->cloneshared->sedgelist;
856: for (i = 0; i < nsvtx; i++) {
857: sedgelist[4 * Nsvtx] = anetnum;
858: sedgelist[4 * Nsvtx + 1] = asvtx[i];
859: sedgelist[4 * Nsvtx + 2] = bnetnum;
860: sedgelist[4 * Nsvtx + 3] = bsvtx[i];
861: Nsvtx++;
862: }
863: PetscCheck(Nsvtx <= 2 * nsubnet, PETSC_COMM_SELF, PETSC_ERR_SUP, "allocate more space for coupling edgelist");
864: network->cloneshared->Nsvtx = Nsvtx;
865: PetscFunctionReturn(PETSC_SUCCESS);
866: }
868: /*@C
869: DMNetworkGetSharedVertices - Returns the info for the shared vertices
871: Not collective
873: Input Parameter:
874: . dm - the DM object
876: Output Parameters:
877: + nsv - number of local shared vertices
878: - svtx - local shared vertices
880: Notes:
881: Cannot call this routine before DMNetworkLayoutSetup()
883: Level: intermediate
885: .seealso: `DMNetworkGetSubnetwork()`, `DMNetworkLayoutSetUp()`, `DMNetworkAddSharedVertices()`
886: @*/
887: PetscErrorCode DMNetworkGetSharedVertices(DM dm, PetscInt *nsv, const PetscInt **svtx)
888: {
889: DM_Network *net = (DM_Network *)dm->data;
891: PetscFunctionBegin;
893: if (nsv) *nsv = net->cloneshared->nsvtx;
894: if (svtx) *svtx = net->cloneshared->svertices;
895: PetscFunctionReturn(PETSC_SUCCESS);
896: }
898: /*@C
899: DMNetworkRegisterComponent - Registers the network component
901: Logically collective on dm
903: Input Parameters:
904: + dm - the network object
905: . name - the component name
906: - size - the storage size in bytes for this component data
908: Output Parameters:
909: . key - an integer key that defines the component
911: Note:
912: This routine should be called by all processors before calling `DMNetworkLayoutSetup()`.
914: Level: beginner
916: .seealso: `DMNetworkCreate()`, `DMNetworkLayoutSetUp()`
917: @*/
918: PetscErrorCode DMNetworkRegisterComponent(DM dm, const char *name, size_t size, PetscInt *key)
919: {
920: DM_Network *network = (DM_Network *)dm->data;
921: DMNetworkComponent *component = NULL, *newcomponent = NULL;
922: PetscBool flg = PETSC_FALSE;
923: PetscInt i;
925: PetscFunctionBegin;
926: if (!network->component) PetscCall(PetscCalloc1(network->max_comps_registered, &network->component));
928: for (i = 0; i < network->ncomponent; i++) {
929: PetscCall(PetscStrcmp(network->component[i].name, name, &flg));
930: if (flg) {
931: *key = i;
932: PetscFunctionReturn(PETSC_SUCCESS);
933: }
934: }
936: if (network->ncomponent == network->max_comps_registered) {
937: /* Reached max allowed so resize component */
938: network->max_comps_registered += 2;
939: PetscCall(PetscCalloc1(network->max_comps_registered, &newcomponent));
940: /* Copy over the previous component info */
941: for (i = 0; i < network->ncomponent; i++) {
942: PetscCall(PetscStrcpy(newcomponent[i].name, network->component[i].name));
943: newcomponent[i].size = network->component[i].size;
944: }
945: /* Free old one */
946: PetscCall(PetscFree(network->component));
947: /* Update pointer */
948: network->component = newcomponent;
949: }
951: component = &network->component[network->ncomponent];
953: PetscCall(PetscStrcpy(component->name, name));
954: component->size = size / sizeof(DMNetworkComponentGenericDataType);
955: *key = network->ncomponent;
956: network->ncomponent++;
957: PetscFunctionReturn(PETSC_SUCCESS);
958: }
960: /*@
961: DMNetworkGetNumVertices - Get the local and global number of vertices for the entire network.
963: Not Collective
965: Input Parameter:
966: . dm - the DMNetwork object
968: Output Parameters:
969: + nVertices - the local number of vertices
970: - NVertices - the global number of vertices
972: Level: beginner
974: .seealso: `DMNetworkGetNumEdges()`
975: @*/
976: PetscErrorCode DMNetworkGetNumVertices(DM dm, PetscInt *nVertices, PetscInt *NVertices)
977: {
978: DM_Network *network = (DM_Network *)dm->data;
980: PetscFunctionBegin;
982: if (nVertices) {
984: *nVertices = network->cloneshared->nVertices;
985: }
986: if (NVertices) {
988: *NVertices = network->cloneshared->NVertices;
989: }
990: PetscFunctionReturn(PETSC_SUCCESS);
991: }
993: /*@
994: DMNetworkGetNumEdges - Get the local and global number of edges for the entire network.
996: Not Collective
998: Input Parameter:
999: . dm - the DMNetwork object
1001: Output Parameters:
1002: + nEdges - the local number of edges
1003: - NEdges - the global number of edges
1005: Level: beginner
1007: .seealso: `DMNetworkGetNumVertices()`
1008: @*/
1009: PetscErrorCode DMNetworkGetNumEdges(DM dm, PetscInt *nEdges, PetscInt *NEdges)
1010: {
1011: DM_Network *network = (DM_Network *)dm->data;
1013: PetscFunctionBegin;
1015: if (nEdges) {
1017: *nEdges = network->cloneshared->nEdges;
1018: }
1019: if (NEdges) {
1021: *NEdges = network->cloneshared->NEdges;
1022: }
1023: PetscFunctionReturn(PETSC_SUCCESS);
1024: }
1026: /*@
1027: DMNetworkGetVertexRange - Get the bounds [start, end) for the local vertices
1029: Not Collective
1031: Input Parameter:
1032: . dm - the DMNetwork object
1034: Output Parameters:
1035: + vStart - the first vertex point
1036: - vEnd - one beyond the last vertex point
1038: Level: beginner
1040: .seealso: `DMNetworkGetEdgeRange()`
1041: @*/
1042: PetscErrorCode DMNetworkGetVertexRange(DM dm, PetscInt *vStart, PetscInt *vEnd)
1043: {
1044: DM_Network *network = (DM_Network *)dm->data;
1046: PetscFunctionBegin;
1047: if (vStart) *vStart = network->cloneshared->vStart;
1048: if (vEnd) *vEnd = network->cloneshared->vEnd;
1049: PetscFunctionReturn(PETSC_SUCCESS);
1050: }
1052: /*@
1053: DMNetworkGetEdgeRange - Get the bounds [start, end) for the local edges
1055: Not Collective
1057: Input Parameter:
1058: . dm - the DMNetwork object
1060: Output Parameters:
1061: + eStart - The first edge point
1062: - eEnd - One beyond the last edge point
1064: Level: beginner
1066: .seealso: `DMNetworkGetVertexRange()`
1067: @*/
1068: PetscErrorCode DMNetworkGetEdgeRange(DM dm, PetscInt *eStart, PetscInt *eEnd)
1069: {
1070: DM_Network *network = (DM_Network *)dm->data;
1072: PetscFunctionBegin;
1074: if (eStart) *eStart = network->cloneshared->eStart;
1075: if (eEnd) *eEnd = network->cloneshared->eEnd;
1076: PetscFunctionReturn(PETSC_SUCCESS);
1077: }
1079: PetscErrorCode DMNetworkGetIndex(DM dm, PetscInt p, PetscInt *index)
1080: {
1081: DM_Network *network = (DM_Network *)dm->data;
1083: PetscFunctionBegin;
1084: if (network->header) {
1085: *index = network->header[p].index;
1086: } else {
1087: PetscInt offsetp;
1088: DMNetworkComponentHeader header;
1090: PetscCall(PetscSectionGetOffset(network->DataSection, p, &offsetp));
1091: header = (DMNetworkComponentHeader)(network->componentdataarray + offsetp);
1092: *index = header->index;
1093: }
1094: PetscFunctionReturn(PETSC_SUCCESS);
1095: }
1097: PetscErrorCode DMNetworkGetSubnetID(DM dm, PetscInt p, PetscInt *subnetid)
1098: {
1099: DM_Network *network = (DM_Network *)dm->data;
1101: PetscFunctionBegin;
1102: if (network->header) {
1103: *subnetid = network->header[p].subnetid;
1104: } else {
1105: PetscInt offsetp;
1106: DMNetworkComponentHeader header;
1108: PetscCall(PetscSectionGetOffset(network->DataSection, p, &offsetp));
1109: header = (DMNetworkComponentHeader)(network->componentdataarray + offsetp);
1110: *subnetid = header->subnetid;
1111: }
1112: PetscFunctionReturn(PETSC_SUCCESS);
1113: }
1115: /*@
1116: DMNetworkGetGlobalEdgeIndex - Get the global numbering for the edge on the network
1118: Not Collective
1120: Input Parameters:
1121: + dm - DMNetwork object
1122: - p - edge point
1124: Output Parameters:
1125: . index - the global numbering for the edge
1127: Level: intermediate
1129: .seealso: `DMNetworkGetGlobalVertexIndex()`
1130: @*/
1131: PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm, PetscInt p, PetscInt *index)
1132: {
1133: PetscFunctionBegin;
1134: PetscCall(DMNetworkGetIndex(dm, p, index));
1135: PetscFunctionReturn(PETSC_SUCCESS);
1136: }
1138: /*@
1139: DMNetworkGetGlobalVertexIndex - Get the global numbering for the vertex on the network
1141: Not Collective
1143: Input Parameters:
1144: + dm - DMNetwork object
1145: - p - vertex point
1147: Output Parameters:
1148: . index - the global numbering for the vertex
1150: Level: intermediate
1152: .seealso: `DMNetworkGetGlobalEdgeIndex()`, `DMNetworkGetLocalVertexIndex()`
1153: @*/
1154: PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm, PetscInt p, PetscInt *index)
1155: {
1156: PetscFunctionBegin;
1157: PetscCall(DMNetworkGetIndex(dm, p, index));
1158: PetscFunctionReturn(PETSC_SUCCESS);
1159: }
1161: /*@
1162: DMNetworkGetNumComponents - Get the number of components at a vertex/edge
1164: Not Collective
1166: Input Parameters:
1167: + dm - the DMNetwork object
1168: - p - vertex/edge point
1170: Output Parameters:
1171: . numcomponents - Number of components at the vertex/edge
1173: Level: beginner
1175: .seealso: `DMNetworkRegisterComponent()`, `DMNetworkAddComponent()`
1176: @*/
1177: PetscErrorCode DMNetworkGetNumComponents(DM dm, PetscInt p, PetscInt *numcomponents)
1178: {
1179: PetscInt offset;
1180: DM_Network *network = (DM_Network *)dm->data;
1182: PetscFunctionBegin;
1183: PetscCall(PetscSectionGetOffset(network->DataSection, p, &offset));
1184: *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray + offset))->ndata;
1185: PetscFunctionReturn(PETSC_SUCCESS);
1186: }
1188: /*@
1189: DMNetworkGetLocalVecOffset - Get the offset for accessing the variables associated with a component at the given vertex/edge from the local vector
1191: Not Collective
1193: Input Parameters:
1194: + dm - the DMNetwork object
1195: . p - the edge or vertex point
1196: - compnum - component number; use ALL_COMPONENTS if no specific component is requested
1198: Output Parameters:
1199: . offset - the local offset
1201: Notes:
1202: These offsets can be passed to MatSetValuesLocal() for matrices obtained with DMCreateMatrix().
1204: For vectors obtained with DMCreateLocalVector() the offsets can be used with VecSetValues().
1206: For vectors obtained with DMCreateLocalVector() and the array obtained with VecGetArray(vec,&array) you can access or set
1207: the vector values with array[offset].
1209: For vectors obtained with DMCreateGlobalVector() the offsets can be used with VecSetValuesLocal().
1211: Level: intermediate
1213: .seealso: `DMGetLocalVector()`, `DMNetworkGetComponent()`, `DMNetworkGetGlobalVecOffset()`, `DMCreateGlobalVector()`, `VecGetArray()`, `VecSetValuesLocal()`, `MatSetValuesLocal()`
1214: @*/
1215: PetscErrorCode DMNetworkGetLocalVecOffset(DM dm, PetscInt p, PetscInt compnum, PetscInt *offset)
1216: {
1217: DM_Network *network = (DM_Network *)dm->data;
1218: PetscInt offsetp, offsetd;
1219: DMNetworkComponentHeader header;
1221: PetscFunctionBegin;
1222: PetscCall(PetscSectionGetOffset(network->plex->localSection, p, &offsetp));
1223: if (compnum == ALL_COMPONENTS) {
1224: *offset = offsetp;
1225: PetscFunctionReturn(PETSC_SUCCESS);
1226: }
1228: PetscCall(PetscSectionGetOffset(network->DataSection, p, &offsetd));
1229: header = (DMNetworkComponentHeader)(network->componentdataarray + offsetd);
1230: *offset = offsetp + header->offsetvarrel[compnum];
1231: PetscFunctionReturn(PETSC_SUCCESS);
1232: }
1234: /*@
1235: DMNetworkGetGlobalVecOffset - Get the global offset for accessing the variables associated with a component for the given vertex/edge from the global vector
1237: Not Collective
1239: Input Parameters:
1240: + dm - the DMNetwork object
1241: . p - the edge or vertex point
1242: - compnum - component number; use ALL_COMPONENTS if no specific component is requested
1244: Output Parameters:
1245: . offsetg - the global offset
1247: Notes:
1248: These offsets can be passed to MatSetValues() for matrices obtained with DMCreateMatrix().
1250: For vectors obtained with DMCreateGlobalVector() the offsets can be used with VecSetValues().
1252: For vectors obtained with DMCreateGlobalVector() and the array obtained with VecGetArray(vec,&array) you can access or set
1253: the vector values with array[offset - rstart] where restart is obtained with VecGetOwnershipRange(v,&rstart,NULL);
1255: Level: intermediate
1257: .seealso: `DMNetworkGetLocalVecOffset()`, `DMGetGlobalVector()`, `DMNetworkGetComponent()`, `DMCreateGlobalVector()`, `VecGetArray()`, `VecSetValues()`, `MatSetValues()`
1258: @*/
1259: PetscErrorCode DMNetworkGetGlobalVecOffset(DM dm, PetscInt p, PetscInt compnum, PetscInt *offsetg)
1260: {
1261: DM_Network *network = (DM_Network *)dm->data;
1262: PetscInt offsetp, offsetd;
1263: DMNetworkComponentHeader header;
1265: PetscFunctionBegin;
1266: PetscCall(PetscSectionGetOffset(network->plex->globalSection, p, &offsetp));
1267: if (offsetp < 0) offsetp = -(offsetp + 1); /* Convert to actual global offset for ghost vertex */
1269: if (compnum == ALL_COMPONENTS) {
1270: *offsetg = offsetp;
1271: PetscFunctionReturn(PETSC_SUCCESS);
1272: }
1273: PetscCall(PetscSectionGetOffset(network->DataSection, p, &offsetd));
1274: header = (DMNetworkComponentHeader)(network->componentdataarray + offsetd);
1275: *offsetg = offsetp + header->offsetvarrel[compnum];
1276: PetscFunctionReturn(PETSC_SUCCESS);
1277: }
1279: /*@
1280: DMNetworkGetEdgeOffset - Get the offset for accessing the variables associated with the given edge from the local subvector
1282: Not Collective
1284: Input Parameters:
1285: + dm - the DMNetwork object
1286: - p - the edge point
1288: Output Parameters:
1289: . offset - the offset
1291: Level: intermediate
1293: .seealso: `DMNetworkGetLocalVecOffset()`, `DMGetLocalVector()`
1294: @*/
1295: PetscErrorCode DMNetworkGetEdgeOffset(DM dm, PetscInt p, PetscInt *offset)
1296: {
1297: DM_Network *network = (DM_Network *)dm->data;
1299: PetscFunctionBegin;
1300: PetscCall(PetscSectionGetOffset(network->edge.DofSection, p, offset));
1301: PetscFunctionReturn(PETSC_SUCCESS);
1302: }
1304: /*@
1305: DMNetworkGetVertexOffset - Get the offset for accessing the variables associated with the given vertex from the local subvector
1307: Not Collective
1309: Input Parameters:
1310: + dm - the DMNetwork object
1311: - p - the vertex point
1313: Output Parameters:
1314: . offset - the offset
1316: Level: intermediate
1318: .seealso: `DMNetworkGetEdgeOffset()`, `DMGetLocalVector()`
1319: @*/
1320: PetscErrorCode DMNetworkGetVertexOffset(DM dm, PetscInt p, PetscInt *offset)
1321: {
1322: DM_Network *network = (DM_Network *)dm->data;
1324: PetscFunctionBegin;
1325: p -= network->cloneshared->vStart;
1326: PetscCall(PetscSectionGetOffset(network->vertex.DofSection, p, offset));
1327: PetscFunctionReturn(PETSC_SUCCESS);
1328: }
1330: /*@
1331: DMNetworkAddComponent - Adds a network component and number of variables at the given point (vertex/edge)
1333: Collective on dm
1335: Input Parameters:
1336: + dm - the DMNetwork
1337: . p - the vertex/edge point. These points are local indices provided by DMNetworkGetSubnetwork()
1338: . componentkey - component key returned while registering the component with DMNetworkRegisterComponent()
1339: . 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
1340: free this space until after DMSetUp() is called.
1341: - 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
1343: Notes:
1344: 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.
1346: DMNetworkLayoutSetUp() must be called before this routine.
1348: Developer Notes:
1349: 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.
1350: Level: beginner
1352: .seealso: `DMNetworkGetComponent()`, `DMNetworkGetSubnetwork()`, `DMNetworkIsGhostVertex()`, `DMNetworkLayoutSetUp()`
1353: @*/
1354: PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p, PetscInt componentkey, void *compvalue, PetscInt nvar)
1355: {
1356: DM_Network *network = (DM_Network *)dm->data;
1357: DMNetworkComponent *component = &network->component[componentkey];
1358: DMNetworkComponentHeader header;
1359: DMNetworkComponentValue cvalue;
1360: PetscInt compnum;
1361: PetscInt *compsize, *compkey, *compoffset, *compnvar, *compoffsetvarrel;
1362: void **compdata;
1364: PetscFunctionBegin;
1365: 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);
1366: 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.");
1367: /* The owning rank and all ghost ranks add nvar */
1368: PetscCall(PetscSectionAddDof(network->DofSection, p, nvar));
1370: /* The owning rank and all ghost ranks add a component, including compvalue=NULL */
1371: header = &network->header[p];
1372: cvalue = &network->cvalue[p];
1373: if (header->ndata == header->maxcomps) {
1374: PetscInt additional_size;
1376: /* Reached limit so resize header component arrays */
1377: header->maxcomps += 2;
1379: /* Allocate arrays for component information and value */
1380: PetscCall(PetscCalloc5(header->maxcomps, &compsize, header->maxcomps, &compkey, header->maxcomps, &compoffset, header->maxcomps, &compnvar, header->maxcomps, &compoffsetvarrel));
1381: PetscCall(PetscMalloc1(header->maxcomps, &compdata));
1383: /* Recalculate header size */
1384: header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5 * header->maxcomps * sizeof(PetscInt);
1386: header->hsize /= sizeof(DMNetworkComponentGenericDataType);
1388: /* Copy over component info */
1389: PetscCall(PetscMemcpy(compsize, header->size, header->ndata * sizeof(PetscInt)));
1390: PetscCall(PetscMemcpy(compkey, header->key, header->ndata * sizeof(PetscInt)));
1391: PetscCall(PetscMemcpy(compoffset, header->offset, header->ndata * sizeof(PetscInt)));
1392: PetscCall(PetscMemcpy(compnvar, header->nvar, header->ndata * sizeof(PetscInt)));
1393: PetscCall(PetscMemcpy(compoffsetvarrel, header->offsetvarrel, header->ndata * sizeof(PetscInt)));
1395: /* Copy over component data pointers */
1396: PetscCall(PetscMemcpy(compdata, cvalue->data, header->ndata * sizeof(void *)));
1398: /* Free old arrays */
1399: PetscCall(PetscFree5(header->size, header->key, header->offset, header->nvar, header->offsetvarrel));
1400: PetscCall(PetscFree(cvalue->data));
1402: /* Update pointers */
1403: header->size = compsize;
1404: header->key = compkey;
1405: header->offset = compoffset;
1406: header->nvar = compnvar;
1407: header->offsetvarrel = compoffsetvarrel;
1409: cvalue->data = compdata;
1411: /* Update DataSection Dofs */
1412: /* 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 */
1413: additional_size = (5 * (header->maxcomps - header->ndata) * sizeof(PetscInt)) / sizeof(DMNetworkComponentGenericDataType);
1414: PetscCall(PetscSectionAddDof(network->DataSection, p, additional_size));
1415: }
1416: header = &network->header[p];
1417: cvalue = &network->cvalue[p];
1419: compnum = header->ndata;
1421: header->size[compnum] = component->size;
1422: PetscCall(PetscSectionAddDof(network->DataSection, p, component->size));
1423: header->key[compnum] = componentkey;
1424: if (compnum != 0) header->offset[compnum] = header->offset[compnum - 1] + header->size[compnum - 1];
1425: cvalue->data[compnum] = (void *)compvalue;
1427: /* variables */
1428: header->nvar[compnum] += nvar;
1429: if (compnum != 0) header->offsetvarrel[compnum] = header->offsetvarrel[compnum - 1] + header->nvar[compnum - 1];
1431: header->ndata++;
1432: PetscFunctionReturn(PETSC_SUCCESS);
1433: }
1435: /*@
1436: DMNetworkGetComponent - Gets the component key, the component data, and the number of variables at a given network point
1438: Not Collective
1440: Input Parameters:
1441: + dm - the DMNetwork object
1442: . p - vertex/edge point
1443: - compnum - component number; use ALL_COMPONENTS if sum up all the components
1445: Output Parameters:
1446: + compkey - the key obtained when registering the component (use NULL if not required)
1447: . component - the component data (use NULL if not required)
1448: - nvar - number of variables (use NULL if not required)
1450: Level: beginner
1452: .seealso: `DMNetworkAddComponent()`, `DMNetworkGetNumComponents()`
1453: @*/
1454: PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *compkey, void **component, PetscInt *nvar)
1455: {
1456: DM_Network *network = (DM_Network *)dm->data;
1457: PetscInt offset = 0;
1458: DMNetworkComponentHeader header;
1460: PetscFunctionBegin;
1461: if (compnum == ALL_COMPONENTS) {
1462: PetscCall(PetscSectionGetDof(network->DofSection, p, nvar));
1463: PetscFunctionReturn(PETSC_SUCCESS);
1464: }
1466: PetscCall(PetscSectionGetOffset(network->DataSection, p, &offset));
1467: header = (DMNetworkComponentHeader)(network->componentdataarray + offset);
1469: if (compnum >= 0) {
1470: if (compkey) *compkey = header->key[compnum];
1471: if (component) {
1472: offset += header->hsize + header->offset[compnum];
1473: *component = network->componentdataarray + offset;
1474: }
1475: }
1477: if (nvar) *nvar = header->nvar[compnum];
1479: PetscFunctionReturn(PETSC_SUCCESS);
1480: }
1482: /*
1483: Sets up the array that holds the data for all components and its associated section.
1484: 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.
1485: */
1486: PetscErrorCode DMNetworkComponentSetUp(DM dm)
1487: {
1488: DM_Network *network = (DM_Network *)dm->data;
1489: PetscInt arr_size, p, offset, offsetp, ncomp, i, *headerarr;
1490: DMNetworkComponentHeader header;
1491: DMNetworkComponentValue cvalue;
1492: DMNetworkComponentHeader headerinfo;
1493: DMNetworkComponentGenericDataType *componentdataarray;
1495: PetscFunctionBegin;
1496: PetscCall(PetscSectionSetUp(network->DataSection));
1497: PetscCall(PetscSectionGetStorageSize(network->DataSection, &arr_size));
1498: /* arr_size+1 fixes pipeline test of opensolaris-misc for src/dm/tests/ex10.c -- Do not know why */
1499: PetscCall(PetscCalloc1(arr_size + 1, &network->componentdataarray));
1500: componentdataarray = network->componentdataarray;
1501: for (p = network->cloneshared->pStart; p < network->cloneshared->pEnd; p++) {
1502: PetscCall(PetscSectionGetOffset(network->DataSection, p, &offsetp));
1503: /* Copy header */
1504: header = &network->header[p];
1505: headerinfo = (DMNetworkComponentHeader)(componentdataarray + offsetp);
1506: PetscCall(PetscMemcpy(headerinfo, header, sizeof(struct _p_DMNetworkComponentHeader)));
1507: headerarr = (PetscInt *)(headerinfo + 1);
1508: PetscCall(PetscMemcpy(headerarr, header->size, header->maxcomps * sizeof(PetscInt)));
1509: headerinfo->size = headerarr;
1510: headerarr += header->maxcomps;
1511: PetscCall(PetscMemcpy(headerarr, header->key, header->maxcomps * sizeof(PetscInt)));
1512: headerinfo->key = headerarr;
1513: headerarr += header->maxcomps;
1514: PetscCall(PetscMemcpy(headerarr, header->offset, header->maxcomps * sizeof(PetscInt)));
1515: headerinfo->offset = headerarr;
1516: headerarr += header->maxcomps;
1517: PetscCall(PetscMemcpy(headerarr, header->nvar, header->maxcomps * sizeof(PetscInt)));
1518: headerinfo->nvar = headerarr;
1519: headerarr += header->maxcomps;
1520: PetscCall(PetscMemcpy(headerarr, header->offsetvarrel, header->maxcomps * sizeof(PetscInt)));
1521: headerinfo->offsetvarrel = headerarr;
1523: /* Copy data */
1524: cvalue = &network->cvalue[p];
1525: ncomp = header->ndata;
1527: for (i = 0; i < ncomp; i++) {
1528: offset = offsetp + header->hsize + header->offset[i];
1529: PetscCall(PetscMemcpy(componentdataarray + offset, cvalue->data[i], header->size[i] * sizeof(DMNetworkComponentGenericDataType)));
1530: }
1531: }
1533: for (i = network->cloneshared->pStart; i < network->cloneshared->pEnd; i++) {
1534: PetscCall(PetscFree5(network->header[i].size, network->header[i].key, network->header[i].offset, network->header[i].nvar, network->header[i].offsetvarrel));
1535: PetscCall(PetscFree(network->cvalue[i].data));
1536: }
1537: PetscCall(PetscFree2(network->header, network->cvalue));
1538: PetscFunctionReturn(PETSC_SUCCESS);
1539: }
1541: /* Sets up the section for dofs. This routine is called during DMSetUp() */
1542: static PetscErrorCode DMNetworkVariablesSetUp(DM dm)
1543: {
1544: DM_Network *network = (DM_Network *)dm->data;
1546: PetscFunctionBegin;
1547: PetscCall(PetscSectionSetUp(network->DofSection));
1548: PetscFunctionReturn(PETSC_SUCCESS);
1549: }
1551: /* Get a subsection from a range of points */
1552: static PetscErrorCode DMNetworkGetSubSection_private(PetscSection main, PetscInt pstart, PetscInt pend, PetscSection *subsection)
1553: {
1554: PetscInt i, nvar;
1556: PetscFunctionBegin;
1557: PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)main), subsection));
1558: PetscCall(PetscSectionSetChart(*subsection, 0, pend - pstart));
1559: for (i = pstart; i < pend; i++) {
1560: PetscCall(PetscSectionGetDof(main, i, &nvar));
1561: PetscCall(PetscSectionSetDof(*subsection, i - pstart, nvar));
1562: }
1564: PetscCall(PetscSectionSetUp(*subsection));
1565: PetscFunctionReturn(PETSC_SUCCESS);
1566: }
1568: /* Create a submap of points with a GlobalToLocal structure */
1569: static PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
1570: {
1571: PetscInt i, *subpoints;
1573: PetscFunctionBegin;
1574: /* Create index sets to map from "points" to "subpoints" */
1575: PetscCall(PetscMalloc1(pend - pstart, &subpoints));
1576: for (i = pstart; i < pend; i++) subpoints[i - pstart] = i;
1577: PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, pend - pstart, subpoints, PETSC_COPY_VALUES, map));
1578: PetscCall(PetscFree(subpoints));
1579: PetscFunctionReturn(PETSC_SUCCESS);
1580: }
1582: /*@
1583: DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute
1585: Collective on dm
1587: Input Parameters:
1588: . dm - the DMNetworkObject
1590: Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
1592: points = [0 1 2 3 4 5 6]
1594: 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]).
1596: With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
1598: Level: intermediate
1600: @*/
1601: PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
1602: {
1603: MPI_Comm comm;
1604: PetscMPIInt size;
1605: DM_Network *network = (DM_Network *)dm->data;
1607: PetscFunctionBegin;
1608: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1609: PetscCallMPI(MPI_Comm_size(comm, &size));
1611: /* Create maps for vertices and edges */
1612: PetscCall(DMNetworkSetSubMap_private(network->cloneshared->vStart, network->cloneshared->vEnd, &network->vertex.mapping));
1613: PetscCall(DMNetworkSetSubMap_private(network->cloneshared->eStart, network->cloneshared->eEnd, &network->edge.mapping));
1615: /* Create local sub-sections */
1616: PetscCall(DMNetworkGetSubSection_private(network->DofSection, network->cloneshared->vStart, network->cloneshared->vEnd, &network->vertex.DofSection));
1617: PetscCall(DMNetworkGetSubSection_private(network->DofSection, network->cloneshared->eStart, network->cloneshared->eEnd, &network->edge.DofSection));
1619: if (size > 1) {
1620: PetscCall(PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf));
1622: PetscCall(PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection));
1623: PetscCall(PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf));
1624: PetscCall(PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection));
1625: } else {
1626: /* create structures for vertex */
1627: PetscCall(PetscSectionClone(network->vertex.DofSection, &network->vertex.GlobalDofSection));
1628: /* create structures for edge */
1629: PetscCall(PetscSectionClone(network->edge.DofSection, &network->edge.GlobalDofSection));
1630: }
1632: /* Add viewers */
1633: PetscCall(PetscObjectSetName((PetscObject)network->edge.GlobalDofSection, "Global edge dof section"));
1634: PetscCall(PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection, "Global vertex dof section"));
1635: PetscCall(PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view"));
1636: PetscCall(PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view"));
1637: PetscFunctionReturn(PETSC_SUCCESS);
1638: }
1640: /*
1641: Setup a lookup btable for the input v's owning subnetworks
1642: - add all owing subnetworks that connect to this v to the btable
1643: vertex_subnetid = supportingedge_subnetid
1644: */
1645: static inline PetscErrorCode SetSubnetIdLookupBT(DM dm, PetscInt v, PetscInt Nsubnet, PetscBT btable)
1646: {
1647: PetscInt e, nedges, offset;
1648: const PetscInt *edges;
1649: DM_Network *newDMnetwork = (DM_Network *)dm->data;
1650: DMNetworkComponentHeader header;
1652: PetscFunctionBegin;
1653: PetscCall(PetscBTMemzero(Nsubnet, btable));
1654: PetscCall(DMNetworkGetSupportingEdges(dm, v, &nedges, &edges));
1655: for (e = 0; e < nedges; e++) {
1656: PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection, edges[e], &offset));
1657: header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);
1658: PetscCall(PetscBTSet(btable, header->subnetid));
1659: }
1660: PetscFunctionReturn(PETSC_SUCCESS);
1661: }
1663: /*
1664: DMNetworkDistributeCoordinates - Internal function to distribute the coordinate network and coordinates.
1666: Collective
1668: Input Parameters:
1669: + dm - The original DMNetwork object
1670: - migrationSF - The PetscSF discribing the migrgration from dm to dmnew
1671: - newDM - The new distributed dmnetwork object.
1672: */
1674: static PetscErrorCode DMNetworkDistributeCoordinates(DM dm, PetscSF migrationSF, DM newDM)
1675: {
1676: DM_Network *newDMnetwork = (DM_Network *)((newDM)->data), *newCoordnetwork, *oldCoordnetwork;
1677: DM cdm, newcdm;
1678: PetscInt cdim, bs, p, pStart, pEnd, offset;
1679: Vec oldCoord, newCoord;
1680: DMNetworkComponentHeader header;
1681: const char *name;
1683: PetscFunctionBegin;
1684: /* Distribute the coordinate network and coordinates */
1685: PetscCall(DMGetCoordinateDim(dm, &cdim));
1686: PetscCall(DMSetCoordinateDim(newDM, cdim));
1688: /* Migrate only if original network had coordinates */
1689: PetscCall(DMGetCoordinatesLocal(dm, &oldCoord));
1690: if (oldCoord) {
1691: PetscCall(DMGetCoordinateDM(dm, &cdm));
1692: PetscCall(DMGetCoordinateDM(newDM, &newcdm));
1693: newCoordnetwork = (DM_Network *)newcdm->data;
1694: oldCoordnetwork = (DM_Network *)cdm->data;
1696: PetscCall(VecCreate(PETSC_COMM_SELF, &newCoord));
1697: PetscCall(PetscObjectGetName((PetscObject)oldCoord, &name));
1698: PetscCall(PetscObjectSetName((PetscObject)newCoord, name));
1699: PetscCall(VecGetBlockSize(oldCoord, &bs));
1700: PetscCall(VecSetBlockSize(newCoord, bs));
1702: PetscCall(DMPlexDistributeField(newDMnetwork->plex, migrationSF, oldCoordnetwork->DofSection, oldCoord, newCoordnetwork->DofSection, newCoord));
1703: PetscCall(DMSetCoordinatesLocal(newDM, newCoord));
1705: PetscCall(VecDestroy(&newCoord));
1706: /* Migrate the components from the orignal coordinate network to the new coordinate network */
1707: PetscCall(DMPlexDistributeData(newDMnetwork->plex, migrationSF, oldCoordnetwork->DataSection, MPIU_INT, (void *)oldCoordnetwork->componentdataarray, newCoordnetwork->DataSection, (void **)&newCoordnetwork->componentdataarray));
1708: /* update the header pointers in the new coordinate network components */
1709: PetscCall(PetscSectionGetChart(newCoordnetwork->DataSection, &pStart, &pEnd));
1710: for (p = pStart; p < pEnd; p++) {
1711: PetscCall(PetscSectionGetOffset(newCoordnetwork->DataSection, p, &offset));
1712: header = (DMNetworkComponentHeader)(newCoordnetwork->componentdataarray + offset);
1713: /* Update pointers */
1714: header->size = (PetscInt *)(header + 1);
1715: header->key = header->size + header->maxcomps;
1716: header->offset = header->key + header->maxcomps;
1717: header->nvar = header->offset + header->maxcomps;
1718: header->offsetvarrel = header->nvar + header->maxcomps;
1719: }
1721: PetscCall(DMSetLocalSection(newCoordnetwork->plex, newCoordnetwork->DofSection));
1722: PetscCall(DMGetGlobalSection(newCoordnetwork->plex, &newCoordnetwork->GlobalDofSection));
1723: newCoordnetwork->componentsetup = PETSC_TRUE;
1724: }
1725: PetscFunctionReturn(PETSC_SUCCESS);
1726: }
1728: /*@
1729: DMNetworkDistribute - Distributes the network and moves associated component data
1731: Collective
1733: Input Parameters:
1734: + DM - the DMNetwork object
1735: - overlap - the overlap of partitions, 0 is the default
1737: Options Database Key:
1738: + -dmnetwork_view - Calls DMView() at the conclusion of DMSetUp()
1739: - -dmnetwork_view_distributed - Calls DMView() at the conclusion of DMNetworkDistribute()
1741: Notes:
1742: Distributes the network with <overlap>-overlapping partitioning of the edges.
1744: Level: intermediate
1746: .seealso: `DMNetworkCreate()`
1747: @*/
1748: PetscErrorCode DMNetworkDistribute(DM *dm, PetscInt overlap)
1749: {
1750: MPI_Comm comm;
1751: PetscMPIInt size;
1752: DM_Network *oldDMnetwork = (DM_Network *)((*dm)->data), *newDMnetwork;
1753: PetscSF pointsf = NULL;
1754: DM newDM;
1755: PetscInt j, e, v, offset, *subnetvtx, *subnetedge, Nsubnet, gidx, svtx_idx, nv, net;
1756: PetscInt *sv;
1757: PetscBT btable;
1758: PetscPartitioner part;
1759: DMNetworkComponentHeader header;
1761: PetscFunctionBegin;
1764: PetscCall(PetscObjectGetComm((PetscObject)*dm, &comm));
1765: PetscCallMPI(MPI_Comm_size(comm, &size));
1766: if (size == 1) {
1767: oldDMnetwork->cloneshared->distributecalled = PETSC_TRUE;
1768: PetscFunctionReturn(PETSC_SUCCESS);
1769: }
1771: PetscCheck(!overlap, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "overlap %" PetscInt_FMT " != 0 is not supported yet", overlap);
1773: /* 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. */
1774: PetscCall(PetscLogEventBegin(DMNetwork_Distribute, dm, 0, 0, 0));
1775: PetscCall(DMNetworkCreate(PetscObjectComm((PetscObject)*dm), &newDM));
1776: newDMnetwork = (DM_Network *)newDM->data;
1777: newDMnetwork->max_comps_registered = oldDMnetwork->max_comps_registered;
1778: PetscCall(PetscMalloc1(newDMnetwork->max_comps_registered, &newDMnetwork->component));
1780: /* Enable runtime options for petscpartitioner */
1781: PetscCall(DMPlexGetPartitioner(oldDMnetwork->plex, &part));
1782: PetscCall(PetscPartitionerSetFromOptions(part));
1784: /* Distribute plex dm */
1785: PetscCall(DMPlexDistribute(oldDMnetwork->plex, overlap, &pointsf, &newDMnetwork->plex));
1787: /* Distribute dof section */
1788: PetscCall(PetscSectionCreate(comm, &newDMnetwork->DofSection));
1789: PetscCall(PetscSFDistributeSection(pointsf, oldDMnetwork->DofSection, NULL, newDMnetwork->DofSection));
1791: /* Distribute data and associated section */
1792: PetscCall(PetscSectionCreate(comm, &newDMnetwork->DataSection));
1793: PetscCall(DMPlexDistributeData(newDMnetwork->plex, pointsf, oldDMnetwork->DataSection, MPIU_INT, (void *)oldDMnetwork->componentdataarray, newDMnetwork->DataSection, (void **)&newDMnetwork->componentdataarray));
1795: PetscCall(PetscSectionGetChart(newDMnetwork->DataSection, &newDMnetwork->cloneshared->pStart, &newDMnetwork->cloneshared->pEnd));
1796: PetscCall(DMPlexGetHeightStratum(newDMnetwork->plex, 0, &newDMnetwork->cloneshared->eStart, &newDMnetwork->cloneshared->eEnd));
1797: PetscCall(DMPlexGetHeightStratum(newDMnetwork->plex, 1, &newDMnetwork->cloneshared->vStart, &newDMnetwork->cloneshared->vEnd));
1798: newDMnetwork->cloneshared->nEdges = newDMnetwork->cloneshared->eEnd - newDMnetwork->cloneshared->eStart;
1799: newDMnetwork->cloneshared->nVertices = newDMnetwork->cloneshared->vEnd - newDMnetwork->cloneshared->vStart;
1800: newDMnetwork->cloneshared->NVertices = oldDMnetwork->cloneshared->NVertices;
1801: newDMnetwork->cloneshared->NEdges = oldDMnetwork->cloneshared->NEdges;
1802: newDMnetwork->cloneshared->svtable = oldDMnetwork->cloneshared->svtable; /* global table! */
1803: oldDMnetwork->cloneshared->svtable = NULL;
1805: /* Set Dof section as the section for dm */
1806: PetscCall(DMSetLocalSection(newDMnetwork->plex, newDMnetwork->DofSection));
1807: PetscCall(DMGetGlobalSection(newDMnetwork->plex, &newDMnetwork->GlobalDofSection));
1809: /* Setup subnetwork info in the newDM */
1810: newDMnetwork->cloneshared->Nsubnet = oldDMnetwork->cloneshared->Nsubnet;
1811: newDMnetwork->cloneshared->Nsvtx = oldDMnetwork->cloneshared->Nsvtx;
1812: oldDMnetwork->cloneshared->Nsvtx = 0;
1813: newDMnetwork->cloneshared->svtx = oldDMnetwork->cloneshared->svtx; /* global vertices! */
1814: oldDMnetwork->cloneshared->svtx = NULL;
1815: PetscCall(PetscCalloc1(newDMnetwork->cloneshared->Nsubnet, &newDMnetwork->cloneshared->subnet));
1817: /* Copy over the global number of vertices and edges in each subnetwork.
1818: Note: these are calculated in DMNetworkLayoutSetUp()
1819: */
1820: Nsubnet = newDMnetwork->cloneshared->Nsubnet;
1821: for (j = 0; j < Nsubnet; j++) {
1822: newDMnetwork->cloneshared->subnet[j].Nvtx = oldDMnetwork->cloneshared->subnet[j].Nvtx;
1823: newDMnetwork->cloneshared->subnet[j].Nedge = oldDMnetwork->cloneshared->subnet[j].Nedge;
1824: }
1826: /* Count local nedges for subnetworks */
1827: for (e = newDMnetwork->cloneshared->eStart; e < newDMnetwork->cloneshared->eEnd; e++) {
1828: PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection, e, &offset));
1829: header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);
1831: /* Update pointers */
1832: header->size = (PetscInt *)(header + 1);
1833: header->key = header->size + header->maxcomps;
1834: header->offset = header->key + header->maxcomps;
1835: header->nvar = header->offset + header->maxcomps;
1836: header->offsetvarrel = header->nvar + header->maxcomps;
1838: newDMnetwork->cloneshared->subnet[header->subnetid].nedge++;
1839: }
1841: /* Setup a btable to keep track subnetworks owned by this process at a shared vertex */
1842: if (newDMnetwork->cloneshared->Nsvtx) PetscCall(PetscBTCreate(Nsubnet, &btable));
1844: /* Count local nvtx for subnetworks */
1845: for (v = newDMnetwork->cloneshared->vStart; v < newDMnetwork->cloneshared->vEnd; v++) {
1846: PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection, v, &offset));
1847: header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);
1849: /* Update pointers */
1850: header->size = (PetscInt *)(header + 1);
1851: header->key = header->size + header->maxcomps;
1852: header->offset = header->key + header->maxcomps;
1853: header->nvar = header->offset + header->maxcomps;
1854: header->offsetvarrel = header->nvar + header->maxcomps;
1856: /* shared vertices: use gidx=header->index to check if v is a shared vertex */
1857: gidx = header->index;
1858: PetscCall(PetscHMapIGetWithDefault(newDMnetwork->cloneshared->svtable, gidx + 1, 0, &svtx_idx));
1859: svtx_idx--;
1861: if (svtx_idx < 0) { /* not a shared vertex */
1862: newDMnetwork->cloneshared->subnet[header->subnetid].nvtx++;
1863: } else { /* a shared vertex belongs to more than one subnetworks, it is being counted by multiple subnets */
1864: /* Setup a lookup btable for this v's owning subnetworks */
1865: PetscCall(SetSubnetIdLookupBT(newDM, v, Nsubnet, btable));
1867: for (j = 0; j < newDMnetwork->cloneshared->svtx[svtx_idx].n; j++) {
1868: sv = newDMnetwork->cloneshared->svtx[svtx_idx].sv + 2 * j;
1869: net = sv[0];
1870: if (PetscBTLookup(btable, net)) newDMnetwork->cloneshared->subnet[net].nvtx++; /* sv is on net owned by this process */
1871: }
1872: }
1873: }
1875: /* Get total local nvtx for subnetworks */
1876: nv = 0;
1877: for (j = 0; j < Nsubnet; j++) nv += newDMnetwork->cloneshared->subnet[j].nvtx;
1878: nv += newDMnetwork->cloneshared->Nsvtx;
1880: /* Now create the vertices and edge arrays for the subnetworks */
1881: PetscCall(PetscCalloc2(newDMnetwork->cloneshared->nEdges, &subnetedge, nv, &subnetvtx)); /* Maps local vertex to local subnetwork's vertex */
1882: newDMnetwork->cloneshared->subnetedge = subnetedge;
1883: newDMnetwork->cloneshared->subnetvtx = subnetvtx;
1884: for (j = 0; j < newDMnetwork->cloneshared->Nsubnet; j++) {
1885: newDMnetwork->cloneshared->subnet[j].edges = subnetedge;
1886: subnetedge += newDMnetwork->cloneshared->subnet[j].nedge;
1888: newDMnetwork->cloneshared->subnet[j].vertices = subnetvtx;
1889: subnetvtx += newDMnetwork->cloneshared->subnet[j].nvtx;
1891: /* 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. */
1892: newDMnetwork->cloneshared->subnet[j].nvtx = newDMnetwork->cloneshared->subnet[j].nedge = 0;
1893: }
1894: newDMnetwork->cloneshared->svertices = subnetvtx;
1896: /* Set the edges and vertices in each subnetwork */
1897: for (e = newDMnetwork->cloneshared->eStart; e < newDMnetwork->cloneshared->eEnd; e++) {
1898: PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection, e, &offset));
1899: header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);
1900: newDMnetwork->cloneshared->subnet[header->subnetid].edges[newDMnetwork->cloneshared->subnet[header->subnetid].nedge++] = e;
1901: }
1903: nv = 0;
1904: for (v = newDMnetwork->cloneshared->vStart; v < newDMnetwork->cloneshared->vEnd; v++) {
1905: PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection, v, &offset));
1906: header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray + offset);
1908: /* coupling vertices: use gidx = header->index to check if v is a coupling vertex */
1909: PetscCall(PetscHMapIGetWithDefault(newDMnetwork->cloneshared->svtable, header->index + 1, 0, &svtx_idx));
1910: svtx_idx--;
1911: if (svtx_idx < 0) {
1912: newDMnetwork->cloneshared->subnet[header->subnetid].vertices[newDMnetwork->cloneshared->subnet[header->subnetid].nvtx++] = v;
1913: } else { /* a shared vertex */
1914: newDMnetwork->cloneshared->svertices[nv++] = v;
1916: /* Setup a lookup btable for this v's owning subnetworks */
1917: PetscCall(SetSubnetIdLookupBT(newDM, v, Nsubnet, btable));
1919: for (j = 0; j < newDMnetwork->cloneshared->svtx[svtx_idx].n; j++) {
1920: sv = newDMnetwork->cloneshared->svtx[svtx_idx].sv + 2 * j;
1921: net = sv[0];
1922: if (PetscBTLookup(btable, net)) newDMnetwork->cloneshared->subnet[net].vertices[newDMnetwork->cloneshared->subnet[net].nvtx++] = v;
1923: }
1924: }
1925: }
1926: newDMnetwork->cloneshared->nsvtx = nv; /* num of local shared vertices */
1928: PetscCall(DMNetworkDistributeCoordinates(*dm, pointsf, newDM));
1929: newDM->setupcalled = (*dm)->setupcalled;
1930: newDMnetwork->cloneshared->distributecalled = PETSC_TRUE;
1932: /* Free spaces */
1933: PetscCall(PetscSFDestroy(&pointsf));
1934: PetscCall(DMDestroy(dm));
1935: if (newDMnetwork->cloneshared->Nsvtx) PetscCall(PetscBTDestroy(&btable));
1937: /* View distributed dmnetwork */
1938: PetscCall(DMViewFromOptions(newDM, NULL, "-dmnetwork_view_distributed"));
1940: *dm = newDM;
1941: PetscCall(PetscLogEventEnd(DMNetwork_Distribute, dm, 0, 0, 0));
1942: PetscFunctionReturn(PETSC_SUCCESS);
1943: }
1945: /*@C
1946: PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering
1948: Collective
1950: Input Parameters:
1951: + mainSF - the original SF structure
1952: - map - a ISLocalToGlobal mapping that contains the subset of points
1954: Output Parameter:
1955: . subSF - a subset of the mainSF for the desired subset.
1957: Level: intermediate
1958: @*/
1959: PetscErrorCode PetscSFGetSubSF(PetscSF mainsf, ISLocalToGlobalMapping map, PetscSF *subSF)
1960: {
1961: PetscInt nroots, nleaves, *ilocal_sub;
1962: PetscInt i, *ilocal_map, nroots_sub, nleaves_sub = 0;
1963: PetscInt *local_points, *remote_points;
1964: PetscSFNode *iremote_sub;
1965: const PetscInt *ilocal;
1966: const PetscSFNode *iremote;
1968: PetscFunctionBegin;
1969: PetscCall(PetscSFGetGraph(mainsf, &nroots, &nleaves, &ilocal, &iremote));
1971: /* Look for leaves that pertain to the subset of points. Get the local ordering */
1972: PetscCall(PetscMalloc1(nleaves, &ilocal_map));
1973: PetscCall(ISGlobalToLocalMappingApply(map, IS_GTOLM_MASK, nleaves, ilocal, NULL, ilocal_map));
1974: for (i = 0; i < nleaves; i++) {
1975: if (ilocal_map[i] != -1) nleaves_sub += 1;
1976: }
1977: /* Re-number ilocal with subset numbering. Need information from roots */
1978: PetscCall(PetscMalloc2(nroots, &local_points, nroots, &remote_points));
1979: for (i = 0; i < nroots; i++) local_points[i] = i;
1980: PetscCall(ISGlobalToLocalMappingApply(map, IS_GTOLM_MASK, nroots, local_points, NULL, local_points));
1981: PetscCall(PetscSFBcastBegin(mainsf, MPIU_INT, local_points, remote_points, MPI_REPLACE));
1982: PetscCall(PetscSFBcastEnd(mainsf, MPIU_INT, local_points, remote_points, MPI_REPLACE));
1983: /* Fill up graph using local (that is, local to the subset) numbering. */
1984: PetscCall(PetscMalloc1(nleaves_sub, &ilocal_sub));
1985: PetscCall(PetscMalloc1(nleaves_sub, &iremote_sub));
1986: nleaves_sub = 0;
1987: for (i = 0; i < nleaves; i++) {
1988: if (ilocal_map[i] != -1) {
1989: ilocal_sub[nleaves_sub] = ilocal_map[i];
1990: iremote_sub[nleaves_sub].rank = iremote[i].rank;
1991: iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
1992: nleaves_sub += 1;
1993: }
1994: }
1995: PetscCall(PetscFree2(local_points, remote_points));
1996: PetscCall(ISLocalToGlobalMappingGetSize(map, &nroots_sub));
1998: /* Create new subSF */
1999: PetscCall(PetscSFCreate(PETSC_COMM_WORLD, subSF));
2000: PetscCall(PetscSFSetFromOptions(*subSF));
2001: PetscCall(PetscSFSetGraph(*subSF, nroots_sub, nleaves_sub, ilocal_sub, PETSC_OWN_POINTER, iremote_sub, PETSC_COPY_VALUES));
2002: PetscCall(PetscFree(ilocal_map));
2003: PetscCall(PetscFree(iremote_sub));
2004: PetscFunctionReturn(PETSC_SUCCESS);
2005: }
2007: /*@C
2008: DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
2010: Not Collective
2012: Input Parameters:
2013: + dm - the DMNetwork object
2014: - p - the vertex point
2016: Output Parameters:
2017: + nedges - number of edges connected to this vertex point
2018: - edges - list of edge points
2020: Level: beginner
2022: Fortran Notes:
2023: Since it returns an array, this routine is only available in Fortran 90, and you must
2024: include petsc.h90 in your code.
2026: .seealso: `DMNetworkCreate()`, `DMNetworkGetConnectedVertices()`
2027: @*/
2028: PetscErrorCode DMNetworkGetSupportingEdges(DM dm, PetscInt vertex, PetscInt *nedges, const PetscInt *edges[])
2029: {
2030: DM_Network *network = (DM_Network *)dm->data;
2032: PetscFunctionBegin;
2033: PetscCall(DMPlexGetSupportSize(network->plex, vertex, nedges));
2034: if (edges) PetscCall(DMPlexGetSupport(network->plex, vertex, edges));
2035: PetscFunctionReturn(PETSC_SUCCESS);
2036: }
2038: /*@C
2039: DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
2041: Not Collective
2043: Input Parameters:
2044: + dm - the DMNetwork object
2045: - p - the edge point
2047: Output Parameters:
2048: . vertices - vertices connected to this edge
2050: Level: beginner
2052: Fortran Notes:
2053: Since it returns an array, this routine is only available in Fortran 90, and you must
2054: include petsc.h90 in your code.
2056: .seealso: `DMNetworkCreate()`, `DMNetworkGetSupportingEdges()`
2057: @*/
2058: PetscErrorCode DMNetworkGetConnectedVertices(DM dm, PetscInt edge, const PetscInt *vertices[])
2059: {
2060: DM_Network *network = (DM_Network *)dm->data;
2062: PetscFunctionBegin;
2063: PetscCall(DMPlexGetCone(network->plex, edge, vertices));
2064: PetscFunctionReturn(PETSC_SUCCESS);
2065: }
2067: /*@
2068: DMNetworkIsSharedVertex - Returns TRUE if the vertex is shared by subnetworks
2070: Not Collective
2072: Input Parameters:
2073: + dm - the DMNetwork object
2074: - p - the vertex point
2076: Output Parameter:
2077: . flag - TRUE if the vertex is shared by subnetworks
2079: Level: beginner
2081: .seealso: `DMNetworkAddSharedVertices()`, `DMNetworkIsGhostVertex()`
2082: @*/
2083: PetscErrorCode DMNetworkIsSharedVertex(DM dm, PetscInt p, PetscBool *flag)
2084: {
2085: PetscFunctionBegin;
2088: if (dm->setupcalled) { /* DMNetworkGetGlobalVertexIndex() requires DMSetUp() be called */
2089: DM_Network *network = (DM_Network *)dm->data;
2090: PetscInt gidx;
2092: PetscCall(DMNetworkGetGlobalVertexIndex(dm, p, &gidx));
2093: PetscCall(PetscHMapIHas(network->cloneshared->svtable, gidx + 1, flag));
2094: } else { /* would be removed? */
2095: PetscInt nv;
2096: const PetscInt *vtx;
2098: PetscCall(DMNetworkGetSharedVertices(dm, &nv, &vtx));
2099: for (PetscInt i = 0; i < nv; i++) {
2100: if (p == vtx[i]) {
2101: *flag = PETSC_TRUE;
2102: PetscFunctionReturn(PETSC_SUCCESS);
2103: }
2104: }
2105: *flag = PETSC_FALSE;
2106: }
2107: PetscFunctionReturn(PETSC_SUCCESS);
2108: }
2110: /*@
2111: DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
2113: Not Collective
2115: Input Parameters:
2116: + dm - the DMNetwork object
2117: - p - the vertex point
2119: Output Parameter:
2120: . isghost - TRUE if the vertex is a ghost point
2122: Level: beginner
2124: .seealso: `DMNetworkGetConnectedVertices()`, `DMNetworkGetVertexRange()`, `DMNetworkIsSharedVertex()`
2125: @*/
2126: PetscErrorCode DMNetworkIsGhostVertex(DM dm, PetscInt p, PetscBool *isghost)
2127: {
2128: DM_Network *network = (DM_Network *)dm->data;
2129: PetscInt offsetg;
2130: PetscSection sectiong;
2132: PetscFunctionBegin;
2133: *isghost = PETSC_FALSE;
2134: PetscCall(DMGetGlobalSection(network->plex, §iong));
2135: PetscCall(PetscSectionGetOffset(sectiong, p, &offsetg));
2136: if (offsetg < 0) *isghost = PETSC_TRUE;
2137: PetscFunctionReturn(PETSC_SUCCESS);
2138: }
2140: PetscErrorCode DMSetUp_Network(DM dm)
2141: {
2142: PetscFunctionBegin;
2143: PetscCall(PetscLogEventBegin(DMNetwork_SetUpNetwork, dm, 0, 0, 0));
2144: PetscCall(DMNetworkFinalizeComponents(dm));
2145: /* View dmnetwork */
2146: PetscCall(DMViewFromOptions(dm, NULL, "-dmnetwork_view"));
2147: PetscCall(PetscLogEventEnd(DMNetwork_SetUpNetwork, dm, 0, 0, 0));
2148: PetscFunctionReturn(PETSC_SUCCESS);
2149: }
2151: /*@
2152: DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
2153: -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
2155: Collective
2157: Input Parameters:
2158: + dm - the DMNetwork object
2159: . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
2160: - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
2162: Level: intermediate
2164: @*/
2165: PetscErrorCode DMNetworkHasJacobian(DM dm, PetscBool eflg, PetscBool vflg)
2166: {
2167: DM_Network *network = (DM_Network *)dm->data;
2168: PetscInt nVertices = network->cloneshared->nVertices;
2170: PetscFunctionBegin;
2171: network->userEdgeJacobian = eflg;
2172: network->userVertexJacobian = vflg;
2174: if (eflg && !network->Je) PetscCall(PetscCalloc1(3 * network->cloneshared->nEdges, &network->Je));
2176: if (vflg && !network->Jv && nVertices) {
2177: PetscInt i, *vptr, nedges, vStart = network->cloneshared->vStart;
2178: PetscInt nedges_total;
2179: const PetscInt *edges;
2181: /* count nvertex_total */
2182: nedges_total = 0;
2183: PetscCall(PetscMalloc1(nVertices + 1, &vptr));
2185: vptr[0] = 0;
2186: for (i = 0; i < nVertices; i++) {
2187: PetscCall(DMNetworkGetSupportingEdges(dm, i + vStart, &nedges, &edges));
2188: nedges_total += nedges;
2189: vptr[i + 1] = vptr[i] + 2 * nedges + 1;
2190: }
2192: PetscCall(PetscCalloc1(2 * nedges_total + nVertices, &network->Jv));
2193: network->Jvptr = vptr;
2194: }
2195: PetscFunctionReturn(PETSC_SUCCESS);
2196: }
2198: /*@
2199: DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
2201: Not Collective
2203: Input Parameters:
2204: + dm - the DMNetwork object
2205: . p - the edge point
2206: - J - array (size = 3) of Jacobian submatrices for this edge point:
2207: J[0]: this edge
2208: J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()
2210: Level: advanced
2212: .seealso: `DMNetworkVertexSetMatrix()`
2213: @*/
2214: PetscErrorCode DMNetworkEdgeSetMatrix(DM dm, PetscInt p, Mat J[])
2215: {
2216: DM_Network *network = (DM_Network *)dm->data;
2218: PetscFunctionBegin;
2219: PetscCheck(network->Je, PetscObjectComm((PetscObject)dm), PETSC_ERR_ORDER, "Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
2221: if (J) {
2222: network->Je[3 * p] = J[0];
2223: network->Je[3 * p + 1] = J[1];
2224: network->Je[3 * p + 2] = J[2];
2225: }
2226: PetscFunctionReturn(PETSC_SUCCESS);
2227: }
2229: /*@
2230: DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
2232: Not Collective
2234: Input Parameters:
2235: + dm - The DMNetwork object
2236: . p - the vertex point
2237: - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
2238: J[0]: this vertex
2239: J[1+2*i]: i-th supporting edge
2240: J[1+2*i+1]: i-th connected vertex
2242: Level: advanced
2244: .seealso: `DMNetworkEdgeSetMatrix()`
2245: @*/
2246: PetscErrorCode DMNetworkVertexSetMatrix(DM dm, PetscInt p, Mat J[])
2247: {
2248: DM_Network *network = (DM_Network *)dm->data;
2249: PetscInt i, *vptr, nedges, vStart = network->cloneshared->vStart;
2250: const PetscInt *edges;
2252: PetscFunctionBegin;
2253: PetscCheck(network->Jv, PetscObjectComm((PetscObject)dm), PETSC_ERR_ORDER, "Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
2255: if (J) {
2256: vptr = network->Jvptr;
2257: network->Jv[vptr[p - vStart]] = J[0]; /* Set Jacobian for this vertex */
2259: /* Set Jacobian for each supporting edge and connected vertex */
2260: PetscCall(DMNetworkGetSupportingEdges(dm, p, &nedges, &edges));
2261: for (i = 1; i <= 2 * nedges; i++) network->Jv[vptr[p - vStart] + i] = J[i];
2262: }
2263: PetscFunctionReturn(PETSC_SUCCESS);
2264: }
2266: static inline PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscBool ghost, Vec vdnz, Vec vonz)
2267: {
2268: PetscInt j;
2269: PetscScalar val = (PetscScalar)ncols;
2271: PetscFunctionBegin;
2272: if (!ghost) {
2273: for (j = 0; j < nrows; j++) PetscCall(VecSetValues(vdnz, 1, &rows[j], &val, ADD_VALUES));
2274: } else {
2275: for (j = 0; j < nrows; j++) PetscCall(VecSetValues(vonz, 1, &rows[j], &val, ADD_VALUES));
2276: }
2277: PetscFunctionReturn(PETSC_SUCCESS);
2278: }
2280: static inline PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju, PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscBool ghost, Vec vdnz, Vec vonz)
2281: {
2282: PetscInt j, ncols_u;
2283: PetscScalar val;
2285: PetscFunctionBegin;
2286: if (!ghost) {
2287: for (j = 0; j < nrows; j++) {
2288: PetscCall(MatGetRow(Ju, j, &ncols_u, NULL, NULL));
2289: val = (PetscScalar)ncols_u;
2290: PetscCall(VecSetValues(vdnz, 1, &rows[j], &val, ADD_VALUES));
2291: PetscCall(MatRestoreRow(Ju, j, &ncols_u, NULL, NULL));
2292: }
2293: } else {
2294: for (j = 0; j < nrows; j++) {
2295: PetscCall(MatGetRow(Ju, j, &ncols_u, NULL, NULL));
2296: val = (PetscScalar)ncols_u;
2297: PetscCall(VecSetValues(vonz, 1, &rows[j], &val, ADD_VALUES));
2298: PetscCall(MatRestoreRow(Ju, j, &ncols_u, NULL, NULL));
2299: }
2300: }
2301: PetscFunctionReturn(PETSC_SUCCESS);
2302: }
2304: static inline PetscErrorCode MatSetPreallocationblock_private(Mat Ju, PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscBool ghost, Vec vdnz, Vec vonz)
2305: {
2306: PetscFunctionBegin;
2307: if (Ju) {
2308: PetscCall(MatSetPreallocationUserblock_private(Ju, nrows, rows, ncols, ghost, vdnz, vonz));
2309: } else {
2310: PetscCall(MatSetPreallocationDenseblock_private(nrows, rows, ncols, ghost, vdnz, vonz));
2311: }
2312: PetscFunctionReturn(PETSC_SUCCESS);
2313: }
2315: static inline PetscErrorCode MatSetDenseblock_private(PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscInt cstart, Mat *J)
2316: {
2317: PetscInt j, *cols;
2318: PetscScalar *zeros;
2320: PetscFunctionBegin;
2321: PetscCall(PetscCalloc2(ncols, &cols, nrows * ncols, &zeros));
2322: for (j = 0; j < ncols; j++) cols[j] = j + cstart;
2323: PetscCall(MatSetValues(*J, nrows, rows, ncols, cols, zeros, INSERT_VALUES));
2324: PetscCall(PetscFree2(cols, zeros));
2325: PetscFunctionReturn(PETSC_SUCCESS);
2326: }
2328: static inline PetscErrorCode MatSetUserblock_private(Mat Ju, PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscInt cstart, Mat *J)
2329: {
2330: PetscInt j, M, N, row, col, ncols_u;
2331: const PetscInt *cols;
2332: PetscScalar zero = 0.0;
2334: PetscFunctionBegin;
2335: PetscCall(MatGetSize(Ju, &M, &N));
2336: 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);
2338: for (row = 0; row < nrows; row++) {
2339: PetscCall(MatGetRow(Ju, row, &ncols_u, &cols, NULL));
2340: for (j = 0; j < ncols_u; j++) {
2341: col = cols[j] + cstart;
2342: PetscCall(MatSetValues(*J, 1, &rows[row], 1, &col, &zero, INSERT_VALUES));
2343: }
2344: PetscCall(MatRestoreRow(Ju, row, &ncols_u, &cols, NULL));
2345: }
2346: PetscFunctionReturn(PETSC_SUCCESS);
2347: }
2349: static inline PetscErrorCode MatSetblock_private(Mat Ju, PetscInt nrows, PetscInt *rows, PetscInt ncols, PetscInt cstart, Mat *J)
2350: {
2351: PetscFunctionBegin;
2352: if (Ju) {
2353: PetscCall(MatSetUserblock_private(Ju, nrows, rows, ncols, cstart, J));
2354: } else {
2355: PetscCall(MatSetDenseblock_private(nrows, rows, ncols, cstart, J));
2356: }
2357: PetscFunctionReturn(PETSC_SUCCESS);
2358: }
2360: /* 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.
2361: */
2362: PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
2363: {
2364: PetscInt i, size, dof;
2365: PetscInt *glob2loc;
2367: PetscFunctionBegin;
2368: PetscCall(PetscSectionGetStorageSize(localsec, &size));
2369: PetscCall(PetscMalloc1(size, &glob2loc));
2371: for (i = 0; i < size; i++) {
2372: PetscCall(PetscSectionGetOffset(globalsec, i, &dof));
2373: dof = (dof >= 0) ? dof : -(dof + 1);
2374: glob2loc[i] = dof;
2375: }
2377: PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, size, glob2loc, PETSC_OWN_POINTER, ltog));
2378: #if 0
2379: PetscCall(PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD));
2380: #endif
2381: PetscFunctionReturn(PETSC_SUCCESS);
2382: }
2384: #include <petsc/private/matimpl.h>
2386: PetscErrorCode DMCreateMatrix_Network_Nest(DM dm, Mat *J)
2387: {
2388: DM_Network *network = (DM_Network *)dm->data;
2389: PetscInt eDof, vDof;
2390: Mat j11, j12, j21, j22, bA[2][2];
2391: MPI_Comm comm;
2392: ISLocalToGlobalMapping eISMap, vISMap;
2394: PetscFunctionBegin;
2395: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2397: PetscCall(PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection, &eDof));
2398: PetscCall(PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection, &vDof));
2400: PetscCall(MatCreate(comm, &j11));
2401: PetscCall(MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE));
2402: PetscCall(MatSetType(j11, MATMPIAIJ));
2404: PetscCall(MatCreate(comm, &j12));
2405: PetscCall(MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE));
2406: PetscCall(MatSetType(j12, MATMPIAIJ));
2408: PetscCall(MatCreate(comm, &j21));
2409: PetscCall(MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE));
2410: PetscCall(MatSetType(j21, MATMPIAIJ));
2412: PetscCall(MatCreate(comm, &j22));
2413: PetscCall(MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE));
2414: PetscCall(MatSetType(j22, MATMPIAIJ));
2416: bA[0][0] = j11;
2417: bA[0][1] = j12;
2418: bA[1][0] = j21;
2419: bA[1][1] = j22;
2421: PetscCall(CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection, network->edge.DofSection, &eISMap));
2422: PetscCall(CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection, network->vertex.DofSection, &vISMap));
2424: PetscCall(MatSetLocalToGlobalMapping(j11, eISMap, eISMap));
2425: PetscCall(MatSetLocalToGlobalMapping(j12, eISMap, vISMap));
2426: PetscCall(MatSetLocalToGlobalMapping(j21, vISMap, eISMap));
2427: PetscCall(MatSetLocalToGlobalMapping(j22, vISMap, vISMap));
2429: PetscCall(MatSetUp(j11));
2430: PetscCall(MatSetUp(j12));
2431: PetscCall(MatSetUp(j21));
2432: PetscCall(MatSetUp(j22));
2434: PetscCall(MatCreateNest(comm, 2, NULL, 2, NULL, &bA[0][0], J));
2435: PetscCall(MatSetUp(*J));
2436: PetscCall(MatNestSetVecType(*J, VECNEST));
2437: PetscCall(MatDestroy(&j11));
2438: PetscCall(MatDestroy(&j12));
2439: PetscCall(MatDestroy(&j21));
2440: PetscCall(MatDestroy(&j22));
2442: PetscCall(MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY));
2443: PetscCall(MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY));
2444: PetscCall(MatSetOption(*J, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
2446: /* Free structures */
2447: PetscCall(ISLocalToGlobalMappingDestroy(&eISMap));
2448: PetscCall(ISLocalToGlobalMappingDestroy(&vISMap));
2449: PetscFunctionReturn(PETSC_SUCCESS);
2450: }
2452: PetscErrorCode DMCreateMatrix_Network(DM dm, Mat *J)
2453: {
2454: DM_Network *network = (DM_Network *)dm->data;
2455: PetscInt eStart, eEnd, vStart, vEnd, rstart, nrows, *rows, localSize;
2456: PetscInt cstart, ncols, j, e, v;
2457: PetscBool ghost, ghost_vc, ghost2, isNest;
2458: Mat Juser;
2459: PetscSection sectionGlobal;
2460: PetscInt nedges, *vptr = NULL, vc, *rows_v; /* suppress maybe-uninitialized warning */
2461: const PetscInt *edges, *cone;
2462: MPI_Comm comm;
2463: MatType mtype;
2464: Vec vd_nz, vo_nz;
2465: PetscInt *dnnz, *onnz;
2466: PetscScalar *vdnz, *vonz;
2468: PetscFunctionBegin;
2469: mtype = dm->mattype;
2470: PetscCall(PetscStrcmp(mtype, MATNEST, &isNest));
2471: if (isNest) {
2472: PetscCall(DMCreateMatrix_Network_Nest(dm, J));
2473: PetscCall(MatSetDM(*J, dm));
2474: PetscFunctionReturn(PETSC_SUCCESS);
2475: }
2477: if (!network->userEdgeJacobian && !network->userVertexJacobian) {
2478: /* user does not provide Jacobian blocks */
2479: PetscCall(DMCreateMatrix_Plex(network->plex, J));
2480: PetscCall(MatSetDM(*J, dm));
2481: PetscFunctionReturn(PETSC_SUCCESS);
2482: }
2484: PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J));
2485: PetscCall(DMGetGlobalSection(network->plex, §ionGlobal));
2486: PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize));
2487: PetscCall(MatSetSizes(*J, localSize, localSize, PETSC_DETERMINE, PETSC_DETERMINE));
2489: PetscCall(MatSetType(*J, MATAIJ));
2490: PetscCall(MatSetFromOptions(*J));
2492: /* (1) Set matrix preallocation */
2493: /*------------------------------*/
2494: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2495: PetscCall(VecCreate(comm, &vd_nz));
2496: PetscCall(VecSetSizes(vd_nz, localSize, PETSC_DECIDE));
2497: PetscCall(VecSetFromOptions(vd_nz));
2498: PetscCall(VecSet(vd_nz, 0.0));
2499: PetscCall(VecDuplicate(vd_nz, &vo_nz));
2501: /* Set preallocation for edges */
2502: /*-----------------------------*/
2503: PetscCall(DMNetworkGetEdgeRange(dm, &eStart, &eEnd));
2505: PetscCall(PetscMalloc1(localSize, &rows));
2506: for (e = eStart; e < eEnd; e++) {
2507: /* Get row indices */
2508: PetscCall(DMNetworkGetGlobalVecOffset(dm, e, ALL_COMPONENTS, &rstart));
2509: PetscCall(PetscSectionGetDof(network->DofSection, e, &nrows));
2510: if (nrows) {
2511: for (j = 0; j < nrows; j++) rows[j] = j + rstart;
2513: /* Set preallocation for connected vertices */
2514: PetscCall(DMNetworkGetConnectedVertices(dm, e, &cone));
2515: for (v = 0; v < 2; v++) {
2516: PetscCall(PetscSectionGetDof(network->DofSection, cone[v], &ncols));
2518: if (network->Je) {
2519: Juser = network->Je[3 * e + 1 + v]; /* Jacobian(e,v) */
2520: } else Juser = NULL;
2521: PetscCall(DMNetworkIsGhostVertex(dm, cone[v], &ghost));
2522: PetscCall(MatSetPreallocationblock_private(Juser, nrows, rows, ncols, ghost, vd_nz, vo_nz));
2523: }
2525: /* Set preallocation for edge self */
2526: cstart = rstart;
2527: if (network->Je) {
2528: Juser = network->Je[3 * e]; /* Jacobian(e,e) */
2529: } else Juser = NULL;
2530: PetscCall(MatSetPreallocationblock_private(Juser, nrows, rows, nrows, PETSC_FALSE, vd_nz, vo_nz));
2531: }
2532: }
2534: /* Set preallocation for vertices */
2535: /*--------------------------------*/
2536: PetscCall(DMNetworkGetVertexRange(dm, &vStart, &vEnd));
2537: if (vEnd - vStart) vptr = network->Jvptr;
2539: for (v = vStart; v < vEnd; v++) {
2540: /* Get row indices */
2541: PetscCall(DMNetworkGetGlobalVecOffset(dm, v, ALL_COMPONENTS, &rstart));
2542: PetscCall(PetscSectionGetDof(network->DofSection, v, &nrows));
2543: if (!nrows) continue;
2545: PetscCall(DMNetworkIsGhostVertex(dm, v, &ghost));
2546: if (ghost) {
2547: PetscCall(PetscMalloc1(nrows, &rows_v));
2548: } else {
2549: rows_v = rows;
2550: }
2552: for (j = 0; j < nrows; j++) rows_v[j] = j + rstart;
2554: /* Get supporting edges and connected vertices */
2555: PetscCall(DMNetworkGetSupportingEdges(dm, v, &nedges, &edges));
2557: for (e = 0; e < nedges; e++) {
2558: /* Supporting edges */
2559: PetscCall(DMNetworkGetGlobalVecOffset(dm, edges[e], ALL_COMPONENTS, &cstart));
2560: PetscCall(PetscSectionGetDof(network->DofSection, edges[e], &ncols));
2562: if (network->Jv) {
2563: Juser = network->Jv[vptr[v - vStart] + 2 * e + 1]; /* Jacobian(v,e) */
2564: } else Juser = NULL;
2565: PetscCall(MatSetPreallocationblock_private(Juser, nrows, rows_v, ncols, ghost, vd_nz, vo_nz));
2567: /* Connected vertices */
2568: PetscCall(DMNetworkGetConnectedVertices(dm, edges[e], &cone));
2569: vc = (v == cone[0]) ? cone[1] : cone[0];
2570: PetscCall(DMNetworkIsGhostVertex(dm, vc, &ghost_vc));
2572: PetscCall(PetscSectionGetDof(network->DofSection, vc, &ncols));
2574: if (network->Jv) {
2575: Juser = network->Jv[vptr[v - vStart] + 2 * e + 2]; /* Jacobian(v,vc) */
2576: } else Juser = NULL;
2577: if (ghost_vc || ghost) {
2578: ghost2 = PETSC_TRUE;
2579: } else {
2580: ghost2 = PETSC_FALSE;
2581: }
2582: PetscCall(MatSetPreallocationblock_private(Juser, nrows, rows_v, ncols, ghost2, vd_nz, vo_nz));
2583: }
2585: /* Set preallocation for vertex self */
2586: PetscCall(DMNetworkIsGhostVertex(dm, v, &ghost));
2587: if (!ghost) {
2588: PetscCall(DMNetworkGetGlobalVecOffset(dm, v, ALL_COMPONENTS, &cstart));
2589: if (network->Jv) {
2590: Juser = network->Jv[vptr[v - vStart]]; /* Jacobian(v,v) */
2591: } else Juser = NULL;
2592: PetscCall(MatSetPreallocationblock_private(Juser, nrows, rows_v, nrows, PETSC_FALSE, vd_nz, vo_nz));
2593: }
2594: if (ghost) PetscCall(PetscFree(rows_v));
2595: }
2597: PetscCall(VecAssemblyBegin(vd_nz));
2598: PetscCall(VecAssemblyBegin(vo_nz));
2600: PetscCall(PetscMalloc2(localSize, &dnnz, localSize, &onnz));
2602: PetscCall(VecAssemblyEnd(vd_nz));
2603: PetscCall(VecAssemblyEnd(vo_nz));
2605: PetscCall(VecGetArray(vd_nz, &vdnz));
2606: PetscCall(VecGetArray(vo_nz, &vonz));
2607: for (j = 0; j < localSize; j++) {
2608: dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
2609: onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
2610: }
2611: PetscCall(VecRestoreArray(vd_nz, &vdnz));
2612: PetscCall(VecRestoreArray(vo_nz, &vonz));
2613: PetscCall(VecDestroy(&vd_nz));
2614: PetscCall(VecDestroy(&vo_nz));
2616: PetscCall(MatSeqAIJSetPreallocation(*J, 0, dnnz));
2617: PetscCall(MatMPIAIJSetPreallocation(*J, 0, dnnz, 0, onnz));
2618: PetscCall(MatSetOption(*J, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
2620: PetscCall(PetscFree2(dnnz, onnz));
2622: /* (2) Set matrix entries for edges */
2623: /*----------------------------------*/
2624: for (e = eStart; e < eEnd; e++) {
2625: /* Get row indices */
2626: PetscCall(DMNetworkGetGlobalVecOffset(dm, e, ALL_COMPONENTS, &rstart));
2627: PetscCall(PetscSectionGetDof(network->DofSection, e, &nrows));
2628: if (nrows) {
2629: for (j = 0; j < nrows; j++) rows[j] = j + rstart;
2631: /* Set matrix entries for connected vertices */
2632: PetscCall(DMNetworkGetConnectedVertices(dm, e, &cone));
2633: for (v = 0; v < 2; v++) {
2634: PetscCall(DMNetworkGetGlobalVecOffset(dm, cone[v], ALL_COMPONENTS, &cstart));
2635: PetscCall(PetscSectionGetDof(network->DofSection, cone[v], &ncols));
2637: if (network->Je) {
2638: Juser = network->Je[3 * e + 1 + v]; /* Jacobian(e,v) */
2639: } else Juser = NULL;
2640: PetscCall(MatSetblock_private(Juser, nrows, rows, ncols, cstart, J));
2641: }
2643: /* Set matrix entries for edge self */
2644: cstart = rstart;
2645: if (network->Je) {
2646: Juser = network->Je[3 * e]; /* Jacobian(e,e) */
2647: } else Juser = NULL;
2648: PetscCall(MatSetblock_private(Juser, nrows, rows, nrows, cstart, J));
2649: }
2650: }
2652: /* Set matrix entries for vertices */
2653: /*---------------------------------*/
2654: for (v = vStart; v < vEnd; v++) {
2655: /* Get row indices */
2656: PetscCall(DMNetworkGetGlobalVecOffset(dm, v, ALL_COMPONENTS, &rstart));
2657: PetscCall(PetscSectionGetDof(network->DofSection, v, &nrows));
2658: if (!nrows) continue;
2660: PetscCall(DMNetworkIsGhostVertex(dm, v, &ghost));
2661: if (ghost) {
2662: PetscCall(PetscMalloc1(nrows, &rows_v));
2663: } else {
2664: rows_v = rows;
2665: }
2666: for (j = 0; j < nrows; j++) rows_v[j] = j + rstart;
2668: /* Get supporting edges and connected vertices */
2669: PetscCall(DMNetworkGetSupportingEdges(dm, v, &nedges, &edges));
2671: for (e = 0; e < nedges; e++) {
2672: /* Supporting edges */
2673: PetscCall(DMNetworkGetGlobalVecOffset(dm, edges[e], ALL_COMPONENTS, &cstart));
2674: PetscCall(PetscSectionGetDof(network->DofSection, edges[e], &ncols));
2676: if (network->Jv) {
2677: Juser = network->Jv[vptr[v - vStart] + 2 * e + 1]; /* Jacobian(v,e) */
2678: } else Juser = NULL;
2679: PetscCall(MatSetblock_private(Juser, nrows, rows_v, ncols, cstart, J));
2681: /* Connected vertices */
2682: PetscCall(DMNetworkGetConnectedVertices(dm, edges[e], &cone));
2683: vc = (v == cone[0]) ? cone[1] : cone[0];
2685: PetscCall(DMNetworkGetGlobalVecOffset(dm, vc, ALL_COMPONENTS, &cstart));
2686: PetscCall(PetscSectionGetDof(network->DofSection, vc, &ncols));
2688: if (network->Jv) {
2689: Juser = network->Jv[vptr[v - vStart] + 2 * e + 2]; /* Jacobian(v,vc) */
2690: } else Juser = NULL;
2691: PetscCall(MatSetblock_private(Juser, nrows, rows_v, ncols, cstart, J));
2692: }
2694: /* Set matrix entries for vertex self */
2695: if (!ghost) {
2696: PetscCall(DMNetworkGetGlobalVecOffset(dm, v, ALL_COMPONENTS, &cstart));
2697: if (network->Jv) {
2698: Juser = network->Jv[vptr[v - vStart]]; /* Jacobian(v,v) */
2699: } else Juser = NULL;
2700: PetscCall(MatSetblock_private(Juser, nrows, rows_v, nrows, cstart, J));
2701: }
2702: if (ghost) PetscCall(PetscFree(rows_v));
2703: }
2704: PetscCall(PetscFree(rows));
2706: PetscCall(MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY));
2707: PetscCall(MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY));
2709: PetscCall(MatSetDM(*J, dm));
2710: PetscFunctionReturn(PETSC_SUCCESS);
2711: }
2713: static PetscErrorCode DMNetworkDestroyComponentData(DM dm)
2714: {
2715: DM_Network *network = (DM_Network *)dm->data;
2716: PetscInt j, np;
2718: PetscFunctionBegin;
2719: if (network->header) {
2720: np = network->cloneshared->pEnd - network->cloneshared->pStart;
2721: for (j = 0; j < np; j++) {
2722: PetscCall(PetscFree5(network->header[j].size, network->header[j].key, network->header[j].offset, network->header[j].nvar, network->header[j].offsetvarrel));
2723: PetscCall(PetscFree(network->cvalue[j].data));
2724: }
2725: PetscCall(PetscFree2(network->header, network->cvalue));
2726: }
2727: PetscFunctionReturn(PETSC_SUCCESS);
2728: }
2730: PetscErrorCode DMDestroy_Network(DM dm)
2731: {
2732: DM_Network *network = (DM_Network *)dm->data;
2733: PetscInt j;
2735: PetscFunctionBegin;
2736: /*
2737: Developers Note: Due to the mixed topological definition of DMNetwork and data defined ON the
2738: network like DofSection, DataSection, *componentdataarray, and friends, when cloning, we share
2739: only the true topological data, and make our own data ON the network. Thus refct only refers
2740: to the number of references to topological data, and data ON the network is always destroyed.
2741: It is understood this is atypical for a DM, but is very intentional.
2742: */
2744: /* Always destroy data ON the network */
2745: PetscCall(PetscFree(network->Je));
2746: if (network->Jv) {
2747: PetscCall(PetscFree(network->Jvptr));
2748: PetscCall(PetscFree(network->Jv));
2749: }
2750: PetscCall(PetscSectionDestroy(&network->DataSection));
2751: PetscCall(PetscSectionDestroy(&network->DofSection));
2752: PetscCall(PetscFree(network->component));
2753: PetscCall(PetscFree(network->componentdataarray));
2754: PetscCall(DMNetworkDestroyComponentData(dm));
2756: PetscCall(DMDestroy(&network->plex)); /* this is cloned in DMClone_Network, so safe to destroy */
2758: /*
2759: Developers Note: The DMNetworkVertexInfo and DMNetworkEdgeInfo data structures are completely
2760: destroyed as they are again a mix of topological data:
2761: ISLocalToGlobalMapping mapping;
2762: PetscSF sf;
2763: and data ON the network
2764: PetscSection DofSection;
2765: PetscSection GlobalDofSection;
2766: And the only way to access them currently is through DMNetworkAssembleGraphStructures which assembles
2767: everything. So we must destroy everything and require DMNetworkAssembleGraphStructures is called again
2768: for any clone.
2769: */
2770: PetscCall(ISLocalToGlobalMappingDestroy(&network->vertex.mapping));
2771: PetscCall(PetscSectionDestroy(&network->vertex.DofSection));
2772: PetscCall(PetscSectionDestroy(&network->vertex.GlobalDofSection));
2773: PetscCall(PetscSFDestroy(&network->vertex.sf));
2774: /* edge */
2775: PetscCall(ISLocalToGlobalMappingDestroy(&network->edge.mapping));
2776: PetscCall(PetscSectionDestroy(&network->edge.DofSection));
2777: PetscCall(PetscSectionDestroy(&network->edge.GlobalDofSection));
2778: PetscCall(PetscSFDestroy(&network->edge.sf));
2779: /* Destroy the potentially cloneshared data */
2780: if (--network->cloneshared->refct <= 0) {
2781: /* Developer Note: I'm not sure if vltog can be reused or not, as I'm not sure what it's purpose is. I
2782: naively think it can be reused. */
2783: PetscCall(PetscFree(network->cloneshared->vltog));
2784: for (j = 0; j < network->cloneshared->Nsvtx; j++) PetscCall(PetscFree(network->cloneshared->svtx[j].sv));
2785: PetscCall(PetscFree(network->cloneshared->svtx));
2786: PetscCall(PetscFree2(network->cloneshared->subnetedge, network->cloneshared->subnetvtx));
2787: PetscCall(PetscHMapIDestroy(&network->cloneshared->svtable));
2788: PetscCall(PetscFree(network->cloneshared->subnet));
2789: PetscCall(PetscFree(network->cloneshared));
2790: }
2791: PetscCall(PetscFree(network)); /* Always freed as this structure is copied in a clone, not cloneshared */
2792: PetscFunctionReturn(PETSC_SUCCESS);
2793: }
2795: PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
2796: {
2797: DM_Network *network = (DM_Network *)dm->data;
2799: PetscFunctionBegin;
2800: PetscCall(DMGlobalToLocalBegin(network->plex, g, mode, l));
2801: PetscFunctionReturn(PETSC_SUCCESS);
2802: }
2804: PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
2805: {
2806: DM_Network *network = (DM_Network *)dm->data;
2808: PetscFunctionBegin;
2809: PetscCall(DMGlobalToLocalEnd(network->plex, g, mode, l));
2810: PetscFunctionReturn(PETSC_SUCCESS);
2811: }
2813: PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
2814: {
2815: DM_Network *network = (DM_Network *)dm->data;
2817: PetscFunctionBegin;
2818: PetscCall(DMLocalToGlobalBegin(network->plex, l, mode, g));
2819: PetscFunctionReturn(PETSC_SUCCESS);
2820: }
2822: PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
2823: {
2824: DM_Network *network = (DM_Network *)dm->data;
2826: PetscFunctionBegin;
2827: PetscCall(DMLocalToGlobalEnd(network->plex, l, mode, g));
2828: PetscFunctionReturn(PETSC_SUCCESS);
2829: }
2831: /*@
2832: DMNetworkGetVertexLocalToGlobalOrdering - Get vertex global index
2834: Not collective
2836: Input Parameters:
2837: + dm - the dm object
2838: - vloc - local vertex ordering, start from 0
2840: Output Parameters:
2841: . vg - global vertex ordering, start from 0
2843: Level: advanced
2845: .seealso: `DMNetworkSetVertexLocalToGlobalOrdering()`
2846: @*/
2847: PetscErrorCode DMNetworkGetVertexLocalToGlobalOrdering(DM dm, PetscInt vloc, PetscInt *vg)
2848: {
2849: DM_Network *network = (DM_Network *)dm->data;
2850: PetscInt *vltog = network->cloneshared->vltog;
2852: PetscFunctionBegin;
2853: PetscCheck(vltog, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Must call DMNetworkSetVertexLocalToGlobalOrdering() first");
2854: *vg = vltog[vloc];
2855: PetscFunctionReturn(PETSC_SUCCESS);
2856: }
2858: /*@
2859: DMNetworkSetVertexLocalToGlobalOrdering - Create and setup vertex local to global map
2861: Collective
2863: Input Parameters:
2864: . dm - the dm object
2866: Level: advanced
2868: .seealso: `DMNetworkGetGlobalVertexIndex()`
2869: @*/
2870: PetscErrorCode DMNetworkSetVertexLocalToGlobalOrdering(DM dm)
2871: {
2872: DM_Network *network = (DM_Network *)dm->data;
2873: MPI_Comm comm;
2874: PetscMPIInt rank, size, *displs = NULL, *recvcounts = NULL, remoterank;
2875: PetscBool ghost;
2876: PetscInt *vltog, nroots, nleaves, i, *vrange, k, N, lidx;
2877: const PetscSFNode *iremote;
2878: PetscSF vsf;
2879: Vec Vleaves, Vleaves_seq;
2880: VecScatter ctx;
2881: PetscScalar *varr, val;
2882: const PetscScalar *varr_read;
2884: PetscFunctionBegin;
2885: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2886: PetscCallMPI(MPI_Comm_size(comm, &size));
2887: PetscCallMPI(MPI_Comm_rank(comm, &rank));
2889: if (size == 1) {
2890: nroots = network->cloneshared->vEnd - network->cloneshared->vStart;
2891: PetscCall(PetscMalloc1(nroots, &vltog));
2892: for (i = 0; i < nroots; i++) vltog[i] = i;
2893: network->cloneshared->vltog = vltog;
2894: PetscFunctionReturn(PETSC_SUCCESS);
2895: }
2897: PetscCheck(network->cloneshared->distributecalled, comm, PETSC_ERR_ARG_WRONGSTATE, "Must call DMNetworkDistribute() first");
2898: if (network->cloneshared->vltog) PetscCall(PetscFree(network->cloneshared->vltog));
2900: PetscCall(DMNetworkSetSubMap_private(network->cloneshared->vStart, network->cloneshared->vEnd, &network->vertex.mapping));
2901: PetscCall(PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf));
2902: vsf = network->vertex.sf;
2904: PetscCall(PetscMalloc3(size + 1, &vrange, size, &displs, size, &recvcounts));
2905: PetscCall(PetscSFGetGraph(vsf, &nroots, &nleaves, NULL, &iremote));
2907: for (i = 0; i < size; i++) {
2908: displs[i] = i;
2909: recvcounts[i] = 1;
2910: }
2912: i = nroots - nleaves; /* local number of vertices, excluding ghosts */
2913: vrange[0] = 0;
2914: PetscCallMPI(MPI_Allgatherv(&i, 1, MPIU_INT, vrange + 1, recvcounts, displs, MPIU_INT, comm));
2915: for (i = 2; i < size + 1; i++) vrange[i] += vrange[i - 1];
2917: PetscCall(PetscMalloc1(nroots, &vltog));
2918: network->cloneshared->vltog = vltog;
2920: /* Set vltog for non-ghost vertices */
2921: k = 0;
2922: for (i = 0; i < nroots; i++) {
2923: PetscCall(DMNetworkIsGhostVertex(dm, i + network->cloneshared->vStart, &ghost));
2924: if (ghost) continue;
2925: vltog[i] = vrange[rank] + k++;
2926: }
2927: PetscCall(PetscFree3(vrange, displs, recvcounts));
2929: /* Set vltog for ghost vertices */
2930: /* (a) create parallel Vleaves and sequential Vleaves_seq to convert local iremote[*].index to global index */
2931: PetscCall(VecCreate(comm, &Vleaves));
2932: PetscCall(VecSetSizes(Vleaves, 2 * nleaves, PETSC_DETERMINE));
2933: PetscCall(VecSetFromOptions(Vleaves));
2934: PetscCall(VecGetArray(Vleaves, &varr));
2935: for (i = 0; i < nleaves; i++) {
2936: varr[2 * i] = (PetscScalar)(iremote[i].rank); /* rank of remote process */
2937: varr[2 * i + 1] = (PetscScalar)(iremote[i].index); /* local index in remote process */
2938: }
2939: PetscCall(VecRestoreArray(Vleaves, &varr));
2941: /* (b) scatter local info to remote processes via VecScatter() */
2942: PetscCall(VecScatterCreateToAll(Vleaves, &ctx, &Vleaves_seq));
2943: PetscCall(VecScatterBegin(ctx, Vleaves, Vleaves_seq, INSERT_VALUES, SCATTER_FORWARD));
2944: PetscCall(VecScatterEnd(ctx, Vleaves, Vleaves_seq, INSERT_VALUES, SCATTER_FORWARD));
2946: /* (c) convert local indices to global indices in parallel vector Vleaves */
2947: PetscCall(VecGetSize(Vleaves_seq, &N));
2948: PetscCall(VecGetArrayRead(Vleaves_seq, &varr_read));
2949: for (i = 0; i < N; i += 2) {
2950: remoterank = (PetscMPIInt)PetscRealPart(varr_read[i]);
2951: if (remoterank == rank) {
2952: k = i + 1; /* row number */
2953: lidx = (PetscInt)PetscRealPart(varr_read[i + 1]);
2954: val = (PetscScalar)vltog[lidx]; /* global index for non-ghost vertex computed above */
2955: PetscCall(VecSetValues(Vleaves, 1, &k, &val, INSERT_VALUES));
2956: }
2957: }
2958: PetscCall(VecRestoreArrayRead(Vleaves_seq, &varr_read));
2959: PetscCall(VecAssemblyBegin(Vleaves));
2960: PetscCall(VecAssemblyEnd(Vleaves));
2962: /* (d) Set vltog for ghost vertices by copying local values of Vleaves */
2963: PetscCall(VecGetArrayRead(Vleaves, &varr_read));
2964: k = 0;
2965: for (i = 0; i < nroots; i++) {
2966: PetscCall(DMNetworkIsGhostVertex(dm, i + network->cloneshared->vStart, &ghost));
2967: if (!ghost) continue;
2968: vltog[i] = (PetscInt)PetscRealPart(varr_read[2 * k + 1]);
2969: k++;
2970: }
2971: PetscCall(VecRestoreArrayRead(Vleaves, &varr_read));
2973: PetscCall(VecDestroy(&Vleaves));
2974: PetscCall(VecDestroy(&Vleaves_seq));
2975: PetscCall(VecScatterDestroy(&ctx));
2976: PetscFunctionReturn(PETSC_SUCCESS);
2977: }
2979: static inline PetscErrorCode DMISAddSize_private(DM_Network *network, PetscInt p, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *nidx)
2980: {
2981: PetscInt i, j, ncomps, nvar, key, offset = 0;
2982: DMNetworkComponentHeader header;
2984: PetscFunctionBegin;
2985: PetscCall(PetscSectionGetOffset(network->DataSection, p, &offset));
2986: ncomps = ((DMNetworkComponentHeader)(network->componentdataarray + offset))->ndata;
2987: header = (DMNetworkComponentHeader)(network->componentdataarray + offset);
2989: for (i = 0; i < ncomps; i++) {
2990: key = header->key[i];
2991: nvar = header->nvar[i];
2992: for (j = 0; j < numkeys; j++) {
2993: if (key == keys[j]) {
2994: if (!blocksize || blocksize[j] == -1) {
2995: *nidx += nvar;
2996: } else {
2997: *nidx += nselectedvar[j] * nvar / blocksize[j];
2998: }
2999: }
3000: }
3001: }
3002: PetscFunctionReturn(PETSC_SUCCESS);
3003: }
3005: static inline PetscErrorCode DMISComputeIdx_private(DM dm, PetscInt p, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *selectedvar[], PetscInt *ii, PetscInt *idx)
3006: {
3007: PetscInt i, j, ncomps, nvar, key, offsetg, k, k1, offset = 0;
3008: DM_Network *network = (DM_Network *)dm->data;
3009: DMNetworkComponentHeader header;
3011: PetscFunctionBegin;
3012: PetscCall(PetscSectionGetOffset(network->DataSection, p, &offset));
3013: ncomps = ((DMNetworkComponentHeader)(network->componentdataarray + offset))->ndata;
3014: header = (DMNetworkComponentHeader)(network->componentdataarray + offset);
3016: for (i = 0; i < ncomps; i++) {
3017: key = header->key[i];
3018: nvar = header->nvar[i];
3019: for (j = 0; j < numkeys; j++) {
3020: if (key != keys[j]) continue;
3022: PetscCall(DMNetworkGetGlobalVecOffset(dm, p, i, &offsetg));
3023: if (!blocksize || blocksize[j] == -1) {
3024: for (k = 0; k < nvar; k++) idx[(*ii)++] = offsetg + k;
3025: } else {
3026: for (k = 0; k < nvar; k += blocksize[j]) {
3027: for (k1 = 0; k1 < nselectedvar[j]; k1++) idx[(*ii)++] = offsetg + k + selectedvar[j][k1];
3028: }
3029: }
3030: }
3031: }
3032: PetscFunctionReturn(PETSC_SUCCESS);
3033: }
3035: /*@
3036: DMNetworkCreateIS - Create an index set object from the global vector of the network
3038: Collective
3040: Input Parameters:
3041: + dm - DMNetwork object
3042: . numkeys - number of keys
3043: . keys - array of keys that define the components of the variables you wish to extract
3044: . blocksize - block size of the variables associated to the component
3045: . nselectedvar - number of variables in each block to select
3046: - selectedvar - the offset into the block of each variable in each block to select
3048: Output Parameters:
3049: . is - the index set
3051: Level: Advanced
3053: Notes:
3054: 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.
3056: .seealso: `DMNetworkCreate()`, `ISCreateGeneral()`, `DMNetworkCreateLocalIS()`
3057: @*/
3058: PetscErrorCode DMNetworkCreateIS(DM dm, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *selectedvar[], IS *is)
3059: {
3060: MPI_Comm comm;
3061: DM_Network *network = (DM_Network *)dm->data;
3062: PetscInt i, p, estart, eend, vstart, vend, nidx, *idx;
3063: PetscBool ghost;
3065: PetscFunctionBegin;
3066: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3068: /* Check input parameters */
3069: for (i = 0; i < numkeys; i++) {
3070: if (!blocksize || blocksize[i] == -1) continue;
3071: 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]);
3072: }
3074: PetscCall(DMNetworkGetEdgeRange(dm, &estart, &eend));
3075: PetscCall(DMNetworkGetVertexRange(dm, &vstart, &vend));
3077: /* Get local number of idx */
3078: nidx = 0;
3079: for (p = estart; p < eend; p++) PetscCall(DMISAddSize_private(network, p, numkeys, keys, blocksize, nselectedvar, &nidx));
3080: for (p = vstart; p < vend; p++) {
3081: PetscCall(DMNetworkIsGhostVertex(dm, p, &ghost));
3082: if (ghost) continue;
3083: PetscCall(DMISAddSize_private(network, p, numkeys, keys, blocksize, nselectedvar, &nidx));
3084: }
3086: /* Compute idx */
3087: PetscCall(PetscMalloc1(nidx, &idx));
3088: i = 0;
3089: for (p = estart; p < eend; p++) PetscCall(DMISComputeIdx_private(dm, p, numkeys, keys, blocksize, nselectedvar, selectedvar, &i, idx));
3090: for (p = vstart; p < vend; p++) {
3091: PetscCall(DMNetworkIsGhostVertex(dm, p, &ghost));
3092: if (ghost) continue;
3093: PetscCall(DMISComputeIdx_private(dm, p, numkeys, keys, blocksize, nselectedvar, selectedvar, &i, idx));
3094: }
3096: /* Create is */
3097: PetscCall(ISCreateGeneral(comm, nidx, idx, PETSC_COPY_VALUES, is));
3098: PetscCall(PetscFree(idx));
3099: PetscFunctionReturn(PETSC_SUCCESS);
3100: }
3102: static inline PetscErrorCode DMISComputeLocalIdx_private(DM dm, PetscInt p, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *selectedvar[], PetscInt *ii, PetscInt *idx)
3103: {
3104: PetscInt i, j, ncomps, nvar, key, offsetl, k, k1, offset = 0;
3105: DM_Network *network = (DM_Network *)dm->data;
3106: DMNetworkComponentHeader header;
3108: PetscFunctionBegin;
3109: PetscCall(PetscSectionGetOffset(network->DataSection, p, &offset));
3110: ncomps = ((DMNetworkComponentHeader)(network->componentdataarray + offset))->ndata;
3111: header = (DMNetworkComponentHeader)(network->componentdataarray + offset);
3113: for (i = 0; i < ncomps; i++) {
3114: key = header->key[i];
3115: nvar = header->nvar[i];
3116: for (j = 0; j < numkeys; j++) {
3117: if (key != keys[j]) continue;
3119: PetscCall(DMNetworkGetLocalVecOffset(dm, p, i, &offsetl));
3120: if (!blocksize || blocksize[j] == -1) {
3121: for (k = 0; k < nvar; k++) idx[(*ii)++] = offsetl + k;
3122: } else {
3123: for (k = 0; k < nvar; k += blocksize[j]) {
3124: for (k1 = 0; k1 < nselectedvar[j]; k1++) idx[(*ii)++] = offsetl + k + selectedvar[j][k1];
3125: }
3126: }
3127: }
3128: }
3129: PetscFunctionReturn(PETSC_SUCCESS);
3130: }
3132: /*@
3133: DMNetworkCreateLocalIS - Create an index set object from the local vector of the network
3135: Not collective
3137: Input Parameters:
3138: + dm - DMNetwork object
3139: . numkeys - number of keys
3140: . keys - array of keys that define the components of the variables you wish to extract
3141: . blocksize - block size of the variables associated to the component
3142: . nselectedvar - number of variables in each block to select
3143: - selectedvar - the offset into the block of each variable in each block to select
3145: Output Parameters:
3146: . is - the index set
3148: Level: Advanced
3150: Notes:
3151: 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.
3153: .seealso: `DMNetworkCreate()`, `DMNetworkCreateIS`, `ISCreateGeneral()`
3154: @*/
3155: PetscErrorCode DMNetworkCreateLocalIS(DM dm, PetscInt numkeys, PetscInt keys[], PetscInt blocksize[], PetscInt nselectedvar[], PetscInt *selectedvar[], IS *is)
3156: {
3157: DM_Network *network = (DM_Network *)dm->data;
3158: PetscInt i, p, pstart, pend, nidx, *idx;
3160: PetscFunctionBegin;
3161: /* Check input parameters */
3162: for (i = 0; i < numkeys; i++) {
3163: if (!blocksize || blocksize[i] == -1) continue;
3164: 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]);
3165: }
3167: pstart = network->cloneshared->pStart;
3168: pend = network->cloneshared->pEnd;
3170: /* Get local number of idx */
3171: nidx = 0;
3172: for (p = pstart; p < pend; p++) PetscCall(DMISAddSize_private(network, p, numkeys, keys, blocksize, nselectedvar, &nidx));
3174: /* Compute local idx */
3175: PetscCall(PetscMalloc1(nidx, &idx));
3176: i = 0;
3177: for (p = pstart; p < pend; p++) PetscCall(DMISComputeLocalIdx_private(dm, p, numkeys, keys, blocksize, nselectedvar, selectedvar, &i, idx));
3179: /* Create is */
3180: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nidx, idx, PETSC_COPY_VALUES, is));
3181: PetscCall(PetscFree(idx));
3182: PetscFunctionReturn(PETSC_SUCCESS);
3183: }
3184: /*@
3185: DMNetworkFinalizeComponents - Sets up internal data structures for the sections and components. It is called after registering new components and adding all components
3186: to the cloned network. After calling this subroutine, no new components can be added to the network.
3188: Collective
3190: Input Parameters:
3191: . dm - the dmnetwork object
3193: Level: beginner
3195: .seealso: `DMNetworkAddComponent()`, `DMNetworkRegisterComponent()`, `DMSetUp()`
3196: @*/
3197: PetscErrorCode DMNetworkFinalizeComponents(DM dm)
3198: {
3199: DM_Network *network = (DM_Network *)dm->data;
3201: PetscFunctionBegin;
3202: if (network->componentsetup) PetscFunctionReturn(PETSC_SUCCESS);
3203: PetscCall(DMNetworkComponentSetUp(dm));
3204: PetscCall(DMNetworkVariablesSetUp(dm));
3205: PetscCall(DMSetLocalSection(network->plex, network->DofSection));
3206: PetscCall(DMGetGlobalSection(network->plex, &network->GlobalDofSection));
3207: network->componentsetup = PETSC_TRUE;
3208: PetscFunctionReturn(PETSC_SUCCESS);
3209: }