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