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