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