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