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