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