Actual source code: network.c
1: #include <petsc/private/dmnetworkimpl.h>
3: /*
4: Creates the component header and value objects for a network point
5: */
6: static PetscErrorCode SetUpNetworkHeaderComponentValue(DM dm,DMNetworkComponentHeader header,DMNetworkComponentValue cvalue)
7: {
8: /* Allocate arrays for component information */
9: PetscCalloc5(header->maxcomps,&header->size,header->maxcomps,&header->key,header->maxcomps,&header->offset,header->maxcomps,&header->nvar,header->maxcomps,&header->offsetvarrel);
10: PetscCalloc1(header->maxcomps,&cvalue->data);
12: /* 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.
13: If the data header struct changes then this header size calculation needs to be updated. */
14: header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5*header->maxcomps*sizeof(PetscInt);
15: header->hsize /= sizeof(DMNetworkComponentGenericDataType);
16: return 0;
17: }
19: /*@
20: DMNetworkGetPlex - Gets the Plex DM associated with this network DM
22: Not collective
24: Input Parameters:
25: . dm - the dm object
27: Output Parameters:
28: . plexdm - the plex dm object
30: Level: Advanced
32: .seealso: DMNetworkCreate()
33: @*/
34: PetscErrorCode DMNetworkGetPlex(DM dm,DM *plexdm)
35: {
36: DM_Network *network = (DM_Network*)dm->data;
38: *plexdm = network->plex;
39: return 0;
40: }
42: /*@
43: DMNetworkGetNumSubNetworks - Gets the the number of subnetworks
45: Not collective
47: Input Parameter:
48: . dm - the dm object
50: Output Parameters:
51: + nsubnet - local number of subnetworks
52: - Nsubnet - global number of subnetworks
54: Level: beginner
56: .seealso: DMNetworkCreate(), DMNetworkSetNumSubNetworks()
57: @*/
58: PetscErrorCode DMNetworkGetNumSubNetworks(DM dm,PetscInt *nsubnet,PetscInt *Nsubnet)
59: {
60: DM_Network *network = (DM_Network*)dm->data;
62: if (nsubnet) *nsubnet = network->nsubnet;
63: if (Nsubnet) *Nsubnet = network->Nsubnet;
64: return 0;
65: }
67: /*@
68: DMNetworkSetNumSubNetworks - Sets the number of subnetworks
70: Collective on dm
72: Input Parameters:
73: + dm - the dm object
74: . nsubnet - local number of subnetworks
75: - Nsubnet - global number of subnetworks
77: Level: beginner
79: .seealso: DMNetworkCreate(), DMNetworkGetNumSubNetworks()
80: @*/
81: PetscErrorCode DMNetworkSetNumSubNetworks(DM dm,PetscInt nsubnet,PetscInt Nsubnet)
82: {
83: DM_Network *network = (DM_Network*)dm->data;
91: if (Nsubnet == PETSC_DECIDE) {
93: MPIU_Allreduce(&nsubnet,&Nsubnet,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));
94: }
97: network->Nsubnet = Nsubnet;
98: network->nsubnet = 0; /* initia value; will be determind by DMNetworkAddSubnetwork() */
99: PetscCalloc1(Nsubnet,&network->subnet);
101: /* num of shared vertices */
102: network->nsvtx = 0;
103: network->Nsvtx = 0;
104: return 0;
105: }
107: /*@
108: DMNetworkAddSubnetwork - Add a subnetwork
110: Collective on dm
112: Input Parameters:
113: + dm - the dm object
114: . name - name of the subnetwork
115: . ne - number of local edges of this subnetwork
116: - 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,
117: $ [first vertex of first edge, second vertex of first edge, first vertex of second edge, second vertex of second edge, etc]
119: Output Parameters:
120: . netnum - global index of the subnetwork
122: Notes:
123: There is no copy involved in this operation, only the pointer is referenced. The edgelist should
124: not be destroyed before the call to DMNetworkLayoutSetUp()
126: 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,
127: 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.
129: Level: beginner
131: Example usage:
132: Consider the following networks:
133: 1) A sigle subnetwork:
134: .vb
135: network 0:
136: rank[0]:
137: v0 --> v2; v1 --> v2
138: rank[1]:
139: v3 --> v5; v4 --> v5
140: .ve
142: The resulting input
143: network 0:
144: rank[0]:
145: ne = 2
146: edgelist = [0 2 | 1 2]
147: rank[1]:
148: ne = 2
149: edgelist = [3 5 | 4 5]
151: 2) Two subnetworks:
152: .vb
153: subnetwork 0:
154: rank[0]:
155: v0 --> v2; v2 --> v1; v1 --> v3;
156: subnetwork 1:
157: rank[1]:
158: v0 --> v3; v3 --> v2; v2 --> v1;
159: .ve
161: The resulting input
162: subnetwork 0:
163: rank[0]:
164: ne = 3
165: edgelist = [0 2 | 2 1 | 1 3]
166: rank[1]:
167: ne = 0
168: edgelist = NULL
170: subnetwork 1:
171: rank[0]:
172: ne = 0
173: edgelist = NULL
174: rank[1]:
175: edgelist = [0 3 | 3 2 | 2 1]
177: .seealso: DMNetworkCreate(), DMNetworkSetNumSubnetworks()
178: @*/
179: PetscErrorCode DMNetworkAddSubnetwork(DM dm,const char* name,PetscInt ne,PetscInt edgelist[],PetscInt *netnum)
180: {
181: DM_Network *network = (DM_Network*)dm->data;
182: PetscInt i,Nedge,j,Nvtx,nvtx,nvtx_min=-1,nvtx_max=0;
183: PetscBT table;
185: for (i=0; i<ne; i++) {
187: }
189: i = 0;
190: if (ne) nvtx_min = nvtx_max = edgelist[0];
191: for (j=0; j<ne; j++) {
192: nvtx_min = PetscMin(nvtx_min, edgelist[i]);
193: nvtx_max = PetscMax(nvtx_max, edgelist[i]);
194: i++;
195: nvtx_min = PetscMin(nvtx_min, edgelist[i]);
196: nvtx_max = PetscMax(nvtx_max, edgelist[i]);
197: i++;
198: }
199: Nvtx = nvtx_max - nvtx_min + 1; /* approximated total local nvtx for this subnet */
201: /* Get exact local nvtx for this subnet: counting local values between nvtx_min and nvtx_max */
202: PetscBTCreate(Nvtx,&table);
203: PetscBTMemzero(Nvtx,table);
204: i = 0;
205: for (j=0; j<ne; j++) {
206: PetscBTSet(table,edgelist[i++]-nvtx_min);
207: PetscBTSet(table,edgelist[i++]-nvtx_min);
208: }
209: nvtx = 0;
210: for (j=0; j<Nvtx; j++) {
211: if (PetscBTLookup(table,j)) nvtx++;
212: }
214: /* Get global total Nvtx = max(edgelist[])+1 for this subnet */
215: MPIU_Allreduce(&nvtx_max,&Nvtx,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)dm));
216: Nvtx++;
217: PetscBTDestroy(&table);
219: /* Get global total Nedge for this subnet */
220: MPIU_Allreduce(&ne,&Nedge,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));
222: i = network->nsubnet;
223: if (name) {
224: PetscStrcpy(network->subnet[i].name,name);
225: }
226: network->subnet[i].nvtx = nvtx; /* include ghost vertices */
227: network->subnet[i].nedge = ne;
228: network->subnet[i].edgelist = edgelist;
229: network->subnet[i].Nvtx = Nvtx;
230: network->subnet[i].Nedge = Nedge;
232: /* ----------------------------------------------------------
233: p=v or e;
234: subnet[0].pStart = 0
235: subnet[i+1].pStart = subnet[i].pEnd = subnet[i].pStart + (nE[i] or NV[i])
236: ----------------------------------------------------------------------- */
237: /* GLOBAL subnet[].vStart and vEnd, used by DMNetworkLayoutSetUp() */
238: network->subnet[i].vStart = network->NVertices;
239: network->subnet[i].vEnd = network->subnet[i].vStart + network->subnet[i].Nvtx; /* global vEnd of subnet[i] */
241: network->nVertices += nvtx; /* include ghost vertices */
242: network->NVertices += network->subnet[i].Nvtx;
244: /* LOCAL subnet[].eStart and eEnd, used by DMNetworkLayoutSetUp() */
245: network->subnet[i].eStart = network->nEdges;
246: network->subnet[i].eEnd = network->subnet[i].eStart + ne;
247: network->nEdges += ne;
248: network->NEdges += network->subnet[i].Nedge;
250: PetscStrcpy(network->subnet[i].name,name);
251: if (netnum) *netnum = network->nsubnet;
252: network->nsubnet++;
253: return 0;
254: }
256: /*@C
257: DMNetworkSharedVertexGetInfo - Get info of a shared vertex struct, see petsc/private/dmnetworkimpl.h
259: Not collective
261: Input Parameters:
262: + dm - the DM object
263: - v - vertex point
265: Output Parameters:
266: + gidx - global number of this shared vertex in the internal dmplex
267: . n - number of subnetworks that share this vertex
268: - sv - array of size n: sv[2*i,2*i+1]=(net[i], idx[i]), i=0,...,n-1
270: Level: intermediate
272: .seealso: DMNetworkGetSharedVertices()
273: @*/
274: PetscErrorCode DMNetworkSharedVertexGetInfo(DM dm,PetscInt v,PetscInt *gidx,PetscInt *n,const PetscInt **sv)
275: {
276: DM_Network *network = (DM_Network*)dm->data;
277: SVtx *svtx = network->svtx;
278: PetscInt i,gidx_tmp;
280: DMNetworkGetGlobalVertexIndex(dm,v,&gidx_tmp);
281: PetscTableFind(network->svtable,gidx_tmp+1,&i);
284: i--;
285: if (gidx) *gidx = gidx_tmp;
286: if (n) *n = svtx[i].n;
287: if (sv) *sv = svtx[i].sv;
288: return 0;
289: }
291: /*
292: VtxGetInfo - Get info of an input vertex=(net,idx)
294: Input Parameters:
295: + Nsvtx - global num of shared vertices
296: . svtx - array of shared vertices (global)
297: - (net,idx) - subnet number and local index for a vertex
299: Output Parameters:
300: + gidx - global index of (net,idx)
301: . svtype - see petsc/private/dmnetworkimpl.h
302: - svtx_idx - ordering in the svtx array
303: */
304: static inline PetscErrorCode VtxGetInfo(PetscInt Nsvtx,SVtx *svtx,PetscInt net,PetscInt idx,PetscInt *gidx,SVtxType *svtype,PetscInt *svtx_idx)
305: {
306: PetscInt i,j,*svto,g_idx;
307: SVtxType vtype;
309: if (!Nsvtx) return 0;
311: g_idx = -1;
312: vtype = SVNONE;
314: for (i=0; i<Nsvtx; i++) {
315: if (net == svtx[i].sv[0] && idx == svtx[i].sv[1]) {
316: g_idx = svtx[i].gidx;
317: vtype = SVFROM;
318: } else { /* loop over svtx[i].n */
319: for (j=1; j<svtx[i].n; j++) {
320: svto = svtx[i].sv + 2*j;
321: if (net == svto[0] && idx == svto[1]) {
322: /* input vertex net.idx is a shared to_vertex, output its global index and its svtype */
323: g_idx = svtx[i].gidx; /* output gidx for to_vertex */
324: vtype = SVTO;
325: }
326: }
327: }
328: if (vtype != SVNONE) break;
329: }
330: if (gidx) *gidx = g_idx;
331: if (svtype) *svtype = vtype;
332: if (svtx_idx) *svtx_idx = i;
333: return 0;
334: }
336: /*
337: TableAddSVtx - Add a new shared vertice from sedgelist[k] to a ctable svta
339: Input: network, sedgelist, k, svta
340: Output: svta, tdata, ta2sv
341: */
342: static inline PetscErrorCode TableAddSVtx(DM_Network *network,PetscInt *sedgelist,PetscInt k,PetscTable svta,PetscInt* tdata,PetscInt *ta2sv)
343: {
344: PetscInt net,idx,gidx;
346: net = sedgelist[k];
347: idx = sedgelist[k+1];
348: gidx = network->subnet[net].vStart + idx;
349: PetscTableAdd(svta,gidx+1,*tdata+1,INSERT_VALUES);
351: ta2sv[*tdata] = k; /* maps tdata to index of sedgelist */
352: (*tdata)++;
353: return 0;
354: }
356: /*
357: SharedVtxCreate - Create an array of global shared vertices. See SVtx and SVtxType in dmnetworkimpl.h
359: Input: dm, Nsedgelist, sedgelist
361: Note: Output svtx is organized as
362: sv(net[0],idx[0]) --> sv(net[1],idx[1])
363: --> sv(net[1],idx[1])
364: ...
365: --> sv(net[n-1],idx[n-1])
366: and net[0] < net[1] < ... < net[n-1]
367: where sv[0] has SVFROM type, sv[i], i>0, has SVTO type.
368: */
369: static PetscErrorCode SharedVtxCreate(DM dm,PetscInt Nsedgelist,PetscInt *sedgelist)
370: {
371: SVtx *svtx = NULL;
372: PetscInt *sv,k,j,nsv,*tdata,**ta2sv;
373: PetscTable *svtas;
374: PetscInt gidx,net,idx,i,nta,ita,idx_from,idx_to,n;
375: DM_Network *network = (DM_Network*)dm->data;
376: PetscTablePosition ppos;
378: /* (1) Crete an array of ctables svtas to map (net,idx) -> gidx; a svtas[] for a shared/merged vertex */
379: PetscCalloc3(Nsedgelist,&svtas,Nsedgelist,&tdata,2*Nsedgelist,&ta2sv);
381: k = 0; /* sedgelist vertex counter j = 4*k */
382: nta = 0; /* num of svta tables created */
384: /* for j=0 */
385: PetscTableCreate(2*Nsedgelist,network->NVertices+1,&svtas[nta]);
386: PetscMalloc1(2*Nsedgelist,&ta2sv[nta]);
388: TableAddSVtx(network,sedgelist,k,svtas[nta],&tdata[nta],ta2sv[nta]);
389: TableAddSVtx(network,sedgelist,k+2,svtas[nta],&tdata[nta],ta2sv[nta]);
390: nta++; k += 4;
392: for (j = 1; j < Nsedgelist; j++) { /* j: sedgelist counter */
393: for (ita = 0; ita < nta; ita++) {
394: /* vfrom */
395: net = sedgelist[k]; idx = sedgelist[k+1];
396: gidx = network->subnet[net].vStart + idx; /* global index of the vertex net.idx before merging shared vertices */
397: PetscTableFind(svtas[ita],gidx+1,&idx_from);
399: /* vto */
400: net = sedgelist[k+2]; idx = sedgelist[k+3];
401: gidx = network->subnet[net].vStart + idx;
402: PetscTableFind(svtas[ita],gidx+1,&idx_to);
404: if (idx_from || idx_to) { /* vfrom or vto is on table svtas[ita] */
405: idx_from--; idx_to--;
406: if (idx_from < 0) { /* vto is on svtas[ita] */
407: TableAddSVtx(network,sedgelist,k,svtas[ita],&tdata[ita],ta2sv[ita]);
408: break;
409: } else if (idx_to < 0) {
410: TableAddSVtx(network,sedgelist,k+2,svtas[ita],&tdata[ita],ta2sv[ita]);
411: break;
412: }
413: }
414: }
416: if (ita == nta) {
417: PetscTableCreate(2*Nsedgelist,network->NVertices+1,&svtas[nta]);
418: PetscMalloc1(2*Nsedgelist, &ta2sv[nta]);
420: TableAddSVtx(network,sedgelist,k,svtas[nta],&tdata[nta],ta2sv[nta]);
421: TableAddSVtx(network,sedgelist,k+2,svtas[nta],&tdata[nta],ta2sv[nta]);
422: nta++;
423: }
424: k += 4;
425: }
427: /* (2) Create svtable for query shared vertices using gidx */
428: PetscTableCreate(nta,network->NVertices+1,&network->svtable);
430: /* (3) Construct svtx from svtas
431: svtx: array of SVtx: sv[0]=(net[0],idx[0]) to vertices sv[k], k=1,...,n-1;
432: net[k], k=0, ...,n-1, are in ascending order */
433: PetscMalloc1(nta,&svtx);
434: for (nsv = 0; nsv < nta; nsv++) {
435: /* for a single svtx, put shared vertices in ascending order of gidx */
436: PetscTableGetCount(svtas[nsv],&n);
437: PetscCalloc1(2*n,&sv);
438: svtx[nsv].sv = sv;
439: svtx[nsv].n = n;
440: svtx[nsv].gidx = network->NVertices; /* initialization */
442: PetscTableGetHeadPosition(svtas[nsv],&ppos);
443: for (k=0; k<n; k++) { /* gidx is sorted in ascending order */
444: PetscTableGetNext(svtas[nsv],&ppos,&gidx,&i);
445: gidx--; i--;
447: if (svtx[nsv].gidx > gidx) svtx[nsv].gidx = gidx; /*svtx[nsv].gidx = min(gidx) */
449: j = ta2sv[nsv][i]; /* maps i to index of sedgelist */
450: sv[2*k] = sedgelist[j]; /* subnet number */
451: sv[2*k+1] = sedgelist[j+1]; /* index on the subnet */
452: }
454: /* Setup svtable for query shared vertices */
455: PetscTableAdd(network->svtable,svtx[nsv].gidx+1,nsv+1,INSERT_VALUES);
456: }
458: for (j=0; j<nta; j++) {
459: PetscTableDestroy(&svtas[j]);
460: PetscFree(ta2sv[j]);
461: }
462: PetscFree3(svtas,tdata,ta2sv);
464: network->Nsvtx = nta;
465: network->svtx = svtx;
466: return 0;
467: }
469: /*
470: GetEdgelist_Coupling - Get an integrated edgelist for dmplex from user-provided subnet[].edgelist when subnets are coupled by shared vertices
472: Input Parameters:
473: . dm - the dmnetwork object
475: Output Parameters:
476: + edges - the integrated edgelist for dmplex
477: - nmerged_ptr - num of vertices being merged
478: */
479: static PetscErrorCode GetEdgelist_Coupling(DM dm,PetscInt *edges,PetscInt *nmerged_ptr)
480: {
481: MPI_Comm comm;
482: PetscMPIInt size,rank,*recvcounts=NULL,*displs=NULL;
483: DM_Network *network = (DM_Network*)dm->data;
484: PetscInt i,j,ctr,np;
485: PetscInt *vidxlTog,Nsv,Nsubnet=network->Nsubnet;
486: PetscInt *sedgelist=network->sedgelist;
487: PetscInt net,idx,gidx,nmerged,*vrange,gidx_from,net_from,sv_idx;
488: SVtxType svtype = SVNONE;
489: SVtx *svtx;
491: PetscObjectGetComm((PetscObject)dm,&comm);
492: MPI_Comm_rank(comm,&rank);
493: MPI_Comm_size(comm,&size);
495: /* (1) Create global svtx[] from sedgelist */
496: /* --------------------------------------- */
497: SharedVtxCreate(dm,network->Nsvtx,sedgelist);
498: Nsv = network->Nsvtx;
499: svtx = network->svtx;
501: /* (2) Merge shared vto vertices to their vfrom vertex with same global vetex index (gidx) */
502: /* --------------------------------------------------------------------------------------- */
503: /* (2.1) compute vrage[rank]: global index of 1st local vertex in proc[rank] */
504: PetscMalloc4(size+1,&vrange,size,&displs,size,&recvcounts,network->nVertices,&vidxlTog);
505: for (i=0; i<size; i++) {displs[i] = i; recvcounts[i] = 1;}
507: vrange[0] = 0;
508: MPI_Allgatherv(&network->nVertices,1,MPIU_INT,vrange+1,recvcounts,displs,MPIU_INT,comm);
509: for (i=2; i<size+1; i++) vrange[i] += vrange[i-1];
511: /* (2.2) Create vidxlTog: maps UN-MERGED local vertex index i to global index gidx (plex, excluding ghost vertices) */
512: i = 0; gidx = 0;
513: nmerged = 0; /* local num of merged vertices */
514: network->nsvtx = 0; /* local num of SVtx structs, including ghosts */
515: for (net=0; net<Nsubnet; net++) {
516: for (idx=0; idx<network->subnet[net].Nvtx; idx++) { /* Note: global subnet[net].Nvtx */
517: VtxGetInfo(Nsv,svtx,net,idx,&gidx_from,&svtype,&sv_idx);
518: if (svtype == SVTO) {
519: if (network->subnet[net].nvtx) {/* this proc owns sv_to */
520: net_from = svtx[sv_idx].sv[0]; /* subnet number of its shared vertex */
521: if (network->subnet[net_from].nvtx == 0) {
522: /* this proc does not own v_from, thus a ghost local vertex */
523: network->nsvtx++;
524: }
525: vidxlTog[i++] = gidx_from; /* gidx before merging! Bug??? */
526: nmerged++; /* a shared vertex -- merged */
527: }
528: } else {
529: if (svtype == SVFROM && network->subnet[net].nvtx) {
530: /* this proc owns this v_from, a new local shared vertex */
531: network->nsvtx++;
532: }
533: if (network->subnet[net].nvtx) vidxlTog[i++] = gidx;
534: gidx++;
535: }
536: }
537: }
538: #if defined(PETSC_USE_DEBUG)
540: #endif
542: /* (2.3) Shared vertices in the subnetworks are merged, update global NVertices: np = sum(local nmerged) */
543: MPI_Allreduce(&nmerged,&np,1,MPIU_INT,MPI_SUM,comm);
544: network->NVertices -= np;
546: ctr = 0;
547: for (net=0; net<Nsubnet; net++) {
548: for (j = 0; j < network->subnet[net].nedge; j++) {
549: /* vfrom: */
550: i = network->subnet[net].edgelist[2*j] + (network->subnet[net].vStart - vrange[rank]);
551: edges[2*ctr] = vidxlTog[i];
553: /* vto */
554: i = network->subnet[net].edgelist[2*j+1] + (network->subnet[net].vStart - vrange[rank]);
555: edges[2*ctr+1] = vidxlTog[i];
556: ctr++;
557: }
558: }
559: PetscFree4(vrange,displs,recvcounts,vidxlTog);
560: PetscFree(sedgelist); /* created in DMNetworkAddSharedVertices() */
562: *nmerged_ptr = nmerged;
563: return 0;
564: }
566: /*@
567: DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
569: Not Collective
571: Input Parameters:
572: . dm - the dmnetwork object
574: Notes:
575: This routine should be called after the network sizes and edgelists have been provided. It creates
576: the bare layout of the network and sets up the network to begin insertion of components.
578: All the components should be registered before calling this routine.
580: Level: beginner
582: .seealso: DMNetworkSetNumSubNetworks(), DMNetworkAddSubnetwork()
583: @*/
584: PetscErrorCode DMNetworkLayoutSetUp(DM dm)
585: {
586: DM_Network *network = (DM_Network*)dm->data;
587: PetscInt i,j,ctr,Nsubnet=network->Nsubnet,*eowners,np,*edges,*subnetvtx,*subnetedge,e,v,vfrom,vto,net;
588: const PetscInt *cone;
589: MPI_Comm comm;
590: PetscMPIInt size,rank;
591: PetscSection sectiong;
592: PetscInt nmerged=0;
596: /* This implementation requires user input each subnet by a single processor when Nsubnet>1, thus subnet[net].nvtx=subnet[net].Nvtx when net>0 */
597: for (net=1; net<Nsubnet; net++) {
599: }
601: PetscObjectGetComm((PetscObject)dm,&comm);
602: MPI_Comm_rank(comm,&rank);
603: MPI_Comm_size(comm,&size);
605: /* Create LOCAL edgelist in global vertex ordering for the network by concatenating local input edgelists of the subnetworks */
606: PetscCalloc2(2*network->nEdges,&edges,size+1,&eowners);
608: if (network->Nsvtx) { /* subnetworks are coupled via shared vertices */
609: GetEdgelist_Coupling(dm,edges,&nmerged);
610: } else { /* subnetworks are not coupled */
611: /* Create a 0-size svtable for query shared vertices */
612: PetscTableCreate(0,network->NVertices+1,&network->svtable);
613: ctr = 0;
614: for (i=0; i < Nsubnet; i++) {
615: for (j = 0; j < network->subnet[i].nedge; j++) {
616: edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
617: edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
618: ctr++;
619: }
620: }
621: }
623: /* Create network->plex; One dimensional network, numCorners=2 */
624: DMCreate(comm,&network->plex);
625: DMSetType(network->plex,DMPLEX);
626: DMSetDimension(network->plex,1);
628: if (size == 1) {
629: DMPlexBuildFromCellList(network->plex,network->nEdges,PETSC_DECIDE,2,edges);
630: } else {
631: DMPlexBuildFromCellListParallel(network->plex,network->nEdges,PETSC_DECIDE,PETSC_DECIDE,2,edges,NULL, NULL);
632: }
634: DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);
635: DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);
636: DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);
638: PetscSectionCreate(comm,&network->DataSection);
639: PetscSectionCreate(comm,&network->DofSection);
640: PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);
641: PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);
643: np = network->pEnd - network->pStart;
644: PetscCalloc2(np,&network->header,np,&network->cvalue);
645: for (i=0; i < np; i++) {
646: network->header[i].maxcomps = 1;
647: SetUpNetworkHeaderComponentValue(dm,&network->header[i],&network->cvalue[i]);
648: }
650: /* Create edge and vertex arrays for the subnetworks
651: This implementation assumes that DMNetwork reads
652: (1) a single subnetwork in parallel; or
653: (2) n subnetworks using n processors, one subnetwork/processor.
654: */
655: PetscCalloc2(network->nEdges,&subnetedge,network->nVertices+network->nsvtx,&subnetvtx); /* Maps local edge/vertex to local subnetwork's edge/vertex */
656: network->subnetedge = subnetedge;
657: network->subnetvtx = subnetvtx;
658: for (j=0; j < Nsubnet; j++) {
659: network->subnet[j].edges = subnetedge;
660: subnetedge += network->subnet[j].nedge;
662: network->subnet[j].vertices = subnetvtx;
663: subnetvtx += network->subnet[j].nvtx;
664: }
665: network->svertices = subnetvtx;
667: /* Get edge ownership */
668: np = network->eEnd - network->eStart;
669: MPI_Allgather(&np,1,MPIU_INT,eowners+1,1,MPIU_INT,comm);
670: eowners[0] = 0;
671: for (i=2; i<=size; i++) eowners[i] += eowners[i-1];
673: /* Setup local edge and vertex arrays for subnetworks */
674: e = 0;
675: for (i=0; i < Nsubnet; i++) {
676: ctr = 0;
677: for (j = 0; j < network->subnet[i].nedge; j++) {
678: /* edge e */
679: network->header[e].index = e + eowners[rank]; /* Global edge index */
680: network->header[e].subnetid = i;
681: network->subnet[i].edges[j] = e;
683: network->header[e].ndata = 0;
684: network->header[e].offset[0] = 0;
685: network->header[e].offsetvarrel[0] = 0;
686: PetscSectionAddDof(network->DataSection,e,network->header[e].hsize);
688: /* connected vertices */
689: DMPlexGetCone(network->plex,e,&cone);
691: /* vertex cone[0] */
692: v = cone[0];
693: network->header[v].index = edges[2*e]; /* Global vertex index */
694: network->header[v].subnetid = i; /* Subnetwork id */
695: if (Nsubnet == 1) {
696: network->subnet[i].vertices[v - network->vStart] = v; /* user's subnet[].idx = petsc's v */
697: } else {
698: vfrom = network->subnet[i].edgelist[2*ctr]; /* =subnet[i].idx, Global index! */
699: network->subnet[i].vertices[vfrom] = v; /* user's subnet[].dix = petsc's v */
700: }
702: /* vertex cone[1] */
703: v = cone[1];
704: network->header[v].index = edges[2*e+1]; /* Global vertex index */
705: network->header[v].subnetid = i; /* Subnetwork id */
706: if (Nsubnet == 1) {
707: network->subnet[i].vertices[v - network->vStart] = v; /* user's subnet[].idx = petsc's v */
708: } else {
709: vto = network->subnet[i].edgelist[2*ctr+1]; /* =subnet[i].idx, Global index! */
710: network->subnet[i].vertices[vto] = v; /* user's subnet[].dix = petsc's v */
711: }
713: e++; ctr++;
714: }
715: }
716: PetscFree2(edges,eowners);
718: /* Set local vertex array for the subnetworks */
719: j = 0;
720: for (v = network->vStart; v < network->vEnd; v++) {
721: network->header[v].ndata = 0;
722: network->header[v].offset[0] = 0;
723: network->header[v].offsetvarrel[0] = 0;
724: PetscSectionAddDof(network->DataSection,v,network->header[v].hsize);
726: /* local shared vertex */
727: PetscTableFind(network->svtable,network->header[v].index+1,&i);
728: if (i) network->svertices[j++] = v;
729: }
731: /* Create a global section to be used by DMNetworkIsGhostVertex() which is a non-collective routine */
732: /* see snes_tutorials_network-ex1_4 */
733: DMGetGlobalSection(network->plex,§iong);
734: return 0;
735: }
737: /*@C
738: DMNetworkGetSubnetwork - Returns the information about a requested subnetwork
740: Not collective
742: Input Parameters:
743: + dm - the DM object
744: - netnum - the global index of the subnetwork
746: Output Parameters:
747: + nv - number of vertices (local)
748: . ne - number of edges (local)
749: . vtx - local vertices of the subnetwork
750: - edge - local edges of the subnetwork
752: Notes:
753: Cannot call this routine before DMNetworkLayoutSetup()
755: The local vertices returned on each rank are determined by DMNetwork. The user does not have any control over what vertices are local.
757: Level: intermediate
759: .seealso: DMNetworkCreate(), DMNetworkAddSubnetwork(), DMNetworkLayoutSetUp()
760: @*/
761: PetscErrorCode DMNetworkGetSubnetwork(DM dm,PetscInt netnum,PetscInt *nv,PetscInt *ne,const PetscInt **vtx,const PetscInt **edge)
762: {
763: DM_Network *network = (DM_Network*)dm->data;
766: if (nv) *nv = network->subnet[netnum].nvtx;
767: if (ne) *ne = network->subnet[netnum].nedge;
768: if (vtx) *vtx = network->subnet[netnum].vertices;
769: if (edge) *edge = network->subnet[netnum].edges;
770: return 0;
771: }
773: /*@
774: DMNetworkAddSharedVertices - Add shared vertices that connect two given subnetworks
776: Collective on dm
778: Input Parameters:
779: + dm - the dm object
780: . anetnum - first subnetwork global numbering returned by DMNetworkAddSubnetwork()
781: . bnetnum - second subnetwork global numbering returned by DMNetworkAddSubnetwork()
782: . nsvtx - number of vertices that are shared by the two subnetworks
783: . asvtx - vertex index in the first subnetwork
784: - bsvtx - vertex index in the second subnetwork
786: Level: beginner
788: .seealso: DMNetworkCreate(), DMNetworkAddSubnetwork(), DMNetworkGetSharedVertices()
789: @*/
790: PetscErrorCode DMNetworkAddSharedVertices(DM dm,PetscInt anetnum,PetscInt bnetnum,PetscInt nsvtx,PetscInt asvtx[],PetscInt bsvtx[])
791: {
792: DM_Network *network = (DM_Network*)dm->data;
793: PetscInt i,nsubnet = network->Nsubnet,*sedgelist,Nsvtx=network->Nsvtx;
797: if (!Nsvtx) {
798: /* allocate network->sedgelist to hold at most 2*nsubnet pairs of shared vertices */
799: PetscMalloc1(2*4*nsubnet,&network->sedgelist);
800: }
802: sedgelist = network->sedgelist;
803: for (i=0; i<nsvtx; i++) {
804: sedgelist[4*Nsvtx] = anetnum; sedgelist[4*Nsvtx+1] = asvtx[i];
805: sedgelist[4*Nsvtx+2] = bnetnum; sedgelist[4*Nsvtx+3] = bsvtx[i];
806: Nsvtx++;
807: }
809: network->Nsvtx = Nsvtx;
810: return 0;
811: }
813: /*@C
814: DMNetworkGetSharedVertices - Returns the info for the shared vertices
816: Not collective
818: Input Parameter:
819: . dm - the DM object
821: Output Parameters:
822: + nsv - number of local shared vertices
823: - svtx - local shared vertices
825: Notes:
826: Cannot call this routine before DMNetworkLayoutSetup()
828: Level: intermediate
830: .seealso: DMNetworkGetSubnetwork(), DMNetworkLayoutSetUp(), DMNetworkAddSharedVertices()
831: @*/
832: PetscErrorCode DMNetworkGetSharedVertices(DM dm,PetscInt *nsv,const PetscInt **svtx)
833: {
834: DM_Network *net = (DM_Network*)dm->data;
837: if (nsv) *nsv = net->nsvtx;
838: if (svtx) *svtx = net->svertices;
839: return 0;
840: }
842: /*@C
843: DMNetworkRegisterComponent - Registers the network component
845: Logically collective on dm
847: Input Parameters:
848: + dm - the network object
849: . name - the component name
850: - size - the storage size in bytes for this component data
852: Output Parameters:
853: . key - an integer key that defines the component
855: Notes
856: This routine should be called by all processors before calling DMNetworkLayoutSetup().
858: Level: beginner
860: .seealso: DMNetworkCreate(), DMNetworkLayoutSetUp()
861: @*/
862: PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,size_t size,PetscInt *key)
863: {
864: DM_Network *network = (DM_Network*) dm->data;
865: DMNetworkComponent *component=NULL,*newcomponent=NULL;
866: PetscBool flg=PETSC_FALSE;
867: PetscInt i;
869: if (!network->component) {
870: PetscCalloc1(network->max_comps_registered,&network->component);
871: }
873: for (i=0; i < network->ncomponent; i++) {
874: PetscStrcmp(network->component[i].name,name,&flg);
875: if (flg) {
876: *key = i;
877: return 0;
878: }
879: }
881: if (network->ncomponent == network->max_comps_registered) {
882: /* Reached max allowed so resize component */
883: network->max_comps_registered += 2;
884: PetscCalloc1(network->max_comps_registered,&newcomponent);
885: /* Copy over the previous component info */
886: for (i=0; i < network->ncomponent; i++) {
887: PetscStrcpy(newcomponent[i].name,network->component[i].name);
888: newcomponent[i].size = network->component[i].size;
889: }
890: /* Free old one */
891: PetscFree(network->component);
892: /* Update pointer */
893: network->component = newcomponent;
894: }
896: component = &network->component[network->ncomponent];
898: PetscStrcpy(component->name,name);
899: component->size = size/sizeof(DMNetworkComponentGenericDataType);
900: *key = network->ncomponent;
901: network->ncomponent++;
902: return 0;
903: }
905: /*@
906: DMNetworkGetVertexRange - Get the bounds [start, end) for the local vertices
908: Not Collective
910: Input Parameter:
911: . dm - the DMNetwork object
913: Output Parameters:
914: + vStart - the first vertex point
915: - vEnd - one beyond the last vertex point
917: Level: beginner
919: .seealso: DMNetworkGetEdgeRange()
920: @*/
921: PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
922: {
923: DM_Network *network = (DM_Network*)dm->data;
925: if (vStart) *vStart = network->vStart;
926: if (vEnd) *vEnd = network->vEnd;
927: return 0;
928: }
930: /*@
931: DMNetworkGetEdgeRange - Get the bounds [start, end) for the local edges
933: Not Collective
935: Input Parameter:
936: . dm - the DMNetwork object
938: Output Parameters:
939: + eStart - The first edge point
940: - eEnd - One beyond the last edge point
942: Level: beginner
944: .seealso: DMNetworkGetVertexRange()
945: @*/
946: PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
947: {
948: DM_Network *network = (DM_Network*)dm->data;
951: if (eStart) *eStart = network->eStart;
952: if (eEnd) *eEnd = network->eEnd;
953: return 0;
954: }
956: static PetscErrorCode DMNetworkGetIndex(DM dm,PetscInt p,PetscInt *index)
957: {
958: DM_Network *network = (DM_Network*)dm->data;
960: if (network->header) {
961: *index = network->header[p].index;
962: } else {
963: PetscInt offsetp;
964: DMNetworkComponentHeader header;
966: PetscSectionGetOffset(network->DataSection,p,&offsetp);
967: header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
968: *index = header->index;
969: }
970: return 0;
971: }
973: /*@
974: DMNetworkGetGlobalEdgeIndex - Get the global numbering for the edge on the network
976: Not Collective
978: Input Parameters:
979: + dm - DMNetwork object
980: - p - edge point
982: Output Parameters:
983: . index - the global numbering for the edge
985: Level: intermediate
987: .seealso: DMNetworkGetGlobalVertexIndex()
988: @*/
989: PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
990: {
991: DMNetworkGetIndex(dm,p,index);
992: return 0;
993: }
995: /*@
996: DMNetworkGetGlobalVertexIndex - Get the global numbering for the vertex on the network
998: Not Collective
1000: Input Parameters:
1001: + dm - DMNetwork object
1002: - p - vertex point
1004: Output Parameters:
1005: . index - the global numbering for the vertex
1007: Level: intermediate
1009: .seealso: DMNetworkGetGlobalEdgeIndex(), DMNetworkGetLocalVertexIndex()
1010: @*/
1011: PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
1012: {
1013: DMNetworkGetIndex(dm,p,index);
1014: return 0;
1015: }
1017: /*@
1018: DMNetworkGetNumComponents - Get the number of components at a vertex/edge
1020: Not Collective
1022: Input Parameters:
1023: + dm - the DMNetwork object
1024: - p - vertex/edge point
1026: Output Parameters:
1027: . numcomponents - Number of components at the vertex/edge
1029: Level: beginner
1031: .seealso: DMNetworkRegisterComponent(), DMNetworkAddComponent()
1032: @*/
1033: PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
1034: {
1035: PetscInt offset;
1036: DM_Network *network = (DM_Network*)dm->data;
1038: PetscSectionGetOffset(network->DataSection,p,&offset);
1039: *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
1040: return 0;
1041: }
1043: /*@
1044: DMNetworkGetLocalVecOffset - Get the offset for accessing the variables associated with a component at the given vertex/edge from the local vector
1046: Not Collective
1048: Input Parameters:
1049: + dm - the DMNetwork object
1050: . p - the edge or vertex point
1051: - compnum - component number; use ALL_COMPONENTS if no specific component is requested
1053: Output Parameters:
1054: . offset - the local offset
1056: Notes:
1057: These offsets can be passed to MatSetValuesLocal() for matrices obtained with DMCreateMatrix().
1059: For vectors obtained with DMCreateLocalVector() the offsets can be used with VecSetValues().
1061: For vectors obtained with DMCreateLocalVector() and the array obtained with VecGetArray(vec,&array) you can access or set
1062: the vector values with array[offset].
1064: For vectors obtained with DMCreateGlobalVector() the offsets can be used with VecSetValuesLocal().
1066: Level: intermediate
1068: .seealso: DMGetLocalVector(), DMNetworkGetComponent(), DMNetworkGetGlobalVecOffset(), DMCreateGlobalVector(), VecGetArray(), VecSetValuesLocal(), MatSetValuesLocal()
1069: @*/
1070: PetscErrorCode DMNetworkGetLocalVecOffset(DM dm,PetscInt p,PetscInt compnum,PetscInt *offset)
1071: {
1072: DM_Network *network = (DM_Network*)dm->data;
1073: PetscInt offsetp,offsetd;
1074: DMNetworkComponentHeader header;
1076: PetscSectionGetOffset(network->plex->localSection,p,&offsetp);
1077: if (compnum == ALL_COMPONENTS) {
1078: *offset = offsetp;
1079: return 0;
1080: }
1082: PetscSectionGetOffset(network->DataSection,p,&offsetd);
1083: header = (DMNetworkComponentHeader)(network->componentdataarray+offsetd);
1084: *offset = offsetp + header->offsetvarrel[compnum];
1085: return 0;
1086: }
1088: /*@
1089: DMNetworkGetGlobalVecOffset - Get the global offset for accessing the variables associated with a component for the given vertex/edge from the global vector
1091: Not Collective
1093: Input Parameters:
1094: + dm - the DMNetwork object
1095: . p - the edge or vertex point
1096: - compnum - component number; use ALL_COMPONENTS if no specific component is requested
1098: Output Parameters:
1099: . offsetg - the global offset
1101: Notes:
1102: These offsets can be passed to MatSetValues() for matrices obtained with DMCreateMatrix().
1104: For vectors obtained with DMCreateGlobalVector() the offsets can be used with VecSetValues().
1106: For vectors obtained with DMCreateGlobalVector() and the array obtained with VecGetArray(vec,&array) you can access or set
1107: the vector values with array[offset - rstart] where restart is obtained with VecGetOwnershipRange(v,&rstart,NULL);
1109: Level: intermediate
1111: .seealso: DMNetworkGetLocalVecOffset(), DMGetGlobalVector(), DMNetworkGetComponent(), DMCreateGlobalVector(), VecGetArray(), VecSetValues(), MatSetValues()
1112: @*/
1113: PetscErrorCode DMNetworkGetGlobalVecOffset(DM dm,PetscInt p,PetscInt compnum,PetscInt *offsetg)
1114: {
1115: DM_Network *network = (DM_Network*)dm->data;
1116: PetscInt offsetp,offsetd;
1117: DMNetworkComponentHeader header;
1119: PetscSectionGetOffset(network->plex->globalSection,p,&offsetp);
1120: if (offsetp < 0) offsetp = -(offsetp + 1); /* Convert to actual global offset for ghost vertex */
1122: if (compnum == ALL_COMPONENTS) {
1123: *offsetg = offsetp;
1124: return 0;
1125: }
1126: PetscSectionGetOffset(network->DataSection,p,&offsetd);
1127: header = (DMNetworkComponentHeader)(network->componentdataarray+offsetd);
1128: *offsetg = offsetp + header->offsetvarrel[compnum];
1129: return 0;
1130: }
1132: /*@
1133: DMNetworkGetEdgeOffset - Get the offset for accessing the variables associated with the given edge from the local subvector
1135: Not Collective
1137: Input Parameters:
1138: + dm - the DMNetwork object
1139: - p - the edge point
1141: Output Parameters:
1142: . offset - the offset
1144: Level: intermediate
1146: .seealso: DMNetworkGetLocalVecOffset(), DMGetLocalVector()
1147: @*/
1148: PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
1149: {
1150: DM_Network *network = (DM_Network*)dm->data;
1152: PetscSectionGetOffset(network->edge.DofSection,p,offset);
1153: return 0;
1154: }
1156: /*@
1157: DMNetworkGetVertexOffset - Get the offset for accessing the variables associated with the given vertex from the local subvector
1159: Not Collective
1161: Input Parameters:
1162: + dm - the DMNetwork object
1163: - p - the vertex point
1165: Output Parameters:
1166: . offset - the offset
1168: Level: intermediate
1170: .seealso: DMNetworkGetEdgeOffset(), DMGetLocalVector()
1171: @*/
1172: PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
1173: {
1174: DM_Network *network = (DM_Network*)dm->data;
1176: p -= network->vStart;
1177: PetscSectionGetOffset(network->vertex.DofSection,p,offset);
1178: return 0;
1179: }
1181: /*@
1182: DMNetworkAddComponent - Adds a network component and number of variables at the given point (vertex/edge)
1184: Collective on dm
1186: Input Parameters:
1187: + dm - the DMNetwork
1188: . p - the vertex/edge point. These points are local indices provided by DMNetworkGetSubnetwork()
1189: . componentkey - component key returned while registering the component with DMNetworkRegisterComponent()
1190: . 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
1191: free this space until after DMSetUp() is called.
1192: - 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
1194: Notes:
1195: 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.
1197: DMNetworkLayoutSetUp() must be called before this routine.
1199: Developer Notes:
1200: 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.
1201: Level: beginner
1203: .seealso: DMNetworkGetComponent(), DMNetworkGetSubnetwork(), DMNetworkIsGhostVertex(), DMNetworkLayoutSetUp()
1204: @*/
1205: PetscErrorCode DMNetworkAddComponent(DM dm,PetscInt p,PetscInt componentkey,void* compvalue,PetscInt nvar)
1206: {
1207: DM_Network *network = (DM_Network*)dm->data;
1208: DMNetworkComponent *component = &network->component[componentkey];
1209: DMNetworkComponentHeader header;
1210: DMNetworkComponentValue cvalue;
1211: PetscInt compnum;
1212: PetscInt *compsize,*compkey,*compoffset,*compnvar,*compoffsetvarrel;
1213: void* *compdata;
1217: /* The owning rank and all ghost ranks add nvar */
1218: PetscSectionAddDof(network->DofSection,p,nvar);
1220: /* The owning rank and all ghost ranks add a component, including compvalue=NULL */
1221: header = &network->header[p];
1222: cvalue = &network->cvalue[p];
1223: if (header->ndata == header->maxcomps) {
1224: PetscInt additional_size;
1226: /* Reached limit so resize header component arrays */
1227: header->maxcomps += 2;
1229: /* Allocate arrays for component information and value */
1230: PetscCalloc5(header->maxcomps,&compsize,header->maxcomps,&compkey,header->maxcomps,&compoffset,header->maxcomps,&compnvar,header->maxcomps,&compoffsetvarrel);
1231: PetscMalloc1(header->maxcomps,&compdata);
1233: /* Recalculate header size */
1234: header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5*header->maxcomps*sizeof(PetscInt);
1236: header->hsize /= sizeof(DMNetworkComponentGenericDataType);
1238: /* Copy over component info */
1239: PetscMemcpy(compsize,header->size,header->ndata*sizeof(PetscInt));
1240: PetscMemcpy(compkey,header->key,header->ndata*sizeof(PetscInt));
1241: PetscMemcpy(compoffset,header->offset,header->ndata*sizeof(PetscInt));
1242: PetscMemcpy(compnvar,header->nvar,header->ndata*sizeof(PetscInt));
1243: PetscMemcpy(compoffsetvarrel,header->offsetvarrel,header->ndata*sizeof(PetscInt));
1245: /* Copy over component data pointers */
1246: PetscMemcpy(compdata,cvalue->data,header->ndata*sizeof(void*));
1248: /* Free old arrays */
1249: PetscFree5(header->size,header->key,header->offset,header->nvar,header->offsetvarrel);
1250: PetscFree(cvalue->data);
1252: /* Update pointers */
1253: header->size = compsize;
1254: header->key = compkey;
1255: header->offset = compoffset;
1256: header->nvar = compnvar;
1257: header->offsetvarrel = compoffsetvarrel;
1259: cvalue->data = compdata;
1261: /* Update DataSection Dofs */
1262: /* 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 */
1263: additional_size = (5*(header->maxcomps - header->ndata)*sizeof(PetscInt))/sizeof(DMNetworkComponentGenericDataType);
1264: PetscSectionAddDof(network->DataSection,p,additional_size);
1265: }
1266: header = &network->header[p];
1267: cvalue = &network->cvalue[p];
1269: compnum = header->ndata;
1271: header->size[compnum] = component->size;
1272: PetscSectionAddDof(network->DataSection,p,component->size);
1273: header->key[compnum] = componentkey;
1274: if (compnum != 0) header->offset[compnum] = header->offset[compnum-1] + header->size[compnum-1];
1275: cvalue->data[compnum] = (void*)compvalue;
1277: /* variables */
1278: header->nvar[compnum] += nvar;
1279: if (compnum != 0) header->offsetvarrel[compnum] = header->offsetvarrel[compnum-1] + header->nvar[compnum-1];
1281: header->ndata++;
1282: return 0;
1283: }
1285: /*@
1286: DMNetworkGetComponent - Gets the component key, the component data, and the number of variables at a given network point
1288: Not Collective
1290: Input Parameters:
1291: + dm - the DMNetwork object
1292: . p - vertex/edge point
1293: - compnum - component number; use ALL_COMPONENTS if sum up all the components
1295: Output Parameters:
1296: + compkey - the key obtained when registering the component (use NULL if not required)
1297: . component - the component data (use NULL if not required)
1298: - nvar - number of variables (use NULL if not required)
1300: Level: beginner
1302: .seealso: DMNetworkAddComponent(), DMNetworkGetNumComponents()
1303: @*/
1304: PetscErrorCode DMNetworkGetComponent(DM dm,PetscInt p,PetscInt compnum,PetscInt *compkey,void **component,PetscInt *nvar)
1305: {
1306: DM_Network *network = (DM_Network*)dm->data;
1307: PetscInt offset = 0;
1308: DMNetworkComponentHeader header;
1310: if (compnum == ALL_COMPONENTS) {
1311: PetscSectionGetDof(network->DofSection,p,nvar);
1312: return 0;
1313: }
1315: PetscSectionGetOffset(network->DataSection,p,&offset);
1316: header = (DMNetworkComponentHeader)(network->componentdataarray+offset);
1318: if (compnum >= 0) {
1319: if (compkey) *compkey = header->key[compnum];
1320: if (component) {
1321: offset += header->hsize+header->offset[compnum];
1322: *component = network->componentdataarray+offset;
1323: }
1324: }
1326: if (nvar) *nvar = header->nvar[compnum];
1328: return 0;
1329: }
1331: /*
1332: Sets up the array that holds the data for all components and its associated section.
1333: 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.
1334: */
1335: PetscErrorCode DMNetworkComponentSetUp(DM dm)
1336: {
1337: DM_Network *network = (DM_Network*)dm->data;
1338: PetscInt arr_size,p,offset,offsetp,ncomp,i,*headerarr;
1339: DMNetworkComponentHeader header;
1340: DMNetworkComponentValue cvalue;
1341: DMNetworkComponentHeader headerinfo;
1342: DMNetworkComponentGenericDataType *componentdataarray;
1344: PetscSectionSetUp(network->DataSection);
1345: PetscSectionGetStorageSize(network->DataSection,&arr_size);
1346: /* arr_size+1 fixes pipeline test of opensolaris-misc for src/dm/tests/ex10.c -- Do not know why */
1347: PetscCalloc1(arr_size+1,&network->componentdataarray);
1349: componentdataarray = network->componentdataarray;
1350: for (p = network->pStart; p < network->pEnd; p++) {
1351: PetscSectionGetOffset(network->DataSection,p,&offsetp);
1352: /* Copy header */
1353: header = &network->header[p];
1354: headerinfo = (DMNetworkComponentHeader)(componentdataarray+offsetp);
1355: PetscMemcpy(headerinfo,header,sizeof(struct _p_DMNetworkComponentHeader));
1356: headerarr = (PetscInt*)(headerinfo+1);
1357: PetscMemcpy(headerarr,header->size,header->maxcomps*sizeof(PetscInt));
1358: headerarr += header->maxcomps;
1359: PetscMemcpy(headerarr,header->key,header->maxcomps*sizeof(PetscInt));
1360: headerarr += header->maxcomps;
1361: PetscMemcpy(headerarr,header->offset,header->maxcomps*sizeof(PetscInt));
1362: headerarr += header->maxcomps;
1363: PetscMemcpy(headerarr,header->nvar,header->maxcomps*sizeof(PetscInt));
1364: headerarr += header->maxcomps;
1365: PetscMemcpy(headerarr,header->offsetvarrel,header->maxcomps*sizeof(PetscInt));
1367: /* Copy data */
1368: cvalue = &network->cvalue[p];
1369: ncomp = header->ndata;
1371: for (i = 0; i < ncomp; i++) {
1372: offset = offsetp + header->hsize + header->offset[i];
1373: PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));
1374: }
1375: }
1376: return 0;
1377: }
1379: /* Sets up the section for dofs. This routine is called during DMSetUp() */
1380: static PetscErrorCode DMNetworkVariablesSetUp(DM dm)
1381: {
1382: DM_Network *network = (DM_Network*)dm->data;
1384: PetscSectionSetUp(network->DofSection);
1385: return 0;
1386: }
1388: /* Get a subsection from a range of points */
1389: static PetscErrorCode DMNetworkGetSubSection_private(PetscSection main,PetscInt pstart,PetscInt pend,PetscSection *subsection)
1390: {
1391: PetscInt i, nvar;
1393: PetscSectionCreate(PetscObjectComm((PetscObject)main), subsection);
1394: PetscSectionSetChart(*subsection, 0, pend - pstart);
1395: for (i = pstart; i < pend; i++) {
1396: PetscSectionGetDof(main,i,&nvar);
1397: PetscSectionSetDof(*subsection, i - pstart, nvar);
1398: }
1400: PetscSectionSetUp(*subsection);
1401: return 0;
1402: }
1404: /* Create a submap of points with a GlobalToLocal structure */
1405: static PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
1406: {
1407: PetscInt i, *subpoints;
1409: /* Create index sets to map from "points" to "subpoints" */
1410: PetscMalloc1(pend - pstart, &subpoints);
1411: for (i = pstart; i < pend; i++) {
1412: subpoints[i - pstart] = i;
1413: }
1414: ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);
1415: PetscFree(subpoints);
1416: return 0;
1417: }
1419: /*@
1420: DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute
1422: Collective on dm
1424: Input Parameters:
1425: . dm - the DMNetworkObject
1427: Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
1429: points = [0 1 2 3 4 5 6]
1431: 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]).
1433: With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
1435: Level: intermediate
1437: @*/
1438: PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
1439: {
1440: MPI_Comm comm;
1441: PetscMPIInt size;
1442: DM_Network *network = (DM_Network*)dm->data;
1444: PetscObjectGetComm((PetscObject)dm,&comm);
1445: MPI_Comm_size(comm, &size);
1447: /* Create maps for vertices and edges */
1448: DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);
1449: DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);
1451: /* Create local sub-sections */
1452: DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);
1453: DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);
1455: if (size > 1) {
1456: PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);
1458: PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);
1459: PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);
1460: PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);
1461: } else {
1462: /* create structures for vertex */
1463: PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);
1464: /* create structures for edge */
1465: PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);
1466: }
1468: /* Add viewers */
1469: PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");
1470: PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");
1471: PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");
1472: PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");
1473: return 0;
1474: }
1476: /*
1477: Setup a lookup btable for the input v's owning subnetworks
1478: - add all owing subnetworks that connect to this v to the btable
1479: vertex_subnetid = supportingedge_subnetid
1480: */
1481: static inline PetscErrorCode SetSubnetIdLookupBT(DM dm,PetscInt v,PetscInt Nsubnet,PetscBT btable)
1482: {
1483: PetscInt e,nedges,offset;
1484: const PetscInt *edges;
1485: DM_Network *newDMnetwork = (DM_Network*)dm->data;
1486: DMNetworkComponentHeader header;
1488: PetscBTMemzero(Nsubnet,btable);
1489: DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);
1490: for (e=0; e<nedges; e++) {
1491: PetscSectionGetOffset(newDMnetwork->DataSection,edges[e],&offset);
1492: header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1493: PetscBTSet(btable,header->subnetid);
1494: }
1495: return 0;
1496: }
1498: /*@
1499: DMNetworkDistribute - Distributes the network and moves associated component data
1501: Collective
1503: Input Parameters:
1504: + DM - the DMNetwork object
1505: - overlap - the overlap of partitions, 0 is the default
1507: Options Database Key:
1508: + -dmnetwork_view - Calls DMView() at the conclusion of DMSetUp()
1509: - -dmnetwork_view_distributed - Calls DMView() at the conclusion of DMNetworkDistribute()
1511: Notes:
1512: Distributes the network with <overlap>-overlapping partitioning of the edges.
1514: Level: intermediate
1516: .seealso: DMNetworkCreate()
1517: @*/
1518: PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
1519: {
1520: MPI_Comm comm;
1521: PetscMPIInt size;
1522: DM_Network *oldDMnetwork = (DM_Network*)((*dm)->data);
1523: DM_Network *newDMnetwork;
1524: PetscSF pointsf=NULL;
1525: DM newDM;
1526: PetscInt j,e,v,offset,*subnetvtx,*subnetedge,Nsubnet,gidx,svtx_idx,nv;
1527: PetscInt net,*sv;
1528: PetscBT btable;
1529: PetscPartitioner part;
1530: DMNetworkComponentHeader header;
1534: PetscObjectGetComm((PetscObject)*dm,&comm);
1535: MPI_Comm_size(comm, &size);
1536: if (size == 1) return 0;
1540: /* 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. */
1541: DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);
1542: newDMnetwork = (DM_Network*)newDM->data;
1543: newDMnetwork->max_comps_registered = oldDMnetwork->max_comps_registered;
1544: PetscMalloc1(newDMnetwork->max_comps_registered,&newDMnetwork->component);
1546: /* Enable runtime options for petscpartitioner */
1547: DMPlexGetPartitioner(oldDMnetwork->plex,&part);
1548: PetscPartitionerSetFromOptions(part);
1550: /* Distribute plex dm */
1551: DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);
1553: /* Distribute dof section */
1554: PetscSectionCreate(comm,&newDMnetwork->DofSection);
1555: PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);
1557: /* Distribute data and associated section */
1558: PetscSectionCreate(comm,&newDMnetwork->DataSection);
1559: DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);
1561: PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);
1562: DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);
1563: DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);
1564: newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
1565: newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
1566: newDMnetwork->NVertices = oldDMnetwork->NVertices;
1567: newDMnetwork->NEdges = oldDMnetwork->NEdges;
1568: newDMnetwork->svtable = oldDMnetwork->svtable; /* global table! */
1569: oldDMnetwork->svtable = NULL;
1571: /* Set Dof section as the section for dm */
1572: DMSetLocalSection(newDMnetwork->plex,newDMnetwork->DofSection);
1573: DMGetGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);
1575: /* Setup subnetwork info in the newDM */
1576: newDMnetwork->Nsubnet = oldDMnetwork->Nsubnet;
1577: newDMnetwork->Nsvtx = oldDMnetwork->Nsvtx;
1578: oldDMnetwork->Nsvtx = 0;
1579: newDMnetwork->svtx = oldDMnetwork->svtx; /* global vertices! */
1580: oldDMnetwork->svtx = NULL;
1581: PetscCalloc1(newDMnetwork->Nsubnet,&newDMnetwork->subnet);
1583: /* Copy over the global number of vertices and edges in each subnetwork.
1584: Note: these are calculated in DMNetworkLayoutSetUp()
1585: */
1586: Nsubnet = newDMnetwork->Nsubnet;
1587: for (j = 0; j < Nsubnet; j++) {
1588: newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx;
1589: newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge;
1590: }
1592: /* Count local nedges for subnetworks */
1593: for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++) {
1594: PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);
1595: header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1597: /* Update pointers */
1598: header->size = (PetscInt*)(header + 1);
1599: header->key = header->size + header->maxcomps;
1600: header->offset = header->key + header->maxcomps;
1601: header->nvar = header->offset + header->maxcomps;
1602: header->offsetvarrel = header->nvar + header->maxcomps;
1604: newDMnetwork->subnet[header->subnetid].nedge++;
1605: }
1607: /* Setup a btable to keep track subnetworks owned by this process at a shared vertex */
1608: if (newDMnetwork->Nsvtx) {
1609: PetscBTCreate(Nsubnet,&btable);
1610: }
1612: /* Count local nvtx for subnetworks */
1613: for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++) {
1614: PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);
1615: header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1617: /* Update pointers */
1618: header->size = (PetscInt*)(header + 1);
1619: header->key = header->size + header->maxcomps;
1620: header->offset = header->key + header->maxcomps;
1621: header->nvar = header->offset + header->maxcomps;
1622: header->offsetvarrel = header->nvar + header->maxcomps;
1624: /* shared vertices: use gidx=header->index to check if v is a shared vertex */
1625: gidx = header->index;
1626: PetscTableFind(newDMnetwork->svtable,gidx+1,&svtx_idx);
1627: svtx_idx--;
1629: if (svtx_idx < 0) { /* not a shared vertex */
1630: newDMnetwork->subnet[header->subnetid].nvtx++;
1631: } else { /* a shared vertex belongs to more than one subnetworks, it is being counted by multiple subnets */
1632: /* Setup a lookup btable for this v's owning subnetworks */
1633: SetSubnetIdLookupBT(newDM,v,Nsubnet,btable);
1635: for (j=0; j<newDMnetwork->svtx[svtx_idx].n; j++) {
1636: sv = newDMnetwork->svtx[svtx_idx].sv + 2*j;
1637: net = sv[0];
1638: if (PetscBTLookup(btable,net)) newDMnetwork->subnet[net].nvtx++; /* sv is on net owned by this proces */
1639: }
1640: }
1641: }
1643: /* Get total local nvtx for subnetworks */
1644: nv = 0;
1645: for (j=0; j<Nsubnet; j++) nv += newDMnetwork->subnet[j].nvtx;
1646: nv += newDMnetwork->Nsvtx;
1648: /* Now create the vertices and edge arrays for the subnetworks */
1649: PetscCalloc2(newDMnetwork->nEdges,&subnetedge,nv,&subnetvtx); /* Maps local vertex to local subnetwork's vertex */
1650: newDMnetwork->subnetedge = subnetedge;
1651: newDMnetwork->subnetvtx = subnetvtx;
1652: for (j=0; j < newDMnetwork->Nsubnet; j++) {
1653: newDMnetwork->subnet[j].edges = subnetedge;
1654: subnetedge += newDMnetwork->subnet[j].nedge;
1656: newDMnetwork->subnet[j].vertices = subnetvtx;
1657: subnetvtx += newDMnetwork->subnet[j].nvtx;
1659: /* 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. */
1660: newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0;
1661: }
1662: newDMnetwork->svertices = subnetvtx;
1664: /* Set the edges and vertices in each subnetwork */
1665: for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++) {
1666: PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);
1667: header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1668: newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++] = e;
1669: }
1671: nv = 0;
1672: for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++) {
1673: PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);
1674: header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1676: /* coupling vertices: use gidx = header->index to check if v is a coupling vertex */
1677: PetscTableFind(newDMnetwork->svtable,header->index+1,&svtx_idx);
1678: svtx_idx--;
1679: if (svtx_idx < 0) {
1680: newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++] = v;
1681: } else { /* a shared vertex */
1682: newDMnetwork->svertices[nv++] = v;
1684: /* Setup a lookup btable for this v's owning subnetworks */
1685: SetSubnetIdLookupBT(newDM,v,Nsubnet,btable);
1687: for (j=0; j<newDMnetwork->svtx[svtx_idx].n; j++) {
1688: sv = newDMnetwork->svtx[svtx_idx].sv + 2*j;
1689: net = sv[0];
1690: if (PetscBTLookup(btable,net))
1691: newDMnetwork->subnet[net].vertices[newDMnetwork->subnet[net].nvtx++] = v;
1692: }
1693: }
1694: }
1695: newDMnetwork->nsvtx = nv; /* num of local shared vertices */
1697: newDM->setupcalled = (*dm)->setupcalled;
1698: newDMnetwork->distributecalled = PETSC_TRUE;
1700: /* Free spaces */
1701: PetscSFDestroy(&pointsf);
1702: DMDestroy(dm);
1703: if (newDMnetwork->Nsvtx) {
1704: PetscBTDestroy(&btable);
1705: }
1707: /* View distributed dmnetwork */
1708: DMViewFromOptions(newDM,NULL,"-dmnetwork_view_distributed");
1710: *dm = newDM;
1711: return 0;
1712: }
1714: /*@C
1715: PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering
1717: Collective
1719: Input Parameters:
1720: + mainSF - the original SF structure
1721: - map - a ISLocalToGlobal mapping that contains the subset of points
1723: Output Parameter:
1724: . subSF - a subset of the mainSF for the desired subset.
1726: Level: intermediate
1727: @*/
1728: PetscErrorCode PetscSFGetSubSF(PetscSF mainsf,ISLocalToGlobalMapping map,PetscSF *subSF)
1729: {
1730: PetscInt nroots, nleaves, *ilocal_sub;
1731: PetscInt i, *ilocal_map, nroots_sub, nleaves_sub = 0;
1732: PetscInt *local_points, *remote_points;
1733: PetscSFNode *iremote_sub;
1734: const PetscInt *ilocal;
1735: const PetscSFNode *iremote;
1737: PetscSFGetGraph(mainsf,&nroots,&nleaves,&ilocal,&iremote);
1739: /* Look for leaves that pertain to the subset of points. Get the local ordering */
1740: PetscMalloc1(nleaves,&ilocal_map);
1741: ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);
1742: for (i = 0; i < nleaves; i++) {
1743: if (ilocal_map[i] != -1) nleaves_sub += 1;
1744: }
1745: /* Re-number ilocal with subset numbering. Need information from roots */
1746: PetscMalloc2(nroots,&local_points,nroots,&remote_points);
1747: for (i = 0; i < nroots; i++) local_points[i] = i;
1748: ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);
1749: PetscSFBcastBegin(mainsf, MPIU_INT, local_points, remote_points,MPI_REPLACE);
1750: PetscSFBcastEnd(mainsf, MPIU_INT, local_points, remote_points,MPI_REPLACE);
1751: /* Fill up graph using local (that is, local to the subset) numbering. */
1752: PetscMalloc1(nleaves_sub,&ilocal_sub);
1753: PetscMalloc1(nleaves_sub,&iremote_sub);
1754: nleaves_sub = 0;
1755: for (i = 0; i < nleaves; i++) {
1756: if (ilocal_map[i] != -1) {
1757: ilocal_sub[nleaves_sub] = ilocal_map[i];
1758: iremote_sub[nleaves_sub].rank = iremote[i].rank;
1759: iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
1760: nleaves_sub += 1;
1761: }
1762: }
1763: PetscFree2(local_points,remote_points);
1764: ISLocalToGlobalMappingGetSize(map,&nroots_sub);
1766: /* Create new subSF */
1767: PetscSFCreate(PETSC_COMM_WORLD,subSF);
1768: PetscSFSetFromOptions(*subSF);
1769: PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);
1770: PetscFree(ilocal_map);
1771: PetscFree(iremote_sub);
1772: return 0;
1773: }
1775: /*@C
1776: DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
1778: Not Collective
1780: Input Parameters:
1781: + dm - the DMNetwork object
1782: - p - the vertex point
1784: Output Parameters:
1785: + nedges - number of edges connected to this vertex point
1786: - edges - list of edge points
1788: Level: beginner
1790: Fortran Notes:
1791: Since it returns an array, this routine is only available in Fortran 90, and you must
1792: include petsc.h90 in your code.
1794: .seealso: DMNetworkCreate(), DMNetworkGetConnectedVertices()
1795: @*/
1796: PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
1797: {
1798: DM_Network *network = (DM_Network*)dm->data;
1800: DMPlexGetSupportSize(network->plex,vertex,nedges);
1801: if (edges) DMPlexGetSupport(network->plex,vertex,edges);
1802: return 0;
1803: }
1805: /*@C
1806: DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
1808: Not Collective
1810: Input Parameters:
1811: + dm - the DMNetwork object
1812: - p - the edge point
1814: Output Parameters:
1815: . vertices - vertices connected to this edge
1817: Level: beginner
1819: Fortran Notes:
1820: Since it returns an array, this routine is only available in Fortran 90, and you must
1821: include petsc.h90 in your code.
1823: .seealso: DMNetworkCreate(), DMNetworkGetSupportingEdges()
1824: @*/
1825: PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
1826: {
1827: DM_Network *network = (DM_Network*)dm->data;
1829: DMPlexGetCone(network->plex,edge,vertices);
1830: return 0;
1831: }
1833: /*@
1834: DMNetworkIsSharedVertex - Returns TRUE if the vertex is shared by subnetworks
1836: Not Collective
1838: Input Parameters:
1839: + dm - the DMNetwork object
1840: - p - the vertex point
1842: Output Parameter:
1843: . flag - TRUE if the vertex is shared by subnetworks
1845: Level: beginner
1847: .seealso: DMNetworkAddSharedVertices(), DMNetworkIsGhostVertex()
1848: @*/
1849: PetscErrorCode DMNetworkIsSharedVertex(DM dm,PetscInt p,PetscBool *flag)
1850: {
1851: PetscInt i;
1853: *flag = PETSC_FALSE;
1855: if (dm->setupcalled) { /* DMNetworkGetGlobalVertexIndex() requires DMSetUp() be called */
1856: DM_Network *network = (DM_Network*)dm->data;
1857: PetscInt gidx;
1858: DMNetworkGetGlobalVertexIndex(dm,p,&gidx);
1859: PetscTableFind(network->svtable,gidx+1,&i);
1860: if (i) *flag = PETSC_TRUE;
1861: } else { /* would be removed? */
1862: PetscInt nv;
1863: const PetscInt *vtx;
1864: DMNetworkGetSharedVertices(dm,&nv,&vtx);
1865: for (i=0; i<nv; i++) {
1866: if (p == vtx[i]) {
1867: *flag = PETSC_TRUE;
1868: break;
1869: }
1870: }
1871: }
1872: return 0;
1873: }
1875: /*@
1876: DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
1878: Not Collective
1880: Input Parameters:
1881: + dm - the DMNetwork object
1882: - p - the vertex point
1884: Output Parameter:
1885: . isghost - TRUE if the vertex is a ghost point
1887: Level: beginner
1889: .seealso: DMNetworkGetConnectedVertices(), DMNetworkGetVertexRange(), DMNetworkIsSharedVertex()
1890: @*/
1891: PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
1892: {
1893: DM_Network *network = (DM_Network*)dm->data;
1894: PetscInt offsetg;
1895: PetscSection sectiong;
1897: *isghost = PETSC_FALSE;
1898: DMGetGlobalSection(network->plex,§iong);
1899: PetscSectionGetOffset(sectiong,p,&offsetg);
1900: if (offsetg < 0) *isghost = PETSC_TRUE;
1901: return 0;
1902: }
1904: PetscErrorCode DMSetUp_Network(DM dm)
1905: {
1906: DM_Network *network=(DM_Network*)dm->data;
1908: DMNetworkComponentSetUp(dm);
1909: DMNetworkVariablesSetUp(dm);
1911: DMSetLocalSection(network->plex,network->DofSection);
1912: DMGetGlobalSection(network->plex,&network->GlobalDofSection);
1914: dm->setupcalled = PETSC_TRUE;
1916: /* View dmnetwork */
1917: DMViewFromOptions(dm,NULL,"-dmnetwork_view");
1918: return 0;
1919: }
1921: /*@
1922: DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
1923: -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
1925: Collective
1927: Input Parameters:
1928: + dm - the DMNetwork object
1929: . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
1930: - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
1932: Level: intermediate
1934: @*/
1935: PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
1936: {
1937: DM_Network *network=(DM_Network*)dm->data;
1938: PetscInt nVertices = network->nVertices;
1940: network->userEdgeJacobian = eflg;
1941: network->userVertexJacobian = vflg;
1943: if (eflg && !network->Je) {
1944: PetscCalloc1(3*network->nEdges,&network->Je);
1945: }
1947: if (vflg && !network->Jv && nVertices) {
1948: PetscInt i,*vptr,nedges,vStart=network->vStart;
1949: PetscInt nedges_total;
1950: const PetscInt *edges;
1952: /* count nvertex_total */
1953: nedges_total = 0;
1954: PetscMalloc1(nVertices+1,&vptr);
1956: vptr[0] = 0;
1957: for (i=0; i<nVertices; i++) {
1958: DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);
1959: nedges_total += nedges;
1960: vptr[i+1] = vptr[i] + 2*nedges + 1;
1961: }
1963: PetscCalloc1(2*nedges_total+nVertices,&network->Jv);
1964: network->Jvptr = vptr;
1965: }
1966: return 0;
1967: }
1969: /*@
1970: DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
1972: Not Collective
1974: Input Parameters:
1975: + dm - the DMNetwork object
1976: . p - the edge point
1977: - J - array (size = 3) of Jacobian submatrices for this edge point:
1978: J[0]: this edge
1979: J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()
1981: Level: advanced
1983: .seealso: DMNetworkVertexSetMatrix()
1984: @*/
1985: PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
1986: {
1987: DM_Network *network=(DM_Network*)dm->data;
1991: if (J) {
1992: network->Je[3*p] = J[0];
1993: network->Je[3*p+1] = J[1];
1994: network->Je[3*p+2] = J[2];
1995: }
1996: return 0;
1997: }
1999: /*@
2000: DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
2002: Not Collective
2004: Input Parameters:
2005: + dm - The DMNetwork object
2006: . p - the vertex point
2007: - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
2008: J[0]: this vertex
2009: J[1+2*i]: i-th supporting edge
2010: J[1+2*i+1]: i-th connected vertex
2012: Level: advanced
2014: .seealso: DMNetworkEdgeSetMatrix()
2015: @*/
2016: PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
2017: {
2018: DM_Network *network=(DM_Network*)dm->data;
2019: PetscInt i,*vptr,nedges,vStart=network->vStart;
2020: const PetscInt *edges;
2024: if (J) {
2025: vptr = network->Jvptr;
2026: network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
2028: /* Set Jacobian for each supporting edge and connected vertex */
2029: DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);
2030: for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
2031: }
2032: return 0;
2033: }
2035: static inline PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
2036: {
2037: PetscInt j;
2038: PetscScalar val=(PetscScalar)ncols;
2040: if (!ghost) {
2041: for (j=0; j<nrows; j++) {
2042: VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);
2043: }
2044: } else {
2045: for (j=0; j<nrows; j++) {
2046: VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);
2047: }
2048: }
2049: return 0;
2050: }
2052: static inline PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
2053: {
2054: PetscInt j,ncols_u;
2055: PetscScalar val;
2057: if (!ghost) {
2058: for (j=0; j<nrows; j++) {
2059: MatGetRow(Ju,j,&ncols_u,NULL,NULL);
2060: val = (PetscScalar)ncols_u;
2061: VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);
2062: MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);
2063: }
2064: } else {
2065: for (j=0; j<nrows; j++) {
2066: MatGetRow(Ju,j,&ncols_u,NULL,NULL);
2067: val = (PetscScalar)ncols_u;
2068: VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);
2069: MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);
2070: }
2071: }
2072: return 0;
2073: }
2075: static inline PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
2076: {
2077: if (Ju) {
2078: MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);
2079: } else {
2080: MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);
2081: }
2082: return 0;
2083: }
2085: static inline PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
2086: {
2087: PetscInt j,*cols;
2088: PetscScalar *zeros;
2090: PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);
2091: for (j=0; j<ncols; j++) cols[j] = j+ cstart;
2092: MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);
2093: PetscFree2(cols,zeros);
2094: return 0;
2095: }
2097: static inline PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
2098: {
2099: PetscInt j,M,N,row,col,ncols_u;
2100: const PetscInt *cols;
2101: PetscScalar zero=0.0;
2103: MatGetSize(Ju,&M,&N);
2106: for (row=0; row<nrows; row++) {
2107: MatGetRow(Ju,row,&ncols_u,&cols,NULL);
2108: for (j=0; j<ncols_u; j++) {
2109: col = cols[j] + cstart;
2110: MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);
2111: }
2112: MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);
2113: }
2114: return 0;
2115: }
2117: static inline PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
2118: {
2119: if (Ju) {
2120: MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);
2121: } else {
2122: MatSetDenseblock_private(nrows,rows,ncols,cstart,J);
2123: }
2124: return 0;
2125: }
2127: /* 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.
2128: */
2129: PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
2130: {
2131: PetscInt i,size,dof;
2132: PetscInt *glob2loc;
2134: PetscSectionGetStorageSize(localsec,&size);
2135: PetscMalloc1(size,&glob2loc);
2137: for (i = 0; i < size; i++) {
2138: PetscSectionGetOffset(globalsec,i,&dof);
2139: dof = (dof >= 0) ? dof : -(dof + 1);
2140: glob2loc[i] = dof;
2141: }
2143: ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);
2144: #if 0
2145: PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);
2146: #endif
2147: return 0;
2148: }
2150: #include <petsc/private/matimpl.h>
2152: PetscErrorCode DMCreateMatrix_Network_Nest(DM dm,Mat *J)
2153: {
2154: DM_Network *network = (DM_Network*)dm->data;
2155: PetscInt eDof,vDof;
2156: Mat j11,j12,j21,j22,bA[2][2];
2157: MPI_Comm comm;
2158: ISLocalToGlobalMapping eISMap,vISMap;
2160: PetscObjectGetComm((PetscObject)dm,&comm);
2162: PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);
2163: PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);
2165: MatCreate(comm, &j11);
2166: MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);
2167: MatSetType(j11, MATMPIAIJ);
2169: MatCreate(comm, &j12);
2170: MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);
2171: MatSetType(j12, MATMPIAIJ);
2173: MatCreate(comm, &j21);
2174: MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);
2175: MatSetType(j21, MATMPIAIJ);
2177: MatCreate(comm, &j22);
2178: MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);
2179: MatSetType(j22, MATMPIAIJ);
2181: bA[0][0] = j11;
2182: bA[0][1] = j12;
2183: bA[1][0] = j21;
2184: bA[1][1] = j22;
2186: CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);
2187: CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);
2189: MatSetLocalToGlobalMapping(j11,eISMap,eISMap);
2190: MatSetLocalToGlobalMapping(j12,eISMap,vISMap);
2191: MatSetLocalToGlobalMapping(j21,vISMap,eISMap);
2192: MatSetLocalToGlobalMapping(j22,vISMap,vISMap);
2194: MatSetUp(j11);
2195: MatSetUp(j12);
2196: MatSetUp(j21);
2197: MatSetUp(j22);
2199: MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);
2200: MatSetUp(*J);
2201: MatNestSetVecType(*J,VECNEST);
2202: MatDestroy(&j11);
2203: MatDestroy(&j12);
2204: MatDestroy(&j21);
2205: MatDestroy(&j22);
2207: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
2208: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
2209: MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
2211: /* Free structures */
2212: ISLocalToGlobalMappingDestroy(&eISMap);
2213: ISLocalToGlobalMappingDestroy(&vISMap);
2214: return 0;
2215: }
2217: PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
2218: {
2219: DM_Network *network = (DM_Network*)dm->data;
2220: PetscInt eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
2221: PetscInt cstart,ncols,j,e,v;
2222: PetscBool ghost,ghost_vc,ghost2,isNest;
2223: Mat Juser;
2224: PetscSection sectionGlobal;
2225: PetscInt nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
2226: const PetscInt *edges,*cone;
2227: MPI_Comm comm;
2228: MatType mtype;
2229: Vec vd_nz,vo_nz;
2230: PetscInt *dnnz,*onnz;
2231: PetscScalar *vdnz,*vonz;
2233: mtype = dm->mattype;
2234: PetscStrcmp(mtype,MATNEST,&isNest);
2235: if (isNest) {
2236: DMCreateMatrix_Network_Nest(dm,J);
2237: MatSetDM(*J,dm);
2238: return 0;
2239: }
2241: if (!network->userEdgeJacobian && !network->userVertexJacobian) {
2242: /* user does not provide Jacobian blocks */
2243: DMCreateMatrix_Plex(network->plex,J);
2244: MatSetDM(*J,dm);
2245: return 0;
2246: }
2248: MatCreate(PetscObjectComm((PetscObject)dm),J);
2249: DMGetGlobalSection(network->plex,§ionGlobal);
2250: PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);
2251: MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);
2253: MatSetType(*J,MATAIJ);
2254: MatSetFromOptions(*J);
2256: /* (1) Set matrix preallocation */
2257: /*------------------------------*/
2258: PetscObjectGetComm((PetscObject)dm,&comm);
2259: VecCreate(comm,&vd_nz);
2260: VecSetSizes(vd_nz,localSize,PETSC_DECIDE);
2261: VecSetFromOptions(vd_nz);
2262: VecSet(vd_nz,0.0);
2263: VecDuplicate(vd_nz,&vo_nz);
2265: /* Set preallocation for edges */
2266: /*-----------------------------*/
2267: DMNetworkGetEdgeRange(dm,&eStart,&eEnd);
2269: PetscMalloc1(localSize,&rows);
2270: for (e=eStart; e<eEnd; e++) {
2271: /* Get row indices */
2272: DMNetworkGetGlobalVecOffset(dm,e,ALL_COMPONENTS,&rstart);
2273: PetscSectionGetDof(network->DofSection,e,&nrows);
2274: if (nrows) {
2275: for (j=0; j<nrows; j++) rows[j] = j + rstart;
2277: /* Set preallocation for connected vertices */
2278: DMNetworkGetConnectedVertices(dm,e,&cone);
2279: for (v=0; v<2; v++) {
2280: PetscSectionGetDof(network->DofSection,cone[v],&ncols);
2282: if (network->Je) {
2283: Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
2284: } else Juser = NULL;
2285: DMNetworkIsGhostVertex(dm,cone[v],&ghost);
2286: MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);
2287: }
2289: /* Set preallocation for edge self */
2290: cstart = rstart;
2291: if (network->Je) {
2292: Juser = network->Je[3*e]; /* Jacobian(e,e) */
2293: } else Juser = NULL;
2294: MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);
2295: }
2296: }
2298: /* Set preallocation for vertices */
2299: /*--------------------------------*/
2300: DMNetworkGetVertexRange(dm,&vStart,&vEnd);
2301: if (vEnd - vStart) vptr = network->Jvptr;
2303: for (v=vStart; v<vEnd; v++) {
2304: /* Get row indices */
2305: DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&rstart);
2306: PetscSectionGetDof(network->DofSection,v,&nrows);
2307: if (!nrows) continue;
2309: DMNetworkIsGhostVertex(dm,v,&ghost);
2310: if (ghost) {
2311: PetscMalloc1(nrows,&rows_v);
2312: } else {
2313: rows_v = rows;
2314: }
2316: for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
2318: /* Get supporting edges and connected vertices */
2319: DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);
2321: for (e=0; e<nedges; e++) {
2322: /* Supporting edges */
2323: DMNetworkGetGlobalVecOffset(dm,edges[e],ALL_COMPONENTS,&cstart);
2324: PetscSectionGetDof(network->DofSection,edges[e],&ncols);
2326: if (network->Jv) {
2327: Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
2328: } else Juser = NULL;
2329: MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);
2331: /* Connected vertices */
2332: DMNetworkGetConnectedVertices(dm,edges[e],&cone);
2333: vc = (v == cone[0]) ? cone[1]:cone[0];
2334: DMNetworkIsGhostVertex(dm,vc,&ghost_vc);
2336: PetscSectionGetDof(network->DofSection,vc,&ncols);
2338: if (network->Jv) {
2339: Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
2340: } else Juser = NULL;
2341: if (ghost_vc||ghost) {
2342: ghost2 = PETSC_TRUE;
2343: } else {
2344: ghost2 = PETSC_FALSE;
2345: }
2346: MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);
2347: }
2349: /* Set preallocation for vertex self */
2350: DMNetworkIsGhostVertex(dm,v,&ghost);
2351: if (!ghost) {
2352: DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&cstart);
2353: if (network->Jv) {
2354: Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
2355: } else Juser = NULL;
2356: MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);
2357: }
2358: if (ghost) {
2359: PetscFree(rows_v);
2360: }
2361: }
2363: VecAssemblyBegin(vd_nz);
2364: VecAssemblyBegin(vo_nz);
2366: PetscMalloc2(localSize,&dnnz,localSize,&onnz);
2368: VecAssemblyEnd(vd_nz);
2369: VecAssemblyEnd(vo_nz);
2371: VecGetArray(vd_nz,&vdnz);
2372: VecGetArray(vo_nz,&vonz);
2373: for (j=0; j<localSize; j++) {
2374: dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
2375: onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
2376: }
2377: VecRestoreArray(vd_nz,&vdnz);
2378: VecRestoreArray(vo_nz,&vonz);
2379: VecDestroy(&vd_nz);
2380: VecDestroy(&vo_nz);
2382: MatSeqAIJSetPreallocation(*J,0,dnnz);
2383: MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);
2384: MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
2386: PetscFree2(dnnz,onnz);
2388: /* (2) Set matrix entries for edges */
2389: /*----------------------------------*/
2390: for (e=eStart; e<eEnd; e++) {
2391: /* Get row indices */
2392: DMNetworkGetGlobalVecOffset(dm,e,ALL_COMPONENTS,&rstart);
2393: PetscSectionGetDof(network->DofSection,e,&nrows);
2394: if (nrows) {
2395: for (j=0; j<nrows; j++) rows[j] = j + rstart;
2397: /* Set matrix entries for connected vertices */
2398: DMNetworkGetConnectedVertices(dm,e,&cone);
2399: for (v=0; v<2; v++) {
2400: DMNetworkGetGlobalVecOffset(dm,cone[v],ALL_COMPONENTS,&cstart);
2401: PetscSectionGetDof(network->DofSection,cone[v],&ncols);
2403: if (network->Je) {
2404: Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
2405: } else Juser = NULL;
2406: MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);
2407: }
2409: /* Set matrix entries for edge self */
2410: cstart = rstart;
2411: if (network->Je) {
2412: Juser = network->Je[3*e]; /* Jacobian(e,e) */
2413: } else Juser = NULL;
2414: MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);
2415: }
2416: }
2418: /* Set matrix entries for vertices */
2419: /*---------------------------------*/
2420: for (v=vStart; v<vEnd; v++) {
2421: /* Get row indices */
2422: DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&rstart);
2423: PetscSectionGetDof(network->DofSection,v,&nrows);
2424: if (!nrows) continue;
2426: DMNetworkIsGhostVertex(dm,v,&ghost);
2427: if (ghost) {
2428: PetscMalloc1(nrows,&rows_v);
2429: } else {
2430: rows_v = rows;
2431: }
2432: for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
2434: /* Get supporting edges and connected vertices */
2435: DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);
2437: for (e=0; e<nedges; e++) {
2438: /* Supporting edges */
2439: DMNetworkGetGlobalVecOffset(dm,edges[e],ALL_COMPONENTS,&cstart);
2440: PetscSectionGetDof(network->DofSection,edges[e],&ncols);
2442: if (network->Jv) {
2443: Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
2444: } else Juser = NULL;
2445: MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);
2447: /* Connected vertices */
2448: DMNetworkGetConnectedVertices(dm,edges[e],&cone);
2449: vc = (v == cone[0]) ? cone[1]:cone[0];
2451: DMNetworkGetGlobalVecOffset(dm,vc,ALL_COMPONENTS,&cstart);
2452: PetscSectionGetDof(network->DofSection,vc,&ncols);
2454: if (network->Jv) {
2455: Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
2456: } else Juser = NULL;
2457: MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);
2458: }
2460: /* Set matrix entries for vertex self */
2461: if (!ghost) {
2462: DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&cstart);
2463: if (network->Jv) {
2464: Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
2465: } else Juser = NULL;
2466: MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);
2467: }
2468: if (ghost) {
2469: PetscFree(rows_v);
2470: }
2471: }
2472: PetscFree(rows);
2474: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
2475: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
2477: MatSetDM(*J,dm);
2478: return 0;
2479: }
2481: PetscErrorCode DMDestroy_Network(DM dm)
2482: {
2483: DM_Network *network = (DM_Network*)dm->data;
2484: PetscInt j,np;
2486: if (--network->refct > 0) return 0;
2487: PetscFree(network->Je);
2488: if (network->Jv) {
2489: PetscFree(network->Jvptr);
2490: PetscFree(network->Jv);
2491: }
2493: ISLocalToGlobalMappingDestroy(&network->vertex.mapping);
2494: PetscSectionDestroy(&network->vertex.DofSection);
2495: PetscSectionDestroy(&network->vertex.GlobalDofSection);
2496: PetscFree(network->vltog);
2497: PetscSFDestroy(&network->vertex.sf);
2498: /* edge */
2499: ISLocalToGlobalMappingDestroy(&network->edge.mapping);
2500: PetscSectionDestroy(&network->edge.DofSection);
2501: PetscSectionDestroy(&network->edge.GlobalDofSection);
2502: PetscSFDestroy(&network->edge.sf);
2503: DMDestroy(&network->plex);
2504: PetscSectionDestroy(&network->DataSection);
2505: PetscSectionDestroy(&network->DofSection);
2507: for (j=0; j<network->Nsvtx; j++) PetscFree(network->svtx[j].sv);
2508: PetscFree(network->svtx);
2509: PetscFree2(network->subnetedge,network->subnetvtx);
2511: PetscTableDestroy(&network->svtable);
2512: PetscFree(network->subnet);
2513: PetscFree(network->component);
2514: PetscFree(network->componentdataarray);
2516: if (network->header) {
2517: np = network->pEnd - network->pStart;
2518: for (j=0; j < np; j++) {
2519: PetscFree5(network->header[j].size,network->header[j].key,network->header[j].offset,network->header[j].nvar,network->header[j].offsetvarrel);
2520: PetscFree(network->cvalue[j].data);
2521: }
2522: PetscFree2(network->header,network->cvalue);
2523: }
2524: PetscFree(network);
2525: return 0;
2526: }
2528: PetscErrorCode DMView_Network(DM dm,PetscViewer viewer)
2529: {
2530: PetscBool iascii;
2531: PetscMPIInt rank;
2535: MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
2536: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
2537: if (iascii) {
2538: const PetscInt *cone,*vtx,*edges;
2539: PetscInt vfrom,vto,i,j,nv,ne,nsv,p,nsubnet;
2540: DM_Network *network = (DM_Network*)dm->data;
2542: nsubnet = network->Nsubnet; /* num of subnetworks */
2543: if (rank == 0) {
2544: PetscPrintf(PETSC_COMM_SELF," NSubnets: %" PetscInt_FMT "; NEdges: %" PetscInt_FMT "; NVertices: %" PetscInt_FMT "; NSharedVertices: %" PetscInt_FMT ".\n",nsubnet,network->NEdges,network->NVertices,network->Nsvtx);
2545: }
2547: DMNetworkGetSharedVertices(dm,&nsv,NULL);
2548: PetscViewerASCIIPushSynchronized(viewer);
2549: PetscViewerASCIISynchronizedPrintf(viewer, " [%d] nEdges: %" PetscInt_FMT "; nVertices: %" PetscInt_FMT "; nSharedVertices: %" PetscInt_FMT "\n",rank,network->nEdges,network->nVertices,nsv);
2551: for (i=0; i<nsubnet; i++) {
2552: DMNetworkGetSubnetwork(dm,i,&nv,&ne,&vtx,&edges);
2553: if (ne) {
2554: PetscViewerASCIISynchronizedPrintf(viewer, " Subnet %" PetscInt_FMT ": nEdges %" PetscInt_FMT ", nVertices(include shared vertices) %" PetscInt_FMT "\n",i,ne,nv);
2555: for (j=0; j<ne; j++) {
2556: p = edges[j];
2557: DMNetworkGetConnectedVertices(dm,p,&cone);
2558: DMNetworkGetGlobalVertexIndex(dm,cone[0],&vfrom);
2559: DMNetworkGetGlobalVertexIndex(dm,cone[1],&vto);
2560: DMNetworkGetGlobalEdgeIndex(dm,edges[j],&p);
2561: PetscViewerASCIISynchronizedPrintf(viewer, " edge %" PetscInt_FMT ": %" PetscInt_FMT " ----> %" PetscInt_FMT "\n",p,vfrom,vto);
2562: }
2563: }
2564: }
2566: /* Shared vertices */
2567: DMNetworkGetSharedVertices(dm,NULL,&vtx);
2568: if (nsv) {
2569: PetscInt gidx;
2570: PetscBool ghost;
2571: const PetscInt *sv=NULL;
2573: PetscViewerASCIISynchronizedPrintf(viewer, " SharedVertices:\n");
2574: for (i=0; i<nsv; i++) {
2575: DMNetworkIsGhostVertex(dm,vtx[i],&ghost);
2576: if (ghost) continue;
2578: DMNetworkSharedVertexGetInfo(dm,vtx[i],&gidx,&nv,&sv);
2579: PetscViewerASCIISynchronizedPrintf(viewer, " svtx %" PetscInt_FMT ": global index %" PetscInt_FMT ", subnet[%" PetscInt_FMT "].%" PetscInt_FMT " ---->\n",i,gidx,sv[0],sv[1]);
2580: for (j=1; j<nv; j++) {
2581: PetscViewerASCIISynchronizedPrintf(viewer, " ----> subnet[%" PetscInt_FMT "].%" PetscInt_FMT "\n",sv[2*j],sv[2*j+1]);
2582: }
2583: }
2584: }
2585: PetscViewerFlush(viewer);
2586: PetscViewerASCIIPopSynchronized(viewer);
2588: return 0;
2589: }
2591: PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
2592: {
2593: DM_Network *network = (DM_Network*)dm->data;
2595: DMGlobalToLocalBegin(network->plex,g,mode,l);
2596: return 0;
2597: }
2599: PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
2600: {
2601: DM_Network *network = (DM_Network*)dm->data;
2603: DMGlobalToLocalEnd(network->plex,g,mode,l);
2604: return 0;
2605: }
2607: PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
2608: {
2609: DM_Network *network = (DM_Network*)dm->data;
2611: DMLocalToGlobalBegin(network->plex,l,mode,g);
2612: return 0;
2613: }
2615: PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
2616: {
2617: DM_Network *network = (DM_Network*)dm->data;
2619: DMLocalToGlobalEnd(network->plex,l,mode,g);
2620: return 0;
2621: }
2623: /*@
2624: DMNetworkGetVertexLocalToGlobalOrdering - Get vertex global index
2626: Not collective
2628: Input Parameters:
2629: + dm - the dm object
2630: - vloc - local vertex ordering, start from 0
2632: Output Parameters:
2633: . vg - global vertex ordering, start from 0
2635: Level: advanced
2637: .seealso: DMNetworkSetVertexLocalToGlobalOrdering()
2638: @*/
2639: PetscErrorCode DMNetworkGetVertexLocalToGlobalOrdering(DM dm,PetscInt vloc,PetscInt *vg)
2640: {
2641: DM_Network *network = (DM_Network*)dm->data;
2642: PetscInt *vltog = network->vltog;
2645: *vg = vltog[vloc];
2646: return 0;
2647: }
2649: /*@
2650: DMNetworkSetVertexLocalToGlobalOrdering - Create and setup vertex local to global map
2652: Collective
2654: Input Parameters:
2655: . dm - the dm object
2657: Level: advanced
2659: .seealso: DMNetworkGetGlobalVertexIndex()
2660: @*/
2661: PetscErrorCode DMNetworkSetVertexLocalToGlobalOrdering(DM dm)
2662: {
2663: DM_Network *network=(DM_Network*)dm->data;
2664: MPI_Comm comm;
2665: PetscMPIInt rank,size,*displs=NULL,*recvcounts=NULL,remoterank;
2666: PetscBool ghost;
2667: PetscInt *vltog,nroots,nleaves,i,*vrange,k,N,lidx;
2668: const PetscSFNode *iremote;
2669: PetscSF vsf;
2670: Vec Vleaves,Vleaves_seq;
2671: VecScatter ctx;
2672: PetscScalar *varr,val;
2673: const PetscScalar *varr_read;
2675: PetscObjectGetComm((PetscObject)dm,&comm);
2676: MPI_Comm_size(comm,&size);
2677: MPI_Comm_rank(comm,&rank);
2679: if (size == 1) {
2680: nroots = network->vEnd - network->vStart;
2681: PetscMalloc1(nroots, &vltog);
2682: for (i=0; i<nroots; i++) vltog[i] = i;
2683: network->vltog = vltog;
2684: return 0;
2685: }
2688: if (network->vltog) {
2689: PetscFree(network->vltog);
2690: }
2692: DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);
2693: PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);
2694: vsf = network->vertex.sf;
2696: PetscMalloc3(size+1,&vrange,size,&displs,size,&recvcounts);
2697: PetscSFGetGraph(vsf,&nroots,&nleaves,NULL,&iremote);
2699: for (i=0; i<size; i++) { displs[i] = i; recvcounts[i] = 1;}
2701: i = nroots - nleaves; /* local number of vertices, excluding ghosts */
2702: vrange[0] = 0;
2703: MPI_Allgatherv(&i,1,MPIU_INT,vrange+1,recvcounts,displs,MPIU_INT,comm);
2704: for (i=2; i<size+1; i++) {vrange[i] += vrange[i-1];}
2706: PetscMalloc1(nroots, &vltog);
2707: network->vltog = vltog;
2709: /* Set vltog for non-ghost vertices */
2710: k = 0;
2711: for (i=0; i<nroots; i++) {
2712: DMNetworkIsGhostVertex(dm,i+network->vStart,&ghost);
2713: if (ghost) continue;
2714: vltog[i] = vrange[rank] + k++;
2715: }
2716: PetscFree3(vrange,displs,recvcounts);
2718: /* Set vltog for ghost vertices */
2719: /* (a) create parallel Vleaves and sequential Vleaves_seq to convert local iremote[*].index to global index */
2720: VecCreate(comm,&Vleaves);
2721: VecSetSizes(Vleaves,2*nleaves,PETSC_DETERMINE);
2722: VecSetFromOptions(Vleaves);
2723: VecGetArray(Vleaves,&varr);
2724: for (i=0; i<nleaves; i++) {
2725: varr[2*i] = (PetscScalar)(iremote[i].rank); /* rank of remote process */
2726: varr[2*i+1] = (PetscScalar)(iremote[i].index); /* local index in remote process */
2727: }
2728: VecRestoreArray(Vleaves,&varr);
2730: /* (b) scatter local info to remote processes via VecScatter() */
2731: VecScatterCreateToAll(Vleaves,&ctx,&Vleaves_seq);
2732: VecScatterBegin(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD);
2733: VecScatterEnd(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD);
2735: /* (c) convert local indices to global indices in parallel vector Vleaves */
2736: VecGetSize(Vleaves_seq,&N);
2737: VecGetArrayRead(Vleaves_seq,&varr_read);
2738: for (i=0; i<N; i+=2) {
2739: remoterank = (PetscMPIInt)PetscRealPart(varr_read[i]);
2740: if (remoterank == rank) {
2741: k = i+1; /* row number */
2742: lidx = (PetscInt)PetscRealPart(varr_read[i+1]);
2743: val = (PetscScalar)vltog[lidx]; /* global index for non-ghost vertex computed above */
2744: VecSetValues(Vleaves,1,&k,&val,INSERT_VALUES);
2745: }
2746: }
2747: VecRestoreArrayRead(Vleaves_seq,&varr_read);
2748: VecAssemblyBegin(Vleaves);
2749: VecAssemblyEnd(Vleaves);
2751: /* (d) Set vltog for ghost vertices by copying local values of Vleaves */
2752: VecGetArrayRead(Vleaves,&varr_read);
2753: k = 0;
2754: for (i=0; i<nroots; i++) {
2755: DMNetworkIsGhostVertex(dm,i+network->vStart,&ghost);
2756: if (!ghost) continue;
2757: vltog[i] = (PetscInt)PetscRealPart(varr_read[2*k+1]); k++;
2758: }
2759: VecRestoreArrayRead(Vleaves,&varr_read);
2761: VecDestroy(&Vleaves);
2762: VecDestroy(&Vleaves_seq);
2763: VecScatterDestroy(&ctx);
2764: return 0;
2765: }
2767: static inline PetscErrorCode DMISAddSize_private(DM_Network *network,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *nidx)
2768: {
2769: PetscInt i,j,ncomps,nvar,key,offset=0;
2770: DMNetworkComponentHeader header;
2772: PetscSectionGetOffset(network->DataSection,p,&offset);
2773: ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
2774: header = (DMNetworkComponentHeader)(network->componentdataarray+offset);
2776: for (i=0; i<ncomps; i++) {
2777: key = header->key[i];
2778: nvar = header->nvar[i];
2779: for (j=0; j<numkeys; j++) {
2780: if (key == keys[j]) {
2781: if (!blocksize || blocksize[j] == -1) {
2782: *nidx += nvar;
2783: } else {
2784: *nidx += nselectedvar[j]*nvar/blocksize[j];
2785: }
2786: }
2787: }
2788: }
2789: return 0;
2790: }
2792: static inline PetscErrorCode DMISComputeIdx_private(DM dm,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],PetscInt *ii,PetscInt *idx)
2793: {
2794: PetscInt i,j,ncomps,nvar,key,offsetg,k,k1,offset=0;
2795: DM_Network *network = (DM_Network*)dm->data;
2796: DMNetworkComponentHeader header;
2798: PetscSectionGetOffset(network->DataSection,p,&offset);
2799: ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
2800: header = (DMNetworkComponentHeader)(network->componentdataarray+offset);
2802: for (i=0; i<ncomps; i++) {
2803: key = header->key[i];
2804: nvar = header->nvar[i];
2805: for (j=0; j<numkeys; j++) {
2806: if (key != keys[j]) continue;
2808: DMNetworkGetGlobalVecOffset(dm,p,i,&offsetg);
2809: if (!blocksize || blocksize[j] == -1) {
2810: for (k=0; k<nvar; k++) idx[(*ii)++] = offsetg + k;
2811: } else {
2812: for (k=0; k<nvar; k+=blocksize[j]) {
2813: for (k1=0; k1<nselectedvar[j]; k1++) idx[(*ii)++] = offsetg + k + selectedvar[j][k1];
2814: }
2815: }
2816: }
2817: }
2818: return 0;
2819: }
2821: /*@
2822: DMNetworkCreateIS - Create an index set object from the global vector of the network
2824: Collective
2826: Input Parameters:
2827: + dm - DMNetwork object
2828: . numkeys - number of keys
2829: . keys - array of keys that define the components of the variables you wish to extract
2830: . blocksize - block size of the variables associated to the component
2831: . nselectedvar - number of variables in each block to select
2832: - selectedvar - the offset into the block of each variable in each block to select
2834: Output Parameters:
2835: . is - the index set
2837: Level: Advanced
2839: Notes:
2840: 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.
2842: .seealso: DMNetworkCreate(), ISCreateGeneral(), DMNetworkCreateLocalIS()
2843: @*/
2844: PetscErrorCode DMNetworkCreateIS(DM dm,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],IS *is)
2845: {
2846: MPI_Comm comm;
2847: DM_Network *network = (DM_Network*)dm->data;
2848: PetscInt i,p,estart,eend,vstart,vend,nidx,*idx;
2849: PetscBool ghost;
2851: PetscObjectGetComm((PetscObject)dm,&comm);
2853: /* Check input parameters */
2854: for (i=0; i<numkeys; i++) {
2855: if (!blocksize || blocksize[i] == -1) continue;
2857: }
2859: DMNetworkGetEdgeRange(dm,&estart,&eend);
2860: DMNetworkGetVertexRange(dm,&vstart,&vend);
2862: /* Get local number of idx */
2863: nidx = 0;
2864: for (p=estart; p<eend; p++) {
2865: DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx);
2866: }
2867: for (p=vstart; p<vend; p++) {
2868: DMNetworkIsGhostVertex(dm,p,&ghost);
2869: if (ghost) continue;
2870: DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx);
2871: }
2873: /* Compute idx */
2874: PetscMalloc1(nidx,&idx);
2875: i = 0;
2876: for (p=estart; p<eend; p++) {
2877: DMISComputeIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx);
2878: }
2879: for (p=vstart; p<vend; p++) {
2880: DMNetworkIsGhostVertex(dm,p,&ghost);
2881: if (ghost) continue;
2882: DMISComputeIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx);
2883: }
2885: /* Create is */
2886: ISCreateGeneral(comm,nidx,idx,PETSC_COPY_VALUES,is);
2887: PetscFree(idx);
2888: return 0;
2889: }
2891: static inline PetscErrorCode DMISComputeLocalIdx_private(DM dm,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],PetscInt *ii,PetscInt *idx)
2892: {
2893: PetscInt i,j,ncomps,nvar,key,offsetl,k,k1,offset=0;
2894: DM_Network *network = (DM_Network*)dm->data;
2895: DMNetworkComponentHeader header;
2897: PetscSectionGetOffset(network->DataSection,p,&offset);
2898: ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
2899: header = (DMNetworkComponentHeader)(network->componentdataarray+offset);
2901: for (i=0; i<ncomps; i++) {
2902: key = header->key[i];
2903: nvar = header->nvar[i];
2904: for (j=0; j<numkeys; j++) {
2905: if (key != keys[j]) continue;
2907: DMNetworkGetLocalVecOffset(dm,p,i,&offsetl);
2908: if (!blocksize || blocksize[j] == -1) {
2909: for (k=0; k<nvar; k++) idx[(*ii)++] = offsetl + k;
2910: } else {
2911: for (k=0; k<nvar; k+=blocksize[j]) {
2912: for (k1=0; k1<nselectedvar[j]; k1++) idx[(*ii)++] = offsetl + k + selectedvar[j][k1];
2913: }
2914: }
2915: }
2916: }
2917: return 0;
2918: }
2920: /*@
2921: DMNetworkCreateLocalIS - Create an index set object from the local vector of the network
2923: Not collective
2925: Input Parameters:
2926: + dm - DMNetwork object
2927: . numkeys - number of keys
2928: . keys - array of keys that define the components of the variables you wish to extract
2929: . blocksize - block size of the variables associated to the component
2930: . nselectedvar - number of variables in each block to select
2931: - selectedvar - the offset into the block of each variable in each block to select
2933: Output Parameters:
2934: . is - the index set
2936: Level: Advanced
2938: Notes:
2939: 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.
2941: .seealso: DMNetworkCreate(), DMNetworkCreateIS, ISCreateGeneral()
2942: @*/
2943: PetscErrorCode DMNetworkCreateLocalIS(DM dm,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],IS *is)
2944: {
2945: DM_Network *network = (DM_Network*)dm->data;
2946: PetscInt i,p,pstart,pend,nidx,*idx;
2948: /* Check input parameters */
2949: for (i=0; i<numkeys; i++) {
2950: if (!blocksize || blocksize[i] == -1) continue;
2952: }
2954: pstart = network->pStart;
2955: pend = network->pEnd;
2957: /* Get local number of idx */
2958: nidx = 0;
2959: for (p=pstart; p<pend; p++) {
2960: DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx);
2961: }
2963: /* Compute local idx */
2964: PetscMalloc1(nidx,&idx);
2965: i = 0;
2966: for (p=pstart; p<pend; p++) {
2967: DMISComputeLocalIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx);
2968: }
2970: /* Create is */
2971: ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,is);
2972: PetscFree(idx);
2973: return 0;
2974: }