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