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