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