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