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 headerarr += header->maxcomps; 1382 PetscCall(PetscMemcpy(headerarr,header->key,header->maxcomps*sizeof(PetscInt))); 1383 headerarr += header->maxcomps; 1384 PetscCall(PetscMemcpy(headerarr,header->offset,header->maxcomps*sizeof(PetscInt))); 1385 headerarr += header->maxcomps; 1386 PetscCall(PetscMemcpy(headerarr,header->nvar,header->maxcomps*sizeof(PetscInt))); 1387 headerarr += header->maxcomps; 1388 PetscCall(PetscMemcpy(headerarr,header->offsetvarrel,header->maxcomps*sizeof(PetscInt))); 1389 1390 /* Copy data */ 1391 cvalue = &network->cvalue[p]; 1392 ncomp = header->ndata; 1393 1394 for (i = 0; i < ncomp; i++) { 1395 offset = offsetp + header->hsize + header->offset[i]; 1396 PetscCall(PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType))); 1397 } 1398 } 1399 PetscFunctionReturn(0); 1400 } 1401 1402 /* Sets up the section for dofs. This routine is called during DMSetUp() */ 1403 static PetscErrorCode DMNetworkVariablesSetUp(DM dm) 1404 { 1405 DM_Network *network = (DM_Network*)dm->data; 1406 1407 PetscFunctionBegin; 1408 PetscCall(PetscSectionSetUp(network->DofSection)); 1409 PetscFunctionReturn(0); 1410 } 1411 1412 /* Get a subsection from a range of points */ 1413 static PetscErrorCode DMNetworkGetSubSection_private(PetscSection main,PetscInt pstart,PetscInt pend,PetscSection *subsection) 1414 { 1415 PetscInt i, nvar; 1416 1417 PetscFunctionBegin; 1418 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)main), subsection)); 1419 PetscCall(PetscSectionSetChart(*subsection, 0, pend - pstart)); 1420 for (i = pstart; i < pend; i++) { 1421 PetscCall(PetscSectionGetDof(main,i,&nvar)); 1422 PetscCall(PetscSectionSetDof(*subsection, i - pstart, nvar)); 1423 } 1424 1425 PetscCall(PetscSectionSetUp(*subsection)); 1426 PetscFunctionReturn(0); 1427 } 1428 1429 /* Create a submap of points with a GlobalToLocal structure */ 1430 static PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map) 1431 { 1432 PetscInt i, *subpoints; 1433 1434 PetscFunctionBegin; 1435 /* Create index sets to map from "points" to "subpoints" */ 1436 PetscCall(PetscMalloc1(pend - pstart, &subpoints)); 1437 for (i = pstart; i < pend; i++) { 1438 subpoints[i - pstart] = i; 1439 } 1440 PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map)); 1441 PetscCall(PetscFree(subpoints)); 1442 PetscFunctionReturn(0); 1443 } 1444 1445 /*@ 1446 DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute 1447 1448 Collective on dm 1449 1450 Input Parameters: 1451 . dm - the DMNetworkObject 1452 1453 Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are: 1454 1455 points = [0 1 2 3 4 5 6] 1456 1457 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]). 1458 1459 With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset. 1460 1461 Level: intermediate 1462 1463 @*/ 1464 PetscErrorCode DMNetworkAssembleGraphStructures(DM dm) 1465 { 1466 MPI_Comm comm; 1467 PetscMPIInt size; 1468 DM_Network *network = (DM_Network*)dm->data; 1469 1470 PetscFunctionBegin; 1471 PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 1472 PetscCallMPI(MPI_Comm_size(comm, &size)); 1473 1474 /* Create maps for vertices and edges */ 1475 PetscCall(DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping)); 1476 PetscCall(DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping)); 1477 1478 /* Create local sub-sections */ 1479 PetscCall(DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection)); 1480 PetscCall(DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection)); 1481 1482 if (size > 1) { 1483 PetscCall(PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf)); 1484 1485 PetscCall(PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection)); 1486 PetscCall(PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf)); 1487 PetscCall(PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection)); 1488 } else { 1489 /* create structures for vertex */ 1490 PetscCall(PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection)); 1491 /* create structures for edge */ 1492 PetscCall(PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection)); 1493 } 1494 1495 /* Add viewers */ 1496 PetscCall(PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section")); 1497 PetscCall(PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section")); 1498 PetscCall(PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view")); 1499 PetscCall(PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view")); 1500 PetscFunctionReturn(0); 1501 } 1502 1503 /* 1504 Setup a lookup btable for the input v's owning subnetworks 1505 - add all owing subnetworks that connect to this v to the btable 1506 vertex_subnetid = supportingedge_subnetid 1507 */ 1508 static inline PetscErrorCode SetSubnetIdLookupBT(DM dm,PetscInt v,PetscInt Nsubnet,PetscBT btable) 1509 { 1510 PetscInt e,nedges,offset; 1511 const PetscInt *edges; 1512 DM_Network *newDMnetwork = (DM_Network*)dm->data; 1513 DMNetworkComponentHeader header; 1514 1515 PetscFunctionBegin; 1516 PetscCall(PetscBTMemzero(Nsubnet,btable)); 1517 PetscCall(DMNetworkGetSupportingEdges(dm,v,&nedges,&edges)); 1518 for (e=0; e<nedges; e++) { 1519 PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,edges[e],&offset)); 1520 header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset); 1521 PetscCall(PetscBTSet(btable,header->subnetid)); 1522 } 1523 PetscFunctionReturn(0); 1524 } 1525 1526 /*@ 1527 DMNetworkDistribute - Distributes the network and moves associated component data 1528 1529 Collective 1530 1531 Input Parameters: 1532 + DM - the DMNetwork object 1533 - overlap - the overlap of partitions, 0 is the default 1534 1535 Options Database Key: 1536 + -dmnetwork_view - Calls DMView() at the conclusion of DMSetUp() 1537 - -dmnetwork_view_distributed - Calls DMView() at the conclusion of DMNetworkDistribute() 1538 1539 Notes: 1540 Distributes the network with <overlap>-overlapping partitioning of the edges. 1541 1542 Level: intermediate 1543 1544 .seealso: `DMNetworkCreate()` 1545 @*/ 1546 PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap) 1547 { 1548 MPI_Comm comm; 1549 PetscMPIInt size; 1550 DM_Network *oldDMnetwork = (DM_Network*)((*dm)->data); 1551 DM_Network *newDMnetwork; 1552 PetscSF pointsf=NULL; 1553 DM newDM; 1554 PetscInt j,e,v,offset,*subnetvtx,*subnetedge,Nsubnet,gidx,svtx_idx,nv; 1555 PetscInt net,*sv; 1556 PetscBT btable; 1557 PetscPartitioner part; 1558 DMNetworkComponentHeader header; 1559 1560 PetscFunctionBegin; 1561 PetscValidPointer(dm,1); 1562 PetscValidHeaderSpecific(*dm,DM_CLASSID,1); 1563 PetscCall(PetscObjectGetComm((PetscObject)*dm,&comm)); 1564 PetscCallMPI(MPI_Comm_size(comm, &size)); 1565 if (size == 1) PetscFunctionReturn(0); 1566 1567 PetscCheck(!overlap,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"overlap %" PetscInt_FMT " != 0 is not supported yet",overlap); 1568 1569 /* 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. */ 1570 PetscCall(DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM)); 1571 newDMnetwork = (DM_Network*)newDM->data; 1572 newDMnetwork->max_comps_registered = oldDMnetwork->max_comps_registered; 1573 PetscCall(PetscMalloc1(newDMnetwork->max_comps_registered,&newDMnetwork->component)); 1574 1575 /* Enable runtime options for petscpartitioner */ 1576 PetscCall(DMPlexGetPartitioner(oldDMnetwork->plex,&part)); 1577 PetscCall(PetscPartitionerSetFromOptions(part)); 1578 1579 /* Distribute plex dm */ 1580 PetscCall(DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex)); 1581 1582 /* Distribute dof section */ 1583 PetscCall(PetscSectionCreate(comm,&newDMnetwork->DofSection)); 1584 PetscCall(PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection)); 1585 1586 /* Distribute data and associated section */ 1587 PetscCall(PetscSectionCreate(comm,&newDMnetwork->DataSection)); 1588 PetscCall(DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray)); 1589 1590 PetscCall(PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd)); 1591 PetscCall(DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd)); 1592 PetscCall(DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd)); 1593 newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart; 1594 newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart; 1595 newDMnetwork->NVertices = oldDMnetwork->NVertices; 1596 newDMnetwork->NEdges = oldDMnetwork->NEdges; 1597 newDMnetwork->svtable = oldDMnetwork->svtable; /* global table! */ 1598 oldDMnetwork->svtable = NULL; 1599 1600 /* Set Dof section as the section for dm */ 1601 PetscCall(DMSetLocalSection(newDMnetwork->plex,newDMnetwork->DofSection)); 1602 PetscCall(DMGetGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection)); 1603 1604 /* Setup subnetwork info in the newDM */ 1605 newDMnetwork->Nsubnet = oldDMnetwork->Nsubnet; 1606 newDMnetwork->Nsvtx = oldDMnetwork->Nsvtx; 1607 oldDMnetwork->Nsvtx = 0; 1608 newDMnetwork->svtx = oldDMnetwork->svtx; /* global vertices! */ 1609 oldDMnetwork->svtx = NULL; 1610 PetscCall(PetscCalloc1(newDMnetwork->Nsubnet,&newDMnetwork->subnet)); 1611 1612 /* Copy over the global number of vertices and edges in each subnetwork. 1613 Note: these are calculated in DMNetworkLayoutSetUp() 1614 */ 1615 Nsubnet = newDMnetwork->Nsubnet; 1616 for (j = 0; j < Nsubnet; j++) { 1617 newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx; 1618 newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge; 1619 } 1620 1621 /* Count local nedges for subnetworks */ 1622 for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++) { 1623 PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset)); 1624 header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset); 1625 1626 /* Update pointers */ 1627 header->size = (PetscInt*)(header + 1); 1628 header->key = header->size + header->maxcomps; 1629 header->offset = header->key + header->maxcomps; 1630 header->nvar = header->offset + header->maxcomps; 1631 header->offsetvarrel = header->nvar + header->maxcomps; 1632 1633 newDMnetwork->subnet[header->subnetid].nedge++; 1634 } 1635 1636 /* Setup a btable to keep track subnetworks owned by this process at a shared vertex */ 1637 if (newDMnetwork->Nsvtx) { 1638 PetscCall(PetscBTCreate(Nsubnet,&btable)); 1639 } 1640 1641 /* Count local nvtx for subnetworks */ 1642 for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++) { 1643 PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset)); 1644 header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset); 1645 1646 /* Update pointers */ 1647 header->size = (PetscInt*)(header + 1); 1648 header->key = header->size + header->maxcomps; 1649 header->offset = header->key + header->maxcomps; 1650 header->nvar = header->offset + header->maxcomps; 1651 header->offsetvarrel = header->nvar + header->maxcomps; 1652 1653 /* shared vertices: use gidx=header->index to check if v is a shared vertex */ 1654 gidx = header->index; 1655 PetscCall(PetscTableFind(newDMnetwork->svtable,gidx+1,&svtx_idx)); 1656 svtx_idx--; 1657 1658 if (svtx_idx < 0) { /* not a shared vertex */ 1659 newDMnetwork->subnet[header->subnetid].nvtx++; 1660 } else { /* a shared vertex belongs to more than one subnetworks, it is being counted by multiple subnets */ 1661 /* Setup a lookup btable for this v's owning subnetworks */ 1662 PetscCall(SetSubnetIdLookupBT(newDM,v,Nsubnet,btable)); 1663 1664 for (j=0; j<newDMnetwork->svtx[svtx_idx].n; j++) { 1665 sv = newDMnetwork->svtx[svtx_idx].sv + 2*j; 1666 net = sv[0]; 1667 if (PetscBTLookup(btable,net)) newDMnetwork->subnet[net].nvtx++; /* sv is on net owned by this proces */ 1668 } 1669 } 1670 } 1671 1672 /* Get total local nvtx for subnetworks */ 1673 nv = 0; 1674 for (j=0; j<Nsubnet; j++) nv += newDMnetwork->subnet[j].nvtx; 1675 nv += newDMnetwork->Nsvtx; 1676 1677 /* Now create the vertices and edge arrays for the subnetworks */ 1678 PetscCall(PetscCalloc2(newDMnetwork->nEdges,&subnetedge,nv,&subnetvtx)); /* Maps local vertex to local subnetwork's vertex */ 1679 newDMnetwork->subnetedge = subnetedge; 1680 newDMnetwork->subnetvtx = subnetvtx; 1681 for (j=0; j < newDMnetwork->Nsubnet; j++) { 1682 newDMnetwork->subnet[j].edges = subnetedge; 1683 subnetedge += newDMnetwork->subnet[j].nedge; 1684 1685 newDMnetwork->subnet[j].vertices = subnetvtx; 1686 subnetvtx += newDMnetwork->subnet[j].nvtx; 1687 1688 /* 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. */ 1689 newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0; 1690 } 1691 newDMnetwork->svertices = subnetvtx; 1692 1693 /* Set the edges and vertices in each subnetwork */ 1694 for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++) { 1695 PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset)); 1696 header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset); 1697 newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++] = e; 1698 } 1699 1700 nv = 0; 1701 for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++) { 1702 PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset)); 1703 header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset); 1704 1705 /* coupling vertices: use gidx = header->index to check if v is a coupling vertex */ 1706 PetscCall(PetscTableFind(newDMnetwork->svtable,header->index+1,&svtx_idx)); 1707 svtx_idx--; 1708 if (svtx_idx < 0) { 1709 newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++] = v; 1710 } else { /* a shared vertex */ 1711 newDMnetwork->svertices[nv++] = v; 1712 1713 /* Setup a lookup btable for this v's owning subnetworks */ 1714 PetscCall(SetSubnetIdLookupBT(newDM,v,Nsubnet,btable)); 1715 1716 for (j=0; j<newDMnetwork->svtx[svtx_idx].n; j++) { 1717 sv = newDMnetwork->svtx[svtx_idx].sv + 2*j; 1718 net = sv[0]; 1719 if (PetscBTLookup(btable,net)) 1720 newDMnetwork->subnet[net].vertices[newDMnetwork->subnet[net].nvtx++] = v; 1721 } 1722 } 1723 } 1724 newDMnetwork->nsvtx = nv; /* num of local shared vertices */ 1725 1726 newDM->setupcalled = (*dm)->setupcalled; 1727 newDMnetwork->distributecalled = PETSC_TRUE; 1728 1729 /* Free spaces */ 1730 PetscCall(PetscSFDestroy(&pointsf)); 1731 PetscCall(DMDestroy(dm)); 1732 if (newDMnetwork->Nsvtx) { 1733 PetscCall(PetscBTDestroy(&btable)); 1734 } 1735 1736 /* View distributed dmnetwork */ 1737 PetscCall(DMViewFromOptions(newDM,NULL,"-dmnetwork_view_distributed")); 1738 1739 *dm = newDM; 1740 PetscFunctionReturn(0); 1741 } 1742 1743 /*@C 1744 PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering 1745 1746 Collective 1747 1748 Input Parameters: 1749 + mainSF - the original SF structure 1750 - map - a ISLocalToGlobal mapping that contains the subset of points 1751 1752 Output Parameter: 1753 . subSF - a subset of the mainSF for the desired subset. 1754 1755 Level: intermediate 1756 @*/ 1757 PetscErrorCode PetscSFGetSubSF(PetscSF mainsf,ISLocalToGlobalMapping map,PetscSF *subSF) 1758 { 1759 PetscInt nroots, nleaves, *ilocal_sub; 1760 PetscInt i, *ilocal_map, nroots_sub, nleaves_sub = 0; 1761 PetscInt *local_points, *remote_points; 1762 PetscSFNode *iremote_sub; 1763 const PetscInt *ilocal; 1764 const PetscSFNode *iremote; 1765 1766 PetscFunctionBegin; 1767 PetscCall(PetscSFGetGraph(mainsf,&nroots,&nleaves,&ilocal,&iremote)); 1768 1769 /* Look for leaves that pertain to the subset of points. Get the local ordering */ 1770 PetscCall(PetscMalloc1(nleaves,&ilocal_map)); 1771 PetscCall(ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map)); 1772 for (i = 0; i < nleaves; i++) { 1773 if (ilocal_map[i] != -1) nleaves_sub += 1; 1774 } 1775 /* Re-number ilocal with subset numbering. Need information from roots */ 1776 PetscCall(PetscMalloc2(nroots,&local_points,nroots,&remote_points)); 1777 for (i = 0; i < nroots; i++) local_points[i] = i; 1778 PetscCall(ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points)); 1779 PetscCall(PetscSFBcastBegin(mainsf, MPIU_INT, local_points, remote_points,MPI_REPLACE)); 1780 PetscCall(PetscSFBcastEnd(mainsf, MPIU_INT, local_points, remote_points,MPI_REPLACE)); 1781 /* Fill up graph using local (that is, local to the subset) numbering. */ 1782 PetscCall(PetscMalloc1(nleaves_sub,&ilocal_sub)); 1783 PetscCall(PetscMalloc1(nleaves_sub,&iremote_sub)); 1784 nleaves_sub = 0; 1785 for (i = 0; i < nleaves; i++) { 1786 if (ilocal_map[i] != -1) { 1787 ilocal_sub[nleaves_sub] = ilocal_map[i]; 1788 iremote_sub[nleaves_sub].rank = iremote[i].rank; 1789 iremote_sub[nleaves_sub].index = remote_points[ilocal[i]]; 1790 nleaves_sub += 1; 1791 } 1792 } 1793 PetscCall(PetscFree2(local_points,remote_points)); 1794 PetscCall(ISLocalToGlobalMappingGetSize(map,&nroots_sub)); 1795 1796 /* Create new subSF */ 1797 PetscCall(PetscSFCreate(PETSC_COMM_WORLD,subSF)); 1798 PetscCall(PetscSFSetFromOptions(*subSF)); 1799 PetscCall(PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES)); 1800 PetscCall(PetscFree(ilocal_map)); 1801 PetscCall(PetscFree(iremote_sub)); 1802 PetscFunctionReturn(0); 1803 } 1804 1805 /*@C 1806 DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point 1807 1808 Not Collective 1809 1810 Input Parameters: 1811 + dm - the DMNetwork object 1812 - p - the vertex point 1813 1814 Output Parameters: 1815 + nedges - number of edges connected to this vertex point 1816 - edges - list of edge points 1817 1818 Level: beginner 1819 1820 Fortran Notes: 1821 Since it returns an array, this routine is only available in Fortran 90, and you must 1822 include petsc.h90 in your code. 1823 1824 .seealso: `DMNetworkCreate()`, `DMNetworkGetConnectedVertices()` 1825 @*/ 1826 PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[]) 1827 { 1828 DM_Network *network = (DM_Network*)dm->data; 1829 1830 PetscFunctionBegin; 1831 PetscCall(DMPlexGetSupportSize(network->plex,vertex,nedges)); 1832 if (edges) PetscCall(DMPlexGetSupport(network->plex,vertex,edges)); 1833 PetscFunctionReturn(0); 1834 } 1835 1836 /*@C 1837 DMNetworkGetConnectedVertices - Return the connected vertices for this edge point 1838 1839 Not Collective 1840 1841 Input Parameters: 1842 + dm - the DMNetwork object 1843 - p - the edge point 1844 1845 Output Parameters: 1846 . vertices - vertices connected to this edge 1847 1848 Level: beginner 1849 1850 Fortran Notes: 1851 Since it returns an array, this routine is only available in Fortran 90, and you must 1852 include petsc.h90 in your code. 1853 1854 .seealso: `DMNetworkCreate()`, `DMNetworkGetSupportingEdges()` 1855 @*/ 1856 PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[]) 1857 { 1858 DM_Network *network = (DM_Network*)dm->data; 1859 1860 PetscFunctionBegin; 1861 PetscCall(DMPlexGetCone(network->plex,edge,vertices)); 1862 PetscFunctionReturn(0); 1863 } 1864 1865 /*@ 1866 DMNetworkIsSharedVertex - Returns TRUE if the vertex is shared by subnetworks 1867 1868 Not Collective 1869 1870 Input Parameters: 1871 + dm - the DMNetwork object 1872 - p - the vertex point 1873 1874 Output Parameter: 1875 . flag - TRUE if the vertex is shared by subnetworks 1876 1877 Level: beginner 1878 1879 .seealso: `DMNetworkAddSharedVertices()`, `DMNetworkIsGhostVertex()` 1880 @*/ 1881 PetscErrorCode DMNetworkIsSharedVertex(DM dm,PetscInt p,PetscBool *flag) 1882 { 1883 PetscInt i; 1884 1885 PetscFunctionBegin; 1886 *flag = PETSC_FALSE; 1887 1888 if (dm->setupcalled) { /* DMNetworkGetGlobalVertexIndex() requires DMSetUp() be called */ 1889 DM_Network *network = (DM_Network*)dm->data; 1890 PetscInt gidx; 1891 PetscCall(DMNetworkGetGlobalVertexIndex(dm,p,&gidx)); 1892 PetscCall(PetscTableFind(network->svtable,gidx+1,&i)); 1893 if (i) *flag = PETSC_TRUE; 1894 } else { /* would be removed? */ 1895 PetscInt nv; 1896 const PetscInt *vtx; 1897 PetscCall(DMNetworkGetSharedVertices(dm,&nv,&vtx)); 1898 for (i=0; i<nv; i++) { 1899 if (p == vtx[i]) { 1900 *flag = PETSC_TRUE; 1901 break; 1902 } 1903 } 1904 } 1905 PetscFunctionReturn(0); 1906 } 1907 1908 /*@ 1909 DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex 1910 1911 Not Collective 1912 1913 Input Parameters: 1914 + dm - the DMNetwork object 1915 - p - the vertex point 1916 1917 Output Parameter: 1918 . isghost - TRUE if the vertex is a ghost point 1919 1920 Level: beginner 1921 1922 .seealso: `DMNetworkGetConnectedVertices()`, `DMNetworkGetVertexRange()`, `DMNetworkIsSharedVertex()` 1923 @*/ 1924 PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost) 1925 { 1926 DM_Network *network = (DM_Network*)dm->data; 1927 PetscInt offsetg; 1928 PetscSection sectiong; 1929 1930 PetscFunctionBegin; 1931 *isghost = PETSC_FALSE; 1932 PetscCall(DMGetGlobalSection(network->plex,§iong)); 1933 PetscCall(PetscSectionGetOffset(sectiong,p,&offsetg)); 1934 if (offsetg < 0) *isghost = PETSC_TRUE; 1935 PetscFunctionReturn(0); 1936 } 1937 1938 PetscErrorCode DMSetUp_Network(DM dm) 1939 { 1940 DM_Network *network=(DM_Network*)dm->data; 1941 1942 PetscFunctionBegin; 1943 PetscCall(DMNetworkComponentSetUp(dm)); 1944 PetscCall(DMNetworkVariablesSetUp(dm)); 1945 1946 PetscCall(DMSetLocalSection(network->plex,network->DofSection)); 1947 PetscCall(DMGetGlobalSection(network->plex,&network->GlobalDofSection)); 1948 1949 dm->setupcalled = PETSC_TRUE; 1950 1951 /* View dmnetwork */ 1952 PetscCall(DMViewFromOptions(dm,NULL,"-dmnetwork_view")); 1953 PetscFunctionReturn(0); 1954 } 1955 1956 /*@ 1957 DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices 1958 -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)? 1959 1960 Collective 1961 1962 Input Parameters: 1963 + dm - the DMNetwork object 1964 . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges 1965 - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices 1966 1967 Level: intermediate 1968 1969 @*/ 1970 PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg) 1971 { 1972 DM_Network *network=(DM_Network*)dm->data; 1973 PetscInt nVertices = network->nVertices; 1974 1975 PetscFunctionBegin; 1976 network->userEdgeJacobian = eflg; 1977 network->userVertexJacobian = vflg; 1978 1979 if (eflg && !network->Je) { 1980 PetscCall(PetscCalloc1(3*network->nEdges,&network->Je)); 1981 } 1982 1983 if (vflg && !network->Jv && nVertices) { 1984 PetscInt i,*vptr,nedges,vStart=network->vStart; 1985 PetscInt nedges_total; 1986 const PetscInt *edges; 1987 1988 /* count nvertex_total */ 1989 nedges_total = 0; 1990 PetscCall(PetscMalloc1(nVertices+1,&vptr)); 1991 1992 vptr[0] = 0; 1993 for (i=0; i<nVertices; i++) { 1994 PetscCall(DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges)); 1995 nedges_total += nedges; 1996 vptr[i+1] = vptr[i] + 2*nedges + 1; 1997 } 1998 1999 PetscCall(PetscCalloc1(2*nedges_total+nVertices,&network->Jv)); 2000 network->Jvptr = vptr; 2001 } 2002 PetscFunctionReturn(0); 2003 } 2004 2005 /*@ 2006 DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network 2007 2008 Not Collective 2009 2010 Input Parameters: 2011 + dm - the DMNetwork object 2012 . p - the edge point 2013 - J - array (size = 3) of Jacobian submatrices for this edge point: 2014 J[0]: this edge 2015 J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices() 2016 2017 Level: advanced 2018 2019 .seealso: `DMNetworkVertexSetMatrix()` 2020 @*/ 2021 PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[]) 2022 { 2023 DM_Network *network=(DM_Network*)dm->data; 2024 2025 PetscFunctionBegin; 2026 PetscCheck(network->Je,PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix"); 2027 2028 if (J) { 2029 network->Je[3*p] = J[0]; 2030 network->Je[3*p+1] = J[1]; 2031 network->Je[3*p+2] = J[2]; 2032 } 2033 PetscFunctionReturn(0); 2034 } 2035 2036 /*@ 2037 DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network 2038 2039 Not Collective 2040 2041 Input Parameters: 2042 + dm - The DMNetwork object 2043 . p - the vertex point 2044 - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point: 2045 J[0]: this vertex 2046 J[1+2*i]: i-th supporting edge 2047 J[1+2*i+1]: i-th connected vertex 2048 2049 Level: advanced 2050 2051 .seealso: `DMNetworkEdgeSetMatrix()` 2052 @*/ 2053 PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[]) 2054 { 2055 DM_Network *network=(DM_Network*)dm->data; 2056 PetscInt i,*vptr,nedges,vStart=network->vStart; 2057 const PetscInt *edges; 2058 2059 PetscFunctionBegin; 2060 PetscCheck(network->Jv,PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix"); 2061 2062 if (J) { 2063 vptr = network->Jvptr; 2064 network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */ 2065 2066 /* Set Jacobian for each supporting edge and connected vertex */ 2067 PetscCall(DMNetworkGetSupportingEdges(dm,p,&nedges,&edges)); 2068 for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i]; 2069 } 2070 PetscFunctionReturn(0); 2071 } 2072 2073 static inline PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 2074 { 2075 PetscInt j; 2076 PetscScalar val=(PetscScalar)ncols; 2077 2078 PetscFunctionBegin; 2079 if (!ghost) { 2080 for (j=0; j<nrows; j++) { 2081 PetscCall(VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES)); 2082 } 2083 } else { 2084 for (j=0; j<nrows; j++) { 2085 PetscCall(VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES)); 2086 } 2087 } 2088 PetscFunctionReturn(0); 2089 } 2090 2091 static inline PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 2092 { 2093 PetscInt j,ncols_u; 2094 PetscScalar val; 2095 2096 PetscFunctionBegin; 2097 if (!ghost) { 2098 for (j=0; j<nrows; j++) { 2099 PetscCall(MatGetRow(Ju,j,&ncols_u,NULL,NULL)); 2100 val = (PetscScalar)ncols_u; 2101 PetscCall(VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES)); 2102 PetscCall(MatRestoreRow(Ju,j,&ncols_u,NULL,NULL)); 2103 } 2104 } else { 2105 for (j=0; j<nrows; j++) { 2106 PetscCall(MatGetRow(Ju,j,&ncols_u,NULL,NULL)); 2107 val = (PetscScalar)ncols_u; 2108 PetscCall(VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES)); 2109 PetscCall(MatRestoreRow(Ju,j,&ncols_u,NULL,NULL)); 2110 } 2111 } 2112 PetscFunctionReturn(0); 2113 } 2114 2115 static inline PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 2116 { 2117 PetscFunctionBegin; 2118 if (Ju) { 2119 PetscCall(MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz)); 2120 } else { 2121 PetscCall(MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz)); 2122 } 2123 PetscFunctionReturn(0); 2124 } 2125 2126 static inline PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 2127 { 2128 PetscInt j,*cols; 2129 PetscScalar *zeros; 2130 2131 PetscFunctionBegin; 2132 PetscCall(PetscCalloc2(ncols,&cols,nrows*ncols,&zeros)); 2133 for (j=0; j<ncols; j++) cols[j] = j+ cstart; 2134 PetscCall(MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES)); 2135 PetscCall(PetscFree2(cols,zeros)); 2136 PetscFunctionReturn(0); 2137 } 2138 2139 static inline PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 2140 { 2141 PetscInt j,M,N,row,col,ncols_u; 2142 const PetscInt *cols; 2143 PetscScalar zero=0.0; 2144 2145 PetscFunctionBegin; 2146 PetscCall(MatGetSize(Ju,&M,&N)); 2147 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); 2148 2149 for (row=0; row<nrows; row++) { 2150 PetscCall(MatGetRow(Ju,row,&ncols_u,&cols,NULL)); 2151 for (j=0; j<ncols_u; j++) { 2152 col = cols[j] + cstart; 2153 PetscCall(MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES)); 2154 } 2155 PetscCall(MatRestoreRow(Ju,row,&ncols_u,&cols,NULL)); 2156 } 2157 PetscFunctionReturn(0); 2158 } 2159 2160 static inline PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 2161 { 2162 PetscFunctionBegin; 2163 if (Ju) { 2164 PetscCall(MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J)); 2165 } else { 2166 PetscCall(MatSetDenseblock_private(nrows,rows,ncols,cstart,J)); 2167 } 2168 PetscFunctionReturn(0); 2169 } 2170 2171 /* 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. 2172 */ 2173 PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog) 2174 { 2175 PetscInt i,size,dof; 2176 PetscInt *glob2loc; 2177 2178 PetscFunctionBegin; 2179 PetscCall(PetscSectionGetStorageSize(localsec,&size)); 2180 PetscCall(PetscMalloc1(size,&glob2loc)); 2181 2182 for (i = 0; i < size; i++) { 2183 PetscCall(PetscSectionGetOffset(globalsec,i,&dof)); 2184 dof = (dof >= 0) ? dof : -(dof + 1); 2185 glob2loc[i] = dof; 2186 } 2187 2188 PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog)); 2189 #if 0 2190 PetscCall(PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD)); 2191 #endif 2192 PetscFunctionReturn(0); 2193 } 2194 2195 #include <petsc/private/matimpl.h> 2196 2197 PetscErrorCode DMCreateMatrix_Network_Nest(DM dm,Mat *J) 2198 { 2199 DM_Network *network = (DM_Network*)dm->data; 2200 PetscInt eDof,vDof; 2201 Mat j11,j12,j21,j22,bA[2][2]; 2202 MPI_Comm comm; 2203 ISLocalToGlobalMapping eISMap,vISMap; 2204 2205 PetscFunctionBegin; 2206 PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 2207 2208 PetscCall(PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof)); 2209 PetscCall(PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof)); 2210 2211 PetscCall(MatCreate(comm, &j11)); 2212 PetscCall(MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE)); 2213 PetscCall(MatSetType(j11, MATMPIAIJ)); 2214 2215 PetscCall(MatCreate(comm, &j12)); 2216 PetscCall(MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE)); 2217 PetscCall(MatSetType(j12, MATMPIAIJ)); 2218 2219 PetscCall(MatCreate(comm, &j21)); 2220 PetscCall(MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE)); 2221 PetscCall(MatSetType(j21, MATMPIAIJ)); 2222 2223 PetscCall(MatCreate(comm, &j22)); 2224 PetscCall(MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE)); 2225 PetscCall(MatSetType(j22, MATMPIAIJ)); 2226 2227 bA[0][0] = j11; 2228 bA[0][1] = j12; 2229 bA[1][0] = j21; 2230 bA[1][1] = j22; 2231 2232 PetscCall(CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap)); 2233 PetscCall(CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap)); 2234 2235 PetscCall(MatSetLocalToGlobalMapping(j11,eISMap,eISMap)); 2236 PetscCall(MatSetLocalToGlobalMapping(j12,eISMap,vISMap)); 2237 PetscCall(MatSetLocalToGlobalMapping(j21,vISMap,eISMap)); 2238 PetscCall(MatSetLocalToGlobalMapping(j22,vISMap,vISMap)); 2239 2240 PetscCall(MatSetUp(j11)); 2241 PetscCall(MatSetUp(j12)); 2242 PetscCall(MatSetUp(j21)); 2243 PetscCall(MatSetUp(j22)); 2244 2245 PetscCall(MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J)); 2246 PetscCall(MatSetUp(*J)); 2247 PetscCall(MatNestSetVecType(*J,VECNEST)); 2248 PetscCall(MatDestroy(&j11)); 2249 PetscCall(MatDestroy(&j12)); 2250 PetscCall(MatDestroy(&j21)); 2251 PetscCall(MatDestroy(&j22)); 2252 2253 PetscCall(MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY)); 2254 PetscCall(MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY)); 2255 PetscCall(MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE)); 2256 2257 /* Free structures */ 2258 PetscCall(ISLocalToGlobalMappingDestroy(&eISMap)); 2259 PetscCall(ISLocalToGlobalMappingDestroy(&vISMap)); 2260 PetscFunctionReturn(0); 2261 } 2262 2263 PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J) 2264 { 2265 DM_Network *network = (DM_Network*)dm->data; 2266 PetscInt eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize; 2267 PetscInt cstart,ncols,j,e,v; 2268 PetscBool ghost,ghost_vc,ghost2,isNest; 2269 Mat Juser; 2270 PetscSection sectionGlobal; 2271 PetscInt nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */ 2272 const PetscInt *edges,*cone; 2273 MPI_Comm comm; 2274 MatType mtype; 2275 Vec vd_nz,vo_nz; 2276 PetscInt *dnnz,*onnz; 2277 PetscScalar *vdnz,*vonz; 2278 2279 PetscFunctionBegin; 2280 mtype = dm->mattype; 2281 PetscCall(PetscStrcmp(mtype,MATNEST,&isNest)); 2282 if (isNest) { 2283 PetscCall(DMCreateMatrix_Network_Nest(dm,J)); 2284 PetscCall(MatSetDM(*J,dm)); 2285 PetscFunctionReturn(0); 2286 } 2287 2288 if (!network->userEdgeJacobian && !network->userVertexJacobian) { 2289 /* user does not provide Jacobian blocks */ 2290 PetscCall(DMCreateMatrix_Plex(network->plex,J)); 2291 PetscCall(MatSetDM(*J,dm)); 2292 PetscFunctionReturn(0); 2293 } 2294 2295 PetscCall(MatCreate(PetscObjectComm((PetscObject)dm),J)); 2296 PetscCall(DMGetGlobalSection(network->plex,§ionGlobal)); 2297 PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize)); 2298 PetscCall(MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE)); 2299 2300 PetscCall(MatSetType(*J,MATAIJ)); 2301 PetscCall(MatSetFromOptions(*J)); 2302 2303 /* (1) Set matrix preallocation */ 2304 /*------------------------------*/ 2305 PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 2306 PetscCall(VecCreate(comm,&vd_nz)); 2307 PetscCall(VecSetSizes(vd_nz,localSize,PETSC_DECIDE)); 2308 PetscCall(VecSetFromOptions(vd_nz)); 2309 PetscCall(VecSet(vd_nz,0.0)); 2310 PetscCall(VecDuplicate(vd_nz,&vo_nz)); 2311 2312 /* Set preallocation for edges */ 2313 /*-----------------------------*/ 2314 PetscCall(DMNetworkGetEdgeRange(dm,&eStart,&eEnd)); 2315 2316 PetscCall(PetscMalloc1(localSize,&rows)); 2317 for (e=eStart; e<eEnd; e++) { 2318 /* Get row indices */ 2319 PetscCall(DMNetworkGetGlobalVecOffset(dm,e,ALL_COMPONENTS,&rstart)); 2320 PetscCall(PetscSectionGetDof(network->DofSection,e,&nrows)); 2321 if (nrows) { 2322 for (j=0; j<nrows; j++) rows[j] = j + rstart; 2323 2324 /* Set preallocation for connected vertices */ 2325 PetscCall(DMNetworkGetConnectedVertices(dm,e,&cone)); 2326 for (v=0; v<2; v++) { 2327 PetscCall(PetscSectionGetDof(network->DofSection,cone[v],&ncols)); 2328 2329 if (network->Je) { 2330 Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 2331 } else Juser = NULL; 2332 PetscCall(DMNetworkIsGhostVertex(dm,cone[v],&ghost)); 2333 PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz)); 2334 } 2335 2336 /* Set preallocation for edge self */ 2337 cstart = rstart; 2338 if (network->Je) { 2339 Juser = network->Je[3*e]; /* Jacobian(e,e) */ 2340 } else Juser = NULL; 2341 PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz)); 2342 } 2343 } 2344 2345 /* Set preallocation for vertices */ 2346 /*--------------------------------*/ 2347 PetscCall(DMNetworkGetVertexRange(dm,&vStart,&vEnd)); 2348 if (vEnd - vStart) vptr = network->Jvptr; 2349 2350 for (v=vStart; v<vEnd; v++) { 2351 /* Get row indices */ 2352 PetscCall(DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&rstart)); 2353 PetscCall(PetscSectionGetDof(network->DofSection,v,&nrows)); 2354 if (!nrows) continue; 2355 2356 PetscCall(DMNetworkIsGhostVertex(dm,v,&ghost)); 2357 if (ghost) { 2358 PetscCall(PetscMalloc1(nrows,&rows_v)); 2359 } else { 2360 rows_v = rows; 2361 } 2362 2363 for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 2364 2365 /* Get supporting edges and connected vertices */ 2366 PetscCall(DMNetworkGetSupportingEdges(dm,v,&nedges,&edges)); 2367 2368 for (e=0; e<nedges; e++) { 2369 /* Supporting edges */ 2370 PetscCall(DMNetworkGetGlobalVecOffset(dm,edges[e],ALL_COMPONENTS,&cstart)); 2371 PetscCall(PetscSectionGetDof(network->DofSection,edges[e],&ncols)); 2372 2373 if (network->Jv) { 2374 Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 2375 } else Juser = NULL; 2376 PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz)); 2377 2378 /* Connected vertices */ 2379 PetscCall(DMNetworkGetConnectedVertices(dm,edges[e],&cone)); 2380 vc = (v == cone[0]) ? cone[1]:cone[0]; 2381 PetscCall(DMNetworkIsGhostVertex(dm,vc,&ghost_vc)); 2382 2383 PetscCall(PetscSectionGetDof(network->DofSection,vc,&ncols)); 2384 2385 if (network->Jv) { 2386 Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 2387 } else Juser = NULL; 2388 if (ghost_vc||ghost) { 2389 ghost2 = PETSC_TRUE; 2390 } else { 2391 ghost2 = PETSC_FALSE; 2392 } 2393 PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz)); 2394 } 2395 2396 /* Set preallocation for vertex self */ 2397 PetscCall(DMNetworkIsGhostVertex(dm,v,&ghost)); 2398 if (!ghost) { 2399 PetscCall(DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&cstart)); 2400 if (network->Jv) { 2401 Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 2402 } else Juser = NULL; 2403 PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz)); 2404 } 2405 if (ghost) PetscCall(PetscFree(rows_v)); 2406 } 2407 2408 PetscCall(VecAssemblyBegin(vd_nz)); 2409 PetscCall(VecAssemblyBegin(vo_nz)); 2410 2411 PetscCall(PetscMalloc2(localSize,&dnnz,localSize,&onnz)); 2412 2413 PetscCall(VecAssemblyEnd(vd_nz)); 2414 PetscCall(VecAssemblyEnd(vo_nz)); 2415 2416 PetscCall(VecGetArray(vd_nz,&vdnz)); 2417 PetscCall(VecGetArray(vo_nz,&vonz)); 2418 for (j=0; j<localSize; j++) { 2419 dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]); 2420 onnz[j] = (PetscInt)PetscRealPart(vonz[j]); 2421 } 2422 PetscCall(VecRestoreArray(vd_nz,&vdnz)); 2423 PetscCall(VecRestoreArray(vo_nz,&vonz)); 2424 PetscCall(VecDestroy(&vd_nz)); 2425 PetscCall(VecDestroy(&vo_nz)); 2426 2427 PetscCall(MatSeqAIJSetPreallocation(*J,0,dnnz)); 2428 PetscCall(MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz)); 2429 PetscCall(MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE)); 2430 2431 PetscCall(PetscFree2(dnnz,onnz)); 2432 2433 /* (2) Set matrix entries for edges */ 2434 /*----------------------------------*/ 2435 for (e=eStart; e<eEnd; e++) { 2436 /* Get row indices */ 2437 PetscCall(DMNetworkGetGlobalVecOffset(dm,e,ALL_COMPONENTS,&rstart)); 2438 PetscCall(PetscSectionGetDof(network->DofSection,e,&nrows)); 2439 if (nrows) { 2440 for (j=0; j<nrows; j++) rows[j] = j + rstart; 2441 2442 /* Set matrix entries for connected vertices */ 2443 PetscCall(DMNetworkGetConnectedVertices(dm,e,&cone)); 2444 for (v=0; v<2; v++) { 2445 PetscCall(DMNetworkGetGlobalVecOffset(dm,cone[v],ALL_COMPONENTS,&cstart)); 2446 PetscCall(PetscSectionGetDof(network->DofSection,cone[v],&ncols)); 2447 2448 if (network->Je) { 2449 Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 2450 } else Juser = NULL; 2451 PetscCall(MatSetblock_private(Juser,nrows,rows,ncols,cstart,J)); 2452 } 2453 2454 /* Set matrix entries for edge self */ 2455 cstart = rstart; 2456 if (network->Je) { 2457 Juser = network->Je[3*e]; /* Jacobian(e,e) */ 2458 } else Juser = NULL; 2459 PetscCall(MatSetblock_private(Juser,nrows,rows,nrows,cstart,J)); 2460 } 2461 } 2462 2463 /* Set matrix entries for vertices */ 2464 /*---------------------------------*/ 2465 for (v=vStart; v<vEnd; v++) { 2466 /* Get row indices */ 2467 PetscCall(DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&rstart)); 2468 PetscCall(PetscSectionGetDof(network->DofSection,v,&nrows)); 2469 if (!nrows) continue; 2470 2471 PetscCall(DMNetworkIsGhostVertex(dm,v,&ghost)); 2472 if (ghost) { 2473 PetscCall(PetscMalloc1(nrows,&rows_v)); 2474 } else { 2475 rows_v = rows; 2476 } 2477 for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 2478 2479 /* Get supporting edges and connected vertices */ 2480 PetscCall(DMNetworkGetSupportingEdges(dm,v,&nedges,&edges)); 2481 2482 for (e=0; e<nedges; e++) { 2483 /* Supporting edges */ 2484 PetscCall(DMNetworkGetGlobalVecOffset(dm,edges[e],ALL_COMPONENTS,&cstart)); 2485 PetscCall(PetscSectionGetDof(network->DofSection,edges[e],&ncols)); 2486 2487 if (network->Jv) { 2488 Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 2489 } else Juser = NULL; 2490 PetscCall(MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J)); 2491 2492 /* Connected vertices */ 2493 PetscCall(DMNetworkGetConnectedVertices(dm,edges[e],&cone)); 2494 vc = (v == cone[0]) ? cone[1]:cone[0]; 2495 2496 PetscCall(DMNetworkGetGlobalVecOffset(dm,vc,ALL_COMPONENTS,&cstart)); 2497 PetscCall(PetscSectionGetDof(network->DofSection,vc,&ncols)); 2498 2499 if (network->Jv) { 2500 Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 2501 } else Juser = NULL; 2502 PetscCall(MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J)); 2503 } 2504 2505 /* Set matrix entries for vertex self */ 2506 if (!ghost) { 2507 PetscCall(DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&cstart)); 2508 if (network->Jv) { 2509 Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 2510 } else Juser = NULL; 2511 PetscCall(MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J)); 2512 } 2513 if (ghost) PetscCall(PetscFree(rows_v)); 2514 } 2515 PetscCall(PetscFree(rows)); 2516 2517 PetscCall(MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY)); 2518 PetscCall(MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY)); 2519 2520 PetscCall(MatSetDM(*J,dm)); 2521 PetscFunctionReturn(0); 2522 } 2523 2524 PetscErrorCode DMDestroy_Network(DM dm) 2525 { 2526 DM_Network *network = (DM_Network*)dm->data; 2527 PetscInt j,np; 2528 2529 PetscFunctionBegin; 2530 if (--network->refct > 0) PetscFunctionReturn(0); 2531 PetscCall(PetscFree(network->Je)); 2532 if (network->Jv) { 2533 PetscCall(PetscFree(network->Jvptr)); 2534 PetscCall(PetscFree(network->Jv)); 2535 } 2536 2537 PetscCall(ISLocalToGlobalMappingDestroy(&network->vertex.mapping)); 2538 PetscCall(PetscSectionDestroy(&network->vertex.DofSection)); 2539 PetscCall(PetscSectionDestroy(&network->vertex.GlobalDofSection)); 2540 PetscCall(PetscFree(network->vltog)); 2541 PetscCall(PetscSFDestroy(&network->vertex.sf)); 2542 /* edge */ 2543 PetscCall(ISLocalToGlobalMappingDestroy(&network->edge.mapping)); 2544 PetscCall(PetscSectionDestroy(&network->edge.DofSection)); 2545 PetscCall(PetscSectionDestroy(&network->edge.GlobalDofSection)); 2546 PetscCall(PetscSFDestroy(&network->edge.sf)); 2547 PetscCall(DMDestroy(&network->plex)); 2548 PetscCall(PetscSectionDestroy(&network->DataSection)); 2549 PetscCall(PetscSectionDestroy(&network->DofSection)); 2550 2551 for (j=0; j<network->Nsvtx; j++) PetscCall(PetscFree(network->svtx[j].sv)); 2552 PetscCall(PetscFree(network->svtx)); 2553 PetscCall(PetscFree2(network->subnetedge,network->subnetvtx)); 2554 2555 PetscCall(PetscTableDestroy(&network->svtable)); 2556 PetscCall(PetscFree(network->subnet)); 2557 PetscCall(PetscFree(network->component)); 2558 PetscCall(PetscFree(network->componentdataarray)); 2559 2560 if (network->header) { 2561 np = network->pEnd - network->pStart; 2562 for (j=0; j < np; j++) { 2563 PetscCall(PetscFree5(network->header[j].size,network->header[j].key,network->header[j].offset,network->header[j].nvar,network->header[j].offsetvarrel)); 2564 PetscCall(PetscFree(network->cvalue[j].data)); 2565 } 2566 PetscCall(PetscFree2(network->header,network->cvalue)); 2567 } 2568 PetscCall(PetscFree(network)); 2569 PetscFunctionReturn(0); 2570 } 2571 2572 PetscErrorCode DMView_Network(DM dm,PetscViewer viewer) 2573 { 2574 PetscBool iascii; 2575 PetscMPIInt rank; 2576 2577 PetscFunctionBegin; 2578 PetscValidHeaderSpecific(dm,DM_CLASSID, 1); 2579 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2580 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank)); 2581 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 2582 if (iascii) { 2583 const PetscInt *cone,*vtx,*edges; 2584 PetscInt vfrom,vto,i,j,nv,ne,nsv,p,nsubnet; 2585 DM_Network *network = (DM_Network*)dm->data; 2586 2587 nsubnet = network->Nsubnet; /* num of subnetworks */ 2588 if (rank == 0) { 2589 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)); 2590 } 2591 2592 PetscCall(DMNetworkGetSharedVertices(dm,&nsv,NULL)); 2593 PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 2594 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d] nEdges: %" PetscInt_FMT "; nVertices: %" PetscInt_FMT "; nSharedVertices: %" PetscInt_FMT "\n",rank,network->nEdges,network->nVertices,nsv)); 2595 2596 for (i=0; i<nsubnet; i++) { 2597 PetscCall(DMNetworkGetSubnetwork(dm,i,&nv,&ne,&vtx,&edges)); 2598 if (ne) { 2599 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " Subnet %" PetscInt_FMT ": nEdges %" PetscInt_FMT ", nVertices(include shared vertices) %" PetscInt_FMT "\n",i,ne,nv)); 2600 for (j=0; j<ne; j++) { 2601 p = edges[j]; 2602 PetscCall(DMNetworkGetConnectedVertices(dm,p,&cone)); 2603 PetscCall(DMNetworkGetGlobalVertexIndex(dm,cone[0],&vfrom)); 2604 PetscCall(DMNetworkGetGlobalVertexIndex(dm,cone[1],&vto)); 2605 PetscCall(DMNetworkGetGlobalEdgeIndex(dm,edges[j],&p)); 2606 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " edge %" PetscInt_FMT ": %" PetscInt_FMT " ----> %" PetscInt_FMT "\n",p,vfrom,vto)); 2607 } 2608 } 2609 } 2610 2611 /* Shared vertices */ 2612 PetscCall(DMNetworkGetSharedVertices(dm,NULL,&vtx)); 2613 if (nsv) { 2614 PetscInt gidx; 2615 PetscBool ghost; 2616 const PetscInt *sv=NULL; 2617 2618 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " SharedVertices:\n")); 2619 for (i=0; i<nsv; i++) { 2620 PetscCall(DMNetworkIsGhostVertex(dm,vtx[i],&ghost)); 2621 if (ghost) continue; 2622 2623 PetscCall(DMNetworkSharedVertexGetInfo(dm,vtx[i],&gidx,&nv,&sv)); 2624 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " svtx %" PetscInt_FMT ": global index %" PetscInt_FMT ", subnet[%" PetscInt_FMT "].%" PetscInt_FMT " ---->\n",i,gidx,sv[0],sv[1])); 2625 for (j=1; j<nv; j++) { 2626 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " ----> subnet[%" PetscInt_FMT "].%" PetscInt_FMT "\n",sv[2*j],sv[2*j+1])); 2627 } 2628 } 2629 } 2630 PetscCall(PetscViewerFlush(viewer)); 2631 PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2632 } else PetscCheck(iascii,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Viewer type %s not yet supported for DMNetwork writing",((PetscObject)viewer)->type_name); 2633 PetscFunctionReturn(0); 2634 } 2635 2636 PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l) 2637 { 2638 DM_Network *network = (DM_Network*)dm->data; 2639 2640 PetscFunctionBegin; 2641 PetscCall(DMGlobalToLocalBegin(network->plex,g,mode,l)); 2642 PetscFunctionReturn(0); 2643 } 2644 2645 PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l) 2646 { 2647 DM_Network *network = (DM_Network*)dm->data; 2648 2649 PetscFunctionBegin; 2650 PetscCall(DMGlobalToLocalEnd(network->plex,g,mode,l)); 2651 PetscFunctionReturn(0); 2652 } 2653 2654 PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g) 2655 { 2656 DM_Network *network = (DM_Network*)dm->data; 2657 2658 PetscFunctionBegin; 2659 PetscCall(DMLocalToGlobalBegin(network->plex,l,mode,g)); 2660 PetscFunctionReturn(0); 2661 } 2662 2663 PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g) 2664 { 2665 DM_Network *network = (DM_Network*)dm->data; 2666 2667 PetscFunctionBegin; 2668 PetscCall(DMLocalToGlobalEnd(network->plex,l,mode,g)); 2669 PetscFunctionReturn(0); 2670 } 2671 2672 /*@ 2673 DMNetworkGetVertexLocalToGlobalOrdering - Get vertex global index 2674 2675 Not collective 2676 2677 Input Parameters: 2678 + dm - the dm object 2679 - vloc - local vertex ordering, start from 0 2680 2681 Output Parameters: 2682 . vg - global vertex ordering, start from 0 2683 2684 Level: advanced 2685 2686 .seealso: `DMNetworkSetVertexLocalToGlobalOrdering()` 2687 @*/ 2688 PetscErrorCode DMNetworkGetVertexLocalToGlobalOrdering(DM dm,PetscInt vloc,PetscInt *vg) 2689 { 2690 DM_Network *network = (DM_Network*)dm->data; 2691 PetscInt *vltog = network->vltog; 2692 2693 PetscFunctionBegin; 2694 PetscCheck(vltog,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Must call DMNetworkSetVertexLocalToGlobalOrdering() first"); 2695 *vg = vltog[vloc]; 2696 PetscFunctionReturn(0); 2697 } 2698 2699 /*@ 2700 DMNetworkSetVertexLocalToGlobalOrdering - Create and setup vertex local to global map 2701 2702 Collective 2703 2704 Input Parameters: 2705 . dm - the dm object 2706 2707 Level: advanced 2708 2709 .seealso: `DMNetworkGetGlobalVertexIndex()` 2710 @*/ 2711 PetscErrorCode DMNetworkSetVertexLocalToGlobalOrdering(DM dm) 2712 { 2713 DM_Network *network=(DM_Network*)dm->data; 2714 MPI_Comm comm; 2715 PetscMPIInt rank,size,*displs=NULL,*recvcounts=NULL,remoterank; 2716 PetscBool ghost; 2717 PetscInt *vltog,nroots,nleaves,i,*vrange,k,N,lidx; 2718 const PetscSFNode *iremote; 2719 PetscSF vsf; 2720 Vec Vleaves,Vleaves_seq; 2721 VecScatter ctx; 2722 PetscScalar *varr,val; 2723 const PetscScalar *varr_read; 2724 2725 PetscFunctionBegin; 2726 PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 2727 PetscCallMPI(MPI_Comm_size(comm,&size)); 2728 PetscCallMPI(MPI_Comm_rank(comm,&rank)); 2729 2730 if (size == 1) { 2731 nroots = network->vEnd - network->vStart; 2732 PetscCall(PetscMalloc1(nroots, &vltog)); 2733 for (i=0; i<nroots; i++) vltog[i] = i; 2734 network->vltog = vltog; 2735 PetscFunctionReturn(0); 2736 } 2737 2738 PetscCheck(network->distributecalled,comm, PETSC_ERR_ARG_WRONGSTATE,"Must call DMNetworkDistribute() first"); 2739 if (network->vltog) PetscCall(PetscFree(network->vltog)); 2740 2741 PetscCall(DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping)); 2742 PetscCall(PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf)); 2743 vsf = network->vertex.sf; 2744 2745 PetscCall(PetscMalloc3(size+1,&vrange,size,&displs,size,&recvcounts)); 2746 PetscCall(PetscSFGetGraph(vsf,&nroots,&nleaves,NULL,&iremote)); 2747 2748 for (i=0; i<size; i++) { displs[i] = i; recvcounts[i] = 1;} 2749 2750 i = nroots - nleaves; /* local number of vertices, excluding ghosts */ 2751 vrange[0] = 0; 2752 PetscCallMPI(MPI_Allgatherv(&i,1,MPIU_INT,vrange+1,recvcounts,displs,MPIU_INT,comm)); 2753 for (i=2; i<size+1; i++) {vrange[i] += vrange[i-1];} 2754 2755 PetscCall(PetscMalloc1(nroots, &vltog)); 2756 network->vltog = vltog; 2757 2758 /* Set vltog for non-ghost vertices */ 2759 k = 0; 2760 for (i=0; i<nroots; i++) { 2761 PetscCall(DMNetworkIsGhostVertex(dm,i+network->vStart,&ghost)); 2762 if (ghost) continue; 2763 vltog[i] = vrange[rank] + k++; 2764 } 2765 PetscCall(PetscFree3(vrange,displs,recvcounts)); 2766 2767 /* Set vltog for ghost vertices */ 2768 /* (a) create parallel Vleaves and sequential Vleaves_seq to convert local iremote[*].index to global index */ 2769 PetscCall(VecCreate(comm,&Vleaves)); 2770 PetscCall(VecSetSizes(Vleaves,2*nleaves,PETSC_DETERMINE)); 2771 PetscCall(VecSetFromOptions(Vleaves)); 2772 PetscCall(VecGetArray(Vleaves,&varr)); 2773 for (i=0; i<nleaves; i++) { 2774 varr[2*i] = (PetscScalar)(iremote[i].rank); /* rank of remote process */ 2775 varr[2*i+1] = (PetscScalar)(iremote[i].index); /* local index in remote process */ 2776 } 2777 PetscCall(VecRestoreArray(Vleaves,&varr)); 2778 2779 /* (b) scatter local info to remote processes via VecScatter() */ 2780 PetscCall(VecScatterCreateToAll(Vleaves,&ctx,&Vleaves_seq)); 2781 PetscCall(VecScatterBegin(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD)); 2782 PetscCall(VecScatterEnd(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD)); 2783 2784 /* (c) convert local indices to global indices in parallel vector Vleaves */ 2785 PetscCall(VecGetSize(Vleaves_seq,&N)); 2786 PetscCall(VecGetArrayRead(Vleaves_seq,&varr_read)); 2787 for (i=0; i<N; i+=2) { 2788 remoterank = (PetscMPIInt)PetscRealPart(varr_read[i]); 2789 if (remoterank == rank) { 2790 k = i+1; /* row number */ 2791 lidx = (PetscInt)PetscRealPart(varr_read[i+1]); 2792 val = (PetscScalar)vltog[lidx]; /* global index for non-ghost vertex computed above */ 2793 PetscCall(VecSetValues(Vleaves,1,&k,&val,INSERT_VALUES)); 2794 } 2795 } 2796 PetscCall(VecRestoreArrayRead(Vleaves_seq,&varr_read)); 2797 PetscCall(VecAssemblyBegin(Vleaves)); 2798 PetscCall(VecAssemblyEnd(Vleaves)); 2799 2800 /* (d) Set vltog for ghost vertices by copying local values of Vleaves */ 2801 PetscCall(VecGetArrayRead(Vleaves,&varr_read)); 2802 k = 0; 2803 for (i=0; i<nroots; i++) { 2804 PetscCall(DMNetworkIsGhostVertex(dm,i+network->vStart,&ghost)); 2805 if (!ghost) continue; 2806 vltog[i] = (PetscInt)PetscRealPart(varr_read[2*k+1]); k++; 2807 } 2808 PetscCall(VecRestoreArrayRead(Vleaves,&varr_read)); 2809 2810 PetscCall(VecDestroy(&Vleaves)); 2811 PetscCall(VecDestroy(&Vleaves_seq)); 2812 PetscCall(VecScatterDestroy(&ctx)); 2813 PetscFunctionReturn(0); 2814 } 2815 2816 static inline PetscErrorCode DMISAddSize_private(DM_Network *network,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *nidx) 2817 { 2818 PetscInt i,j,ncomps,nvar,key,offset=0; 2819 DMNetworkComponentHeader header; 2820 2821 PetscFunctionBegin; 2822 PetscCall(PetscSectionGetOffset(network->DataSection,p,&offset)); 2823 ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 2824 header = (DMNetworkComponentHeader)(network->componentdataarray+offset); 2825 2826 for (i=0; i<ncomps; i++) { 2827 key = header->key[i]; 2828 nvar = header->nvar[i]; 2829 for (j=0; j<numkeys; j++) { 2830 if (key == keys[j]) { 2831 if (!blocksize || blocksize[j] == -1) { 2832 *nidx += nvar; 2833 } else { 2834 *nidx += nselectedvar[j]*nvar/blocksize[j]; 2835 } 2836 } 2837 } 2838 } 2839 PetscFunctionReturn(0); 2840 } 2841 2842 static inline PetscErrorCode DMISComputeIdx_private(DM dm,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],PetscInt *ii,PetscInt *idx) 2843 { 2844 PetscInt i,j,ncomps,nvar,key,offsetg,k,k1,offset=0; 2845 DM_Network *network = (DM_Network*)dm->data; 2846 DMNetworkComponentHeader header; 2847 2848 PetscFunctionBegin; 2849 PetscCall(PetscSectionGetOffset(network->DataSection,p,&offset)); 2850 ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 2851 header = (DMNetworkComponentHeader)(network->componentdataarray+offset); 2852 2853 for (i=0; i<ncomps; i++) { 2854 key = header->key[i]; 2855 nvar = header->nvar[i]; 2856 for (j=0; j<numkeys; j++) { 2857 if (key != keys[j]) continue; 2858 2859 PetscCall(DMNetworkGetGlobalVecOffset(dm,p,i,&offsetg)); 2860 if (!blocksize || blocksize[j] == -1) { 2861 for (k=0; k<nvar; k++) idx[(*ii)++] = offsetg + k; 2862 } else { 2863 for (k=0; k<nvar; k+=blocksize[j]) { 2864 for (k1=0; k1<nselectedvar[j]; k1++) idx[(*ii)++] = offsetg + k + selectedvar[j][k1]; 2865 } 2866 } 2867 } 2868 } 2869 PetscFunctionReturn(0); 2870 } 2871 2872 /*@ 2873 DMNetworkCreateIS - Create an index set object from the global vector of the network 2874 2875 Collective 2876 2877 Input Parameters: 2878 + dm - DMNetwork object 2879 . numkeys - number of keys 2880 . keys - array of keys that define the components of the variables you wish to extract 2881 . blocksize - block size of the variables associated to the component 2882 . nselectedvar - number of variables in each block to select 2883 - selectedvar - the offset into the block of each variable in each block to select 2884 2885 Output Parameters: 2886 . is - the index set 2887 2888 Level: Advanced 2889 2890 Notes: 2891 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. 2892 2893 .seealso: `DMNetworkCreate()`, `ISCreateGeneral()`, `DMNetworkCreateLocalIS()` 2894 @*/ 2895 PetscErrorCode DMNetworkCreateIS(DM dm,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],IS *is) 2896 { 2897 MPI_Comm comm; 2898 DM_Network *network = (DM_Network*)dm->data; 2899 PetscInt i,p,estart,eend,vstart,vend,nidx,*idx; 2900 PetscBool ghost; 2901 2902 PetscFunctionBegin; 2903 PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 2904 2905 /* Check input parameters */ 2906 for (i=0; i<numkeys; i++) { 2907 if (!blocksize || blocksize[i] == -1) continue; 2908 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]); 2909 } 2910 2911 PetscCall(DMNetworkGetEdgeRange(dm,&estart,&eend)); 2912 PetscCall(DMNetworkGetVertexRange(dm,&vstart,&vend)); 2913 2914 /* Get local number of idx */ 2915 nidx = 0; 2916 for (p=estart; p<eend; p++) { 2917 PetscCall(DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx)); 2918 } 2919 for (p=vstart; p<vend; p++) { 2920 PetscCall(DMNetworkIsGhostVertex(dm,p,&ghost)); 2921 if (ghost) continue; 2922 PetscCall(DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx)); 2923 } 2924 2925 /* Compute idx */ 2926 PetscCall(PetscMalloc1(nidx,&idx)); 2927 i = 0; 2928 for (p=estart; p<eend; p++) { 2929 PetscCall(DMISComputeIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx)); 2930 } 2931 for (p=vstart; p<vend; p++) { 2932 PetscCall(DMNetworkIsGhostVertex(dm,p,&ghost)); 2933 if (ghost) continue; 2934 PetscCall(DMISComputeIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx)); 2935 } 2936 2937 /* Create is */ 2938 PetscCall(ISCreateGeneral(comm,nidx,idx,PETSC_COPY_VALUES,is)); 2939 PetscCall(PetscFree(idx)); 2940 PetscFunctionReturn(0); 2941 } 2942 2943 static inline PetscErrorCode DMISComputeLocalIdx_private(DM dm,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],PetscInt *ii,PetscInt *idx) 2944 { 2945 PetscInt i,j,ncomps,nvar,key,offsetl,k,k1,offset=0; 2946 DM_Network *network = (DM_Network*)dm->data; 2947 DMNetworkComponentHeader header; 2948 2949 PetscFunctionBegin; 2950 PetscCall(PetscSectionGetOffset(network->DataSection,p,&offset)); 2951 ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 2952 header = (DMNetworkComponentHeader)(network->componentdataarray+offset); 2953 2954 for (i=0; i<ncomps; i++) { 2955 key = header->key[i]; 2956 nvar = header->nvar[i]; 2957 for (j=0; j<numkeys; j++) { 2958 if (key != keys[j]) continue; 2959 2960 PetscCall(DMNetworkGetLocalVecOffset(dm,p,i,&offsetl)); 2961 if (!blocksize || blocksize[j] == -1) { 2962 for (k=0; k<nvar; k++) idx[(*ii)++] = offsetl + k; 2963 } else { 2964 for (k=0; k<nvar; k+=blocksize[j]) { 2965 for (k1=0; k1<nselectedvar[j]; k1++) idx[(*ii)++] = offsetl + k + selectedvar[j][k1]; 2966 } 2967 } 2968 } 2969 } 2970 PetscFunctionReturn(0); 2971 } 2972 2973 /*@ 2974 DMNetworkCreateLocalIS - Create an index set object from the local vector of the network 2975 2976 Not collective 2977 2978 Input Parameters: 2979 + dm - DMNetwork object 2980 . numkeys - number of keys 2981 . keys - array of keys that define the components of the variables you wish to extract 2982 . blocksize - block size of the variables associated to the component 2983 . nselectedvar - number of variables in each block to select 2984 - selectedvar - the offset into the block of each variable in each block to select 2985 2986 Output Parameters: 2987 . is - the index set 2988 2989 Level: Advanced 2990 2991 Notes: 2992 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. 2993 2994 .seealso: `DMNetworkCreate()`, `DMNetworkCreateIS`, `ISCreateGeneral()` 2995 @*/ 2996 PetscErrorCode DMNetworkCreateLocalIS(DM dm,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],IS *is) 2997 { 2998 DM_Network *network = (DM_Network*)dm->data; 2999 PetscInt i,p,pstart,pend,nidx,*idx; 3000 3001 PetscFunctionBegin; 3002 /* Check input parameters */ 3003 for (i=0; i<numkeys; i++) { 3004 if (!blocksize || blocksize[i] == -1) continue; 3005 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]); 3006 } 3007 3008 pstart = network->pStart; 3009 pend = network->pEnd; 3010 3011 /* Get local number of idx */ 3012 nidx = 0; 3013 for (p=pstart; p<pend; p++) { 3014 PetscCall(DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx)); 3015 } 3016 3017 /* Compute local idx */ 3018 PetscCall(PetscMalloc1(nidx,&idx)); 3019 i = 0; 3020 for (p=pstart; p<pend; p++) { 3021 PetscCall(DMISComputeLocalIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx)); 3022 } 3023 3024 /* Create is */ 3025 PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,is)); 3026 PetscCall(PetscFree(idx)); 3027 PetscFunctionReturn(0); 3028 } 3029