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