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