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