1 #include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/ 2 #include <petscsf.h> 3 #include <petscds.h> 4 5 PetscClassId DM_CLASSID; 6 PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal, DM_LocalToLocal; 7 8 const char *const DMBoundaryTypes[] = {"NONE","GHOSTED","MIRROR","PERIODIC","TWIST","DM_BOUNDARY_",0}; 9 10 #undef __FUNCT__ 11 #define __FUNCT__ "DMCreate" 12 /*@ 13 DMCreate - Creates an empty DM object. The type can then be set with DMSetType(). 14 15 If you never call DMSetType() it will generate an 16 error when you try to use the vector. 17 18 Collective on MPI_Comm 19 20 Input Parameter: 21 . comm - The communicator for the DM object 22 23 Output Parameter: 24 . dm - The DM object 25 26 Level: beginner 27 28 .seealso: DMSetType(), DMDA, DMSLICED, DMCOMPOSITE 29 @*/ 30 PetscErrorCode DMCreate(MPI_Comm comm,DM *dm) 31 { 32 DM v; 33 PetscErrorCode ierr; 34 35 PetscFunctionBegin; 36 PetscValidPointer(dm,2); 37 *dm = NULL; 38 ierr = PetscSysInitializePackage();CHKERRQ(ierr); 39 ierr = VecInitializePackage();CHKERRQ(ierr); 40 ierr = MatInitializePackage();CHKERRQ(ierr); 41 ierr = DMInitializePackage();CHKERRQ(ierr); 42 43 ierr = PetscHeaderCreate(v, _p_DM, struct _DMOps, DM_CLASSID, "DM", "Distribution Manager", "DM", comm, DMDestroy, DMView);CHKERRQ(ierr); 44 ierr = PetscMemzero(v->ops, sizeof(struct _DMOps));CHKERRQ(ierr); 45 46 47 v->ltogmap = NULL; 48 v->bs = 1; 49 v->coloringtype = IS_COLORING_GLOBAL; 50 ierr = PetscSFCreate(comm, &v->sf);CHKERRQ(ierr); 51 ierr = PetscSFCreate(comm, &v->defaultSF);CHKERRQ(ierr); 52 v->defaultSection = NULL; 53 v->defaultGlobalSection = NULL; 54 v->L = NULL; 55 v->maxCell = NULL; 56 { 57 PetscInt i; 58 for (i = 0; i < 10; ++i) { 59 v->nullspaceConstructors[i] = NULL; 60 } 61 } 62 ierr = PetscDSCreate(comm, &v->prob);CHKERRQ(ierr); 63 v->dmBC = NULL; 64 v->outputSequenceNum = -1; 65 v->outputSequenceVal = 0.0; 66 ierr = DMSetVecType(v,VECSTANDARD);CHKERRQ(ierr); 67 ierr = DMSetMatType(v,MATAIJ);CHKERRQ(ierr); 68 *dm = v; 69 PetscFunctionReturn(0); 70 } 71 72 #undef __FUNCT__ 73 #define __FUNCT__ "DMClone" 74 /*@ 75 DMClone - Creates a DM object with the same topology as the original. 76 77 Collective on MPI_Comm 78 79 Input Parameter: 80 . dm - The original DM object 81 82 Output Parameter: 83 . newdm - The new DM object 84 85 Level: beginner 86 87 .keywords: DM, topology, create 88 @*/ 89 PetscErrorCode DMClone(DM dm, DM *newdm) 90 { 91 PetscSF sf; 92 Vec coords; 93 void *ctx; 94 PetscErrorCode ierr; 95 96 PetscFunctionBegin; 97 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 98 PetscValidPointer(newdm,2); 99 ierr = DMCreate(PetscObjectComm((PetscObject)dm), newdm);CHKERRQ(ierr); 100 if (dm->ops->clone) { 101 ierr = (*dm->ops->clone)(dm, newdm);CHKERRQ(ierr); 102 } 103 (*newdm)->setupcalled = PETSC_TRUE; 104 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 105 ierr = DMSetPointSF(*newdm, sf);CHKERRQ(ierr); 106 ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); 107 ierr = DMSetApplicationContext(*newdm, ctx);CHKERRQ(ierr); 108 ierr = DMGetCoordinatesLocal(dm, &coords);CHKERRQ(ierr); 109 if (coords) { 110 ierr = DMSetCoordinatesLocal(*newdm, coords);CHKERRQ(ierr); 111 } else { 112 ierr = DMGetCoordinates(dm, &coords);CHKERRQ(ierr); 113 if (coords) {ierr = DMSetCoordinates(*newdm, coords);CHKERRQ(ierr);} 114 } 115 if (dm->maxCell) { 116 const PetscReal *maxCell, *L; 117 ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr); 118 ierr = DMSetPeriodicity(*newdm, maxCell, L);CHKERRQ(ierr); 119 } 120 PetscFunctionReturn(0); 121 } 122 123 #undef __FUNCT__ 124 #define __FUNCT__ "DMSetVecType" 125 /*@C 126 DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() 127 128 Logically Collective on DMDA 129 130 Input Parameter: 131 + da - initial distributed array 132 . ctype - the vector type, currently either VECSTANDARD or VECCUSP 133 134 Options Database: 135 . -dm_vec_type ctype 136 137 Level: intermediate 138 139 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType, DMGetVecType() 140 @*/ 141 PetscErrorCode DMSetVecType(DM da,VecType ctype) 142 { 143 PetscErrorCode ierr; 144 145 PetscFunctionBegin; 146 PetscValidHeaderSpecific(da,DM_CLASSID,1); 147 ierr = PetscFree(da->vectype);CHKERRQ(ierr); 148 ierr = PetscStrallocpy(ctype,(char**)&da->vectype);CHKERRQ(ierr); 149 PetscFunctionReturn(0); 150 } 151 152 #undef __FUNCT__ 153 #define __FUNCT__ "DMGetVecType" 154 /*@C 155 DMGetVecType - Gets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() 156 157 Logically Collective on DMDA 158 159 Input Parameter: 160 . da - initial distributed array 161 162 Output Parameter: 163 . ctype - the vector type 164 165 Level: intermediate 166 167 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType 168 @*/ 169 PetscErrorCode DMGetVecType(DM da,VecType *ctype) 170 { 171 PetscFunctionBegin; 172 PetscValidHeaderSpecific(da,DM_CLASSID,1); 173 *ctype = da->vectype; 174 PetscFunctionReturn(0); 175 } 176 177 #undef __FUNCT__ 178 #define __FUNCT__ "VecGetDM" 179 /*@ 180 VecGetDM - Gets the DM defining the data layout of the vector 181 182 Not collective 183 184 Input Parameter: 185 . v - The Vec 186 187 Output Parameter: 188 . dm - The DM 189 190 Level: intermediate 191 192 .seealso: VecSetDM(), DMGetLocalVector(), DMGetGlobalVector(), DMSetVecType() 193 @*/ 194 PetscErrorCode VecGetDM(Vec v, DM *dm) 195 { 196 PetscErrorCode ierr; 197 198 PetscFunctionBegin; 199 PetscValidHeaderSpecific(v,VEC_CLASSID,1); 200 PetscValidPointer(dm,2); 201 ierr = PetscObjectQuery((PetscObject) v, "__PETSc_dm", (PetscObject*) dm);CHKERRQ(ierr); 202 PetscFunctionReturn(0); 203 } 204 205 #undef __FUNCT__ 206 #define __FUNCT__ "VecSetDM" 207 /*@ 208 VecSetDM - Sets the DM defining the data layout of the vector 209 210 Not collective 211 212 Input Parameters: 213 + v - The Vec 214 - dm - The DM 215 216 Level: intermediate 217 218 .seealso: VecGetDM(), DMGetLocalVector(), DMGetGlobalVector(), DMSetVecType() 219 @*/ 220 PetscErrorCode VecSetDM(Vec v, DM dm) 221 { 222 PetscErrorCode ierr; 223 224 PetscFunctionBegin; 225 PetscValidHeaderSpecific(v,VEC_CLASSID,1); 226 if (dm) PetscValidHeaderSpecific(dm,DM_CLASSID,2); 227 ierr = PetscObjectCompose((PetscObject) v, "__PETSc_dm", (PetscObject) dm);CHKERRQ(ierr); 228 PetscFunctionReturn(0); 229 } 230 231 #undef __FUNCT__ 232 #define __FUNCT__ "DMSetMatType" 233 /*@C 234 DMSetMatType - Sets the type of matrix created with DMCreateMatrix() 235 236 Logically Collective on DM 237 238 Input Parameter: 239 + dm - the DM context 240 . ctype - the matrix type 241 242 Options Database: 243 . -dm_mat_type ctype 244 245 Level: intermediate 246 247 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType, DMGetMatType() 248 @*/ 249 PetscErrorCode DMSetMatType(DM dm,MatType ctype) 250 { 251 PetscErrorCode ierr; 252 253 PetscFunctionBegin; 254 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 255 ierr = PetscFree(dm->mattype);CHKERRQ(ierr); 256 ierr = PetscStrallocpy(ctype,(char**)&dm->mattype);CHKERRQ(ierr); 257 PetscFunctionReturn(0); 258 } 259 260 #undef __FUNCT__ 261 #define __FUNCT__ "DMGetMatType" 262 /*@C 263 DMGetMatType - Gets the type of matrix created with DMCreateMatrix() 264 265 Logically Collective on DM 266 267 Input Parameter: 268 . dm - the DM context 269 270 Output Parameter: 271 . ctype - the matrix type 272 273 Options Database: 274 . -dm_mat_type ctype 275 276 Level: intermediate 277 278 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType, DMSetMatType() 279 @*/ 280 PetscErrorCode DMGetMatType(DM dm,MatType *ctype) 281 { 282 PetscFunctionBegin; 283 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 284 *ctype = dm->mattype; 285 PetscFunctionReturn(0); 286 } 287 288 #undef __FUNCT__ 289 #define __FUNCT__ "MatGetDM" 290 /*@ 291 MatGetDM - Gets the DM defining the data layout of the matrix 292 293 Not collective 294 295 Input Parameter: 296 . A - The Mat 297 298 Output Parameter: 299 . dm - The DM 300 301 Level: intermediate 302 303 .seealso: MatSetDM(), DMCreateMatrix(), DMSetMatType() 304 @*/ 305 PetscErrorCode MatGetDM(Mat A, DM *dm) 306 { 307 PetscErrorCode ierr; 308 309 PetscFunctionBegin; 310 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 311 PetscValidPointer(dm,2); 312 ierr = PetscObjectQuery((PetscObject) A, "__PETSc_dm", (PetscObject*) dm);CHKERRQ(ierr); 313 PetscFunctionReturn(0); 314 } 315 316 #undef __FUNCT__ 317 #define __FUNCT__ "MatSetDM" 318 /*@ 319 MatSetDM - Sets the DM defining the data layout of the matrix 320 321 Not collective 322 323 Input Parameters: 324 + A - The Mat 325 - dm - The DM 326 327 Level: intermediate 328 329 .seealso: MatGetDM(), DMCreateMatrix(), DMSetMatType() 330 @*/ 331 PetscErrorCode MatSetDM(Mat A, DM dm) 332 { 333 PetscErrorCode ierr; 334 335 PetscFunctionBegin; 336 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 337 if (dm) PetscValidHeaderSpecific(dm,DM_CLASSID,2); 338 ierr = PetscObjectCompose((PetscObject) A, "__PETSc_dm", (PetscObject) dm);CHKERRQ(ierr); 339 PetscFunctionReturn(0); 340 } 341 342 #undef __FUNCT__ 343 #define __FUNCT__ "DMSetOptionsPrefix" 344 /*@C 345 DMSetOptionsPrefix - Sets the prefix used for searching for all 346 DMDA options in the database. 347 348 Logically Collective on DMDA 349 350 Input Parameter: 351 + da - the DMDA context 352 - prefix - the prefix to prepend to all option names 353 354 Notes: 355 A hyphen (-) must NOT be given at the beginning of the prefix name. 356 The first character of all runtime options is AUTOMATICALLY the hyphen. 357 358 Level: advanced 359 360 .keywords: DMDA, set, options, prefix, database 361 362 .seealso: DMSetFromOptions() 363 @*/ 364 PetscErrorCode DMSetOptionsPrefix(DM dm,const char prefix[]) 365 { 366 PetscErrorCode ierr; 367 368 PetscFunctionBegin; 369 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 370 ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr); 371 PetscFunctionReturn(0); 372 } 373 374 #undef __FUNCT__ 375 #define __FUNCT__ "DMDestroy" 376 /*@ 377 DMDestroy - Destroys a vector packer or DMDA. 378 379 Collective on DM 380 381 Input Parameter: 382 . dm - the DM object to destroy 383 384 Level: developer 385 386 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 387 388 @*/ 389 PetscErrorCode DMDestroy(DM *dm) 390 { 391 PetscInt i, cnt = 0; 392 DMNamedVecLink nlink,nnext; 393 PetscErrorCode ierr; 394 395 PetscFunctionBegin; 396 if (!*dm) PetscFunctionReturn(0); 397 PetscValidHeaderSpecific((*dm),DM_CLASSID,1); 398 399 /* count all the circular references of DM and its contained Vecs */ 400 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 401 if ((*dm)->localin[i]) cnt++; 402 if ((*dm)->globalin[i]) cnt++; 403 } 404 for (nlink=(*dm)->namedglobal; nlink; nlink=nlink->next) cnt++; 405 for (nlink=(*dm)->namedlocal; nlink; nlink=nlink->next) cnt++; 406 if ((*dm)->x) { 407 DM obj; 408 ierr = VecGetDM((*dm)->x, &obj);CHKERRQ(ierr); 409 if (obj == *dm) cnt++; 410 } 411 412 if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);} 413 /* 414 Need this test because the dm references the vectors that 415 reference the dm, so destroying the dm calls destroy on the 416 vectors that cause another destroy on the dm 417 */ 418 if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0); 419 ((PetscObject) (*dm))->refct = 0; 420 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 421 if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()"); 422 ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr); 423 } 424 nnext=(*dm)->namedglobal; 425 (*dm)->namedglobal = NULL; 426 for (nlink=nnext; nlink; nlink=nnext) { /* Destroy the named vectors */ 427 nnext = nlink->next; 428 if (nlink->status != DMVEC_STATUS_IN) SETERRQ1(((PetscObject)*dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM still has Vec named '%s' checked out",nlink->name); 429 ierr = PetscFree(nlink->name);CHKERRQ(ierr); 430 ierr = VecDestroy(&nlink->X);CHKERRQ(ierr); 431 ierr = PetscFree(nlink);CHKERRQ(ierr); 432 } 433 nnext=(*dm)->namedlocal; 434 (*dm)->namedlocal = NULL; 435 for (nlink=nnext; nlink; nlink=nnext) { /* Destroy the named local vectors */ 436 nnext = nlink->next; 437 if (nlink->status != DMVEC_STATUS_IN) SETERRQ1(((PetscObject)*dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM still has Vec named '%s' checked out",nlink->name); 438 ierr = PetscFree(nlink->name);CHKERRQ(ierr); 439 ierr = VecDestroy(&nlink->X);CHKERRQ(ierr); 440 ierr = PetscFree(nlink);CHKERRQ(ierr); 441 } 442 443 /* Destroy the list of hooks */ 444 { 445 DMCoarsenHookLink link,next; 446 for (link=(*dm)->coarsenhook; link; link=next) { 447 next = link->next; 448 ierr = PetscFree(link);CHKERRQ(ierr); 449 } 450 (*dm)->coarsenhook = NULL; 451 } 452 { 453 DMRefineHookLink link,next; 454 for (link=(*dm)->refinehook; link; link=next) { 455 next = link->next; 456 ierr = PetscFree(link);CHKERRQ(ierr); 457 } 458 (*dm)->refinehook = NULL; 459 } 460 { 461 DMSubDomainHookLink link,next; 462 for (link=(*dm)->subdomainhook; link; link=next) { 463 next = link->next; 464 ierr = PetscFree(link);CHKERRQ(ierr); 465 } 466 (*dm)->subdomainhook = NULL; 467 } 468 { 469 DMGlobalToLocalHookLink link,next; 470 for (link=(*dm)->gtolhook; link; link=next) { 471 next = link->next; 472 ierr = PetscFree(link);CHKERRQ(ierr); 473 } 474 (*dm)->gtolhook = NULL; 475 } 476 /* Destroy the work arrays */ 477 { 478 DMWorkLink link,next; 479 if ((*dm)->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Work array still checked out"); 480 for (link=(*dm)->workin; link; link=next) { 481 next = link->next; 482 ierr = PetscFree(link->mem);CHKERRQ(ierr); 483 ierr = PetscFree(link);CHKERRQ(ierr); 484 } 485 (*dm)->workin = NULL; 486 } 487 488 ierr = PetscObjectDestroy(&(*dm)->dmksp);CHKERRQ(ierr); 489 ierr = PetscObjectDestroy(&(*dm)->dmsnes);CHKERRQ(ierr); 490 ierr = PetscObjectDestroy(&(*dm)->dmts);CHKERRQ(ierr); 491 492 if ((*dm)->ctx && (*dm)->ctxdestroy) { 493 ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr); 494 } 495 ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr); 496 ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr); 497 ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr); 498 ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr); 499 ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr); 500 ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr); 501 502 ierr = PetscSectionDestroy(&(*dm)->defaultSection);CHKERRQ(ierr); 503 ierr = PetscSectionDestroy(&(*dm)->defaultGlobalSection);CHKERRQ(ierr); 504 ierr = PetscLayoutDestroy(&(*dm)->map);CHKERRQ(ierr); 505 ierr = PetscSFDestroy(&(*dm)->sf);CHKERRQ(ierr); 506 ierr = PetscSFDestroy(&(*dm)->defaultSF);CHKERRQ(ierr); 507 508 ierr = DMDestroy(&(*dm)->coordinateDM);CHKERRQ(ierr); 509 ierr = VecDestroy(&(*dm)->coordinates);CHKERRQ(ierr); 510 ierr = VecDestroy(&(*dm)->coordinatesLocal);CHKERRQ(ierr); 511 ierr = PetscFree((*dm)->maxCell);CHKERRQ(ierr); 512 ierr = PetscFree((*dm)->L);CHKERRQ(ierr); 513 514 ierr = PetscDSDestroy(&(*dm)->prob);CHKERRQ(ierr); 515 ierr = DMDestroy(&(*dm)->dmBC);CHKERRQ(ierr); 516 /* if memory was published with SAWs then destroy it */ 517 ierr = PetscObjectSAWsViewOff((PetscObject)*dm);CHKERRQ(ierr); 518 519 ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr); 520 /* We do not destroy (*dm)->data here so that we can reference count backend objects */ 521 ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr); 522 PetscFunctionReturn(0); 523 } 524 525 #undef __FUNCT__ 526 #define __FUNCT__ "DMSetUp" 527 /*@ 528 DMSetUp - sets up the data structures inside a DM object 529 530 Collective on DM 531 532 Input Parameter: 533 . dm - the DM object to setup 534 535 Level: developer 536 537 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 538 539 @*/ 540 PetscErrorCode DMSetUp(DM dm) 541 { 542 PetscErrorCode ierr; 543 544 PetscFunctionBegin; 545 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 546 if (dm->setupcalled) PetscFunctionReturn(0); 547 if (dm->ops->setup) { 548 ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr); 549 } 550 dm->setupcalled = PETSC_TRUE; 551 PetscFunctionReturn(0); 552 } 553 554 #undef __FUNCT__ 555 #define __FUNCT__ "DMSetFromOptions" 556 /*@ 557 DMSetFromOptions - sets parameters in a DM from the options database 558 559 Collective on DM 560 561 Input Parameter: 562 . dm - the DM object to set options for 563 564 Options Database: 565 + -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros 566 . -dm_vec_type <type> type of vector to create inside DM 567 . -dm_mat_type <type> type of matrix to create inside DM 568 - -dm_coloring_type <global or ghosted> 569 570 Level: developer 571 572 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 573 574 @*/ 575 PetscErrorCode DMSetFromOptions(DM dm) 576 { 577 char typeName[256]; 578 PetscBool flg; 579 PetscErrorCode ierr; 580 581 PetscFunctionBegin; 582 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 583 ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr); 584 ierr = PetscOptionsBool("-dm_preallocate_only","only preallocate matrix, but do not set column indices","DMSetMatrixPreallocateOnly",dm->prealloc_only,&dm->prealloc_only,NULL);CHKERRQ(ierr); 585 ierr = PetscOptionsFList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr); 586 if (flg) { 587 ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr); 588 } 589 ierr = PetscOptionsFList("-dm_mat_type","Matrix type used for created matrices","DMSetMatType",MatList,dm->mattype ? dm->mattype : typeName,typeName,sizeof(typeName),&flg);CHKERRQ(ierr); 590 if (flg) { 591 ierr = DMSetMatType(dm,typeName);CHKERRQ(ierr); 592 } 593 ierr = PetscOptionsEnum("-dm_is_coloring_type","Global or local coloring of Jacobian","ISColoringType",ISColoringTypes,(PetscEnum)dm->coloringtype,(PetscEnum*)&dm->coloringtype,NULL);CHKERRQ(ierr); 594 if (dm->ops->setfromoptions) { 595 ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr); 596 } 597 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 598 ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr); 599 ierr = PetscOptionsEnd();CHKERRQ(ierr); 600 PetscFunctionReturn(0); 601 } 602 603 #undef __FUNCT__ 604 #define __FUNCT__ "DMView" 605 /*@C 606 DMView - Views a vector packer or DMDA. 607 608 Collective on DM 609 610 Input Parameter: 611 + dm - the DM object to view 612 - v - the viewer 613 614 Level: developer 615 616 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 617 618 @*/ 619 PetscErrorCode DMView(DM dm,PetscViewer v) 620 { 621 PetscErrorCode ierr; 622 PetscBool isbinary; 623 624 PetscFunctionBegin; 625 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 626 if (!v) { 627 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)dm),&v);CHKERRQ(ierr); 628 } 629 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)dm,v);CHKERRQ(ierr); 630 ierr = PetscObjectTypeCompare((PetscObject)v,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 631 if (isbinary) { 632 PetscInt classid = DM_FILE_CLASSID; 633 char type[256]; 634 635 ierr = PetscViewerBinaryWrite(v,&classid,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 636 ierr = PetscStrncpy(type,((PetscObject)dm)->type_name,256);CHKERRQ(ierr); 637 ierr = PetscViewerBinaryWrite(v,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 638 } 639 if (dm->ops->view) { 640 ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr); 641 } 642 PetscFunctionReturn(0); 643 } 644 645 #undef __FUNCT__ 646 #define __FUNCT__ "DMCreateGlobalVector" 647 /*@ 648 DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object 649 650 Collective on DM 651 652 Input Parameter: 653 . dm - the DM object 654 655 Output Parameter: 656 . vec - the global vector 657 658 Level: beginner 659 660 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 661 662 @*/ 663 PetscErrorCode DMCreateGlobalVector(DM dm,Vec *vec) 664 { 665 PetscErrorCode ierr; 666 667 PetscFunctionBegin; 668 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 669 ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr); 670 PetscFunctionReturn(0); 671 } 672 673 #undef __FUNCT__ 674 #define __FUNCT__ "DMCreateLocalVector" 675 /*@ 676 DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object 677 678 Not Collective 679 680 Input Parameter: 681 . dm - the DM object 682 683 Output Parameter: 684 . vec - the local vector 685 686 Level: beginner 687 688 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 689 690 @*/ 691 PetscErrorCode DMCreateLocalVector(DM dm,Vec *vec) 692 { 693 PetscErrorCode ierr; 694 695 PetscFunctionBegin; 696 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 697 ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr); 698 PetscFunctionReturn(0); 699 } 700 701 #undef __FUNCT__ 702 #define __FUNCT__ "DMGetLocalToGlobalMapping" 703 /*@ 704 DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM. 705 706 Collective on DM 707 708 Input Parameter: 709 . dm - the DM that provides the mapping 710 711 Output Parameter: 712 . ltog - the mapping 713 714 Level: intermediate 715 716 Notes: 717 This mapping can then be used by VecSetLocalToGlobalMapping() or 718 MatSetLocalToGlobalMapping(). 719 720 .seealso: DMCreateLocalVector() 721 @*/ 722 PetscErrorCode DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog) 723 { 724 PetscErrorCode ierr; 725 726 PetscFunctionBegin; 727 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 728 PetscValidPointer(ltog,2); 729 if (!dm->ltogmap) { 730 PetscSection section, sectionGlobal; 731 732 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 733 if (section) { 734 PetscInt *ltog; 735 PetscInt pStart, pEnd, size, p, l; 736 737 ierr = DMGetDefaultGlobalSection(dm, §ionGlobal);CHKERRQ(ierr); 738 ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 739 ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr); 740 ierr = PetscMalloc1(size, <og);CHKERRQ(ierr); /* We want the local+overlap size */ 741 for (p = pStart, l = 0; p < pEnd; ++p) { 742 PetscInt dof, off, c; 743 744 /* Should probably use constrained dofs */ 745 ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 746 ierr = PetscSectionGetOffset(sectionGlobal, p, &off);CHKERRQ(ierr); 747 for (c = 0; c < dof; ++c, ++l) { 748 ltog[l] = off+c; 749 } 750 } 751 ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF, 1,size, ltog, PETSC_OWN_POINTER, &dm->ltogmap);CHKERRQ(ierr); 752 ierr = PetscLogObjectParent((PetscObject)dm, (PetscObject)dm->ltogmap);CHKERRQ(ierr); 753 } else { 754 if (!dm->ops->getlocaltoglobalmapping) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping"); 755 ierr = (*dm->ops->getlocaltoglobalmapping)(dm);CHKERRQ(ierr); 756 } 757 } 758 *ltog = dm->ltogmap; 759 PetscFunctionReturn(0); 760 } 761 762 #undef __FUNCT__ 763 #define __FUNCT__ "DMGetBlockSize" 764 /*@ 765 DMGetBlockSize - Gets the inherent block size associated with a DM 766 767 Not Collective 768 769 Input Parameter: 770 . dm - the DM with block structure 771 772 Output Parameter: 773 . bs - the block size, 1 implies no exploitable block structure 774 775 Level: intermediate 776 777 .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMapping() 778 @*/ 779 PetscErrorCode DMGetBlockSize(DM dm,PetscInt *bs) 780 { 781 PetscFunctionBegin; 782 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 783 PetscValidPointer(bs,2); 784 if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet"); 785 *bs = dm->bs; 786 PetscFunctionReturn(0); 787 } 788 789 #undef __FUNCT__ 790 #define __FUNCT__ "DMCreateInterpolation" 791 /*@ 792 DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects 793 794 Collective on DM 795 796 Input Parameter: 797 + dm1 - the DM object 798 - dm2 - the second, finer DM object 799 800 Output Parameter: 801 + mat - the interpolation 802 - vec - the scaling (optional) 803 804 Level: developer 805 806 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 807 DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation. 808 809 For DMDA objects you can use this interpolation (more precisely the interpolation from the DMGetCoordinateDM()) to interpolate the mesh coordinate vectors 810 EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic. 811 812 813 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen() 814 815 @*/ 816 PetscErrorCode DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec) 817 { 818 PetscErrorCode ierr; 819 820 PetscFunctionBegin; 821 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 822 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 823 ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr); 824 PetscFunctionReturn(0); 825 } 826 827 #undef __FUNCT__ 828 #define __FUNCT__ "DMCreateInjection" 829 /*@ 830 DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects 831 832 Collective on DM 833 834 Input Parameter: 835 + dm1 - the DM object 836 - dm2 - the second, finer DM object 837 838 Output Parameter: 839 . ctx - the injection 840 841 Level: developer 842 843 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 844 DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection. 845 846 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation() 847 848 @*/ 849 PetscErrorCode DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx) 850 { 851 PetscErrorCode ierr; 852 853 PetscFunctionBegin; 854 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 855 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 856 ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr); 857 PetscFunctionReturn(0); 858 } 859 860 #undef __FUNCT__ 861 #define __FUNCT__ "DMCreateColoring" 862 /*@ 863 DMCreateColoring - Gets coloring for a DM 864 865 Collective on DM 866 867 Input Parameter: 868 + dm - the DM object 869 - ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL 870 871 Output Parameter: 872 . coloring - the coloring 873 874 Level: developer 875 876 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix(), DMSetMatType() 877 878 @*/ 879 PetscErrorCode DMCreateColoring(DM dm,ISColoringType ctype,ISColoring *coloring) 880 { 881 PetscErrorCode ierr; 882 883 PetscFunctionBegin; 884 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 885 if (!dm->ops->getcoloring) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No coloring for this type of DM yet"); 886 ierr = (*dm->ops->getcoloring)(dm,ctype,coloring);CHKERRQ(ierr); 887 PetscFunctionReturn(0); 888 } 889 890 #undef __FUNCT__ 891 #define __FUNCT__ "DMCreateMatrix" 892 /*@ 893 DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite 894 895 Collective on DM 896 897 Input Parameter: 898 . dm - the DM object 899 900 Output Parameter: 901 . mat - the empty Jacobian 902 903 Level: beginner 904 905 Notes: This properly preallocates the number of nonzeros in the sparse matrix so you 906 do not need to do it yourself. 907 908 By default it also sets the nonzero structure and puts in the zero entries. To prevent setting 909 the nonzero pattern call DMDASetMatPreallocateOnly() 910 911 For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used 912 internally by PETSc. 913 914 For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires 915 the indices for the global numbering for DMDAs which is complicated. 916 917 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMSetMatType() 918 919 @*/ 920 PetscErrorCode DMCreateMatrix(DM dm,Mat *mat) 921 { 922 PetscErrorCode ierr; 923 924 PetscFunctionBegin; 925 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 926 ierr = MatInitializePackage();CHKERRQ(ierr); 927 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 928 PetscValidPointer(mat,3); 929 ierr = (*dm->ops->creatematrix)(dm,mat);CHKERRQ(ierr); 930 PetscFunctionReturn(0); 931 } 932 933 #undef __FUNCT__ 934 #define __FUNCT__ "DMSetMatrixPreallocateOnly" 935 /*@ 936 DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly 937 preallocated but the nonzero structure and zero values will not be set. 938 939 Logically Collective on DMDA 940 941 Input Parameter: 942 + dm - the DM 943 - only - PETSC_TRUE if only want preallocation 944 945 Level: developer 946 .seealso DMCreateMatrix() 947 @*/ 948 PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only) 949 { 950 PetscFunctionBegin; 951 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 952 dm->prealloc_only = only; 953 PetscFunctionReturn(0); 954 } 955 956 #undef __FUNCT__ 957 #define __FUNCT__ "DMGetWorkArray" 958 /*@C 959 DMGetWorkArray - Gets a work array guaranteed to be at least the input size, restore with DMRestoreWorkArray() 960 961 Not Collective 962 963 Input Parameters: 964 + dm - the DM object 965 . count - The minium size 966 - dtype - data type (PETSC_REAL, PETSC_SCALAR, PETSC_INT) 967 968 Output Parameter: 969 . array - the work array 970 971 Level: developer 972 973 .seealso DMDestroy(), DMCreate() 974 @*/ 975 PetscErrorCode DMGetWorkArray(DM dm,PetscInt count,PetscDataType dtype,void *mem) 976 { 977 PetscErrorCode ierr; 978 DMWorkLink link; 979 size_t size; 980 981 PetscFunctionBegin; 982 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 983 PetscValidPointer(mem,4); 984 if (dm->workin) { 985 link = dm->workin; 986 dm->workin = dm->workin->next; 987 } else { 988 ierr = PetscNewLog(dm,&link);CHKERRQ(ierr); 989 } 990 ierr = PetscDataTypeGetSize(dtype,&size);CHKERRQ(ierr); 991 if (size*count > link->bytes) { 992 ierr = PetscFree(link->mem);CHKERRQ(ierr); 993 ierr = PetscMalloc(size*count,&link->mem);CHKERRQ(ierr); 994 link->bytes = size*count; 995 } 996 link->next = dm->workout; 997 dm->workout = link; 998 *(void**)mem = link->mem; 999 PetscFunctionReturn(0); 1000 } 1001 1002 #undef __FUNCT__ 1003 #define __FUNCT__ "DMRestoreWorkArray" 1004 /*@C 1005 DMRestoreWorkArray - Restores a work array guaranteed to be at least the input size, restore with DMRestoreWorkArray() 1006 1007 Not Collective 1008 1009 Input Parameters: 1010 + dm - the DM object 1011 . count - The minium size 1012 - dtype - data type (PETSC_REAL, PETSC_SCALAR, PETSC_INT) 1013 1014 Output Parameter: 1015 . array - the work array 1016 1017 Level: developer 1018 1019 .seealso DMDestroy(), DMCreate() 1020 @*/ 1021 PetscErrorCode DMRestoreWorkArray(DM dm,PetscInt count,PetscDataType dtype,void *mem) 1022 { 1023 DMWorkLink *p,link; 1024 1025 PetscFunctionBegin; 1026 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1027 PetscValidPointer(mem,4); 1028 for (p=&dm->workout; (link=*p); p=&link->next) { 1029 if (link->mem == *(void**)mem) { 1030 *p = link->next; 1031 link->next = dm->workin; 1032 dm->workin = link; 1033 *(void**)mem = NULL; 1034 PetscFunctionReturn(0); 1035 } 1036 } 1037 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Array was not checked out"); 1038 PetscFunctionReturn(0); 1039 } 1040 1041 #undef __FUNCT__ 1042 #define __FUNCT__ "DMSetNullSpaceConstructor" 1043 PetscErrorCode DMSetNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (*nullsp)(DM dm, PetscInt field, MatNullSpace *nullSpace)) 1044 { 1045 PetscFunctionBegin; 1046 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1047 if (field >= 10) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %d >= 10 fields", field); 1048 dm->nullspaceConstructors[field] = nullsp; 1049 PetscFunctionReturn(0); 1050 } 1051 1052 #undef __FUNCT__ 1053 #define __FUNCT__ "DMCreateFieldIS" 1054 /*@C 1055 DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field 1056 1057 Not collective 1058 1059 Input Parameter: 1060 . dm - the DM object 1061 1062 Output Parameters: 1063 + numFields - The number of fields (or NULL if not requested) 1064 . fieldNames - The name for each field (or NULL if not requested) 1065 - fields - The global indices for each field (or NULL if not requested) 1066 1067 Level: intermediate 1068 1069 Notes: 1070 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 1071 PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with 1072 PetscFree(). 1073 1074 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 1075 @*/ 1076 PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields) 1077 { 1078 PetscSection section, sectionGlobal; 1079 PetscErrorCode ierr; 1080 1081 PetscFunctionBegin; 1082 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1083 if (numFields) { 1084 PetscValidPointer(numFields,2); 1085 *numFields = 0; 1086 } 1087 if (fieldNames) { 1088 PetscValidPointer(fieldNames,3); 1089 *fieldNames = NULL; 1090 } 1091 if (fields) { 1092 PetscValidPointer(fields,4); 1093 *fields = NULL; 1094 } 1095 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1096 if (section) { 1097 PetscInt *fieldSizes, **fieldIndices; 1098 PetscInt nF, f, pStart, pEnd, p; 1099 1100 ierr = DMGetDefaultGlobalSection(dm, §ionGlobal);CHKERRQ(ierr); 1101 ierr = PetscSectionGetNumFields(section, &nF);CHKERRQ(ierr); 1102 ierr = PetscMalloc2(nF,&fieldSizes,nF,&fieldIndices);CHKERRQ(ierr); 1103 ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr); 1104 for (f = 0; f < nF; ++f) { 1105 fieldSizes[f] = 0; 1106 } 1107 for (p = pStart; p < pEnd; ++p) { 1108 PetscInt gdof; 1109 1110 ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr); 1111 if (gdof > 0) { 1112 for (f = 0; f < nF; ++f) { 1113 PetscInt fdof, fcdof; 1114 1115 ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr); 1116 ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr); 1117 fieldSizes[f] += fdof-fcdof; 1118 } 1119 } 1120 } 1121 for (f = 0; f < nF; ++f) { 1122 ierr = PetscMalloc1(fieldSizes[f], &fieldIndices[f]);CHKERRQ(ierr); 1123 fieldSizes[f] = 0; 1124 } 1125 for (p = pStart; p < pEnd; ++p) { 1126 PetscInt gdof, goff; 1127 1128 ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr); 1129 if (gdof > 0) { 1130 ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 1131 for (f = 0; f < nF; ++f) { 1132 PetscInt fdof, fcdof, fc; 1133 1134 ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr); 1135 ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr); 1136 for (fc = 0; fc < fdof-fcdof; ++fc, ++fieldSizes[f]) { 1137 fieldIndices[f][fieldSizes[f]] = goff++; 1138 } 1139 } 1140 } 1141 } 1142 if (numFields) *numFields = nF; 1143 if (fieldNames) { 1144 ierr = PetscMalloc1(nF, fieldNames);CHKERRQ(ierr); 1145 for (f = 0; f < nF; ++f) { 1146 const char *fieldName; 1147 1148 ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr); 1149 ierr = PetscStrallocpy(fieldName, (char**) &(*fieldNames)[f]);CHKERRQ(ierr); 1150 } 1151 } 1152 if (fields) { 1153 ierr = PetscMalloc1(nF, fields);CHKERRQ(ierr); 1154 for (f = 0; f < nF; ++f) { 1155 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), fieldSizes[f], fieldIndices[f], PETSC_OWN_POINTER, &(*fields)[f]);CHKERRQ(ierr); 1156 } 1157 } 1158 ierr = PetscFree2(fieldSizes,fieldIndices);CHKERRQ(ierr); 1159 } else if (dm->ops->createfieldis) { 1160 ierr = (*dm->ops->createfieldis)(dm, numFields, fieldNames, fields);CHKERRQ(ierr); 1161 } 1162 PetscFunctionReturn(0); 1163 } 1164 1165 1166 #undef __FUNCT__ 1167 #define __FUNCT__ "DMCreateFieldDecomposition" 1168 /*@C 1169 DMCreateFieldDecomposition - Returns a list of IS objects defining a decomposition of a problem into subproblems 1170 corresponding to different fields: each IS contains the global indices of the dofs of the 1171 corresponding field. The optional list of DMs define the DM for each subproblem. 1172 Generalizes DMCreateFieldIS(). 1173 1174 Not collective 1175 1176 Input Parameter: 1177 . dm - the DM object 1178 1179 Output Parameters: 1180 + len - The number of subproblems in the field decomposition (or NULL if not requested) 1181 . namelist - The name for each field (or NULL if not requested) 1182 . islist - The global indices for each field (or NULL if not requested) 1183 - dmlist - The DMs for each field subproblem (or NULL, if not requested; if NULL is returned, no DMs are defined) 1184 1185 Level: intermediate 1186 1187 Notes: 1188 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 1189 PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(), 1190 and all of the arrays should be freed with PetscFree(). 1191 1192 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS() 1193 @*/ 1194 PetscErrorCode DMCreateFieldDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 1195 { 1196 PetscErrorCode ierr; 1197 1198 PetscFunctionBegin; 1199 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1200 if (len) { 1201 PetscValidPointer(len,2); 1202 *len = 0; 1203 } 1204 if (namelist) { 1205 PetscValidPointer(namelist,3); 1206 *namelist = 0; 1207 } 1208 if (islist) { 1209 PetscValidPointer(islist,4); 1210 *islist = 0; 1211 } 1212 if (dmlist) { 1213 PetscValidPointer(dmlist,5); 1214 *dmlist = 0; 1215 } 1216 /* 1217 Is it a good idea to apply the following check across all impls? 1218 Perhaps some impls can have a well-defined decomposition before DMSetUp? 1219 This, however, follows the general principle that accessors are not well-behaved until the object is set up. 1220 */ 1221 if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE, "Decomposition defined only after DMSetUp"); 1222 if (!dm->ops->createfielddecomposition) { 1223 PetscSection section; 1224 PetscInt numFields, f; 1225 1226 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1227 if (section) {ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);} 1228 if (section && numFields && dm->ops->createsubdm) { 1229 *len = numFields; 1230 if (namelist) {ierr = PetscMalloc1(numFields,namelist);CHKERRQ(ierr);} 1231 if (islist) {ierr = PetscMalloc1(numFields,islist);CHKERRQ(ierr);} 1232 if (dmlist) {ierr = PetscMalloc1(numFields,dmlist);CHKERRQ(ierr);} 1233 for (f = 0; f < numFields; ++f) { 1234 const char *fieldName; 1235 1236 ierr = DMCreateSubDM(dm, 1, &f, islist ? &(*islist)[f] : NULL, dmlist ? &(*dmlist)[f] : NULL);CHKERRQ(ierr); 1237 if (namelist) { 1238 ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr); 1239 ierr = PetscStrallocpy(fieldName, (char**) &(*namelist)[f]);CHKERRQ(ierr); 1240 } 1241 } 1242 } else { 1243 ierr = DMCreateFieldIS(dm, len, namelist, islist);CHKERRQ(ierr); 1244 /* By default there are no DMs associated with subproblems. */ 1245 if (dmlist) *dmlist = NULL; 1246 } 1247 } else { 1248 ierr = (*dm->ops->createfielddecomposition)(dm,len,namelist,islist,dmlist);CHKERRQ(ierr); 1249 } 1250 PetscFunctionReturn(0); 1251 } 1252 1253 #undef __FUNCT__ 1254 #define __FUNCT__ "DMCreateSubDM" 1255 /*@C 1256 DMCreateSubDM - Returns an IS and DM encapsulating a subproblem defined by the fields passed in. 1257 The fields are defined by DMCreateFieldIS(). 1258 1259 Not collective 1260 1261 Input Parameters: 1262 + dm - the DM object 1263 . numFields - number of fields in this subproblem 1264 - len - The number of subproblems in the decomposition (or NULL if not requested) 1265 1266 Output Parameters: 1267 . is - The global indices for the subproblem 1268 - dm - The DM for the subproblem 1269 1270 Level: intermediate 1271 1272 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS() 1273 @*/ 1274 PetscErrorCode DMCreateSubDM(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm) 1275 { 1276 PetscErrorCode ierr; 1277 1278 PetscFunctionBegin; 1279 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1280 PetscValidPointer(fields,3); 1281 if (is) PetscValidPointer(is,4); 1282 if (subdm) PetscValidPointer(subdm,5); 1283 if (dm->ops->createsubdm) { 1284 ierr = (*dm->ops->createsubdm)(dm, numFields, fields, is, subdm);CHKERRQ(ierr); 1285 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This type has no DMCreateSubDM implementation defined"); 1286 PetscFunctionReturn(0); 1287 } 1288 1289 1290 #undef __FUNCT__ 1291 #define __FUNCT__ "DMCreateDomainDecomposition" 1292 /*@C 1293 DMCreateDomainDecomposition - Returns lists of IS objects defining a decomposition of a problem into subproblems 1294 corresponding to restrictions to pairs nested subdomains: each IS contains the global 1295 indices of the dofs of the corresponding subdomains. The inner subdomains conceptually 1296 define a nonoverlapping covering, while outer subdomains can overlap. 1297 The optional list of DMs define the DM for each subproblem. 1298 1299 Not collective 1300 1301 Input Parameter: 1302 . dm - the DM object 1303 1304 Output Parameters: 1305 + len - The number of subproblems in the domain decomposition (or NULL if not requested) 1306 . namelist - The name for each subdomain (or NULL if not requested) 1307 . innerislist - The global indices for each inner subdomain (or NULL, if not requested) 1308 . outerislist - The global indices for each outer subdomain (or NULL, if not requested) 1309 - dmlist - The DMs for each subdomain subproblem (or NULL, if not requested; if NULL is returned, no DMs are defined) 1310 1311 Level: intermediate 1312 1313 Notes: 1314 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 1315 PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(), 1316 and all of the arrays should be freed with PetscFree(). 1317 1318 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateDomainDecompositionDM(), DMCreateFieldDecomposition() 1319 @*/ 1320 PetscErrorCode DMCreateDomainDecomposition(DM dm, PetscInt *len, char ***namelist, IS **innerislist, IS **outerislist, DM **dmlist) 1321 { 1322 PetscErrorCode ierr; 1323 DMSubDomainHookLink link; 1324 PetscInt i,l; 1325 1326 PetscFunctionBegin; 1327 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1328 if (len) {PetscValidPointer(len,2); *len = 0;} 1329 if (namelist) {PetscValidPointer(namelist,3); *namelist = NULL;} 1330 if (innerislist) {PetscValidPointer(innerislist,4); *innerislist = NULL;} 1331 if (outerislist) {PetscValidPointer(outerislist,5); *outerislist = NULL;} 1332 if (dmlist) {PetscValidPointer(dmlist,6); *dmlist = NULL;} 1333 /* 1334 Is it a good idea to apply the following check across all impls? 1335 Perhaps some impls can have a well-defined decomposition before DMSetUp? 1336 This, however, follows the general principle that accessors are not well-behaved until the object is set up. 1337 */ 1338 if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE, "Decomposition defined only after DMSetUp"); 1339 if (dm->ops->createdomaindecomposition) { 1340 ierr = (*dm->ops->createdomaindecomposition)(dm,&l,namelist,innerislist,outerislist,dmlist);CHKERRQ(ierr); 1341 /* copy subdomain hooks and context over to the subdomain DMs */ 1342 if (dmlist) { 1343 for (i = 0; i < l; i++) { 1344 for (link=dm->subdomainhook; link; link=link->next) { 1345 if (link->ddhook) {ierr = (*link->ddhook)(dm,(*dmlist)[i],link->ctx);CHKERRQ(ierr);} 1346 } 1347 (*dmlist)[i]->ctx = dm->ctx; 1348 } 1349 } 1350 if (len) *len = l; 1351 } 1352 PetscFunctionReturn(0); 1353 } 1354 1355 1356 #undef __FUNCT__ 1357 #define __FUNCT__ "DMCreateDomainDecompositionScatters" 1358 /*@C 1359 DMCreateDomainDecompositionScatters - Returns scatters to the subdomain vectors from the global vector 1360 1361 Not collective 1362 1363 Input Parameters: 1364 + dm - the DM object 1365 . n - the number of subdomain scatters 1366 - subdms - the local subdomains 1367 1368 Output Parameters: 1369 + n - the number of scatters returned 1370 . iscat - scatter from global vector to nonoverlapping global vector entries on subdomain 1371 . oscat - scatter from global vector to overlapping global vector entries on subdomain 1372 - gscat - scatter from global vector to local vector on subdomain (fills in ghosts) 1373 1374 Notes: This is an alternative to the iis and ois arguments in DMCreateDomainDecomposition that allow for the solution 1375 of general nonlinear problems with overlapping subdomain methods. While merely having index sets that enable subsets 1376 of the residual equations to be created is fine for linear problems, nonlinear problems require local assembly of 1377 solution and residual data. 1378 1379 Level: developer 1380 1381 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS() 1382 @*/ 1383 PetscErrorCode DMCreateDomainDecompositionScatters(DM dm,PetscInt n,DM *subdms,VecScatter **iscat,VecScatter **oscat,VecScatter **gscat) 1384 { 1385 PetscErrorCode ierr; 1386 1387 PetscFunctionBegin; 1388 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1389 PetscValidPointer(subdms,3); 1390 if (dm->ops->createddscatters) { 1391 ierr = (*dm->ops->createddscatters)(dm,n,subdms,iscat,oscat,gscat);CHKERRQ(ierr); 1392 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This type has no DMCreateDomainDecompositionLocalScatter implementation defined"); 1393 PetscFunctionReturn(0); 1394 } 1395 1396 #undef __FUNCT__ 1397 #define __FUNCT__ "DMRefine" 1398 /*@ 1399 DMRefine - Refines a DM object 1400 1401 Collective on DM 1402 1403 Input Parameter: 1404 + dm - the DM object 1405 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 1406 1407 Output Parameter: 1408 . dmf - the refined DM, or NULL 1409 1410 Note: If no refinement was done, the return value is NULL 1411 1412 Level: developer 1413 1414 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1415 @*/ 1416 PetscErrorCode DMRefine(DM dm,MPI_Comm comm,DM *dmf) 1417 { 1418 PetscErrorCode ierr; 1419 DMRefineHookLink link; 1420 1421 PetscFunctionBegin; 1422 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1423 ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr); 1424 if (*dmf) { 1425 (*dmf)->ops->creatematrix = dm->ops->creatematrix; 1426 1427 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr); 1428 1429 (*dmf)->ctx = dm->ctx; 1430 (*dmf)->leveldown = dm->leveldown; 1431 (*dmf)->levelup = dm->levelup + 1; 1432 1433 ierr = DMSetMatType(*dmf,dm->mattype);CHKERRQ(ierr); 1434 for (link=dm->refinehook; link; link=link->next) { 1435 if (link->refinehook) { 1436 ierr = (*link->refinehook)(dm,*dmf,link->ctx);CHKERRQ(ierr); 1437 } 1438 } 1439 } 1440 PetscFunctionReturn(0); 1441 } 1442 1443 #undef __FUNCT__ 1444 #define __FUNCT__ "DMRefineHookAdd" 1445 /*@C 1446 DMRefineHookAdd - adds a callback to be run when interpolating a nonlinear problem to a finer grid 1447 1448 Logically Collective 1449 1450 Input Arguments: 1451 + coarse - nonlinear solver context on which to run a hook when restricting to a coarser level 1452 . refinehook - function to run when setting up a coarser level 1453 . interphook - function to run to update data on finer levels (once per SNESSolve()) 1454 - ctx - [optional] user-defined context for provide data for the hooks (may be NULL) 1455 1456 Calling sequence of refinehook: 1457 $ refinehook(DM coarse,DM fine,void *ctx); 1458 1459 + coarse - coarse level DM 1460 . fine - fine level DM to interpolate problem to 1461 - ctx - optional user-defined function context 1462 1463 Calling sequence for interphook: 1464 $ interphook(DM coarse,Mat interp,DM fine,void *ctx) 1465 1466 + coarse - coarse level DM 1467 . interp - matrix interpolating a coarse-level solution to the finer grid 1468 . fine - fine level DM to update 1469 - ctx - optional user-defined function context 1470 1471 Level: advanced 1472 1473 Notes: 1474 This function is only needed if auxiliary data needs to be passed to fine grids while grid sequencing 1475 1476 If this function is called multiple times, the hooks will be run in the order they are added. 1477 1478 This function is currently not available from Fortran. 1479 1480 .seealso: DMCoarsenHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 1481 @*/ 1482 PetscErrorCode DMRefineHookAdd(DM coarse,PetscErrorCode (*refinehook)(DM,DM,void*),PetscErrorCode (*interphook)(DM,Mat,DM,void*),void *ctx) 1483 { 1484 PetscErrorCode ierr; 1485 DMRefineHookLink link,*p; 1486 1487 PetscFunctionBegin; 1488 PetscValidHeaderSpecific(coarse,DM_CLASSID,1); 1489 for (p=&coarse->refinehook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */ 1490 ierr = PetscMalloc(sizeof(struct _DMRefineHookLink),&link);CHKERRQ(ierr); 1491 link->refinehook = refinehook; 1492 link->interphook = interphook; 1493 link->ctx = ctx; 1494 link->next = NULL; 1495 *p = link; 1496 PetscFunctionReturn(0); 1497 } 1498 1499 #undef __FUNCT__ 1500 #define __FUNCT__ "DMInterpolate" 1501 /*@ 1502 DMInterpolate - interpolates user-defined problem data to a finer DM by running hooks registered by DMRefineHookAdd() 1503 1504 Collective if any hooks are 1505 1506 Input Arguments: 1507 + coarse - coarser DM to use as a base 1508 . restrct - interpolation matrix, apply using MatInterpolate() 1509 - fine - finer DM to update 1510 1511 Level: developer 1512 1513 .seealso: DMRefineHookAdd(), MatInterpolate() 1514 @*/ 1515 PetscErrorCode DMInterpolate(DM coarse,Mat interp,DM fine) 1516 { 1517 PetscErrorCode ierr; 1518 DMRefineHookLink link; 1519 1520 PetscFunctionBegin; 1521 for (link=fine->refinehook; link; link=link->next) { 1522 if (link->interphook) { 1523 ierr = (*link->interphook)(coarse,interp,fine,link->ctx);CHKERRQ(ierr); 1524 } 1525 } 1526 PetscFunctionReturn(0); 1527 } 1528 1529 #undef __FUNCT__ 1530 #define __FUNCT__ "DMGetRefineLevel" 1531 /*@ 1532 DMGetRefineLevel - Get's the number of refinements that have generated this DM. 1533 1534 Not Collective 1535 1536 Input Parameter: 1537 . dm - the DM object 1538 1539 Output Parameter: 1540 . level - number of refinements 1541 1542 Level: developer 1543 1544 .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1545 1546 @*/ 1547 PetscErrorCode DMGetRefineLevel(DM dm,PetscInt *level) 1548 { 1549 PetscFunctionBegin; 1550 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1551 *level = dm->levelup; 1552 PetscFunctionReturn(0); 1553 } 1554 1555 #undef __FUNCT__ 1556 #define __FUNCT__ "DMGlobalToLocalHookAdd" 1557 /*@C 1558 DMGlobalToLocalHookAdd - adds a callback to be run when global to local is called 1559 1560 Logically Collective 1561 1562 Input Arguments: 1563 + dm - the DM 1564 . beginhook - function to run at the beginning of DMGlobalToLocalBegin() 1565 . endhook - function to run after DMGlobalToLocalEnd() has completed 1566 - ctx - [optional] user-defined context for provide data for the hooks (may be NULL) 1567 1568 Calling sequence for beginhook: 1569 $ beginhook(DM fine,VecScatter out,VecScatter in,DM coarse,void *ctx) 1570 1571 + dm - global DM 1572 . g - global vector 1573 . mode - mode 1574 . l - local vector 1575 - ctx - optional user-defined function context 1576 1577 1578 Calling sequence for endhook: 1579 $ endhook(DM fine,VecScatter out,VecScatter in,DM coarse,void *ctx) 1580 1581 + global - global DM 1582 - ctx - optional user-defined function context 1583 1584 Level: advanced 1585 1586 .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 1587 @*/ 1588 PetscErrorCode DMGlobalToLocalHookAdd(DM dm,PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*),void *ctx) 1589 { 1590 PetscErrorCode ierr; 1591 DMGlobalToLocalHookLink link,*p; 1592 1593 PetscFunctionBegin; 1594 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1595 for (p=&dm->gtolhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */ 1596 ierr = PetscMalloc(sizeof(struct _DMGlobalToLocalHookLink),&link);CHKERRQ(ierr); 1597 link->beginhook = beginhook; 1598 link->endhook = endhook; 1599 link->ctx = ctx; 1600 link->next = NULL; 1601 *p = link; 1602 PetscFunctionReturn(0); 1603 } 1604 1605 #undef __FUNCT__ 1606 #define __FUNCT__ "DMGlobalToLocalBegin" 1607 /*@ 1608 DMGlobalToLocalBegin - Begins updating local vectors from global vector 1609 1610 Neighbor-wise Collective on DM 1611 1612 Input Parameters: 1613 + dm - the DM object 1614 . g - the global vector 1615 . mode - INSERT_VALUES or ADD_VALUES 1616 - l - the local vector 1617 1618 1619 Level: beginner 1620 1621 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 1622 1623 @*/ 1624 PetscErrorCode DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) 1625 { 1626 PetscSF sf; 1627 PetscErrorCode ierr; 1628 DMGlobalToLocalHookLink link; 1629 1630 PetscFunctionBegin; 1631 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1632 for (link=dm->gtolhook; link; link=link->next) { 1633 if (link->beginhook) { 1634 ierr = (*link->beginhook)(dm,g,mode,l,link->ctx);CHKERRQ(ierr); 1635 } 1636 } 1637 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 1638 if (sf) { 1639 PetscScalar *lArray, *gArray; 1640 1641 if (mode == ADD_VALUES) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 1642 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 1643 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 1644 ierr = PetscSFBcastBegin(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr); 1645 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 1646 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 1647 } else { 1648 ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 1649 } 1650 PetscFunctionReturn(0); 1651 } 1652 1653 #undef __FUNCT__ 1654 #define __FUNCT__ "DMGlobalToLocalEnd" 1655 /*@ 1656 DMGlobalToLocalEnd - Ends updating local vectors from global vector 1657 1658 Neighbor-wise Collective on DM 1659 1660 Input Parameters: 1661 + dm - the DM object 1662 . g - the global vector 1663 . mode - INSERT_VALUES or ADD_VALUES 1664 - l - the local vector 1665 1666 1667 Level: beginner 1668 1669 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 1670 1671 @*/ 1672 PetscErrorCode DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) 1673 { 1674 PetscSF sf; 1675 PetscErrorCode ierr; 1676 PetscScalar *lArray, *gArray; 1677 DMGlobalToLocalHookLink link; 1678 1679 PetscFunctionBegin; 1680 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1681 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 1682 if (sf) { 1683 if (mode == ADD_VALUES) SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 1684 1685 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 1686 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 1687 ierr = PetscSFBcastEnd(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr); 1688 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 1689 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 1690 } else { 1691 ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 1692 } 1693 for (link=dm->gtolhook; link; link=link->next) { 1694 if (link->endhook) {ierr = (*link->endhook)(dm,g,mode,l,link->ctx);CHKERRQ(ierr);} 1695 } 1696 PetscFunctionReturn(0); 1697 } 1698 1699 #undef __FUNCT__ 1700 #define __FUNCT__ "DMLocalToGlobalBegin" 1701 /*@ 1702 DMLocalToGlobalBegin - updates global vectors from local vectors 1703 1704 Neighbor-wise Collective on DM 1705 1706 Input Parameters: 1707 + dm - the DM object 1708 . l - the local vector 1709 . 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 1710 base point. 1711 - - the global vector 1712 1713 Notes: In the ADD_VALUES case you normally would zero the receiving vector before beginning this operation. If you would like to simply add the non-ghosted values in the local 1714 array into the global array you need to either (1) zero the ghosted locations and use ADD_VALUES or (2) use INSERT_VALUES into a work global array and then add the work 1715 global array to the final global array with VecAXPY(). 1716 1717 Level: beginner 1718 1719 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin() 1720 1721 @*/ 1722 PetscErrorCode DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g) 1723 { 1724 PetscSF sf; 1725 PetscErrorCode ierr; 1726 1727 PetscFunctionBegin; 1728 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1729 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 1730 if (sf) { 1731 MPI_Op op; 1732 PetscScalar *lArray, *gArray; 1733 1734 switch (mode) { 1735 case INSERT_VALUES: 1736 case INSERT_ALL_VALUES: 1737 op = MPIU_REPLACE; break; 1738 case ADD_VALUES: 1739 case ADD_ALL_VALUES: 1740 op = MPI_SUM; break; 1741 default: 1742 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 1743 } 1744 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 1745 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 1746 ierr = PetscSFReduceBegin(sf, MPIU_SCALAR, lArray, gArray, op);CHKERRQ(ierr); 1747 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 1748 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 1749 } else { 1750 ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 1751 } 1752 PetscFunctionReturn(0); 1753 } 1754 1755 #undef __FUNCT__ 1756 #define __FUNCT__ "DMLocalToGlobalEnd" 1757 /*@ 1758 DMLocalToGlobalEnd - updates global vectors from local vectors 1759 1760 Neighbor-wise Collective on DM 1761 1762 Input Parameters: 1763 + dm - the DM object 1764 . l - the local vector 1765 . mode - INSERT_VALUES or ADD_VALUES 1766 - g - the global vector 1767 1768 1769 Level: beginner 1770 1771 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd() 1772 1773 @*/ 1774 PetscErrorCode DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g) 1775 { 1776 PetscSF sf; 1777 PetscErrorCode ierr; 1778 1779 PetscFunctionBegin; 1780 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1781 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 1782 if (sf) { 1783 MPI_Op op; 1784 PetscScalar *lArray, *gArray; 1785 1786 switch (mode) { 1787 case INSERT_VALUES: 1788 case INSERT_ALL_VALUES: 1789 op = MPIU_REPLACE; break; 1790 case ADD_VALUES: 1791 case ADD_ALL_VALUES: 1792 op = MPI_SUM; break; 1793 default: 1794 SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 1795 } 1796 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 1797 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 1798 ierr = PetscSFReduceEnd(sf, MPIU_SCALAR, lArray, gArray, op);CHKERRQ(ierr); 1799 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 1800 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 1801 } else { 1802 ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 1803 } 1804 PetscFunctionReturn(0); 1805 } 1806 1807 #undef __FUNCT__ 1808 #define __FUNCT__ "DMLocalToLocalBegin" 1809 /*@ 1810 DMLocalToLocalBegin - Maps from a local vector (including ghost points 1811 that contain irrelevant values) to another local vector where the ghost 1812 points in the second are set correctly. Must be followed by DMLocalToLocalEnd(). 1813 1814 Neighbor-wise Collective on DM and Vec 1815 1816 Input Parameters: 1817 + dm - the DM object 1818 . g - the original local vector 1819 - mode - one of INSERT_VALUES or ADD_VALUES 1820 1821 Output Parameter: 1822 . l - the local vector with correct ghost values 1823 1824 Level: intermediate 1825 1826 Notes: 1827 The local vectors used here need not be the same as those 1828 obtained from DMCreateLocalVector(), BUT they 1829 must have the same parallel data layout; they could, for example, be 1830 obtained with VecDuplicate() from the DM originating vectors. 1831 1832 .keywords: DM, local-to-local, begin 1833 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateLocalVector(), DMCreateGlobalVector(), DMCreateInterpolation(), DMLocalToLocalEnd(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 1834 1835 @*/ 1836 PetscErrorCode DMLocalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) 1837 { 1838 PetscErrorCode ierr; 1839 1840 PetscFunctionBegin; 1841 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1842 ierr = (*dm->ops->localtolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 1843 PetscFunctionReturn(0); 1844 } 1845 1846 #undef __FUNCT__ 1847 #define __FUNCT__ "DMLocalToLocalEnd" 1848 /*@ 1849 DMLocalToLocalEnd - Maps from a local vector (including ghost points 1850 that contain irrelevant values) to another local vector where the ghost 1851 points in the second are set correctly. Must be preceded by DMLocalToLocalBegin(). 1852 1853 Neighbor-wise Collective on DM and Vec 1854 1855 Input Parameters: 1856 + da - the DM object 1857 . g - the original local vector 1858 - mode - one of INSERT_VALUES or ADD_VALUES 1859 1860 Output Parameter: 1861 . l - the local vector with correct ghost values 1862 1863 Level: intermediate 1864 1865 Notes: 1866 The local vectors used here need not be the same as those 1867 obtained from DMCreateLocalVector(), BUT they 1868 must have the same parallel data layout; they could, for example, be 1869 obtained with VecDuplicate() from the DM originating vectors. 1870 1871 .keywords: DM, local-to-local, end 1872 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateLocalVector(), DMCreateGlobalVector(), DMCreateInterpolation(), DMLocalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 1873 1874 @*/ 1875 PetscErrorCode DMLocalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) 1876 { 1877 PetscErrorCode ierr; 1878 1879 PetscFunctionBegin; 1880 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1881 ierr = (*dm->ops->localtolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 1882 PetscFunctionReturn(0); 1883 } 1884 1885 1886 #undef __FUNCT__ 1887 #define __FUNCT__ "DMCoarsen" 1888 /*@ 1889 DMCoarsen - Coarsens a DM object 1890 1891 Collective on DM 1892 1893 Input Parameter: 1894 + dm - the DM object 1895 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 1896 1897 Output Parameter: 1898 . dmc - the coarsened DM 1899 1900 Level: developer 1901 1902 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1903 1904 @*/ 1905 PetscErrorCode DMCoarsen(DM dm, MPI_Comm comm, DM *dmc) 1906 { 1907 PetscErrorCode ierr; 1908 DMCoarsenHookLink link; 1909 1910 PetscFunctionBegin; 1911 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1912 ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr); 1913 (*dmc)->ops->creatematrix = dm->ops->creatematrix; 1914 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr); 1915 (*dmc)->ctx = dm->ctx; 1916 (*dmc)->levelup = dm->levelup; 1917 (*dmc)->leveldown = dm->leveldown + 1; 1918 ierr = DMSetMatType(*dmc,dm->mattype);CHKERRQ(ierr); 1919 for (link=dm->coarsenhook; link; link=link->next) { 1920 if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);} 1921 } 1922 PetscFunctionReturn(0); 1923 } 1924 1925 #undef __FUNCT__ 1926 #define __FUNCT__ "DMCoarsenHookAdd" 1927 /*@C 1928 DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid 1929 1930 Logically Collective 1931 1932 Input Arguments: 1933 + fine - nonlinear solver context on which to run a hook when restricting to a coarser level 1934 . coarsenhook - function to run when setting up a coarser level 1935 . restricthook - function to run to update data on coarser levels (once per SNESSolve()) 1936 - ctx - [optional] user-defined context for provide data for the hooks (may be NULL) 1937 1938 Calling sequence of coarsenhook: 1939 $ coarsenhook(DM fine,DM coarse,void *ctx); 1940 1941 + fine - fine level DM 1942 . coarse - coarse level DM to restrict problem to 1943 - ctx - optional user-defined function context 1944 1945 Calling sequence for restricthook: 1946 $ restricthook(DM fine,Mat mrestrict,Vec rscale,Mat inject,DM coarse,void *ctx) 1947 1948 + fine - fine level DM 1949 . mrestrict - matrix restricting a fine-level solution to the coarse grid 1950 . rscale - scaling vector for restriction 1951 . inject - matrix restricting by injection 1952 . coarse - coarse level DM to update 1953 - ctx - optional user-defined function context 1954 1955 Level: advanced 1956 1957 Notes: 1958 This function is only needed if auxiliary data needs to be set up on coarse grids. 1959 1960 If this function is called multiple times, the hooks will be run in the order they are added. 1961 1962 In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to 1963 extract the finest level information from its context (instead of from the SNES). 1964 1965 This function is currently not available from Fortran. 1966 1967 .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 1968 @*/ 1969 PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx) 1970 { 1971 PetscErrorCode ierr; 1972 DMCoarsenHookLink link,*p; 1973 1974 PetscFunctionBegin; 1975 PetscValidHeaderSpecific(fine,DM_CLASSID,1); 1976 for (p=&fine->coarsenhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */ 1977 ierr = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr); 1978 link->coarsenhook = coarsenhook; 1979 link->restricthook = restricthook; 1980 link->ctx = ctx; 1981 link->next = NULL; 1982 *p = link; 1983 PetscFunctionReturn(0); 1984 } 1985 1986 #undef __FUNCT__ 1987 #define __FUNCT__ "DMRestrict" 1988 /*@ 1989 DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd() 1990 1991 Collective if any hooks are 1992 1993 Input Arguments: 1994 + fine - finer DM to use as a base 1995 . restrct - restriction matrix, apply using MatRestrict() 1996 . inject - injection matrix, also use MatRestrict() 1997 - coarse - coarer DM to update 1998 1999 Level: developer 2000 2001 .seealso: DMCoarsenHookAdd(), MatRestrict() 2002 @*/ 2003 PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse) 2004 { 2005 PetscErrorCode ierr; 2006 DMCoarsenHookLink link; 2007 2008 PetscFunctionBegin; 2009 for (link=fine->coarsenhook; link; link=link->next) { 2010 if (link->restricthook) { 2011 ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr); 2012 } 2013 } 2014 PetscFunctionReturn(0); 2015 } 2016 2017 #undef __FUNCT__ 2018 #define __FUNCT__ "DMSubDomainHookAdd" 2019 /*@C 2020 DMSubDomainHookAdd - adds a callback to be run when restricting a problem to the coarse grid 2021 2022 Logically Collective 2023 2024 Input Arguments: 2025 + global - global DM 2026 . ddhook - function to run to pass data to the decomposition DM upon its creation 2027 . restricthook - function to run to update data on block solve (at the beginning of the block solve) 2028 - ctx - [optional] user-defined context for provide data for the hooks (may be NULL) 2029 2030 2031 Calling sequence for ddhook: 2032 $ ddhook(DM global,DM block,void *ctx) 2033 2034 + global - global DM 2035 . block - block DM 2036 - ctx - optional user-defined function context 2037 2038 Calling sequence for restricthook: 2039 $ restricthook(DM global,VecScatter out,VecScatter in,DM block,void *ctx) 2040 2041 + global - global DM 2042 . out - scatter to the outer (with ghost and overlap points) block vector 2043 . in - scatter to block vector values only owned locally 2044 . block - block DM 2045 - ctx - optional user-defined function context 2046 2047 Level: advanced 2048 2049 Notes: 2050 This function is only needed if auxiliary data needs to be set up on subdomain DMs. 2051 2052 If this function is called multiple times, the hooks will be run in the order they are added. 2053 2054 In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to 2055 extract the global information from its context (instead of from the SNES). 2056 2057 This function is currently not available from Fortran. 2058 2059 .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 2060 @*/ 2061 PetscErrorCode DMSubDomainHookAdd(DM global,PetscErrorCode (*ddhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*),void *ctx) 2062 { 2063 PetscErrorCode ierr; 2064 DMSubDomainHookLink link,*p; 2065 2066 PetscFunctionBegin; 2067 PetscValidHeaderSpecific(global,DM_CLASSID,1); 2068 for (p=&global->subdomainhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */ 2069 ierr = PetscMalloc(sizeof(struct _DMSubDomainHookLink),&link);CHKERRQ(ierr); 2070 link->restricthook = restricthook; 2071 link->ddhook = ddhook; 2072 link->ctx = ctx; 2073 link->next = NULL; 2074 *p = link; 2075 PetscFunctionReturn(0); 2076 } 2077 2078 #undef __FUNCT__ 2079 #define __FUNCT__ "DMSubDomainRestrict" 2080 /*@ 2081 DMSubDomainRestrict - restricts user-defined problem data to a block DM by running hooks registered by DMSubDomainHookAdd() 2082 2083 Collective if any hooks are 2084 2085 Input Arguments: 2086 + fine - finer DM to use as a base 2087 . oscatter - scatter from domain global vector filling subdomain global vector with overlap 2088 . gscatter - scatter from domain global vector filling subdomain local vector with ghosts 2089 - coarse - coarer DM to update 2090 2091 Level: developer 2092 2093 .seealso: DMCoarsenHookAdd(), MatRestrict() 2094 @*/ 2095 PetscErrorCode DMSubDomainRestrict(DM global,VecScatter oscatter,VecScatter gscatter,DM subdm) 2096 { 2097 PetscErrorCode ierr; 2098 DMSubDomainHookLink link; 2099 2100 PetscFunctionBegin; 2101 for (link=global->subdomainhook; link; link=link->next) { 2102 if (link->restricthook) { 2103 ierr = (*link->restricthook)(global,oscatter,gscatter,subdm,link->ctx);CHKERRQ(ierr); 2104 } 2105 } 2106 PetscFunctionReturn(0); 2107 } 2108 2109 #undef __FUNCT__ 2110 #define __FUNCT__ "DMGetCoarsenLevel" 2111 /*@ 2112 DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM. 2113 2114 Not Collective 2115 2116 Input Parameter: 2117 . dm - the DM object 2118 2119 Output Parameter: 2120 . level - number of coarsenings 2121 2122 Level: developer 2123 2124 .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 2125 2126 @*/ 2127 PetscErrorCode DMGetCoarsenLevel(DM dm,PetscInt *level) 2128 { 2129 PetscFunctionBegin; 2130 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2131 *level = dm->leveldown; 2132 PetscFunctionReturn(0); 2133 } 2134 2135 2136 2137 #undef __FUNCT__ 2138 #define __FUNCT__ "DMRefineHierarchy" 2139 /*@C 2140 DMRefineHierarchy - Refines a DM object, all levels at once 2141 2142 Collective on DM 2143 2144 Input Parameter: 2145 + dm - the DM object 2146 - nlevels - the number of levels of refinement 2147 2148 Output Parameter: 2149 . dmf - the refined DM hierarchy 2150 2151 Level: developer 2152 2153 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 2154 2155 @*/ 2156 PetscErrorCode DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[]) 2157 { 2158 PetscErrorCode ierr; 2159 2160 PetscFunctionBegin; 2161 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2162 if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 2163 if (nlevels == 0) PetscFunctionReturn(0); 2164 if (dm->ops->refinehierarchy) { 2165 ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr); 2166 } else if (dm->ops->refine) { 2167 PetscInt i; 2168 2169 ierr = DMRefine(dm,PetscObjectComm((PetscObject)dm),&dmf[0]);CHKERRQ(ierr); 2170 for (i=1; i<nlevels; i++) { 2171 ierr = DMRefine(dmf[i-1],PetscObjectComm((PetscObject)dm),&dmf[i]);CHKERRQ(ierr); 2172 } 2173 } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No RefineHierarchy for this DM yet"); 2174 PetscFunctionReturn(0); 2175 } 2176 2177 #undef __FUNCT__ 2178 #define __FUNCT__ "DMCoarsenHierarchy" 2179 /*@C 2180 DMCoarsenHierarchy - Coarsens a DM object, all levels at once 2181 2182 Collective on DM 2183 2184 Input Parameter: 2185 + dm - the DM object 2186 - nlevels - the number of levels of coarsening 2187 2188 Output Parameter: 2189 . dmc - the coarsened DM hierarchy 2190 2191 Level: developer 2192 2193 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 2194 2195 @*/ 2196 PetscErrorCode DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[]) 2197 { 2198 PetscErrorCode ierr; 2199 2200 PetscFunctionBegin; 2201 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2202 if (nlevels < 0) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 2203 if (nlevels == 0) PetscFunctionReturn(0); 2204 PetscValidPointer(dmc,3); 2205 if (dm->ops->coarsenhierarchy) { 2206 ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr); 2207 } else if (dm->ops->coarsen) { 2208 PetscInt i; 2209 2210 ierr = DMCoarsen(dm,PetscObjectComm((PetscObject)dm),&dmc[0]);CHKERRQ(ierr); 2211 for (i=1; i<nlevels; i++) { 2212 ierr = DMCoarsen(dmc[i-1],PetscObjectComm((PetscObject)dm),&dmc[i]);CHKERRQ(ierr); 2213 } 2214 } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet"); 2215 PetscFunctionReturn(0); 2216 } 2217 2218 #undef __FUNCT__ 2219 #define __FUNCT__ "DMCreateAggregates" 2220 /*@ 2221 DMCreateAggregates - Gets the aggregates that map between 2222 grids associated with two DMs. 2223 2224 Collective on DM 2225 2226 Input Parameters: 2227 + dmc - the coarse grid DM 2228 - dmf - the fine grid DM 2229 2230 Output Parameters: 2231 . rest - the restriction matrix (transpose of the projection matrix) 2232 2233 Level: intermediate 2234 2235 .keywords: interpolation, restriction, multigrid 2236 2237 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation() 2238 @*/ 2239 PetscErrorCode DMCreateAggregates(DM dmc, DM dmf, Mat *rest) 2240 { 2241 PetscErrorCode ierr; 2242 2243 PetscFunctionBegin; 2244 PetscValidHeaderSpecific(dmc,DM_CLASSID,1); 2245 PetscValidHeaderSpecific(dmf,DM_CLASSID,2); 2246 ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr); 2247 PetscFunctionReturn(0); 2248 } 2249 2250 #undef __FUNCT__ 2251 #define __FUNCT__ "DMSetApplicationContextDestroy" 2252 /*@C 2253 DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed 2254 2255 Not Collective 2256 2257 Input Parameters: 2258 + dm - the DM object 2259 - destroy - the destroy function 2260 2261 Level: intermediate 2262 2263 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 2264 2265 @*/ 2266 PetscErrorCode DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**)) 2267 { 2268 PetscFunctionBegin; 2269 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2270 dm->ctxdestroy = destroy; 2271 PetscFunctionReturn(0); 2272 } 2273 2274 #undef __FUNCT__ 2275 #define __FUNCT__ "DMSetApplicationContext" 2276 /*@ 2277 DMSetApplicationContext - Set a user context into a DM object 2278 2279 Not Collective 2280 2281 Input Parameters: 2282 + dm - the DM object 2283 - ctx - the user context 2284 2285 Level: intermediate 2286 2287 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 2288 2289 @*/ 2290 PetscErrorCode DMSetApplicationContext(DM dm,void *ctx) 2291 { 2292 PetscFunctionBegin; 2293 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2294 dm->ctx = ctx; 2295 PetscFunctionReturn(0); 2296 } 2297 2298 #undef __FUNCT__ 2299 #define __FUNCT__ "DMGetApplicationContext" 2300 /*@ 2301 DMGetApplicationContext - Gets a user context from a DM object 2302 2303 Not Collective 2304 2305 Input Parameter: 2306 . dm - the DM object 2307 2308 Output Parameter: 2309 . ctx - the user context 2310 2311 Level: intermediate 2312 2313 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 2314 2315 @*/ 2316 PetscErrorCode DMGetApplicationContext(DM dm,void *ctx) 2317 { 2318 PetscFunctionBegin; 2319 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2320 *(void**)ctx = dm->ctx; 2321 PetscFunctionReturn(0); 2322 } 2323 2324 #undef __FUNCT__ 2325 #define __FUNCT__ "DMSetVariableBounds" 2326 /*@C 2327 DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI. 2328 2329 Logically Collective on DM 2330 2331 Input Parameter: 2332 + dm - the DM object 2333 - f - the function that computes variable bounds used by SNESVI (use NULL to cancel a previous function that was set) 2334 2335 Level: intermediate 2336 2337 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), 2338 DMSetJacobian() 2339 2340 @*/ 2341 PetscErrorCode DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 2342 { 2343 PetscFunctionBegin; 2344 dm->ops->computevariablebounds = f; 2345 PetscFunctionReturn(0); 2346 } 2347 2348 #undef __FUNCT__ 2349 #define __FUNCT__ "DMHasVariableBounds" 2350 /*@ 2351 DMHasVariableBounds - does the DM object have a variable bounds function? 2352 2353 Not Collective 2354 2355 Input Parameter: 2356 . dm - the DM object to destroy 2357 2358 Output Parameter: 2359 . flg - PETSC_TRUE if the variable bounds function exists 2360 2361 Level: developer 2362 2363 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 2364 2365 @*/ 2366 PetscErrorCode DMHasVariableBounds(DM dm,PetscBool *flg) 2367 { 2368 PetscFunctionBegin; 2369 *flg = (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE; 2370 PetscFunctionReturn(0); 2371 } 2372 2373 #undef __FUNCT__ 2374 #define __FUNCT__ "DMComputeVariableBounds" 2375 /*@C 2376 DMComputeVariableBounds - compute variable bounds used by SNESVI. 2377 2378 Logically Collective on DM 2379 2380 Input Parameters: 2381 + dm - the DM object to destroy 2382 - x - current solution at which the bounds are computed 2383 2384 Output parameters: 2385 + xl - lower bound 2386 - xu - upper bound 2387 2388 Level: intermediate 2389 2390 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 2391 2392 @*/ 2393 PetscErrorCode DMComputeVariableBounds(DM dm, Vec xl, Vec xu) 2394 { 2395 PetscErrorCode ierr; 2396 2397 PetscFunctionBegin; 2398 PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 2399 PetscValidHeaderSpecific(xu,VEC_CLASSID,2); 2400 if (dm->ops->computevariablebounds) { 2401 ierr = (*dm->ops->computevariablebounds)(dm, xl,xu);CHKERRQ(ierr); 2402 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "This DM is incapable of computing variable bounds."); 2403 PetscFunctionReturn(0); 2404 } 2405 2406 #undef __FUNCT__ 2407 #define __FUNCT__ "DMHasColoring" 2408 /*@ 2409 DMHasColoring - does the DM object have a method of providing a coloring? 2410 2411 Not Collective 2412 2413 Input Parameter: 2414 . dm - the DM object 2415 2416 Output Parameter: 2417 . flg - PETSC_TRUE if the DM has facilities for DMCreateColoring(). 2418 2419 Level: developer 2420 2421 .seealso DMHasFunction(), DMCreateColoring() 2422 2423 @*/ 2424 PetscErrorCode DMHasColoring(DM dm,PetscBool *flg) 2425 { 2426 PetscFunctionBegin; 2427 *flg = (dm->ops->getcoloring) ? PETSC_TRUE : PETSC_FALSE; 2428 PetscFunctionReturn(0); 2429 } 2430 2431 #undef __FUNCT__ 2432 #define __FUNCT__ "DMSetVec" 2433 /*@C 2434 DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear. 2435 2436 Collective on DM 2437 2438 Input Parameter: 2439 + dm - the DM object 2440 - x - location to compute residual and Jacobian, if NULL is passed to those routines; will be NULL for linear problems. 2441 2442 Level: developer 2443 2444 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 2445 2446 @*/ 2447 PetscErrorCode DMSetVec(DM dm,Vec x) 2448 { 2449 PetscErrorCode ierr; 2450 2451 PetscFunctionBegin; 2452 if (x) { 2453 if (!dm->x) { 2454 ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 2455 } 2456 ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 2457 } else if (dm->x) { 2458 ierr = VecDestroy(&dm->x);CHKERRQ(ierr); 2459 } 2460 PetscFunctionReturn(0); 2461 } 2462 2463 PetscFunctionList DMList = NULL; 2464 PetscBool DMRegisterAllCalled = PETSC_FALSE; 2465 2466 #undef __FUNCT__ 2467 #define __FUNCT__ "DMSetType" 2468 /*@C 2469 DMSetType - Builds a DM, for a particular DM implementation. 2470 2471 Collective on DM 2472 2473 Input Parameters: 2474 + dm - The DM object 2475 - method - The name of the DM type 2476 2477 Options Database Key: 2478 . -dm_type <type> - Sets the DM type; use -help for a list of available types 2479 2480 Notes: 2481 See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D). 2482 2483 Level: intermediate 2484 2485 .keywords: DM, set, type 2486 .seealso: DMGetType(), DMCreate() 2487 @*/ 2488 PetscErrorCode DMSetType(DM dm, DMType method) 2489 { 2490 PetscErrorCode (*r)(DM); 2491 PetscBool match; 2492 PetscErrorCode ierr; 2493 2494 PetscFunctionBegin; 2495 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 2496 ierr = PetscObjectTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr); 2497 if (match) PetscFunctionReturn(0); 2498 2499 if (!DMRegisterAllCalled) {ierr = DMRegisterAll();CHKERRQ(ierr);} 2500 ierr = PetscFunctionListFind(DMList,method,&r);CHKERRQ(ierr); 2501 if (!r) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method); 2502 2503 if (dm->ops->destroy) { 2504 ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr); 2505 dm->ops->destroy = NULL; 2506 } 2507 ierr = (*r)(dm);CHKERRQ(ierr); 2508 ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr); 2509 PetscFunctionReturn(0); 2510 } 2511 2512 #undef __FUNCT__ 2513 #define __FUNCT__ "DMGetType" 2514 /*@C 2515 DMGetType - Gets the DM type name (as a string) from the DM. 2516 2517 Not Collective 2518 2519 Input Parameter: 2520 . dm - The DM 2521 2522 Output Parameter: 2523 . type - The DM type name 2524 2525 Level: intermediate 2526 2527 .keywords: DM, get, type, name 2528 .seealso: DMSetType(), DMCreate() 2529 @*/ 2530 PetscErrorCode DMGetType(DM dm, DMType *type) 2531 { 2532 PetscErrorCode ierr; 2533 2534 PetscFunctionBegin; 2535 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 2536 PetscValidCharPointer(type,2); 2537 if (!DMRegisterAllCalled) { 2538 ierr = DMRegisterAll();CHKERRQ(ierr); 2539 } 2540 *type = ((PetscObject)dm)->type_name; 2541 PetscFunctionReturn(0); 2542 } 2543 2544 #undef __FUNCT__ 2545 #define __FUNCT__ "DMConvert" 2546 /*@C 2547 DMConvert - Converts a DM to another DM, either of the same or different type. 2548 2549 Collective on DM 2550 2551 Input Parameters: 2552 + dm - the DM 2553 - newtype - new DM type (use "same" for the same type) 2554 2555 Output Parameter: 2556 . M - pointer to new DM 2557 2558 Notes: 2559 Cannot be used to convert a sequential DM to parallel or parallel to sequential, 2560 the MPI communicator of the generated DM is always the same as the communicator 2561 of the input DM. 2562 2563 Level: intermediate 2564 2565 .seealso: DMCreate() 2566 @*/ 2567 PetscErrorCode DMConvert(DM dm, DMType newtype, DM *M) 2568 { 2569 DM B; 2570 char convname[256]; 2571 PetscBool sametype, issame; 2572 PetscErrorCode ierr; 2573 2574 PetscFunctionBegin; 2575 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2576 PetscValidType(dm,1); 2577 PetscValidPointer(M,3); 2578 ierr = PetscObjectTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr); 2579 ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr); 2580 { 2581 PetscErrorCode (*conv)(DM, DMType, DM*) = NULL; 2582 2583 /* 2584 Order of precedence: 2585 1) See if a specialized converter is known to the current DM. 2586 2) See if a specialized converter is known to the desired DM class. 2587 3) See if a good general converter is registered for the desired class 2588 4) See if a good general converter is known for the current matrix. 2589 5) Use a really basic converter. 2590 */ 2591 2592 /* 1) See if a specialized converter is known to the current DM and the desired class */ 2593 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 2594 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 2595 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 2596 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 2597 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 2598 ierr = PetscObjectQueryFunction((PetscObject)dm,convname,&conv);CHKERRQ(ierr); 2599 if (conv) goto foundconv; 2600 2601 /* 2) See if a specialized converter is known to the desired DM class. */ 2602 ierr = DMCreate(PetscObjectComm((PetscObject)dm), &B);CHKERRQ(ierr); 2603 ierr = DMSetType(B, newtype);CHKERRQ(ierr); 2604 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 2605 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 2606 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 2607 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 2608 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 2609 ierr = PetscObjectQueryFunction((PetscObject)B,convname,&conv);CHKERRQ(ierr); 2610 if (conv) { 2611 ierr = DMDestroy(&B);CHKERRQ(ierr); 2612 goto foundconv; 2613 } 2614 2615 #if 0 2616 /* 3) See if a good general converter is registered for the desired class */ 2617 conv = B->ops->convertfrom; 2618 ierr = DMDestroy(&B);CHKERRQ(ierr); 2619 if (conv) goto foundconv; 2620 2621 /* 4) See if a good general converter is known for the current matrix */ 2622 if (dm->ops->convert) { 2623 conv = dm->ops->convert; 2624 } 2625 if (conv) goto foundconv; 2626 #endif 2627 2628 /* 5) Use a really basic converter. */ 2629 SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype); 2630 2631 foundconv: 2632 ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 2633 ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr); 2634 ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 2635 } 2636 ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr); 2637 PetscFunctionReturn(0); 2638 } 2639 2640 /*--------------------------------------------------------------------------------------------------------------------*/ 2641 2642 #undef __FUNCT__ 2643 #define __FUNCT__ "DMRegister" 2644 /*@C 2645 DMRegister - Adds a new DM component implementation 2646 2647 Not Collective 2648 2649 Input Parameters: 2650 + name - The name of a new user-defined creation routine 2651 - create_func - The creation routine itself 2652 2653 Notes: 2654 DMRegister() may be called multiple times to add several user-defined DMs 2655 2656 2657 Sample usage: 2658 .vb 2659 DMRegister("my_da", MyDMCreate); 2660 .ve 2661 2662 Then, your DM type can be chosen with the procedural interface via 2663 .vb 2664 DMCreate(MPI_Comm, DM *); 2665 DMSetType(DM,"my_da"); 2666 .ve 2667 or at runtime via the option 2668 .vb 2669 -da_type my_da 2670 .ve 2671 2672 Level: advanced 2673 2674 .keywords: DM, register 2675 .seealso: DMRegisterAll(), DMRegisterDestroy() 2676 2677 @*/ 2678 PetscErrorCode DMRegister(const char sname[],PetscErrorCode (*function)(DM)) 2679 { 2680 PetscErrorCode ierr; 2681 2682 PetscFunctionBegin; 2683 ierr = PetscFunctionListAdd(&DMList,sname,function);CHKERRQ(ierr); 2684 PetscFunctionReturn(0); 2685 } 2686 2687 #undef __FUNCT__ 2688 #define __FUNCT__ "DMLoad" 2689 /*@C 2690 DMLoad - Loads a DM that has been stored in binary with DMView(). 2691 2692 Collective on PetscViewer 2693 2694 Input Parameters: 2695 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or 2696 some related function before a call to DMLoad(). 2697 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 2698 HDF5 file viewer, obtained from PetscViewerHDF5Open() 2699 2700 Level: intermediate 2701 2702 Notes: 2703 The type is determined by the data in the file, any type set into the DM before this call is ignored. 2704 2705 Notes for advanced users: 2706 Most users should not need to know the details of the binary storage 2707 format, since DMLoad() and DMView() completely hide these details. 2708 But for anyone who's interested, the standard binary matrix storage 2709 format is 2710 .vb 2711 has not yet been determined 2712 .ve 2713 2714 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad() 2715 @*/ 2716 PetscErrorCode DMLoad(DM newdm, PetscViewer viewer) 2717 { 2718 PetscBool isbinary, ishdf5; 2719 PetscErrorCode ierr; 2720 2721 PetscFunctionBegin; 2722 PetscValidHeaderSpecific(newdm,DM_CLASSID,1); 2723 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2724 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 2725 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 2726 if (isbinary) { 2727 PetscInt classid; 2728 char type[256]; 2729 2730 ierr = PetscViewerBinaryRead(viewer,&classid,1,PETSC_INT);CHKERRQ(ierr); 2731 if (classid != DM_FILE_CLASSID) SETERRQ1(PetscObjectComm((PetscObject)newdm),PETSC_ERR_ARG_WRONG,"Not DM next in file, classid found %d",(int)classid); 2732 ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr); 2733 ierr = DMSetType(newdm, type);CHKERRQ(ierr); 2734 if (newdm->ops->load) {ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);} 2735 } else if (ishdf5) { 2736 if (newdm->ops->load) {ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr);} 2737 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen() or PetscViewerHDF5Open()"); 2738 PetscFunctionReturn(0); 2739 } 2740 2741 /******************************** FEM Support **********************************/ 2742 2743 #undef __FUNCT__ 2744 #define __FUNCT__ "DMPrintCellVector" 2745 PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[]) 2746 { 2747 PetscInt f; 2748 PetscErrorCode ierr; 2749 2750 PetscFunctionBegin; 2751 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 2752 for (f = 0; f < len; ++f) { 2753 ierr = PetscPrintf(PETSC_COMM_SELF, " | %g |\n", (double)PetscRealPart(x[f]));CHKERRQ(ierr); 2754 } 2755 PetscFunctionReturn(0); 2756 } 2757 2758 #undef __FUNCT__ 2759 #define __FUNCT__ "DMPrintCellMatrix" 2760 PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[]) 2761 { 2762 PetscInt f, g; 2763 PetscErrorCode ierr; 2764 2765 PetscFunctionBegin; 2766 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 2767 for (f = 0; f < rows; ++f) { 2768 ierr = PetscPrintf(PETSC_COMM_SELF, " |");CHKERRQ(ierr); 2769 for (g = 0; g < cols; ++g) { 2770 ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5g", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr); 2771 } 2772 ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr); 2773 } 2774 PetscFunctionReturn(0); 2775 } 2776 2777 #undef __FUNCT__ 2778 #define __FUNCT__ "DMPrintLocalVec" 2779 PetscErrorCode DMPrintLocalVec(DM dm, const char name[], PetscReal tol, Vec X) 2780 { 2781 PetscMPIInt rank, numProcs; 2782 PetscInt p; 2783 PetscErrorCode ierr; 2784 2785 PetscFunctionBegin; 2786 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 2787 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr); 2788 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "%s:\n", name);CHKERRQ(ierr); 2789 for (p = 0; p < numProcs; ++p) { 2790 if (p == rank) { 2791 Vec x; 2792 2793 ierr = VecDuplicate(X, &x);CHKERRQ(ierr); 2794 ierr = VecCopy(X, x);CHKERRQ(ierr); 2795 ierr = VecChop(x, tol);CHKERRQ(ierr); 2796 ierr = VecView(x, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 2797 ierr = VecDestroy(&x);CHKERRQ(ierr); 2798 ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 2799 } 2800 ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 2801 } 2802 PetscFunctionReturn(0); 2803 } 2804 2805 #undef __FUNCT__ 2806 #define __FUNCT__ "DMGetDefaultSection" 2807 /*@ 2808 DMGetDefaultSection - Get the PetscSection encoding the local data layout for the DM. 2809 2810 Input Parameter: 2811 . dm - The DM 2812 2813 Output Parameter: 2814 . section - The PetscSection 2815 2816 Level: intermediate 2817 2818 Note: This gets a borrowed reference, so the user should not destroy this PetscSection. 2819 2820 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection() 2821 @*/ 2822 PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *section) 2823 { 2824 PetscErrorCode ierr; 2825 2826 PetscFunctionBegin; 2827 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2828 PetscValidPointer(section, 2); 2829 if (!dm->defaultSection && dm->ops->createdefaultsection) {ierr = (*dm->ops->createdefaultsection)(dm);CHKERRQ(ierr);} 2830 *section = dm->defaultSection; 2831 PetscFunctionReturn(0); 2832 } 2833 2834 #undef __FUNCT__ 2835 #define __FUNCT__ "DMSetDefaultSection" 2836 /*@ 2837 DMSetDefaultSection - Set the PetscSection encoding the local data layout for the DM. 2838 2839 Input Parameters: 2840 + dm - The DM 2841 - section - The PetscSection 2842 2843 Level: intermediate 2844 2845 Note: Any existing Section will be destroyed 2846 2847 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection() 2848 @*/ 2849 PetscErrorCode DMSetDefaultSection(DM dm, PetscSection section) 2850 { 2851 PetscInt numFields; 2852 PetscInt f; 2853 PetscErrorCode ierr; 2854 2855 PetscFunctionBegin; 2856 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2857 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2); 2858 ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr); 2859 ierr = PetscSectionDestroy(&dm->defaultSection);CHKERRQ(ierr); 2860 dm->defaultSection = section; 2861 ierr = PetscSectionGetNumFields(dm->defaultSection, &numFields);CHKERRQ(ierr); 2862 if (numFields) { 2863 ierr = DMSetNumFields(dm, numFields);CHKERRQ(ierr); 2864 for (f = 0; f < numFields; ++f) { 2865 PetscObject disc; 2866 const char *name; 2867 2868 ierr = PetscSectionGetFieldName(dm->defaultSection, f, &name);CHKERRQ(ierr); 2869 ierr = DMGetField(dm, f, &disc);CHKERRQ(ierr); 2870 ierr = PetscObjectSetName(disc, name);CHKERRQ(ierr); 2871 } 2872 } 2873 /* The global section will be rebuilt in the next call to DMGetDefaultGlobalSection(). */ 2874 ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr); 2875 PetscFunctionReturn(0); 2876 } 2877 2878 #undef __FUNCT__ 2879 #define __FUNCT__ "DMGetDefaultGlobalSection" 2880 /*@ 2881 DMGetDefaultGlobalSection - Get the PetscSection encoding the global data layout for the DM. 2882 2883 Collective on DM 2884 2885 Input Parameter: 2886 . dm - The DM 2887 2888 Output Parameter: 2889 . section - The PetscSection 2890 2891 Level: intermediate 2892 2893 Note: This gets a borrowed reference, so the user should not destroy this PetscSection. 2894 2895 .seealso: DMSetDefaultSection(), DMGetDefaultSection() 2896 @*/ 2897 PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *section) 2898 { 2899 PetscErrorCode ierr; 2900 2901 PetscFunctionBegin; 2902 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2903 PetscValidPointer(section, 2); 2904 if (!dm->defaultGlobalSection) { 2905 PetscSection s; 2906 2907 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 2908 if (!s) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONGSTATE, "DM must have a default PetscSection in order to create a global PetscSection"); 2909 if (!dm->sf) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM must have a default PetscSF in order to create a global PetscSection"); 2910 ierr = PetscSectionCreateGlobalSection(s, dm->sf, PETSC_FALSE, &dm->defaultGlobalSection);CHKERRQ(ierr); 2911 ierr = PetscLayoutDestroy(&dm->map);CHKERRQ(ierr); 2912 ierr = PetscSectionGetValueLayout(PetscObjectComm((PetscObject)dm), dm->defaultGlobalSection, &dm->map);CHKERRQ(ierr); 2913 } 2914 *section = dm->defaultGlobalSection; 2915 PetscFunctionReturn(0); 2916 } 2917 2918 #undef __FUNCT__ 2919 #define __FUNCT__ "DMSetDefaultGlobalSection" 2920 /*@ 2921 DMSetDefaultGlobalSection - Set the PetscSection encoding the global data layout for the DM. 2922 2923 Input Parameters: 2924 + dm - The DM 2925 - section - The PetscSection, or NULL 2926 2927 Level: intermediate 2928 2929 Note: Any existing Section will be destroyed 2930 2931 .seealso: DMGetDefaultGlobalSection(), DMSetDefaultSection() 2932 @*/ 2933 PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection section) 2934 { 2935 PetscErrorCode ierr; 2936 2937 PetscFunctionBegin; 2938 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2939 if (section) PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2); 2940 ierr = PetscObjectReference((PetscObject)section);CHKERRQ(ierr); 2941 ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr); 2942 dm->defaultGlobalSection = section; 2943 PetscFunctionReturn(0); 2944 } 2945 2946 #undef __FUNCT__ 2947 #define __FUNCT__ "DMGetDefaultSF" 2948 /*@ 2949 DMGetDefaultSF - Get the PetscSF encoding the parallel dof overlap for the DM. If it has not been set, 2950 it is created from the default PetscSection layouts in the DM. 2951 2952 Input Parameter: 2953 . dm - The DM 2954 2955 Output Parameter: 2956 . sf - The PetscSF 2957 2958 Level: intermediate 2959 2960 Note: This gets a borrowed reference, so the user should not destroy this PetscSF. 2961 2962 .seealso: DMSetDefaultSF(), DMCreateDefaultSF() 2963 @*/ 2964 PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *sf) 2965 { 2966 PetscInt nroots; 2967 PetscErrorCode ierr; 2968 2969 PetscFunctionBegin; 2970 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2971 PetscValidPointer(sf, 2); 2972 ierr = PetscSFGetGraph(dm->defaultSF, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 2973 if (nroots < 0) { 2974 PetscSection section, gSection; 2975 2976 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 2977 if (section) { 2978 ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr); 2979 ierr = DMCreateDefaultSF(dm, section, gSection);CHKERRQ(ierr); 2980 } else { 2981 *sf = NULL; 2982 PetscFunctionReturn(0); 2983 } 2984 } 2985 *sf = dm->defaultSF; 2986 PetscFunctionReturn(0); 2987 } 2988 2989 #undef __FUNCT__ 2990 #define __FUNCT__ "DMSetDefaultSF" 2991 /*@ 2992 DMSetDefaultSF - Set the PetscSF encoding the parallel dof overlap for the DM 2993 2994 Input Parameters: 2995 + dm - The DM 2996 - sf - The PetscSF 2997 2998 Level: intermediate 2999 3000 Note: Any previous SF is destroyed 3001 3002 .seealso: DMGetDefaultSF(), DMCreateDefaultSF() 3003 @*/ 3004 PetscErrorCode DMSetDefaultSF(DM dm, PetscSF sf) 3005 { 3006 PetscErrorCode ierr; 3007 3008 PetscFunctionBegin; 3009 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3010 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 3011 ierr = PetscSFDestroy(&dm->defaultSF);CHKERRQ(ierr); 3012 dm->defaultSF = sf; 3013 PetscFunctionReturn(0); 3014 } 3015 3016 #undef __FUNCT__ 3017 #define __FUNCT__ "DMCreateDefaultSF" 3018 /*@C 3019 DMCreateDefaultSF - Create the PetscSF encoding the parallel dof overlap for the DM based upon the PetscSections 3020 describing the data layout. 3021 3022 Input Parameters: 3023 + dm - The DM 3024 . localSection - PetscSection describing the local data layout 3025 - globalSection - PetscSection describing the global data layout 3026 3027 Level: intermediate 3028 3029 .seealso: DMGetDefaultSF(), DMSetDefaultSF() 3030 @*/ 3031 PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection localSection, PetscSection globalSection) 3032 { 3033 MPI_Comm comm; 3034 PetscLayout layout; 3035 const PetscInt *ranges; 3036 PetscInt *local; 3037 PetscSFNode *remote; 3038 PetscInt pStart, pEnd, p, nroots, nleaves = 0, l; 3039 PetscMPIInt size, rank; 3040 PetscErrorCode ierr; 3041 3042 PetscFunctionBegin; 3043 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 3044 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3045 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 3046 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 3047 ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr); 3048 ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr); 3049 ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 3050 ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr); 3051 ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr); 3052 ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 3053 ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr); 3054 for (p = pStart; p < pEnd; ++p) { 3055 PetscInt gdof, gcdof; 3056 3057 ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr); 3058 ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr); 3059 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)); 3060 nleaves += gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof; 3061 } 3062 ierr = PetscMalloc1(nleaves, &local);CHKERRQ(ierr); 3063 ierr = PetscMalloc1(nleaves, &remote);CHKERRQ(ierr); 3064 for (p = pStart, l = 0; p < pEnd; ++p) { 3065 const PetscInt *cind; 3066 PetscInt dof, cdof, off, gdof, gcdof, goff, gsize, d, c; 3067 3068 ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr); 3069 ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr); 3070 ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr); 3071 ierr = PetscSectionGetConstraintIndices(localSection, p, &cind);CHKERRQ(ierr); 3072 ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr); 3073 ierr = PetscSectionGetConstraintDof(globalSection, p, &gcdof);CHKERRQ(ierr); 3074 ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr); 3075 if (!gdof) continue; /* Censored point */ 3076 gsize = gdof < 0 ? -(gdof+1)-gcdof : gdof-gcdof; 3077 if (gsize != dof-cdof) { 3078 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); 3079 cdof = 0; /* Ignore constraints */ 3080 } 3081 for (d = 0, c = 0; d < dof; ++d) { 3082 if ((c < cdof) && (cind[c] == d)) {++c; continue;} 3083 local[l+d-c] = off+d; 3084 } 3085 if (gdof < 0) { 3086 for (d = 0; d < gsize; ++d, ++l) { 3087 PetscInt offset = -(goff+1) + d, r; 3088 3089 ierr = PetscFindInt(offset,size+1,ranges,&r);CHKERRQ(ierr); 3090 if (r < 0) r = -(r+2); 3091 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); 3092 remote[l].rank = r; 3093 remote[l].index = offset - ranges[r]; 3094 } 3095 } else { 3096 for (d = 0; d < gsize; ++d, ++l) { 3097 remote[l].rank = rank; 3098 remote[l].index = goff+d - ranges[rank]; 3099 } 3100 } 3101 } 3102 if (l != nleaves) SETERRQ2(comm, PETSC_ERR_PLIB, "Iteration error, l %d != nleaves %d", l, nleaves); 3103 ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 3104 ierr = PetscSFSetGraph(dm->defaultSF, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr); 3105 PetscFunctionReturn(0); 3106 } 3107 3108 #undef __FUNCT__ 3109 #define __FUNCT__ "DMGetPointSF" 3110 /*@ 3111 DMGetPointSF - Get the PetscSF encoding the parallel section point overlap for the DM. 3112 3113 Input Parameter: 3114 . dm - The DM 3115 3116 Output Parameter: 3117 . sf - The PetscSF 3118 3119 Level: intermediate 3120 3121 Note: This gets a borrowed reference, so the user should not destroy this PetscSF. 3122 3123 .seealso: DMSetPointSF(), DMGetDefaultSF(), DMSetDefaultSF(), DMCreateDefaultSF() 3124 @*/ 3125 PetscErrorCode DMGetPointSF(DM dm, PetscSF *sf) 3126 { 3127 PetscFunctionBegin; 3128 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3129 PetscValidPointer(sf, 2); 3130 *sf = dm->sf; 3131 PetscFunctionReturn(0); 3132 } 3133 3134 #undef __FUNCT__ 3135 #define __FUNCT__ "DMSetPointSF" 3136 /*@ 3137 DMSetPointSF - Set the PetscSF encoding the parallel section point overlap for the DM. 3138 3139 Input Parameters: 3140 + dm - The DM 3141 - sf - The PetscSF 3142 3143 Level: intermediate 3144 3145 .seealso: DMGetPointSF(), DMGetDefaultSF(), DMSetDefaultSF(), DMCreateDefaultSF() 3146 @*/ 3147 PetscErrorCode DMSetPointSF(DM dm, PetscSF sf) 3148 { 3149 PetscErrorCode ierr; 3150 3151 PetscFunctionBegin; 3152 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3153 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 1); 3154 ierr = PetscSFDestroy(&dm->sf);CHKERRQ(ierr); 3155 ierr = PetscObjectReference((PetscObject) sf);CHKERRQ(ierr); 3156 dm->sf = sf; 3157 PetscFunctionReturn(0); 3158 } 3159 3160 #undef __FUNCT__ 3161 #define __FUNCT__ "DMGetDS" 3162 /*@ 3163 DMGetDS - Get the PetscDS 3164 3165 Input Parameter: 3166 . dm - The DM 3167 3168 Output Parameter: 3169 . prob - The PetscDS 3170 3171 Level: developer 3172 3173 .seealso: DMSetDS() 3174 @*/ 3175 PetscErrorCode DMGetDS(DM dm, PetscDS *prob) 3176 { 3177 PetscFunctionBegin; 3178 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3179 PetscValidPointer(prob, 2); 3180 *prob = dm->prob; 3181 PetscFunctionReturn(0); 3182 } 3183 3184 #undef __FUNCT__ 3185 #define __FUNCT__ "DMSetDS" 3186 /*@ 3187 DMSetDS - Set the PetscDS 3188 3189 Input Parameters: 3190 + dm - The DM 3191 - prob - The PetscDS 3192 3193 Level: developer 3194 3195 .seealso: DMGetDS() 3196 @*/ 3197 PetscErrorCode DMSetDS(DM dm, PetscDS prob) 3198 { 3199 PetscErrorCode ierr; 3200 3201 PetscFunctionBegin; 3202 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3203 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2); 3204 ierr = PetscDSDestroy(&dm->prob);CHKERRQ(ierr); 3205 dm->prob = prob; 3206 ierr = PetscObjectReference((PetscObject) dm->prob);CHKERRQ(ierr); 3207 PetscFunctionReturn(0); 3208 } 3209 3210 #undef __FUNCT__ 3211 #define __FUNCT__ "DMGetNumFields" 3212 PetscErrorCode DMGetNumFields(DM dm, PetscInt *numFields) 3213 { 3214 PetscErrorCode ierr; 3215 3216 PetscFunctionBegin; 3217 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3218 ierr = PetscDSGetNumFields(dm->prob, numFields);CHKERRQ(ierr); 3219 PetscFunctionReturn(0); 3220 } 3221 3222 #undef __FUNCT__ 3223 #define __FUNCT__ "DMSetNumFields" 3224 PetscErrorCode DMSetNumFields(DM dm, PetscInt numFields) 3225 { 3226 PetscInt Nf, f; 3227 PetscErrorCode ierr; 3228 3229 PetscFunctionBegin; 3230 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3231 ierr = PetscDSGetNumFields(dm->prob, &Nf);CHKERRQ(ierr); 3232 for (f = Nf; f < numFields; ++f) { 3233 PetscContainer obj; 3234 3235 ierr = PetscContainerCreate(PetscObjectComm((PetscObject) dm), &obj);CHKERRQ(ierr); 3236 ierr = PetscDSSetDiscretization(dm->prob, f, (PetscObject) obj);CHKERRQ(ierr); 3237 ierr = PetscContainerDestroy(&obj);CHKERRQ(ierr); 3238 } 3239 PetscFunctionReturn(0); 3240 } 3241 3242 #undef __FUNCT__ 3243 #define __FUNCT__ "DMGetField" 3244 /*@ 3245 DMGetField - Return the discretization object for a given DM field 3246 3247 Not collective 3248 3249 Input Parameters: 3250 + dm - The DM 3251 - f - The field number 3252 3253 Output Parameter: 3254 . field - The discretization object 3255 3256 Level: developer 3257 3258 .seealso: DMSetField() 3259 @*/ 3260 PetscErrorCode DMGetField(DM dm, PetscInt f, PetscObject *field) 3261 { 3262 PetscErrorCode ierr; 3263 3264 PetscFunctionBegin; 3265 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3266 ierr = PetscDSGetDiscretization(dm->prob, f, field);CHKERRQ(ierr); 3267 PetscFunctionReturn(0); 3268 } 3269 3270 #undef __FUNCT__ 3271 #define __FUNCT__ "DMSetField" 3272 /*@ 3273 DMSetField - Set the discretization object for a given DM field 3274 3275 Logically collective on DM 3276 3277 Input Parameters: 3278 + dm - The DM 3279 . f - The field number 3280 - field - The discretization object 3281 3282 Level: developer 3283 3284 .seealso: DMGetField() 3285 @*/ 3286 PetscErrorCode DMSetField(DM dm, PetscInt f, PetscObject field) 3287 { 3288 PetscErrorCode ierr; 3289 3290 PetscFunctionBegin; 3291 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3292 ierr = PetscDSSetDiscretization(dm->prob, f, field);CHKERRQ(ierr); 3293 PetscFunctionReturn(0); 3294 } 3295 3296 #undef __FUNCT__ 3297 #define __FUNCT__ "DMRestrictHook_Coordinates" 3298 PetscErrorCode DMRestrictHook_Coordinates(DM dm,DM dmc,void *ctx) 3299 { 3300 DM dm_coord,dmc_coord; 3301 PetscErrorCode ierr; 3302 Vec coords,ccoords; 3303 VecScatter scat; 3304 PetscFunctionBegin; 3305 ierr = DMGetCoordinateDM(dm,&dm_coord);CHKERRQ(ierr); 3306 ierr = DMGetCoordinateDM(dmc,&dmc_coord);CHKERRQ(ierr); 3307 ierr = DMGetCoordinates(dm,&coords);CHKERRQ(ierr); 3308 ierr = DMGetCoordinates(dmc,&ccoords);CHKERRQ(ierr); 3309 if (coords && !ccoords) { 3310 ierr = DMCreateGlobalVector(dmc_coord,&ccoords);CHKERRQ(ierr); 3311 ierr = DMCreateInjection(dmc_coord,dm_coord,&scat);CHKERRQ(ierr); 3312 ierr = VecScatterBegin(scat,coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3313 ierr = VecScatterEnd(scat,coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3314 ierr = DMSetCoordinates(dmc,ccoords);CHKERRQ(ierr); 3315 ierr = VecScatterDestroy(&scat);CHKERRQ(ierr); 3316 ierr = VecDestroy(&ccoords);CHKERRQ(ierr); 3317 } 3318 PetscFunctionReturn(0); 3319 } 3320 3321 #undef __FUNCT__ 3322 #define __FUNCT__ "DMSubDomainHook_Coordinates" 3323 static PetscErrorCode DMSubDomainHook_Coordinates(DM dm,DM subdm,void *ctx) 3324 { 3325 DM dm_coord,subdm_coord; 3326 PetscErrorCode ierr; 3327 Vec coords,ccoords,clcoords; 3328 VecScatter *scat_i,*scat_g; 3329 PetscFunctionBegin; 3330 ierr = DMGetCoordinateDM(dm,&dm_coord);CHKERRQ(ierr); 3331 ierr = DMGetCoordinateDM(subdm,&subdm_coord);CHKERRQ(ierr); 3332 ierr = DMGetCoordinates(dm,&coords);CHKERRQ(ierr); 3333 ierr = DMGetCoordinates(subdm,&ccoords);CHKERRQ(ierr); 3334 if (coords && !ccoords) { 3335 ierr = DMCreateGlobalVector(subdm_coord,&ccoords);CHKERRQ(ierr); 3336 ierr = DMCreateLocalVector(subdm_coord,&clcoords);CHKERRQ(ierr); 3337 ierr = DMCreateDomainDecompositionScatters(dm_coord,1,&subdm_coord,NULL,&scat_i,&scat_g);CHKERRQ(ierr); 3338 ierr = VecScatterBegin(scat_i[0],coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3339 ierr = VecScatterBegin(scat_g[0],coords,clcoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3340 ierr = VecScatterEnd(scat_i[0],coords,ccoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3341 ierr = VecScatterEnd(scat_g[0],coords,clcoords,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3342 ierr = DMSetCoordinates(subdm,ccoords);CHKERRQ(ierr); 3343 ierr = DMSetCoordinatesLocal(subdm,clcoords);CHKERRQ(ierr); 3344 ierr = VecScatterDestroy(&scat_i[0]);CHKERRQ(ierr); 3345 ierr = VecScatterDestroy(&scat_g[0]);CHKERRQ(ierr); 3346 ierr = VecDestroy(&ccoords);CHKERRQ(ierr); 3347 ierr = VecDestroy(&clcoords);CHKERRQ(ierr); 3348 ierr = PetscFree(scat_i);CHKERRQ(ierr); 3349 ierr = PetscFree(scat_g);CHKERRQ(ierr); 3350 } 3351 PetscFunctionReturn(0); 3352 } 3353 3354 #undef __FUNCT__ 3355 #define __FUNCT__ "DMSetCoordinates" 3356 /*@ 3357 DMSetCoordinates - Sets into the DM a global vector that holds the coordinates 3358 3359 Collective on DM 3360 3361 Input Parameters: 3362 + dm - the DM 3363 - c - coordinate vector 3364 3365 Note: 3366 The coordinates do include those for ghost points, which are in the local vector 3367 3368 Level: intermediate 3369 3370 .keywords: distributed array, get, corners, nodes, local indices, coordinates 3371 .seealso: DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLoca(), DMGetCoordinateDM() 3372 @*/ 3373 PetscErrorCode DMSetCoordinates(DM dm, Vec c) 3374 { 3375 PetscErrorCode ierr; 3376 3377 PetscFunctionBegin; 3378 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3379 PetscValidHeaderSpecific(c,VEC_CLASSID,2); 3380 ierr = PetscObjectReference((PetscObject) c);CHKERRQ(ierr); 3381 ierr = VecDestroy(&dm->coordinates);CHKERRQ(ierr); 3382 dm->coordinates = c; 3383 ierr = VecDestroy(&dm->coordinatesLocal);CHKERRQ(ierr); 3384 ierr = DMCoarsenHookAdd(dm,DMRestrictHook_Coordinates,NULL,NULL);CHKERRQ(ierr); 3385 ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_Coordinates,NULL,NULL);CHKERRQ(ierr); 3386 PetscFunctionReturn(0); 3387 } 3388 3389 #undef __FUNCT__ 3390 #define __FUNCT__ "DMSetCoordinatesLocal" 3391 /*@ 3392 DMSetCoordinatesLocal - Sets into the DM a local vector that holds the coordinates 3393 3394 Collective on DM 3395 3396 Input Parameters: 3397 + dm - the DM 3398 - c - coordinate vector 3399 3400 Note: 3401 The coordinates of ghost points can be set using DMSetCoordinates() 3402 followed by DMGetCoordinatesLocal(). This is intended to enable the 3403 setting of ghost coordinates outside of the domain. 3404 3405 Level: intermediate 3406 3407 .keywords: distributed array, get, corners, nodes, local indices, coordinates 3408 .seealso: DMGetCoordinatesLocal(), DMSetCoordinates(), DMGetCoordinates(), DMGetCoordinateDM() 3409 @*/ 3410 PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c) 3411 { 3412 PetscErrorCode ierr; 3413 3414 PetscFunctionBegin; 3415 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3416 PetscValidHeaderSpecific(c,VEC_CLASSID,2); 3417 ierr = PetscObjectReference((PetscObject) c);CHKERRQ(ierr); 3418 ierr = VecDestroy(&dm->coordinatesLocal);CHKERRQ(ierr); 3419 3420 dm->coordinatesLocal = c; 3421 3422 ierr = VecDestroy(&dm->coordinates);CHKERRQ(ierr); 3423 PetscFunctionReturn(0); 3424 } 3425 3426 #undef __FUNCT__ 3427 #define __FUNCT__ "DMGetCoordinates" 3428 /*@ 3429 DMGetCoordinates - Gets a global vector with the coordinates associated with the DM. 3430 3431 Not Collective 3432 3433 Input Parameter: 3434 . dm - the DM 3435 3436 Output Parameter: 3437 . c - global coordinate vector 3438 3439 Note: 3440 This is a borrowed reference, so the user should NOT destroy this vector 3441 3442 Each process has only the local coordinates (does NOT have the ghost coordinates). 3443 3444 For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 3445 and (x_0,y_0,z_0,x_1,y_1,z_1...) 3446 3447 Level: intermediate 3448 3449 .keywords: distributed array, get, corners, nodes, local indices, coordinates 3450 .seealso: DMSetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM() 3451 @*/ 3452 PetscErrorCode DMGetCoordinates(DM dm, Vec *c) 3453 { 3454 PetscErrorCode ierr; 3455 3456 PetscFunctionBegin; 3457 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3458 PetscValidPointer(c,2); 3459 if (!dm->coordinates && dm->coordinatesLocal) { 3460 DM cdm = NULL; 3461 3462 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 3463 ierr = DMCreateGlobalVector(cdm, &dm->coordinates);CHKERRQ(ierr); 3464 ierr = PetscObjectSetName((PetscObject) dm->coordinates, "coordinates");CHKERRQ(ierr); 3465 ierr = DMLocalToGlobalBegin(cdm, dm->coordinatesLocal, INSERT_VALUES, dm->coordinates);CHKERRQ(ierr); 3466 ierr = DMLocalToGlobalEnd(cdm, dm->coordinatesLocal, INSERT_VALUES, dm->coordinates);CHKERRQ(ierr); 3467 } 3468 *c = dm->coordinates; 3469 PetscFunctionReturn(0); 3470 } 3471 3472 #undef __FUNCT__ 3473 #define __FUNCT__ "DMGetCoordinatesLocal" 3474 /*@ 3475 DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the DM. 3476 3477 Collective on DM 3478 3479 Input Parameter: 3480 . dm - the DM 3481 3482 Output Parameter: 3483 . c - coordinate vector 3484 3485 Note: 3486 This is a borrowed reference, so the user should NOT destroy this vector 3487 3488 Each process has the local and ghost coordinates 3489 3490 For DMDA, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...) 3491 and (x_0,y_0,z_0,x_1,y_1,z_1...) 3492 3493 Level: intermediate 3494 3495 .keywords: distributed array, get, corners, nodes, local indices, coordinates 3496 .seealso: DMSetCoordinatesLocal(), DMGetCoordinates(), DMSetCoordinates(), DMGetCoordinateDM() 3497 @*/ 3498 PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c) 3499 { 3500 PetscErrorCode ierr; 3501 3502 PetscFunctionBegin; 3503 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3504 PetscValidPointer(c,2); 3505 if (!dm->coordinatesLocal && dm->coordinates) { 3506 DM cdm = NULL; 3507 3508 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 3509 ierr = DMCreateLocalVector(cdm, &dm->coordinatesLocal);CHKERRQ(ierr); 3510 ierr = PetscObjectSetName((PetscObject) dm->coordinatesLocal, "coordinates");CHKERRQ(ierr); 3511 ierr = DMGlobalToLocalBegin(cdm, dm->coordinates, INSERT_VALUES, dm->coordinatesLocal);CHKERRQ(ierr); 3512 ierr = DMGlobalToLocalEnd(cdm, dm->coordinates, INSERT_VALUES, dm->coordinatesLocal);CHKERRQ(ierr); 3513 } 3514 *c = dm->coordinatesLocal; 3515 PetscFunctionReturn(0); 3516 } 3517 3518 #undef __FUNCT__ 3519 #define __FUNCT__ "DMGetCoordinateDM" 3520 /*@ 3521 DMGetCoordinateDM - Gets the DM that prescribes coordinate layout and scatters between global and local coordinates 3522 3523 Collective on DM 3524 3525 Input Parameter: 3526 . dm - the DM 3527 3528 Output Parameter: 3529 . cdm - coordinate DM 3530 3531 Level: intermediate 3532 3533 .keywords: distributed array, get, corners, nodes, local indices, coordinates 3534 .seealso: DMSetCoordinateDM(), DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal() 3535 @*/ 3536 PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm) 3537 { 3538 PetscErrorCode ierr; 3539 3540 PetscFunctionBegin; 3541 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3542 PetscValidPointer(cdm,2); 3543 if (!dm->coordinateDM) { 3544 if (!dm->ops->createcoordinatedm) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unable to create coordinates for this DM"); 3545 ierr = (*dm->ops->createcoordinatedm)(dm, &dm->coordinateDM);CHKERRQ(ierr); 3546 } 3547 *cdm = dm->coordinateDM; 3548 PetscFunctionReturn(0); 3549 } 3550 3551 #undef __FUNCT__ 3552 #define __FUNCT__ "DMSetCoordinateDM" 3553 /*@ 3554 DMSetCoordinateDM - Sets the DM that prescribes coordinate layout and scatters between global and local coordinates 3555 3556 Logically Collective on DM 3557 3558 Input Parameters: 3559 + dm - the DM 3560 - cdm - coordinate DM 3561 3562 Level: intermediate 3563 3564 .keywords: distributed array, get, corners, nodes, local indices, coordinates 3565 .seealso: DMGetCoordinateDM(), DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal() 3566 @*/ 3567 PetscErrorCode DMSetCoordinateDM(DM dm, DM cdm) 3568 { 3569 PetscErrorCode ierr; 3570 3571 PetscFunctionBegin; 3572 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3573 PetscValidHeaderSpecific(cdm,DM_CLASSID,2); 3574 ierr = DMDestroy(&dm->coordinateDM);CHKERRQ(ierr); 3575 dm->coordinateDM = cdm; 3576 ierr = PetscObjectReference((PetscObject) dm->coordinateDM);CHKERRQ(ierr); 3577 PetscFunctionReturn(0); 3578 } 3579 3580 #undef __FUNCT__ 3581 #define __FUNCT__ "DMGetCoordinateSection" 3582 /*@ 3583 DMGetCoordinateSection - Retrieve the layout of coordinate values over the mesh. 3584 3585 Not Collective 3586 3587 Input Parameter: 3588 . dm - The DM object 3589 3590 Output Parameter: 3591 . section - The PetscSection object 3592 3593 Level: intermediate 3594 3595 .keywords: mesh, coordinates 3596 .seealso: DMGetCoordinateDM(), DMGetDefaultSection(), DMSetDefaultSection() 3597 @*/ 3598 PetscErrorCode DMGetCoordinateSection(DM dm, PetscSection *section) 3599 { 3600 DM cdm; 3601 PetscErrorCode ierr; 3602 3603 PetscFunctionBegin; 3604 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3605 PetscValidPointer(section, 2); 3606 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 3607 ierr = DMGetDefaultSection(cdm, section);CHKERRQ(ierr); 3608 PetscFunctionReturn(0); 3609 } 3610 3611 #undef __FUNCT__ 3612 #define __FUNCT__ "DMSetCoordinateSection" 3613 /*@ 3614 DMSetCoordinateSection - Set the layout of coordinate values over the mesh. 3615 3616 Not Collective 3617 3618 Input Parameters: 3619 + dm - The DM object 3620 - section - The PetscSection object 3621 3622 Level: intermediate 3623 3624 .keywords: mesh, coordinates 3625 .seealso: DMGetCoordinateSection(), DMGetDefaultSection(), DMSetDefaultSection() 3626 @*/ 3627 PetscErrorCode DMSetCoordinateSection(DM dm, PetscSection section) 3628 { 3629 DM cdm; 3630 PetscErrorCode ierr; 3631 3632 PetscFunctionBegin; 3633 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3634 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,2); 3635 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 3636 ierr = DMSetDefaultSection(cdm, section);CHKERRQ(ierr); 3637 PetscFunctionReturn(0); 3638 } 3639 3640 #undef __FUNCT__ 3641 #define __FUNCT__ "DMGetPeriodicity" 3642 PetscErrorCode DMGetPeriodicity(DM dm, const PetscReal **maxCell, const PetscReal **L) 3643 { 3644 PetscFunctionBegin; 3645 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3646 if (L) *L = dm->L; 3647 if (maxCell) *maxCell = dm->maxCell; 3648 PetscFunctionReturn(0); 3649 } 3650 3651 #undef __FUNCT__ 3652 #define __FUNCT__ "DMSetPeriodicity" 3653 PetscErrorCode DMSetPeriodicity(DM dm, const PetscReal maxCell[], const PetscReal L[]) 3654 { 3655 Vec coordinates; 3656 PetscInt dim, d; 3657 PetscErrorCode ierr; 3658 3659 PetscFunctionBegin; 3660 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3661 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 3662 ierr = VecGetBlockSize(coordinates, &dim);CHKERRQ(ierr); 3663 ierr = PetscMalloc1(dim,&dm->L);CHKERRQ(ierr); 3664 ierr = PetscMalloc1(dim,&dm->maxCell);CHKERRQ(ierr); 3665 for (d = 0; d < dim; ++d) {dm->L[d] = L[d]; dm->maxCell[d] = maxCell[d];} 3666 PetscFunctionReturn(0); 3667 } 3668 3669 #undef __FUNCT__ 3670 #define __FUNCT__ "DMLocatePoints" 3671 /*@ 3672 DMLocatePoints - Locate the points in v in the mesh and return an IS of the containing cells 3673 3674 Not collective 3675 3676 Input Parameters: 3677 + dm - The DM 3678 - v - The Vec of points 3679 3680 Output Parameter: 3681 . cells - The local cell numbers for cells which contain the points 3682 3683 Level: developer 3684 3685 .keywords: point location, mesh 3686 .seealso: DMSetCoordinates(), DMSetCoordinatesLocal(), DMGetCoordinates(), DMGetCoordinatesLocal() 3687 @*/ 3688 PetscErrorCode DMLocatePoints(DM dm, Vec v, IS *cells) 3689 { 3690 PetscErrorCode ierr; 3691 3692 PetscFunctionBegin; 3693 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3694 PetscValidHeaderSpecific(v,VEC_CLASSID,2); 3695 PetscValidPointer(cells,3); 3696 if (dm->ops->locatepoints) { 3697 ierr = (*dm->ops->locatepoints)(dm,v,cells);CHKERRQ(ierr); 3698 } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Point location not available for this DM"); 3699 PetscFunctionReturn(0); 3700 } 3701 3702 #undef __FUNCT__ 3703 #define __FUNCT__ "DMGetOutputDM" 3704 /*@ 3705 DMGetOutputDM - Retrieve the DM associated with the layout for output 3706 3707 Input Parameter: 3708 . dm - The original DM 3709 3710 Output Parameter: 3711 . odm - The DM which provides the layout for output 3712 3713 Level: intermediate 3714 3715 .seealso: VecView(), DMGetDefaultGlobalSection() 3716 @*/ 3717 PetscErrorCode DMGetOutputDM(DM dm, DM *odm) 3718 { 3719 PetscSection section; 3720 PetscBool hasConstraints; 3721 PetscErrorCode ierr; 3722 3723 PetscFunctionBegin; 3724 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3725 PetscValidPointer(odm,2); 3726 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 3727 ierr = PetscSectionHasConstraints(section, &hasConstraints);CHKERRQ(ierr); 3728 if (!hasConstraints) { 3729 *odm = dm; 3730 PetscFunctionReturn(0); 3731 } 3732 if (!dm->dmBC) { 3733 PetscSection newSection, gsection; 3734 PetscSF sf; 3735 3736 ierr = DMClone(dm, &dm->dmBC);CHKERRQ(ierr); 3737 ierr = PetscSectionClone(section, &newSection);CHKERRQ(ierr); 3738 ierr = DMSetDefaultSection(dm->dmBC, newSection);CHKERRQ(ierr); 3739 ierr = PetscSectionDestroy(&newSection);CHKERRQ(ierr); 3740 ierr = DMGetPointSF(dm->dmBC, &sf);CHKERRQ(ierr); 3741 ierr = PetscSectionCreateGlobalSection(section, sf, PETSC_TRUE, &gsection);CHKERRQ(ierr); 3742 ierr = DMSetDefaultGlobalSection(dm->dmBC, gsection);CHKERRQ(ierr); 3743 ierr = PetscSectionDestroy(&gsection);CHKERRQ(ierr); 3744 } 3745 *odm = dm->dmBC; 3746 PetscFunctionReturn(0); 3747 } 3748 3749 #undef __FUNCT__ 3750 #define __FUNCT__ "DMGetOutputSequenceNumber" 3751 /*@ 3752 DMGetOutputSequenceNumber - Retrieve the sequence number/value for output 3753 3754 Input Parameter: 3755 . dm - The original DM 3756 3757 Output Parameters: 3758 + num - The output sequence number 3759 - val - The output sequence value 3760 3761 Level: intermediate 3762 3763 Note: This is intended for output that should appear in sequence, for instance 3764 a set of timesteps in an HDF5 file, or a set of realizations of a stochastic system. 3765 3766 .seealso: VecView() 3767 @*/ 3768 PetscErrorCode DMGetOutputSequenceNumber(DM dm, PetscInt *num, PetscReal *val) 3769 { 3770 PetscFunctionBegin; 3771 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3772 if (num) {PetscValidPointer(num,2); *num = dm->outputSequenceNum;} 3773 if (val) {PetscValidPointer(val,3);*val = dm->outputSequenceVal;} 3774 PetscFunctionReturn(0); 3775 } 3776 3777 #undef __FUNCT__ 3778 #define __FUNCT__ "DMSetOutputSequenceNumber" 3779 /*@ 3780 DMSetOutputSequenceNumber - Set the sequence number/value for output 3781 3782 Input Parameters: 3783 + dm - The original DM 3784 . num - The output sequence number 3785 - val - The output sequence value 3786 3787 Level: intermediate 3788 3789 Note: This is intended for output that should appear in sequence, for instance 3790 a set of timesteps in an HDF5 file, or a set of realizations of a stochastic system. 3791 3792 .seealso: VecView() 3793 @*/ 3794 PetscErrorCode DMSetOutputSequenceNumber(DM dm, PetscInt num, PetscReal val) 3795 { 3796 PetscFunctionBegin; 3797 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3798 dm->outputSequenceNum = num; 3799 dm->outputSequenceVal = val; 3800 PetscFunctionReturn(0); 3801 } 3802 3803 #undef __FUNCT__ 3804 #define __FUNCT__ "DMOutputSequenceLoad" 3805 /*@C 3806 DMOutputSequenceLoad - Retrieve the sequence value from a Viewer 3807 3808 Input Parameters: 3809 + dm - The original DM 3810 . name - The sequence name 3811 - num - The output sequence number 3812 3813 Output Parameter: 3814 . val - The output sequence value 3815 3816 Level: intermediate 3817 3818 Note: This is intended for output that should appear in sequence, for instance 3819 a set of timesteps in an HDF5 file, or a set of realizations of a stochastic system. 3820 3821 .seealso: DMGetOutputSequenceNumber(), DMSetOutputSequenceNumber(), VecView() 3822 @*/ 3823 PetscErrorCode DMOutputSequenceLoad(DM dm, PetscViewer viewer, const char *name, PetscInt num, PetscReal *val) 3824 { 3825 PetscBool ishdf5; 3826 PetscErrorCode ierr; 3827 3828 PetscFunctionBegin; 3829 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3830 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 3831 PetscValidPointer(val,4); 3832 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5);CHKERRQ(ierr); 3833 if (ishdf5) { 3834 #if defined(PETSC_HAVE_HDF5) 3835 PetscScalar value; 3836 3837 ierr = DMSequenceLoad_HDF5(dm, name, num, &value, viewer);CHKERRQ(ierr); 3838 *val = PetscRealPart(value); 3839 #endif 3840 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerHDF5Open()"); 3841 PetscFunctionReturn(0); 3842 } 3843