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