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 of Jacobian submatrices: 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 DMNetworkCreateJacobian() collectively before calling DMNetworkElementSetMatrix"); 795 if (!network->Je) { 796 ierr = PetscCalloc1(network->nEdges,&network->Je);CHKERRQ(ierr); 797 } 798 network->Je[p] = J[0]; 799 PetscFunctionReturn(0); 800 } 801 802 #undef __FUNCT__ 803 #define __FUNCT__ "DMNetworkVetexSetMatrix" 804 /*@ 805 DMNetworkVetexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network 806 807 Not Collective 808 809 Input Parameters: 810 + dm - The DMNetwork object 811 . p - the vertex point 812 - J - array of Jacobian submatrices: 813 814 Level: intermediate 815 816 .seealso: DMNetworkEdgeSetMatrix 817 @*/ 818 PetscErrorCode DMNetworkVetexSetMatrix(DM dm,PetscInt p,Mat J[]) 819 { 820 PetscErrorCode ierr; 821 DM_Network *network=(DM_Network*)dm->data; 822 823 PetscFunctionBegin; 824 if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkCreateJacobian() collectively before calling DMNetworkElementSetMatrix"); 825 if (!network->Jv) { 826 ierr = PetscCalloc1(network->nNodes,&network->Jv);CHKERRQ(ierr); 827 } 828 network->Jv[p] = J[0]; 829 PetscFunctionReturn(0); 830 } 831 832 #undef __FUNCT__ 833 #define __FUNCT__ "DMCreateMatrix_Network" 834 PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J) 835 { 836 PetscErrorCode ierr; 837 DM_Network *network = (DM_Network*) dm->data; 838 PetscInt eStart,eEnd,vStart,vEnd,rstart,rend,row,row_e,nrows,localSize; 839 PetscInt cstart,ncols,col,j,e,v,*dnz,*onz,*dnzu,*onzu; 840 const PetscInt *cols; 841 PetscScalar *zeros,zero = 0.0; 842 PetscBool ghost; 843 Mat Je; 844 PetscSection sectionGlobal; 845 PetscInt nedges; 846 const PetscInt *edges; 847 848 PetscFunctionBegin; 849 if (!network->userEdgeJacobian && !network->userVertexJacobian) { /* user does not provide Jacobian blocks */ 850 ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr); 851 ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 852 PetscFunctionReturn(0); 853 } 854 855 ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr); 856 ierr = DMGetDefaultGlobalSection(network->plex,§ionGlobal);CHKERRQ(ierr); 857 ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr); 858 ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 859 860 ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr); 861 ierr = MatSetFromOptions(*J);CHKERRQ(ierr); 862 863 /* Preallocation - submatrix for an element (edge/vertex) is allocated as a dense block, see DMCreateMatrix_Plex() */ 864 ierr = PetscCalloc4(localSize,&dnz,localSize,&onz,localSize,&dnzu,localSize,&onzu);CHKERRQ(ierr); 865 ierr = DMPlexPreallocateOperator(network->plex,1,dnz,onz,dnzu,onzu,*J,PETSC_FALSE);CHKERRQ(ierr); 866 ierr = PetscFree4(dnz,onz,dnzu,onzu);CHKERRQ(ierr); 867 868 /* Set matrix entries for edges */ 869 /*------------------------------*/ 870 ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr); 871 for (e=eStart; e<eEnd; e++) { 872 /* Get row indices */ 873 ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 874 ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 875 876 PetscInt rows[nrows],*cols_tmp; 877 for (j=0; j<nrows; j++) rows[j] = j + rstart; 878 879 /* Set matrix entries for conntected vertices */ 880 const PetscInt *cone; 881 ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr); 882 883 for (v=0; v<2; v++) { 884 ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr); 885 ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr); 886 if (ghost) cstart = -(cstart + 1); /* Convert to actual global offset for ghost nodes */ 887 ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 888 ierr = PetscCalloc2(ncols,&cols_tmp,nrows*ncols,&zeros);CHKERRQ(ierr); 889 for (j=0; j<ncols; j++) cols_tmp[j] = j+ cstart; 890 891 ierr = MatSetValues(*J,nrows,rows,ncols,cols_tmp,zeros,INSERT_VALUES);CHKERRQ(ierr); 892 ierr = PetscFree2(cols_tmp,zeros);CHKERRQ(ierr); 893 } 894 895 /* Set matrix entries for edge self */ 896 Je = network->Je[e]; 897 if (Je) { /* use user-provided Je */ 898 PetscInt M,N; 899 ierr = MatGetSize(Je,&M,&N);CHKERRQ(ierr); 900 if (nrows != M || nrows != N) SETERRQ3(PetscObjectComm((PetscObject)Je),PETSC_ERR_USER,"%D must equal %D and %D",rend-rstart,M,N); 901 902 rend = rstart + nrows; 903 for (row=rstart; row<rend; row++) { 904 row_e = row - rstart; 905 ierr = MatGetRow(Je,row_e,&ncols,&cols,NULL);CHKERRQ(ierr); 906 for (j=0; j<ncols; j++) { 907 col = cols[j] + rstart; 908 ierr = MatSetValues(*J,1,&row,1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr); 909 } 910 ierr = MatRestoreRow(Je,row_e,&ncols,&cols,NULL);CHKERRQ(ierr); 911 } 912 } else { /* set a dense block */ 913 PetscInt *cols_tmp,cstart; 914 ierr = PetscCalloc2(nrows,&cols_tmp,nrows*nrows,&zeros);CHKERRQ(ierr); 915 ierr = DMNetworkGetVariableGlobalOffset(dm,e,&cstart);CHKERRQ(ierr); 916 for (j=0; j<nrows; j++) cols_tmp[j] = j+ cstart; 917 ierr = MatSetValues(*J,nrows,rows,nrows,cols_tmp,zeros,INSERT_VALUES);CHKERRQ(ierr); 918 ierr = PetscFree2(cols_tmp,zeros);CHKERRQ(ierr); 919 } 920 } 921 922 /* Set matrix entries for vertices */ 923 /*---------------------------------*/ 924 ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 925 for (v=vStart; v<vEnd; v++) { 926 /* Get row indices */ 927 ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 928 ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 929 if (ghost) rstart = -(rstart + 1); /* Convert to actual global offset for ghost nodes */ 930 ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 931 932 PetscInt rows[nrows]; 933 for (j=0; j<nrows; j++) rows[j] = j + rstart; 934 935 /* Get supporting edges and connected vertices */ 936 ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 937 938 for (e=0; e<nedges; e++) { 939 /* Supporting edges */ 940 PetscInt *cols_tmp; 941 ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 942 ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 943 944 ierr = PetscCalloc2(ncols,&cols_tmp,nrows*ncols,&zeros);CHKERRQ(ierr); 945 for (j=0; j<ncols; j++) cols_tmp[j] = j+ cstart; 946 ierr = MatSetValues(*J,nrows,rows,ncols,cols_tmp,zeros,INSERT_VALUES);CHKERRQ(ierr); 947 ierr = PetscFree2(cols_tmp,zeros);CHKERRQ(ierr); 948 949 /* Connected vertices */ 950 const PetscInt *cone; 951 PetscInt vc; 952 ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr); 953 vc = (v == cone[0]) ? cone[1]:cone[0]; 954 955 ierr = DMNetworkIsGhostVertex(dm,vc,&ghost);CHKERRQ(ierr); 956 ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr); 957 if (ghost) cstart = -(cstart + 1); /* Convert to actual global offset for ghost nodes */ 958 ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 959 960 ierr = PetscCalloc2(ncols,&cols_tmp,nrows*ncols,&zeros);CHKERRQ(ierr); 961 for (j=0; j<ncols; j++) cols_tmp[j] = j+ cstart; 962 ierr = MatSetValues(*J,nrows,rows,ncols,cols_tmp,zeros,INSERT_VALUES);CHKERRQ(ierr); 963 ierr = PetscFree2(cols_tmp,zeros);CHKERRQ(ierr); 964 } 965 966 /* Set matrix entries for vertex self */ 967 ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 968 if (!ghost) { 969 PetscInt *cols_tmp,cstart; 970 ierr = PetscCalloc2(nrows,&cols_tmp,nrows*nrows,&zeros);CHKERRQ(ierr); 971 ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 972 for (j=0; j<nrows; j++) cols_tmp[j] = j+ cstart; 973 ierr = MatSetValues(*J,nrows,rows,nrows,cols_tmp,zeros,INSERT_VALUES);CHKERRQ(ierr); 974 ierr = PetscFree2(cols_tmp,zeros);CHKERRQ(ierr); 975 } 976 } 977 ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 978 ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 979 #if 0 980 printf("\nMatrix J:\n"); 981 ierr = MatView(*J,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 982 #endif 983 ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 984 PetscFunctionReturn(0); 985 } 986 987 #undef __FUNCT__ 988 #define __FUNCT__ "DMDestroy_Network" 989 PetscErrorCode DMDestroy_Network(DM dm) 990 { 991 PetscErrorCode ierr; 992 DM_Network *network = (DM_Network*) dm->data; 993 PetscInt i; 994 995 PetscFunctionBegin; 996 if (--network->refct > 0) PetscFunctionReturn(0); 997 if (network->Je) { 998 for (i=0; i<network->nEdges; i++) { 999 ierr = MatDestroy(&network->Je[i]);CHKERRQ(ierr); 1000 } 1001 ierr = PetscFree(network->Je);CHKERRQ(ierr); 1002 } 1003 if (network->Jv) { 1004 for (i=0; i<network->nNodes; i++) { 1005 ierr = MatDestroy(&network->Jv[i]);CHKERRQ(ierr); 1006 } 1007 ierr = PetscFree(network->Jv);CHKERRQ(ierr); 1008 } 1009 ierr = DMDestroy(&network->plex);CHKERRQ(ierr); 1010 network->edges = NULL; 1011 ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr); 1012 ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr); 1013 1014 ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr); 1015 ierr = PetscFree(network->cvalue);CHKERRQ(ierr); 1016 ierr = PetscFree(network->header);CHKERRQ(ierr); 1017 ierr = PetscFree(network);CHKERRQ(ierr); 1018 PetscFunctionReturn(0); 1019 } 1020 1021 #undef __FUNCT__ 1022 #define __FUNCT__ "DMView_Network" 1023 PetscErrorCode DMView_Network(DM dm, PetscViewer viewer) 1024 { 1025 PetscErrorCode ierr; 1026 DM_Network *network = (DM_Network*) dm->data; 1027 1028 PetscFunctionBegin; 1029 ierr = DMView(network->plex,viewer);CHKERRQ(ierr); 1030 PetscFunctionReturn(0); 1031 } 1032 1033 #undef __FUNCT__ 1034 #define __FUNCT__ "DMGlobalToLocalBegin_Network" 1035 PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l) 1036 { 1037 PetscErrorCode ierr; 1038 DM_Network *network = (DM_Network*) dm->data; 1039 1040 PetscFunctionBegin; 1041 ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr); 1042 PetscFunctionReturn(0); 1043 } 1044 1045 #undef __FUNCT__ 1046 #define __FUNCT__ "DMGlobalToLocalEnd_Network" 1047 PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l) 1048 { 1049 PetscErrorCode ierr; 1050 DM_Network *network = (DM_Network*) dm->data; 1051 1052 PetscFunctionBegin; 1053 ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr); 1054 PetscFunctionReturn(0); 1055 } 1056 1057 #undef __FUNCT__ 1058 #define __FUNCT__ "DMLocalToGlobalBegin_Network" 1059 PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g) 1060 { 1061 PetscErrorCode ierr; 1062 DM_Network *network = (DM_Network*) dm->data; 1063 1064 PetscFunctionBegin; 1065 ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr); 1066 PetscFunctionReturn(0); 1067 } 1068 1069 #undef __FUNCT__ 1070 #define __FUNCT__ "DMLocalToGlobalEnd_Network" 1071 PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g) 1072 { 1073 PetscErrorCode ierr; 1074 DM_Network *network = (DM_Network*) dm->data; 1075 1076 PetscFunctionBegin; 1077 ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr); 1078 PetscFunctionReturn(0); 1079 } 1080