1 #include <petsc/private/dmnetworkimpl.h> /*I "petscdmnetwork.h" I*/ 2 #include <petscdmplex.h> 3 #include <petscsf.h> 4 5 /*@ 6 DMNetworkGetPlex - Gets the Plex DM associated with this network DM 7 8 Not collective 9 10 Input Parameters: 11 + netdm - the dm object 12 - plexmdm - the plex dm object 13 14 Level: Advanced 15 16 .seealso: DMNetworkCreate() 17 @*/ 18 PetscErrorCode DMNetworkGetPlex(DM netdm, DM *plexdm) 19 { 20 DM_Network *network = (DM_Network*) netdm->data; 21 22 PetscFunctionBegin; 23 *plexdm = network->plex; 24 PetscFunctionReturn(0); 25 } 26 27 /*@ 28 DMNetworkSetSizes - Sets the number of subnetworks, local and global vertices and edges for each subnetwork. 29 30 Collective on DM 31 32 Input Parameters: 33 + dm - the dm object 34 . Nsubnet - number of subnetworks 35 . NsubnetCouple - number of coupling subnetworks 36 . nV - number of local vertices for each subnetwork and coupling subnetwork 37 . nE - number of local edges for each subnetwork and coupling subnetwork 38 . NV - number of global vertices (or PETSC_DETERMINE) for each subnetwork and coupling subnetwork 39 - NE - number of global edges (or PETSC_DETERMINE) for each subnetwork and coupling subnetwork 40 41 Notes 42 If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang. 43 44 You cannot change the sizes once they have been set. nV, nE, NV and NE are arrays of length Nsubnet+NsubnetCouple. 45 46 Level: intermediate 47 48 .seealso: DMNetworkCreate() 49 @*/ 50 PetscErrorCode DMNetworkSetSizes(DM dm,PetscInt Nsubnet,PetscInt NsubnetCouple,PetscInt nV[], PetscInt nE[], PetscInt NV[], PetscInt NE[]) 51 { 52 PetscErrorCode ierr; 53 DM_Network *network = (DM_Network*) dm->data; 54 PetscInt a[2],b[2],i; 55 56 PetscFunctionBegin; 57 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 58 if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet); 59 if (NsubnetCouple < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of coupling subnetworks %D cannot be less than 0",NsubnetCouple); 60 61 if (Nsubnet > 0) PetscValidLogicalCollectiveInt(dm,Nsubnet,2); 62 if (network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network"); 63 64 network->nsubnet = Nsubnet + NsubnetCouple; 65 Nsubnet += NsubnetCouple; 66 network->ncsubnet = NsubnetCouple; 67 ierr = PetscCalloc1(Nsubnet,&network->subnet);CHKERRQ(ierr); 68 for(i=0; i < Nsubnet; i++) { 69 if (NV) { 70 if (NV[i] > 0) PetscValidLogicalCollectiveInt(dm,NV[i],6); 71 if (NV[i] > 0 && nV[i] > NV[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local vertex size %D cannot be larger than global vertex size %D",i,nV[i],NV[i]); 72 } 73 if (NE) { 74 if (NE[i] > 0) PetscValidLogicalCollectiveInt(dm,NE[i],7); 75 if (NE[i] > 0 && nE[i] > NE[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local edge size %D cannot be larger than global edge size %D",i,nE[i],NE[i]); 76 } 77 78 a[0] = nV[i]; a[1] = nE[i]; 79 ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 80 network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1]; 81 82 network->subnet[i].id = i; 83 84 network->subnet[i].nvtx = nV[i]; 85 network->subnet[i].vStart = network->nVertices; 86 network->subnet[i].vEnd = network->subnet[i].vStart + network->subnet[i].Nvtx; 87 network->nVertices += network->subnet[i].nvtx; 88 network->NVertices += network->subnet[i].Nvtx; 89 90 network->subnet[i].nedge = nE[i]; 91 network->subnet[i].eStart = network->nEdges; 92 network->subnet[i].eEnd = network->subnet[i].eStart + network->subnet[i].Nedge; 93 network->nEdges += network->subnet[i].nedge; 94 network->NEdges += network->subnet[i].Nedge; 95 } 96 PetscFunctionReturn(0); 97 } 98 99 /*@ 100 DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network 101 102 Logically collective on DM 103 104 Input Parameters: 105 + dm - the dm object 106 . edgelist - list of edges for each subnetwork 107 - edgelistCouple - list of edges for each coupling subnetwork 108 109 Notes: 110 There is no copy involved in this operation, only the pointer is referenced. The edgelist should 111 not be destroyed before the call to DMNetworkLayoutSetUp 112 113 Level: intermediate 114 115 Example usage: 116 Consider the following 2 separate networks and a coupling network: 117 118 .vb 119 network 0: v0 -> v1 -> v2 -> v3 120 network 1: v1 -> v2 -> v0 121 coupling network: network 1: v2 -> network 0: v0 122 .ve 123 124 The resulting input 125 edgelist[0] = [0 1 | 1 2 | 2 3]; 126 edgelist[1] = [1 2 | 2 0] 127 edgelistCouple[0] = [(network)1 (v)2 (network)0 (v)0]. 128 129 .seealso: DMNetworkCreate, DMNetworkSetSizes 130 @*/ 131 PetscErrorCode DMNetworkSetEdgeList(DM dm,int *edgelist[],int *edgelistCouple[]) 132 { 133 DM_Network *network = (DM_Network*) dm->data; 134 PetscInt i,nsubnet,ncsubnet=network->ncsubnet; 135 136 PetscFunctionBegin; 137 nsubnet = network->nsubnet - ncsubnet; 138 for(i=0; i < nsubnet; i++) { 139 network->subnet[i].edgelist = edgelist[i]; 140 } 141 if (edgelistCouple) { 142 PetscInt j; 143 j = 0; 144 nsubnet = network->nsubnet; 145 while (i < nsubnet) network->subnet[i++].edgelist = edgelistCouple[j++]; 146 } 147 PetscFunctionReturn(0); 148 } 149 150 /*@ 151 DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network 152 153 Collective on DM 154 155 Input Parameters 156 . DM - the dmnetwork object 157 158 Notes: 159 This routine should be called after the network sizes and edgelists have been provided. It creates 160 the bare layout of the network and sets up the network to begin insertion of components. 161 162 All the components should be registered before calling this routine. 163 164 Level: intermediate 165 166 .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList 167 @*/ 168 PetscErrorCode DMNetworkLayoutSetUp(DM dm) 169 { 170 PetscErrorCode ierr; 171 DM_Network *network = (DM_Network*)dm->data; 172 PetscInt dim = 1; /* One dimensional network */ 173 PetscInt numCorners=2,spacedim=2; 174 double *vertexcoords=NULL; 175 PetscInt i,j,ndata,ctr=0,nsubnet; 176 PetscInt k,netid,vid; 177 int *edgelist_couple=NULL; 178 179 PetscFunctionBegin; 180 if (network->nVertices) { 181 ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr); 182 } 183 184 /* Create the edgelist for the network by concatenating edgelists of the subnetworks */ 185 nsubnet = network->nsubnet - network->ncsubnet; 186 ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr); 187 for (i=0; i < nsubnet; i++) { 188 for (j = 0; j < network->subnet[i].nedge; j++) { 189 network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j]; 190 network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1]; 191 ctr++; 192 } 193 } 194 195 i = nsubnet; /* netid of coupling subnet */ 196 nsubnet = network->nsubnet; 197 while (i < nsubnet) { 198 edgelist_couple = network->subnet[i].edgelist; 199 k = 0; 200 for (j = 0; j < network->subnet[i].nedge; j++) { 201 netid = edgelist_couple[k]; vid = edgelist_couple[k+1]; 202 network->edges[2*ctr] = network->subnet[netid].vStart + vid; k += 2; 203 204 netid = edgelist_couple[k]; vid = edgelist_couple[k+1]; 205 network->edges[2*ctr+1] = network->subnet[netid].vStart + vid; k+=2; 206 ctr++; 207 } 208 i++; 209 } 210 211 #if 0 212 for(i=0; i < network->nEdges; i++) { 213 ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr); 214 printf("\n"); 215 } 216 #endif 217 218 ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr); 219 if (network->nVertices) { 220 ierr = PetscFree(vertexcoords);CHKERRQ(ierr); 221 } 222 ierr = PetscFree(network->edges);CHKERRQ(ierr); 223 224 ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr); 225 ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr); 226 ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr); 227 228 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr); 229 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr); 230 ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr); 231 ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr); 232 233 /* Create vertices and edges array for the subnetworks */ 234 for (j=0; j < network->nsubnet; j++) { 235 ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr); 236 ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr); 237 /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. 238 These get updated when the vertices and edges are added. */ 239 network->subnet[j].nvtx = network->subnet[j].nedge = 0; 240 } 241 242 network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 243 ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr); 244 for (i=network->eStart; i < network->eEnd; i++) { 245 network->header[i].index = i; /* Global edge number */ 246 for (j=0; j < network->nsubnet; j++) { 247 if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) { 248 network->header[i].subnetid = j; /* Subnetwork id */ 249 network->subnet[j].edges[network->subnet[j].nedge++] = i; 250 break; 251 } 252 } 253 network->header[i].ndata = 0; 254 ndata = network->header[i].ndata; 255 ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 256 network->header[i].offset[ndata] = 0; 257 } 258 259 for(i=network->vStart; i < network->vEnd; i++) { 260 network->header[i].index = i - network->vStart; 261 for (j=0; j < network->nsubnet; j++) { 262 if ((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) { 263 network->header[i].subnetid = j; 264 network->subnet[j].vertices[network->subnet[j].nvtx++] = i; 265 break; 266 } 267 } 268 network->header[i].ndata = 0; 269 ndata = network->header[i].ndata; 270 ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 271 network->header[i].offset[ndata] = 0; 272 } 273 274 ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr); 275 PetscFunctionReturn(0); 276 } 277 278 /*@C 279 DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork 280 281 Input Parameters 282 + dm - the number object 283 - id - the ID (integer) of the subnetwork 284 285 Output Parameters 286 + nv - number of vertices (local) 287 . ne - number of edges (local) 288 . vtx - local vertices for this subnetwork 289 . edge - local edges for this subnetwork 290 291 Notes: 292 Cannot call this routine before DMNetworkLayoutSetup() 293 294 .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 295 @*/ 296 PetscErrorCode DMNetworkGetSubnetworkInfo(DM netdm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge) 297 { 298 DM_Network *network = (DM_Network*) netdm->data; 299 300 PetscFunctionBegin; 301 *nv = network->subnet[id].nvtx; 302 *ne = network->subnet[id].nedge; 303 *vtx = network->subnet[id].vertices; 304 *edge = network->subnet[id].edges; 305 PetscFunctionReturn(0); 306 } 307 308 /*@C 309 DMNetworkRegisterComponent - Registers the network component 310 311 Logically collective on DM 312 313 Input Parameters 314 + dm - the network object 315 . name - the component name 316 - size - the storage size in bytes for this component data 317 318 Output Parameters 319 . key - an integer key that defines the component 320 321 Notes 322 This routine should be called by all processors before calling DMNetworkLayoutSetup(). 323 324 Level: intermediate 325 326 .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 327 @*/ 328 PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key) 329 { 330 PetscErrorCode ierr; 331 DM_Network *network = (DM_Network*) dm->data; 332 DMNetworkComponent *component=&network->component[network->ncomponent]; 333 PetscBool flg=PETSC_FALSE; 334 PetscInt i; 335 336 PetscFunctionBegin; 337 for (i=0; i < network->ncomponent; i++) { 338 ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr); 339 if (flg) { 340 *key = i; 341 PetscFunctionReturn(0); 342 } 343 } 344 if(network->ncomponent == MAX_COMPONENTS) { 345 SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components registered exceeds the max %D",MAX_COMPONENTS); 346 } 347 348 ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr); 349 component->size = size/sizeof(DMNetworkComponentGenericDataType); 350 *key = network->ncomponent; 351 network->ncomponent++; 352 PetscFunctionReturn(0); 353 } 354 355 /*@ 356 DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices. 357 358 Not Collective 359 360 Input Parameters: 361 + dm - The DMNetwork object 362 363 Output Paramters: 364 + vStart - The first vertex point 365 - vEnd - One beyond the last vertex point 366 367 Level: intermediate 368 369 .seealso: DMNetworkGetEdgeRange 370 @*/ 371 PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd) 372 { 373 DM_Network *network = (DM_Network*)dm->data; 374 375 PetscFunctionBegin; 376 if (vStart) *vStart = network->vStart; 377 if (vEnd) *vEnd = network->vEnd; 378 PetscFunctionReturn(0); 379 } 380 381 /*@ 382 DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges. 383 384 Not Collective 385 386 Input Parameters: 387 + dm - The DMNetwork object 388 389 Output Paramters: 390 + eStart - The first edge point 391 - eEnd - One beyond the last edge point 392 393 Level: intermediate 394 395 .seealso: DMNetworkGetVertexRange 396 @*/ 397 PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd) 398 { 399 DM_Network *network = (DM_Network*)dm->data; 400 401 PetscFunctionBegin; 402 if (eStart) *eStart = network->eStart; 403 if (eEnd) *eEnd = network->eEnd; 404 PetscFunctionReturn(0); 405 } 406 407 /*@ 408 DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge. 409 410 Not Collective 411 412 Input Parameters: 413 + dm - DMNetwork object 414 - p - edge point 415 416 Output Paramters: 417 . index - user global numbering for the edge 418 419 Level: intermediate 420 421 .seealso: DMNetworkGetGlobalVertexIndex 422 @*/ 423 PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index) 424 { 425 PetscErrorCode ierr; 426 DM_Network *network = (DM_Network*)dm->data; 427 PetscInt offsetp; 428 DMNetworkComponentHeader header; 429 430 PetscFunctionBegin; 431 ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 432 header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 433 *index = header->index; 434 PetscFunctionReturn(0); 435 } 436 437 /*@ 438 DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex. 439 440 Not Collective 441 442 Input Parameters: 443 + dm - DMNetwork object 444 - p - vertex point 445 446 Output Paramters: 447 . index - user global numbering for the vertex 448 449 Level: intermediate 450 451 .seealso: DMNetworkGetGlobalEdgeIndex 452 @*/ 453 PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index) 454 { 455 PetscErrorCode ierr; 456 DM_Network *network = (DM_Network*)dm->data; 457 PetscInt offsetp; 458 DMNetworkComponentHeader header; 459 460 PetscFunctionBegin; 461 ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 462 header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 463 *index = header->index; 464 PetscFunctionReturn(0); 465 } 466 467 /* 468 DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the 469 component value from the component data array 470 471 Not Collective 472 473 Input Parameters: 474 + dm - The DMNetwork object 475 . p - vertex/edge point 476 - compnum - component number 477 478 Output Parameters: 479 + compkey - the key obtained when registering the component 480 - offset - offset into the component data array associated with the vertex/edge point 481 482 Notes: 483 Typical usage: 484 485 DMNetworkGetComponentDataArray(dm, &arr); 486 DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 487 Loop over vertices or edges 488 DMNetworkGetNumComponents(dm,v,&numcomps); 489 Loop over numcomps 490 DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset); 491 compdata = (UserCompDataType)(arr+offset); 492 493 Level: intermediate 494 495 .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray, 496 */ 497 PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset) 498 { 499 PetscErrorCode ierr; 500 PetscInt offsetp; 501 DMNetworkComponentHeader header; 502 DM_Network *network = (DM_Network*)dm->data; 503 504 PetscFunctionBegin; 505 ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 506 header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 507 if (compkey) *compkey = header->key[compnum]; 508 if (offset) *offset = offsetp+network->dataheadersize+header->offset[compnum]; 509 PetscFunctionReturn(0); 510 } 511 512 /*@ 513 DMNetworkGetComponent - Returns the network component and its key 514 515 Not Collective 516 517 Input Parameters 518 + dm - DMNetwork object 519 . p - edge or vertex point 520 - compnum - component number 521 522 Output Parameters: 523 + compkey - the key set for this computing during registration 524 - component - the component data 525 526 Notes: 527 Typical usage: 528 529 DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 530 Loop over vertices or edges 531 DMNetworkGetNumComponents(dm,v,&numcomps); 532 Loop over numcomps 533 DMNetworkGetComponent(dm,v,compnum,&key,&component); 534 535 Level: intermediate 536 537 .seealso: DMNetworkGetNumComponents, DMNetworkGetVariableOffset 538 @*/ 539 PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *key, void **component) 540 { 541 PetscErrorCode ierr; 542 DM_Network *network = (DM_Network*)dm->data; 543 PetscInt offsetd; 544 545 PetscFunctionBegin; 546 547 ierr = DMNetworkGetComponentKeyOffset(dm,p,compnum,key,&offsetd);CHKERRQ(ierr); 548 *component = network->componentdataarray+offsetd; 549 550 PetscFunctionReturn(0); 551 } 552 553 /*@ 554 DMNetworkAddComponent - Adds a network component at the given point (vertex/edge) 555 556 Not Collective 557 558 Input Parameters: 559 + dm - The DMNetwork object 560 . p - vertex/edge point 561 . componentkey - component key returned while registering the component 562 - compvalue - pointer to the data structure for the component 563 564 Level: intermediate 565 566 .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent 567 @*/ 568 PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue) 569 { 570 DM_Network *network = (DM_Network*)dm->data; 571 DMNetworkComponent *component = &network->component[componentkey]; 572 DMNetworkComponentHeader header = &network->header[p]; 573 DMNetworkComponentValue cvalue = &network->cvalue[p]; 574 PetscErrorCode ierr; 575 576 PetscFunctionBegin; 577 if (header->ndata == MAX_DATA_AT_POINT) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components at a point exceeds the max %D",MAX_DATA_AT_POINT); 578 579 header->size[header->ndata] = component->size; 580 ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr); 581 header->key[header->ndata] = componentkey; 582 if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1]; 583 584 cvalue->data[header->ndata] = (void*)compvalue; 585 header->ndata++; 586 PetscFunctionReturn(0); 587 } 588 589 /*@ 590 DMNetworkGetNumComponents - Get the number of components at a vertex/edge 591 592 Not Collective 593 594 Input Parameters: 595 + dm - The DMNetwork object 596 . p - vertex/edge point 597 598 Output Parameters: 599 . numcomponents - Number of components at the vertex/edge 600 601 Level: intermediate 602 603 .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent 604 @*/ 605 PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents) 606 { 607 PetscErrorCode ierr; 608 PetscInt offset; 609 DM_Network *network = (DM_Network*)dm->data; 610 611 PetscFunctionBegin; 612 ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 613 *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 614 PetscFunctionReturn(0); 615 } 616 617 /*@ 618 DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector. 619 620 Not Collective 621 622 Input Parameters: 623 + dm - The DMNetwork object 624 - p - the edge/vertex point 625 626 Output Parameters: 627 . offset - the offset 628 629 Level: intermediate 630 631 .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 632 @*/ 633 PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset) 634 { 635 PetscErrorCode ierr; 636 DM_Network *network = (DM_Network*)dm->data; 637 638 PetscFunctionBegin; 639 ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr); 640 PetscFunctionReturn(0); 641 } 642 643 /*@ 644 DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector. 645 646 Not Collective 647 648 Input Parameters: 649 + dm - The DMNetwork object 650 - p - the edge/vertex point 651 652 Output Parameters: 653 . offsetg - the offset 654 655 Level: intermediate 656 657 .seealso: DMNetworkGetVariableOffset, DMGetLocalVector 658 @*/ 659 PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg) 660 { 661 PetscErrorCode ierr; 662 DM_Network *network = (DM_Network*)dm->data; 663 664 PetscFunctionBegin; 665 ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr); 666 if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */ 667 PetscFunctionReturn(0); 668 } 669 670 /*@ 671 DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector. 672 673 Not Collective 674 675 Input Parameters: 676 + dm - The DMNetwork object 677 - p - the edge point 678 679 Output Parameters: 680 . offset - the offset 681 682 Level: intermediate 683 684 .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 685 @*/ 686 PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset) 687 { 688 PetscErrorCode ierr; 689 DM_Network *network = (DM_Network*)dm->data; 690 691 PetscFunctionBegin; 692 693 ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr); 694 PetscFunctionReturn(0); 695 } 696 697 /*@ 698 DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector. 699 700 Not Collective 701 702 Input Parameters: 703 + dm - The DMNetwork object 704 - p - the vertex point 705 706 Output Parameters: 707 . offset - the offset 708 709 Level: intermediate 710 711 .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 712 @*/ 713 PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset) 714 { 715 PetscErrorCode ierr; 716 DM_Network *network = (DM_Network*)dm->data; 717 718 PetscFunctionBegin; 719 720 p -= network->vStart; 721 722 ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr); 723 PetscFunctionReturn(0); 724 } 725 /*@ 726 DMNetworkAddNumVariables - Add number of variables associated with a given point. 727 728 Not Collective 729 730 Input Parameters: 731 + dm - The DMNetworkObject 732 . p - the vertex/edge point 733 - nvar - number of additional variables 734 735 Level: intermediate 736 737 .seealso: DMNetworkSetNumVariables 738 @*/ 739 PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar) 740 { 741 PetscErrorCode ierr; 742 DM_Network *network = (DM_Network*)dm->data; 743 744 PetscFunctionBegin; 745 ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr); 746 PetscFunctionReturn(0); 747 } 748 749 /*@ 750 DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point. 751 752 Not Collective 753 754 Input Parameters: 755 + dm - The DMNetworkObject 756 - p - the vertex/edge point 757 758 Output Parameters: 759 . nvar - number of variables 760 761 Level: intermediate 762 763 .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables 764 @*/ 765 PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar) 766 { 767 PetscErrorCode ierr; 768 DM_Network *network = (DM_Network*)dm->data; 769 770 PetscFunctionBegin; 771 ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 772 PetscFunctionReturn(0); 773 } 774 775 /*@ 776 DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point. 777 778 Not Collective 779 780 Input Parameters: 781 + dm - The DMNetworkObject 782 . p - the vertex/edge point 783 - nvar - number of variables 784 785 Level: intermediate 786 787 .seealso: DMNetworkAddNumVariables 788 @*/ 789 PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar) 790 { 791 PetscErrorCode ierr; 792 DM_Network *network = (DM_Network*)dm->data; 793 794 PetscFunctionBegin; 795 ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 796 PetscFunctionReturn(0); 797 } 798 799 /* Sets up the array that holds the data for all components and its associated section. This 800 function is called during DMSetUp() */ 801 PetscErrorCode DMNetworkComponentSetUp(DM dm) 802 { 803 PetscErrorCode ierr; 804 DM_Network *network = (DM_Network*)dm->data; 805 PetscInt arr_size,p,offset,offsetp,ncomp,i; 806 DMNetworkComponentHeader header; 807 DMNetworkComponentValue cvalue; 808 DMNetworkComponentGenericDataType *componentdataarray; 809 810 PetscFunctionBegin; 811 ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr); 812 ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr); 813 ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr); 814 componentdataarray = network->componentdataarray; 815 for (p = network->pStart; p < network->pEnd; p++) { 816 ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 817 /* Copy header */ 818 header = &network->header[p]; 819 ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 820 /* Copy data */ 821 cvalue = &network->cvalue[p]; 822 ncomp = header->ndata; 823 for (i = 0; i < ncomp; i++) { 824 offset = offsetp + network->dataheadersize + header->offset[i]; 825 ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 826 } 827 } 828 PetscFunctionReturn(0); 829 } 830 831 /* Sets up the section for dofs. This routine is called during DMSetUp() */ 832 PetscErrorCode DMNetworkVariablesSetUp(DM dm) 833 { 834 PetscErrorCode ierr; 835 DM_Network *network = (DM_Network*)dm->data; 836 837 PetscFunctionBegin; 838 ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr); 839 PetscFunctionReturn(0); 840 } 841 842 /*@C 843 DMNetworkGetComponentDataArray - Returns the component data array 844 845 Not Collective 846 847 Input Parameters: 848 . dm - The DMNetwork Object 849 850 Output Parameters: 851 . componentdataarray - array that holds data for all components 852 853 Level: intermediate 854 855 .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents 856 @*/ 857 PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray) 858 { 859 DM_Network *network = (DM_Network*)dm->data; 860 861 PetscFunctionBegin; 862 *componentdataarray = network->componentdataarray; 863 PetscFunctionReturn(0); 864 } 865 866 /* Get a subsection from a range of points */ 867 PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection) 868 { 869 PetscErrorCode ierr; 870 PetscInt i, nvar; 871 872 PetscFunctionBegin; 873 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr); 874 ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr); 875 for (i = pstart; i < pend; i++) { 876 ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr); 877 ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr); 878 } 879 880 ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr); 881 PetscFunctionReturn(0); 882 } 883 884 /* Create a submap of points with a GlobalToLocal structure */ 885 PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map) 886 { 887 PetscErrorCode ierr; 888 PetscInt i, *subpoints; 889 890 PetscFunctionBegin; 891 /* Create index sets to map from "points" to "subpoints" */ 892 ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr); 893 for (i = pstart; i < pend; i++) { 894 subpoints[i - pstart] = i; 895 } 896 ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr); 897 ierr = PetscFree(subpoints);CHKERRQ(ierr); 898 PetscFunctionReturn(0); 899 } 900 901 /*@ 902 DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute. 903 904 Collective 905 906 Input Parameters: 907 . dm - The DMNetworkObject 908 909 Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are: 910 911 points = [0 1 2 3 4 5 6] 912 913 where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]). 914 915 With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset. 916 917 Level: intermediate 918 919 @*/ 920 PetscErrorCode DMNetworkAssembleGraphStructures(DM dm) 921 { 922 PetscErrorCode ierr; 923 MPI_Comm comm; 924 PetscMPIInt rank, size; 925 DM_Network *network = (DM_Network*)dm->data; 926 927 PetscFunctionBegin; 928 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 929 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 930 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 931 932 /* Create maps for vertices and edges */ 933 ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr); 934 ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr); 935 936 /* Create local sub-sections */ 937 ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr); 938 ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr); 939 940 if (size > 1) { 941 ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr); 942 ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr); 943 ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr); 944 ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr); 945 } else { 946 /* create structures for vertex */ 947 ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr); 948 /* create structures for edge */ 949 ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr); 950 } 951 952 953 /* Add viewers */ 954 ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr); 955 ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr); 956 ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr); 957 ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr); 958 959 PetscFunctionReturn(0); 960 } 961 962 /*@ 963 DMNetworkDistribute - Distributes the network and moves associated component data. 964 965 Collective 966 967 Input Parameter: 968 + DM - the DMNetwork object 969 - overlap - The overlap of partitions, 0 is the default 970 971 Notes: 972 Distributes the network with <overlap>-overlapping partitioning of the edges. 973 974 Level: intermediate 975 976 .seealso: DMNetworkCreate 977 @*/ 978 PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap) 979 { 980 MPI_Comm comm; 981 PetscErrorCode ierr; 982 PetscMPIInt size; 983 DM_Network *oldDMnetwork = (DM_Network*)((*dm)->data); 984 DM_Network *newDMnetwork; 985 PetscSF pointsf; 986 DM newDM; 987 PetscPartitioner part; 988 PetscInt j,e,v,offset; 989 DMNetworkComponentHeader header; 990 991 PetscFunctionBegin; 992 993 ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr); 994 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 995 if (size == 1) PetscFunctionReturn(0); 996 997 ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr); 998 newDMnetwork = (DM_Network*)newDM->data; 999 newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 1000 1001 /* Enable runtime options for petscpartitioner */ 1002 ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr); 1003 ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 1004 1005 /* Distribute plex dm and dof section */ 1006 ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr); 1007 1008 /* Distribute dof section */ 1009 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr); 1010 ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr); 1011 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr); 1012 1013 /* Distribute data and associated section */ 1014 ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr); 1015 1016 ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr); 1017 ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr); 1018 ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr); 1019 newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart; 1020 newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart; 1021 newDMnetwork->NVertices = oldDMnetwork->NVertices; 1022 newDMnetwork->NEdges = oldDMnetwork->NEdges; 1023 1024 /* Set Dof section as the default section for dm */ 1025 ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr); 1026 ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr); 1027 1028 /* Set up subnetwork info in the newDM */ 1029 newDMnetwork->nsubnet = oldDMnetwork->nsubnet; 1030 ierr = PetscCalloc1(newDMnetwork->nsubnet,&newDMnetwork->subnet);CHKERRQ(ierr); 1031 /* Copy over the global number of vertices and edges in each subnetwork. Note that these are already 1032 calculated in DMNetworkLayoutSetUp() 1033 */ 1034 for(j=0; j < newDMnetwork->nsubnet; j++) { 1035 newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx; 1036 newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge; 1037 } 1038 1039 for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) { 1040 ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr); 1041 header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1042 newDMnetwork->subnet[header->subnetid].nedge++; 1043 } 1044 1045 for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) { 1046 ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr); 1047 header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1048 newDMnetwork->subnet[header->subnetid].nvtx++; 1049 } 1050 1051 /* Now create the vertices and edge arrays for the subnetworks */ 1052 for(j=0; j < newDMnetwork->nsubnet; j++) { 1053 ierr = PetscCalloc1(newDMnetwork->subnet[j].nedge,&newDMnetwork->subnet[j].edges);CHKERRQ(ierr); 1054 ierr = PetscCalloc1(newDMnetwork->subnet[j].nvtx,&newDMnetwork->subnet[j].vertices);CHKERRQ(ierr); 1055 /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. 1056 These get updated when the vertices and edges are added. */ 1057 newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0; 1058 } 1059 1060 /* Set the vertices and edges in each subnetwork */ 1061 for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) { 1062 ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr); 1063 header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1064 newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++] = e; 1065 } 1066 1067 for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) { 1068 ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr); 1069 header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1070 newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++] = v; 1071 } 1072 1073 /* Destroy point SF */ 1074 ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr); 1075 1076 ierr = DMDestroy(dm);CHKERRQ(ierr); 1077 *dm = newDM; 1078 PetscFunctionReturn(0); 1079 } 1080 1081 /*@C 1082 PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering. 1083 1084 Input Parameters: 1085 + masterSF - the original SF structure 1086 - map - a ISLocalToGlobal mapping that contains the subset of points 1087 1088 Output Parameters: 1089 . subSF - a subset of the masterSF for the desired subset. 1090 */ 1091 1092 PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) { 1093 1094 PetscErrorCode ierr; 1095 PetscInt nroots, nleaves, *ilocal_sub; 1096 PetscInt i, *ilocal_map, nroots_sub, nleaves_sub = 0; 1097 PetscInt *local_points, *remote_points; 1098 PetscSFNode *iremote_sub; 1099 const PetscInt *ilocal; 1100 const PetscSFNode *iremote; 1101 1102 PetscFunctionBegin; 1103 ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 1104 1105 /* Look for leaves that pertain to the subset of points. Get the local ordering */ 1106 ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr); 1107 ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr); 1108 for (i = 0; i < nleaves; i++) { 1109 if (ilocal_map[i] != -1) nleaves_sub += 1; 1110 } 1111 /* Re-number ilocal with subset numbering. Need information from roots */ 1112 ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr); 1113 for (i = 0; i < nroots; i++) local_points[i] = i; 1114 ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr); 1115 ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 1116 ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 1117 /* Fill up graph using local (that is, local to the subset) numbering. */ 1118 ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr); 1119 ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr); 1120 nleaves_sub = 0; 1121 for (i = 0; i < nleaves; i++) { 1122 if (ilocal_map[i] != -1) { 1123 ilocal_sub[nleaves_sub] = ilocal_map[i]; 1124 iremote_sub[nleaves_sub].rank = iremote[i].rank; 1125 iremote_sub[nleaves_sub].index = remote_points[ilocal[i]]; 1126 nleaves_sub += 1; 1127 } 1128 } 1129 ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr); 1130 ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr); 1131 1132 /* Create new subSF */ 1133 ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr); 1134 ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr); 1135 ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr); 1136 ierr = PetscFree(ilocal_map);CHKERRQ(ierr); 1137 ierr = PetscFree(iremote_sub);CHKERRQ(ierr); 1138 PetscFunctionReturn(0); 1139 } 1140 1141 /*@C 1142 DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point 1143 1144 Not Collective 1145 1146 Input Parameters: 1147 + dm - The DMNetwork object 1148 - p - the vertex point 1149 1150 Output Paramters: 1151 + nedges - number of edges connected to this vertex point 1152 - edges - List of edge points 1153 1154 Level: intermediate 1155 1156 Fortran Notes: 1157 Since it returns an array, this routine is only available in Fortran 90, and you must 1158 include petsc.h90 in your code. 1159 1160 .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices 1161 @*/ 1162 PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[]) 1163 { 1164 PetscErrorCode ierr; 1165 DM_Network *network = (DM_Network*)dm->data; 1166 1167 PetscFunctionBegin; 1168 ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr); 1169 ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr); 1170 PetscFunctionReturn(0); 1171 } 1172 1173 /*@C 1174 DMNetworkGetConnectedVertices - Return the connected vertices for this edge point 1175 1176 Not Collective 1177 1178 Input Parameters: 1179 + dm - The DMNetwork object 1180 - p - the edge point 1181 1182 Output Paramters: 1183 . vertices - vertices connected to this edge 1184 1185 Level: intermediate 1186 1187 Fortran Notes: 1188 Since it returns an array, this routine is only available in Fortran 90, and you must 1189 include petsc.h90 in your code. 1190 1191 .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges 1192 @*/ 1193 PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[]) 1194 { 1195 PetscErrorCode ierr; 1196 DM_Network *network = (DM_Network*)dm->data; 1197 1198 PetscFunctionBegin; 1199 ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr); 1200 PetscFunctionReturn(0); 1201 } 1202 1203 /*@ 1204 DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex 1205 1206 Not Collective 1207 1208 Input Parameters: 1209 + dm - The DMNetwork object 1210 . p - the vertex point 1211 1212 Output Parameter: 1213 . isghost - TRUE if the vertex is a ghost point 1214 1215 Level: intermediate 1216 1217 .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange 1218 @*/ 1219 PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost) 1220 { 1221 PetscErrorCode ierr; 1222 DM_Network *network = (DM_Network*)dm->data; 1223 PetscInt offsetg; 1224 PetscSection sectiong; 1225 1226 PetscFunctionBegin; 1227 *isghost = PETSC_FALSE; 1228 ierr = DMGetDefaultGlobalSection(network->plex,§iong);CHKERRQ(ierr); 1229 ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr); 1230 if (offsetg < 0) *isghost = PETSC_TRUE; 1231 PetscFunctionReturn(0); 1232 } 1233 1234 PetscErrorCode DMSetUp_Network(DM dm) 1235 { 1236 PetscErrorCode ierr; 1237 DM_Network *network=(DM_Network*)dm->data; 1238 1239 PetscFunctionBegin; 1240 ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr); 1241 ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr); 1242 1243 ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr); 1244 ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr); 1245 PetscFunctionReturn(0); 1246 } 1247 1248 /*@ 1249 DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices 1250 -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)? 1251 1252 Collective 1253 1254 Input Parameters: 1255 + dm - The DMNetwork object 1256 . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges 1257 - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices 1258 1259 Level: intermediate 1260 1261 @*/ 1262 PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg) 1263 { 1264 DM_Network *network=(DM_Network*)dm->data; 1265 PetscErrorCode ierr; 1266 1267 PetscFunctionBegin; 1268 network->userEdgeJacobian = eflg; 1269 network->userVertexJacobian = vflg; 1270 1271 if (eflg && !network->Je) { 1272 ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr); 1273 } 1274 1275 if (vflg && !network->Jv) { 1276 PetscInt i,*vptr,nedges,vStart=network->vStart; 1277 PetscInt nVertices = network->nVertices,nedges_total; 1278 const PetscInt *edges; 1279 1280 /* count nvertex_total */ 1281 nedges_total = 0; 1282 ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr); 1283 1284 vptr[0] = 0; 1285 for (i=0; i<nVertices; i++) { 1286 ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr); 1287 nedges_total += nedges; 1288 vptr[i+1] = vptr[i] + 2*nedges + 1; 1289 } 1290 1291 ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr); 1292 network->Jvptr = vptr; 1293 } 1294 PetscFunctionReturn(0); 1295 } 1296 1297 /*@ 1298 DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network 1299 1300 Not Collective 1301 1302 Input Parameters: 1303 + dm - The DMNetwork object 1304 . p - the edge point 1305 - J - array (size = 3) of Jacobian submatrices for this edge point: 1306 J[0]: this edge 1307 J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices() 1308 1309 Level: intermediate 1310 1311 .seealso: DMNetworkVertexSetMatrix 1312 @*/ 1313 PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[]) 1314 { 1315 DM_Network *network=(DM_Network*)dm->data; 1316 1317 PetscFunctionBegin; 1318 if (!network->Je) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix"); 1319 1320 if (J) { 1321 network->Je[3*p] = J[0]; 1322 network->Je[3*p+1] = J[1]; 1323 network->Je[3*p+2] = J[2]; 1324 } 1325 PetscFunctionReturn(0); 1326 } 1327 1328 /*@ 1329 DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network 1330 1331 Not Collective 1332 1333 Input Parameters: 1334 + dm - The DMNetwork object 1335 . p - the vertex point 1336 - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point: 1337 J[0]: this vertex 1338 J[1+2*i]: i-th supporting edge 1339 J[1+2*i+1]: i-th connected vertex 1340 1341 Level: intermediate 1342 1343 .seealso: DMNetworkEdgeSetMatrix 1344 @*/ 1345 PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[]) 1346 { 1347 PetscErrorCode ierr; 1348 DM_Network *network=(DM_Network*)dm->data; 1349 PetscInt i,*vptr,nedges,vStart=network->vStart; 1350 const PetscInt *edges; 1351 1352 PetscFunctionBegin; 1353 if (!network->Jv) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix"); 1354 1355 if (J) { 1356 vptr = network->Jvptr; 1357 network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */ 1358 1359 /* Set Jacobian for each supporting edge and connected vertex */ 1360 ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr); 1361 for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i]; 1362 } 1363 PetscFunctionReturn(0); 1364 } 1365 1366 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 1367 { 1368 PetscErrorCode ierr; 1369 PetscInt j; 1370 PetscScalar val=(PetscScalar)ncols; 1371 1372 PetscFunctionBegin; 1373 if (!ghost) { 1374 for (j=0; j<nrows; j++) { 1375 ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 1376 } 1377 } else { 1378 for (j=0; j<nrows; j++) { 1379 ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 1380 } 1381 } 1382 PetscFunctionReturn(0); 1383 } 1384 1385 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 1386 { 1387 PetscErrorCode ierr; 1388 PetscInt j,ncols_u; 1389 PetscScalar val; 1390 1391 PetscFunctionBegin; 1392 if (!ghost) { 1393 for (j=0; j<nrows; j++) { 1394 ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 1395 val = (PetscScalar)ncols_u; 1396 ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 1397 ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 1398 } 1399 } else { 1400 for (j=0; j<nrows; j++) { 1401 ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 1402 val = (PetscScalar)ncols_u; 1403 ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 1404 ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 1405 } 1406 } 1407 PetscFunctionReturn(0); 1408 } 1409 1410 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 1411 { 1412 PetscErrorCode ierr; 1413 1414 PetscFunctionBegin; 1415 if (Ju) { 1416 ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 1417 } else { 1418 ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 1419 } 1420 PetscFunctionReturn(0); 1421 } 1422 1423 PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1424 { 1425 PetscErrorCode ierr; 1426 PetscInt j,*cols; 1427 PetscScalar *zeros; 1428 1429 PetscFunctionBegin; 1430 ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr); 1431 for (j=0; j<ncols; j++) cols[j] = j+ cstart; 1432 ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr); 1433 ierr = PetscFree2(cols,zeros);CHKERRQ(ierr); 1434 PetscFunctionReturn(0); 1435 } 1436 1437 PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1438 { 1439 PetscErrorCode ierr; 1440 PetscInt j,M,N,row,col,ncols_u; 1441 const PetscInt *cols; 1442 PetscScalar zero=0.0; 1443 1444 PetscFunctionBegin; 1445 ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr); 1446 if (nrows != M || ncols != N) SETERRQ4(PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%D by %D must equal %D by %D",nrows,ncols,M,N); 1447 1448 for (row=0; row<nrows; row++) { 1449 ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 1450 for (j=0; j<ncols_u; j++) { 1451 col = cols[j] + cstart; 1452 ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr); 1453 } 1454 ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 1455 } 1456 PetscFunctionReturn(0); 1457 } 1458 1459 PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1460 { 1461 PetscErrorCode ierr; 1462 1463 PetscFunctionBegin; 1464 if (Ju) { 1465 ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1466 } else { 1467 ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1468 } 1469 PetscFunctionReturn(0); 1470 } 1471 1472 /* 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. 1473 */ 1474 PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog) 1475 { 1476 PetscErrorCode ierr; 1477 PetscInt i, size, dof; 1478 PetscInt *glob2loc; 1479 1480 PetscFunctionBegin; 1481 ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr); 1482 ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr); 1483 1484 for (i = 0; i < size; i++) { 1485 ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr); 1486 dof = (dof >= 0) ? dof : -(dof + 1); 1487 glob2loc[i] = dof; 1488 } 1489 1490 ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 1491 #if 0 1492 ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1493 #endif 1494 PetscFunctionReturn(0); 1495 } 1496 1497 #include <petsc/private/matimpl.h> 1498 PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J) 1499 { 1500 PetscErrorCode ierr; 1501 PetscMPIInt rank, size; 1502 DM_Network *network = (DM_Network*)dm->data; 1503 PetscInt eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize; 1504 PetscInt cstart,ncols,j,e,v; 1505 PetscBool ghost,ghost_vc,ghost2,isNest; 1506 Mat Juser; 1507 PetscSection sectionGlobal; 1508 PetscInt nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */ 1509 const PetscInt *edges,*cone; 1510 MPI_Comm comm; 1511 MatType mtype; 1512 Vec vd_nz,vo_nz; 1513 PetscInt *dnnz,*onnz; 1514 PetscScalar *vdnz,*vonz; 1515 1516 PetscFunctionBegin; 1517 mtype = dm->mattype; 1518 ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr); 1519 1520 if (isNest) { 1521 /* ierr = DMCreateMatrix_Network_Nest(); */ 1522 PetscInt eDof, vDof; 1523 Mat j11, j12, j21, j22, bA[2][2]; 1524 ISLocalToGlobalMapping eISMap, vISMap; 1525 1526 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1527 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1528 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1529 1530 ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr); 1531 ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr); 1532 1533 ierr = MatCreate(comm, &j11);CHKERRQ(ierr); 1534 ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 1535 ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr); 1536 1537 ierr = MatCreate(comm, &j12);CHKERRQ(ierr); 1538 ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr); 1539 ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr); 1540 1541 ierr = MatCreate(comm, &j21);CHKERRQ(ierr); 1542 ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 1543 ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr); 1544 1545 ierr = MatCreate(comm, &j22);CHKERRQ(ierr); 1546 ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 1547 ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr); 1548 1549 bA[0][0] = j11; 1550 bA[0][1] = j12; 1551 bA[1][0] = j21; 1552 bA[1][1] = j22; 1553 1554 ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr); 1555 ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr); 1556 1557 ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr); 1558 ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr); 1559 ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr); 1560 ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr); 1561 1562 ierr = MatSetUp(j11);CHKERRQ(ierr); 1563 ierr = MatSetUp(j12);CHKERRQ(ierr); 1564 ierr = MatSetUp(j21);CHKERRQ(ierr); 1565 ierr = MatSetUp(j22);CHKERRQ(ierr); 1566 1567 ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr); 1568 ierr = MatSetUp(*J);CHKERRQ(ierr); 1569 ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr); 1570 ierr = MatDestroy(&j11);CHKERRQ(ierr); 1571 ierr = MatDestroy(&j12);CHKERRQ(ierr); 1572 ierr = MatDestroy(&j21);CHKERRQ(ierr); 1573 ierr = MatDestroy(&j22);CHKERRQ(ierr); 1574 1575 ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1576 ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1577 1578 /* Free structures */ 1579 ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr); 1580 ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr); 1581 1582 PetscFunctionReturn(0); 1583 } else if (!network->userEdgeJacobian && !network->userVertexJacobian) { 1584 /* user does not provide Jacobian blocks */ 1585 ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr); 1586 ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 1587 PetscFunctionReturn(0); 1588 } 1589 1590 ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr); 1591 ierr = DMGetDefaultGlobalSection(network->plex,§ionGlobal);CHKERRQ(ierr); 1592 ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr); 1593 ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1594 1595 ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr); 1596 ierr = MatSetFromOptions(*J);CHKERRQ(ierr); 1597 1598 /* (1) Set matrix preallocation */ 1599 /*------------------------------*/ 1600 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1601 ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr); 1602 ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr); 1603 ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr); 1604 ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr); 1605 ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr); 1606 1607 /* Set preallocation for edges */ 1608 /*-----------------------------*/ 1609 ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr); 1610 1611 ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr); 1612 for (e=eStart; e<eEnd; e++) { 1613 /* Get row indices */ 1614 ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 1615 ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 1616 if (nrows) { 1617 for (j=0; j<nrows; j++) rows[j] = j + rstart; 1618 1619 /* Set preallocation for conntected vertices */ 1620 ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr); 1621 for (v=0; v<2; v++) { 1622 ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 1623 1624 if (network->Je) { 1625 Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 1626 } else Juser = NULL; 1627 ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr); 1628 ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1629 } 1630 1631 /* Set preallocation for edge self */ 1632 cstart = rstart; 1633 if (network->Je) { 1634 Juser = network->Je[3*e]; /* Jacobian(e,e) */ 1635 } else Juser = NULL; 1636 ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1637 } 1638 } 1639 1640 /* Set preallocation for vertices */ 1641 /*--------------------------------*/ 1642 ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 1643 if (vEnd - vStart) vptr = network->Jvptr; 1644 1645 for (v=vStart; v<vEnd; v++) { 1646 /* Get row indices */ 1647 ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1648 ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 1649 if (!nrows) continue; 1650 1651 ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1652 if (ghost) { 1653 ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1654 } else { 1655 rows_v = rows; 1656 } 1657 1658 for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1659 1660 /* Get supporting edges and connected vertices */ 1661 ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1662 1663 for (e=0; e<nedges; e++) { 1664 /* Supporting edges */ 1665 ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1666 ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1667 1668 if (network->Jv) { 1669 Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 1670 } else Juser = NULL; 1671 ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1672 1673 /* Connected vertices */ 1674 ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 1675 vc = (v == cone[0]) ? cone[1]:cone[0]; 1676 ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr); 1677 1678 ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1679 1680 if (network->Jv) { 1681 Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 1682 } else Juser = NULL; 1683 if (ghost_vc||ghost) { 1684 ghost2 = PETSC_TRUE; 1685 } else { 1686 ghost2 = PETSC_FALSE; 1687 } 1688 ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr); 1689 } 1690 1691 /* Set preallocation for vertex self */ 1692 ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1693 if (!ghost) { 1694 ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 1695 if (network->Jv) { 1696 Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 1697 } else Juser = NULL; 1698 ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1699 } 1700 if (ghost) { 1701 ierr = PetscFree(rows_v);CHKERRQ(ierr); 1702 } 1703 } 1704 1705 ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr); 1706 ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr); 1707 1708 ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr); 1709 1710 ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr); 1711 ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr); 1712 1713 ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr); 1714 ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr); 1715 for (j=0; j<localSize; j++) { 1716 dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]); 1717 onnz[j] = (PetscInt)PetscRealPart(vonz[j]); 1718 } 1719 ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr); 1720 ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr); 1721 ierr = VecDestroy(&vd_nz);CHKERRQ(ierr); 1722 ierr = VecDestroy(&vo_nz);CHKERRQ(ierr); 1723 1724 ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr); 1725 ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr); 1726 ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 1727 1728 ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr); 1729 1730 /* (2) Set matrix entries for edges */ 1731 /*----------------------------------*/ 1732 for (e=eStart; e<eEnd; e++) { 1733 /* Get row indices */ 1734 ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 1735 ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 1736 if (nrows) { 1737 for (j=0; j<nrows; j++) rows[j] = j + rstart; 1738 1739 /* Set matrix entries for conntected vertices */ 1740 ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr); 1741 for (v=0; v<2; v++) { 1742 ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr); 1743 ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 1744 1745 if (network->Je) { 1746 Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 1747 } else Juser = NULL; 1748 ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1749 } 1750 1751 /* Set matrix entries for edge self */ 1752 cstart = rstart; 1753 if (network->Je) { 1754 Juser = network->Je[3*e]; /* Jacobian(e,e) */ 1755 } else Juser = NULL; 1756 ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr); 1757 } 1758 } 1759 1760 /* Set matrix entries for vertices */ 1761 /*---------------------------------*/ 1762 for (v=vStart; v<vEnd; v++) { 1763 /* Get row indices */ 1764 ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1765 ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 1766 if (!nrows) continue; 1767 1768 ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1769 if (ghost) { 1770 ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1771 } else { 1772 rows_v = rows; 1773 } 1774 for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1775 1776 /* Get supporting edges and connected vertices */ 1777 ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1778 1779 for (e=0; e<nedges; e++) { 1780 /* Supporting edges */ 1781 ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1782 ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1783 1784 if (network->Jv) { 1785 Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 1786 } else Juser = NULL; 1787 ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1788 1789 /* Connected vertices */ 1790 ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 1791 vc = (v == cone[0]) ? cone[1]:cone[0]; 1792 1793 ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr); 1794 ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1795 1796 if (network->Jv) { 1797 Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 1798 } else Juser = NULL; 1799 ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1800 } 1801 1802 /* Set matrix entries for vertex self */ 1803 if (!ghost) { 1804 ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 1805 if (network->Jv) { 1806 Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 1807 } else Juser = NULL; 1808 ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr); 1809 } 1810 if (ghost) { 1811 ierr = PetscFree(rows_v);CHKERRQ(ierr); 1812 } 1813 } 1814 ierr = PetscFree(rows);CHKERRQ(ierr); 1815 1816 ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1817 ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1818 1819 ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 1820 PetscFunctionReturn(0); 1821 } 1822 1823 PetscErrorCode DMDestroy_Network(DM dm) 1824 { 1825 PetscErrorCode ierr; 1826 DM_Network *network = (DM_Network*)dm->data; 1827 PetscInt j; 1828 1829 PetscFunctionBegin; 1830 if (--network->refct > 0) PetscFunctionReturn(0); 1831 if (network->Je) { 1832 ierr = PetscFree(network->Je);CHKERRQ(ierr); 1833 } 1834 if (network->Jv) { 1835 ierr = PetscFree(network->Jvptr);CHKERRQ(ierr); 1836 ierr = PetscFree(network->Jv);CHKERRQ(ierr); 1837 } 1838 1839 ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr); 1840 ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr); 1841 ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr); 1842 if (network->vertex.sf) { 1843 ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr); 1844 } 1845 /* edge */ 1846 ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr); 1847 ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr); 1848 ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr); 1849 if (network->edge.sf) { 1850 ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr); 1851 } 1852 ierr = DMDestroy(&network->plex);CHKERRQ(ierr); 1853 network->edges = NULL; 1854 ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr); 1855 ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr); 1856 1857 for(j=0; j < network->nsubnet; j++) { 1858 ierr = PetscFree(network->subnet[j].edges);CHKERRQ(ierr); 1859 ierr = PetscFree(network->subnet[j].vertices);CHKERRQ(ierr); 1860 } 1861 ierr = PetscFree(network->subnet);CHKERRQ(ierr); 1862 ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr); 1863 ierr = PetscFree(network->cvalue);CHKERRQ(ierr); 1864 ierr = PetscFree(network->header);CHKERRQ(ierr); 1865 ierr = PetscFree(network);CHKERRQ(ierr); 1866 PetscFunctionReturn(0); 1867 } 1868 1869 PetscErrorCode DMView_Network(DM dm,PetscViewer viewer) 1870 { 1871 PetscErrorCode ierr; 1872 DM_Network *network = (DM_Network*)dm->data; 1873 1874 PetscFunctionBegin; 1875 ierr = DMView(network->plex,viewer);CHKERRQ(ierr); 1876 PetscFunctionReturn(0); 1877 } 1878 1879 PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l) 1880 { 1881 PetscErrorCode ierr; 1882 DM_Network *network = (DM_Network*)dm->data; 1883 1884 PetscFunctionBegin; 1885 ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr); 1886 PetscFunctionReturn(0); 1887 } 1888 1889 PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l) 1890 { 1891 PetscErrorCode ierr; 1892 DM_Network *network = (DM_Network*)dm->data; 1893 1894 PetscFunctionBegin; 1895 ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr); 1896 PetscFunctionReturn(0); 1897 } 1898 1899 PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g) 1900 { 1901 PetscErrorCode ierr; 1902 DM_Network *network = (DM_Network*)dm->data; 1903 1904 PetscFunctionBegin; 1905 ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr); 1906 PetscFunctionReturn(0); 1907 } 1908 1909 PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g) 1910 { 1911 PetscErrorCode ierr; 1912 DM_Network *network = (DM_Network*)dm->data; 1913 1914 PetscFunctionBegin; 1915 ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr); 1916 PetscFunctionReturn(0); 1917 } 1918