1 2 #include <private/dmimpl.h> /*I "petscdm.h" I*/ 3 4 PetscClassId DM_CLASSID; 5 PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal; 6 7 #undef __FUNCT__ 8 #define __FUNCT__ "DMCreate" 9 /*@ 10 DMCreate - Creates an empty vector object. The type can then be set with DMetType(). 11 12 If you never call DMSetType() it will generate an 13 error when you try to use the vector. 14 15 Collective on MPI_Comm 16 17 Input Parameter: 18 . comm - The communicator for the DM object 19 20 Output Parameter: 21 . dm - The DM object 22 23 Level: beginner 24 25 .seealso: DMSetType(), DMDA, DMSLICED, DMCOMPOSITE 26 @*/ 27 PetscErrorCode DMCreate(MPI_Comm comm,DM *dm) 28 { 29 DM v; 30 PetscErrorCode ierr; 31 32 PetscFunctionBegin; 33 PetscValidPointer(dm,2); 34 *dm = PETSC_NULL; 35 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 36 ierr = DMInitializePackage(PETSC_NULL);CHKERRQ(ierr); 37 #endif 38 39 ierr = PetscHeaderCreate(v, _p_DM, struct _DMOps, DM_CLASSID, -1, "DM", comm, DMDestroy, DMView);CHKERRQ(ierr); 40 ierr = PetscMemzero(v->ops, sizeof(struct _DMOps));CHKERRQ(ierr); 41 42 v->ltogmap = PETSC_NULL; 43 v->ltogmapb = PETSC_NULL; 44 v->bs = 1; 45 46 *dm = v; 47 PetscFunctionReturn(0); 48 } 49 50 51 #undef __FUNCT__ 52 #define __FUNCT__ "DMSetVecType" 53 /*@C 54 DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() 55 56 Logically Collective on DMDA 57 58 Input Parameter: 59 + da - initial distributed array 60 . ctype - the vector type, currently either VECSTANDARD or VECCUSP 61 62 Options Database: 63 . -da_vec_type ctype 64 65 Level: intermediate 66 67 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType 68 @*/ 69 PetscErrorCode DMSetVecType(DM da,const VecType ctype) 70 { 71 PetscErrorCode ierr; 72 73 PetscFunctionBegin; 74 PetscValidHeaderSpecific(da,DM_CLASSID,1); 75 ierr = PetscFree(da->vectype);CHKERRQ(ierr); 76 ierr = PetscStrallocpy(ctype,&da->vectype);CHKERRQ(ierr); 77 PetscFunctionReturn(0); 78 } 79 80 #undef __FUNCT__ 81 #define __FUNCT__ "DMSetOptionsPrefix" 82 /*@C 83 DMSetOptionsPrefix - Sets the prefix used for searching for all 84 DMDA options in the database. 85 86 Logically Collective on DMDA 87 88 Input Parameter: 89 + da - the DMDA context 90 - prefix - the prefix to prepend to all option names 91 92 Notes: 93 A hyphen (-) must NOT be given at the beginning of the prefix name. 94 The first character of all runtime options is AUTOMATICALLY the hyphen. 95 96 Level: advanced 97 98 .keywords: DMDA, set, options, prefix, database 99 100 .seealso: DMSetFromOptions() 101 @*/ 102 PetscErrorCode DMSetOptionsPrefix(DM dm,const char prefix[]) 103 { 104 PetscErrorCode ierr; 105 106 PetscFunctionBegin; 107 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 108 ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr); 109 PetscFunctionReturn(0); 110 } 111 112 #undef __FUNCT__ 113 #define __FUNCT__ "DMDestroy" 114 /*@ 115 DMDestroy - Destroys a vector packer or DMDA. 116 117 Collective on DM 118 119 Input Parameter: 120 . dm - the DM object to destroy 121 122 Level: developer 123 124 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix() 125 126 @*/ 127 PetscErrorCode DMDestroy(DM *dm) 128 { 129 PetscInt i, cnt = 0; 130 PetscErrorCode ierr; 131 132 PetscFunctionBegin; 133 if (!*dm) PetscFunctionReturn(0); 134 PetscValidHeaderSpecific((*dm),DM_CLASSID,1); 135 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 136 if ((*dm)->localin[i]) {cnt++;} 137 if ((*dm)->globalin[i]) {cnt++;} 138 } 139 140 if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);} 141 /* 142 Need this test because the dm references the vectors that 143 reference the dm, so destroying the dm calls destroy on the 144 vectors that cause another destroy on the dm 145 */ 146 if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0); 147 ((PetscObject) (*dm))->refct = 0; 148 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 149 if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()"); 150 ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr); 151 } 152 ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr); 153 ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr); 154 ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr); 155 ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr); 156 157 /* if memory was published with AMS then destroy it */ 158 ierr = PetscObjectDepublish(*dm);CHKERRQ(ierr); 159 160 ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr); 161 ierr = PetscFree((*dm)->data);CHKERRQ(ierr); 162 ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr); 163 PetscFunctionReturn(0); 164 } 165 166 #undef __FUNCT__ 167 #define __FUNCT__ "DMSetUp" 168 /*@ 169 DMSetUp - sets up the data structures inside a DM object 170 171 Collective on DM 172 173 Input Parameter: 174 . dm - the DM object to setup 175 176 Level: developer 177 178 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix() 179 180 @*/ 181 PetscErrorCode DMSetUp(DM dm) 182 { 183 PetscErrorCode ierr; 184 185 PetscFunctionBegin; 186 if (dm->ops->setup) { 187 ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr); 188 } 189 PetscFunctionReturn(0); 190 } 191 192 #undef __FUNCT__ 193 #define __FUNCT__ "DMSetFromOptions" 194 /*@ 195 DMSetFromOptions - sets parameters in a DM from the options database 196 197 Collective on DM 198 199 Input Parameter: 200 . dm - the DM object to set options for 201 202 Options Database: 203 . -dm_preallocate_only: Only preallocate the matrix for DMGetMatrix(), but do not fill it with zeros 204 205 Level: developer 206 207 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix() 208 209 @*/ 210 PetscErrorCode DMSetFromOptions(DM dm) 211 { 212 PetscErrorCode ierr; 213 214 PetscFunctionBegin; 215 ierr = PetscOptionsGetBool(((PetscObject) dm)->prefix, "-dm_preallocate_only", &dm->prealloc_only, PETSC_NULL);CHKERRQ(ierr); 216 if (dm->ops->setfromoptions) { 217 ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr); 218 } 219 PetscFunctionReturn(0); 220 } 221 222 #undef __FUNCT__ 223 #define __FUNCT__ "DMView" 224 /*@C 225 DMView - Views a vector packer or DMDA. 226 227 Collective on DM 228 229 Input Parameter: 230 + dm - the DM object to view 231 - v - the viewer 232 233 Level: developer 234 235 .seealso DMDestroy(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix() 236 237 @*/ 238 PetscErrorCode DMView(DM dm,PetscViewer v) 239 { 240 PetscErrorCode ierr; 241 242 PetscFunctionBegin; 243 if (!v) { 244 ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr); 245 } 246 if (dm->ops->view) { 247 ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr); 248 } 249 PetscFunctionReturn(0); 250 } 251 252 #undef __FUNCT__ 253 #define __FUNCT__ "DMCreateGlobalVector" 254 /*@ 255 DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object 256 257 Collective on DM 258 259 Input Parameter: 260 . dm - the DM object 261 262 Output Parameter: 263 . vec - the global vector 264 265 Level: developer 266 267 .seealso DMDestroy(), DMView(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix() 268 269 @*/ 270 PetscErrorCode DMCreateGlobalVector(DM dm,Vec *vec) 271 { 272 PetscErrorCode ierr; 273 274 PetscFunctionBegin; 275 ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr); 276 PetscFunctionReturn(0); 277 } 278 279 #undef __FUNCT__ 280 #define __FUNCT__ "DMCreateLocalVector" 281 /*@ 282 DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object 283 284 Not Collective 285 286 Input Parameter: 287 . dm - the DM object 288 289 Output Parameter: 290 . vec - the local vector 291 292 Level: developer 293 294 .seealso DMDestroy(), DMView(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix() 295 296 @*/ 297 PetscErrorCode DMCreateLocalVector(DM dm,Vec *vec) 298 { 299 PetscErrorCode ierr; 300 301 PetscFunctionBegin; 302 ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr); 303 PetscFunctionReturn(0); 304 } 305 306 #undef __FUNCT__ 307 #define __FUNCT__ "DMGetLocalToGlobalMapping" 308 /*@ 309 DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM. 310 311 Collective on DM 312 313 Input Parameter: 314 . dm - the DM that provides the mapping 315 316 Output Parameter: 317 . ltog - the mapping 318 319 Level: intermediate 320 321 Notes: 322 This mapping can then be used by VecSetLocalToGlobalMapping() or 323 MatSetLocalToGlobalMapping(). 324 325 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock() 326 @*/ 327 PetscErrorCode DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog) 328 { 329 PetscErrorCode ierr; 330 331 PetscFunctionBegin; 332 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 333 PetscValidPointer(ltog,2); 334 if (!dm->ltogmap) { 335 if (!dm->ops->createlocaltoglobalmapping) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping"); 336 ierr = (*dm->ops->createlocaltoglobalmapping)(dm);CHKERRQ(ierr); 337 } 338 *ltog = dm->ltogmap; 339 PetscFunctionReturn(0); 340 } 341 342 #undef __FUNCT__ 343 #define __FUNCT__ "DMGetLocalToGlobalMappingBlock" 344 /*@ 345 DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM. 346 347 Collective on DM 348 349 Input Parameter: 350 . da - the distributed array that provides the mapping 351 352 Output Parameter: 353 . ltog - the block mapping 354 355 Level: intermediate 356 357 Notes: 358 This mapping can then be used by VecSetLocalToGlobalMappingBlock() or 359 MatSetLocalToGlobalMappingBlock(). 360 361 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize() 362 @*/ 363 PetscErrorCode DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog) 364 { 365 PetscErrorCode ierr; 366 367 PetscFunctionBegin; 368 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 369 PetscValidPointer(ltog,2); 370 if (!dm->ltogmapb) { 371 PetscInt bs; 372 ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr); 373 if (bs > 1) { 374 if (!dm->ops->createlocaltoglobalmappingblock) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock"); 375 ierr = (*dm->ops->createlocaltoglobalmappingblock)(dm);CHKERRQ(ierr); 376 } else { 377 ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr); 378 ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr); 379 } 380 } 381 *ltog = dm->ltogmapb; 382 PetscFunctionReturn(0); 383 } 384 385 #undef __FUNCT__ 386 #define __FUNCT__ "DMGetBlockSize" 387 /*@ 388 DMGetBlockSize - Gets the inherent block size associated with a DM 389 390 Not Collective 391 392 Input Parameter: 393 . dm - the DM with block structure 394 395 Output Parameter: 396 . bs - the block size, 1 implies no exploitable block structure 397 398 Level: intermediate 399 400 .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock() 401 @*/ 402 PetscErrorCode DMGetBlockSize(DM dm,PetscInt *bs) 403 { 404 PetscFunctionBegin; 405 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 406 PetscValidPointer(bs,2); 407 if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet"); 408 *bs = dm->bs; 409 PetscFunctionReturn(0); 410 } 411 412 #undef __FUNCT__ 413 #define __FUNCT__ "DMGetInterpolation" 414 /*@ 415 DMGetInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects 416 417 Collective on DM 418 419 Input Parameter: 420 + dm1 - the DM object 421 - dm2 - the second, finer DM object 422 423 Output Parameter: 424 + mat - the interpolation 425 - vec - the scaling (optional) 426 427 Level: developer 428 429 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix() 430 431 @*/ 432 PetscErrorCode DMGetInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec) 433 { 434 PetscErrorCode ierr; 435 436 PetscFunctionBegin; 437 ierr = (*dm1->ops->getinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr); 438 PetscFunctionReturn(0); 439 } 440 441 #undef __FUNCT__ 442 #define __FUNCT__ "DMGetInjection" 443 /*@ 444 DMGetInjection - Gets injection matrix between two DMDA or DMComposite objects 445 446 Collective on DM 447 448 Input Parameter: 449 + dm1 - the DM object 450 - dm2 - the second, finer DM object 451 452 Output Parameter: 453 . ctx - the injection 454 455 Level: developer 456 457 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix(), DMGetInterpolation() 458 459 @*/ 460 PetscErrorCode DMGetInjection(DM dm1,DM dm2,VecScatter *ctx) 461 { 462 PetscErrorCode ierr; 463 464 PetscFunctionBegin; 465 ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr); 466 PetscFunctionReturn(0); 467 } 468 469 #undef __FUNCT__ 470 #define __FUNCT__ "DMGetColoring" 471 /*@C 472 DMGetColoring - Gets coloring for a DMDA or DMComposite 473 474 Collective on DM 475 476 Input Parameter: 477 + dm - the DM object 478 . ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL 479 - matype - either MATAIJ or MATBAIJ 480 481 Output Parameter: 482 . coloring - the coloring 483 484 Level: developer 485 486 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetMatrix() 487 488 @*/ 489 PetscErrorCode DMGetColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring) 490 { 491 PetscErrorCode ierr; 492 493 PetscFunctionBegin; 494 if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet"); 495 ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr); 496 PetscFunctionReturn(0); 497 } 498 499 #undef __FUNCT__ 500 #define __FUNCT__ "DMGetMatrix" 501 /*@C 502 DMGetMatrix - Gets empty Jacobian for a DMDA or DMComposite 503 504 Collective on DM 505 506 Input Parameter: 507 + dm - the DM object 508 - mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or 509 any type which inherits from one of these (such as MATAIJ) 510 511 Output Parameter: 512 . mat - the empty Jacobian 513 514 Level: developer 515 516 Notes: This properly preallocates the number of nonzeros in the sparse matrix so you 517 do not need to do it yourself. 518 519 By default it also sets the nonzero structure and puts in the zero entries. To prevent setting 520 the nonzero pattern call DMDASetMatPreallocateOnly() 521 522 For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used 523 internally by PETSc. 524 525 For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires 526 the indices for the global numbering for DMDAs which is complicated. 527 528 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 529 530 @*/ 531 PetscErrorCode DMGetMatrix(DM dm,const MatType mtype,Mat *mat) 532 { 533 PetscErrorCode ierr; 534 char ttype[256]; 535 PetscBool flg; 536 537 PetscFunctionBegin; 538 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 539 PetscValidPointer(mat,3); 540 ierr = PetscStrncpy(ttype,mtype,sizeof ttype);CHKERRQ(ierr); 541 ttype[sizeof ttype-1] = 0; 542 ierr = PetscOptionsBegin(((PetscObject)dm)->comm,((PetscObject)dm)->prefix,"DM options","Mat");CHKERRQ(ierr); 543 ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,ttype,ttype,sizeof ttype,&flg);CHKERRQ(ierr); 544 ierr = PetscOptionsEnd(); 545 if (flg || mtype) { 546 ierr = (*dm->ops->getmatrix)(dm,ttype,mat);CHKERRQ(ierr); 547 } else { /* Let the implementation decide */ 548 ierr = (*dm->ops->getmatrix)(dm,PETSC_NULL,mat);CHKERRQ(ierr); 549 } 550 PetscFunctionReturn(0); 551 } 552 553 #undef __FUNCT__ 554 #define __FUNCT__ "DMSetMatrixPreallocateOnly" 555 /*@ 556 DMSetMatrixPreallocateOnly - When DMGetMatrix() is called the matrix will be properly 557 preallocated but the nonzero structure and zero values will not be set. 558 559 Logically Collective on DMDA 560 561 Input Parameter: 562 + dm - the DM 563 - only - PETSC_TRUE if only want preallocation 564 565 Level: developer 566 .seealso DMGetMatrix() 567 @*/ 568 PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only) 569 { 570 PetscFunctionBegin; 571 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 572 dm->prealloc_only = only; 573 PetscFunctionReturn(0); 574 } 575 576 #undef __FUNCT__ 577 #define __FUNCT__ "DMRefine" 578 /*@ 579 DMRefine - Refines a DM object 580 581 Collective on DM 582 583 Input Parameter: 584 + dm - the DM object 585 - comm - the communicator to contain the new DM object (or PETSC_NULL) 586 587 Output Parameter: 588 . dmf - the refined DM 589 590 Level: developer 591 592 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 593 594 @*/ 595 PetscErrorCode DMRefine(DM dm,MPI_Comm comm,DM *dmf) 596 { 597 PetscErrorCode ierr; 598 599 PetscFunctionBegin; 600 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 601 ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr); 602 (*dmf)->levelup = dm->levelup + 1; 603 PetscFunctionReturn(0); 604 } 605 606 #undef __FUNCT__ 607 #define __FUNCT__ "DMGetRefineLevel" 608 /*@ 609 DMGetRefineLevel - Get's the number of refinements that have generated this DM. 610 611 Not Collective 612 613 Input Parameter: 614 . dm - the DM object 615 616 Output Parameter: 617 . level - number of refinements 618 619 Level: developer 620 621 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 622 623 @*/ 624 PetscErrorCode DMGetRefineLevel(DM dm,PetscInt *level) 625 { 626 PetscFunctionBegin; 627 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 628 *level = dm->levelup; 629 PetscFunctionReturn(0); 630 } 631 632 #undef __FUNCT__ 633 #define __FUNCT__ "DMGlobalToLocalBegin" 634 /*@ 635 DMGlobalToLocalBegin - Begins updating local vectors from global vector 636 637 Neighbor-wise Collective on DM 638 639 Input Parameters: 640 + dm - the DM object 641 . g - the global vector 642 . mode - INSERT_VALUES or ADD_VALUES 643 - l - the local vector 644 645 646 Level: beginner 647 648 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 649 650 @*/ 651 PetscErrorCode DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) 652 { 653 PetscErrorCode ierr; 654 655 PetscFunctionBegin; 656 ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode,l);CHKERRQ(ierr); 657 PetscFunctionReturn(0); 658 } 659 660 #undef __FUNCT__ 661 #define __FUNCT__ "DMGlobalToLocalEnd" 662 /*@ 663 DMGlobalToLocalEnd - Ends updating local vectors from global vector 664 665 Neighbor-wise Collective on DM 666 667 Input Parameters: 668 + dm - the DM object 669 . g - the global vector 670 . mode - INSERT_VALUES or ADD_VALUES 671 - l - the local vector 672 673 674 Level: beginner 675 676 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 677 678 @*/ 679 PetscErrorCode DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) 680 { 681 PetscErrorCode ierr; 682 683 PetscFunctionBegin; 684 ierr = (*dm->ops->globaltolocalend)(dm,g,mode,l);CHKERRQ(ierr); 685 PetscFunctionReturn(0); 686 } 687 688 #undef __FUNCT__ 689 #define __FUNCT__ "DMLocalToGlobalBegin" 690 /*@ 691 DMLocalToGlobalBegin - updates global vectors from local vectors 692 693 Neighbor-wise Collective on DM 694 695 Input Parameters: 696 + dm - the DM object 697 . l - the local vector 698 . 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 699 base point. 700 - - the global vector 701 702 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 703 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 704 global array to the final global array with VecAXPY(). 705 706 Level: beginner 707 708 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin() 709 710 @*/ 711 PetscErrorCode DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g) 712 { 713 PetscErrorCode ierr; 714 715 PetscFunctionBegin; 716 ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode,g);CHKERRQ(ierr); 717 PetscFunctionReturn(0); 718 } 719 720 #undef __FUNCT__ 721 #define __FUNCT__ "DMLocalToGlobalEnd" 722 /*@ 723 DMLocalToGlobalEnd - updates global vectors from local vectors 724 725 Neighbor-wise Collective on DM 726 727 Input Parameters: 728 + dm - the DM object 729 . l - the local vector 730 . mode - INSERT_VALUES or ADD_VALUES 731 - g - the global vector 732 733 734 Level: beginner 735 736 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd() 737 738 @*/ 739 PetscErrorCode DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g) 740 { 741 PetscErrorCode ierr; 742 743 PetscFunctionBegin; 744 ierr = (*dm->ops->localtoglobalend)(dm,l,mode,g);CHKERRQ(ierr); 745 PetscFunctionReturn(0); 746 } 747 748 #undef __FUNCT__ 749 #define __FUNCT__ "DMComputeJacobianDefault" 750 /*@ 751 DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided 752 753 Collective on DM 754 755 Input Parameter: 756 + dm - the DM object 757 . x - location to compute Jacobian at; may be ignored for linear problems 758 . A - matrix that defines the operator for the linear solve 759 - B - the matrix used to construct the preconditioner 760 761 Level: developer 762 763 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 764 DMSetFunction() 765 766 @*/ 767 PetscErrorCode DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 768 { 769 PetscErrorCode ierr; 770 PetscFunctionBegin; 771 *stflag = SAME_NONZERO_PATTERN; 772 ierr = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr); 773 if (A != B) { 774 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 775 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 776 } 777 PetscFunctionReturn(0); 778 } 779 780 #undef __FUNCT__ 781 #define __FUNCT__ "DMCoarsen" 782 /*@ 783 DMCoarsen - Coarsens a DM object 784 785 Collective on DM 786 787 Input Parameter: 788 + dm - the DM object 789 - comm - the communicator to contain the new DM object (or PETSC_NULL) 790 791 Output Parameter: 792 . dmc - the coarsened DM 793 794 Level: developer 795 796 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 797 798 @*/ 799 PetscErrorCode DMCoarsen(DM dm, MPI_Comm comm, DM *dmc) 800 { 801 PetscErrorCode ierr; 802 803 PetscFunctionBegin; 804 ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr); 805 (*dmc)->ops->initialguess = dm->ops->initialguess; 806 (*dmc)->ops->function = dm->ops->function; 807 (*dmc)->ops->functionj = dm->ops->functionj; 808 if (dm->ops->jacobian != DMComputeJacobianDefault) { 809 (*dmc)->ops->jacobian = dm->ops->jacobian; 810 } 811 (*dmc)->leveldown = dm->leveldown + 1; 812 PetscFunctionReturn(0); 813 } 814 815 #undef __FUNCT__ 816 #define __FUNCT__ "DMRefineHierarchy" 817 /*@C 818 DMRefineHierarchy - Refines a DM object, all levels at once 819 820 Collective on DM 821 822 Input Parameter: 823 + dm - the DM object 824 - nlevels - the number of levels of refinement 825 826 Output Parameter: 827 . dmf - the refined DM hierarchy 828 829 Level: developer 830 831 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 832 833 @*/ 834 PetscErrorCode DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[]) 835 { 836 PetscErrorCode ierr; 837 838 PetscFunctionBegin; 839 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 840 if (nlevels == 0) PetscFunctionReturn(0); 841 if (dm->ops->refinehierarchy) { 842 ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr); 843 } else if (dm->ops->refine) { 844 PetscInt i; 845 846 ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr); 847 for (i=1; i<nlevels; i++) { 848 ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr); 849 } 850 } else { 851 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet"); 852 } 853 PetscFunctionReturn(0); 854 } 855 856 #undef __FUNCT__ 857 #define __FUNCT__ "DMCoarsenHierarchy" 858 /*@C 859 DMCoarsenHierarchy - Coarsens a DM object, all levels at once 860 861 Collective on DM 862 863 Input Parameter: 864 + dm - the DM object 865 - nlevels - the number of levels of coarsening 866 867 Output Parameter: 868 . dmc - the coarsened DM hierarchy 869 870 Level: developer 871 872 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 873 874 @*/ 875 PetscErrorCode DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[]) 876 { 877 PetscErrorCode ierr; 878 879 PetscFunctionBegin; 880 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 881 if (nlevels == 0) PetscFunctionReturn(0); 882 PetscValidPointer(dmc,3); 883 if (dm->ops->coarsenhierarchy) { 884 ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr); 885 } else if (dm->ops->coarsen) { 886 PetscInt i; 887 888 ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr); 889 for (i=1; i<nlevels; i++) { 890 ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr); 891 } 892 } else { 893 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet"); 894 } 895 PetscFunctionReturn(0); 896 } 897 898 #undef __FUNCT__ 899 #define __FUNCT__ "DMGetAggregates" 900 /*@ 901 DMGetAggregates - Gets the aggregates that map between 902 grids associated with two DMs. 903 904 Collective on DM 905 906 Input Parameters: 907 + dmc - the coarse grid DM 908 - dmf - the fine grid DM 909 910 Output Parameters: 911 . rest - the restriction matrix (transpose of the projection matrix) 912 913 Level: intermediate 914 915 .keywords: interpolation, restriction, multigrid 916 917 .seealso: DMRefine(), DMGetInjection(), DMGetInterpolation() 918 @*/ 919 PetscErrorCode DMGetAggregates(DM dmc, DM dmf, Mat *rest) 920 { 921 PetscErrorCode ierr; 922 923 PetscFunctionBegin; 924 ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr); 925 PetscFunctionReturn(0); 926 } 927 928 #undef __FUNCT__ 929 #define __FUNCT__ "DMSetApplicationContext" 930 /*@C 931 DMSetApplicationContext - Set a user context into a DM object 932 933 Not Collective 934 935 Input Parameters: 936 + dm - the DM object 937 - ctx - the user context 938 939 Level: intermediate 940 941 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext() 942 943 @*/ 944 PetscErrorCode DMSetApplicationContext(DM dm,void *ctx) 945 { 946 PetscFunctionBegin; 947 dm->ctx = ctx; 948 PetscFunctionReturn(0); 949 } 950 951 #undef __FUNCT__ 952 #define __FUNCT__ "DMGetApplicationContext" 953 /*@ 954 DMGetApplicationContext - Gets a user context from a DM object 955 956 Not Collective 957 958 Input Parameter: 959 . dm - the DM object 960 961 Output Parameter: 962 . ctx - the user context 963 964 Level: intermediate 965 966 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext() 967 968 @*/ 969 PetscErrorCode DMGetApplicationContext(DM dm,void *ctx) 970 { 971 PetscFunctionBegin; 972 *(void**)ctx = dm->ctx; 973 PetscFunctionReturn(0); 974 } 975 976 #undef __FUNCT__ 977 #define __FUNCT__ "DMSetInitialGuess" 978 /*@ 979 DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers 980 981 Logically Collective on DM 982 983 Input Parameter: 984 + dm - the DM object to destroy 985 - f - the function to compute the initial guess 986 987 Level: intermediate 988 989 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 990 991 @*/ 992 PetscErrorCode DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec)) 993 { 994 PetscFunctionBegin; 995 dm->ops->initialguess = f; 996 PetscFunctionReturn(0); 997 } 998 999 #undef __FUNCT__ 1000 #define __FUNCT__ "DMSetFunction" 1001 /*@ 1002 DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES 1003 1004 Logically Collective on DM 1005 1006 Input Parameter: 1007 + dm - the DM object 1008 - f - the function to compute (use PETSC_NULL to cancel a previous function that was set) 1009 1010 Level: intermediate 1011 1012 Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian 1013 computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian. 1014 1015 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1016 DMSetJacobian() 1017 1018 @*/ 1019 PetscErrorCode DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 1020 { 1021 PetscFunctionBegin; 1022 dm->ops->function = f; 1023 if (f) { 1024 dm->ops->functionj = f; 1025 } 1026 PetscFunctionReturn(0); 1027 } 1028 1029 #undef __FUNCT__ 1030 #define __FUNCT__ "DMSetJacobian" 1031 /*@ 1032 DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES 1033 1034 Logically Collective on DM 1035 1036 Input Parameter: 1037 + dm - the DM object to destroy 1038 - f - the function to compute the matrix entries 1039 1040 Level: intermediate 1041 1042 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1043 DMSetFunction() 1044 1045 @*/ 1046 PetscErrorCode DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*)) 1047 { 1048 PetscFunctionBegin; 1049 dm->ops->jacobian = f; 1050 PetscFunctionReturn(0); 1051 } 1052 1053 #undef __FUNCT__ 1054 #define __FUNCT__ "DMComputeInitialGuess" 1055 /*@ 1056 DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers 1057 1058 Collective on DM 1059 1060 Input Parameter: 1061 + dm - the DM object to destroy 1062 - x - the vector to hold the initial guess values 1063 1064 Level: developer 1065 1066 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat() 1067 1068 @*/ 1069 PetscErrorCode DMComputeInitialGuess(DM dm,Vec x) 1070 { 1071 PetscErrorCode ierr; 1072 PetscFunctionBegin; 1073 if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()"); 1074 ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr); 1075 PetscFunctionReturn(0); 1076 } 1077 1078 #undef __FUNCT__ 1079 #define __FUNCT__ "DMHasInitialGuess" 1080 /*@ 1081 DMHasInitialGuess - does the DM object have an initial guess function 1082 1083 Not Collective 1084 1085 Input Parameter: 1086 . dm - the DM object to destroy 1087 1088 Output Parameter: 1089 . flg - PETSC_TRUE if function exists 1090 1091 Level: developer 1092 1093 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1094 1095 @*/ 1096 PetscErrorCode DMHasInitialGuess(DM dm,PetscBool *flg) 1097 { 1098 PetscFunctionBegin; 1099 *flg = (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE; 1100 PetscFunctionReturn(0); 1101 } 1102 1103 #undef __FUNCT__ 1104 #define __FUNCT__ "DMHasFunction" 1105 /*@ 1106 DMHasFunction - does the DM object have a function 1107 1108 Not Collective 1109 1110 Input Parameter: 1111 . dm - the DM object to destroy 1112 1113 Output Parameter: 1114 . flg - PETSC_TRUE if function exists 1115 1116 Level: developer 1117 1118 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1119 1120 @*/ 1121 PetscErrorCode DMHasFunction(DM dm,PetscBool *flg) 1122 { 1123 PetscFunctionBegin; 1124 *flg = (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE; 1125 PetscFunctionReturn(0); 1126 } 1127 1128 #undef __FUNCT__ 1129 #define __FUNCT__ "DMHasJacobian" 1130 /*@ 1131 DMHasJacobian - does the DM object have a matrix function 1132 1133 Not Collective 1134 1135 Input Parameter: 1136 . dm - the DM object to destroy 1137 1138 Output Parameter: 1139 . flg - PETSC_TRUE if function exists 1140 1141 Level: developer 1142 1143 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1144 1145 @*/ 1146 PetscErrorCode DMHasJacobian(DM dm,PetscBool *flg) 1147 { 1148 PetscFunctionBegin; 1149 *flg = (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE; 1150 PetscFunctionReturn(0); 1151 } 1152 1153 #undef __FUNCT__ 1154 #define __FUNCT__ "DMComputeFunction" 1155 /*@ 1156 DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES 1157 1158 Collective on DM 1159 1160 Input Parameter: 1161 + dm - the DM object to destroy 1162 . x - the location where the function is evaluationed, may be ignored for linear problems 1163 - b - the vector to hold the right hand side entries 1164 1165 Level: developer 1166 1167 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1168 DMSetJacobian() 1169 1170 @*/ 1171 PetscErrorCode DMComputeFunction(DM dm,Vec x,Vec b) 1172 { 1173 PetscErrorCode ierr; 1174 PetscFunctionBegin; 1175 if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()"); 1176 if (!x) x = dm->x; 1177 ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr); 1178 PetscFunctionReturn(0); 1179 } 1180 1181 1182 #undef __FUNCT__ 1183 #define __FUNCT__ "DMComputeJacobian" 1184 /*@ 1185 DMComputeJacobian - compute the matrix entries for the solver 1186 1187 Collective on DM 1188 1189 Input Parameter: 1190 + dm - the DM object 1191 . x - location to compute Jacobian at; may be ignored for linear problems 1192 . A - matrix that defines the operator for the linear solve 1193 - B - the matrix used to construct the preconditioner 1194 1195 Level: developer 1196 1197 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1198 DMSetFunction() 1199 1200 @*/ 1201 PetscErrorCode DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 1202 { 1203 PetscErrorCode ierr; 1204 1205 PetscFunctionBegin; 1206 if (!dm->ops->jacobian) { 1207 ISColoring coloring; 1208 MatFDColoring fd; 1209 1210 ierr = DMGetColoring(dm,IS_COLORING_GHOSTED,MATAIJ,&coloring);CHKERRQ(ierr); 1211 ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr); 1212 ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 1213 ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr); 1214 dm->fd = fd; 1215 dm->ops->jacobian = DMComputeJacobianDefault; 1216 1217 if (!dm->x) { 1218 ierr = MatGetVecs(B,&dm->x,PETSC_NULL);CHKERRQ(ierr); 1219 } 1220 } 1221 if (!x) x = dm->x; 1222 ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr); 1223 PetscFunctionReturn(0); 1224 } 1225 1226 PetscFList DMList = PETSC_NULL; 1227 PetscBool DMRegisterAllCalled = PETSC_FALSE; 1228 1229 #undef __FUNCT__ 1230 #define __FUNCT__ "DMSetType" 1231 /*@C 1232 DMSetType - Builds a DM, for a particular DM implementation. 1233 1234 Collective on DM 1235 1236 Input Parameters: 1237 + dm - The DM object 1238 - method - The name of the DM type 1239 1240 Options Database Key: 1241 . -dm_type <type> - Sets the DM type; use -help for a list of available types 1242 1243 Notes: 1244 See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D). 1245 1246 Level: intermediate 1247 1248 .keywords: DM, set, type 1249 .seealso: DMGetType(), DMCreate() 1250 @*/ 1251 PetscErrorCode DMSetType(DM dm, const DMType method) 1252 { 1253 PetscErrorCode (*r)(DM); 1254 PetscBool match; 1255 PetscErrorCode ierr; 1256 1257 PetscFunctionBegin; 1258 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1259 ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr); 1260 if (match) PetscFunctionReturn(0); 1261 1262 if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1263 ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 1264 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method); 1265 1266 if (dm->ops->destroy) { 1267 ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr); 1268 } 1269 ierr = (*r)(dm);CHKERRQ(ierr); 1270 ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr); 1271 PetscFunctionReturn(0); 1272 } 1273 1274 #undef __FUNCT__ 1275 #define __FUNCT__ "DMGetType" 1276 /*@C 1277 DMGetType - Gets the DM type name (as a string) from the DM. 1278 1279 Not Collective 1280 1281 Input Parameter: 1282 . dm - The DM 1283 1284 Output Parameter: 1285 . type - The DM type name 1286 1287 Level: intermediate 1288 1289 .keywords: DM, get, type, name 1290 .seealso: DMSetType(), DMCreate() 1291 @*/ 1292 PetscErrorCode DMGetType(DM dm, const DMType *type) 1293 { 1294 PetscErrorCode ierr; 1295 1296 PetscFunctionBegin; 1297 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1298 PetscValidCharPointer(type,2); 1299 if (!DMRegisterAllCalled) { 1300 ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1301 } 1302 *type = ((PetscObject)dm)->type_name; 1303 PetscFunctionReturn(0); 1304 } 1305 1306 #undef __FUNCT__ 1307 #define __FUNCT__ "DMConvert" 1308 /*@C 1309 DMConvert - Converts a DM to another DM, either of the same or different type. 1310 1311 Collective on DM 1312 1313 Input Parameters: 1314 + dm - the DM 1315 - newtype - new DM type (use "same" for the same type) 1316 1317 Output Parameter: 1318 . M - pointer to new DM 1319 1320 Notes: 1321 Cannot be used to convert a sequential DM to parallel or parallel to sequential, 1322 the MPI communicator of the generated DM is always the same as the communicator 1323 of the input DM. 1324 1325 Level: intermediate 1326 1327 .seealso: DMCreate() 1328 @*/ 1329 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M) 1330 { 1331 DM B; 1332 char convname[256]; 1333 PetscBool sametype, issame; 1334 PetscErrorCode ierr; 1335 1336 PetscFunctionBegin; 1337 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1338 PetscValidType(dm,1); 1339 PetscValidPointer(M,3); 1340 ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr); 1341 ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr); 1342 { 1343 PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL; 1344 1345 /* 1346 Order of precedence: 1347 1) See if a specialized converter is known to the current DM. 1348 2) See if a specialized converter is known to the desired DM class. 1349 3) See if a good general converter is registered for the desired class 1350 4) See if a good general converter is known for the current matrix. 1351 5) Use a really basic converter. 1352 */ 1353 1354 /* 1) See if a specialized converter is known to the current DM and the desired class */ 1355 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 1356 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 1357 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 1358 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 1359 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 1360 ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr); 1361 if (conv) goto foundconv; 1362 1363 /* 2) See if a specialized converter is known to the desired DM class. */ 1364 ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr); 1365 ierr = DMSetType(B, newtype);CHKERRQ(ierr); 1366 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 1367 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 1368 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 1369 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 1370 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 1371 ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr); 1372 if (conv) { 1373 ierr = DMDestroy(&B);CHKERRQ(ierr); 1374 goto foundconv; 1375 } 1376 1377 #if 0 1378 /* 3) See if a good general converter is registered for the desired class */ 1379 conv = B->ops->convertfrom; 1380 ierr = DMDestroy(&B);CHKERRQ(ierr); 1381 if (conv) goto foundconv; 1382 1383 /* 4) See if a good general converter is known for the current matrix */ 1384 if (dm->ops->convert) { 1385 conv = dm->ops->convert; 1386 } 1387 if (conv) goto foundconv; 1388 #endif 1389 1390 /* 5) Use a really basic converter. */ 1391 SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype); 1392 1393 foundconv: 1394 ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 1395 ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr); 1396 ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 1397 } 1398 ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr); 1399 PetscFunctionReturn(0); 1400 } 1401 1402 /*--------------------------------------------------------------------------------------------------------------------*/ 1403 1404 #undef __FUNCT__ 1405 #define __FUNCT__ "DMRegister" 1406 /*@C 1407 DMRegister - See DMRegisterDynamic() 1408 1409 Level: advanced 1410 @*/ 1411 PetscErrorCode DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM)) 1412 { 1413 char fullname[PETSC_MAX_PATH_LEN]; 1414 PetscErrorCode ierr; 1415 1416 PetscFunctionBegin; 1417 ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr); 1418 ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr); 1419 ierr = PetscStrcat(fullname, name);CHKERRQ(ierr); 1420 ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr); 1421 PetscFunctionReturn(0); 1422 } 1423 1424 1425 /*--------------------------------------------------------------------------------------------------------------------*/ 1426 #undef __FUNCT__ 1427 #define __FUNCT__ "DMRegisterDestroy" 1428 /*@C 1429 DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic(). 1430 1431 Not Collective 1432 1433 Level: advanced 1434 1435 .keywords: DM, register, destroy 1436 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic() 1437 @*/ 1438 PetscErrorCode DMRegisterDestroy(void) 1439 { 1440 PetscErrorCode ierr; 1441 1442 PetscFunctionBegin; 1443 ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr); 1444 DMRegisterAllCalled = PETSC_FALSE; 1445 PetscFunctionReturn(0); 1446 } 1447 1448 #if defined(PETSC_HAVE_MATLAB_ENGINE) 1449 #include <mex.h> 1450 1451 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext; 1452 1453 #undef __FUNCT__ 1454 #define __FUNCT__ "DMComputeFunction_Matlab" 1455 /* 1456 DMComputeFunction_Matlab - Calls the function that has been set with 1457 DMSetFunctionMatlab(). 1458 1459 For linear problems x is null 1460 1461 .seealso: DMSetFunction(), DMGetFunction() 1462 */ 1463 PetscErrorCode DMComputeFunction_Matlab(DM dm,Vec x,Vec y) 1464 { 1465 PetscErrorCode ierr; 1466 DMMatlabContext *sctx; 1467 int nlhs = 1,nrhs = 4; 1468 mxArray *plhs[1],*prhs[4]; 1469 long long int lx = 0,ly = 0,ls = 0; 1470 1471 PetscFunctionBegin; 1472 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1473 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 1474 PetscCheckSameComm(dm,1,y,3); 1475 1476 /* call Matlab function in ctx with arguments x and y */ 1477 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1478 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 1479 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 1480 ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr); 1481 prhs[0] = mxCreateDoubleScalar((double)ls); 1482 prhs[1] = mxCreateDoubleScalar((double)lx); 1483 prhs[2] = mxCreateDoubleScalar((double)ly); 1484 prhs[3] = mxCreateString(sctx->funcname); 1485 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr); 1486 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 1487 mxDestroyArray(prhs[0]); 1488 mxDestroyArray(prhs[1]); 1489 mxDestroyArray(prhs[2]); 1490 mxDestroyArray(prhs[3]); 1491 mxDestroyArray(plhs[0]); 1492 PetscFunctionReturn(0); 1493 } 1494 1495 1496 #undef __FUNCT__ 1497 #define __FUNCT__ "DMSetFunctionMatlab" 1498 /* 1499 DMSetFunctionMatlab - Sets the function evaluation routine 1500 1501 */ 1502 PetscErrorCode DMSetFunctionMatlab(DM dm,const char *func) 1503 { 1504 PetscErrorCode ierr; 1505 DMMatlabContext *sctx; 1506 1507 PetscFunctionBegin; 1508 /* currently sctx is memory bleed */ 1509 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1510 if (!sctx) { 1511 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 1512 } 1513 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 1514 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 1515 ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr); 1516 PetscFunctionReturn(0); 1517 } 1518 1519 #undef __FUNCT__ 1520 #define __FUNCT__ "DMComputeJacobian_Matlab" 1521 /* 1522 DMComputeJacobian_Matlab - Calls the function that has been set with 1523 DMSetJacobianMatlab(). 1524 1525 For linear problems x is null 1526 1527 .seealso: DMSetFunction(), DMGetFunction() 1528 */ 1529 PetscErrorCode DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str) 1530 { 1531 PetscErrorCode ierr; 1532 DMMatlabContext *sctx; 1533 int nlhs = 2,nrhs = 5; 1534 mxArray *plhs[2],*prhs[5]; 1535 long long int lx = 0,lA = 0,lB = 0,ls = 0; 1536 1537 PetscFunctionBegin; 1538 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1539 PetscValidHeaderSpecific(A,MAT_CLASSID,3); 1540 1541 /* call MATLAB function in ctx with arguments x, A, and B */ 1542 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1543 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 1544 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 1545 ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr); 1546 ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr); 1547 prhs[0] = mxCreateDoubleScalar((double)ls); 1548 prhs[1] = mxCreateDoubleScalar((double)lx); 1549 prhs[2] = mxCreateDoubleScalar((double)lA); 1550 prhs[3] = mxCreateDoubleScalar((double)lB); 1551 prhs[4] = mxCreateString(sctx->jacname); 1552 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr); 1553 *str = (MatStructure) mxGetScalar(plhs[0]); 1554 ierr = (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr); 1555 mxDestroyArray(prhs[0]); 1556 mxDestroyArray(prhs[1]); 1557 mxDestroyArray(prhs[2]); 1558 mxDestroyArray(prhs[3]); 1559 mxDestroyArray(prhs[4]); 1560 mxDestroyArray(plhs[0]); 1561 mxDestroyArray(plhs[1]); 1562 PetscFunctionReturn(0); 1563 } 1564 1565 1566 #undef __FUNCT__ 1567 #define __FUNCT__ "DMSetJacobianMatlab" 1568 /* 1569 DMSetJacobianMatlab - Sets the Jacobian function evaluation routine 1570 1571 */ 1572 PetscErrorCode DMSetJacobianMatlab(DM dm,const char *func) 1573 { 1574 PetscErrorCode ierr; 1575 DMMatlabContext *sctx; 1576 1577 PetscFunctionBegin; 1578 /* currently sctx is memory bleed */ 1579 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1580 if (!sctx) { 1581 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 1582 } 1583 ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr); 1584 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 1585 ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr); 1586 PetscFunctionReturn(0); 1587 } 1588 #endif 1589