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