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