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