1 #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 2 #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabel.h" I*/ 3 #include <petsc/private/petscdsimpl.h> /*I "petscds.h" I*/ 4 #include <petscdmplex.h> 5 #include <petscsf.h> 6 #include <petscds.h> 7 8 PetscClassId DM_CLASSID; 9 PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal, DM_LocalToLocal, DM_LocatePoints, DM_Coarsen, DM_Refine, DM_CreateInterpolation, DM_CreateRestriction; 10 11 const char *const DMBoundaryTypes[] = {"NONE","GHOSTED","MIRROR","PERIODIC","TWIST","DM_BOUNDARY_",0}; 12 13 /*@ 14 DMCreate - Creates an empty DM object. The type can then be set with DMSetType(). 15 16 If you never call DMSetType() it will generate an 17 error when you try to use the vector. 18 19 Collective on MPI_Comm 20 21 Input Parameter: 22 . comm - The communicator for the DM object 23 24 Output Parameter: 25 . dm - The DM object 26 27 Level: beginner 28 29 .seealso: DMSetType(), DMDA, DMSLICED, DMCOMPOSITE, DMPLEX, DMMOAB, DMNETWORK 30 @*/ 31 PetscErrorCode DMCreate(MPI_Comm comm,DM *dm) 32 { 33 DM v; 34 PetscErrorCode ierr; 35 36 PetscFunctionBegin; 37 PetscValidPointer(dm,2); 38 *dm = NULL; 39 ierr = PetscSysInitializePackage();CHKERRQ(ierr); 40 ierr = VecInitializePackage();CHKERRQ(ierr); 41 ierr = MatInitializePackage();CHKERRQ(ierr); 42 ierr = DMInitializePackage();CHKERRQ(ierr); 43 44 ierr = PetscHeaderCreate(v, DM_CLASSID, "DM", "Distribution Manager", "DM", comm, DMDestroy, DMView);CHKERRQ(ierr); 45 46 v->ltogmap = NULL; 47 v->bs = 1; 48 v->coloringtype = IS_COLORING_GLOBAL; 49 ierr = PetscSFCreate(comm, &v->sf);CHKERRQ(ierr); 50 ierr = PetscSFCreate(comm, &v->defaultSF);CHKERRQ(ierr); 51 v->labels = NULL; 52 v->depthLabel = NULL; 53 v->defaultSection = NULL; 54 v->defaultGlobalSection = NULL; 55 v->defaultConstraintSection = NULL; 56 v->defaultConstraintMat = NULL; 57 v->L = NULL; 58 v->maxCell = NULL; 59 v->bdtype = NULL; 60 v->dimEmbed = PETSC_DEFAULT; 61 { 62 PetscInt i; 63 for (i = 0; i < 10; ++i) { 64 v->nullspaceConstructors[i] = NULL; 65 } 66 } 67 ierr = PetscDSCreate(comm, &v->prob);CHKERRQ(ierr); 68 v->dmBC = NULL; 69 v->coarseMesh = NULL; 70 v->outputSequenceNum = -1; 71 v->outputSequenceVal = 0.0; 72 ierr = DMSetVecType(v,VECSTANDARD);CHKERRQ(ierr); 73 ierr = DMSetMatType(v,MATAIJ);CHKERRQ(ierr); 74 ierr = PetscNew(&(v->labels));CHKERRQ(ierr); 75 v->labels->refct = 1; 76 *dm = v; 77 PetscFunctionReturn(0); 78 } 79 80 /*@ 81 DMClone - Creates a DM object with the same topology as the original. 82 83 Collective on MPI_Comm 84 85 Input Parameter: 86 . dm - The original DM object 87 88 Output Parameter: 89 . newdm - The new DM object 90 91 Level: beginner 92 93 .keywords: DM, topology, create 94 @*/ 95 PetscErrorCode DMClone(DM dm, DM *newdm) 96 { 97 PetscSF sf; 98 Vec coords; 99 void *ctx; 100 PetscInt dim, cdim; 101 PetscErrorCode ierr; 102 103 PetscFunctionBegin; 104 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 105 PetscValidPointer(newdm,2); 106 ierr = DMCreate(PetscObjectComm((PetscObject) dm), newdm);CHKERRQ(ierr); 107 ierr = PetscFree((*newdm)->labels);CHKERRQ(ierr); 108 dm->labels->refct++; 109 (*newdm)->labels = dm->labels; 110 (*newdm)->depthLabel = dm->depthLabel; 111 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 112 ierr = DMSetDimension(*newdm, dim);CHKERRQ(ierr); 113 if (dm->ops->clone) { 114 ierr = (*dm->ops->clone)(dm, newdm);CHKERRQ(ierr); 115 } 116 (*newdm)->setupcalled = dm->setupcalled; 117 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 118 ierr = DMSetPointSF(*newdm, sf);CHKERRQ(ierr); 119 ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); 120 ierr = DMSetApplicationContext(*newdm, ctx);CHKERRQ(ierr); 121 if (dm->coordinateDM) { 122 DM ncdm; 123 PetscSection cs; 124 PetscInt pEnd = -1, pEndMax = -1; 125 126 ierr = DMGetDefaultSection(dm->coordinateDM, &cs);CHKERRQ(ierr); 127 if (cs) {ierr = PetscSectionGetChart(cs, NULL, &pEnd);CHKERRQ(ierr);} 128 ierr = MPI_Allreduce(&pEnd,&pEndMax,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 129 if (pEndMax >= 0) { 130 ierr = DMClone(dm->coordinateDM, &ncdm);CHKERRQ(ierr); 131 ierr = DMSetDefaultSection(ncdm, cs);CHKERRQ(ierr); 132 ierr = DMSetCoordinateDM(*newdm, ncdm);CHKERRQ(ierr); 133 ierr = DMDestroy(&ncdm);CHKERRQ(ierr); 134 } 135 } 136 ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr); 137 ierr = DMSetCoordinateDim(*newdm, cdim);CHKERRQ(ierr); 138 ierr = DMGetCoordinatesLocal(dm, &coords);CHKERRQ(ierr); 139 if (coords) { 140 ierr = DMSetCoordinatesLocal(*newdm, coords);CHKERRQ(ierr); 141 } else { 142 ierr = DMGetCoordinates(dm, &coords);CHKERRQ(ierr); 143 if (coords) {ierr = DMSetCoordinates(*newdm, coords);CHKERRQ(ierr);} 144 } 145 { 146 PetscBool isper; 147 const PetscReal *maxCell, *L; 148 const DMBoundaryType *bd; 149 ierr = DMGetPeriodicity(dm, &isper, &maxCell, &L, &bd);CHKERRQ(ierr); 150 ierr = DMSetPeriodicity(*newdm, isper, maxCell, L, bd);CHKERRQ(ierr); 151 } 152 PetscFunctionReturn(0); 153 } 154 155 /*@C 156 DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() 157 158 Logically Collective on DM 159 160 Input Parameter: 161 + da - initial distributed array 162 . ctype - the vector type, currently either VECSTANDARD or VECCUSP 163 164 Options Database: 165 . -dm_vec_type ctype 166 167 Level: intermediate 168 169 .seealso: DMCreate(), DMDestroy(), DM, DMDAInterpolationType, VecType, DMGetVecType() 170 @*/ 171 PetscErrorCode DMSetVecType(DM da,VecType ctype) 172 { 173 PetscErrorCode ierr; 174 175 PetscFunctionBegin; 176 PetscValidHeaderSpecific(da,DM_CLASSID,1); 177 ierr = PetscFree(da->vectype);CHKERRQ(ierr); 178 ierr = PetscStrallocpy(ctype,(char**)&da->vectype);CHKERRQ(ierr); 179 PetscFunctionReturn(0); 180 } 181 182 /*@C 183 DMGetVecType - Gets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() 184 185 Logically Collective on DM 186 187 Input Parameter: 188 . da - initial distributed array 189 190 Output Parameter: 191 . ctype - the vector type 192 193 Level: intermediate 194 195 .seealso: DMCreate(), DMDestroy(), DM, DMDAInterpolationType, VecType 196 @*/ 197 PetscErrorCode DMGetVecType(DM da,VecType *ctype) 198 { 199 PetscFunctionBegin; 200 PetscValidHeaderSpecific(da,DM_CLASSID,1); 201 *ctype = da->vectype; 202 PetscFunctionReturn(0); 203 } 204 205 /*@ 206 VecGetDM - Gets the DM defining the data layout of the vector 207 208 Not collective 209 210 Input Parameter: 211 . v - The Vec 212 213 Output Parameter: 214 . dm - The DM 215 216 Level: intermediate 217 218 .seealso: VecSetDM(), DMGetLocalVector(), DMGetGlobalVector(), DMSetVecType() 219 @*/ 220 PetscErrorCode VecGetDM(Vec v, DM *dm) 221 { 222 PetscErrorCode ierr; 223 224 PetscFunctionBegin; 225 PetscValidHeaderSpecific(v,VEC_CLASSID,1); 226 PetscValidPointer(dm,2); 227 ierr = PetscObjectQuery((PetscObject) v, "__PETSc_dm", (PetscObject*) dm);CHKERRQ(ierr); 228 PetscFunctionReturn(0); 229 } 230 231 /*@ 232 VecSetDM - Sets the DM defining the data layout of the vector. 233 234 Not collective 235 236 Input Parameters: 237 + v - The Vec 238 - dm - The DM 239 240 Note: This is NOT the same as DMCreateGlobalVector() since it does not change the view methods or perform other customization, but merely sets the DM member. 241 242 Level: intermediate 243 244 .seealso: VecGetDM(), DMGetLocalVector(), DMGetGlobalVector(), DMSetVecType() 245 @*/ 246 PetscErrorCode VecSetDM(Vec v, DM dm) 247 { 248 PetscErrorCode ierr; 249 250 PetscFunctionBegin; 251 PetscValidHeaderSpecific(v,VEC_CLASSID,1); 252 if (dm) PetscValidHeaderSpecific(dm,DM_CLASSID,2); 253 ierr = PetscObjectCompose((PetscObject) v, "__PETSc_dm", (PetscObject) dm);CHKERRQ(ierr); 254 PetscFunctionReturn(0); 255 } 256 257 /*@C 258 DMSetMatType - Sets the type of matrix created with DMCreateMatrix() 259 260 Logically Collective on DM 261 262 Input Parameter: 263 + dm - the DM context 264 - ctype - the matrix type 265 266 Options Database: 267 . -dm_mat_type ctype 268 269 Level: intermediate 270 271 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType, DMGetMatType() 272 @*/ 273 PetscErrorCode DMSetMatType(DM dm,MatType ctype) 274 { 275 PetscErrorCode ierr; 276 277 PetscFunctionBegin; 278 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 279 ierr = PetscFree(dm->mattype);CHKERRQ(ierr); 280 ierr = PetscStrallocpy(ctype,(char**)&dm->mattype);CHKERRQ(ierr); 281 PetscFunctionReturn(0); 282 } 283 284 /*@C 285 DMGetMatType - Gets the type of matrix created with DMCreateMatrix() 286 287 Logically Collective on DM 288 289 Input Parameter: 290 . dm - the DM context 291 292 Output Parameter: 293 . ctype - the matrix type 294 295 Options Database: 296 . -dm_mat_type ctype 297 298 Level: intermediate 299 300 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType, DMSetMatType() 301 @*/ 302 PetscErrorCode DMGetMatType(DM dm,MatType *ctype) 303 { 304 PetscFunctionBegin; 305 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 306 *ctype = dm->mattype; 307 PetscFunctionReturn(0); 308 } 309 310 /*@ 311 MatGetDM - Gets the DM defining the data layout of the matrix 312 313 Not collective 314 315 Input Parameter: 316 . A - The Mat 317 318 Output Parameter: 319 . dm - The DM 320 321 Level: intermediate 322 323 .seealso: MatSetDM(), DMCreateMatrix(), DMSetMatType() 324 @*/ 325 PetscErrorCode MatGetDM(Mat A, DM *dm) 326 { 327 PetscErrorCode ierr; 328 329 PetscFunctionBegin; 330 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 331 PetscValidPointer(dm,2); 332 ierr = PetscObjectQuery((PetscObject) A, "__PETSc_dm", (PetscObject*) dm);CHKERRQ(ierr); 333 PetscFunctionReturn(0); 334 } 335 336 /*@ 337 MatSetDM - Sets the DM defining the data layout of the matrix 338 339 Not collective 340 341 Input Parameters: 342 + A - The Mat 343 - dm - The DM 344 345 Level: intermediate 346 347 .seealso: MatGetDM(), DMCreateMatrix(), DMSetMatType() 348 @*/ 349 PetscErrorCode MatSetDM(Mat A, DM dm) 350 { 351 PetscErrorCode ierr; 352 353 PetscFunctionBegin; 354 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 355 if (dm) PetscValidHeaderSpecific(dm,DM_CLASSID,2); 356 ierr = PetscObjectCompose((PetscObject) A, "__PETSc_dm", (PetscObject) dm);CHKERRQ(ierr); 357 PetscFunctionReturn(0); 358 } 359 360 /*@C 361 DMSetOptionsPrefix - Sets the prefix used for searching for all 362 DM options in the database. 363 364 Logically Collective on DM 365 366 Input Parameter: 367 + da - the DM context 368 - prefix - the prefix to prepend to all option names 369 370 Notes: 371 A hyphen (-) must NOT be given at the beginning of the prefix name. 372 The first character of all runtime options is AUTOMATICALLY the hyphen. 373 374 Level: advanced 375 376 .keywords: DM, set, options, prefix, database 377 378 .seealso: DMSetFromOptions() 379 @*/ 380 PetscErrorCode DMSetOptionsPrefix(DM dm,const char prefix[]) 381 { 382 PetscErrorCode ierr; 383 384 PetscFunctionBegin; 385 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 386 ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr); 387 if (dm->sf) { 388 ierr = PetscObjectSetOptionsPrefix((PetscObject)dm->sf,prefix);CHKERRQ(ierr); 389 } 390 if (dm->defaultSF) { 391 ierr = PetscObjectSetOptionsPrefix((PetscObject)dm->defaultSF,prefix);CHKERRQ(ierr); 392 } 393 PetscFunctionReturn(0); 394 } 395 396 /*@C 397 DMAppendOptionsPrefix - Appends to the prefix used for searching for all 398 DM options in the database. 399 400 Logically Collective on DM 401 402 Input Parameters: 403 + dm - the DM context 404 - prefix - the prefix string to prepend to all DM option requests 405 406 Notes: 407 A hyphen (-) must NOT be given at the beginning of the prefix name. 408 The first character of all runtime options is AUTOMATICALLY the hyphen. 409 410 Level: advanced 411 412 .keywords: DM, append, options, prefix, database 413 414 .seealso: DMSetOptionsPrefix(), DMGetOptionsPrefix() 415 @*/ 416 PetscErrorCode DMAppendOptionsPrefix(DM dm,const char prefix[]) 417 { 418 PetscErrorCode ierr; 419 420 PetscFunctionBegin; 421 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 422 ierr = PetscObjectAppendOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr); 423 PetscFunctionReturn(0); 424 } 425 426 /*@C 427 DMGetOptionsPrefix - Gets the prefix used for searching for all 428 DM options in the database. 429 430 Not Collective 431 432 Input Parameters: 433 . dm - the DM context 434 435 Output Parameters: 436 . prefix - pointer to the prefix string used is returned 437 438 Notes: On the fortran side, the user should pass in a string 'prefix' of 439 sufficient length to hold the prefix. 440 441 Level: advanced 442 443 .keywords: DM, set, options, prefix, database 444 445 .seealso: DMSetOptionsPrefix(), DMAppendOptionsPrefix() 446 @*/ 447 PetscErrorCode DMGetOptionsPrefix(DM dm,const char *prefix[]) 448 { 449 PetscErrorCode ierr; 450 451 PetscFunctionBegin; 452 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 453 ierr = PetscObjectGetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr); 454 PetscFunctionReturn(0); 455 } 456 457 static PetscErrorCode DMCountNonCyclicReferences(DM dm, PetscBool recurseCoarse, PetscBool recurseFine, PetscInt *ncrefct) 458 { 459 PetscInt i, refct = ((PetscObject) dm)->refct; 460 DMNamedVecLink nlink; 461 PetscErrorCode ierr; 462 463 PetscFunctionBegin; 464 *ncrefct = 0; 465 /* count all the circular references of DM and its contained Vecs */ 466 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 467 if (dm->localin[i]) refct--; 468 if (dm->globalin[i]) refct--; 469 } 470 for (nlink=dm->namedglobal; nlink; nlink=nlink->next) refct--; 471 for (nlink=dm->namedlocal; nlink; nlink=nlink->next) refct--; 472 if (dm->x) { 473 DM obj; 474 ierr = VecGetDM(dm->x, &obj);CHKERRQ(ierr); 475 if (obj == dm) refct--; 476 } 477 if (dm->coarseMesh && dm->coarseMesh->fineMesh == dm) { 478 refct--; 479 if (recurseCoarse) { 480 PetscInt coarseCount; 481 482 ierr = DMCountNonCyclicReferences(dm->coarseMesh, PETSC_TRUE, PETSC_FALSE,&coarseCount);CHKERRQ(ierr); 483 refct += coarseCount; 484 } 485 } 486 if (dm->fineMesh && dm->fineMesh->coarseMesh == dm) { 487 refct--; 488 if (recurseFine) { 489 PetscInt fineCount; 490 491 ierr = DMCountNonCyclicReferences(dm->fineMesh, PETSC_FALSE, PETSC_TRUE,&fineCount);CHKERRQ(ierr); 492 refct += fineCount; 493 } 494 } 495 *ncrefct = refct; 496 PetscFunctionReturn(0); 497 } 498 499 PetscErrorCode DMDestroyLabelLinkList(DM dm) 500 { 501 PetscErrorCode ierr; 502 503 PetscFunctionBegin; 504 if (!--(dm->labels->refct)) { 505 DMLabelLink next = dm->labels->next; 506 507 /* destroy the labels */ 508 while (next) { 509 DMLabelLink tmp = next->next; 510 511 ierr = DMLabelDestroy(&next->label);CHKERRQ(ierr); 512 ierr = PetscFree(next);CHKERRQ(ierr); 513 next = tmp; 514 } 515 ierr = PetscFree(dm->labels);CHKERRQ(ierr); 516 } 517 PetscFunctionReturn(0); 518 } 519 520 /*@ 521 DMDestroy - Destroys a vector packer or DM. 522 523 Collective on DM 524 525 Input Parameter: 526 . dm - the DM object to destroy 527 528 Level: developer 529 530 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 531 532 @*/ 533 PetscErrorCode DMDestroy(DM *dm) 534 { 535 PetscInt i, cnt; 536 DMNamedVecLink nlink,nnext; 537 PetscErrorCode ierr; 538 539 PetscFunctionBegin; 540 if (!*dm) PetscFunctionReturn(0); 541 PetscValidHeaderSpecific((*dm),DM_CLASSID,1); 542 543 /* count all non-cyclic references in the doubly-linked list of coarse<->fine meshes */ 544 ierr = DMCountNonCyclicReferences(*dm,PETSC_TRUE,PETSC_TRUE,&cnt);CHKERRQ(ierr); 545 --((PetscObject)(*dm))->refct; 546 if (--cnt > 0) {*dm = 0; PetscFunctionReturn(0);} 547 /* 548 Need this test because the dm references the vectors that 549 reference the dm, so destroying the dm calls destroy on the 550 vectors that cause another destroy on the dm 551 */ 552 if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0); 553 ((PetscObject) (*dm))->refct = 0; 554 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 555 if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()"); 556 ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr); 557 } 558 nnext=(*dm)->namedglobal; 559 (*dm)->namedglobal = NULL; 560 for (nlink=nnext; nlink; nlink=nnext) { /* Destroy the named vectors */ 561 nnext = nlink->next; 562 if (nlink->status != DMVEC_STATUS_IN) SETERRQ1(((PetscObject)*dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM still has Vec named '%s' checked out",nlink->name); 563 ierr = PetscFree(nlink->name);CHKERRQ(ierr); 564 ierr = VecDestroy(&nlink->X);CHKERRQ(ierr); 565 ierr = PetscFree(nlink);CHKERRQ(ierr); 566 } 567 nnext=(*dm)->namedlocal; 568 (*dm)->namedlocal = NULL; 569 for (nlink=nnext; nlink; nlink=nnext) { /* Destroy the named local vectors */ 570 nnext = nlink->next; 571 if (nlink->status != DMVEC_STATUS_IN) SETERRQ1(((PetscObject)*dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM still has Vec named '%s' checked out",nlink->name); 572 ierr = PetscFree(nlink->name);CHKERRQ(ierr); 573 ierr = VecDestroy(&nlink->X);CHKERRQ(ierr); 574 ierr = PetscFree(nlink);CHKERRQ(ierr); 575 } 576 577 /* Destroy the list of hooks */ 578 { 579 DMCoarsenHookLink link,next; 580 for (link=(*dm)->coarsenhook; link; link=next) { 581 next = link->next; 582 ierr = PetscFree(link);CHKERRQ(ierr); 583 } 584 (*dm)->coarsenhook = NULL; 585 } 586 { 587 DMRefineHookLink link,next; 588 for (link=(*dm)->refinehook; link; link=next) { 589 next = link->next; 590 ierr = PetscFree(link);CHKERRQ(ierr); 591 } 592 (*dm)->refinehook = NULL; 593 } 594 { 595 DMSubDomainHookLink link,next; 596 for (link=(*dm)->subdomainhook; link; link=next) { 597 next = link->next; 598 ierr = PetscFree(link);CHKERRQ(ierr); 599 } 600 (*dm)->subdomainhook = NULL; 601 } 602 { 603 DMGlobalToLocalHookLink link,next; 604 for (link=(*dm)->gtolhook; link; link=next) { 605 next = link->next; 606 ierr = PetscFree(link);CHKERRQ(ierr); 607 } 608 (*dm)->gtolhook = NULL; 609 } 610 { 611 DMLocalToGlobalHookLink link,next; 612 for (link=(*dm)->ltoghook; link; link=next) { 613 next = link->next; 614 ierr = PetscFree(link);CHKERRQ(ierr); 615 } 616 (*dm)->ltoghook = NULL; 617 } 618 /* Destroy the work arrays */ 619 { 620 DMWorkLink link,next; 621 if ((*dm)->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Work array still checked out"); 622 for (link=(*dm)->workin; link; link=next) { 623 next = link->next; 624 ierr = PetscFree(link->mem);CHKERRQ(ierr); 625 ierr = PetscFree(link);CHKERRQ(ierr); 626 } 627 (*dm)->workin = NULL; 628 } 629 if (!--((*dm)->labels->refct)) { 630 DMLabelLink next = (*dm)->labels->next; 631 632 /* destroy the labels */ 633 while (next) { 634 DMLabelLink tmp = next->next; 635 636 ierr = DMLabelDestroy(&next->label);CHKERRQ(ierr); 637 ierr = PetscFree(next);CHKERRQ(ierr); 638 next = tmp; 639 } 640 ierr = PetscFree((*dm)->labels);CHKERRQ(ierr); 641 } 642 { 643 DMBoundary next = (*dm)->boundary; 644 while (next) { 645 DMBoundary b = next; 646 647 next = b->next; 648 ierr = PetscFree(b);CHKERRQ(ierr); 649 } 650 } 651 652 ierr = PetscObjectDestroy(&(*dm)->dmksp);CHKERRQ(ierr); 653 ierr = PetscObjectDestroy(&(*dm)->dmsnes);CHKERRQ(ierr); 654 ierr = PetscObjectDestroy(&(*dm)->dmts);CHKERRQ(ierr); 655 656 if ((*dm)->ctx && (*dm)->ctxdestroy) { 657 ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr); 658 } 659 ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr); 660 ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr); 661 ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr); 662 ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr); 663 ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr); 664 ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr); 665 666 ierr = PetscSectionDestroy(&(*dm)->defaultSection);CHKERRQ(ierr); 667 ierr = PetscSectionDestroy(&(*dm)->defaultGlobalSection);CHKERRQ(ierr); 668 ierr = PetscLayoutDestroy(&(*dm)->map);CHKERRQ(ierr); 669 ierr = PetscSectionDestroy(&(*dm)->defaultConstraintSection);CHKERRQ(ierr); 670 ierr = MatDestroy(&(*dm)->defaultConstraintMat);CHKERRQ(ierr); 671 ierr = PetscSFDestroy(&(*dm)->sf);CHKERRQ(ierr); 672 ierr = PetscSFDestroy(&(*dm)->defaultSF);CHKERRQ(ierr); 673 ierr = PetscSFDestroy(&(*dm)->sfNatural);CHKERRQ(ierr); 674 675 if ((*dm)->coarseMesh && (*dm)->coarseMesh->fineMesh == *dm) { 676 ierr = DMSetFineDM((*dm)->coarseMesh,NULL);CHKERRQ(ierr); 677 } 678 ierr = DMDestroy(&(*dm)->coarseMesh);CHKERRQ(ierr); 679 if ((*dm)->fineMesh && (*dm)->fineMesh->coarseMesh == *dm) { 680 ierr = DMSetCoarseDM((*dm)->fineMesh,NULL);CHKERRQ(ierr); 681 } 682 ierr = DMDestroy(&(*dm)->fineMesh);CHKERRQ(ierr); 683 ierr = DMDestroy(&(*dm)->coordinateDM);CHKERRQ(ierr); 684 ierr = VecDestroy(&(*dm)->coordinates);CHKERRQ(ierr); 685 ierr = VecDestroy(&(*dm)->coordinatesLocal);CHKERRQ(ierr); 686 ierr = PetscFree3((*dm)->L,(*dm)->maxCell,(*dm)->bdtype);CHKERRQ(ierr); 687 688 ierr = PetscDSDestroy(&(*dm)->prob);CHKERRQ(ierr); 689 ierr = DMDestroy(&(*dm)->dmBC);CHKERRQ(ierr); 690 /* if memory was published with SAWs then destroy it */ 691 ierr = PetscObjectSAWsViewOff((PetscObject)*dm);CHKERRQ(ierr); 692 693 ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr); 694 /* We do not destroy (*dm)->data here so that we can reference count backend objects */ 695 ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr); 696 PetscFunctionReturn(0); 697 } 698 699 /*@ 700 DMSetUp - sets up the data structures inside a DM object 701 702 Collective on DM 703 704 Input Parameter: 705 . dm - the DM object to setup 706 707 Level: developer 708 709 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 710 711 @*/ 712 PetscErrorCode DMSetUp(DM dm) 713 { 714 PetscErrorCode ierr; 715 716 PetscFunctionBegin; 717 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 718 if (dm->setupcalled) PetscFunctionReturn(0); 719 if (dm->ops->setup) { 720 ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr); 721 } 722 dm->setupcalled = PETSC_TRUE; 723 PetscFunctionReturn(0); 724 } 725 726 /*@ 727 DMSetFromOptions - sets parameters in a DM from the options database 728 729 Collective on DM 730 731 Input Parameter: 732 . dm - the DM object to set options for 733 734 Options Database: 735 + -dm_preallocate_only - Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros 736 . -dm_vec_type <type> - type of vector to create inside DM 737 . -dm_mat_type <type> - type of matrix to create inside DM 738 - -dm_coloring_type - <global or ghosted> 739 740 Level: developer 741 742 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 743 744 @*/ 745 PetscErrorCode DMSetFromOptions(DM dm) 746 { 747 char typeName[256]; 748 PetscBool flg; 749 PetscErrorCode ierr; 750 751 PetscFunctionBegin; 752 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 753 if (dm->prob) { 754 ierr = PetscDSSetFromOptions(dm->prob);CHKERRQ(ierr); 755 } 756 if (dm->sf) { 757 ierr = PetscSFSetFromOptions(dm->sf);CHKERRQ(ierr); 758 } 759 if (dm->defaultSF) { 760 ierr = PetscSFSetFromOptions(dm->defaultSF);CHKERRQ(ierr); 761 } 762 ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr); 763 ierr = PetscOptionsBool("-dm_preallocate_only","only preallocate matrix, but do not set column indices","DMSetMatrixPreallocateOnly",dm->prealloc_only,&dm->prealloc_only,NULL);CHKERRQ(ierr); 764 ierr = PetscOptionsFList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr); 765 if (flg) { 766 ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr); 767 } 768 ierr = PetscOptionsFList("-dm_mat_type","Matrix type used for created matrices","DMSetMatType",MatList,dm->mattype ? dm->mattype : typeName,typeName,sizeof(typeName),&flg);CHKERRQ(ierr); 769 if (flg) { 770 ierr = DMSetMatType(dm,typeName);CHKERRQ(ierr); 771 } 772 ierr = PetscOptionsEnum("-dm_is_coloring_type","Global or local coloring of Jacobian","ISColoringType",ISColoringTypes,(PetscEnum)dm->coloringtype,(PetscEnum*)&dm->coloringtype,NULL);CHKERRQ(ierr); 773 if (dm->ops->setfromoptions) { 774 ierr = (*dm->ops->setfromoptions)(PetscOptionsObject,dm);CHKERRQ(ierr); 775 } 776 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 777 ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) dm);CHKERRQ(ierr); 778 ierr = PetscOptionsEnd();CHKERRQ(ierr); 779 PetscFunctionReturn(0); 780 } 781 782 /*@C 783 DMView - Views a DM 784 785 Collective on DM 786 787 Input Parameter: 788 + dm - the DM object to view 789 - v - the viewer 790 791 Level: beginner 792 793 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 794 795 @*/ 796 PetscErrorCode DMView(DM dm,PetscViewer v) 797 { 798 PetscErrorCode ierr; 799 PetscBool isbinary; 800 801 PetscFunctionBegin; 802 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 803 if (!v) { 804 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)dm),&v);CHKERRQ(ierr); 805 } 806 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)dm,v);CHKERRQ(ierr); 807 ierr = PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 808 if (isbinary) { 809 PetscInt classid = DM_FILE_CLASSID; 810 char type[256]; 811 812 ierr = PetscViewerBinaryWrite(v,&classid,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 813 ierr = PetscStrncpy(type,((PetscObject)dm)->type_name,256);CHKERRQ(ierr); 814 ierr = PetscViewerBinaryWrite(v,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 815 } 816 if (dm->ops->view) { 817 ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr); 818 } 819 PetscFunctionReturn(0); 820 } 821 822 /*@ 823 DMCreateGlobalVector - Creates a global vector from a DM object 824 825 Collective on DM 826 827 Input Parameter: 828 . dm - the DM object 829 830 Output Parameter: 831 . vec - the global vector 832 833 Level: beginner 834 835 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 836 837 @*/ 838 PetscErrorCode DMCreateGlobalVector(DM dm,Vec *vec) 839 { 840 PetscErrorCode ierr; 841 842 PetscFunctionBegin; 843 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 844 ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr); 845 PetscFunctionReturn(0); 846 } 847 848 /*@ 849 DMCreateLocalVector - Creates a local vector from a DM object 850 851 Not Collective 852 853 Input Parameter: 854 . dm - the DM object 855 856 Output Parameter: 857 . vec - the local vector 858 859 Level: beginner 860 861 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 862 863 @*/ 864 PetscErrorCode DMCreateLocalVector(DM dm,Vec *vec) 865 { 866 PetscErrorCode ierr; 867 868 PetscFunctionBegin; 869 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 870 ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr); 871 PetscFunctionReturn(0); 872 } 873 874 /*@ 875 DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM. 876 877 Collective on DM 878 879 Input Parameter: 880 . dm - the DM that provides the mapping 881 882 Output Parameter: 883 . ltog - the mapping 884 885 Level: intermediate 886 887 Notes: 888 This mapping can then be used by VecSetLocalToGlobalMapping() or 889 MatSetLocalToGlobalMapping(). 890 891 .seealso: DMCreateLocalVector() 892 @*/ 893 PetscErrorCode DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog) 894 { 895 PetscInt bs = -1, bsLocal[2], bsMinMax[2]; 896 PetscErrorCode ierr; 897 898 PetscFunctionBegin; 899 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 900 PetscValidPointer(ltog,2); 901 if (!dm->ltogmap) { 902 PetscSection section, sectionGlobal; 903 904 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 905 if (section) { 906 const PetscInt *cdofs; 907 PetscInt *ltog; 908 PetscInt pStart, pEnd, n, p, k, l; 909 910 ierr = DMGetDefaultGlobalSection(dm, §ionGlobal);CHKERRQ(ierr); 911 ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 912 ierr = PetscSectionGetStorageSize(section, &n);CHKERRQ(ierr); 913 ierr = PetscMalloc1(n, <og);CHKERRQ(ierr); /* We want the local+overlap size */ 914 for (p = pStart, l = 0; p < pEnd; ++p) { 915 PetscInt bdof, cdof, dof, off, c, cind = 0; 916 917 /* Should probably use constrained dofs */ 918 ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 919 ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 920 ierr = PetscSectionGetConstraintIndices(section, p, &cdofs);CHKERRQ(ierr); 921 ierr = PetscSectionGetOffset(sectionGlobal, p, &off);CHKERRQ(ierr); 922 /* If you have dofs, and constraints, and they are unequal, we set the blocksize to 1 */ 923 bdof = cdof && (dof-cdof) ? 1 : dof; 924 if (dof) { 925 if (bs < 0) {bs = bdof;} 926 else if (bs != bdof) {bs = 1;} 927 } 928 for (c = 0; c < dof; ++c, ++l) { 929 if ((cind < cdof) && (c == cdofs[cind])) ltog[l] = off < 0 ? off-c : off+c; 930 else ltog[l] = (off < 0 ? -(off+1) : off) + c; 931 } 932 } 933 /* Must have same blocksize on all procs (some might have no points) */ 934 bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs; bsLocal[1] = bs; 935 ierr = PetscGlobalMinMaxInt(PetscObjectComm((PetscObject) dm), bsLocal, bsMinMax);CHKERRQ(ierr); 936 if (bsMinMax[0] != bsMinMax[1]) {bs = 1;} 937 else {bs = bsMinMax[0];} 938 bs = bs < 0 ? 1 : bs; 939 /* Must reduce indices by blocksize */ 940 if (bs > 1) { 941 for (l = 0, k = 0; l < n; l += bs, ++k) ltog[k] = ltog[l]/bs; 942 n /= bs; 943 } 944 ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)dm), bs, n, ltog, PETSC_OWN_POINTER, &dm->ltogmap);CHKERRQ(ierr); 945 ierr = PetscLogObjectParent((PetscObject)dm, (PetscObject)dm->ltogmap);CHKERRQ(ierr); 946 } else { 947 if (!dm->ops->getlocaltoglobalmapping) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping"); 948 ierr = (*dm->ops->getlocaltoglobalmapping)(dm);CHKERRQ(ierr); 949 } 950 } 951 *ltog = dm->ltogmap; 952 PetscFunctionReturn(0); 953 } 954 955 /*@ 956 DMGetBlockSize - Gets the inherent block size associated with a DM 957 958 Not Collective 959 960 Input Parameter: 961 . dm - the DM with block structure 962 963 Output Parameter: 964 . bs - the block size, 1 implies no exploitable block structure 965 966 Level: intermediate 967 968 .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMapping() 969 @*/ 970 PetscErrorCode DMGetBlockSize(DM dm,PetscInt *bs) 971 { 972 PetscFunctionBegin; 973 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 974 PetscValidPointer(bs,2); 975 if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet"); 976 *bs = dm->bs; 977 PetscFunctionReturn(0); 978 } 979 980 /*@ 981 DMCreateInterpolation - Gets interpolation matrix between two DM objects 982 983 Collective on DM 984 985 Input Parameter: 986 + dm1 - the DM object 987 - dm2 - the second, finer DM object 988 989 Output Parameter: 990 + mat - the interpolation 991 - vec - the scaling (optional) 992 993 Level: developer 994 995 Notes: For DMDA objects this only works for "uniform refinement", that is the refined mesh was obtained DMRefine() or the coarse mesh was obtained by 996 DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation. 997 998 For DMDA objects you can use this interpolation (more precisely the interpolation from the DMGetCoordinateDM()) to interpolate the mesh coordinate vectors 999 EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic. 1000 1001 1002 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen(), DMCreateRestriction() 1003 1004 @*/ 1005 PetscErrorCode DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec) 1006 { 1007 PetscErrorCode ierr; 1008 1009 PetscFunctionBegin; 1010 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 1011 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 1012 ierr = PetscLogEventBegin(DM_CreateInterpolation,dm1,dm2,0,0);CHKERRQ(ierr); 1013 ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr); 1014 ierr = PetscLogEventEnd(DM_CreateInterpolation,dm1,dm2,0,0);CHKERRQ(ierr); 1015 PetscFunctionReturn(0); 1016 } 1017 1018 /*@ 1019 DMCreateRestriction - Gets restriction matrix between two DM objects 1020 1021 Collective on DM 1022 1023 Input Parameter: 1024 + dm1 - the DM object 1025 - dm2 - the second, finer DM object 1026 1027 Output Parameter: 1028 . mat - the restriction 1029 1030 1031 Level: developer 1032 1033 Notes: For DMDA objects this only works for "uniform refinement", that is the refined mesh was obtained DMRefine() or the coarse mesh was obtained by 1034 DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation. 1035 1036 1037 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen(), DMCreateInterpolation() 1038 1039 @*/ 1040 PetscErrorCode DMCreateRestriction(DM dm1,DM dm2,Mat *mat) 1041 { 1042 PetscErrorCode ierr; 1043 1044 PetscFunctionBegin; 1045 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 1046 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 1047 ierr = PetscLogEventBegin(DM_CreateRestriction,dm1,dm2,0,0);CHKERRQ(ierr); 1048 if (!dm1->ops->createrestriction) SETERRQ(PetscObjectComm((PetscObject)dm1),PETSC_ERR_SUP,"DMCreateRestriction not implemented for this type"); 1049 ierr = (*dm1->ops->createrestriction)(dm1,dm2,mat);CHKERRQ(ierr); 1050 ierr = PetscLogEventEnd(DM_CreateRestriction,dm1,dm2,0,0);CHKERRQ(ierr); 1051 PetscFunctionReturn(0); 1052 } 1053 1054 /*@ 1055 DMCreateInjection - Gets injection matrix between two DM objects 1056 1057 Collective on DM 1058 1059 Input Parameter: 1060 + dm1 - the DM object 1061 - dm2 - the second, finer DM object 1062 1063 Output Parameter: 1064 . mat - the injection 1065 1066 Level: developer 1067 1068 Notes: For DMDA objects this only works for "uniform refinement", that is the refined mesh was obtained DMRefine() or the coarse mesh was obtained by 1069 DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection. 1070 1071 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation() 1072 1073 @*/ 1074 PetscErrorCode DMCreateInjection(DM dm1,DM dm2,Mat *mat) 1075 { 1076 PetscErrorCode ierr; 1077 1078 PetscFunctionBegin; 1079 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 1080 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 1081 if (!dm1->ops->getinjection) SETERRQ(PetscObjectComm((PetscObject)dm1),PETSC_ERR_SUP,"DMCreateInjection not implemented for this type"); 1082 ierr = (*dm1->ops->getinjection)(dm1,dm2,mat);CHKERRQ(ierr); 1083 PetscFunctionReturn(0); 1084 } 1085 1086 /*@ 1087 DMCreateMassMatrix - Gets mass matrix between two DM objects, M_ij = \int \phi_i \psi_j 1088 1089 Collective on DM 1090 1091 Input Parameter: 1092 + dm1 - the DM object 1093 - dm2 - the second, finer DM object 1094 1095 Output Parameter: 1096 . mat - the interpolation 1097 1098 Level: developer 1099 1100 .seealso DMCreateMatrix(), DMRefine(), DMCoarsen(), DMCreateRestriction(), DMCreateInterpolation(), DMCreateInjection() 1101 @*/ 1102 PetscErrorCode DMCreateMassMatrix(DM dm1, DM dm2, Mat *mat) 1103 { 1104 PetscErrorCode ierr; 1105 1106 PetscFunctionBegin; 1107 PetscValidHeaderSpecific(dm1, DM_CLASSID, 1); 1108 PetscValidHeaderSpecific(dm2, DM_CLASSID, 2); 1109 ierr = (*dm1->ops->createmassmatrix)(dm1, dm2, mat);CHKERRQ(ierr); 1110 PetscFunctionReturn(0); 1111 } 1112 1113 /*@ 1114 DMCreateColoring - Gets coloring for a DM 1115 1116 Collective on DM 1117 1118 Input Parameter: 1119 + dm - the DM object 1120 - ctype - IS_COLORING_LOCAL or IS_COLORING_GLOBAL 1121 1122 Output Parameter: 1123 . coloring - the coloring 1124 1125 Level: developer 1126 1127 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix(), DMSetMatType() 1128 1129 @*/ 1130 PetscErrorCode DMCreateColoring(DM dm,ISColoringType ctype,ISColoring *coloring) 1131 { 1132 PetscErrorCode ierr; 1133 1134 PetscFunctionBegin; 1135 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1136 if (!dm->ops->getcoloring) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No coloring for this type of DM yet"); 1137 ierr = (*dm->ops->getcoloring)(dm,ctype,coloring);CHKERRQ(ierr); 1138 PetscFunctionReturn(0); 1139 } 1140 1141 /*@ 1142 DMCreateMatrix - Gets empty Jacobian for a DM 1143 1144 Collective on DM 1145 1146 Input Parameter: 1147 . dm - the DM object 1148 1149 Output Parameter: 1150 . mat - the empty Jacobian 1151 1152 Level: beginner 1153 1154 Notes: This properly preallocates the number of nonzeros in the sparse matrix so you 1155 do not need to do it yourself. 1156 1157 By default it also sets the nonzero structure and puts in the zero entries. To prevent setting 1158 the nonzero pattern call DMSetMatrixPreallocateOnly() 1159 1160 For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used 1161 internally by PETSc. 1162 1163 For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires 1164 the indices for the global numbering for DMDAs which is complicated. 1165 1166 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMSetMatType() 1167 1168 @*/ 1169 PetscErrorCode DMCreateMatrix(DM dm,Mat *mat) 1170 { 1171 PetscErrorCode ierr; 1172 1173 PetscFunctionBegin; 1174 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1175 ierr = MatInitializePackage();CHKERRQ(ierr); 1176 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1177 PetscValidPointer(mat,3); 1178 ierr = (*dm->ops->creatematrix)(dm,mat);CHKERRQ(ierr); 1179 PetscFunctionReturn(0); 1180 } 1181 1182 /*@ 1183 DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly 1184 preallocated but the nonzero structure and zero values will not be set. 1185 1186 Logically Collective on DM 1187 1188 Input Parameter: 1189 + dm - the DM 1190 - only - PETSC_TRUE if only want preallocation 1191 1192 Level: developer 1193 .seealso DMCreateMatrix(), DMSetMatrixStructureOnly() 1194 @*/ 1195 PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only) 1196 { 1197 PetscFunctionBegin; 1198 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1199 dm->prealloc_only = only; 1200 PetscFunctionReturn(0); 1201 } 1202 1203 /*@ 1204 DMSetMatrixStructureOnly - When DMCreateMatrix() is called, the matrix structure will be created 1205 but the array for values will not be allocated. 1206 1207 Logically Collective on DM 1208 1209 Input Parameter: 1210 + dm - the DM 1211 - only - PETSC_TRUE if only want matrix stucture 1212 1213 Level: developer 1214 .seealso DMCreateMatrix(), DMSetMatrixPreallocateOnly() 1215 @*/ 1216 PetscErrorCode DMSetMatrixStructureOnly(DM dm, PetscBool only) 1217 { 1218 PetscFunctionBegin; 1219 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1220 dm->structure_only = only; 1221 PetscFunctionReturn(0); 1222 } 1223 1224 /*@C 1225 DMGetWorkArray - Gets a work array guaranteed to be at least the input size, restore with DMRestoreWorkArray() 1226 1227 Not Collective 1228 1229 Input Parameters: 1230 + dm - the DM object 1231 . count - The minium size 1232 - dtype - data type (PETSC_REAL, PETSC_SCALAR, PETSC_INT) 1233 1234 Output Parameter: 1235 . array - the work array 1236 1237 Level: developer 1238 1239 .seealso DMDestroy(), DMCreate() 1240 @*/ 1241 PetscErrorCode DMGetWorkArray(DM dm,PetscInt count,PetscDataType dtype,void *mem) 1242 { 1243 PetscErrorCode ierr; 1244 DMWorkLink link; 1245 size_t dsize; 1246 1247 PetscFunctionBegin; 1248 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1249 PetscValidPointer(mem,4); 1250 if (dm->workin) { 1251 link = dm->workin; 1252 dm->workin = dm->workin->next; 1253 } else { 1254 ierr = PetscNewLog(dm,&link);CHKERRQ(ierr); 1255 } 1256 ierr = PetscDataTypeGetSize(dtype,&dsize);CHKERRQ(ierr); 1257 if (dsize*count > link->bytes) { 1258 ierr = PetscFree(link->mem);CHKERRQ(ierr); 1259 ierr = PetscMalloc(dsize*count,&link->mem);CHKERRQ(ierr); 1260 link->bytes = dsize*count; 1261 } 1262 link->next = dm->workout; 1263 dm->workout = link; 1264 *(void**)mem = link->mem; 1265 PetscFunctionReturn(0); 1266 } 1267 1268 /*@C 1269 DMRestoreWorkArray - Restores a work array guaranteed to be at least the input size, restore with DMRestoreWorkArray() 1270 1271 Not Collective 1272 1273 Input Parameters: 1274 + dm - the DM object 1275 . count - The minium size 1276 - dtype - data type (PETSC_REAL, PETSC_SCALAR, PETSC_INT) 1277 1278 Output Parameter: 1279 . array - the work array 1280 1281 Level: developer 1282 1283 .seealso DMDestroy(), DMCreate() 1284 @*/ 1285 PetscErrorCode DMRestoreWorkArray(DM dm,PetscInt count,PetscDataType dtype,void *mem) 1286 { 1287 DMWorkLink *p,link; 1288 1289 PetscFunctionBegin; 1290 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1291 PetscValidPointer(mem,4); 1292 for (p=&dm->workout; (link=*p); p=&link->next) { 1293 if (link->mem == *(void**)mem) { 1294 *p = link->next; 1295 link->next = dm->workin; 1296 dm->workin = link; 1297 *(void**)mem = NULL; 1298 PetscFunctionReturn(0); 1299 } 1300 } 1301 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Array was not checked out"); 1302 } 1303 1304 PetscErrorCode DMSetNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (*nullsp)(DM dm, PetscInt field, MatNullSpace *nullSpace)) 1305 { 1306 PetscFunctionBegin; 1307 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1308 if (field >= 10) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %d >= 10 fields", field); 1309 dm->nullspaceConstructors[field] = nullsp; 1310 PetscFunctionReturn(0); 1311 } 1312 1313 /*@C 1314 DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field 1315 1316 Not collective 1317 1318 Input Parameter: 1319 . dm - the DM object 1320 1321 Output Parameters: 1322 + numFields - The number of fields (or NULL if not requested) 1323 . fieldNames - The name for each field (or NULL if not requested) 1324 - fields - The global indices for each field (or NULL if not requested) 1325 1326 Level: intermediate 1327 1328 Notes: 1329 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 1330 PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with 1331 PetscFree(). 1332 1333 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 1334 @*/ 1335 PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields) 1336 { 1337 PetscSection section, sectionGlobal; 1338 PetscErrorCode ierr; 1339 1340 PetscFunctionBegin; 1341 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1342 if (numFields) { 1343 PetscValidPointer(numFields,2); 1344 *numFields = 0; 1345 } 1346 if (fieldNames) { 1347 PetscValidPointer(fieldNames,3); 1348 *fieldNames = NULL; 1349 } 1350 if (fields) { 1351 PetscValidPointer(fields,4); 1352 *fields = NULL; 1353 } 1354 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1355 if (section) { 1356 PetscInt *fieldSizes, **fieldIndices; 1357 PetscInt nF, f, pStart, pEnd, p; 1358 1359 ierr = DMGetDefaultGlobalSection(dm, §ionGlobal);CHKERRQ(ierr); 1360 ierr = PetscSectionGetNumFields(section, &nF);CHKERRQ(ierr); 1361 ierr = PetscMalloc2(nF,&fieldSizes,nF,&fieldIndices);CHKERRQ(ierr); 1362 ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr); 1363 for (f = 0; f < nF; ++f) { 1364 fieldSizes[f] = 0; 1365 } 1366 for (p = pStart; p < pEnd; ++p) { 1367 PetscInt gdof; 1368 1369 ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr); 1370 if (gdof > 0) { 1371 for (f = 0; f < nF; ++f) { 1372 PetscInt fdof, fcdof; 1373 1374 ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr); 1375 ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr); 1376 fieldSizes[f] += fdof-fcdof; 1377 } 1378 } 1379 } 1380 for (f = 0; f < nF; ++f) { 1381 ierr = PetscMalloc1(fieldSizes[f], &fieldIndices[f]);CHKERRQ(ierr); 1382 fieldSizes[f] = 0; 1383 } 1384 for (p = pStart; p < pEnd; ++p) { 1385 PetscInt gdof, goff; 1386 1387 ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr); 1388 if (gdof > 0) { 1389 ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 1390 for (f = 0; f < nF; ++f) { 1391 PetscInt fdof, fcdof, fc; 1392 1393 ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr); 1394 ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr); 1395 for (fc = 0; fc < fdof-fcdof; ++fc, ++fieldSizes[f]) { 1396 fieldIndices[f][fieldSizes[f]] = goff++; 1397 } 1398 } 1399 } 1400 } 1401 if (numFields) *numFields = nF; 1402 if (fieldNames) { 1403 ierr = PetscMalloc1(nF, fieldNames);CHKERRQ(ierr); 1404 for (f = 0; f < nF; ++f) { 1405 const char *fieldName; 1406 1407 ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr); 1408 ierr = PetscStrallocpy(fieldName, (char**) &(*fieldNames)[f]);CHKERRQ(ierr); 1409 } 1410 } 1411 if (fields) { 1412 ierr = PetscMalloc1(nF, fields);CHKERRQ(ierr); 1413 for (f = 0; f < nF; ++f) { 1414 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), fieldSizes[f], fieldIndices[f], PETSC_OWN_POINTER, &(*fields)[f]);CHKERRQ(ierr); 1415 } 1416 } 1417 ierr = PetscFree2(fieldSizes,fieldIndices);CHKERRQ(ierr); 1418 } else if (dm->ops->createfieldis) { 1419 ierr = (*dm->ops->createfieldis)(dm, numFields, fieldNames, fields);CHKERRQ(ierr); 1420 } 1421 PetscFunctionReturn(0); 1422 } 1423 1424 1425 /*@C 1426 DMCreateFieldDecomposition - Returns a list of IS objects defining a decomposition of a problem into subproblems 1427 corresponding to different fields: each IS contains the global indices of the dofs of the 1428 corresponding field. The optional list of DMs define the DM for each subproblem. 1429 Generalizes DMCreateFieldIS(). 1430 1431 Not collective 1432 1433 Input Parameter: 1434 . dm - the DM object 1435 1436 Output Parameters: 1437 + len - The number of subproblems in the field decomposition (or NULL if not requested) 1438 . namelist - The name for each field (or NULL if not requested) 1439 . islist - The global indices for each field (or NULL if not requested) 1440 - dmlist - The DMs for each field subproblem (or NULL, if not requested; if NULL is returned, no DMs are defined) 1441 1442 Level: intermediate 1443 1444 Notes: 1445 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 1446 PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(), 1447 and all of the arrays should be freed with PetscFree(). 1448 1449 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS() 1450 @*/ 1451 PetscErrorCode DMCreateFieldDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 1452 { 1453 PetscErrorCode ierr; 1454 1455 PetscFunctionBegin; 1456 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1457 if (len) { 1458 PetscValidPointer(len,2); 1459 *len = 0; 1460 } 1461 if (namelist) { 1462 PetscValidPointer(namelist,3); 1463 *namelist = 0; 1464 } 1465 if (islist) { 1466 PetscValidPointer(islist,4); 1467 *islist = 0; 1468 } 1469 if (dmlist) { 1470 PetscValidPointer(dmlist,5); 1471 *dmlist = 0; 1472 } 1473 /* 1474 Is it a good idea to apply the following check across all impls? 1475 Perhaps some impls can have a well-defined decomposition before DMSetUp? 1476 This, however, follows the general principle that accessors are not well-behaved until the object is set up. 1477 */ 1478 if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE, "Decomposition defined only after DMSetUp"); 1479 if (!dm->ops->createfielddecomposition) { 1480 PetscSection section; 1481 PetscInt numFields, f; 1482 1483 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1484 if (section) {ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);} 1485 if (section && numFields && dm->ops->createsubdm) { 1486 if (len) *len = numFields; 1487 if (namelist) {ierr = PetscMalloc1(numFields,namelist);CHKERRQ(ierr);} 1488 if (islist) {ierr = PetscMalloc1(numFields,islist);CHKERRQ(ierr);} 1489 if (dmlist) {ierr = PetscMalloc1(numFields,dmlist);CHKERRQ(ierr);} 1490 for (f = 0; f < numFields; ++f) { 1491 const char *fieldName; 1492 1493 ierr = DMCreateSubDM(dm, 1, &f, islist ? &(*islist)[f] : NULL, dmlist ? &(*dmlist)[f] : NULL);CHKERRQ(ierr); 1494 if (namelist) { 1495 ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr); 1496 ierr = PetscStrallocpy(fieldName, (char**) &(*namelist)[f]);CHKERRQ(ierr); 1497 } 1498 } 1499 } else { 1500 ierr = DMCreateFieldIS(dm, len, namelist, islist);CHKERRQ(ierr); 1501 /* By default there are no DMs associated with subproblems. */ 1502 if (dmlist) *dmlist = NULL; 1503 } 1504 } else { 1505 ierr = (*dm->ops->createfielddecomposition)(dm,len,namelist,islist,dmlist);CHKERRQ(ierr); 1506 } 1507 PetscFunctionReturn(0); 1508 } 1509 1510 /*@ 1511 DMCreateSubDM - Returns an IS and DM encapsulating a subproblem defined by the fields passed in. 1512 The fields are defined by DMCreateFieldIS(). 1513 1514 Not collective 1515 1516 Input Parameters: 1517 + dm - The DM object 1518 . numFields - The number of fields in this subproblem 1519 - fields - The field numbers of the selected fields 1520 1521 Output Parameters: 1522 + is - The global indices for the subproblem 1523 - subdm - The DM for the subproblem 1524 1525 Level: intermediate 1526 1527 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS() 1528 @*/ 1529 PetscErrorCode DMCreateSubDM(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm) 1530 { 1531 PetscErrorCode ierr; 1532 1533 PetscFunctionBegin; 1534 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1535 PetscValidPointer(fields,3); 1536 if (is) PetscValidPointer(is,4); 1537 if (subdm) PetscValidPointer(subdm,5); 1538 if (dm->ops->createsubdm) { 1539 ierr = (*dm->ops->createsubdm)(dm, numFields, fields, is, subdm);CHKERRQ(ierr); 1540 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This type has no DMCreateSubDM implementation defined"); 1541 PetscFunctionReturn(0); 1542 } 1543 1544 /*@C 1545 DMCreateSuperDM - Returns an arrays of ISes and DM encapsulating a superproblem defined by the DMs passed in. 1546 1547 Not collective 1548 1549 Input Parameter: 1550 + dms - The DM objects 1551 - len - The number of DMs 1552 1553 Output Parameters: 1554 + is - The global indices for the subproblem 1555 - superdm - The DM for the superproblem 1556 1557 Level: intermediate 1558 1559 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS() 1560 @*/ 1561 PetscErrorCode DMCreateSuperDM(DM dms[], PetscInt len, IS **is, DM *superdm) 1562 { 1563 PetscInt i; 1564 PetscErrorCode ierr; 1565 1566 PetscFunctionBegin; 1567 PetscValidPointer(dms,1); 1568 for (i = 0; i < len; ++i) {PetscValidHeaderSpecific(dms[i],DM_CLASSID,1);} 1569 if (is) PetscValidPointer(is,3); 1570 if (superdm) PetscValidPointer(superdm,4); 1571 if (len < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of DMs must be nonnegative: %D", len); 1572 if (len) { 1573 if (dms[0]->ops->createsuperdm) { 1574 ierr = (*dms[0]->ops->createsuperdm)(dms, len, is, superdm);CHKERRQ(ierr); 1575 } else SETERRQ(PetscObjectComm((PetscObject) dms[0]), PETSC_ERR_SUP, "This type has no DMCreateSuperDM implementation defined"); 1576 } 1577 PetscFunctionReturn(0); 1578 } 1579 1580 1581 /*@C 1582 DMCreateDomainDecomposition - Returns lists of IS objects defining a decomposition of a problem into subproblems 1583 corresponding to restrictions to pairs nested subdomains: each IS contains the global 1584 indices of the dofs of the corresponding subdomains. The inner subdomains conceptually 1585 define a nonoverlapping covering, while outer subdomains can overlap. 1586 The optional list of DMs define the DM for each subproblem. 1587 1588 Not collective 1589 1590 Input Parameter: 1591 . dm - the DM object 1592 1593 Output Parameters: 1594 + len - The number of subproblems in the domain decomposition (or NULL if not requested) 1595 . namelist - The name for each subdomain (or NULL if not requested) 1596 . innerislist - The global indices for each inner subdomain (or NULL, if not requested) 1597 . outerislist - The global indices for each outer subdomain (or NULL, if not requested) 1598 - dmlist - The DMs for each subdomain subproblem (or NULL, if not requested; if NULL is returned, no DMs are defined) 1599 1600 Level: intermediate 1601 1602 Notes: 1603 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 1604 PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(), 1605 and all of the arrays should be freed with PetscFree(). 1606 1607 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateDomainDecompositionDM(), DMCreateFieldDecomposition() 1608 @*/ 1609 PetscErrorCode DMCreateDomainDecomposition(DM dm, PetscInt *len, char ***namelist, IS **innerislist, IS **outerislist, DM **dmlist) 1610 { 1611 PetscErrorCode ierr; 1612 DMSubDomainHookLink link; 1613 PetscInt i,l; 1614 1615 PetscFunctionBegin; 1616 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1617 if (len) {PetscValidPointer(len,2); *len = 0;} 1618 if (namelist) {PetscValidPointer(namelist,3); *namelist = NULL;} 1619 if (innerislist) {PetscValidPointer(innerislist,4); *innerislist = NULL;} 1620 if (outerislist) {PetscValidPointer(outerislist,5); *outerislist = NULL;} 1621 if (dmlist) {PetscValidPointer(dmlist,6); *dmlist = NULL;} 1622 /* 1623 Is it a good idea to apply the following check across all impls? 1624 Perhaps some impls can have a well-defined decomposition before DMSetUp? 1625 This, however, follows the general principle that accessors are not well-behaved until the object is set up. 1626 */ 1627 if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE, "Decomposition defined only after DMSetUp"); 1628 if (dm->ops->createdomaindecomposition) { 1629 ierr = (*dm->ops->createdomaindecomposition)(dm,&l,namelist,innerislist,outerislist,dmlist);CHKERRQ(ierr); 1630 /* copy subdomain hooks and context over to the subdomain DMs */ 1631 if (dmlist && *dmlist) { 1632 for (i = 0; i < l; i++) { 1633 for (link=dm->subdomainhook; link; link=link->next) { 1634 if (link->ddhook) {ierr = (*link->ddhook)(dm,(*dmlist)[i],link->ctx);CHKERRQ(ierr);} 1635 } 1636 if (dm->ctx) (*dmlist)[i]->ctx = dm->ctx; 1637 } 1638 } 1639 if (len) *len = l; 1640 } 1641 PetscFunctionReturn(0); 1642 } 1643 1644 1645 /*@C 1646 DMCreateDomainDecompositionScatters - Returns scatters to the subdomain vectors from the global vector 1647 1648 Not collective 1649 1650 Input Parameters: 1651 + dm - the DM object 1652 . n - the number of subdomain scatters 1653 - subdms - the local subdomains 1654 1655 Output Parameters: 1656 + n - the number of scatters returned 1657 . iscat - scatter from global vector to nonoverlapping global vector entries on subdomain 1658 . oscat - scatter from global vector to overlapping global vector entries on subdomain 1659 - gscat - scatter from global vector to local vector on subdomain (fills in ghosts) 1660 1661 Notes: This is an alternative to the iis and ois arguments in DMCreateDomainDecomposition that allow for the solution 1662 of general nonlinear problems with overlapping subdomain methods. While merely having index sets that enable subsets 1663 of the residual equations to be created is fine for linear problems, nonlinear problems require local assembly of 1664 solution and residual data. 1665 1666 Level: developer 1667 1668 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS() 1669 @*/ 1670 PetscErrorCode DMCreateDomainDecompositionScatters(DM dm,PetscInt n,DM *subdms,VecScatter **iscat,VecScatter **oscat,VecScatter **gscat) 1671 { 1672 PetscErrorCode ierr; 1673 1674 PetscFunctionBegin; 1675 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1676 PetscValidPointer(subdms,3); 1677 if (dm->ops->createddscatters) { 1678 ierr = (*dm->ops->createddscatters)(dm,n,subdms,iscat,oscat,gscat);CHKERRQ(ierr); 1679 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This type has no DMCreateDomainDecompositionScatter implementation defined"); 1680 PetscFunctionReturn(0); 1681 } 1682 1683 /*@ 1684 DMRefine - Refines a DM object 1685 1686 Collective on DM 1687 1688 Input Parameter: 1689 + dm - the DM object 1690 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 1691 1692 Output Parameter: 1693 . dmf - the refined DM, or NULL 1694 1695 Note: If no refinement was done, the return value is NULL 1696 1697 Level: developer 1698 1699 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1700 @*/ 1701 PetscErrorCode DMRefine(DM dm,MPI_Comm comm,DM *dmf) 1702 { 1703 PetscErrorCode ierr; 1704 DMRefineHookLink link; 1705 1706 PetscFunctionBegin; 1707 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1708 ierr = PetscLogEventBegin(DM_Refine,dm,0,0,0);CHKERRQ(ierr); 1709 if (!dm->ops->refine) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"This DM cannot refine"); 1710 ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr); 1711 if (*dmf) { 1712 (*dmf)->ops->creatematrix = dm->ops->creatematrix; 1713 1714 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr); 1715 1716 (*dmf)->ctx = dm->ctx; 1717 (*dmf)->leveldown = dm->leveldown; 1718 (*dmf)->levelup = dm->levelup + 1; 1719 1720 ierr = DMSetMatType(*dmf,dm->mattype);CHKERRQ(ierr); 1721 for (link=dm->refinehook; link; link=link->next) { 1722 if (link->refinehook) { 1723 ierr = (*link->refinehook)(dm,*dmf,link->ctx);CHKERRQ(ierr); 1724 } 1725 } 1726 } 1727 ierr = PetscLogEventEnd(DM_Refine,dm,0,0,0);CHKERRQ(ierr); 1728 PetscFunctionReturn(0); 1729 } 1730 1731 /*@C 1732 DMRefineHookAdd - adds a callback to be run when interpolating a nonlinear problem to a finer grid 1733 1734 Logically Collective 1735 1736 Input Arguments: 1737 + coarse - nonlinear solver context on which to run a hook when restricting to a coarser level 1738 . refinehook - function to run when setting up a coarser level 1739 . interphook - function to run to update data on finer levels (once per SNESSolve()) 1740 - ctx - [optional] user-defined context for provide data for the hooks (may be NULL) 1741 1742 Calling sequence of refinehook: 1743 $ refinehook(DM coarse,DM fine,void *ctx); 1744 1745 + coarse - coarse level DM 1746 . fine - fine level DM to interpolate problem to 1747 - ctx - optional user-defined function context 1748 1749 Calling sequence for interphook: 1750 $ interphook(DM coarse,Mat interp,DM fine,void *ctx) 1751 1752 + coarse - coarse level DM 1753 . interp - matrix interpolating a coarse-level solution to the finer grid 1754 . fine - fine level DM to update 1755 - ctx - optional user-defined function context 1756 1757 Level: advanced 1758 1759 Notes: 1760 This function is only needed if auxiliary data needs to be passed to fine grids while grid sequencing 1761 1762 If this function is called multiple times, the hooks will be run in the order they are added. 1763 1764 This function is currently not available from Fortran. 1765 1766 .seealso: DMCoarsenHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 1767 @*/ 1768 PetscErrorCode DMRefineHookAdd(DM coarse,PetscErrorCode (*refinehook)(DM,DM,void*),PetscErrorCode (*interphook)(DM,Mat,DM,void*),void *ctx) 1769 { 1770 PetscErrorCode ierr; 1771 DMRefineHookLink link,*p; 1772 1773 PetscFunctionBegin; 1774 PetscValidHeaderSpecific(coarse,DM_CLASSID,1); 1775 for (p=&coarse->refinehook; *p; p=&(*p)->next) { /* Scan to the end of the current list of hooks */ 1776 if ((*p)->refinehook == refinehook && (*p)->interphook == interphook && (*p)->ctx == ctx) PetscFunctionReturn(0); 1777 } 1778 ierr = PetscNew(&link);CHKERRQ(ierr); 1779 link->refinehook = refinehook; 1780 link->interphook = interphook; 1781 link->ctx = ctx; 1782 link->next = NULL; 1783 *p = link; 1784 PetscFunctionReturn(0); 1785 } 1786 1787 /*@C 1788 DMRefineHookRemove - remove a callback from the list of hooks to be run when interpolating a nonlinear problem to a finer grid 1789 1790 Logically Collective 1791 1792 Input Arguments: 1793 + coarse - nonlinear solver context on which to run a hook when restricting to a coarser level 1794 . refinehook - function to run when setting up a coarser level 1795 . interphook - function to run to update data on finer levels (once per SNESSolve()) 1796 - ctx - [optional] user-defined context for provide data for the hooks (may be NULL) 1797 1798 Level: advanced 1799 1800 Notes: 1801 This function does nothing if the hook is not in the list. 1802 1803 This function is currently not available from Fortran. 1804 1805 .seealso: DMCoarsenHookRemove(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 1806 @*/ 1807 PetscErrorCode DMRefineHookRemove(DM coarse,PetscErrorCode (*refinehook)(DM,DM,void*),PetscErrorCode (*interphook)(DM,Mat,DM,void*),void *ctx) 1808 { 1809 PetscErrorCode ierr; 1810 DMRefineHookLink link,*p; 1811 1812 PetscFunctionBegin; 1813 PetscValidHeaderSpecific(coarse,DM_CLASSID,1); 1814 for (p=&coarse->refinehook; *p; p=&(*p)->next) { /* Search the list of current hooks */ 1815 if ((*p)->refinehook == refinehook && (*p)->interphook == interphook && (*p)->ctx == ctx) { 1816 link = *p; 1817 *p = link->next; 1818 ierr = PetscFree(link);CHKERRQ(ierr); 1819 break; 1820 } 1821 } 1822 PetscFunctionReturn(0); 1823 } 1824 1825 /*@ 1826 DMInterpolate - interpolates user-defined problem data to a finer DM by running hooks registered by DMRefineHookAdd() 1827 1828 Collective if any hooks are 1829 1830 Input Arguments: 1831 + coarse - coarser DM to use as a base 1832 . restrct - interpolation matrix, apply using MatInterpolate() 1833 - fine - finer DM to update 1834 1835 Level: developer 1836 1837 .seealso: DMRefineHookAdd(), MatInterpolate() 1838 @*/ 1839 PetscErrorCode DMInterpolate(DM coarse,Mat interp,DM fine) 1840 { 1841 PetscErrorCode ierr; 1842 DMRefineHookLink link; 1843 1844 PetscFunctionBegin; 1845 for (link=fine->refinehook; link; link=link->next) { 1846 if (link->interphook) { 1847 ierr = (*link->interphook)(coarse,interp,fine,link->ctx);CHKERRQ(ierr); 1848 } 1849 } 1850 PetscFunctionReturn(0); 1851 } 1852 1853 /*@ 1854 DMGetRefineLevel - Get's the number of refinements that have generated this DM. 1855 1856 Not Collective 1857 1858 Input Parameter: 1859 . dm - the DM object 1860 1861 Output Parameter: 1862 . level - number of refinements 1863 1864 Level: developer 1865 1866 .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1867 1868 @*/ 1869 PetscErrorCode DMGetRefineLevel(DM dm,PetscInt *level) 1870 { 1871 PetscFunctionBegin; 1872 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1873 *level = dm->levelup; 1874 PetscFunctionReturn(0); 1875 } 1876 1877 /*@ 1878 DMSetRefineLevel - Set's the number of refinements that have generated this DM. 1879 1880 Not Collective 1881 1882 Input Parameter: 1883 + dm - the DM object 1884 - level - number of refinements 1885 1886 Level: advanced 1887 1888 Notes: This value is used by PCMG to determine how many multigrid levels to use 1889 1890 .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1891 1892 @*/ 1893 PetscErrorCode DMSetRefineLevel(DM dm,PetscInt level) 1894 { 1895 PetscFunctionBegin; 1896 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1897 dm->levelup = level; 1898 PetscFunctionReturn(0); 1899 } 1900 1901 /*@C 1902 DMGlobalToLocalHookAdd - adds a callback to be run when global to local is called 1903 1904 Logically Collective 1905 1906 Input Arguments: 1907 + dm - the DM 1908 . beginhook - function to run at the beginning of DMGlobalToLocalBegin() 1909 . endhook - function to run after DMGlobalToLocalEnd() has completed 1910 - ctx - [optional] user-defined context for provide data for the hooks (may be NULL) 1911 1912 Calling sequence for beginhook: 1913 $ beginhook(DM fine,VecScatter out,VecScatter in,DM coarse,void *ctx) 1914 1915 + dm - global DM 1916 . g - global vector 1917 . mode - mode 1918 . l - local vector 1919 - ctx - optional user-defined function context 1920 1921 1922 Calling sequence for endhook: 1923 $ endhook(DM fine,VecScatter out,VecScatter in,DM coarse,void *ctx) 1924 1925 + global - global DM 1926 - ctx - optional user-defined function context 1927 1928 Level: advanced 1929 1930 .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 1931 @*/ 1932 PetscErrorCode DMGlobalToLocalHookAdd(DM dm,PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*),void *ctx) 1933 { 1934 PetscErrorCode ierr; 1935 DMGlobalToLocalHookLink link,*p; 1936 1937 PetscFunctionBegin; 1938 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1939 for (p=&dm->gtolhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */ 1940 ierr = PetscNew(&link);CHKERRQ(ierr); 1941 link->beginhook = beginhook; 1942 link->endhook = endhook; 1943 link->ctx = ctx; 1944 link->next = NULL; 1945 *p = link; 1946 PetscFunctionReturn(0); 1947 } 1948 1949 static PetscErrorCode DMGlobalToLocalHook_Constraints(DM dm, Vec g, InsertMode mode, Vec l, void *ctx) 1950 { 1951 Mat cMat; 1952 Vec cVec; 1953 PetscSection section, cSec; 1954 PetscInt pStart, pEnd, p, dof; 1955 PetscErrorCode ierr; 1956 1957 PetscFunctionBegin; 1958 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1959 ierr = DMGetDefaultConstraints(dm,&cSec,&cMat);CHKERRQ(ierr); 1960 if (cMat && (mode == INSERT_VALUES || mode == INSERT_ALL_VALUES || mode == INSERT_BC_VALUES)) { 1961 PetscInt nRows; 1962 1963 ierr = MatGetSize(cMat,&nRows,NULL);CHKERRQ(ierr); 1964 if (nRows <= 0) PetscFunctionReturn(0); 1965 ierr = DMGetDefaultSection(dm,§ion);CHKERRQ(ierr); 1966 ierr = MatCreateVecs(cMat,NULL,&cVec);CHKERRQ(ierr); 1967 ierr = MatMult(cMat,l,cVec);CHKERRQ(ierr); 1968 ierr = PetscSectionGetChart(cSec,&pStart,&pEnd);CHKERRQ(ierr); 1969 for (p = pStart; p < pEnd; p++) { 1970 ierr = PetscSectionGetDof(cSec,p,&dof);CHKERRQ(ierr); 1971 if (dof) { 1972 PetscScalar *vals; 1973 ierr = VecGetValuesSection(cVec,cSec,p,&vals);CHKERRQ(ierr); 1974 ierr = VecSetValuesSection(l,section,p,vals,INSERT_ALL_VALUES);CHKERRQ(ierr); 1975 } 1976 } 1977 ierr = VecDestroy(&cVec);CHKERRQ(ierr); 1978 } 1979 PetscFunctionReturn(0); 1980 } 1981 1982 /*@ 1983 DMGlobalToLocalBegin - Begins updating local vectors from global vector 1984 1985 Neighbor-wise Collective on DM 1986 1987 Input Parameters: 1988 + dm - the DM object 1989 . g - the global vector 1990 . mode - INSERT_VALUES or ADD_VALUES 1991 - l - the local vector 1992 1993 1994 Level: beginner 1995 1996 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 1997 1998 @*/ 1999 PetscErrorCode DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) 2000 { 2001 PetscSF sf; 2002 PetscErrorCode ierr; 2003 DMGlobalToLocalHookLink link; 2004 2005 PetscFunctionBegin; 2006 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2007 for (link=dm->gtolhook; link; link=link->next) { 2008 if (link->beginhook) { 2009 ierr = (*link->beginhook)(dm,g,mode,l,link->ctx);CHKERRQ(ierr); 2010 } 2011 } 2012 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 2013 if (sf) { 2014 const PetscScalar *gArray; 2015 PetscScalar *lArray; 2016 2017 if (mode == ADD_VALUES) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 2018 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 2019 ierr = VecGetArrayRead(g, &gArray);CHKERRQ(ierr); 2020 ierr = PetscSFBcastBegin(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr); 2021 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 2022 ierr = VecRestoreArrayRead(g, &gArray);CHKERRQ(ierr); 2023 } else { 2024 ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 2025 } 2026 PetscFunctionReturn(0); 2027 } 2028 2029 /*@ 2030 DMGlobalToLocalEnd - Ends updating local vectors from global vector 2031 2032 Neighbor-wise Collective on DM 2033 2034 Input Parameters: 2035 + dm - the DM object 2036 . g - the global vector 2037 . mode - INSERT_VALUES or ADD_VALUES 2038 - l - the local vector 2039 2040 2041 Level: beginner 2042 2043 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 2044 2045 @*/ 2046 PetscErrorCode DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) 2047 { 2048 PetscSF sf; 2049 PetscErrorCode ierr; 2050 const PetscScalar *gArray; 2051 PetscScalar *lArray; 2052 DMGlobalToLocalHookLink link; 2053 2054 PetscFunctionBegin; 2055 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2056 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 2057 if (sf) { 2058 if (mode == ADD_VALUES) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 2059 2060 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 2061 ierr = VecGetArrayRead(g, &gArray);CHKERRQ(ierr); 2062 ierr = PetscSFBcastEnd(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr); 2063 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 2064 ierr = VecRestoreArrayRead(g, &gArray);CHKERRQ(ierr); 2065 } else { 2066 ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 2067 } 2068 ierr = DMGlobalToLocalHook_Constraints(dm,g,mode,l,NULL);CHKERRQ(ierr); 2069 for (link=dm->gtolhook; link; link=link->next) { 2070 if (link->endhook) {ierr = (*link->endhook)(dm,g,mode,l,link->ctx);CHKERRQ(ierr);} 2071 } 2072 PetscFunctionReturn(0); 2073 } 2074 2075 /*@C 2076 DMLocalToGlobalHookAdd - adds a callback to be run when a local to global is called 2077 2078 Logically Collective 2079 2080 Input Arguments: 2081 + dm - the DM 2082 . beginhook - function to run at the beginning of DMLocalToGlobalBegin() 2083 . endhook - function to run after DMLocalToGlobalEnd() has completed 2084 - ctx - [optional] user-defined context for provide data for the hooks (may be NULL) 2085 2086 Calling sequence for beginhook: 2087 $ beginhook(DM fine,Vec l,InsertMode mode,Vec g,void *ctx) 2088 2089 + dm - global DM 2090 . l - local vector 2091 . mode - mode 2092 . g - global vector 2093 - ctx - optional user-defined function context 2094 2095 2096 Calling sequence for endhook: 2097 $ endhook(DM fine,Vec l,InsertMode mode,Vec g,void *ctx) 2098 2099 + global - global DM 2100 . l - local vector 2101 . mode - mode 2102 . g - global vector 2103 - ctx - optional user-defined function context 2104 2105 Level: advanced 2106 2107 .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 2108 @*/ 2109 PetscErrorCode DMLocalToGlobalHookAdd(DM dm,PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*),void *ctx) 2110 { 2111 PetscErrorCode ierr; 2112 DMLocalToGlobalHookLink link,*p; 2113 2114 PetscFunctionBegin; 2115 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2116 for (p=&dm->ltoghook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */ 2117 ierr = PetscNew(&link);CHKERRQ(ierr); 2118 link->beginhook = beginhook; 2119 link->endhook = endhook; 2120 link->ctx = ctx; 2121 link->next = NULL; 2122 *p = link; 2123 PetscFunctionReturn(0); 2124 } 2125 2126 static PetscErrorCode DMLocalToGlobalHook_Constraints(DM dm, Vec l, InsertMode mode, Vec g, void *ctx) 2127 { 2128 Mat cMat; 2129 Vec cVec; 2130 PetscSection section, cSec; 2131 PetscInt pStart, pEnd, p, dof; 2132 PetscErrorCode ierr; 2133 2134 PetscFunctionBegin; 2135 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2136 ierr = DMGetDefaultConstraints(dm,&cSec,&cMat);CHKERRQ(ierr); 2137 if (cMat && (mode == ADD_VALUES || mode == ADD_ALL_VALUES || mode == ADD_BC_VALUES)) { 2138 PetscInt nRows; 2139 2140 ierr = MatGetSize(cMat,&nRows,NULL);CHKERRQ(ierr); 2141 if (nRows <= 0) PetscFunctionReturn(0); 2142 ierr = DMGetDefaultSection(dm,§ion);CHKERRQ(ierr); 2143 ierr = MatCreateVecs(cMat,NULL,&cVec);CHKERRQ(ierr); 2144 ierr = PetscSectionGetChart(cSec,&pStart,&pEnd);CHKERRQ(ierr); 2145 for (p = pStart; p < pEnd; p++) { 2146 ierr = PetscSectionGetDof(cSec,p,&dof);CHKERRQ(ierr); 2147 if (dof) { 2148 PetscInt d; 2149 PetscScalar *vals; 2150 ierr = VecGetValuesSection(l,section,p,&vals);CHKERRQ(ierr); 2151 ierr = VecSetValuesSection(cVec,cSec,p,vals,mode);CHKERRQ(ierr); 2152 /* for this to be the true transpose, we have to zero the values that 2153 * we just extracted */ 2154 for (d = 0; d < dof; d++) { 2155 vals[d] = 0.; 2156 } 2157 } 2158 } 2159 ierr = MatMultTransposeAdd(cMat,cVec,l,l);CHKERRQ(ierr); 2160 ierr = VecDestroy(&cVec);CHKERRQ(ierr); 2161 } 2162 PetscFunctionReturn(0); 2163 } 2164 2165 /*@ 2166 DMLocalToGlobalBegin - updates global vectors from local vectors 2167 2168 Neighbor-wise Collective on DM 2169 2170 Input Parameters: 2171 + dm - the DM object 2172 . l - the local vector 2173 . mode - if INSERT_VALUES then no parallel communication is used, if ADD_VALUES then all ghost points from the same base point accumulate into that base point. 2174 - g - the global vector 2175 2176 Notes: In the ADD_VALUES case you normally would zero the receiving vector before beginning this operation. 2177 INSERT_VALUES is not supported for DMDA, in that case simply compute the values directly into a global vector instead of a local one. 2178 2179 Level: beginner 2180 2181 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin() 2182 2183 @*/ 2184 PetscErrorCode DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g) 2185 { 2186 PetscSF sf; 2187 PetscSection s, gs; 2188 DMLocalToGlobalHookLink link; 2189 const PetscScalar *lArray; 2190 PetscScalar *gArray; 2191 PetscBool isInsert; 2192 PetscErrorCode ierr; 2193 2194 PetscFunctionBegin; 2195 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2196 for (link=dm->ltoghook; link; link=link->next) { 2197 if (link->beginhook) { 2198 ierr = (*link->beginhook)(dm,l,mode,g,link->ctx);CHKERRQ(ierr); 2199 } 2200 } 2201 ierr = DMLocalToGlobalHook_Constraints(dm,l,mode,g,NULL);CHKERRQ(ierr); 2202 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 2203 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 2204 switch (mode) { 2205 case INSERT_VALUES: 2206 case INSERT_ALL_VALUES: 2207 case INSERT_BC_VALUES: 2208 isInsert = PETSC_TRUE; break; 2209 case ADD_VALUES: 2210 case ADD_ALL_VALUES: 2211 case ADD_BC_VALUES: 2212 isInsert = PETSC_FALSE; break; 2213 default: 2214 SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 2215 } 2216 if (sf && !isInsert) { 2217 ierr = VecGetArrayRead(l, &lArray);CHKERRQ(ierr); 2218 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 2219 ierr = PetscSFReduceBegin(sf, MPIU_SCALAR, lArray, gArray, MPIU_SUM);CHKERRQ(ierr); 2220 ierr = VecRestoreArrayRead(l, &lArray);CHKERRQ(ierr); 2221 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 2222 } else if (s && isInsert) { 2223 PetscInt gStart, pStart, pEnd, p; 2224 2225 ierr = DMGetDefaultGlobalSection(dm, &gs);CHKERRQ(ierr); 2226 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 2227 ierr = VecGetOwnershipRange(g, &gStart, NULL);CHKERRQ(ierr); 2228 ierr = VecGetArrayRead(l, &lArray);CHKERRQ(ierr); 2229 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 2230 for (p = pStart; p < pEnd; ++p) { 2231 PetscInt dof, gdof, cdof, gcdof, off, goff, d, e; 2232 2233 ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr); 2234 ierr = PetscSectionGetDof(gs, p, &gdof);CHKERRQ(ierr); 2235 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 2236 ierr = PetscSectionGetConstraintDof(gs, p, &gcdof);CHKERRQ(ierr); 2237 ierr = PetscSectionGetOffset(s, p, &off);CHKERRQ(ierr); 2238 ierr = PetscSectionGetOffset(gs, p, &goff);CHKERRQ(ierr); 2239 /* Ignore off-process data and points with no global data */ 2240 if (!gdof || goff < 0) continue; 2241 if (dof != gdof) SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Inconsistent sizes, p: %d dof: %d gdof: %d cdof: %d gcdof: %d", p, dof, gdof, cdof, gcdof); 2242 /* If no constraints are enforced in the global vector */ 2243 if (!gcdof) { 2244 for (d = 0; d < dof; ++d) gArray[goff-gStart+d] = lArray[off+d]; 2245 /* If constraints are enforced in the global vector */ 2246 } else if (cdof == gcdof) { 2247 const PetscInt *cdofs; 2248 PetscInt cind = 0; 2249 2250 ierr = PetscSectionGetConstraintIndices(s, p, &cdofs);CHKERRQ(ierr); 2251 for (d = 0, e = 0; d < dof; ++d) { 2252 if ((cind < cdof) && (d == cdofs[cind])) {++cind; continue;} 2253 gArray[goff-gStart+e++] = lArray[off+d]; 2254 } 2255 } else SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Inconsistent sizes, p: %d dof: %d gdof: %d cdof: %d gcdof: %d", p, dof, gdof, cdof, gcdof); 2256 } 2257 ierr = VecRestoreArrayRead(l, &lArray);CHKERRQ(ierr); 2258 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 2259 } else { 2260 ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 2261 } 2262 PetscFunctionReturn(0); 2263 } 2264 2265 /*@ 2266 DMLocalToGlobalEnd - updates global vectors from local vectors 2267 2268 Neighbor-wise Collective on DM 2269 2270 Input Parameters: 2271 + dm - the DM object 2272 . l - the local vector 2273 . mode - INSERT_VALUES or ADD_VALUES 2274 - g - the global vector 2275 2276 2277 Level: beginner 2278 2279 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd() 2280 2281 @*/ 2282 PetscErrorCode DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g) 2283 { 2284 PetscSF sf; 2285 PetscSection s; 2286 DMLocalToGlobalHookLink link; 2287 PetscBool isInsert; 2288 PetscErrorCode ierr; 2289 2290 PetscFunctionBegin; 2291 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2292 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 2293 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 2294 switch (mode) { 2295 case INSERT_VALUES: 2296 case INSERT_ALL_VALUES: 2297 isInsert = PETSC_TRUE; break; 2298 case ADD_VALUES: 2299 case ADD_ALL_VALUES: 2300 isInsert = PETSC_FALSE; break; 2301 default: 2302 SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 2303 } 2304 if (sf && !isInsert) { 2305 const PetscScalar *lArray; 2306 PetscScalar *gArray; 2307 2308 ierr = VecGetArrayRead(l, &lArray);CHKERRQ(ierr); 2309 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 2310 ierr = PetscSFReduceEnd(sf, MPIU_SCALAR, lArray, gArray, MPIU_SUM);CHKERRQ(ierr); 2311 ierr = VecRestoreArrayRead(l, &lArray);CHKERRQ(ierr); 2312 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 2313 } else if (s && isInsert) { 2314 } else { 2315 ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 2316 } 2317 for (link=dm->ltoghook; link; link=link->next) { 2318 if (link->endhook) {ierr = (*link->endhook)(dm,g,mode,l,link->ctx);CHKERRQ(ierr);} 2319 } 2320 PetscFunctionReturn(0); 2321 } 2322 2323 /*@ 2324 DMLocalToLocalBegin - Maps from a local vector (including ghost points 2325 that contain irrelevant values) to another local vector where the ghost 2326 points in the second are set correctly. Must be followed by DMLocalToLocalEnd(). 2327 2328 Neighbor-wise Collective on DM and Vec 2329 2330 Input Parameters: 2331 + dm - the DM object 2332 . g - the original local vector 2333 - mode - one of INSERT_VALUES or ADD_VALUES 2334 2335 Output Parameter: 2336 . l - the local vector with correct ghost values 2337 2338 Level: intermediate 2339 2340 Notes: 2341 The local vectors used here need not be the same as those 2342 obtained from DMCreateLocalVector(), BUT they 2343 must have the same parallel data layout; they could, for example, be 2344 obtained with VecDuplicate() from the DM originating vectors. 2345 2346 .keywords: DM, local-to-local, begin 2347 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateLocalVector(), DMCreateGlobalVector(), DMCreateInterpolation(), DMLocalToLocalEnd(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 2348 2349 @*/ 2350 PetscErrorCode DMLocalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) 2351 { 2352 PetscErrorCode ierr; 2353 2354 PetscFunctionBegin; 2355 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2356 ierr = (*dm->ops->localtolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 2357 PetscFunctionReturn(0); 2358 } 2359 2360 /*@ 2361 DMLocalToLocalEnd - Maps from a local vector (including ghost points 2362 that contain irrelevant values) to another local vector where the ghost 2363 points in the second are set correctly. Must be preceded by DMLocalToLocalBegin(). 2364 2365 Neighbor-wise Collective on DM and Vec 2366 2367 Input Parameters: 2368 + da - the DM object 2369 . g - the original local vector 2370 - mode - one of INSERT_VALUES or ADD_VALUES 2371 2372 Output Parameter: 2373 . l - the local vector with correct ghost values 2374 2375 Level: intermediate 2376 2377 Notes: 2378 The local vectors used here need not be the same as those 2379 obtained from DMCreateLocalVector(), BUT they 2380 must have the same parallel data layout; they could, for example, be 2381 obtained with VecDuplicate() from the DM originating vectors. 2382 2383 .keywords: DM, local-to-local, end 2384 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateLocalVector(), DMCreateGlobalVector(), DMCreateInterpolation(), DMLocalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 2385 2386 @*/ 2387 PetscErrorCode DMLocalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) 2388 { 2389 PetscErrorCode ierr; 2390 2391 PetscFunctionBegin; 2392 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2393 ierr = (*dm->ops->localtolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 2394 PetscFunctionReturn(0); 2395 } 2396 2397 2398 /*@ 2399 DMCoarsen - Coarsens a DM object 2400 2401 Collective on DM 2402 2403 Input Parameter: 2404 + dm - the DM object 2405 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 2406 2407 Output Parameter: 2408 . dmc - the coarsened DM 2409 2410 Level: developer 2411 2412 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 2413 2414 @*/ 2415 PetscErrorCode DMCoarsen(DM dm, MPI_Comm comm, DM *dmc) 2416 { 2417 PetscErrorCode ierr; 2418 DMCoarsenHookLink link; 2419 2420 PetscFunctionBegin; 2421 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2422 if (!dm->ops->coarsen) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"This DM cannot coarsen"); 2423 ierr = PetscLogEventBegin(DM_Coarsen,dm,0,0,0);CHKERRQ(ierr); 2424 ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr); 2425 if (!(*dmc)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "NULL coarse mesh produced"); 2426 ierr = DMSetCoarseDM(dm,*dmc);CHKERRQ(ierr); 2427 (*dmc)->ops->creatematrix = dm->ops->creatematrix; 2428 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr); 2429 (*dmc)->ctx = dm->ctx; 2430 (*dmc)->levelup = dm->levelup; 2431 (*dmc)->leveldown = dm->leveldown + 1; 2432 ierr = DMSetMatType(*dmc,dm->mattype);CHKERRQ(ierr); 2433 for (link=dm->coarsenhook; link; link=link->next) { 2434 if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);} 2435 } 2436 ierr = PetscLogEventEnd(DM_Coarsen,dm,0,0,0);CHKERRQ(ierr); 2437 PetscFunctionReturn(0); 2438 } 2439 2440 /*@C 2441 DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid 2442 2443 Logically Collective 2444 2445 Input Arguments: 2446 + fine - nonlinear solver context on which to run a hook when restricting to a coarser level 2447 . coarsenhook - function to run when setting up a coarser level 2448 . restricthook - function to run to update data on coarser levels (once per SNESSolve()) 2449 - ctx - [optional] user-defined context for provide data for the hooks (may be NULL) 2450 2451 Calling sequence of coarsenhook: 2452 $ coarsenhook(DM fine,DM coarse,void *ctx); 2453 2454 + fine - fine level DM 2455 . coarse - coarse level DM to restrict problem to 2456 - ctx - optional user-defined function context 2457 2458 Calling sequence for restricthook: 2459 $ restricthook(DM fine,Mat mrestrict,Vec rscale,Mat inject,DM coarse,void *ctx) 2460 2461 + fine - fine level DM 2462 . mrestrict - matrix restricting a fine-level solution to the coarse grid 2463 . rscale - scaling vector for restriction 2464 . inject - matrix restricting by injection 2465 . coarse - coarse level DM to update 2466 - ctx - optional user-defined function context 2467 2468 Level: advanced 2469 2470 Notes: 2471 This function is only needed if auxiliary data needs to be set up on coarse grids. 2472 2473 If this function is called multiple times, the hooks will be run in the order they are added. 2474 2475 In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to 2476 extract the finest level information from its context (instead of from the SNES). 2477 2478 This function is currently not available from Fortran. 2479 2480 .seealso: DMCoarsenHookRemove(), DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 2481 @*/ 2482 PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx) 2483 { 2484 PetscErrorCode ierr; 2485 DMCoarsenHookLink link,*p; 2486 2487 PetscFunctionBegin; 2488 PetscValidHeaderSpecific(fine,DM_CLASSID,1); 2489 for (p=&fine->coarsenhook; *p; p=&(*p)->next) { /* Scan to the end of the current list of hooks */ 2490 if ((*p)->coarsenhook == coarsenhook && (*p)->restricthook == restricthook && (*p)->ctx == ctx) PetscFunctionReturn(0); 2491 } 2492 ierr = PetscNew(&link);CHKERRQ(ierr); 2493 link->coarsenhook = coarsenhook; 2494 link->restricthook = restricthook; 2495 link->ctx = ctx; 2496 link->next = NULL; 2497 *p = link; 2498 PetscFunctionReturn(0); 2499 } 2500 2501 /*@C 2502 DMCoarsenHookRemove - remove a callback from the list of hooks to be run when restricting a nonlinear problem to the coarse grid 2503 2504 Logically Collective 2505 2506 Input Arguments: 2507 + fine - nonlinear solver context on which to run a hook when restricting to a coarser level 2508 . coarsenhook - function to run when setting up a coarser level 2509 . restricthook - function to run to update data on coarser levels (once per SNESSolve()) 2510 - ctx - [optional] user-defined context for provide data for the hooks (may be NULL) 2511 2512 Level: advanced 2513 2514 Notes: 2515 This function does nothing if the hook is not in the list. 2516 2517 This function is currently not available from Fortran. 2518 2519 .seealso: DMCoarsenHookAdd(), DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 2520 @*/ 2521 PetscErrorCode DMCoarsenHookRemove(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx) 2522 { 2523 PetscErrorCode ierr; 2524 DMCoarsenHookLink link,*p; 2525 2526 PetscFunctionBegin; 2527 PetscValidHeaderSpecific(fine,DM_CLASSID,1); 2528 for (p=&fine->coarsenhook; *p; p=&(*p)->next) { /* Search the list of current hooks */ 2529 if ((*p)->coarsenhook == coarsenhook && (*p)->restricthook == restricthook && (*p)->ctx == ctx) { 2530 link = *p; 2531 *p = link->next; 2532 ierr = PetscFree(link);CHKERRQ(ierr); 2533 break; 2534 } 2535 } 2536 PetscFunctionReturn(0); 2537 } 2538 2539 2540 /*@ 2541 DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd() 2542 2543 Collective if any hooks are 2544 2545 Input Arguments: 2546 + fine - finer DM to use as a base 2547 . restrct - restriction matrix, apply using MatRestrict() 2548 . inject - injection matrix, also use MatRestrict() 2549 - coarse - coarer DM to update 2550 2551 Level: developer 2552 2553 .seealso: DMCoarsenHookAdd(), MatRestrict() 2554 @*/ 2555 PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse) 2556 { 2557 PetscErrorCode ierr; 2558 DMCoarsenHookLink link; 2559 2560 PetscFunctionBegin; 2561 for (link=fine->coarsenhook; link; link=link->next) { 2562 if (link->restricthook) { 2563 ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr); 2564 } 2565 } 2566 PetscFunctionReturn(0); 2567 } 2568 2569 /*@C 2570 DMSubDomainHookAdd - adds a callback to be run when restricting a problem to the coarse grid 2571 2572 Logically Collective 2573 2574 Input Arguments: 2575 + global - global DM 2576 . ddhook - function to run to pass data to the decomposition DM upon its creation 2577 . restricthook - function to run to update data on block solve (at the beginning of the block solve) 2578 - ctx - [optional] user-defined context for provide data for the hooks (may be NULL) 2579 2580 2581 Calling sequence for ddhook: 2582 $ ddhook(DM global,DM block,void *ctx) 2583 2584 + global - global DM 2585 . block - block DM 2586 - ctx - optional user-defined function context 2587 2588 Calling sequence for restricthook: 2589 $ restricthook(DM global,VecScatter out,VecScatter in,DM block,void *ctx) 2590 2591 + global - global DM 2592 . out - scatter to the outer (with ghost and overlap points) block vector 2593 . in - scatter to block vector values only owned locally 2594 . block - block DM 2595 - ctx - optional user-defined function context 2596 2597 Level: advanced 2598 2599 Notes: 2600 This function is only needed if auxiliary data needs to be set up on subdomain DMs. 2601 2602 If this function is called multiple times, the hooks will be run in the order they are added. 2603 2604 In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to 2605 extract the global information from its context (instead of from the SNES). 2606 2607 This function is currently not available from Fortran. 2608 2609 .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 2610 @*/ 2611 PetscErrorCode DMSubDomainHookAdd(DM global,PetscErrorCode (*ddhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*),void *ctx) 2612 { 2613 PetscErrorCode ierr; 2614 DMSubDomainHookLink link,*p; 2615 2616 PetscFunctionBegin; 2617 PetscValidHeaderSpecific(global,DM_CLASSID,1); 2618 for (p=&global->subdomainhook; *p; p=&(*p)->next) { /* Scan to the end of the current list of hooks */ 2619 if ((*p)->ddhook == ddhook && (*p)->restricthook == restricthook && (*p)->ctx == ctx) PetscFunctionReturn(0); 2620 } 2621 ierr = PetscNew(&link);CHKERRQ(ierr); 2622 link->restricthook = restricthook; 2623 link->ddhook = ddhook; 2624 link->ctx = ctx; 2625 link->next = NULL; 2626 *p = link; 2627 PetscFunctionReturn(0); 2628 } 2629 2630 /*@C 2631 DMSubDomainHookRemove - remove a callback from the list to be run when restricting a problem to the coarse grid 2632 2633 Logically Collective 2634 2635 Input Arguments: 2636 + global - global DM 2637 . ddhook - function to run to pass data to the decomposition DM upon its creation 2638 . restricthook - function to run to update data on block solve (at the beginning of the block solve) 2639 - ctx - [optional] user-defined context for provide data for the hooks (may be NULL) 2640 2641 Level: advanced 2642 2643 Notes: 2644 2645 This function is currently not available from Fortran. 2646 2647 .seealso: DMSubDomainHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 2648 @*/ 2649 PetscErrorCode DMSubDomainHookRemove(DM global,PetscErrorCode (*ddhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*),void *ctx) 2650 { 2651 PetscErrorCode ierr; 2652 DMSubDomainHookLink link,*p; 2653 2654 PetscFunctionBegin; 2655 PetscValidHeaderSpecific(global,DM_CLASSID,1); 2656 for (p=&global->subdomainhook; *p; p=&(*p)->next) { /* Search the list of current hooks */ 2657 if ((*p)->ddhook == ddhook && (*p)->restricthook == restricthook && (*p)->ctx == ctx) { 2658 link = *p; 2659 *p = link->next; 2660 ierr = PetscFree(link);CHKERRQ(ierr); 2661 break; 2662 } 2663 } 2664 PetscFunctionReturn(0); 2665 } 2666 2667 /*@ 2668 DMSubDomainRestrict - restricts user-defined problem data to a block DM by running hooks registered by DMSubDomainHookAdd() 2669 2670 Collective if any hooks are 2671 2672 Input Arguments: 2673 + fine - finer DM to use as a base 2674 . oscatter - scatter from domain global vector filling subdomain global vector with overlap 2675 . gscatter - scatter from domain global vector filling subdomain local vector with ghosts 2676 - coarse - coarer DM to update 2677 2678 Level: developer 2679 2680 .seealso: DMCoarsenHookAdd(), MatRestrict() 2681 @*/ 2682 PetscErrorCode DMSubDomainRestrict(DM global,VecScatter oscatter,VecScatter gscatter,DM subdm) 2683 { 2684 PetscErrorCode ierr; 2685 DMSubDomainHookLink link; 2686 2687 PetscFunctionBegin; 2688 for (link=global->subdomainhook; link; link=link->next) { 2689 if (link->restricthook) { 2690 ierr = (*link->restricthook)(global,oscatter,gscatter,subdm,link->ctx);CHKERRQ(ierr); 2691 } 2692 } 2693 PetscFunctionReturn(0); 2694 } 2695 2696 /*@ 2697 DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM. 2698 2699 Not Collective 2700 2701 Input Parameter: 2702 . dm - the DM object 2703 2704 Output Parameter: 2705 . level - number of coarsenings 2706 2707 Level: developer 2708 2709 .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 2710 2711 @*/ 2712 PetscErrorCode DMGetCoarsenLevel(DM dm,PetscInt *level) 2713 { 2714 PetscFunctionBegin; 2715 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2716 *level = dm->leveldown; 2717 PetscFunctionReturn(0); 2718 } 2719 2720 2721 2722 /*@C 2723 DMRefineHierarchy - Refines a DM object, all levels at once 2724 2725 Collective on DM 2726 2727 Input Parameter: 2728 + dm - the DM object 2729 - nlevels - the number of levels of refinement 2730 2731 Output Parameter: 2732 . dmf - the refined DM hierarchy 2733 2734 Level: developer 2735 2736 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 2737 2738 @*/ 2739 PetscErrorCode DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[]) 2740 { 2741 PetscErrorCode ierr; 2742 2743 PetscFunctionBegin; 2744 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2745 if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 2746 if (nlevels == 0) PetscFunctionReturn(0); 2747 if (dm->ops->refinehierarchy) { 2748 ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr); 2749 } else if (dm->ops->refine) { 2750 PetscInt i; 2751 2752 ierr = DMRefine(dm,PetscObjectComm((PetscObject)dm),&dmf[0]);CHKERRQ(ierr); 2753 for (i=1; i<nlevels; i++) { 2754 ierr = DMRefine(dmf[i-1],PetscObjectComm((PetscObject)dm),&dmf[i]);CHKERRQ(ierr); 2755 } 2756 } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No RefineHierarchy for this DM yet"); 2757 PetscFunctionReturn(0); 2758 } 2759 2760 /*@C 2761 DMCoarsenHierarchy - Coarsens a DM object, all levels at once 2762 2763 Collective on DM 2764 2765 Input Parameter: 2766 + dm - the DM object 2767 - nlevels - the number of levels of coarsening 2768 2769 Output Parameter: 2770 . dmc - the coarsened DM hierarchy 2771 2772 Level: developer 2773 2774 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 2775 2776 @*/ 2777 PetscErrorCode DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[]) 2778 { 2779 PetscErrorCode ierr; 2780 2781 PetscFunctionBegin; 2782 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2783 if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 2784 if (nlevels == 0) PetscFunctionReturn(0); 2785 PetscValidPointer(dmc,3); 2786 if (dm->ops->coarsenhierarchy) { 2787 ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr); 2788 } else if (dm->ops->coarsen) { 2789 PetscInt i; 2790 2791 ierr = DMCoarsen(dm,PetscObjectComm((PetscObject)dm),&dmc[0]);CHKERRQ(ierr); 2792 for (i=1; i<nlevels; i++) { 2793 ierr = DMCoarsen(dmc[i-1],PetscObjectComm((PetscObject)dm),&dmc[i]);CHKERRQ(ierr); 2794 } 2795 } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet"); 2796 PetscFunctionReturn(0); 2797 } 2798 2799 /*@ 2800 DMCreateAggregates - Gets the aggregates that map between 2801 grids associated with two DMs. 2802 2803 Collective on DM 2804 2805 Input Parameters: 2806 + dmc - the coarse grid DM 2807 - dmf - the fine grid DM 2808 2809 Output Parameters: 2810 . rest - the restriction matrix (transpose of the projection matrix) 2811 2812 Level: intermediate 2813 2814 .keywords: interpolation, restriction, multigrid 2815 2816 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation() 2817 @*/ 2818 PetscErrorCode DMCreateAggregates(DM dmc, DM dmf, Mat *rest) 2819 { 2820 PetscErrorCode ierr; 2821 2822 PetscFunctionBegin; 2823 PetscValidHeaderSpecific(dmc,DM_CLASSID,1); 2824 PetscValidHeaderSpecific(dmf,DM_CLASSID,2); 2825 ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr); 2826 PetscFunctionReturn(0); 2827 } 2828 2829 /*@C 2830 DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed 2831 2832 Not Collective 2833 2834 Input Parameters: 2835 + dm - the DM object 2836 - destroy - the destroy function 2837 2838 Level: intermediate 2839 2840 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 2841 2842 @*/ 2843 PetscErrorCode DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**)) 2844 { 2845 PetscFunctionBegin; 2846 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2847 dm->ctxdestroy = destroy; 2848 PetscFunctionReturn(0); 2849 } 2850 2851 /*@ 2852 DMSetApplicationContext - Set a user context into a DM object 2853 2854 Not Collective 2855 2856 Input Parameters: 2857 + dm - the DM object 2858 - ctx - the user context 2859 2860 Level: intermediate 2861 2862 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 2863 2864 @*/ 2865 PetscErrorCode DMSetApplicationContext(DM dm,void *ctx) 2866 { 2867 PetscFunctionBegin; 2868 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2869 dm->ctx = ctx; 2870 PetscFunctionReturn(0); 2871 } 2872 2873 /*@ 2874 DMGetApplicationContext - Gets a user context from a DM object 2875 2876 Not Collective 2877 2878 Input Parameter: 2879 . dm - the DM object 2880 2881 Output Parameter: 2882 . ctx - the user context 2883 2884 Level: intermediate 2885 2886 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 2887 2888 @*/ 2889 PetscErrorCode DMGetApplicationContext(DM dm,void *ctx) 2890 { 2891 PetscFunctionBegin; 2892 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2893 *(void**)ctx = dm->ctx; 2894 PetscFunctionReturn(0); 2895 } 2896 2897 /*@C 2898 DMSetVariableBounds - sets a function to compute the lower and upper bound vectors for SNESVI. 2899 2900 Logically Collective on DM 2901 2902 Input Parameter: 2903 + dm - the DM object 2904 - f - the function that computes variable bounds used by SNESVI (use NULL to cancel a previous function that was set) 2905 2906 Level: intermediate 2907 2908 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), 2909 DMSetJacobian() 2910 2911 @*/ 2912 PetscErrorCode DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 2913 { 2914 PetscFunctionBegin; 2915 dm->ops->computevariablebounds = f; 2916 PetscFunctionReturn(0); 2917 } 2918 2919 /*@ 2920 DMHasVariableBounds - does the DM object have a variable bounds function? 2921 2922 Not Collective 2923 2924 Input Parameter: 2925 . dm - the DM object to destroy 2926 2927 Output Parameter: 2928 . flg - PETSC_TRUE if the variable bounds function exists 2929 2930 Level: developer 2931 2932 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 2933 2934 @*/ 2935 PetscErrorCode DMHasVariableBounds(DM dm,PetscBool *flg) 2936 { 2937 PetscFunctionBegin; 2938 *flg = (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE; 2939 PetscFunctionReturn(0); 2940 } 2941 2942 /*@C 2943 DMComputeVariableBounds - compute variable bounds used by SNESVI. 2944 2945 Logically Collective on DM 2946 2947 Input Parameters: 2948 . dm - the DM object 2949 2950 Output parameters: 2951 + xl - lower bound 2952 - xu - upper bound 2953 2954 Level: advanced 2955 2956 Notes: This is generally not called by users. It calls the function provided by the user with DMSetVariableBounds() 2957 2958 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 2959 2960 @*/ 2961 PetscErrorCode DMComputeVariableBounds(DM dm, Vec xl, Vec xu) 2962 { 2963 PetscErrorCode ierr; 2964 2965 PetscFunctionBegin; 2966 PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 2967 PetscValidHeaderSpecific(xu,VEC_CLASSID,2); 2968 if (dm->ops->computevariablebounds) { 2969 ierr = (*dm->ops->computevariablebounds)(dm, xl,xu);CHKERRQ(ierr); 2970 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "This DM is incapable of computing variable bounds."); 2971 PetscFunctionReturn(0); 2972 } 2973 2974 /*@ 2975 DMHasColoring - does the DM object have a method of providing a coloring? 2976 2977 Not Collective 2978 2979 Input Parameter: 2980 . dm - the DM object 2981 2982 Output Parameter: 2983 . flg - PETSC_TRUE if the DM has facilities for DMCreateColoring(). 2984 2985 Level: developer 2986 2987 .seealso DMHasFunction(), DMCreateColoring() 2988 2989 @*/ 2990 PetscErrorCode DMHasColoring(DM dm,PetscBool *flg) 2991 { 2992 PetscFunctionBegin; 2993 *flg = (dm->ops->getcoloring) ? PETSC_TRUE : PETSC_FALSE; 2994 PetscFunctionReturn(0); 2995 } 2996 2997 /*@ 2998 DMHasCreateRestriction - does the DM object have a method of providing a restriction? 2999 3000 Not Collective 3001 3002 Input Parameter: 3003 . dm - the DM object 3004 3005 Output Parameter: 3006 . flg - PETSC_TRUE if the DM has facilities for DMCreateRestriction(). 3007 3008 Level: developer 3009 3010 .seealso DMHasFunction(), DMCreateRestriction() 3011 3012 @*/ 3013 PetscErrorCode DMHasCreateRestriction(DM dm,PetscBool *flg) 3014 { 3015 PetscFunctionBegin; 3016 *flg = (dm->ops->createrestriction) ? PETSC_TRUE : PETSC_FALSE; 3017 PetscFunctionReturn(0); 3018 } 3019 3020 /*@C 3021 DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear. 3022 3023 Collective on DM 3024 3025 Input Parameter: 3026 + dm - the DM object 3027 - x - location to compute residual and Jacobian, if NULL is passed to those routines; will be NULL for linear problems. 3028 3029 Level: developer 3030 3031 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 3032 3033 @*/ 3034 PetscErrorCode DMSetVec(DM dm,Vec x) 3035 { 3036 PetscErrorCode ierr; 3037 3038 PetscFunctionBegin; 3039 if (x) { 3040 if (!dm->x) { 3041 ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 3042 } 3043 ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 3044 } else if (dm->x) { 3045 ierr = VecDestroy(&dm->x);CHKERRQ(ierr); 3046 } 3047 PetscFunctionReturn(0); 3048 } 3049 3050 PetscFunctionList DMList = NULL; 3051 PetscBool DMRegisterAllCalled = PETSC_FALSE; 3052 3053 /*@C 3054 DMSetType - Builds a DM, for a particular DM implementation. 3055 3056 Collective on DM 3057 3058 Input Parameters: 3059 + dm - The DM object 3060 - method - The name of the DM type 3061 3062 Options Database Key: 3063 . -dm_type <type> - Sets the DM type; use -help for a list of available types 3064 3065 Notes: 3066 See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D). 3067 3068 Level: intermediate 3069 3070 .keywords: DM, set, type 3071 .seealso: DMGetType(), DMCreate() 3072 @*/ 3073 PetscErrorCode DMSetType(DM dm, DMType method) 3074 { 3075 PetscErrorCode (*r)(DM); 3076 PetscBool match; 3077 PetscErrorCode ierr; 3078 3079 PetscFunctionBegin; 3080 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 3081 ierr = PetscObjectTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr); 3082 if (match) PetscFunctionReturn(0); 3083 3084 ierr = DMRegisterAll();CHKERRQ(ierr); 3085 ierr = PetscFunctionListFind(DMList,method,&r);CHKERRQ(ierr); 3086 if (!r) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method); 3087 3088 if (dm->ops->destroy) { 3089 ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr); 3090 dm->ops->destroy = NULL; 3091 } 3092 ierr = (*r)(dm);CHKERRQ(ierr); 3093 ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr); 3094 PetscFunctionReturn(0); 3095 } 3096 3097 /*@C 3098 DMGetType - Gets the DM type name (as a string) from the DM. 3099 3100 Not Collective 3101 3102 Input Parameter: 3103 . dm - The DM 3104 3105 Output Parameter: 3106 . type - The DM type name 3107 3108 Level: intermediate 3109 3110 .keywords: DM, get, type, name 3111 .seealso: DMSetType(), DMCreate() 3112 @*/ 3113 PetscErrorCode DMGetType(DM dm, DMType *type) 3114 { 3115 PetscErrorCode ierr; 3116 3117 PetscFunctionBegin; 3118 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 3119 PetscValidPointer(type,2); 3120 ierr = DMRegisterAll();CHKERRQ(ierr); 3121 *type = ((PetscObject)dm)->type_name; 3122 PetscFunctionReturn(0); 3123 } 3124 3125 /*@C 3126 DMConvert - Converts a DM to another DM, either of the same or different type. 3127 3128 Collective on DM 3129 3130 Input Parameters: 3131 + dm - the DM 3132 - newtype - new DM type (use "same" for the same type) 3133 3134 Output Parameter: 3135 . M - pointer to new DM 3136 3137 Notes: 3138 Cannot be used to convert a sequential DM to parallel or parallel to sequential, 3139 the MPI communicator of the generated DM is always the same as the communicator 3140 of the input DM. 3141 3142 Level: intermediate 3143 3144 .seealso: DMCreate() 3145 @*/ 3146 PetscErrorCode DMConvert(DM dm, DMType newtype, DM *M) 3147 { 3148 DM B; 3149 char convname[256]; 3150 PetscBool sametype/*, issame */; 3151 PetscErrorCode ierr; 3152 3153 PetscFunctionBegin; 3154 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3155 PetscValidType(dm,1); 3156 PetscValidPointer(M,3); 3157 ierr = PetscObjectTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr); 3158 /* ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr); */ 3159 if (sametype) { 3160 *M = dm; 3161 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 3162 PetscFunctionReturn(0); 3163 } else { 3164 PetscErrorCode (*conv)(DM, DMType, DM*) = NULL; 3165 3166 /* 3167 Order of precedence: 3168 1) See if a specialized converter is known to the current DM. 3169 2) See if a specialized converter is known to the desired DM class. 3170 3) See if a good general converter is registered for the desired class 3171 4) See if a good general converter is known for the current matrix. 3172 5) Use a really basic converter. 3173 */ 3174 3175 /* 1) See if a specialized converter is known to the current DM and the desired class */ 3176 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 3177 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 3178 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 3179 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 3180 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 3181 ierr = PetscObjectQueryFunction((PetscObject)dm,convname,&conv);CHKERRQ(ierr); 3182 if (conv) goto foundconv; 3183 3184 /* 2) See if a specialized converter is known to the desired DM class. */ 3185 ierr = DMCreate(PetscObjectComm((PetscObject)dm), &B);CHKERRQ(ierr); 3186 ierr = DMSetType(B, newtype);CHKERRQ(ierr); 3187 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 3188 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 3189 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 3190 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 3191 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 3192 ierr = PetscObjectQueryFunction((PetscObject)B,convname,&conv);CHKERRQ(ierr); 3193 if (conv) { 3194 ierr = DMDestroy(&B);CHKERRQ(ierr); 3195 goto foundconv; 3196 } 3197 3198 #if 0 3199 /* 3) See if a good general converter is registered for the desired class */ 3200 conv = B->ops->convertfrom; 3201 ierr = DMDestroy(&B);CHKERRQ(ierr); 3202 if (conv) goto foundconv; 3203 3204 /* 4) See if a good general converter is known for the current matrix */ 3205 if (dm->ops->convert) { 3206 conv = dm->ops->convert; 3207 } 3208 if (conv) goto foundconv; 3209 #endif 3210 3211 /* 5) Use a really basic converter. */ 3212 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype); 3213 3214 foundconv: 3215 ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 3216 ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr); 3217 /* Things that are independent of DM type: We should consult DMClone() here */ 3218 { 3219 PetscBool isper; 3220 const PetscReal *maxCell, *L; 3221 const DMBoundaryType *bd; 3222 ierr = DMGetPeriodicity(dm, &isper, &maxCell, &L, &bd);CHKERRQ(ierr); 3223 ierr = DMSetPeriodicity(*M, isper, maxCell, L, bd);CHKERRQ(ierr); 3224 } 3225 ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 3226 } 3227 ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr); 3228 PetscFunctionReturn(0); 3229 } 3230 3231 /*--------------------------------------------------------------------------------------------------------------------*/ 3232 3233 /*@C 3234 DMRegister - Adds a new DM component implementation 3235 3236 Not Collective 3237 3238 Input Parameters: 3239 + name - The name of a new user-defined creation routine 3240 - create_func - The creation routine itself 3241 3242 Notes: 3243 DMRegister() may be called multiple times to add several user-defined DMs 3244 3245 3246 Sample usage: 3247 .vb 3248 DMRegister("my_da", MyDMCreate); 3249 .ve 3250 3251 Then, your DM type can be chosen with the procedural interface via 3252 .vb 3253 DMCreate(MPI_Comm, DM *); 3254 DMSetType(DM,"my_da"); 3255 .ve 3256 or at runtime via the option 3257 .vb 3258 -da_type my_da 3259 .ve 3260 3261 Level: advanced 3262 3263 .keywords: DM, register 3264 .seealso: DMRegisterAll(), DMRegisterDestroy() 3265 3266 @*/ 3267 PetscErrorCode DMRegister(const char sname[],PetscErrorCode (*function)(DM)) 3268 { 3269 PetscErrorCode ierr; 3270 3271 PetscFunctionBegin; 3272 ierr = PetscFunctionListAdd(&DMList,sname,function);CHKERRQ(ierr); 3273 PetscFunctionReturn(0); 3274 } 3275 3276 /*@C 3277 DMLoad - Loads a DM that has been stored in binary with DMView(). 3278 3279 Collective on PetscViewer 3280 3281 Input Parameters: 3282 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or 3283 some related function before a call to DMLoad(). 3284 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 3285 HDF5 file viewer, obtained from PetscViewerHDF5Open() 3286 3287 Level: intermediate 3288 3289 Notes: 3290 The type is determined by the data in the file, any type set into the DM before this call is ignored. 3291 3292 Notes for advanced users: 3293 Most users should not need to know the details of the binary storage 3294 format, since DMLoad() and DMView() completely hide these details. 3295 But for anyone who's interested, the standard binary matrix storage 3296 format is 3297 .vb 3298 has not yet been determined 3299 .ve 3300 3301 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad() 3302 @*/ 3303 PetscErrorCode DMLoad(DM newdm, PetscViewer viewer) 3304 { 3305 PetscBool isbinary, ishdf5; 3306 PetscErrorCode ierr; 3307 3308 PetscFunctionBegin; 3309 PetscValidHeaderSpecific(newdm,DM_CLASSID,1); 3310 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 3311 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 3312 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 3313 if (isbinary) { 3314 PetscInt classid; 3315 char type[256]; 3316 3317 ierr = PetscViewerBinaryRead(viewer,&classid,1,NULL,PETSC_INT);CHKERRQ(ierr); 3318 if (classid != DM_FILE_CLASSID) SETERRQ1(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not DM next in file, classid found %d",(int)classid); 3319 ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr); 3320 ierr = DMSetType(newdm, type);CHKERRQ(ierr); 3321 if (newdm->ops->load) {ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);} 3322 } else if (ishdf5) { 3323 if (newdm->ops->load) {ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);} 3324 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen() or PetscViewerHDF5Open()"); 3325 PetscFunctionReturn(0); 3326 } 3327 3328 /******************************** FEM Support **********************************/ 3329 3330 PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[]) 3331 { 3332 PetscInt f; 3333 PetscErrorCode ierr; 3334 3335 PetscFunctionBegin; 3336 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 3337 for (f = 0; f < len; ++f) { 3338 ierr = PetscPrintf(PETSC_COMM_SELF, " | %g |\n", (double)PetscRealPart(x[f]));CHKERRQ(ierr); 3339 } 3340 PetscFunctionReturn(0); 3341 } 3342 3343 PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[]) 3344 { 3345 PetscInt f, g; 3346 PetscErrorCode ierr; 3347 3348 PetscFunctionBegin; 3349 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 3350 for (f = 0; f < rows; ++f) { 3351 ierr = PetscPrintf(PETSC_COMM_SELF, " |");CHKERRQ(ierr); 3352 for (g = 0; g < cols; ++g) { 3353 ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5g", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr); 3354 } 3355 ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr); 3356 } 3357 PetscFunctionReturn(0); 3358 } 3359 3360 PetscErrorCode DMPrintLocalVec(DM dm, const char name[], PetscReal tol, Vec X) 3361 { 3362 PetscInt localSize, bs; 3363 PetscMPIInt size; 3364 Vec x, xglob; 3365 const PetscScalar *xarray; 3366 PetscErrorCode ierr; 3367 3368 PetscFunctionBegin; 3369 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm),&size);CHKERRQ(ierr); 3370 ierr = VecDuplicate(X, &x);CHKERRQ(ierr); 3371 ierr = VecCopy(X, x);CHKERRQ(ierr); 3372 ierr = VecChop(x, tol);CHKERRQ(ierr); 3373 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm),"%s:\n",name);CHKERRQ(ierr); 3374 if (size > 1) { 3375 ierr = VecGetLocalSize(x,&localSize);CHKERRQ(ierr); 3376 ierr = VecGetArrayRead(x,&xarray);CHKERRQ(ierr); 3377 ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 3378 ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject) dm),bs,localSize,PETSC_DETERMINE,xarray,&xglob);CHKERRQ(ierr); 3379 } else { 3380 xglob = x; 3381 } 3382 ierr = VecView(xglob,PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject) dm)));CHKERRQ(ierr); 3383 if (size > 1) { 3384 ierr = VecDestroy(&xglob);CHKERRQ(ierr); 3385 ierr = VecRestoreArrayRead(x,&xarray);CHKERRQ(ierr); 3386 } 3387 ierr = VecDestroy(&x);CHKERRQ(ierr); 3388 PetscFunctionReturn(0); 3389 } 3390 3391 /*@ 3392 DMGetDefaultSection - Get the PetscSection encoding the local data layout for the DM. 3393 3394 Input Parameter: 3395 . dm - The DM 3396 3397 Output Parameter: 3398 . section - The PetscSection 3399 3400 Level: intermediate 3401 3402 Note: This gets a borrowed reference, so the user should not destroy this PetscSection. 3403 3404 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection() 3405 @*/ 3406 PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *section) 3407 { 3408 PetscErrorCode ierr; 3409 3410 PetscFunctionBegin; 3411 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3412 PetscValidPointer(section, 2); 3413 if (!dm->defaultSection && dm->ops->createdefaultsection) { 3414 ierr = (*dm->ops->createdefaultsection)(dm);CHKERRQ(ierr); 3415 if (dm->defaultSection) {ierr = PetscObjectViewFromOptions((PetscObject) dm->defaultSection, NULL, "-dm_petscsection_view");CHKERRQ(ierr);} 3416 } 3417 *section = dm->defaultSection; 3418 PetscFunctionReturn(0); 3419 } 3420 3421 /*@ 3422 DMSetDefaultSection - Set the PetscSection encoding the local data layout for the DM. 3423 3424 Input Parameters: 3425 + dm - The DM 3426 - section - The PetscSection 3427 3428 Level: intermediate 3429 3430 Note: Any existing Section will be destroyed 3431 3432 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection() 3433 @*/ 3434 PetscErrorCode DMSetDefaultSection(DM dm, PetscSection section) 3435 { 3436 PetscInt numFields = 0; 3437 PetscInt f; 3438 PetscErrorCode ierr; 3439 3440 PetscFunctionBegin; 3441 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3442 if (section) { 3443 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2); 3444 ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr); 3445 } 3446 ierr = PetscSectionDestroy(&dm->defaultSection);CHKERRQ(ierr); 3447 dm->defaultSection = section; 3448 if (section) {ierr = PetscSectionGetNumFields(dm->defaultSection, &numFields);CHKERRQ(ierr);} 3449 if (numFields) { 3450 ierr = DMSetNumFields(dm, numFields);CHKERRQ(ierr); 3451 for (f = 0; f < numFields; ++f) { 3452 PetscObject disc; 3453 const char *name; 3454 3455 ierr = PetscSectionGetFieldName(dm->defaultSection, f, &name);CHKERRQ(ierr); 3456 ierr = DMGetField(dm, f, &disc);CHKERRQ(ierr); 3457 ierr = PetscObjectSetName(disc, name);CHKERRQ(ierr); 3458 } 3459 } 3460 /* The global section will be rebuilt in the next call to DMGetDefaultGlobalSection(). */ 3461 ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr); 3462 PetscFunctionReturn(0); 3463 } 3464 3465 /*@ 3466 DMGetDefaultConstraints - Get the PetscSection and Mat the specify the local constraint interpolation. See DMSetDefaultConstraints() for a description of the purpose of constraint interpolation. 3467 3468 not collective 3469 3470 Input Parameter: 3471 . dm - The DM 3472 3473 Output Parameter: 3474 + section - The PetscSection describing the range of the constraint matrix: relates rows of the constraint matrix to dofs of the default section. Returns NULL if there are no local constraints. 3475 - mat - The Mat that interpolates local constraints: its width should be the layout size of the default section. Returns NULL if there are no local constraints. 3476 3477 Level: advanced 3478 3479 Note: This gets borrowed references, so the user should not destroy the PetscSection or the Mat. 3480 3481 .seealso: DMSetDefaultConstraints() 3482 @*/ 3483 PetscErrorCode DMGetDefaultConstraints(DM dm, PetscSection *section, Mat *mat) 3484 { 3485 PetscErrorCode ierr; 3486 3487 PetscFunctionBegin; 3488 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3489 if (!dm->defaultConstraintSection && !dm->defaultConstraintMat && dm->ops->createdefaultconstraints) {ierr = (*dm->ops->createdefaultconstraints)(dm);CHKERRQ(ierr);} 3490 if (section) {*section = dm->defaultConstraintSection;} 3491 if (mat) {*mat = dm->defaultConstraintMat;} 3492 PetscFunctionReturn(0); 3493 } 3494 3495 /*@ 3496 DMSetDefaultConstraints - Set the PetscSection and Mat the specify the local constraint interpolation. 3497 3498 If a constraint matrix is specified, then it is applied during DMGlobalToLocalEnd() when mode is INSERT_VALUES, INSERT_BC_VALUES, or INSERT_ALL_VALUES. Without a constraint matrix, the local vector l returned by DMGlobalToLocalEnd() contains values that have been scattered from a global vector without modification; with a constraint matrix A, l is modified by computing c = A * l, l[s[i]] = c[i], where the scatter s is defined by the PetscSection returned by DMGetDefaultConstraintMatrix(). 3499 3500 If a constraint matrix is specified, then its adjoint is applied during DMLocalToGlobalBegin() when mode is ADD_VALUES, ADD_BC_VALUES, or ADD_ALL_VALUES. Without a constraint matrix, the local vector l is accumulated into a global vector without modification; with a constraint matrix A, l is first modified by computing c[i] = l[s[i]], l[s[i]] = 0, l = l + A'*c, which is the adjoint of the operation described above. 3501 3502 collective on dm 3503 3504 Input Parameters: 3505 + dm - The DM 3506 + section - The PetscSection describing the range of the constraint matrix: relates rows of the constraint matrix to dofs of the default section. Must have a local communicator (PETSC_COMM_SELF or derivative). 3507 - mat - The Mat that interpolates local constraints: its width should be the layout size of the default section: NULL indicates no constraints. Must have a local communicator (PETSC_COMM_SELF or derivative). 3508 3509 Level: advanced 3510 3511 Note: This increments the references of the PetscSection and the Mat, so they user can destroy them 3512 3513 .seealso: DMGetDefaultConstraints() 3514 @*/ 3515 PetscErrorCode DMSetDefaultConstraints(DM dm, PetscSection section, Mat mat) 3516 { 3517 PetscMPIInt result; 3518 PetscErrorCode ierr; 3519 3520 PetscFunctionBegin; 3521 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3522 if (section) { 3523 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2); 3524 ierr = MPI_Comm_compare(PETSC_COMM_SELF,PetscObjectComm((PetscObject)section),&result);CHKERRQ(ierr); 3525 if (result != MPI_CONGRUENT && result != MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMECOMM,"constraint section must have local communicator"); 3526 } 3527 if (mat) { 3528 PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 3529 ierr = MPI_Comm_compare(PETSC_COMM_SELF,PetscObjectComm((PetscObject)mat),&result);CHKERRQ(ierr); 3530 if (result != MPI_CONGRUENT && result != MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMECOMM,"constraint matrix must have local communicator"); 3531 } 3532 ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr); 3533 ierr = PetscSectionDestroy(&dm->defaultConstraintSection);CHKERRQ(ierr); 3534 dm->defaultConstraintSection = section; 3535 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 3536 ierr = MatDestroy(&dm->defaultConstraintMat);CHKERRQ(ierr); 3537 dm->defaultConstraintMat = mat; 3538 PetscFunctionReturn(0); 3539 } 3540 3541 #ifdef PETSC_USE_DEBUG 3542 /* 3543 DMDefaultSectionCheckConsistency - Check the consistentcy of the global and local sections. 3544 3545 Input Parameters: 3546 + dm - The DM 3547 . localSection - PetscSection describing the local data layout 3548 - globalSection - PetscSection describing the global data layout 3549 3550 Level: intermediate 3551 3552 .seealso: DMGetDefaultSF(), DMSetDefaultSF() 3553 */ 3554 static PetscErrorCode DMDefaultSectionCheckConsistency_Internal(DM dm, PetscSection localSection, PetscSection globalSection) 3555 { 3556 MPI_Comm comm; 3557 PetscLayout layout; 3558 const PetscInt *ranges; 3559 PetscInt pStart, pEnd, p, nroots; 3560 PetscMPIInt size, rank; 3561 PetscBool valid = PETSC_TRUE, gvalid; 3562 PetscErrorCode ierr; 3563 3564 PetscFunctionBegin; 3565 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 3566 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3567 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 3568 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 3569 ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr); 3570 ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr); 3571 ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 3572 ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr); 3573 ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr); 3574 ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 3575 ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr); 3576 for (p = pStart; p < pEnd; ++p) { 3577 PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d; 3578 3579 ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr); 3580 ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr); 3581 ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr); 3582 ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr); 3583 ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr); 3584 ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr); 3585 if (!gdof) continue; /* Censored point */ 3586 if ((gdof < 0 ? -(gdof+1) : gdof) != dof) {ierr = PetscSynchronizedPrintf(comm, "[%d]Global dof %d for point %d not equal to local dof %d\n", rank, gdof, p, dof);CHKERRQ(ierr); valid = PETSC_FALSE;} 3587 if (gcdof && (gcdof != cdof)) {ierr = PetscSynchronizedPrintf(comm, "[%d]Global constraints %d for point %d not equal to local constraints %d\n", rank, gcdof, p, cdof);CHKERRQ(ierr); valid = PETSC_FALSE;} 3588 if (gdof < 0) { 3589 gsize = gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof; 3590 for (d = 0; d < gsize; ++d) { 3591 PetscInt offset = -(goff+1) + d, r; 3592 3593 ierr = PetscFindInt(offset,size+1,ranges,&r);CHKERRQ(ierr); 3594 if (r < 0) r = -(r+2); 3595 if ((r < 0) || (r >= size)) {ierr = PetscSynchronizedPrintf(comm, "[%d]Point %d mapped to invalid process %d (%d, %d)\n", rank, p, r, gdof, goff);CHKERRQ(ierr); valid = PETSC_FALSE;break;} 3596 } 3597 } 3598 } 3599 ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 3600 ierr = PetscSynchronizedFlush(comm, NULL);CHKERRQ(ierr); 3601 ierr = MPIU_Allreduce(&valid, &gvalid, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRQ(ierr); 3602 if (!gvalid) { 3603 ierr = DMView(dm, NULL);CHKERRQ(ierr); 3604 SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Inconsistent local and global sections"); 3605 } 3606 PetscFunctionReturn(0); 3607 } 3608 #endif 3609 3610 /*@ 3611 DMGetDefaultGlobalSection - Get the PetscSection encoding the global data layout for the DM. 3612 3613 Collective on DM 3614 3615 Input Parameter: 3616 . dm - The DM 3617 3618 Output Parameter: 3619 . section - The PetscSection 3620 3621 Level: intermediate 3622 3623 Note: This gets a borrowed reference, so the user should not destroy this PetscSection. 3624 3625 .seealso: DMSetDefaultSection(), DMGetDefaultSection() 3626 @*/ 3627 PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *section) 3628 { 3629 PetscErrorCode ierr; 3630 3631 PetscFunctionBegin; 3632 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3633 PetscValidPointer(section, 2); 3634 if (!dm->defaultGlobalSection) { 3635 PetscSection s; 3636 3637 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 3638 if (!s) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM must have a default PetscSection in order to create a global PetscSection"); 3639 if (!dm->sf) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM must have a default PetscSF in order to create a global PetscSection"); 3640 ierr = PetscSectionCreateGlobalSection(s, dm->sf, PETSC_FALSE, PETSC_FALSE, &dm->defaultGlobalSection);CHKERRQ(ierr); 3641 ierr = PetscLayoutDestroy(&dm->map);CHKERRQ(ierr); 3642 ierr = PetscSectionGetValueLayout(PetscObjectComm((PetscObject)dm), dm->defaultGlobalSection, &dm->map);CHKERRQ(ierr); 3643 ierr = PetscSectionViewFromOptions(dm->defaultGlobalSection, NULL, "-global_section_view");CHKERRQ(ierr); 3644 } 3645 *section = dm->defaultGlobalSection; 3646 PetscFunctionReturn(0); 3647 } 3648 3649 /*@ 3650 DMSetDefaultGlobalSection - Set the PetscSection encoding the global data layout for the DM. 3651 3652 Input Parameters: 3653 + dm - The DM 3654 - section - The PetscSection, or NULL 3655 3656 Level: intermediate 3657 3658 Note: Any existing Section will be destroyed 3659 3660 .seealso: DMGetDefaultGlobalSection(), DMSetDefaultSection() 3661 @*/ 3662 PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection section) 3663 { 3664 PetscErrorCode ierr; 3665 3666 PetscFunctionBegin; 3667 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3668 if (section) PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2); 3669 ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr); 3670 ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr); 3671 dm->defaultGlobalSection = section; 3672 #ifdef PETSC_USE_DEBUG 3673 if (section) {ierr = DMDefaultSectionCheckConsistency_Internal(dm, dm->defaultSection, section);CHKERRQ(ierr);} 3674 #endif 3675 PetscFunctionReturn(0); 3676 } 3677 3678 /*@ 3679 DMGetDefaultSF - Get the PetscSF encoding the parallel dof overlap for the DM. If it has not been set, 3680 it is created from the default PetscSection layouts in the DM. 3681 3682 Input Parameter: 3683 . dm - The DM 3684 3685 Output Parameter: 3686 . sf - The PetscSF 3687 3688 Level: intermediate 3689 3690 Note: This gets a borrowed reference, so the user should not destroy this PetscSF. 3691 3692 .seealso: DMSetDefaultSF(), DMCreateDefaultSF() 3693 @*/ 3694 PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *sf) 3695 { 3696 PetscInt nroots; 3697 PetscErrorCode ierr; 3698 3699 PetscFunctionBegin; 3700 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3701 PetscValidPointer(sf, 2); 3702 ierr = PetscSFGetGraph(dm->defaultSF, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 3703 if (nroots < 0) { 3704 PetscSection section, gSection; 3705 3706 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 3707 if (section) { 3708 ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr); 3709 ierr = DMCreateDefaultSF(dm, section, gSection);CHKERRQ(ierr); 3710 } else { 3711 *sf = NULL; 3712 PetscFunctionReturn(0); 3713 } 3714 } 3715 *sf = dm->defaultSF; 3716 PetscFunctionReturn(0); 3717 } 3718 3719 /*@ 3720 DMSetDefaultSF - Set the PetscSF encoding the parallel dof overlap for the DM 3721 3722 Input Parameters: 3723 + dm - The DM 3724 - sf - The PetscSF 3725 3726 Level: intermediate 3727 3728 Note: Any previous SF is destroyed 3729 3730 .seealso: DMGetDefaultSF(), DMCreateDefaultSF() 3731 @*/ 3732 PetscErrorCode DMSetDefaultSF(DM dm, PetscSF sf) 3733 { 3734 PetscErrorCode ierr; 3735 3736 PetscFunctionBegin; 3737 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3738 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 3739 ierr = PetscSFDestroy(&dm->defaultSF);CHKERRQ(ierr); 3740 dm->defaultSF = sf; 3741 PetscFunctionReturn(0); 3742 } 3743 3744 /*@C 3745 DMCreateDefaultSF - Create the PetscSF encoding the parallel dof overlap for the DM based upon the PetscSections 3746 describing the data layout. 3747 3748 Input Parameters: 3749 + dm - The DM 3750 . localSection - PetscSection describing the local data layout 3751 - globalSection - PetscSection describing the global data layout 3752 3753 Level: intermediate 3754 3755 .seealso: DMGetDefaultSF(), DMSetDefaultSF() 3756 @*/ 3757 PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection localSection, PetscSection globalSection) 3758 { 3759 MPI_Comm comm; 3760 PetscLayout layout; 3761 const PetscInt *ranges; 3762 PetscInt *local; 3763 PetscSFNode *remote; 3764 PetscInt pStart, pEnd, p, nroots, nleaves = 0, l; 3765 PetscMPIInt size, rank; 3766 PetscErrorCode ierr; 3767 3768 PetscFunctionBegin; 3769 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 3770 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3771 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 3772 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 3773 ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr); 3774 ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr); 3775 ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 3776 ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr); 3777 ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr); 3778 ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 3779 ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr); 3780 for (p = pStart; p < pEnd; ++p) { 3781 PetscInt gdof, gcdof; 3782 3783 ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr); 3784 ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr); 3785 if (gcdof > (gdof < 0 ? -(gdof+1) : gdof)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d has %d constraints > %d dof", p, gcdof, (gdof < 0 ? -(gdof+1) : gdof)); 3786 nleaves += gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof; 3787 } 3788 ierr = PetscMalloc1(nleaves, &local);CHKERRQ(ierr); 3789 ierr = PetscMalloc1(nleaves, &remote);CHKERRQ(ierr); 3790 for (p = pStart, l = 0; p < pEnd; ++p) { 3791 const PetscInt *cind; 3792 PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d, c; 3793 3794 ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr); 3795 ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr); 3796 ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr); 3797 ierr = PetscSectionGetConstraintIndices(localSection, p, &cind);CHKERRQ(ierr); 3798 ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr); 3799 ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr); 3800 ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr); 3801 if (!gdof) continue; /* Censored point */ 3802 gsize = gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof; 3803 if (gsize != dof-cdof) { 3804 if (gsize != dof) SETERRQ4(comm, PETSC_ERR_ARG_WRONG, "Global dof %d for point %d is neither the constrained size %d, nor the unconstrained %d", gsize, p, dof-cdof, dof); 3805 cdof = 0; /* Ignore constraints */ 3806 } 3807 for (d = 0, c = 0; d < dof; ++d) { 3808 if ((c < cdof) && (cind[c] == d)) {++c; continue;} 3809 local[l+d-c] = off+d; 3810 } 3811 if (gdof < 0) { 3812 for (d = 0; d < gsize; ++d, ++l) { 3813 PetscInt offset = -(goff+1) + d, r; 3814 3815 ierr = PetscFindInt(offset,size+1,ranges,&r);CHKERRQ(ierr); 3816 if (r < 0) r = -(r+2); 3817 if ((r < 0) || (r >= size)) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d mapped to invalid process %d (%d, %d)", p, r, gdof, goff); 3818 remote[l].rank = r; 3819 remote[l].index = offset - ranges[r]; 3820 } 3821 } else { 3822 for (d = 0; d < gsize; ++d, ++l) { 3823 remote[l].rank = rank; 3824 remote[l].index = goff+d - ranges[rank]; 3825 } 3826 } 3827 } 3828 if (l != nleaves) SETERRQ2(comm, PETSC_ERR_PLIB, "Iteration error, l %d != nleaves %d", l, nleaves); 3829 ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 3830 ierr = PetscSFSetGraph(dm->defaultSF, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr); 3831 PetscFunctionReturn(0); 3832 } 3833 3834 /*@ 3835 DMGetPointSF - Get the PetscSF encoding the parallel section point overlap for the DM. 3836 3837 Input Parameter: 3838 . dm - The DM 3839 3840 Output Parameter: 3841 . sf - The PetscSF 3842 3843 Level: intermediate 3844 3845 Note: This gets a borrowed reference, so the user should not destroy this PetscSF. 3846 3847 .seealso: DMSetPointSF(), DMGetDefaultSF(), DMSetDefaultSF(), DMCreateDefaultSF() 3848 @*/ 3849 PetscErrorCode DMGetPointSF(DM dm, PetscSF *sf) 3850 { 3851 PetscFunctionBegin; 3852 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3853 PetscValidPointer(sf, 2); 3854 *sf = dm->sf; 3855 PetscFunctionReturn(0); 3856 } 3857 3858 /*@ 3859 DMSetPointSF - Set the PetscSF encoding the parallel section point overlap for the DM. 3860 3861 Input Parameters: 3862 + dm - The DM 3863 - sf - The PetscSF 3864 3865 Level: intermediate 3866 3867 .seealso: DMGetPointSF(), DMGetDefaultSF(), DMSetDefaultSF(), DMCreateDefaultSF() 3868 @*/ 3869 PetscErrorCode DMSetPointSF(DM dm, PetscSF sf) 3870 { 3871 PetscErrorCode ierr; 3872 3873 PetscFunctionBegin; 3874 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3875 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 3876 ierr = PetscSFDestroy(&dm->sf);CHKERRQ(ierr); 3877 ierr = PetscObjectReference((PetscObject) sf);CHKERRQ(ierr); 3878 dm->sf = sf; 3879 PetscFunctionReturn(0); 3880 } 3881 3882 /*@ 3883 DMGetDS - Get the PetscDS 3884 3885 Input Parameter: 3886 . dm - The DM 3887 3888 Output Parameter: 3889 . prob - The PetscDS 3890 3891 Level: developer 3892 3893 .seealso: DMSetDS() 3894 @*/ 3895 PetscErrorCode DMGetDS(DM dm, PetscDS *prob) 3896 { 3897 PetscFunctionBegin; 3898 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3899 PetscValidPointer(prob, 2); 3900 *prob = dm->prob; 3901 PetscFunctionReturn(0); 3902 } 3903 3904 /*@ 3905 DMSetDS - Set the PetscDS 3906 3907 Input Parameters: 3908 + dm - The DM 3909 - prob - The PetscDS 3910 3911 Level: developer 3912 3913 .seealso: DMGetDS() 3914 @*/ 3915 PetscErrorCode DMSetDS(DM dm, PetscDS prob) 3916 { 3917 PetscErrorCode ierr; 3918 3919 PetscFunctionBegin; 3920 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3921 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2); 3922 ierr = PetscObjectReference((PetscObject) prob);CHKERRQ(ierr); 3923 ierr = PetscDSDestroy(&dm->prob);CHKERRQ(ierr); 3924 dm->prob = prob; 3925 PetscFunctionReturn(0); 3926 } 3927 3928 PetscErrorCode DMGetNumFields(DM dm, PetscInt *numFields) 3929 { 3930 PetscErrorCode ierr; 3931 3932 PetscFunctionBegin; 3933 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3934 ierr = PetscDSGetNumFields(dm->prob, numFields);CHKERRQ(ierr); 3935 PetscFunctionReturn(0); 3936 } 3937 3938 PetscErrorCode DMSetNumFields(DM dm, PetscInt numFields) 3939 { 3940 PetscInt Nf, f; 3941 PetscErrorCode ierr; 3942 3943 PetscFunctionBegin; 3944 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3945 ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr); 3946 for (f = Nf; f < numFields; ++f) { 3947 PetscContainer obj; 3948 3949 ierr = PetscContainerCreate(PetscObjectComm((PetscObject) dm), &obj);CHKERRQ(ierr); 3950 ierr = PetscDSSetDiscretization(dm->prob, f, (PetscObject) obj);CHKERRQ(ierr); 3951 ierr = PetscContainerDestroy(&obj);CHKERRQ(ierr); 3952 } 3953 PetscFunctionReturn(0); 3954 } 3955 3956 /*@ 3957 DMGetField - Return the discretization object for a given DM field 3958 3959 Not collective 3960 3961 Input Parameters: 3962 + dm - The DM 3963 - f - The field number 3964 3965 Output Parameter: 3966 . field - The discretization object 3967 3968 Level: developer 3969 3970 .seealso: DMSetField() 3971 @*/ 3972 PetscErrorCode DMGetField(DM dm, PetscInt f, PetscObject *field) 3973 { 3974 PetscErrorCode ierr; 3975 3976 PetscFunctionBegin; 3977 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3978 ierr = PetscDSGetDiscretization(dm->prob, f, field);CHKERRQ(ierr); 3979 PetscFunctionReturn(0); 3980 } 3981 3982 /*@ 3983 DMSetField - Set the discretization object for a given DM field 3984 3985 Logically collective on DM 3986 3987 Input Parameters: 3988 + dm - The DM 3989 . f - The field number 3990 - field - The discretization object 3991 3992 Level: developer 3993 3994 .seealso: DMGetField() 3995 @*/ 3996 PetscErrorCode DMSetField(DM dm, PetscInt f, PetscObject field) 3997 { 3998 PetscErrorCode ierr; 3999 4000 PetscFunctionBegin; 4001 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4002 ierr = PetscDSSetDiscretization(dm->prob, f, field);CHKERRQ(ierr); 4003 PetscFunctionReturn(0); 4004 } 4005 4006 PetscErrorCode DMRestrictHook_Coordinates(DM dm,DM dmc,void *ctx) 4007 { 4008 DM dm_coord,dmc_coord; 4009 PetscErrorCode ierr; 4010 Vec coords,ccoords; 4011 Mat inject; 4012 PetscFunctionBegin; 4013 ierr = DMGetCoordinateDM(dm,&dm_coord);CHKERRQ(ierr); 4014 ierr = DMGetCoordinateDM(dmc,&dmc_coord);CHKERRQ(ierr); 4015 ierr = DMGetCoordinates(dm,&coords);CHKERRQ(ierr); 4016 ierr = DMGetCoordinates(dmc,&ccoords);CHKERRQ(ierr); 4017 if (coords && !ccoords) { 4018 ierr = DMCreateGlobalVector(dmc_coord,&ccoords);CHKERRQ(ierr); 4019 ierr = PetscObjectSetName((PetscObject)ccoords,"coordinates");CHKERRQ(ierr); 4020 ierr = DMCreateInjection(dmc_coord,dm_coord,&inject);CHKERRQ(ierr); 4021 ierr = MatRestrict(inject,coords,ccoords);CHKERRQ(ierr); 4022 ierr = MatDestroy(&inject);CHKERRQ(ierr); 4023 ierr = DMSetCoordinates(dmc,ccoords);CHKERRQ(ierr); 4024 ierr = VecDestroy(&ccoords);CHKERRQ(ierr); 4025 } 4026 PetscFunctionReturn(0); 4027 } 4028 4029 static PetscErrorCode DMSubDomainHook_Coordinates(DM dm,DM subdm,void *ctx) 4030 { 4031 DM dm_coord,subdm_coord; 4032 PetscErrorCode ierr; 4033 Vec coords,ccoords,clcoords; 4034 VecScatter *scat_i,*scat_g; 4035 PetscFunctionBegin; 4036 ierr = DMGetCoordinateDM(dm,&dm_coord);CHKERRQ(ierr); 4037 ierr = DMGetCoordinateDM(subdm,&subdm_coord);CHKERRQ(ierr); 4038 ierr = DMGetCoordinates(dm,&coords);CHKERRQ(ierr); 4039 ierr = DMGetCoordinates(subdm,&ccoords);CHKERRQ(ierr); 4040 if (coords && !ccoords) { 4041 ierr = DMCreateGlobalVector(subdm_coord,&ccoords);CHKERRQ(ierr); 4042 ierr = PetscObjectSetName((PetscObject)ccoords,"coordinates");CHKERRQ(ierr); 4043 ierr = DMCreateLocalVector(subdm_coord,&clcoords);CHKERRQ(ierr); 4044 ierr = PetscObjectSetName((PetscObject)clcoords,"coordinates");CHKERRQ(ierr); 4045 ierr = DMCreateDomainDecompositionScatters(dm_coord,1,&subdm_coord,NULL,&scat_i,&scat_g);CHKERRQ(ierr); 4046 ierr = VecScatterBegin(scat_i[0],coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4047 ierr = VecScatterBegin(scat_g[0],coords,clcoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4048 ierr = VecScatterEnd(scat_i[0],coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4049 ierr = VecScatterEnd(scat_g[0],coords,clcoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4050 ierr = DMSetCoordinates(subdm,ccoords);CHKERRQ(ierr); 4051 ierr = DMSetCoordinatesLocal(subdm,clcoords);CHKERRQ(ierr); 4052 ierr = VecScatterDestroy(&scat_i[0]);CHKERRQ(ierr); 4053 ierr = VecScatterDestroy(&scat_g[0]);CHKERRQ(ierr); 4054 ierr = VecDestroy(&ccoords);CHKERRQ(ierr); 4055 ierr = VecDestroy(&clcoords);CHKERRQ(ierr); 4056 ierr = PetscFree(scat_i);CHKERRQ(ierr); 4057 ierr = PetscFree(scat_g);CHKERRQ(ierr); 4058 } 4059 PetscFunctionReturn(0); 4060 } 4061 4062 /*@ 4063 DMGetDimension - Return the topological dimension of the DM 4064 4065 Not collective 4066 4067 Input Parameter: 4068 . dm - The DM 4069 4070 Output Parameter: 4071 . dim - The topological dimension 4072 4073 Level: beginner 4074 4075 .seealso: DMSetDimension(), DMCreate() 4076 @*/ 4077 PetscErrorCode DMGetDimension(DM dm, PetscInt *dim) 4078 { 4079 PetscFunctionBegin; 4080 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4081 PetscValidPointer(dim, 2); 4082 *dim = dm->dim; 4083 PetscFunctionReturn(0); 4084 } 4085 4086 /*@ 4087 DMSetDimension - Set the topological dimension of the DM 4088 4089 Collective on dm 4090 4091 Input Parameters: 4092 + dm - The DM 4093 - dim - The topological dimension 4094 4095 Level: beginner 4096 4097 .seealso: DMGetDimension(), DMCreate() 4098 @*/ 4099 PetscErrorCode DMSetDimension(DM dm, PetscInt dim) 4100 { 4101 PetscFunctionBegin; 4102 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4103 PetscValidLogicalCollectiveInt(dm, dim, 2); 4104 dm->dim = dim; 4105 PetscFunctionReturn(0); 4106 } 4107 4108 /*@ 4109 DMGetDimPoints - Get the half-open interval for all points of a given dimension 4110 4111 Collective on DM 4112 4113 Input Parameters: 4114 + dm - the DM 4115 - dim - the dimension 4116 4117 Output Parameters: 4118 + pStart - The first point of the given dimension 4119 . pEnd - The first point following points of the given dimension 4120 4121 Note: 4122 The points are vertices in the Hasse diagram encoding the topology. This is explained in 4123 http://arxiv.org/abs/0908.4427. If not points exist of this dimension in the storage scheme, 4124 then the interval is empty. 4125 4126 Level: intermediate 4127 4128 .keywords: point, Hasse Diagram, dimension 4129 .seealso: DMPLEX, DMPlexGetDepthStratum(), DMPlexGetHeightStratum() 4130 @*/ 4131 PetscErrorCode DMGetDimPoints(DM dm, PetscInt dim, PetscInt *pStart, PetscInt *pEnd) 4132 { 4133 PetscInt d; 4134 PetscErrorCode ierr; 4135 4136 PetscFunctionBegin; 4137 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4138 ierr = DMGetDimension(dm, &d);CHKERRQ(ierr); 4139 if ((dim < 0) || (dim > d)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %d 1", dim, d); 4140 ierr = (*dm->ops->getdimpoints)(dm, dim, pStart, pEnd);CHKERRQ(ierr); 4141 PetscFunctionReturn(0); 4142 } 4143 4144 /*@ 4145 DMSetCoordinates - Sets into the DM a global vector that holds the coordinates 4146 4147 Collective on DM 4148 4149 Input Parameters: 4150 + dm - the DM 4151 - c - coordinate vector 4152 4153 Note: 4154 The coordinates do include those for ghost points, which are in the local vector 4155 4156 Level: intermediate 4157 4158 .keywords: distributed array, get, corners, nodes, local indices, coordinates 4159 .seealso: DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLoca(), DMGetCoordinateDM() 4160 @*/ 4161 PetscErrorCode DMSetCoordinates(DM dm, Vec c) 4162 { 4163 PetscErrorCode ierr; 4164 4165 PetscFunctionBegin; 4166 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4167 PetscValidHeaderSpecific(c,VEC_CLASSID,2); 4168 ierr = PetscObjectReference((PetscObject) c);CHKERRQ(ierr); 4169 ierr = VecDestroy(&dm->coordinates);CHKERRQ(ierr); 4170 dm->coordinates = c; 4171 ierr = VecDestroy(&dm->coordinatesLocal);CHKERRQ(ierr); 4172 ierr = DMCoarsenHookAdd(dm,DMRestrictHook_Coordinates,NULL,NULL);CHKERRQ(ierr); 4173 ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_Coordinates,NULL,NULL);CHKERRQ(ierr); 4174 PetscFunctionReturn(0); 4175 } 4176 4177 /*@ 4178 DMSetCoordinatesLocal - Sets into the DM a local vector that holds the coordinates 4179 4180 Collective on DM 4181 4182 Input Parameters: 4183 + dm - the DM 4184 - c - coordinate vector 4185 4186 Note: 4187 The coordinates of ghost points can be set using DMSetCoordinates() 4188 followed by DMGetCoordinatesLocal(). This is intended to enable the 4189 setting of ghost coordinates outside of the domain. 4190 4191 Level: intermediate 4192 4193 .keywords: distributed array, get, corners, nodes, local indices, coordinates 4194 .seealso: DMGetCoordinatesLocal(), DMSetCoordinates(), DMGetCoordinates(), DMGetCoordinateDM() 4195 @*/ 4196 PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c) 4197 { 4198 PetscErrorCode ierr; 4199 4200 PetscFunctionBegin; 4201 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4202 PetscValidHeaderSpecific(c,VEC_CLASSID,2); 4203 ierr = PetscObjectReference((PetscObject) c);CHKERRQ(ierr); 4204 ierr = VecDestroy(&dm->coordinatesLocal);CHKERRQ(ierr); 4205 4206 dm->coordinatesLocal = c; 4207 4208 ierr = VecDestroy(&dm->coordinates);CHKERRQ(ierr); 4209 PetscFunctionReturn(0); 4210 } 4211 4212 /*@ 4213 DMGetCoordinates - Gets a global vector with the coordinates associated with the DM. 4214 4215 Not Collective 4216 4217 Input Parameter: 4218 . dm - the DM 4219 4220 Output Parameter: 4221 . c - global coordinate vector 4222 4223 Note: 4224 This is a borrowed reference, so the user should NOT destroy this vector 4225 4226 Each process has only the local coordinates (does NOT have the ghost coordinates). 4227 4228 For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 4229 and (x_0,y_0,z_0,x_1,y_1,z_1...) 4230 4231 Level: intermediate 4232 4233 .keywords: distributed array, get, corners, nodes, local indices, coordinates 4234 .seealso: DMSetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM() 4235 @*/ 4236 PetscErrorCode DMGetCoordinates(DM dm, Vec *c) 4237 { 4238 PetscErrorCode ierr; 4239 4240 PetscFunctionBegin; 4241 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4242 PetscValidPointer(c,2); 4243 if (!dm->coordinates && dm->coordinatesLocal) { 4244 DM cdm = NULL; 4245 4246 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 4247 ierr = DMCreateGlobalVector(cdm, &dm->coordinates);CHKERRQ(ierr); 4248 ierr = PetscObjectSetName((PetscObject) dm->coordinates, "coordinates");CHKERRQ(ierr); 4249 ierr = DMLocalToGlobalBegin(cdm, dm->coordinatesLocal, INSERT_VALUES, dm->coordinates);CHKERRQ(ierr); 4250 ierr = DMLocalToGlobalEnd(cdm, dm->coordinatesLocal, INSERT_VALUES, dm->coordinates);CHKERRQ(ierr); 4251 } 4252 *c = dm->coordinates; 4253 PetscFunctionReturn(0); 4254 } 4255 4256 /*@ 4257 DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the DM. 4258 4259 Collective on DM 4260 4261 Input Parameter: 4262 . dm - the DM 4263 4264 Output Parameter: 4265 . c - coordinate vector 4266 4267 Note: 4268 This is a borrowed reference, so the user should NOT destroy this vector 4269 4270 Each process has the local and ghost coordinates 4271 4272 For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 4273 and (x_0,y_0,z_0,x_1,y_1,z_1...) 4274 4275 Level: intermediate 4276 4277 .keywords: distributed array, get, corners, nodes, local indices, coordinates 4278 .seealso: DMSetCoordinatesLocal(), DMGetCoordinates(), DMSetCoordinates(), DMGetCoordinateDM() 4279 @*/ 4280 PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c) 4281 { 4282 PetscErrorCode ierr; 4283 4284 PetscFunctionBegin; 4285 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4286 PetscValidPointer(c,2); 4287 if (!dm->coordinatesLocal && dm->coordinates) { 4288 DM cdm = NULL; 4289 4290 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 4291 ierr = DMCreateLocalVector(cdm, &dm->coordinatesLocal);CHKERRQ(ierr); 4292 ierr = PetscObjectSetName((PetscObject) dm->coordinatesLocal, "coordinates");CHKERRQ(ierr); 4293 ierr = DMGlobalToLocalBegin(cdm, dm->coordinates, INSERT_VALUES, dm->coordinatesLocal);CHKERRQ(ierr); 4294 ierr = DMGlobalToLocalEnd(cdm, dm->coordinates, INSERT_VALUES, dm->coordinatesLocal);CHKERRQ(ierr); 4295 } 4296 *c = dm->coordinatesLocal; 4297 PetscFunctionReturn(0); 4298 } 4299 4300 /*@ 4301 DMGetCoordinateDM - Gets the DM that prescribes coordinate layout and scatters between global and local coordinates 4302 4303 Collective on DM 4304 4305 Input Parameter: 4306 . dm - the DM 4307 4308 Output Parameter: 4309 . cdm - coordinate DM 4310 4311 Level: intermediate 4312 4313 .keywords: distributed array, get, corners, nodes, local indices, coordinates 4314 .seealso: DMSetCoordinateDM(), DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal() 4315 @*/ 4316 PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm) 4317 { 4318 PetscErrorCode ierr; 4319 4320 PetscFunctionBegin; 4321 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4322 PetscValidPointer(cdm,2); 4323 if (!dm->coordinateDM) { 4324 if (!dm->ops->createcoordinatedm) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unable to create coordinates for this DM"); 4325 ierr = (*dm->ops->createcoordinatedm)(dm, &dm->coordinateDM);CHKERRQ(ierr); 4326 } 4327 *cdm = dm->coordinateDM; 4328 PetscFunctionReturn(0); 4329 } 4330 4331 /*@ 4332 DMSetCoordinateDM - Sets the DM that prescribes coordinate layout and scatters between global and local coordinates 4333 4334 Logically Collective on DM 4335 4336 Input Parameters: 4337 + dm - the DM 4338 - cdm - coordinate DM 4339 4340 Level: intermediate 4341 4342 .keywords: distributed array, get, corners, nodes, local indices, coordinates 4343 .seealso: DMGetCoordinateDM(), DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal() 4344 @*/ 4345 PetscErrorCode DMSetCoordinateDM(DM dm, DM cdm) 4346 { 4347 PetscErrorCode ierr; 4348 4349 PetscFunctionBegin; 4350 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4351 PetscValidHeaderSpecific(cdm,DM_CLASSID,2); 4352 ierr = PetscObjectReference((PetscObject)cdm);CHKERRQ(ierr); 4353 ierr = DMDestroy(&dm->coordinateDM);CHKERRQ(ierr); 4354 dm->coordinateDM = cdm; 4355 PetscFunctionReturn(0); 4356 } 4357 4358 /*@ 4359 DMGetCoordinateDim - Retrieve the dimension of embedding space for coordinate values. 4360 4361 Not Collective 4362 4363 Input Parameter: 4364 . dm - The DM object 4365 4366 Output Parameter: 4367 . dim - The embedding dimension 4368 4369 Level: intermediate 4370 4371 .keywords: mesh, coordinates 4372 .seealso: DMSetCoordinateDim(), DMGetCoordinateSection(), DMGetCoordinateDM(), DMGetDefaultSection(), DMSetDefaultSection() 4373 @*/ 4374 PetscErrorCode DMGetCoordinateDim(DM dm, PetscInt *dim) 4375 { 4376 PetscFunctionBegin; 4377 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4378 PetscValidPointer(dim, 2); 4379 if (dm->dimEmbed == PETSC_DEFAULT) { 4380 dm->dimEmbed = dm->dim; 4381 } 4382 *dim = dm->dimEmbed; 4383 PetscFunctionReturn(0); 4384 } 4385 4386 /*@ 4387 DMSetCoordinateDim - Set the dimension of the embedding space for coordinate values. 4388 4389 Not Collective 4390 4391 Input Parameters: 4392 + dm - The DM object 4393 - dim - The embedding dimension 4394 4395 Level: intermediate 4396 4397 .keywords: mesh, coordinates 4398 .seealso: DMGetCoordinateDim(), DMSetCoordinateSection(), DMGetCoordinateSection(), DMGetDefaultSection(), DMSetDefaultSection() 4399 @*/ 4400 PetscErrorCode DMSetCoordinateDim(DM dm, PetscInt dim) 4401 { 4402 PetscFunctionBegin; 4403 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4404 dm->dimEmbed = dim; 4405 PetscFunctionReturn(0); 4406 } 4407 4408 /*@ 4409 DMGetCoordinateSection - Retrieve the layout of coordinate values over the mesh. 4410 4411 Not Collective 4412 4413 Input Parameter: 4414 . dm - The DM object 4415 4416 Output Parameter: 4417 . section - The PetscSection object 4418 4419 Level: intermediate 4420 4421 .keywords: mesh, coordinates 4422 .seealso: DMGetCoordinateDM(), DMGetDefaultSection(), DMSetDefaultSection() 4423 @*/ 4424 PetscErrorCode DMGetCoordinateSection(DM dm, PetscSection *section) 4425 { 4426 DM cdm; 4427 PetscErrorCode ierr; 4428 4429 PetscFunctionBegin; 4430 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4431 PetscValidPointer(section, 2); 4432 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 4433 ierr = DMGetDefaultSection(cdm, section);CHKERRQ(ierr); 4434 PetscFunctionReturn(0); 4435 } 4436 4437 /*@ 4438 DMSetCoordinateSection - Set the layout of coordinate values over the mesh. 4439 4440 Not Collective 4441 4442 Input Parameters: 4443 + dm - The DM object 4444 . dim - The embedding dimension, or PETSC_DETERMINE 4445 - section - The PetscSection object 4446 4447 Level: intermediate 4448 4449 .keywords: mesh, coordinates 4450 .seealso: DMGetCoordinateSection(), DMGetDefaultSection(), DMSetDefaultSection() 4451 @*/ 4452 PetscErrorCode DMSetCoordinateSection(DM dm, PetscInt dim, PetscSection section) 4453 { 4454 DM cdm; 4455 PetscErrorCode ierr; 4456 4457 PetscFunctionBegin; 4458 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4459 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,3); 4460 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 4461 ierr = DMSetDefaultSection(cdm, section);CHKERRQ(ierr); 4462 if (dim == PETSC_DETERMINE) { 4463 PetscInt d = PETSC_DEFAULT; 4464 PetscInt pStart, pEnd, vStart, vEnd, v, dd; 4465 4466 ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 4467 ierr = DMGetDimPoints(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 4468 pStart = PetscMax(vStart, pStart); 4469 pEnd = PetscMin(vEnd, pEnd); 4470 for (v = pStart; v < pEnd; ++v) { 4471 ierr = PetscSectionGetDof(section, v, &dd);CHKERRQ(ierr); 4472 if (dd) {d = dd; break;} 4473 } 4474 if (d < 0) d = PETSC_DEFAULT; 4475 ierr = DMSetCoordinateDim(dm, d);CHKERRQ(ierr); 4476 } 4477 PetscFunctionReturn(0); 4478 } 4479 4480 /*@C 4481 DMGetPeriodicity - Get the description of mesh periodicity 4482 4483 Input Parameters: 4484 . dm - The DM object 4485 4486 Output Parameters: 4487 + per - Whether the DM is periodic or not 4488 . maxCell - Over distances greater than this, we can assume a point has crossed over to another sheet, when trying to localize cell coordinates 4489 . L - If we assume the mesh is a torus, this is the length of each coordinate 4490 - bd - This describes the type of periodicity in each topological dimension 4491 4492 Level: developer 4493 4494 .seealso: DMGetPeriodicity() 4495 @*/ 4496 PetscErrorCode DMGetPeriodicity(DM dm, PetscBool *per, const PetscReal **maxCell, const PetscReal **L, const DMBoundaryType **bd) 4497 { 4498 PetscFunctionBegin; 4499 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4500 if (per) *per = dm->periodic; 4501 if (L) *L = dm->L; 4502 if (maxCell) *maxCell = dm->maxCell; 4503 if (bd) *bd = dm->bdtype; 4504 PetscFunctionReturn(0); 4505 } 4506 4507 /*@C 4508 DMSetPeriodicity - Set the description of mesh periodicity 4509 4510 Input Parameters: 4511 + dm - The DM object 4512 . per - Whether the DM is periodic or not. If maxCell is not provided, coordinates need to be localized 4513 . maxCell - Over distances greater than this, we can assume a point has crossed over to another sheet, when trying to localize cell coordinates 4514 . L - If we assume the mesh is a torus, this is the length of each coordinate 4515 - bd - This describes the type of periodicity in each topological dimension 4516 4517 Level: developer 4518 4519 .seealso: DMGetPeriodicity() 4520 @*/ 4521 PetscErrorCode DMSetPeriodicity(DM dm, PetscBool per, const PetscReal maxCell[], const PetscReal L[], const DMBoundaryType bd[]) 4522 { 4523 PetscInt dim, d; 4524 PetscErrorCode ierr; 4525 4526 PetscFunctionBegin; 4527 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4528 PetscValidLogicalCollectiveBool(dm,per,2); 4529 if (maxCell) { 4530 PetscValidPointer(maxCell,3); 4531 PetscValidPointer(L,4); 4532 PetscValidPointer(bd,5); 4533 } 4534 ierr = PetscFree3(dm->L,dm->maxCell,dm->bdtype);CHKERRQ(ierr); 4535 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 4536 if (maxCell) { 4537 ierr = PetscMalloc3(dim,&dm->L,dim,&dm->maxCell,dim,&dm->bdtype);CHKERRQ(ierr); 4538 for (d = 0; d < dim; ++d) {dm->L[d] = L[d]; dm->maxCell[d] = maxCell[d]; dm->bdtype[d] = bd[d];} 4539 dm->periodic = PETSC_TRUE; 4540 } else { 4541 dm->periodic = per; 4542 } 4543 PetscFunctionReturn(0); 4544 } 4545 4546 /*@ 4547 DMLocalizeCoordinate - If a mesh is periodic (a torus with lengths L_i, some of which can be infinite), project the coordinate onto [0, L_i) in each dimension. 4548 4549 Input Parameters: 4550 + dm - The DM 4551 . in - The input coordinate point (dim numbers) 4552 - endpoint - Include the endpoint L_i 4553 4554 Output Parameter: 4555 . out - The localized coordinate point 4556 4557 Level: developer 4558 4559 .seealso: DMLocalizeCoordinates(), DMLocalizeAddCoordinate() 4560 @*/ 4561 PetscErrorCode DMLocalizeCoordinate(DM dm, const PetscScalar in[], PetscBool endpoint, PetscScalar out[]) 4562 { 4563 PetscInt dim, d; 4564 PetscErrorCode ierr; 4565 4566 PetscFunctionBegin; 4567 ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 4568 if (!dm->maxCell) { 4569 for (d = 0; d < dim; ++d) out[d] = in[d]; 4570 } else { 4571 if (endpoint) { 4572 for (d = 0; d < dim; ++d) { 4573 if ((PetscAbsReal(PetscRealPart(in[d])/dm->L[d] - PetscFloorReal(PetscRealPart(in[d])/dm->L[d])) < PETSC_SMALL) && (PetscRealPart(in[d])/dm->L[d] > PETSC_SMALL)) { 4574 out[d] = in[d] - dm->L[d]*(PetscFloorReal(PetscRealPart(in[d])/dm->L[d]) - 1); 4575 } else { 4576 out[d] = in[d] - dm->L[d]*PetscFloorReal(PetscRealPart(in[d])/dm->L[d]); 4577 } 4578 } 4579 } else { 4580 for (d = 0; d < dim; ++d) { 4581 out[d] = in[d] - dm->L[d]*PetscFloorReal(PetscRealPart(in[d])/dm->L[d]); 4582 } 4583 } 4584 } 4585 PetscFunctionReturn(0); 4586 } 4587 4588 /* 4589 DMLocalizeCoordinate_Internal - If a mesh is periodic, and the input point is far from the anchor, pick the coordinate sheet of the torus which moves it closer. 4590 4591 Input Parameters: 4592 + dm - The DM 4593 . dim - The spatial dimension 4594 . anchor - The anchor point, the input point can be no more than maxCell away from it 4595 - in - The input coordinate point (dim numbers) 4596 4597 Output Parameter: 4598 . out - The localized coordinate point 4599 4600 Level: developer 4601 4602 Note: This is meant to get a set of coordinates close to each other, as in a cell. The anchor is usually the one of the vertices on a containing cell 4603 4604 .seealso: DMLocalizeCoordinates(), DMLocalizeAddCoordinate() 4605 */ 4606 PetscErrorCode DMLocalizeCoordinate_Internal(DM dm, PetscInt dim, const PetscScalar anchor[], const PetscScalar in[], PetscScalar out[]) 4607 { 4608 PetscInt d; 4609 4610 PetscFunctionBegin; 4611 if (!dm->maxCell) { 4612 for (d = 0; d < dim; ++d) out[d] = in[d]; 4613 } else { 4614 for (d = 0; d < dim; ++d) { 4615 if ((dm->bdtype[d] != DM_BOUNDARY_NONE) && (PetscAbsScalar(anchor[d] - in[d]) > dm->maxCell[d])) { 4616 out[d] = PetscRealPart(anchor[d]) > PetscRealPart(in[d]) ? dm->L[d] + in[d] : in[d] - dm->L[d]; 4617 } else { 4618 out[d] = in[d]; 4619 } 4620 } 4621 } 4622 PetscFunctionReturn(0); 4623 } 4624 PetscErrorCode DMLocalizeCoordinateReal_Internal(DM dm, PetscInt dim, const PetscReal anchor[], const PetscReal in[], PetscReal out[]) 4625 { 4626 PetscInt d; 4627 4628 PetscFunctionBegin; 4629 if (!dm->maxCell) { 4630 for (d = 0; d < dim; ++d) out[d] = in[d]; 4631 } else { 4632 for (d = 0; d < dim; ++d) { 4633 if ((dm->bdtype[d] != DM_BOUNDARY_NONE) && (PetscAbsReal(anchor[d] - in[d]) > dm->maxCell[d])) { 4634 out[d] = anchor[d] > in[d] ? dm->L[d] + in[d] : in[d] - dm->L[d]; 4635 } else { 4636 out[d] = in[d]; 4637 } 4638 } 4639 } 4640 PetscFunctionReturn(0); 4641 } 4642 4643 /* 4644 DMLocalizeAddCoordinate_Internal - If a mesh is periodic, and the input point is far from the anchor, pick the coordinate sheet of the torus which moves it closer. 4645 4646 Input Parameters: 4647 + dm - The DM 4648 . dim - The spatial dimension 4649 . anchor - The anchor point, the input point can be no more than maxCell away from it 4650 . in - The input coordinate delta (dim numbers) 4651 - out - The input coordinate point (dim numbers) 4652 4653 Output Parameter: 4654 . out - The localized coordinate in + out 4655 4656 Level: developer 4657 4658 Note: This is meant to get a set of coordinates close to each other, as in a cell. The anchor is usually the one of the vertices on a containing cell 4659 4660 .seealso: DMLocalizeCoordinates(), DMLocalizeCoordinate() 4661 */ 4662 PetscErrorCode DMLocalizeAddCoordinate_Internal(DM dm, PetscInt dim, const PetscScalar anchor[], const PetscScalar in[], PetscScalar out[]) 4663 { 4664 PetscInt d; 4665 4666 PetscFunctionBegin; 4667 if (!dm->maxCell) { 4668 for (d = 0; d < dim; ++d) out[d] += in[d]; 4669 } else { 4670 for (d = 0; d < dim; ++d) { 4671 if ((dm->bdtype[d] != DM_BOUNDARY_NONE) && (PetscAbsScalar(anchor[d] - in[d]) > dm->maxCell[d])) { 4672 out[d] += PetscRealPart(anchor[d]) > PetscRealPart(in[d]) ? dm->L[d] + in[d] : in[d] - dm->L[d]; 4673 } else { 4674 out[d] += in[d]; 4675 } 4676 } 4677 } 4678 PetscFunctionReturn(0); 4679 } 4680 4681 /*@ 4682 DMGetCoordinatesLocalized - Check if the DM coordinates have been localized for cells 4683 4684 Input Parameter: 4685 . dm - The DM 4686 4687 Output Parameter: 4688 areLocalized - True if localized 4689 4690 Level: developer 4691 4692 .seealso: DMLocalizeCoordinates() 4693 @*/ 4694 PetscErrorCode DMGetCoordinatesLocalized(DM dm,PetscBool *areLocalized) 4695 { 4696 DM cdm; 4697 PetscSection coordSection; 4698 PetscInt cStart, cEnd, c, sStart, sEnd, dof; 4699 PetscBool alreadyLocalized, alreadyLocalizedGlobal; 4700 PetscErrorCode ierr; 4701 4702 PetscFunctionBegin; 4703 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4704 if (!dm->periodic) { 4705 *areLocalized = PETSC_FALSE; 4706 PetscFunctionReturn(0); 4707 } 4708 /* We need some generic way of refering to cells/vertices */ 4709 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 4710 { 4711 PetscBool isplex; 4712 4713 ierr = PetscObjectTypeCompare((PetscObject) cdm, DMPLEX, &isplex);CHKERRQ(ierr); 4714 if (isplex) { 4715 ierr = DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd);CHKERRQ(ierr); 4716 } else SETERRQ(PetscObjectComm((PetscObject) cdm), PETSC_ERR_ARG_WRONG, "Coordinate localization requires a DMPLEX coordinate DM"); 4717 } 4718 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 4719 ierr = PetscSectionGetChart(coordSection,&sStart,&sEnd);CHKERRQ(ierr); 4720 alreadyLocalized = alreadyLocalizedGlobal = PETSC_FALSE; 4721 for (c = cStart; c < cEnd; ++c) { 4722 if (c < sStart || c >= sEnd) { 4723 alreadyLocalized = PETSC_FALSE; 4724 break; 4725 } 4726 ierr = PetscSectionGetDof(coordSection, c, &dof);CHKERRQ(ierr); 4727 if (dof) { 4728 alreadyLocalized = PETSC_TRUE; 4729 break; 4730 } 4731 } 4732 ierr = MPI_Allreduce(&alreadyLocalized,&alreadyLocalizedGlobal,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 4733 *areLocalized = alreadyLocalizedGlobal; 4734 PetscFunctionReturn(0); 4735 } 4736 4737 4738 /*@ 4739 DMLocalizeCoordinates - If a mesh is periodic, create local coordinates for each cell 4740 4741 Input Parameter: 4742 . dm - The DM 4743 4744 Level: developer 4745 4746 .seealso: DMLocalizeCoordinate(), DMLocalizeAddCoordinate() 4747 @*/ 4748 PetscErrorCode DMLocalizeCoordinates(DM dm) 4749 { 4750 DM cdm; 4751 PetscSection coordSection, cSection; 4752 Vec coordinates, cVec; 4753 PetscScalar *coords, *coords2, *anchor, *localized; 4754 PetscInt Nc, vStart, vEnd, v, sStart, sEnd, newStart = PETSC_MAX_INT, newEnd = PETSC_MIN_INT, dof, d, off, off2, bs, coordSize; 4755 PetscBool alreadyLocalized, alreadyLocalizedGlobal; 4756 PetscInt maxHeight = 0, h; 4757 PetscInt *pStart = NULL, *pEnd = NULL; 4758 PetscErrorCode ierr; 4759 4760 PetscFunctionBegin; 4761 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4762 if (!dm->periodic) PetscFunctionReturn(0); 4763 /* We need some generic way of refering to cells/vertices */ 4764 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 4765 { 4766 PetscBool isplex; 4767 4768 ierr = PetscObjectTypeCompare((PetscObject) cdm, DMPLEX, &isplex);CHKERRQ(ierr); 4769 if (isplex) { 4770 ierr = DMPlexGetDepthStratum(cdm, 0, &vStart, &vEnd);CHKERRQ(ierr); 4771 ierr = DMPlexGetMaxProjectionHeight(cdm,&maxHeight);CHKERRQ(ierr); 4772 ierr = DMGetWorkArray(dm,2*(maxHeight + 1),PETSC_INT,&pStart);CHKERRQ(ierr); 4773 pEnd = &pStart[maxHeight + 1]; 4774 newStart = vStart; 4775 newEnd = vEnd; 4776 for (h = 0; h <= maxHeight; h++) { 4777 ierr = DMPlexGetHeightStratum(cdm, h, &pStart[h], &pEnd[h]);CHKERRQ(ierr); 4778 newStart = PetscMin(newStart,pStart[h]); 4779 newEnd = PetscMax(newEnd,pEnd[h]); 4780 } 4781 } else SETERRQ(PetscObjectComm((PetscObject) cdm), PETSC_ERR_ARG_WRONG, "Coordinate localization requires a DMPLEX coordinate DM"); 4782 } 4783 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 4784 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 4785 ierr = VecGetBlockSize(coordinates, &bs);CHKERRQ(ierr); 4786 ierr = PetscSectionGetChart(coordSection,&sStart,&sEnd);CHKERRQ(ierr); 4787 4788 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &cSection);CHKERRQ(ierr); 4789 ierr = PetscSectionSetNumFields(cSection, 1);CHKERRQ(ierr); 4790 ierr = PetscSectionGetFieldComponents(coordSection, 0, &Nc);CHKERRQ(ierr); 4791 ierr = PetscSectionSetFieldComponents(cSection, 0, Nc);CHKERRQ(ierr); 4792 ierr = PetscSectionSetChart(cSection, newStart, newEnd);CHKERRQ(ierr); 4793 4794 ierr = DMGetWorkArray(dm, 2 * bs, PETSC_SCALAR, &anchor);CHKERRQ(ierr); 4795 localized = &anchor[bs]; 4796 alreadyLocalized = alreadyLocalizedGlobal = PETSC_TRUE; 4797 for (h = 0; h <= maxHeight; h++) { 4798 PetscInt cStart = pStart[h], cEnd = pEnd[h], c; 4799 4800 for (c = cStart; c < cEnd; ++c) { 4801 PetscScalar *cellCoords = NULL; 4802 PetscInt b; 4803 4804 if (c < sStart || c >= sEnd) alreadyLocalized = PETSC_FALSE; 4805 ierr = DMPlexVecGetClosure(cdm, coordSection, coordinates, c, &dof, &cellCoords);CHKERRQ(ierr); 4806 for (b = 0; b < bs; ++b) anchor[b] = cellCoords[b]; 4807 for (d = 0; d < dof/bs; ++d) { 4808 ierr = DMLocalizeCoordinate_Internal(dm, bs, anchor, &cellCoords[d*bs], localized);CHKERRQ(ierr); 4809 for (b = 0; b < bs; b++) { 4810 if (cellCoords[d*bs + b] != localized[b]) break; 4811 } 4812 if (b < bs) break; 4813 } 4814 if (d < dof/bs) { 4815 if (c >= sStart && c < sEnd) { 4816 PetscInt cdof; 4817 4818 ierr = PetscSectionGetDof(coordSection, c, &cdof);CHKERRQ(ierr); 4819 if (cdof != dof) alreadyLocalized = PETSC_FALSE; 4820 } 4821 ierr = PetscSectionSetDof(cSection, c, dof);CHKERRQ(ierr); 4822 ierr = PetscSectionSetFieldDof(cSection, c, 0, dof);CHKERRQ(ierr); 4823 } 4824 ierr = DMPlexVecRestoreClosure(cdm, coordSection, coordinates, c, &dof, &cellCoords);CHKERRQ(ierr); 4825 } 4826 } 4827 ierr = MPI_Allreduce(&alreadyLocalized,&alreadyLocalizedGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 4828 if (alreadyLocalizedGlobal) { 4829 ierr = DMRestoreWorkArray(dm, 2 * bs, PETSC_SCALAR, &anchor);CHKERRQ(ierr); 4830 ierr = PetscSectionDestroy(&cSection);CHKERRQ(ierr); 4831 ierr = DMRestoreWorkArray(dm,2*(maxHeight + 1),PETSC_INT,&pStart);CHKERRQ(ierr); 4832 PetscFunctionReturn(0); 4833 } 4834 for (v = vStart; v < vEnd; ++v) { 4835 ierr = PetscSectionGetDof(coordSection, v, &dof);CHKERRQ(ierr); 4836 ierr = PetscSectionSetDof(cSection, v, dof);CHKERRQ(ierr); 4837 ierr = PetscSectionSetFieldDof(cSection, v, 0, dof);CHKERRQ(ierr); 4838 } 4839 ierr = PetscSectionSetUp(cSection);CHKERRQ(ierr); 4840 ierr = PetscSectionGetStorageSize(cSection, &coordSize);CHKERRQ(ierr); 4841 ierr = VecCreate(PETSC_COMM_SELF, &cVec);CHKERRQ(ierr); 4842 ierr = PetscObjectSetName((PetscObject)cVec,"coordinates");CHKERRQ(ierr); 4843 ierr = VecSetBlockSize(cVec, bs);CHKERRQ(ierr); 4844 ierr = VecSetSizes(cVec, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 4845 ierr = VecSetType(cVec, VECSTANDARD);CHKERRQ(ierr); 4846 ierr = VecGetArrayRead(coordinates, (const PetscScalar**)&coords);CHKERRQ(ierr); 4847 ierr = VecGetArray(cVec, &coords2);CHKERRQ(ierr); 4848 for (v = vStart; v < vEnd; ++v) { 4849 ierr = PetscSectionGetDof(coordSection, v, &dof);CHKERRQ(ierr); 4850 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 4851 ierr = PetscSectionGetOffset(cSection, v, &off2);CHKERRQ(ierr); 4852 for (d = 0; d < dof; ++d) coords2[off2+d] = coords[off+d]; 4853 } 4854 for (h = 0; h <= maxHeight; h++) { 4855 PetscInt cStart = pStart[h], cEnd = pEnd[h], c; 4856 4857 for (c = cStart; c < cEnd; ++c) { 4858 PetscScalar *cellCoords = NULL; 4859 PetscInt b, cdof; 4860 4861 ierr = PetscSectionGetDof(cSection,c,&cdof);CHKERRQ(ierr); 4862 if (!cdof) continue; 4863 ierr = DMPlexVecGetClosure(cdm, coordSection, coordinates, c, &dof, &cellCoords);CHKERRQ(ierr); 4864 ierr = PetscSectionGetOffset(cSection, c, &off2);CHKERRQ(ierr); 4865 for (b = 0; b < bs; ++b) anchor[b] = cellCoords[b]; 4866 for (d = 0; d < dof/bs; ++d) {ierr = DMLocalizeCoordinate_Internal(dm, bs, anchor, &cellCoords[d*bs], &coords2[off2+d*bs]);CHKERRQ(ierr);} 4867 ierr = DMPlexVecRestoreClosure(cdm, coordSection, coordinates, c, &dof, &cellCoords);CHKERRQ(ierr); 4868 } 4869 } 4870 ierr = DMRestoreWorkArray(dm, 2 * bs, PETSC_SCALAR, &anchor);CHKERRQ(ierr); 4871 ierr = DMRestoreWorkArray(dm,2*(maxHeight + 1),PETSC_INT,&pStart);CHKERRQ(ierr); 4872 ierr = VecRestoreArrayRead(coordinates, (const PetscScalar**)&coords);CHKERRQ(ierr); 4873 ierr = VecRestoreArray(cVec, &coords2);CHKERRQ(ierr); 4874 ierr = DMSetCoordinateSection(dm, PETSC_DETERMINE, cSection);CHKERRQ(ierr); 4875 ierr = DMSetCoordinatesLocal(dm, cVec);CHKERRQ(ierr); 4876 ierr = VecDestroy(&cVec);CHKERRQ(ierr); 4877 ierr = PetscSectionDestroy(&cSection);CHKERRQ(ierr); 4878 PetscFunctionReturn(0); 4879 } 4880 4881 /*@ 4882 DMLocatePoints - Locate the points in v in the mesh and return a PetscSF of the containing cells 4883 4884 Collective on Vec v (see explanation below) 4885 4886 Input Parameters: 4887 + dm - The DM 4888 . v - The Vec of points 4889 . ltype - The type of point location, e.g. DM_POINTLOCATION_NONE or DM_POINTLOCATION_NEAREST 4890 - cells - Points to either NULL, or a PetscSF with guesses for which cells contain each point. 4891 4892 Output Parameter: 4893 + v - The Vec of points, which now contains the nearest mesh points to the given points if DM_POINTLOCATION_NEAREST is used 4894 - cells - The PetscSF containing the ranks and local indices of the containing points. 4895 4896 4897 Level: developer 4898 4899 Notes: 4900 To do a search of the local cells of the mesh, v should have PETSC_COMM_SELF as its communicator. 4901 To do a search of all the cells in the distributed mesh, v should have the same communicator as dm. 4902 4903 If *cellSF is NULL on input, a PetscSF will be created. 4904 If *cellSF is not NULL on input, it should point to an existing PetscSF, whose graph will be used as initial guesses. 4905 4906 An array that maps each point to its containing cell can be obtained with 4907 4908 $ const PetscSFNode *cells; 4909 $ PetscInt nFound; 4910 $ const PetscSFNode *found; 4911 $ 4912 $ PetscSFGetGraph(cells,NULL,&nFound,&found,&cells); 4913 4914 Where cells[i].rank is the rank of the cell containing point found[i] (or i if found == NULL), and cells[i].index is 4915 the index of the cell in its rank's local numbering. 4916 4917 .keywords: point location, mesh 4918 .seealso: DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMPointLocationType 4919 @*/ 4920 PetscErrorCode DMLocatePoints(DM dm, Vec v, DMPointLocationType ltype, PetscSF *cellSF) 4921 { 4922 PetscErrorCode ierr; 4923 4924 PetscFunctionBegin; 4925 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4926 PetscValidHeaderSpecific(v,VEC_CLASSID,2); 4927 PetscValidPointer(cellSF,4); 4928 if (*cellSF) { 4929 PetscMPIInt result; 4930 4931 PetscValidHeaderSpecific(*cellSF,PETSCSF_CLASSID,4); 4932 ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)v),PetscObjectComm((PetscObject)*cellSF),&result);CHKERRQ(ierr); 4933 if (result != MPI_IDENT && result != MPI_CONGRUENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"cellSF must have a communicator congruent to v's"); 4934 } else { 4935 ierr = PetscSFCreate(PetscObjectComm((PetscObject)v),cellSF);CHKERRQ(ierr); 4936 } 4937 ierr = PetscLogEventBegin(DM_LocatePoints,dm,0,0,0);CHKERRQ(ierr); 4938 if (dm->ops->locatepoints) { 4939 ierr = (*dm->ops->locatepoints)(dm,v,ltype,*cellSF);CHKERRQ(ierr); 4940 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Point location not available for this DM"); 4941 ierr = PetscLogEventEnd(DM_LocatePoints,dm,0,0,0);CHKERRQ(ierr); 4942 PetscFunctionReturn(0); 4943 } 4944 4945 /*@ 4946 DMGetOutputDM - Retrieve the DM associated with the layout for output 4947 4948 Input Parameter: 4949 . dm - The original DM 4950 4951 Output Parameter: 4952 . odm - The DM which provides the layout for output 4953 4954 Level: intermediate 4955 4956 .seealso: VecView(), DMGetDefaultGlobalSection() 4957 @*/ 4958 PetscErrorCode DMGetOutputDM(DM dm, DM *odm) 4959 { 4960 PetscSection section; 4961 PetscBool hasConstraints, ghasConstraints; 4962 PetscErrorCode ierr; 4963 4964 PetscFunctionBegin; 4965 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4966 PetscValidPointer(odm,2); 4967 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 4968 ierr = PetscSectionHasConstraints(section, &hasConstraints);CHKERRQ(ierr); 4969 ierr = MPI_Allreduce(&hasConstraints, &ghasConstraints, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr); 4970 if (!ghasConstraints) { 4971 *odm = dm; 4972 PetscFunctionReturn(0); 4973 } 4974 if (!dm->dmBC) { 4975 PetscDS ds; 4976 PetscSection newSection, gsection; 4977 PetscSF sf; 4978 4979 ierr = DMClone(dm, &dm->dmBC);CHKERRQ(ierr); 4980 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 4981 ierr = DMSetDS(dm->dmBC, ds);CHKERRQ(ierr); 4982 ierr = PetscSectionClone(section, &newSection);CHKERRQ(ierr); 4983 ierr = DMSetDefaultSection(dm->dmBC, newSection);CHKERRQ(ierr); 4984 ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr); 4985 ierr = DMGetPointSF(dm->dmBC, &sf);CHKERRQ(ierr); 4986 ierr = PetscSectionCreateGlobalSection(section, sf, PETSC_TRUE, PETSC_FALSE, &gsection);CHKERRQ(ierr); 4987 ierr = DMSetDefaultGlobalSection(dm->dmBC, gsection);CHKERRQ(ierr); 4988 ierr = PetscSectionDestroy(&gsection);CHKERRQ(ierr); 4989 } 4990 *odm = dm->dmBC; 4991 PetscFunctionReturn(0); 4992 } 4993 4994 /*@ 4995 DMGetOutputSequenceNumber - Retrieve the sequence number/value for output 4996 4997 Input Parameter: 4998 . dm - The original DM 4999 5000 Output Parameters: 5001 + num - The output sequence number 5002 - val - The output sequence value 5003 5004 Level: intermediate 5005 5006 Note: This is intended for output that should appear in sequence, for instance 5007 a set of timesteps in an HDF5 file, or a set of realizations of a stochastic system. 5008 5009 .seealso: VecView() 5010 @*/ 5011 PetscErrorCode DMGetOutputSequenceNumber(DM dm, PetscInt *num, PetscReal *val) 5012 { 5013 PetscFunctionBegin; 5014 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5015 if (num) {PetscValidPointer(num,2); *num = dm->outputSequenceNum;} 5016 if (val) {PetscValidPointer(val,3);*val = dm->outputSequenceVal;} 5017 PetscFunctionReturn(0); 5018 } 5019 5020 /*@ 5021 DMSetOutputSequenceNumber - Set the sequence number/value for output 5022 5023 Input Parameters: 5024 + dm - The original DM 5025 . num - The output sequence number 5026 - val - The output sequence value 5027 5028 Level: intermediate 5029 5030 Note: This is intended for output that should appear in sequence, for instance 5031 a set of timesteps in an HDF5 file, or a set of realizations of a stochastic system. 5032 5033 .seealso: VecView() 5034 @*/ 5035 PetscErrorCode DMSetOutputSequenceNumber(DM dm, PetscInt num, PetscReal val) 5036 { 5037 PetscFunctionBegin; 5038 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5039 dm->outputSequenceNum = num; 5040 dm->outputSequenceVal = val; 5041 PetscFunctionReturn(0); 5042 } 5043 5044 /*@C 5045 DMOutputSequenceLoad - Retrieve the sequence value from a Viewer 5046 5047 Input Parameters: 5048 + dm - The original DM 5049 . name - The sequence name 5050 - num - The output sequence number 5051 5052 Output Parameter: 5053 . val - The output sequence value 5054 5055 Level: intermediate 5056 5057 Note: This is intended for output that should appear in sequence, for instance 5058 a set of timesteps in an HDF5 file, or a set of realizations of a stochastic system. 5059 5060 .seealso: DMGetOutputSequenceNumber(), DMSetOutputSequenceNumber(), VecView() 5061 @*/ 5062 PetscErrorCode DMOutputSequenceLoad(DM dm, PetscViewer viewer, const char *name, PetscInt num, PetscReal *val) 5063 { 5064 PetscBool ishdf5; 5065 PetscErrorCode ierr; 5066 5067 PetscFunctionBegin; 5068 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5069 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 5070 PetscValidPointer(val,4); 5071 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);CHKERRQ(ierr); 5072 if (ishdf5) { 5073 #if defined(PETSC_HAVE_HDF5) 5074 PetscScalar value; 5075 5076 ierr = DMSequenceLoad_HDF5_Internal(dm, name, num, &value, viewer);CHKERRQ(ierr); 5077 *val = PetscRealPart(value); 5078 #endif 5079 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerHDF5Open()"); 5080 PetscFunctionReturn(0); 5081 } 5082 5083 /*@ 5084 DMGetUseNatural - Get the flag for creating a mapping to the natural order on distribution 5085 5086 Not collective 5087 5088 Input Parameter: 5089 . dm - The DM 5090 5091 Output Parameter: 5092 . useNatural - The flag to build the mapping to a natural order during distribution 5093 5094 Level: beginner 5095 5096 .seealso: DMSetUseNatural(), DMCreate() 5097 @*/ 5098 PetscErrorCode DMGetUseNatural(DM dm, PetscBool *useNatural) 5099 { 5100 PetscFunctionBegin; 5101 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5102 PetscValidPointer(useNatural, 2); 5103 *useNatural = dm->useNatural; 5104 PetscFunctionReturn(0); 5105 } 5106 5107 /*@ 5108 DMSetUseNatural - Set the flag for creating a mapping to the natural order on distribution 5109 5110 Collective on dm 5111 5112 Input Parameters: 5113 + dm - The DM 5114 - useNatural - The flag to build the mapping to a natural order during distribution 5115 5116 Level: beginner 5117 5118 .seealso: DMGetUseNatural(), DMCreate() 5119 @*/ 5120 PetscErrorCode DMSetUseNatural(DM dm, PetscBool useNatural) 5121 { 5122 PetscFunctionBegin; 5123 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5124 PetscValidLogicalCollectiveInt(dm, useNatural, 2); 5125 dm->useNatural = useNatural; 5126 PetscFunctionReturn(0); 5127 } 5128 5129 5130 /*@C 5131 DMCreateLabel - Create a label of the given name if it does not already exist 5132 5133 Not Collective 5134 5135 Input Parameters: 5136 + dm - The DM object 5137 - name - The label name 5138 5139 Level: intermediate 5140 5141 .keywords: mesh 5142 .seealso: DMLabelCreate(), DMHasLabel(), DMGetLabelValue(), DMSetLabelValue(), DMGetStratumIS() 5143 @*/ 5144 PetscErrorCode DMCreateLabel(DM dm, const char name[]) 5145 { 5146 DMLabelLink next = dm->labels->next; 5147 PetscBool flg = PETSC_FALSE; 5148 PetscErrorCode ierr; 5149 5150 PetscFunctionBegin; 5151 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5152 PetscValidCharPointer(name, 2); 5153 while (next) { 5154 ierr = PetscStrcmp(name, next->label->name, &flg);CHKERRQ(ierr); 5155 if (flg) break; 5156 next = next->next; 5157 } 5158 if (!flg) { 5159 DMLabelLink tmpLabel; 5160 5161 ierr = PetscCalloc1(1, &tmpLabel);CHKERRQ(ierr); 5162 ierr = DMLabelCreate(name, &tmpLabel->label);CHKERRQ(ierr); 5163 tmpLabel->output = PETSC_TRUE; 5164 tmpLabel->next = dm->labels->next; 5165 dm->labels->next = tmpLabel; 5166 } 5167 PetscFunctionReturn(0); 5168 } 5169 5170 /*@C 5171 DMGetLabelValue - Get the value in a Sieve Label for the given point, with 0 as the default 5172 5173 Not Collective 5174 5175 Input Parameters: 5176 + dm - The DM object 5177 . name - The label name 5178 - point - The mesh point 5179 5180 Output Parameter: 5181 . value - The label value for this point, or -1 if the point is not in the label 5182 5183 Level: beginner 5184 5185 .keywords: mesh 5186 .seealso: DMLabelGetValue(), DMSetLabelValue(), DMGetStratumIS() 5187 @*/ 5188 PetscErrorCode DMGetLabelValue(DM dm, const char name[], PetscInt point, PetscInt *value) 5189 { 5190 DMLabel label; 5191 PetscErrorCode ierr; 5192 5193 PetscFunctionBegin; 5194 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5195 PetscValidCharPointer(name, 2); 5196 ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 5197 if (!label) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No label named %s was found", name); 5198 ierr = DMLabelGetValue(label, point, value);CHKERRQ(ierr); 5199 PetscFunctionReturn(0); 5200 } 5201 5202 /*@C 5203 DMSetLabelValue - Add a point to a Sieve Label with given value 5204 5205 Not Collective 5206 5207 Input Parameters: 5208 + dm - The DM object 5209 . name - The label name 5210 . point - The mesh point 5211 - value - The label value for this point 5212 5213 Output Parameter: 5214 5215 Level: beginner 5216 5217 .keywords: mesh 5218 .seealso: DMLabelSetValue(), DMGetStratumIS(), DMClearLabelValue() 5219 @*/ 5220 PetscErrorCode DMSetLabelValue(DM dm, const char name[], PetscInt point, PetscInt value) 5221 { 5222 DMLabel label; 5223 PetscErrorCode ierr; 5224 5225 PetscFunctionBegin; 5226 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5227 PetscValidCharPointer(name, 2); 5228 ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 5229 if (!label) { 5230 ierr = DMCreateLabel(dm, name);CHKERRQ(ierr); 5231 ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 5232 } 5233 ierr = DMLabelSetValue(label, point, value);CHKERRQ(ierr); 5234 PetscFunctionReturn(0); 5235 } 5236 5237 /*@C 5238 DMClearLabelValue - Remove a point from a Sieve Label with given value 5239 5240 Not Collective 5241 5242 Input Parameters: 5243 + dm - The DM object 5244 . name - The label name 5245 . point - The mesh point 5246 - value - The label value for this point 5247 5248 Output Parameter: 5249 5250 Level: beginner 5251 5252 .keywords: mesh 5253 .seealso: DMLabelClearValue(), DMSetLabelValue(), DMGetStratumIS() 5254 @*/ 5255 PetscErrorCode DMClearLabelValue(DM dm, const char name[], PetscInt point, PetscInt value) 5256 { 5257 DMLabel label; 5258 PetscErrorCode ierr; 5259 5260 PetscFunctionBegin; 5261 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5262 PetscValidCharPointer(name, 2); 5263 ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 5264 if (!label) PetscFunctionReturn(0); 5265 ierr = DMLabelClearValue(label, point, value);CHKERRQ(ierr); 5266 PetscFunctionReturn(0); 5267 } 5268 5269 /*@C 5270 DMGetLabelSize - Get the number of different integer ids in a Label 5271 5272 Not Collective 5273 5274 Input Parameters: 5275 + dm - The DM object 5276 - name - The label name 5277 5278 Output Parameter: 5279 . size - The number of different integer ids, or 0 if the label does not exist 5280 5281 Level: beginner 5282 5283 .keywords: mesh 5284 .seealso: DMLabeGetNumValues(), DMSetLabelValue() 5285 @*/ 5286 PetscErrorCode DMGetLabelSize(DM dm, const char name[], PetscInt *size) 5287 { 5288 DMLabel label; 5289 PetscErrorCode ierr; 5290 5291 PetscFunctionBegin; 5292 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5293 PetscValidCharPointer(name, 2); 5294 PetscValidPointer(size, 3); 5295 ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 5296 *size = 0; 5297 if (!label) PetscFunctionReturn(0); 5298 ierr = DMLabelGetNumValues(label, size);CHKERRQ(ierr); 5299 PetscFunctionReturn(0); 5300 } 5301 5302 /*@C 5303 DMGetLabelIdIS - Get the integer ids in a label 5304 5305 Not Collective 5306 5307 Input Parameters: 5308 + mesh - The DM object 5309 - name - The label name 5310 5311 Output Parameter: 5312 . ids - The integer ids, or NULL if the label does not exist 5313 5314 Level: beginner 5315 5316 .keywords: mesh 5317 .seealso: DMLabelGetValueIS(), DMGetLabelSize() 5318 @*/ 5319 PetscErrorCode DMGetLabelIdIS(DM dm, const char name[], IS *ids) 5320 { 5321 DMLabel label; 5322 PetscErrorCode ierr; 5323 5324 PetscFunctionBegin; 5325 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5326 PetscValidCharPointer(name, 2); 5327 PetscValidPointer(ids, 3); 5328 ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 5329 *ids = NULL; 5330 if (!label) PetscFunctionReturn(0); 5331 ierr = DMLabelGetValueIS(label, ids);CHKERRQ(ierr); 5332 PetscFunctionReturn(0); 5333 } 5334 5335 /*@C 5336 DMGetStratumSize - Get the number of points in a label stratum 5337 5338 Not Collective 5339 5340 Input Parameters: 5341 + dm - The DM object 5342 . name - The label name 5343 - value - The stratum value 5344 5345 Output Parameter: 5346 . size - The stratum size 5347 5348 Level: beginner 5349 5350 .keywords: mesh 5351 .seealso: DMLabelGetStratumSize(), DMGetLabelSize(), DMGetLabelIds() 5352 @*/ 5353 PetscErrorCode DMGetStratumSize(DM dm, const char name[], PetscInt value, PetscInt *size) 5354 { 5355 DMLabel label; 5356 PetscErrorCode ierr; 5357 5358 PetscFunctionBegin; 5359 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5360 PetscValidCharPointer(name, 2); 5361 PetscValidPointer(size, 4); 5362 ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 5363 *size = 0; 5364 if (!label) PetscFunctionReturn(0); 5365 ierr = DMLabelGetStratumSize(label, value, size);CHKERRQ(ierr); 5366 PetscFunctionReturn(0); 5367 } 5368 5369 /*@C 5370 DMGetStratumIS - Get the points in a label stratum 5371 5372 Not Collective 5373 5374 Input Parameters: 5375 + dm - The DM object 5376 . name - The label name 5377 - value - The stratum value 5378 5379 Output Parameter: 5380 . points - The stratum points, or NULL if the label does not exist or does not have that value 5381 5382 Level: beginner 5383 5384 .keywords: mesh 5385 .seealso: DMLabelGetStratumIS(), DMGetStratumSize() 5386 @*/ 5387 PetscErrorCode DMGetStratumIS(DM dm, const char name[], PetscInt value, IS *points) 5388 { 5389 DMLabel label; 5390 PetscErrorCode ierr; 5391 5392 PetscFunctionBegin; 5393 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5394 PetscValidCharPointer(name, 2); 5395 PetscValidPointer(points, 4); 5396 ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 5397 *points = NULL; 5398 if (!label) PetscFunctionReturn(0); 5399 ierr = DMLabelGetStratumIS(label, value, points);CHKERRQ(ierr); 5400 PetscFunctionReturn(0); 5401 } 5402 5403 /*@C 5404 DMGetStratumIS - Set the points in a label stratum 5405 5406 Not Collective 5407 5408 Input Parameters: 5409 + dm - The DM object 5410 . name - The label name 5411 . value - The stratum value 5412 - points - The stratum points 5413 5414 Level: beginner 5415 5416 .keywords: mesh 5417 .seealso: DMLabelSetStratumIS(), DMGetStratumSize() 5418 @*/ 5419 PetscErrorCode DMSetStratumIS(DM dm, const char name[], PetscInt value, IS points) 5420 { 5421 DMLabel label; 5422 PetscErrorCode ierr; 5423 5424 PetscFunctionBegin; 5425 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5426 PetscValidCharPointer(name, 2); 5427 PetscValidPointer(points, 4); 5428 ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 5429 if (!label) PetscFunctionReturn(0); 5430 ierr = DMLabelSetStratumIS(label, value, points);CHKERRQ(ierr); 5431 PetscFunctionReturn(0); 5432 } 5433 5434 /*@C 5435 DMClearLabelStratum - Remove all points from a stratum from a Sieve Label 5436 5437 Not Collective 5438 5439 Input Parameters: 5440 + dm - The DM object 5441 . name - The label name 5442 - value - The label value for this point 5443 5444 Output Parameter: 5445 5446 Level: beginner 5447 5448 .keywords: mesh 5449 .seealso: DMLabelClearStratum(), DMSetLabelValue(), DMGetStratumIS(), DMClearLabelValue() 5450 @*/ 5451 PetscErrorCode DMClearLabelStratum(DM dm, const char name[], PetscInt value) 5452 { 5453 DMLabel label; 5454 PetscErrorCode ierr; 5455 5456 PetscFunctionBegin; 5457 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5458 PetscValidCharPointer(name, 2); 5459 ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 5460 if (!label) PetscFunctionReturn(0); 5461 ierr = DMLabelClearStratum(label, value);CHKERRQ(ierr); 5462 PetscFunctionReturn(0); 5463 } 5464 5465 /*@ 5466 DMGetNumLabels - Return the number of labels defined by the mesh 5467 5468 Not Collective 5469 5470 Input Parameter: 5471 . dm - The DM object 5472 5473 Output Parameter: 5474 . numLabels - the number of Labels 5475 5476 Level: intermediate 5477 5478 .keywords: mesh 5479 .seealso: DMGetLabelValue(), DMSetLabelValue(), DMGetStratumIS() 5480 @*/ 5481 PetscErrorCode DMGetNumLabels(DM dm, PetscInt *numLabels) 5482 { 5483 DMLabelLink next = dm->labels->next; 5484 PetscInt n = 0; 5485 5486 PetscFunctionBegin; 5487 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5488 PetscValidPointer(numLabels, 2); 5489 while (next) {++n; next = next->next;} 5490 *numLabels = n; 5491 PetscFunctionReturn(0); 5492 } 5493 5494 /*@C 5495 DMGetLabelName - Return the name of nth label 5496 5497 Not Collective 5498 5499 Input Parameters: 5500 + dm - The DM object 5501 - n - the label number 5502 5503 Output Parameter: 5504 . name - the label name 5505 5506 Level: intermediate 5507 5508 .keywords: mesh 5509 .seealso: DMGetLabelValue(), DMSetLabelValue(), DMGetStratumIS() 5510 @*/ 5511 PetscErrorCode DMGetLabelName(DM dm, PetscInt n, const char **name) 5512 { 5513 DMLabelLink next = dm->labels->next; 5514 PetscInt l = 0; 5515 5516 PetscFunctionBegin; 5517 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5518 PetscValidPointer(name, 3); 5519 while (next) { 5520 if (l == n) { 5521 *name = next->label->name; 5522 PetscFunctionReturn(0); 5523 } 5524 ++l; 5525 next = next->next; 5526 } 5527 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %D does not exist in this DM", n); 5528 } 5529 5530 /*@C 5531 DMHasLabel - Determine whether the mesh has a label of a given name 5532 5533 Not Collective 5534 5535 Input Parameters: 5536 + dm - The DM object 5537 - name - The label name 5538 5539 Output Parameter: 5540 . hasLabel - PETSC_TRUE if the label is present 5541 5542 Level: intermediate 5543 5544 .keywords: mesh 5545 .seealso: DMCreateLabel(), DMGetLabelValue(), DMSetLabelValue(), DMGetStratumIS() 5546 @*/ 5547 PetscErrorCode DMHasLabel(DM dm, const char name[], PetscBool *hasLabel) 5548 { 5549 DMLabelLink next = dm->labels->next; 5550 PetscErrorCode ierr; 5551 5552 PetscFunctionBegin; 5553 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5554 PetscValidCharPointer(name, 2); 5555 PetscValidPointer(hasLabel, 3); 5556 *hasLabel = PETSC_FALSE; 5557 while (next) { 5558 ierr = PetscStrcmp(name, next->label->name, hasLabel);CHKERRQ(ierr); 5559 if (*hasLabel) break; 5560 next = next->next; 5561 } 5562 PetscFunctionReturn(0); 5563 } 5564 5565 /*@C 5566 DMGetLabel - Return the label of a given name, or NULL 5567 5568 Not Collective 5569 5570 Input Parameters: 5571 + dm - The DM object 5572 - name - The label name 5573 5574 Output Parameter: 5575 . label - The DMLabel, or NULL if the label is absent 5576 5577 Level: intermediate 5578 5579 .keywords: mesh 5580 .seealso: DMCreateLabel(), DMHasLabel(), DMGetLabelValue(), DMSetLabelValue(), DMGetStratumIS() 5581 @*/ 5582 PetscErrorCode DMGetLabel(DM dm, const char name[], DMLabel *label) 5583 { 5584 DMLabelLink next = dm->labels->next; 5585 PetscBool hasLabel; 5586 PetscErrorCode ierr; 5587 5588 PetscFunctionBegin; 5589 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5590 PetscValidCharPointer(name, 2); 5591 PetscValidPointer(label, 3); 5592 *label = NULL; 5593 while (next) { 5594 ierr = PetscStrcmp(name, next->label->name, &hasLabel);CHKERRQ(ierr); 5595 if (hasLabel) { 5596 *label = next->label; 5597 break; 5598 } 5599 next = next->next; 5600 } 5601 PetscFunctionReturn(0); 5602 } 5603 5604 /*@C 5605 DMGetLabelByNum - Return the nth label 5606 5607 Not Collective 5608 5609 Input Parameters: 5610 + dm - The DM object 5611 - n - the label number 5612 5613 Output Parameter: 5614 . label - the label 5615 5616 Level: intermediate 5617 5618 .keywords: mesh 5619 .seealso: DMGetLabelValue(), DMSetLabelValue(), DMGetStratumIS() 5620 @*/ 5621 PetscErrorCode DMGetLabelByNum(DM dm, PetscInt n, DMLabel *label) 5622 { 5623 DMLabelLink next = dm->labels->next; 5624 PetscInt l = 0; 5625 5626 PetscFunctionBegin; 5627 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5628 PetscValidPointer(label, 3); 5629 while (next) { 5630 if (l == n) { 5631 *label = next->label; 5632 PetscFunctionReturn(0); 5633 } 5634 ++l; 5635 next = next->next; 5636 } 5637 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %D does not exist in this DM", n); 5638 } 5639 5640 /*@C 5641 DMAddLabel - Add the label to this mesh 5642 5643 Not Collective 5644 5645 Input Parameters: 5646 + dm - The DM object 5647 - label - The DMLabel 5648 5649 Level: developer 5650 5651 .keywords: mesh 5652 .seealso: DMCreateLabel(), DMHasLabel(), DMGetLabelValue(), DMSetLabelValue(), DMGetStratumIS() 5653 @*/ 5654 PetscErrorCode DMAddLabel(DM dm, DMLabel label) 5655 { 5656 DMLabelLink tmpLabel; 5657 PetscBool hasLabel; 5658 PetscErrorCode ierr; 5659 5660 PetscFunctionBegin; 5661 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5662 ierr = DMHasLabel(dm, label->name, &hasLabel);CHKERRQ(ierr); 5663 if (hasLabel) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label %s already exists in this DM", label->name); 5664 ierr = PetscCalloc1(1, &tmpLabel);CHKERRQ(ierr); 5665 tmpLabel->label = label; 5666 tmpLabel->output = PETSC_TRUE; 5667 tmpLabel->next = dm->labels->next; 5668 dm->labels->next = tmpLabel; 5669 PetscFunctionReturn(0); 5670 } 5671 5672 /*@C 5673 DMRemoveLabel - Remove the label from this mesh 5674 5675 Not Collective 5676 5677 Input Parameters: 5678 + dm - The DM object 5679 - name - The label name 5680 5681 Output Parameter: 5682 . label - The DMLabel, or NULL if the label is absent 5683 5684 Level: developer 5685 5686 .keywords: mesh 5687 .seealso: DMCreateLabel(), DMHasLabel(), DMGetLabelValue(), DMSetLabelValue(), DMGetStratumIS() 5688 @*/ 5689 PetscErrorCode DMRemoveLabel(DM dm, const char name[], DMLabel *label) 5690 { 5691 DMLabelLink next = dm->labels->next; 5692 DMLabelLink last = NULL; 5693 PetscBool hasLabel; 5694 PetscErrorCode ierr; 5695 5696 PetscFunctionBegin; 5697 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5698 ierr = DMHasLabel(dm, name, &hasLabel);CHKERRQ(ierr); 5699 *label = NULL; 5700 if (!hasLabel) PetscFunctionReturn(0); 5701 while (next) { 5702 ierr = PetscStrcmp(name, next->label->name, &hasLabel);CHKERRQ(ierr); 5703 if (hasLabel) { 5704 if (last) last->next = next->next; 5705 else dm->labels->next = next->next; 5706 next->next = NULL; 5707 *label = next->label; 5708 ierr = PetscStrcmp(name, "depth", &hasLabel);CHKERRQ(ierr); 5709 if (hasLabel) { 5710 dm->depthLabel = NULL; 5711 } 5712 ierr = PetscFree(next);CHKERRQ(ierr); 5713 break; 5714 } 5715 last = next; 5716 next = next->next; 5717 } 5718 PetscFunctionReturn(0); 5719 } 5720 5721 /*@C 5722 DMGetLabelOutput - Get the output flag for a given label 5723 5724 Not Collective 5725 5726 Input Parameters: 5727 + dm - The DM object 5728 - name - The label name 5729 5730 Output Parameter: 5731 . output - The flag for output 5732 5733 Level: developer 5734 5735 .keywords: mesh 5736 .seealso: DMSetLabelOutput(), DMCreateLabel(), DMHasLabel(), DMGetLabelValue(), DMSetLabelValue(), DMGetStratumIS() 5737 @*/ 5738 PetscErrorCode DMGetLabelOutput(DM dm, const char name[], PetscBool *output) 5739 { 5740 DMLabelLink next = dm->labels->next; 5741 PetscErrorCode ierr; 5742 5743 PetscFunctionBegin; 5744 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5745 PetscValidPointer(name, 2); 5746 PetscValidPointer(output, 3); 5747 while (next) { 5748 PetscBool flg; 5749 5750 ierr = PetscStrcmp(name, next->label->name, &flg);CHKERRQ(ierr); 5751 if (flg) {*output = next->output; PetscFunctionReturn(0);} 5752 next = next->next; 5753 } 5754 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No label named %s was present in this dm", name); 5755 } 5756 5757 /*@C 5758 DMSetLabelOutput - Set the output flag for a given label 5759 5760 Not Collective 5761 5762 Input Parameters: 5763 + dm - The DM object 5764 . name - The label name 5765 - output - The flag for output 5766 5767 Level: developer 5768 5769 .keywords: mesh 5770 .seealso: DMGetLabelOutput(), DMCreateLabel(), DMHasLabel(), DMGetLabelValue(), DMSetLabelValue(), DMGetStratumIS() 5771 @*/ 5772 PetscErrorCode DMSetLabelOutput(DM dm, const char name[], PetscBool output) 5773 { 5774 DMLabelLink next = dm->labels->next; 5775 PetscErrorCode ierr; 5776 5777 PetscFunctionBegin; 5778 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5779 PetscValidPointer(name, 2); 5780 while (next) { 5781 PetscBool flg; 5782 5783 ierr = PetscStrcmp(name, next->label->name, &flg);CHKERRQ(ierr); 5784 if (flg) {next->output = output; PetscFunctionReturn(0);} 5785 next = next->next; 5786 } 5787 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No label named %s was present in this dm", name); 5788 } 5789 5790 5791 /*@ 5792 DMCopyLabels - Copy labels from one mesh to another with a superset of the points 5793 5794 Collective on DM 5795 5796 Input Parameter: 5797 . dmA - The DM object with initial labels 5798 5799 Output Parameter: 5800 . dmB - The DM object with copied labels 5801 5802 Level: intermediate 5803 5804 Note: This is typically used when interpolating or otherwise adding to a mesh 5805 5806 .keywords: mesh 5807 .seealso: DMCopyCoordinates(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection() 5808 @*/ 5809 PetscErrorCode DMCopyLabels(DM dmA, DM dmB) 5810 { 5811 PetscInt numLabels, l; 5812 PetscErrorCode ierr; 5813 5814 PetscFunctionBegin; 5815 if (dmA == dmB) PetscFunctionReturn(0); 5816 ierr = DMGetNumLabels(dmA, &numLabels);CHKERRQ(ierr); 5817 for (l = 0; l < numLabels; ++l) { 5818 DMLabel label, labelNew; 5819 const char *name; 5820 PetscBool flg; 5821 5822 ierr = DMGetLabelName(dmA, l, &name);CHKERRQ(ierr); 5823 ierr = PetscStrcmp(name, "depth", &flg);CHKERRQ(ierr); 5824 if (flg) continue; 5825 ierr = DMGetLabel(dmA, name, &label);CHKERRQ(ierr); 5826 ierr = DMLabelDuplicate(label, &labelNew);CHKERRQ(ierr); 5827 ierr = DMAddLabel(dmB, labelNew);CHKERRQ(ierr); 5828 } 5829 PetscFunctionReturn(0); 5830 } 5831 5832 /*@ 5833 DMGetCoarseDM - Get the coarse mesh from which this was obtained by refinement 5834 5835 Input Parameter: 5836 . dm - The DM object 5837 5838 Output Parameter: 5839 . cdm - The coarse DM 5840 5841 Level: intermediate 5842 5843 .seealso: DMSetCoarseDM() 5844 @*/ 5845 PetscErrorCode DMGetCoarseDM(DM dm, DM *cdm) 5846 { 5847 PetscFunctionBegin; 5848 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5849 PetscValidPointer(cdm, 2); 5850 *cdm = dm->coarseMesh; 5851 PetscFunctionReturn(0); 5852 } 5853 5854 /*@ 5855 DMSetCoarseDM - Set the coarse mesh from which this was obtained by refinement 5856 5857 Input Parameters: 5858 + dm - The DM object 5859 - cdm - The coarse DM 5860 5861 Level: intermediate 5862 5863 .seealso: DMGetCoarseDM() 5864 @*/ 5865 PetscErrorCode DMSetCoarseDM(DM dm, DM cdm) 5866 { 5867 PetscErrorCode ierr; 5868 5869 PetscFunctionBegin; 5870 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5871 if (cdm) PetscValidHeaderSpecific(cdm, DM_CLASSID, 2); 5872 ierr = PetscObjectReference((PetscObject)cdm);CHKERRQ(ierr); 5873 ierr = DMDestroy(&dm->coarseMesh);CHKERRQ(ierr); 5874 dm->coarseMesh = cdm; 5875 PetscFunctionReturn(0); 5876 } 5877 5878 /*@ 5879 DMGetFineDM - Get the fine mesh from which this was obtained by refinement 5880 5881 Input Parameter: 5882 . dm - The DM object 5883 5884 Output Parameter: 5885 . fdm - The fine DM 5886 5887 Level: intermediate 5888 5889 .seealso: DMSetFineDM() 5890 @*/ 5891 PetscErrorCode DMGetFineDM(DM dm, DM *fdm) 5892 { 5893 PetscFunctionBegin; 5894 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5895 PetscValidPointer(fdm, 2); 5896 *fdm = dm->fineMesh; 5897 PetscFunctionReturn(0); 5898 } 5899 5900 /*@ 5901 DMSetFineDM - Set the fine mesh from which this was obtained by refinement 5902 5903 Input Parameters: 5904 + dm - The DM object 5905 - fdm - The fine DM 5906 5907 Level: intermediate 5908 5909 .seealso: DMGetFineDM() 5910 @*/ 5911 PetscErrorCode DMSetFineDM(DM dm, DM fdm) 5912 { 5913 PetscErrorCode ierr; 5914 5915 PetscFunctionBegin; 5916 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5917 if (fdm) PetscValidHeaderSpecific(fdm, DM_CLASSID, 2); 5918 ierr = PetscObjectReference((PetscObject)fdm);CHKERRQ(ierr); 5919 ierr = DMDestroy(&dm->fineMesh);CHKERRQ(ierr); 5920 dm->fineMesh = fdm; 5921 PetscFunctionReturn(0); 5922 } 5923 5924 /*=== DMBoundary code ===*/ 5925 5926 PetscErrorCode DMCopyBoundary(DM dm, DM dmNew) 5927 { 5928 PetscErrorCode ierr; 5929 5930 PetscFunctionBegin; 5931 ierr = PetscDSCopyBoundary(dm->prob,dmNew->prob);CHKERRQ(ierr); 5932 PetscFunctionReturn(0); 5933 } 5934 5935 /*@C 5936 DMAddBoundary - Add a boundary condition to the model 5937 5938 Input Parameters: 5939 + dm - The DM, with a PetscDS that matches the problem being constrained 5940 . type - The type of condition, e.g. DM_BC_ESSENTIAL_ANALYTIC/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann) 5941 . name - The BC name 5942 . labelname - The label defining constrained points 5943 . field - The field to constrain 5944 . numcomps - The number of constrained field components 5945 . comps - An array of constrained component numbers 5946 . bcFunc - A pointwise function giving boundary values 5947 . numids - The number of DMLabel ids for constrained points 5948 . ids - An array of ids for constrained points 5949 - ctx - An optional user context for bcFunc 5950 5951 Options Database Keys: 5952 + -bc_<boundary name> <num> - Overrides the boundary ids 5953 - -bc_<boundary name>_comp <num> - Overrides the boundary components 5954 5955 Level: developer 5956 5957 .seealso: DMGetBoundary() 5958 @*/ 5959 PetscErrorCode DMAddBoundary(DM dm, DMBoundaryConditionType type, const char name[], const char labelname[], PetscInt field, PetscInt numcomps, const PetscInt *comps, void (*bcFunc)(void), PetscInt numids, const PetscInt *ids, void *ctx) 5960 { 5961 PetscErrorCode ierr; 5962 5963 PetscFunctionBegin; 5964 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5965 ierr = PetscDSAddBoundary(dm->prob,type,name,labelname,field,numcomps,comps,bcFunc,numids,ids,ctx);CHKERRQ(ierr); 5966 PetscFunctionReturn(0); 5967 } 5968 5969 /*@ 5970 DMGetNumBoundary - Get the number of registered BC 5971 5972 Input Parameters: 5973 . dm - The mesh object 5974 5975 Output Parameters: 5976 . numBd - The number of BC 5977 5978 Level: intermediate 5979 5980 .seealso: DMAddBoundary(), DMGetBoundary() 5981 @*/ 5982 PetscErrorCode DMGetNumBoundary(DM dm, PetscInt *numBd) 5983 { 5984 PetscErrorCode ierr; 5985 5986 PetscFunctionBegin; 5987 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5988 ierr = PetscDSGetNumBoundary(dm->prob,numBd);CHKERRQ(ierr); 5989 PetscFunctionReturn(0); 5990 } 5991 5992 /*@C 5993 DMGetBoundary - Get a model boundary condition 5994 5995 Input Parameters: 5996 + dm - The mesh object 5997 - bd - The BC number 5998 5999 Output Parameters: 6000 + type - The type of condition, e.g. DM_BC_ESSENTIAL_ANALYTIC/DM_BC_ESSENTIAL_FIELD (Dirichlet), or DM_BC_NATURAL (Neumann) 6001 . name - The BC name 6002 . labelname - The label defining constrained points 6003 . field - The field to constrain 6004 . numcomps - The number of constrained field components 6005 . comps - An array of constrained component numbers 6006 . bcFunc - A pointwise function giving boundary values 6007 . numids - The number of DMLabel ids for constrained points 6008 . ids - An array of ids for constrained points 6009 - ctx - An optional user context for bcFunc 6010 6011 Options Database Keys: 6012 + -bc_<boundary name> <num> - Overrides the boundary ids 6013 - -bc_<boundary name>_comp <num> - Overrides the boundary components 6014 6015 Level: developer 6016 6017 .seealso: DMAddBoundary() 6018 @*/ 6019 PetscErrorCode DMGetBoundary(DM dm, PetscInt bd, DMBoundaryConditionType *type, const char **name, const char **labelname, PetscInt *field, PetscInt *numcomps, const PetscInt **comps, void (**func)(void), PetscInt *numids, const PetscInt **ids, void **ctx) 6020 { 6021 PetscErrorCode ierr; 6022 6023 PetscFunctionBegin; 6024 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6025 ierr = PetscDSGetBoundary(dm->prob,bd,type,name,labelname,field,numcomps,comps,func,numids,ids,ctx);CHKERRQ(ierr); 6026 PetscFunctionReturn(0); 6027 } 6028 6029 static PetscErrorCode DMPopulateBoundary(DM dm) 6030 { 6031 DMBoundary *lastnext; 6032 DSBoundary dsbound; 6033 PetscErrorCode ierr; 6034 6035 PetscFunctionBegin; 6036 dsbound = dm->prob->boundary; 6037 if (dm->boundary) { 6038 DMBoundary next = dm->boundary; 6039 6040 /* quick check to see if the PetscDS has changed */ 6041 if (next->dsboundary == dsbound) PetscFunctionReturn(0); 6042 /* the PetscDS has changed: tear down and rebuild */ 6043 while (next) { 6044 DMBoundary b = next; 6045 6046 next = b->next; 6047 ierr = PetscFree(b);CHKERRQ(ierr); 6048 } 6049 dm->boundary = NULL; 6050 } 6051 6052 lastnext = &(dm->boundary); 6053 while (dsbound) { 6054 DMBoundary dmbound; 6055 6056 ierr = PetscNew(&dmbound);CHKERRQ(ierr); 6057 dmbound->dsboundary = dsbound; 6058 ierr = DMGetLabel(dm, dsbound->labelname, &(dmbound->label));CHKERRQ(ierr); 6059 if (!dmbound->label) PetscInfo2(dm, "DSBoundary %s wants label %s, which is not in this dm.\n",dsbound->name,dsbound->labelname);CHKERRQ(ierr); 6060 /* push on the back instead of the front so that it is in the same order as in the PetscDS */ 6061 *lastnext = dmbound; 6062 lastnext = &(dmbound->next); 6063 dsbound = dsbound->next; 6064 } 6065 PetscFunctionReturn(0); 6066 } 6067 6068 PetscErrorCode DMIsBoundaryPoint(DM dm, PetscInt point, PetscBool *isBd) 6069 { 6070 DMBoundary b; 6071 PetscErrorCode ierr; 6072 6073 PetscFunctionBegin; 6074 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6075 PetscValidPointer(isBd, 3); 6076 *isBd = PETSC_FALSE; 6077 ierr = DMPopulateBoundary(dm);CHKERRQ(ierr); 6078 b = dm->boundary; 6079 while (b && !(*isBd)) { 6080 DMLabel label = b->label; 6081 DSBoundary dsb = b->dsboundary; 6082 6083 if (label) { 6084 PetscInt i; 6085 6086 for (i = 0; i < dsb->numids && !(*isBd); ++i) { 6087 ierr = DMLabelStratumHasPoint(label, dsb->ids[i], point, isBd);CHKERRQ(ierr); 6088 } 6089 } 6090 b = b->next; 6091 } 6092 PetscFunctionReturn(0); 6093 } 6094 6095 /*@C 6096 DMProjectFunction - This projects the given function into the function space provided. 6097 6098 Input Parameters: 6099 + dm - The DM 6100 . time - The time 6101 . funcs - The coordinate functions to evaluate, one per field 6102 . ctxs - Optional array of contexts to pass to each coordinate function. ctxs itself may be null. 6103 - mode - The insertion mode for values 6104 6105 Output Parameter: 6106 . X - vector 6107 6108 Calling sequence of func: 6109 $ func(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx); 6110 6111 + dim - The spatial dimension 6112 . x - The coordinates 6113 . Nf - The number of fields 6114 . u - The output field values 6115 - ctx - optional user-defined function context 6116 6117 Level: developer 6118 6119 .seealso: DMComputeL2Diff() 6120 @*/ 6121 PetscErrorCode DMProjectFunction(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X) 6122 { 6123 Vec localX; 6124 PetscErrorCode ierr; 6125 6126 PetscFunctionBegin; 6127 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6128 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 6129 ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, mode, localX);CHKERRQ(ierr); 6130 ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 6131 ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 6132 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 6133 PetscFunctionReturn(0); 6134 } 6135 6136 PetscErrorCode DMProjectFunctionLocal(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 6137 { 6138 PetscErrorCode ierr; 6139 6140 PetscFunctionBegin; 6141 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6142 PetscValidHeaderSpecific(localX,VEC_CLASSID,5); 6143 if (!dm->ops->projectfunctionlocal) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM type %s does not implemnt DMProjectFunctionLocal",((PetscObject)dm)->type_name); 6144 ierr = (dm->ops->projectfunctionlocal) (dm, time, funcs, ctxs, mode, localX);CHKERRQ(ierr); 6145 PetscFunctionReturn(0); 6146 } 6147 6148 PetscErrorCode DMProjectFunctionLabelLocal(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Nc, const PetscInt comps[], PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 6149 { 6150 PetscErrorCode ierr; 6151 6152 PetscFunctionBegin; 6153 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6154 PetscValidHeaderSpecific(localX,VEC_CLASSID,5); 6155 if (!dm->ops->projectfunctionlabellocal) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM type %s does not implemnt DMProjectFunctionLabelLocal",((PetscObject)dm)->type_name); 6156 ierr = (dm->ops->projectfunctionlabellocal) (dm, time, label, numIds, ids, Nc, comps, funcs, ctxs, mode, localX);CHKERRQ(ierr); 6157 PetscFunctionReturn(0); 6158 } 6159 6160 PetscErrorCode DMProjectFieldLocal(DM dm, PetscReal time, Vec localU, 6161 void (**funcs)(PetscInt, PetscInt, PetscInt, 6162 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6163 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6164 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 6165 InsertMode mode, Vec localX) 6166 { 6167 PetscErrorCode ierr; 6168 6169 PetscFunctionBegin; 6170 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6171 PetscValidHeaderSpecific(localU,VEC_CLASSID,3); 6172 PetscValidHeaderSpecific(localX,VEC_CLASSID,6); 6173 if (!dm->ops->projectfieldlocal) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM type %s does not implemnt DMProjectFieldLocal",((PetscObject)dm)->type_name); 6174 ierr = (dm->ops->projectfieldlocal) (dm, time, localU, funcs, mode, localX);CHKERRQ(ierr); 6175 PetscFunctionReturn(0); 6176 } 6177 6178 PetscErrorCode DMProjectFieldLabelLocal(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Nc, const PetscInt comps[], Vec localU, 6179 void (**funcs)(PetscInt, PetscInt, PetscInt, 6180 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6181 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 6182 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 6183 InsertMode mode, Vec localX) 6184 { 6185 PetscErrorCode ierr; 6186 6187 PetscFunctionBegin; 6188 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6189 PetscValidHeaderSpecific(localU,VEC_CLASSID,6); 6190 PetscValidHeaderSpecific(localX,VEC_CLASSID,9); 6191 if (!dm->ops->projectfieldlabellocal) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM type %s does not implemnt DMProjectFieldLocal",((PetscObject)dm)->type_name); 6192 ierr = (dm->ops->projectfieldlabellocal)(dm, time, label, numIds, ids, Nc, comps, localU, funcs, mode, localX);CHKERRQ(ierr); 6193 PetscFunctionReturn(0); 6194 } 6195 6196 /*@C 6197 DMComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 6198 6199 Input Parameters: 6200 + dm - The DM 6201 . time - The time 6202 . funcs - The functions to evaluate for each field component 6203 . ctxs - Optional array of contexts to pass to each function, or NULL. 6204 - X - The coefficient vector u_h 6205 6206 Output Parameter: 6207 . diff - The diff ||u - u_h||_2 6208 6209 Level: developer 6210 6211 .seealso: DMProjectFunction(), DMComputeL2FieldDiff(), DMComputeL2GradientDiff() 6212 @*/ 6213 PetscErrorCode DMComputeL2Diff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 6214 { 6215 PetscErrorCode ierr; 6216 6217 PetscFunctionBegin; 6218 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6219 PetscValidHeaderSpecific(X,VEC_CLASSID,5); 6220 if (!dm->ops->computel2diff) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM type %s does not implemnt DMComputeL2Diff",((PetscObject)dm)->type_name); 6221 ierr = (dm->ops->computel2diff)(dm,time,funcs,ctxs,X,diff);CHKERRQ(ierr); 6222 PetscFunctionReturn(0); 6223 } 6224 6225 /*@C 6226 DMComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h. 6227 6228 Input Parameters: 6229 + dm - The DM 6230 , time - The time 6231 . funcs - The gradient functions to evaluate for each field component 6232 . ctxs - Optional array of contexts to pass to each function, or NULL. 6233 . X - The coefficient vector u_h 6234 - n - The vector to project along 6235 6236 Output Parameter: 6237 . diff - The diff ||(grad u - grad u_h) . n||_2 6238 6239 Level: developer 6240 6241 .seealso: DMProjectFunction(), DMComputeL2Diff() 6242 @*/ 6243 PetscErrorCode DMComputeL2GradientDiff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 6244 { 6245 PetscErrorCode ierr; 6246 6247 PetscFunctionBegin; 6248 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6249 PetscValidHeaderSpecific(X,VEC_CLASSID,5); 6250 if (!dm->ops->computel2gradientdiff) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM type %s does not implement DMComputeL2GradientDiff",((PetscObject)dm)->type_name); 6251 ierr = (dm->ops->computel2gradientdiff)(dm,time,funcs,ctxs,X,n,diff);CHKERRQ(ierr); 6252 PetscFunctionReturn(0); 6253 } 6254 6255 /*@C 6256 DMComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components. 6257 6258 Input Parameters: 6259 + dm - The DM 6260 . time - The time 6261 . funcs - The functions to evaluate for each field component 6262 . ctxs - Optional array of contexts to pass to each function, or NULL. 6263 - X - The coefficient vector u_h 6264 6265 Output Parameter: 6266 . diff - The array of differences, ||u^f - u^f_h||_2 6267 6268 Level: developer 6269 6270 .seealso: DMProjectFunction(), DMComputeL2FieldDiff(), DMComputeL2GradientDiff() 6271 @*/ 6272 PetscErrorCode DMComputeL2FieldDiff(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[]) 6273 { 6274 PetscErrorCode ierr; 6275 6276 PetscFunctionBegin; 6277 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6278 PetscValidHeaderSpecific(X,VEC_CLASSID,5); 6279 if (!dm->ops->computel2fielddiff) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM type %s does not implemnt DMComputeL2FieldDiff",((PetscObject)dm)->type_name); 6280 ierr = (dm->ops->computel2fielddiff)(dm,time,funcs,ctxs,X,diff);CHKERRQ(ierr); 6281 PetscFunctionReturn(0); 6282 } 6283 6284 /*@C 6285 DMAdaptLabel - Adapt a dm based on a label with values interpreted as coarsening and refining flags. Specific implementations of DM maybe have 6286 specialized flags, but all implementations should accept flag values DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, DM_ADAPT_REFINE, and DM_ADAPT_COARSEN. 6287 6288 Collective on dm 6289 6290 Input parameters: 6291 + dm - the pre-adaptation DM object 6292 - label - label with the flags 6293 6294 Output parameters: 6295 . dmAdapt - the adapted DM object: may be NULL if an adapted DM could not be produced. 6296 6297 Level: intermediate 6298 6299 .seealso: DMAdaptMetric(), DMCoarsen(), DMRefine() 6300 @*/ 6301 PetscErrorCode DMAdaptLabel(DM dm, DMLabel label, DM *dmAdapt) 6302 { 6303 PetscErrorCode ierr; 6304 6305 PetscFunctionBegin; 6306 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6307 PetscValidPointer(label,2); 6308 PetscValidPointer(dmAdapt,3); 6309 *dmAdapt = NULL; 6310 if (!dm->ops->adaptlabel) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM type %s does not implement DMAdaptLabel",((PetscObject)dm)->type_name); 6311 ierr = (dm->ops->adaptlabel)(dm, label, dmAdapt);CHKERRQ(ierr); 6312 PetscFunctionReturn(0); 6313 } 6314 6315 /*@C 6316 DMAdaptMetric - Generates a mesh adapted to the specified metric field using the pragmatic library. 6317 6318 Input Parameters: 6319 + dm - The DM object 6320 . metric - The metric to which the mesh is adapted, defined vertex-wise. 6321 - bdLabel - Label for boundary tags, which will be preserved in the output mesh. bdLabel should be NULL if there is no such label, and should be different from "_boundary_". 6322 6323 Output Parameter: 6324 . dmAdapt - Pointer to the DM object containing the adapted mesh 6325 6326 Note: The label in the adapted mesh will be registered under the name of the input DMLabel object 6327 6328 Level: advanced 6329 6330 .seealso: DMAdaptLabel(), DMCoarsen(), DMRefine() 6331 @*/ 6332 PetscErrorCode DMAdaptMetric(DM dm, Vec metric, DMLabel bdLabel, DM *dmAdapt) 6333 { 6334 PetscErrorCode ierr; 6335 6336 PetscFunctionBegin; 6337 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6338 PetscValidHeaderSpecific(metric, VEC_CLASSID, 2); 6339 if (bdLabel) PetscValidPointer(bdLabel, 3); 6340 PetscValidPointer(dmAdapt, 4); 6341 *dmAdapt = NULL; 6342 if (!dm->ops->adaptmetric) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM type %s does not implement DMAdaptMetric",((PetscObject)dm)->type_name); 6343 ierr = (dm->ops->adaptmetric)(dm, metric, bdLabel, dmAdapt);CHKERRQ(ierr); 6344 PetscFunctionReturn(0); 6345 } 6346 6347 /*@C 6348 DMGetNeighbors - Gets an array containing the MPI rank of all the processes neighbors 6349 6350 Not Collective 6351 6352 Input Parameter: 6353 . dm - The DM 6354 6355 Output Parameter: 6356 . nranks - the number of neighbours 6357 . ranks - the neighbors ranks 6358 6359 Notes: 6360 Do not free the array, it is freed when the DM is destroyed. 6361 6362 Level: beginner 6363 6364 .seealso: DMDAGetNeighbors(), PetscSFGetRanks() 6365 @*/ 6366 PetscErrorCode DMGetNeighbors(DM dm,PetscInt *nranks,const PetscMPIInt *ranks[]) 6367 { 6368 PetscErrorCode ierr; 6369 6370 PetscFunctionBegin; 6371 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6372 if (!dm->ops->getneighbors) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM type %s does not implemnt DMGetNeighbors",((PetscObject)dm)->type_name); 6373 ierr = (dm->ops->getneighbors)(dm,nranks,ranks);CHKERRQ(ierr); 6374 PetscFunctionReturn(0); 6375 } 6376 6377 #include <petsc/private/matimpl.h> /* Needed because of coloring->ctype below */ 6378 6379 /* 6380 Converts the input vector to a ghosted vector and then calls the standard coloring code. 6381 This has be a different function because it requires DM which is not defined in the Mat library 6382 */ 6383 PetscErrorCode MatFDColoringApply_AIJDM(Mat J,MatFDColoring coloring,Vec x1,void *sctx) 6384 { 6385 PetscErrorCode ierr; 6386 6387 PetscFunctionBegin; 6388 if (coloring->ctype == IS_COLORING_LOCAL) { 6389 Vec x1local; 6390 DM dm; 6391 ierr = MatGetDM(J,&dm);CHKERRQ(ierr); 6392 if (!dm) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_INCOMP,"IS_COLORING_LOCAL requires a DM"); 6393 ierr = DMGetLocalVector(dm,&x1local);CHKERRQ(ierr); 6394 ierr = DMGlobalToLocalBegin(dm,x1,INSERT_VALUES,x1local);CHKERRQ(ierr); 6395 ierr = DMGlobalToLocalEnd(dm,x1,INSERT_VALUES,x1local);CHKERRQ(ierr); 6396 x1 = x1local; 6397 } 6398 ierr = MatFDColoringApply_AIJ(J,coloring,x1,sctx);CHKERRQ(ierr); 6399 if (coloring->ctype == IS_COLORING_LOCAL) { 6400 DM dm; 6401 ierr = MatGetDM(J,&dm);CHKERRQ(ierr); 6402 ierr = DMRestoreLocalVector(dm,&x1);CHKERRQ(ierr); 6403 } 6404 PetscFunctionReturn(0); 6405 } 6406 6407 /*@ 6408 MatFDColoringUseDM - allows a MatFDColoring object to use the DM associated with the matrix to use a IS_COLORING_LOCAL coloring 6409 6410 Input Parameter: 6411 . coloring - the MatFDColoring object 6412 6413 Developer Notes: this routine exists because the PETSc Mat library does not know about the DM objects 6414 6415 Level: advanced 6416 6417 .seealso: MatFDColoring, MatFDColoringCreate(), ISColoringType 6418 @*/ 6419 PetscErrorCode MatFDColoringUseDM(Mat coloring,MatFDColoring fdcoloring) 6420 { 6421 PetscFunctionBegin; 6422 coloring->ops->fdcoloringapply = MatFDColoringApply_AIJDM; 6423 PetscFunctionReturn(0); 6424 } 6425