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