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