1 #include <petsc/private/dmnetworkimpl.h> /*I "petscdmnetwork.h" I*/ 2 #include <petscdmplex.h> 3 #include <petscsf.h> 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "DMNetworkSetSizes" 7 /*@ 8 DMNetworkSetSizes - Sets the local and global vertices and edges. 9 10 Collective on DM 11 12 Input Parameters: 13 + dm - the dm object 14 . nV - number of local vertices 15 . nE - number of local edges 16 . NV - number of global vertices (or PETSC_DETERMINE) 17 - NE - number of global edges (or PETSC_DETERMINE) 18 19 Notes 20 If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang. 21 22 You cannot change the sizes once they have been set 23 24 Level: intermediate 25 26 .seealso: DMNetworkCreate 27 @*/ 28 PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt nV, PetscInt nE, PetscInt NV, PetscInt NE) 29 { 30 PetscErrorCode ierr; 31 DM_Network *network = (DM_Network*) dm->data; 32 PetscInt a[2],b[2]; 33 34 PetscFunctionBegin; 35 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 36 if (NV > 0) PetscValidLogicalCollectiveInt(dm,NV,4); 37 if (NE > 0) PetscValidLogicalCollectiveInt(dm,NE,5); 38 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); 39 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); 40 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); 41 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); 42 if (NE < 0 || NV < 0) { 43 a[0] = nV; a[1] = nE; 44 ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 45 NV = b[0]; NE = b[1]; 46 } 47 network->nNodes = nV; 48 network->NNodes = NV; 49 network->nEdges = nE; 50 network->NEdges = NE; 51 PetscFunctionReturn(0); 52 } 53 54 #undef __FUNCT__ 55 #define __FUNCT__ "DMNetworkSetEdgeList" 56 /*@ 57 DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network 58 59 Logically collective on DM 60 61 Input Parameters: 62 . edges - list of edges 63 64 Notes: 65 There is no copy involved in this operation, only the pointer is referenced. The edgelist should 66 not be destroyed before the call to DMNetworkLayoutSetUp 67 68 Level: intermediate 69 70 .seealso: DMNetworkCreate, DMNetworkSetSizes 71 @*/ 72 PetscErrorCode DMNetworkSetEdgeList(DM dm, int edgelist[]) 73 { 74 DM_Network *network = (DM_Network*) dm->data; 75 76 PetscFunctionBegin; 77 network->edges = edgelist; 78 PetscFunctionReturn(0); 79 } 80 81 #undef __FUNCT__ 82 #define __FUNCT__ "DMNetworkLayoutSetUp" 83 /*@ 84 DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network 85 86 Collective on DM 87 88 Input Parameters 89 . DM - the dmnetwork object 90 91 Notes: 92 This routine should be called after the network sizes and edgelists have been provided. It creates 93 the bare layout of the network and sets up the network to begin insertion of components. 94 95 All the components should be registered before calling this routine. 96 97 Level: intermediate 98 99 .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList 100 @*/ 101 PetscErrorCode DMNetworkLayoutSetUp(DM dm) 102 { 103 PetscErrorCode ierr; 104 DM_Network *network = (DM_Network*) dm->data; 105 PetscInt dim = 1; /* One dimensional network */ 106 PetscInt numCorners=2; 107 PetscInt spacedim=2; 108 double *vertexcoords=NULL; 109 PetscInt i; 110 PetscInt ndata; 111 112 PetscFunctionBegin; 113 if (network->nNodes) { 114 ierr = PetscCalloc1(numCorners*network->nNodes,&vertexcoords);CHKERRQ(ierr); 115 } 116 ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nNodes,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr); 117 if (network->nNodes) { 118 ierr = PetscFree(vertexcoords);CHKERRQ(ierr); 119 } 120 ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr); 121 ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr); 122 ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr); 123 124 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr); 125 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr); 126 ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr); 127 ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr); 128 129 network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 130 ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr); 131 for (i = network->pStart; i < network->pEnd; i++) { 132 network->header[i].ndata = 0; 133 ndata = network->header[i].ndata; 134 ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 135 network->header[i].offset[ndata] = 0; 136 } 137 ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr); 138 PetscFunctionReturn(0); 139 } 140 141 #undef __FUNCT__ 142 #define __FUNCT__ "DMNetworkRegisterComponent" 143 /*@ 144 DMNetworkRegisterComponent - Registers the network component 145 146 Logically collective on DM 147 148 Input Parameters 149 + dm - the network object 150 . name - the component name 151 - size - the storage size in bytes for this component data 152 153 Output Parameters 154 . key - an integer key that defines the component 155 156 Notes 157 This routine should be called by all processors before calling DMNetworkLayoutSetup(). 158 159 Level: intermediate 160 161 .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 162 @*/ 163 PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key) 164 { 165 PetscErrorCode ierr; 166 DM_Network *network = (DM_Network*) dm->data; 167 DMNetworkComponent *component=&network->component[network->ncomponent]; 168 PetscBool flg=PETSC_FALSE; 169 PetscInt i; 170 171 PetscFunctionBegin; 172 173 for (i=0; i < network->ncomponent; i++) { 174 ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr); 175 if (flg) { 176 *key = i; 177 PetscFunctionReturn(0); 178 } 179 } 180 181 ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr); 182 component->size = size/sizeof(DMNetworkComponentGenericDataType); 183 *key = network->ncomponent; 184 network->ncomponent++; 185 PetscFunctionReturn(0); 186 } 187 188 #undef __FUNCT__ 189 #define __FUNCT__ "DMNetworkGetVertexRange" 190 /*@ 191 DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices. 192 193 Not Collective 194 195 Input Parameters: 196 + dm - The DMNetwork object 197 198 Output Paramters: 199 + vStart - The first vertex point 200 - vEnd - One beyond the last vertex point 201 202 Level: intermediate 203 204 .seealso: DMNetworkGetEdgeRange 205 @*/ 206 PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd) 207 { 208 DM_Network *network = (DM_Network*)dm->data; 209 210 PetscFunctionBegin; 211 if (vStart) *vStart = network->vStart; 212 if (vEnd) *vEnd = network->vEnd; 213 PetscFunctionReturn(0); 214 } 215 216 #undef __FUNCT__ 217 #define __FUNCT__ "DMNetworkGetEdgeRange" 218 /*@ 219 DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges. 220 221 Not Collective 222 223 Input Parameters: 224 + dm - The DMNetwork object 225 226 Output Paramters: 227 + eStart - The first edge point 228 - eEnd - One beyond the last edge point 229 230 Level: intermediate 231 232 .seealso: DMNetworkGetVertexRange 233 @*/ 234 PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd) 235 { 236 DM_Network *network = (DM_Network*)dm->data; 237 238 PetscFunctionBegin; 239 if (eStart) *eStart = network->eStart; 240 if (eEnd) *eEnd = network->eEnd; 241 PetscFunctionReturn(0); 242 } 243 244 #undef __FUNCT__ 245 #define __FUNCT__ "DMNetworkAddComponent" 246 /*@ 247 DMNetworkAddComponent - Adds a network component at the given point (vertex/edge) 248 249 Not Collective 250 251 Input Parameters: 252 + dm - The DMNetwork object 253 . p - vertex/edge point 254 . componentkey - component key returned while registering the component 255 - compvalue - pointer to the data structure for the component 256 257 Level: intermediate 258 259 .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent 260 @*/ 261 PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue) 262 { 263 DM_Network *network = (DM_Network*)dm->data; 264 DMNetworkComponent *component = &network->component[componentkey]; 265 DMNetworkComponentHeader header = &network->header[p]; 266 DMNetworkComponentValue cvalue = &network->cvalue[p]; 267 PetscErrorCode ierr; 268 269 PetscFunctionBegin; 270 header->size[header->ndata] = component->size; 271 ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr); 272 header->key[header->ndata] = componentkey; 273 if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1]; 274 275 cvalue->data[header->ndata] = (void*)compvalue; 276 header->ndata++; 277 PetscFunctionReturn(0); 278 } 279 280 #undef __FUNCT__ 281 #define __FUNCT__ "DMNetworkGetNumComponents" 282 /*@ 283 DMNetworkGetNumComponents - Get the number of components at a vertex/edge 284 285 Not Collective 286 287 Input Parameters: 288 + dm - The DMNetwork object 289 . p - vertex/edge point 290 291 Output Parameters: 292 . numcomponents - Number of components at the vertex/edge 293 294 Level: intermediate 295 296 .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent 297 @*/ 298 PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents) 299 { 300 PetscErrorCode ierr; 301 PetscInt offset; 302 DM_Network *network = (DM_Network*)dm->data; 303 304 PetscFunctionBegin; 305 ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 306 *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 307 PetscFunctionReturn(0); 308 } 309 310 #undef __FUNCT__ 311 #define __FUNCT__ "DMNetworkGetComponentTypeOffset" 312 /*@ 313 DMNetworkGetComponentTypeOffset - Gets the type along with the offset for indexing the 314 component value from the component data array 315 316 Not Collective 317 318 Input Parameters: 319 + dm - The DMNetwork object 320 . p - vertex/edge point 321 - compnum - component number 322 323 Output Parameters: 324 + compkey - the key obtained when registering the component 325 - offset - offset into the component data array associated with the vertex/edge point 326 327 Notes: 328 Typical usage: 329 330 DMNetworkGetComponentDataArray(dm, &arr); 331 DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 332 Loop over vertices or edges 333 DMNetworkGetNumComponents(dm,v,&numcomps); 334 Loop over numcomps 335 DMNetworkGetComponentTypeOffset(dm,v,compnum,&key,&offset); 336 compdata = (UserCompDataType)(arr+offset); 337 338 Level: intermediate 339 340 .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray, 341 @*/ 342 PetscErrorCode DMNetworkGetComponentTypeOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset) 343 { 344 PetscErrorCode ierr; 345 PetscInt offsetp; 346 DMNetworkComponentHeader header; 347 DM_Network *network = (DM_Network*)dm->data; 348 349 PetscFunctionBegin; 350 ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 351 header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 352 *compkey = header->key[compnum]; 353 *offset = offsetp+network->dataheadersize+header->offset[compnum]; 354 PetscFunctionReturn(0); 355 } 356 357 #undef __FUNCT__ 358 #define __FUNCT__ "DMNetworkGetVariableOffset" 359 /*@ 360 DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector. 361 362 Not Collective 363 364 Input Parameters: 365 + dm - The DMNetwork object 366 - p - the edge/vertex point 367 368 Output Parameters: 369 . offset - the offset 370 371 Level: intermediate 372 373 .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 374 @*/ 375 PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset) 376 { 377 PetscErrorCode ierr; 378 DM_Network *network = (DM_Network*)dm->data; 379 380 PetscFunctionBegin; 381 ierr = PetscSectionGetOffset(network->DofSection,p,offset);CHKERRQ(ierr); 382 PetscFunctionReturn(0); 383 } 384 385 #undef __FUNCT__ 386 #define __FUNCT__ "DMNetworkGetVariableGlobalOffset" 387 /*@ 388 DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector. 389 390 Not Collective 391 392 Input Parameters: 393 + dm - The DMNetwork object 394 - p - the edge/vertex point 395 396 Output Parameters: 397 . offsetg - the offset 398 399 Level: intermediate 400 401 .seealso: DMNetworkGetVariableOffset, DMGetLocalVector 402 @*/ 403 PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg) 404 { 405 PetscErrorCode ierr; 406 DM_Network *network = (DM_Network*)dm->data; 407 408 PetscFunctionBegin; 409 ierr = PetscSectionGetOffset(network->GlobalDofSection,p,offsetg);CHKERRQ(ierr); 410 PetscFunctionReturn(0); 411 } 412 413 #undef __FUNCT__ 414 #define __FUNCT__ "DMNetworkAddNumVariables" 415 /*@ 416 DMNetworkAddNumVariables - Add number of variables associated with a given point. 417 418 Not Collective 419 420 Input Parameters: 421 + dm - The DMNetworkObject 422 . p - the vertex/edge point 423 - nvar - number of additional variables 424 425 Level: intermediate 426 427 .seealso: DMNetworkSetNumVariables 428 @*/ 429 PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar) 430 { 431 PetscErrorCode ierr; 432 DM_Network *network = (DM_Network*)dm->data; 433 434 PetscFunctionBegin; 435 ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr); 436 PetscFunctionReturn(0); 437 } 438 439 #undef __FUNCT__ 440 #define __FUNCT__ "DMNetworkGetNumVariables" 441 /*@ 442 DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point. 443 444 Not Collective 445 446 Input Parameters: 447 + dm - The DMNetworkObject 448 - p - the vertex/edge point 449 450 Output Parameters: 451 . nvar - number of variables 452 453 Level: intermediate 454 455 .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables 456 @*/ 457 PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar) 458 { 459 PetscErrorCode ierr; 460 DM_Network *network = (DM_Network*)dm->data; 461 462 PetscFunctionBegin; 463 ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 464 PetscFunctionReturn(0); 465 } 466 467 #undef __FUNCT__ 468 #define __FUNCT__ "DMNetworkSetNumVariables" 469 /*@ 470 DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point. 471 472 Not Collective 473 474 Input Parameters: 475 + dm - The DMNetworkObject 476 . p - the vertex/edge point 477 - nvar - number of variables 478 479 Level: intermediate 480 481 .seealso: DMNetworkAddNumVariables 482 @*/ 483 PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar) 484 { 485 PetscErrorCode ierr; 486 DM_Network *network = (DM_Network*)dm->data; 487 488 PetscFunctionBegin; 489 ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 490 PetscFunctionReturn(0); 491 } 492 493 /* Sets up the array that holds the data for all components and its associated section. This 494 function is called during DMSetUp() */ 495 #undef __FUNCT__ 496 #define __FUNCT__ "DMNetworkComponentSetUp" 497 PetscErrorCode DMNetworkComponentSetUp(DM dm) 498 { 499 PetscErrorCode ierr; 500 DM_Network *network = (DM_Network*)dm->data; 501 PetscInt arr_size; 502 PetscInt p,offset,offsetp; 503 DMNetworkComponentHeader header; 504 DMNetworkComponentValue cvalue; 505 DMNetworkComponentGenericDataType *componentdataarray; 506 PetscInt ncomp, i; 507 508 PetscFunctionBegin; 509 ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr); 510 ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr); 511 ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr); 512 componentdataarray = network->componentdataarray; 513 for (p = network->pStart; p < network->pEnd; p++) { 514 ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 515 /* Copy header */ 516 header = &network->header[p]; 517 ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 518 /* Copy data */ 519 cvalue = &network->cvalue[p]; 520 ncomp = header->ndata; 521 for (i = 0; i < ncomp; i++) { 522 offset = offsetp + network->dataheadersize + header->offset[i]; 523 ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 524 } 525 } 526 PetscFunctionReturn(0); 527 } 528 529 /* Sets up the section for dofs. This routine is called during DMSetUp() */ 530 #undef __FUNCT__ 531 #define __FUNCT__ "DMNetworkVariablesSetUp" 532 PetscErrorCode DMNetworkVariablesSetUp(DM dm) 533 { 534 PetscErrorCode ierr; 535 DM_Network *network = (DM_Network*)dm->data; 536 537 PetscFunctionBegin; 538 ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr); 539 PetscFunctionReturn(0); 540 } 541 542 #undef __FUNCT__ 543 #define __FUNCT__ "DMNetworkGetComponentDataArray" 544 /*@C 545 DMNetworkGetComponentDataArray - Returns the component data array 546 547 Not Collective 548 549 Input Parameters: 550 . dm - The DMNetwork Object 551 552 Output Parameters: 553 . componentdataarray - array that holds data for all components 554 555 Level: intermediate 556 557 .seealso: DMNetworkGetComponentTypeOffset, DMNetworkGetNumComponents 558 @*/ 559 PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray) 560 { 561 DM_Network *network = (DM_Network*)dm->data; 562 563 PetscFunctionBegin; 564 *componentdataarray = network->componentdataarray; 565 PetscFunctionReturn(0); 566 } 567 568 #undef __FUNCT__ 569 #define __FUNCT__ "DMNetworkDistribute" 570 /*@ 571 DMNetworkDistribute - Distributes the network and moves associated component data. 572 573 Collective 574 575 Input Parameter: 576 + oldDM - the original DMNetwork object 577 - overlap - The overlap of partitions, 0 is the default 578 579 Output Parameter: 580 . distDM - the distributed DMNetwork object 581 582 Notes: 583 This routine should be called only when using multiple processors. 584 585 Distributes the network with <overlap>-overlapping partitioning of the edges. 586 587 Level: intermediate 588 589 .seealso: DMNetworkCreate 590 @*/ 591 PetscErrorCode DMNetworkDistribute(DM oldDM, PetscInt overlap,DM *distDM) 592 { 593 PetscErrorCode ierr; 594 DM_Network *oldDMnetwork = (DM_Network*)oldDM->data; 595 PetscSF pointsf; 596 DM newDM; 597 DM_Network *newDMnetwork; 598 599 PetscFunctionBegin; 600 ierr = DMNetworkCreate(PetscObjectComm((PetscObject)oldDM),&newDM);CHKERRQ(ierr); 601 newDMnetwork = (DM_Network*)newDM->data; 602 newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 603 /* Distribute plex dm and dof section */ 604 ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr); 605 /* Distribute dof section */ 606 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DofSection);CHKERRQ(ierr); 607 ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr); 608 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DataSection);CHKERRQ(ierr); 609 /* Distribute data and associated section */ 610 ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr); 611 /* Destroy point SF */ 612 ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr); 613 614 ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr); 615 ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr); 616 ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr); 617 newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart; 618 newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart; 619 newDMnetwork->NNodes = oldDMnetwork->NNodes; 620 newDMnetwork->NEdges = oldDMnetwork->NEdges; 621 /* Set Dof section as the default section for dm */ 622 ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr); 623 ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr); 624 625 *distDM = newDM; 626 PetscFunctionReturn(0); 627 } 628 629 #undef __FUNCT__ 630 #define __FUNCT__ "DMNetworkGetSupportingEdges" 631 /*@C 632 DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point 633 634 Not Collective 635 636 Input Parameters: 637 + dm - The DMNetwork object 638 - p - the vertex point 639 640 Output Paramters: 641 + nedges - number of edges connected to this vertex point 642 - edges - List of edge points 643 644 Level: intermediate 645 646 Fortran Notes: 647 Since it returns an array, this routine is only available in Fortran 90, and you must 648 include petsc.h90 in your code. 649 650 .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes 651 @*/ 652 PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[]) 653 { 654 PetscErrorCode ierr; 655 DM_Network *network = (DM_Network*)dm->data; 656 657 PetscFunctionBegin; 658 ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr); 659 ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr); 660 PetscFunctionReturn(0); 661 } 662 663 #undef __FUNCT__ 664 #define __FUNCT__ "DMNetworkGetConnectedNodes" 665 /*@C 666 DMNetworkGetConnectedNodes - Return the connected vertices for this edge point 667 668 Not Collective 669 670 Input Parameters: 671 + dm - The DMNetwork object 672 - p - the edge point 673 674 Output Paramters: 675 . vertices - vertices connected to this edge 676 677 Level: intermediate 678 679 Fortran Notes: 680 Since it returns an array, this routine is only available in Fortran 90, and you must 681 include petsc.h90 in your code. 682 683 .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges 684 @*/ 685 PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[]) 686 { 687 PetscErrorCode ierr; 688 DM_Network *network = (DM_Network*)dm->data; 689 690 PetscFunctionBegin; 691 ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr); 692 PetscFunctionReturn(0); 693 } 694 695 #undef __FUNCT__ 696 #define __FUNCT__ "DMNetworkIsGhostVertex" 697 /*@ 698 DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex 699 700 Not Collective 701 702 Input Parameters: 703 + dm - The DMNetwork object 704 . p - the vertex point 705 706 Output Parameter: 707 . isghost - TRUE if the vertex is a ghost point 708 709 Level: intermediate 710 711 .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange 712 @*/ 713 PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost) 714 { 715 PetscErrorCode ierr; 716 DM_Network *network = (DM_Network*)dm->data; 717 PetscInt offsetg; 718 PetscSection sectiong; 719 720 PetscFunctionBegin; 721 *isghost = PETSC_FALSE; 722 ierr = DMGetDefaultGlobalSection(network->plex,§iong);CHKERRQ(ierr); 723 ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr); 724 if (offsetg < 0) *isghost = PETSC_TRUE; 725 PetscFunctionReturn(0); 726 } 727 728 #undef __FUNCT__ 729 #define __FUNCT__ "DMSetUp_Network" 730 PetscErrorCode DMSetUp_Network(DM dm) 731 { 732 PetscErrorCode ierr; 733 DM_Network *network=(DM_Network*)dm->data; 734 735 PetscFunctionBegin; 736 ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr); 737 ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr); 738 739 ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr); 740 ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr); 741 PetscFunctionReturn(0); 742 } 743 744 #undef __FUNCT__ 745 #define __FUNCT__ "DMNetworkHasJacobian" 746 /*@ 747 DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices 748 -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)? 749 750 Collective 751 752 Input Parameters: 753 + dm - The DMNetwork object 754 . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges 755 - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices 756 757 Level: intermediate 758 759 @*/ 760 PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg) 761 { 762 DM_Network *network=(DM_Network*)dm->data; 763 764 PetscFunctionBegin; 765 network->userEdgeJacobian = eflg; 766 network->userVertexJacobian = vflg; 767 PetscFunctionReturn(0); 768 } 769 770 #undef __FUNCT__ 771 #define __FUNCT__ "DMNetworkEdgeSetMatrix" 772 /*@ 773 DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network 774 775 Not Collective 776 777 Input Parameters: 778 + dm - The DMNetwork object 779 . p - the edge point 780 - J - array (size = 3) of Jacobian submatrices for this edge point: 781 J[0]: this edge 782 J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedNodes() 783 784 Level: intermediate 785 786 .seealso: DMNetworkVertexSetMatrix 787 @*/ 788 PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[]) 789 { 790 PetscErrorCode ierr; 791 DM_Network *network=(DM_Network*)dm->data; 792 793 PetscFunctionBegin; 794 if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix"); 795 if (!network->Je) { 796 ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr); 797 } 798 network->Je[3*p] = J[0]; 799 network->Je[3*p+1] = J[1]; 800 network->Je[3*p+2] = J[2]; 801 PetscFunctionReturn(0); 802 } 803 804 #undef __FUNCT__ 805 #define __FUNCT__ "DMNetworkVertexSetMatrix" 806 /*@ 807 DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network 808 809 Not Collective 810 811 Input Parameters: 812 + dm - The DMNetwork object 813 . p - the vertex point 814 - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point: 815 J[0]: this vertex 816 J[1+2*i]: i-th supporting edge 817 J[1+2*i+1]: i-th connected vertex 818 819 Level: intermediate 820 821 .seealso: DMNetworkEdgeSetMatrix 822 @*/ 823 PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[]) 824 { 825 PetscErrorCode ierr; 826 DM_Network *network=(DM_Network*)dm->data; 827 PetscInt i,*vptr,nedges,vStart,vEnd; 828 const PetscInt *edges; 829 830 PetscFunctionBegin; 831 if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix"); 832 833 ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 834 835 if (!network->Jv) { 836 PetscInt nNodes = network->nNodes,nedges_total; 837 /* count nvertex_total */ 838 nedges_total = 0; 839 ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 840 ierr = PetscMalloc1(nNodes+1,&vptr);CHKERRQ(ierr); 841 842 vptr[0] = 0; 843 for (i=0; i<nNodes; i++) { 844 ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr); 845 nedges_total += nedges; 846 vptr[i+1] = vptr[i] + 2*nedges + 1; 847 } 848 849 ierr = PetscCalloc1(2*nedges_total+nNodes,&network->Jv);CHKERRQ(ierr); 850 network->Jvptr = vptr; 851 } 852 853 vptr = network->Jvptr; 854 network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */ 855 856 /* Set Jacobian for each supporting edge and connected vertex */ 857 ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr); 858 for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i]; 859 PetscFunctionReturn(0); 860 } 861 862 #undef __FUNCT__ 863 #define __FUNCT__ "MatSetDenseblock_private" 864 PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 865 { 866 PetscErrorCode ierr; 867 PetscInt j,*cols; 868 PetscScalar *zeros; 869 870 PetscFunctionBegin; 871 ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr); 872 for (j=0; j<ncols; j++) cols[j] = j+ cstart; 873 ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr); 874 ierr = PetscFree2(cols,zeros);CHKERRQ(ierr); 875 PetscFunctionReturn(0); 876 } 877 878 #undef __FUNCT__ 879 #define __FUNCT__ "MatSetUserblock_private" 880 PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 881 { 882 PetscErrorCode ierr; 883 PetscInt j,M,N,row,col,ncols_u; 884 const PetscInt *cols; 885 PetscScalar zero=0.0; 886 887 PetscFunctionBegin; 888 ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr); 889 if (nrows != M || ncols != N) SETERRQ4(PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%D by %D must equal %D by %D",nrows,ncols,M,N); 890 891 for (row=0; row<nrows; row++) { 892 ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 893 for (j=0; j<ncols_u; j++) { 894 col = cols[j] + cstart; 895 ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr); 896 } 897 ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 898 } 899 PetscFunctionReturn(0); 900 } 901 902 #undef __FUNCT__ 903 #define __FUNCT__ "MatSetblock_private" 904 PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 905 { 906 PetscErrorCode ierr; 907 908 PetscFunctionBegin; 909 if (Ju) { 910 ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 911 } else { 912 ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 913 } 914 PetscFunctionReturn(0); 915 } 916 917 #undef __FUNCT__ 918 #define __FUNCT__ "DMCreateMatrix_Network" 919 PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J) 920 { 921 PetscErrorCode ierr; 922 DM_Network *network = (DM_Network*) dm->data; 923 PetscInt eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize; 924 PetscInt cstart,ncols,j,e,v,*dnz,*onz,*dnzu,*onzu; 925 PetscBool ghost; 926 Mat Juser; 927 PetscSection sectionGlobal; 928 PetscInt nedges,*vptr,vc; 929 const PetscInt *edges,*cone; 930 931 PetscFunctionBegin; 932 if (!network->userEdgeJacobian && !network->userVertexJacobian) { 933 /* user does not provide Jacobian blocks */ 934 ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr); 935 ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 936 PetscFunctionReturn(0); 937 } 938 939 ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr); 940 ierr = DMGetDefaultGlobalSection(network->plex,§ionGlobal);CHKERRQ(ierr); 941 ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr); 942 ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 943 944 ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr); 945 ierr = MatSetFromOptions(*J);CHKERRQ(ierr); 946 947 /* Preallocation - submatrix for an element (edge/vertex) is allocated as a dense block, 948 see DMCreateMatrix_Plex() */ 949 ierr = PetscCalloc4(localSize,&dnz,localSize,&onz,localSize,&dnzu,localSize,&onzu);CHKERRQ(ierr); 950 ierr = DMPlexPreallocateOperator(network->plex,1,dnz,onz,dnzu,onzu,*J,PETSC_FALSE);CHKERRQ(ierr); 951 ierr = PetscFree4(dnz,onz,dnzu,onzu);CHKERRQ(ierr); 952 953 /* Set matrix entries for edges */ 954 /*------------------------------*/ 955 ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr); 956 for (e=eStart; e<eEnd; e++) { 957 /* Get row indices */ 958 ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 959 ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 960 ierr = PetscMalloc1(nrows,&rows);CHKERRQ(ierr); 961 for (j=0; j<nrows; j++) rows[j] = j + rstart; 962 963 /* Set matrix entries for conntected vertices */ 964 ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr); 965 for (v=0; v<2; v++) { 966 ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr); 967 ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr); 968 if (ghost) cstart = -(cstart + 1); /* Convert to actual global offset for ghost nodes */ 969 ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 970 971 Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 972 ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 973 } 974 975 /* Set matrix entries for edge self */ 976 cstart = rstart; 977 Juser = network->Je[3*e]; /* Jacobian(e,e) */ 978 ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr); 979 980 ierr = PetscFree(rows);CHKERRQ(ierr); 981 } 982 983 /* Set matrix entries for vertices */ 984 /*---------------------------------*/ 985 ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 986 for (v=vStart; v<vEnd; v++) { 987 /* Get row indices */ 988 ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 989 ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 990 if (ghost) rstart = -(rstart + 1); /* Convert to actual global offset for ghost nodes */ 991 ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 992 993 ierr = PetscMalloc1(nrows,&rows);CHKERRQ(ierr); 994 for (j=0; j<nrows; j++) rows[j] = j + rstart; 995 996 /* Get supporting edges and connected vertices */ 997 ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 998 999 for (e=0; e<nedges; e++) { 1000 /* Supporting edges */ 1001 ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1002 ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1003 1004 vptr = network->Jvptr; 1005 Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 1006 ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1007 1008 /* Connected vertices */ 1009 ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr); 1010 vc = (v == cone[0]) ? cone[1]:cone[0]; 1011 1012 ierr = DMNetworkIsGhostVertex(dm,vc,&ghost);CHKERRQ(ierr); 1013 ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr); 1014 if (ghost) cstart = -(cstart + 1); /* Convert to actual global offset for ghost nodes */ 1015 ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1016 1017 Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 1018 ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr); 1019 } 1020 1021 /* Set matrix entries for vertex self */ 1022 ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1023 if (!ghost) { 1024 ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 1025 Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 1026 ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr); 1027 } 1028 ierr = PetscFree(rows);CHKERRQ(ierr); 1029 } 1030 ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1031 ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1032 #if 0 1033 printf("\nMatrix J:\n"); 1034 ierr = MatView(*J,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1035 #endif 1036 ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 1037 PetscFunctionReturn(0); 1038 } 1039 1040 #undef __FUNCT__ 1041 #define __FUNCT__ "DMDestroy_Network" 1042 PetscErrorCode DMDestroy_Network(DM dm) 1043 { 1044 PetscErrorCode ierr; 1045 DM_Network *network = (DM_Network*) dm->data; 1046 PetscInt i; 1047 1048 PetscFunctionBegin; 1049 if (--network->refct > 0) PetscFunctionReturn(0); 1050 if (network->Je) { 1051 for (i=0; i<3*network->nEdges; i++) { 1052 ierr = MatDestroy(&network->Je[i]);CHKERRQ(ierr); 1053 } 1054 ierr = PetscFree(network->Je);CHKERRQ(ierr); 1055 } 1056 if (network->Jv) { 1057 PetscInt v,vStart,vEnd,nedges,*vptr=network->Jvptr; 1058 const PetscInt *edges; 1059 ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 1060 for (v=vStart; v<vEnd; v++) { 1061 ierr = MatDestroy(&network->Jv[vptr[v-vStart]]);CHKERRQ(ierr); 1062 ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1063 for (i=0; i<2*nedges; i++) { 1064 ierr = MatDestroy(&network->Jv[vptr[v-vStart]+i+1]);CHKERRQ(ierr); 1065 } 1066 } 1067 ierr = PetscFree(network->Jvptr);CHKERRQ(ierr); 1068 ierr = PetscFree(network->Jv);CHKERRQ(ierr); 1069 } 1070 ierr = DMDestroy(&network->plex);CHKERRQ(ierr); 1071 network->edges = NULL; 1072 ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr); 1073 ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr); 1074 1075 ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr); 1076 ierr = PetscFree(network->cvalue);CHKERRQ(ierr); 1077 ierr = PetscFree(network->header);CHKERRQ(ierr); 1078 ierr = PetscFree(network);CHKERRQ(ierr); 1079 PetscFunctionReturn(0); 1080 } 1081 1082 #undef __FUNCT__ 1083 #define __FUNCT__ "DMView_Network" 1084 PetscErrorCode DMView_Network(DM dm, PetscViewer viewer) 1085 { 1086 PetscErrorCode ierr; 1087 DM_Network *network = (DM_Network*) dm->data; 1088 1089 PetscFunctionBegin; 1090 ierr = DMView(network->plex,viewer);CHKERRQ(ierr); 1091 PetscFunctionReturn(0); 1092 } 1093 1094 #undef __FUNCT__ 1095 #define __FUNCT__ "DMGlobalToLocalBegin_Network" 1096 PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l) 1097 { 1098 PetscErrorCode ierr; 1099 DM_Network *network = (DM_Network*) dm->data; 1100 1101 PetscFunctionBegin; 1102 ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr); 1103 PetscFunctionReturn(0); 1104 } 1105 1106 #undef __FUNCT__ 1107 #define __FUNCT__ "DMGlobalToLocalEnd_Network" 1108 PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l) 1109 { 1110 PetscErrorCode ierr; 1111 DM_Network *network = (DM_Network*) dm->data; 1112 1113 PetscFunctionBegin; 1114 ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr); 1115 PetscFunctionReturn(0); 1116 } 1117 1118 #undef __FUNCT__ 1119 #define __FUNCT__ "DMLocalToGlobalBegin_Network" 1120 PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g) 1121 { 1122 PetscErrorCode ierr; 1123 DM_Network *network = (DM_Network*) dm->data; 1124 1125 PetscFunctionBegin; 1126 ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr); 1127 PetscFunctionReturn(0); 1128 } 1129 1130 #undef __FUNCT__ 1131 #define __FUNCT__ "DMLocalToGlobalEnd_Network" 1132 PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g) 1133 { 1134 PetscErrorCode ierr; 1135 DM_Network *network = (DM_Network*) dm->data; 1136 1137 PetscFunctionBegin; 1138 ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr); 1139 PetscFunctionReturn(0); 1140 } 1141