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