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