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