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 C@*/ 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(), DMGetMatrix() 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 PetscFunctionReturn(0); 603 } 604 605 #undef __FUNCT__ 606 #define __FUNCT__ "DMGlobalToLocalBegin" 607 /*@ 608 DMGlobalToLocalBegin - Begins updating local vectors from global vector 609 610 Neighbor-wise Collective on DM 611 612 Input Parameters: 613 + dm - the DM object 614 . g - the global vector 615 . mode - INSERT_VALUES or ADD_VALUES 616 - l - the local vector 617 618 619 Level: beginner 620 621 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 622 623 @*/ 624 PetscErrorCode DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) 625 { 626 PetscErrorCode ierr; 627 628 PetscFunctionBegin; 629 ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode,l);CHKERRQ(ierr); 630 PetscFunctionReturn(0); 631 } 632 633 #undef __FUNCT__ 634 #define __FUNCT__ "DMGlobalToLocalEnd" 635 /*@ 636 DMGlobalToLocalEnd - Ends updating local vectors from global vector 637 638 Neighbor-wise Collective on DM 639 640 Input Parameters: 641 + dm - the DM object 642 . g - the global vector 643 . mode - INSERT_VALUES or ADD_VALUES 644 - l - the local vector 645 646 647 Level: beginner 648 649 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 650 651 @*/ 652 PetscErrorCode DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) 653 { 654 PetscErrorCode ierr; 655 656 PetscFunctionBegin; 657 ierr = (*dm->ops->globaltolocalend)(dm,g,mode,l);CHKERRQ(ierr); 658 PetscFunctionReturn(0); 659 } 660 661 #undef __FUNCT__ 662 #define __FUNCT__ "DMLocalToGlobalBegin" 663 /*@ 664 DMLocalToGlobalBegin - updates global vectors from local vectors 665 666 Neighbor-wise Collective on DM 667 668 Input Parameters: 669 + dm - the DM object 670 . l - the local vector 671 . 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 672 base point. 673 - - the global vector 674 675 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 676 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 677 global array to the final global array with VecAXPY(). 678 679 Level: beginner 680 681 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin() 682 683 @*/ 684 PetscErrorCode DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g) 685 { 686 PetscErrorCode ierr; 687 688 PetscFunctionBegin; 689 ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode,g);CHKERRQ(ierr); 690 PetscFunctionReturn(0); 691 } 692 693 #undef __FUNCT__ 694 #define __FUNCT__ "DMLocalToGlobalEnd" 695 /*@ 696 DMLocalToGlobalEnd - updates global vectors from local vectors 697 698 Neighbor-wise Collective on DM 699 700 Input Parameters: 701 + dm - the DM object 702 . l - the local vector 703 . mode - INSERT_VALUES or ADD_VALUES 704 - g - the global vector 705 706 707 Level: beginner 708 709 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd() 710 711 @*/ 712 PetscErrorCode DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g) 713 { 714 PetscErrorCode ierr; 715 716 PetscFunctionBegin; 717 ierr = (*dm->ops->localtoglobalend)(dm,l,mode,g);CHKERRQ(ierr); 718 PetscFunctionReturn(0); 719 } 720 721 #undef __FUNCT__ 722 #define __FUNCT__ "DMComputeJacobianDefault" 723 /*@ 724 DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided 725 726 Collective on DM 727 728 Input Parameter: 729 + dm - the DM object 730 . x - location to compute Jacobian at; may be ignored for linear problems 731 . A - matrix that defines the operator for the linear solve 732 - B - the matrix used to construct the preconditioner 733 734 Level: developer 735 736 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(), 737 DMSetFunction() 738 739 @*/ 740 PetscErrorCode DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 741 { 742 PetscErrorCode ierr; 743 PetscFunctionBegin; 744 *stflag = SAME_NONZERO_PATTERN; 745 ierr = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr); 746 if (A != B) { 747 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 748 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 749 } 750 PetscFunctionReturn(0); 751 } 752 753 #undef __FUNCT__ 754 #define __FUNCT__ "DMCoarsen" 755 /*@ 756 DMCoarsen - Coarsens a DM object 757 758 Collective on DM 759 760 Input Parameter: 761 + dm - the DM object 762 - comm - the communicator to contain the new DM object (or PETSC_NULL) 763 764 Output Parameter: 765 . dmc - the coarsened DM 766 767 Level: developer 768 769 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 770 771 @*/ 772 PetscErrorCode DMCoarsen(DM dm, MPI_Comm comm, DM *dmc) 773 { 774 PetscErrorCode ierr; 775 776 PetscFunctionBegin; 777 ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr); 778 (*dmc)->ops->initialguess = dm->ops->initialguess; 779 (*dmc)->ops->function = dm->ops->function; 780 (*dmc)->ops->functionj = dm->ops->functionj; 781 if (dm->ops->jacobian != DMComputeJacobianDefault) { 782 (*dmc)->ops->jacobian = dm->ops->jacobian; 783 } 784 PetscFunctionReturn(0); 785 } 786 787 #undef __FUNCT__ 788 #define __FUNCT__ "DMRefineHierarchy" 789 /*@C 790 DMRefineHierarchy - Refines a DM object, all levels at once 791 792 Collective on DM 793 794 Input Parameter: 795 + dm - the DM object 796 - nlevels - the number of levels of refinement 797 798 Output Parameter: 799 . dmf - the refined DM hierarchy 800 801 Level: developer 802 803 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 804 805 @*/ 806 PetscErrorCode DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[]) 807 { 808 PetscErrorCode ierr; 809 810 PetscFunctionBegin; 811 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 812 if (nlevels == 0) PetscFunctionReturn(0); 813 if (dm->ops->refinehierarchy) { 814 ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr); 815 } else if (dm->ops->refine) { 816 PetscInt i; 817 818 ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr); 819 for (i=1; i<nlevels; i++) { 820 ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr); 821 } 822 } else { 823 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet"); 824 } 825 PetscFunctionReturn(0); 826 } 827 828 #undef __FUNCT__ 829 #define __FUNCT__ "DMCoarsenHierarchy" 830 /*@C 831 DMCoarsenHierarchy - Coarsens a DM object, all levels at once 832 833 Collective on DM 834 835 Input Parameter: 836 + dm - the DM object 837 - nlevels - the number of levels of coarsening 838 839 Output Parameter: 840 . dmc - the coarsened DM hierarchy 841 842 Level: developer 843 844 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 845 846 @*/ 847 PetscErrorCode DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[]) 848 { 849 PetscErrorCode ierr; 850 851 PetscFunctionBegin; 852 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 853 if (nlevels == 0) PetscFunctionReturn(0); 854 PetscValidPointer(dmc,3); 855 if (dm->ops->coarsenhierarchy) { 856 ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr); 857 } else if (dm->ops->coarsen) { 858 PetscInt i; 859 860 ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr); 861 for (i=1; i<nlevels; i++) { 862 ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr); 863 } 864 } else { 865 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet"); 866 } 867 PetscFunctionReturn(0); 868 } 869 870 #undef __FUNCT__ 871 #define __FUNCT__ "DMGetAggregates" 872 /*@ 873 DMGetAggregates - Gets the aggregates that map between 874 grids associated with two DMs. 875 876 Collective on DM 877 878 Input Parameters: 879 + dmc - the coarse grid DM 880 - dmf - the fine grid DM 881 882 Output Parameters: 883 . rest - the restriction matrix (transpose of the projection matrix) 884 885 Level: intermediate 886 887 .keywords: interpolation, restriction, multigrid 888 889 .seealso: DMRefine(), DMGetInjection(), DMGetInterpolation() 890 @*/ 891 PetscErrorCode DMGetAggregates(DM dmc, DM dmf, Mat *rest) 892 { 893 PetscErrorCode ierr; 894 895 PetscFunctionBegin; 896 ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr); 897 PetscFunctionReturn(0); 898 } 899 900 #undef __FUNCT__ 901 #define __FUNCT__ "DMSetContext" 902 /*@ 903 DMSetContext - Set a user context into a DM object 904 905 Not Collective 906 907 Input Parameters: 908 + dm - the DM object 909 - ctx - the user context 910 911 Level: intermediate 912 913 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext() 914 915 @*/ 916 PetscErrorCode DMSetContext(DM dm,void *ctx) 917 { 918 PetscFunctionBegin; 919 dm->ctx = ctx; 920 PetscFunctionReturn(0); 921 } 922 923 #undef __FUNCT__ 924 #define __FUNCT__ "DMGetContext" 925 /*@ 926 DMGetContext - Gets a user context from a DM object 927 928 Not Collective 929 930 Input Parameter: 931 . dm - the DM object 932 933 Output Parameter: 934 . ctx - the user context 935 936 Level: intermediate 937 938 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext() 939 940 @*/ 941 PetscErrorCode DMGetContext(DM dm,void **ctx) 942 { 943 PetscFunctionBegin; 944 *ctx = dm->ctx; 945 PetscFunctionReturn(0); 946 } 947 948 #undef __FUNCT__ 949 #define __FUNCT__ "DMSetInitialGuess" 950 /*@ 951 DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers 952 953 Logically Collective on DM 954 955 Input Parameter: 956 + dm - the DM object to destroy 957 - f - the function to compute the initial guess 958 959 Level: intermediate 960 961 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian() 962 963 @*/ 964 PetscErrorCode DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec)) 965 { 966 PetscFunctionBegin; 967 dm->ops->initialguess = f; 968 PetscFunctionReturn(0); 969 } 970 971 #undef __FUNCT__ 972 #define __FUNCT__ "DMSetFunction" 973 /*@ 974 DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES 975 976 Logically Collective on DM 977 978 Input Parameter: 979 + dm - the DM object 980 - f - the function to compute (use PETSC_NULL to cancel a previous function that was set) 981 982 Level: intermediate 983 984 Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian 985 computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian. 986 987 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(), 988 DMSetJacobian() 989 990 @*/ 991 PetscErrorCode DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 992 { 993 PetscFunctionBegin; 994 dm->ops->function = f; 995 if (f) { 996 dm->ops->functionj = f; 997 } 998 PetscFunctionReturn(0); 999 } 1000 1001 #undef __FUNCT__ 1002 #define __FUNCT__ "DMSetJacobian" 1003 /*@ 1004 DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES 1005 1006 Logically Collective on DM 1007 1008 Input Parameter: 1009 + dm - the DM object to destroy 1010 - f - the function to compute the matrix entries 1011 1012 Level: intermediate 1013 1014 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(), 1015 DMSetFunction() 1016 1017 @*/ 1018 PetscErrorCode DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*)) 1019 { 1020 PetscFunctionBegin; 1021 dm->ops->jacobian = f; 1022 PetscFunctionReturn(0); 1023 } 1024 1025 #undef __FUNCT__ 1026 #define __FUNCT__ "DMComputeInitialGuess" 1027 /*@ 1028 DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers 1029 1030 Collective on DM 1031 1032 Input Parameter: 1033 + dm - the DM object to destroy 1034 - x - the vector to hold the initial guess values 1035 1036 Level: developer 1037 1038 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetRhs(), DMSetMat() 1039 1040 @*/ 1041 PetscErrorCode DMComputeInitialGuess(DM dm,Vec x) 1042 { 1043 PetscErrorCode ierr; 1044 PetscFunctionBegin; 1045 if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()"); 1046 ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr); 1047 PetscFunctionReturn(0); 1048 } 1049 1050 #undef __FUNCT__ 1051 #define __FUNCT__ "DMHasInitialGuess" 1052 /*@ 1053 DMHasInitialGuess - does the DM object have an initial guess function 1054 1055 Not Collective 1056 1057 Input Parameter: 1058 . dm - the DM object to destroy 1059 1060 Output Parameter: 1061 . flg - PETSC_TRUE if function exists 1062 1063 Level: developer 1064 1065 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian() 1066 1067 @*/ 1068 PetscErrorCode DMHasInitialGuess(DM dm,PetscBool *flg) 1069 { 1070 PetscFunctionBegin; 1071 *flg = (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE; 1072 PetscFunctionReturn(0); 1073 } 1074 1075 #undef __FUNCT__ 1076 #define __FUNCT__ "DMHasFunction" 1077 /*@ 1078 DMHasFunction - does the DM object have a function 1079 1080 Not Collective 1081 1082 Input Parameter: 1083 . dm - the DM object to destroy 1084 1085 Output Parameter: 1086 . flg - PETSC_TRUE if function exists 1087 1088 Level: developer 1089 1090 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian() 1091 1092 @*/ 1093 PetscErrorCode DMHasFunction(DM dm,PetscBool *flg) 1094 { 1095 PetscFunctionBegin; 1096 *flg = (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE; 1097 PetscFunctionReturn(0); 1098 } 1099 1100 #undef __FUNCT__ 1101 #define __FUNCT__ "DMHasJacobian" 1102 /*@ 1103 DMHasJacobian - does the DM object have a matrix function 1104 1105 Not Collective 1106 1107 Input Parameter: 1108 . dm - the DM object to destroy 1109 1110 Output Parameter: 1111 . flg - PETSC_TRUE if function exists 1112 1113 Level: developer 1114 1115 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetFunction(), DMSetJacobian() 1116 1117 @*/ 1118 PetscErrorCode DMHasJacobian(DM dm,PetscBool *flg) 1119 { 1120 PetscFunctionBegin; 1121 *flg = (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE; 1122 PetscFunctionReturn(0); 1123 } 1124 1125 #undef __FUNCT__ 1126 #define __FUNCT__ "DMComputeFunction" 1127 /*@ 1128 DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES 1129 1130 Collective on DM 1131 1132 Input Parameter: 1133 + dm - the DM object to destroy 1134 . x - the location where the function is evaluationed, may be ignored for linear problems 1135 - b - the vector to hold the right hand side entries 1136 1137 Level: developer 1138 1139 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(), 1140 DMSetJacobian() 1141 1142 @*/ 1143 PetscErrorCode DMComputeFunction(DM dm,Vec x,Vec b) 1144 { 1145 PetscErrorCode ierr; 1146 PetscFunctionBegin; 1147 if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()"); 1148 if (!x) x = dm->x; 1149 ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr); 1150 PetscFunctionReturn(0); 1151 } 1152 1153 1154 #undef __FUNCT__ 1155 #define __FUNCT__ "DMComputeJacobian" 1156 /*@ 1157 DMComputeJacobian - compute the matrix entries for the solver 1158 1159 Collective on DM 1160 1161 Input Parameter: 1162 + dm - the DM object 1163 . x - location to compute Jacobian at; may be ignored for linear problems 1164 . A - matrix that defines the operator for the linear solve 1165 - B - the matrix used to construct the preconditioner 1166 1167 Level: developer 1168 1169 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetContext(), DMSetInitialGuess(), 1170 DMSetFunction() 1171 1172 @*/ 1173 PetscErrorCode DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 1174 { 1175 PetscErrorCode ierr; 1176 1177 PetscFunctionBegin; 1178 if (!dm->ops->jacobian) { 1179 ISColoring coloring; 1180 MatFDColoring fd; 1181 1182 ierr = DMGetColoring(dm,IS_COLORING_GHOSTED,MATAIJ,&coloring);CHKERRQ(ierr); 1183 ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr); 1184 ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 1185 ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr); 1186 dm->fd = fd; 1187 dm->ops->jacobian = DMComputeJacobianDefault; 1188 1189 if (!dm->x) { 1190 ierr = MatGetVecs(B,&dm->x,PETSC_NULL);CHKERRQ(ierr); 1191 } 1192 } 1193 if (!x) x = dm->x; 1194 ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr); 1195 PetscFunctionReturn(0); 1196 } 1197 1198 PetscFList DMList = PETSC_NULL; 1199 PetscBool DMRegisterAllCalled = PETSC_FALSE; 1200 1201 #undef __FUNCT__ 1202 #define __FUNCT__ "DMSetType" 1203 /*@C 1204 DMSetType - Builds a DM, for a particular DM implementation. 1205 1206 Collective on DM 1207 1208 Input Parameters: 1209 + dm - The DM object 1210 - method - The name of the DM type 1211 1212 Options Database Key: 1213 . -dm_type <type> - Sets the DM type; use -help for a list of available types 1214 1215 Notes: 1216 See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D). 1217 1218 Level: intermediate 1219 1220 .keywords: DM, set, type 1221 .seealso: DMGetType(), DMCreate() 1222 @*/ 1223 PetscErrorCode DMSetType(DM dm, const DMType method) 1224 { 1225 PetscErrorCode (*r)(DM); 1226 PetscBool match; 1227 PetscErrorCode ierr; 1228 1229 PetscFunctionBegin; 1230 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1231 ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr); 1232 if (match) PetscFunctionReturn(0); 1233 1234 if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1235 ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 1236 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method); 1237 1238 if (dm->ops->destroy) { 1239 ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr); 1240 } 1241 ierr = (*r)(dm);CHKERRQ(ierr); 1242 ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr); 1243 PetscFunctionReturn(0); 1244 } 1245 1246 #undef __FUNCT__ 1247 #define __FUNCT__ "DMGetType" 1248 /*@C 1249 DMGetType - Gets the DM type name (as a string) from the DM. 1250 1251 Not Collective 1252 1253 Input Parameter: 1254 . dm - The DM 1255 1256 Output Parameter: 1257 . type - The DM type name 1258 1259 Level: intermediate 1260 1261 .keywords: DM, get, type, name 1262 .seealso: DMSetType(), DMCreate() 1263 @*/ 1264 PetscErrorCode DMGetType(DM dm, const DMType *type) 1265 { 1266 PetscErrorCode ierr; 1267 1268 PetscFunctionBegin; 1269 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1270 PetscValidCharPointer(type,2); 1271 if (!DMRegisterAllCalled) { 1272 ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1273 } 1274 *type = ((PetscObject)dm)->type_name; 1275 PetscFunctionReturn(0); 1276 } 1277 1278 #undef __FUNCT__ 1279 #define __FUNCT__ "DMConvert" 1280 /*@C 1281 DMConvert - Converts a DM to another DM, either of the same or different type. 1282 1283 Collective on DM 1284 1285 Input Parameters: 1286 + dm - the DM 1287 - newtype - new DM type (use "same" for the same type) 1288 1289 Output Parameter: 1290 . M - pointer to new DM 1291 1292 Notes: 1293 Cannot be used to convert a sequential DM to parallel or parallel to sequential, 1294 the MPI communicator of the generated DM is always the same as the communicator 1295 of the input DM. 1296 1297 Level: intermediate 1298 1299 .seealso: DMCreate() 1300 @*/ 1301 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M) 1302 { 1303 DM B; 1304 char convname[256]; 1305 PetscBool sametype, issame; 1306 PetscErrorCode ierr; 1307 1308 PetscFunctionBegin; 1309 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1310 PetscValidType(dm,1); 1311 PetscValidPointer(M,3); 1312 ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr); 1313 ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr); 1314 { 1315 PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL; 1316 1317 /* 1318 Order of precedence: 1319 1) See if a specialized converter is known to the current DM. 1320 2) See if a specialized converter is known to the desired DM class. 1321 3) See if a good general converter is registered for the desired class 1322 4) See if a good general converter is known for the current matrix. 1323 5) Use a really basic converter. 1324 */ 1325 1326 /* 1) See if a specialized converter is known to the current DM and the desired class */ 1327 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 1328 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 1329 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 1330 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 1331 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 1332 ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr); 1333 if (conv) goto foundconv; 1334 1335 /* 2) See if a specialized converter is known to the desired DM class. */ 1336 ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr); 1337 ierr = DMSetType(B, newtype);CHKERRQ(ierr); 1338 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 1339 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 1340 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 1341 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 1342 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 1343 ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr); 1344 if (conv) { 1345 ierr = DMDestroy(&B);CHKERRQ(ierr); 1346 goto foundconv; 1347 } 1348 1349 #if 0 1350 /* 3) See if a good general converter is registered for the desired class */ 1351 conv = B->ops->convertfrom; 1352 ierr = DMDestroy(&B);CHKERRQ(ierr); 1353 if (conv) goto foundconv; 1354 1355 /* 4) See if a good general converter is known for the current matrix */ 1356 if (dm->ops->convert) { 1357 conv = dm->ops->convert; 1358 } 1359 if (conv) goto foundconv; 1360 #endif 1361 1362 /* 5) Use a really basic converter. */ 1363 SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype); 1364 1365 foundconv: 1366 ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 1367 ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr); 1368 ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 1369 } 1370 ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr); 1371 PetscFunctionReturn(0); 1372 } 1373 1374 /*--------------------------------------------------------------------------------------------------------------------*/ 1375 1376 #undef __FUNCT__ 1377 #define __FUNCT__ "DMRegister" 1378 /*@C 1379 DMRegister - See DMRegisterDynamic() 1380 1381 Level: advanced 1382 @*/ 1383 PetscErrorCode DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM)) 1384 { 1385 char fullname[PETSC_MAX_PATH_LEN]; 1386 PetscErrorCode ierr; 1387 1388 PetscFunctionBegin; 1389 ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr); 1390 ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr); 1391 ierr = PetscStrcat(fullname, name);CHKERRQ(ierr); 1392 ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr); 1393 PetscFunctionReturn(0); 1394 } 1395 1396 1397 /*--------------------------------------------------------------------------------------------------------------------*/ 1398 #undef __FUNCT__ 1399 #define __FUNCT__ "DMRegisterDestroy" 1400 /*@C 1401 DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic(). 1402 1403 Not Collective 1404 1405 Level: advanced 1406 1407 .keywords: DM, register, destroy 1408 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic() 1409 @*/ 1410 PetscErrorCode DMRegisterDestroy(void) 1411 { 1412 PetscErrorCode ierr; 1413 1414 PetscFunctionBegin; 1415 ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr); 1416 DMRegisterAllCalled = PETSC_FALSE; 1417 PetscFunctionReturn(0); 1418 } 1419 1420 #if defined(PETSC_HAVE_MATLAB_ENGINE) 1421 #include <mex.h> 1422 1423 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext; 1424 1425 #undef __FUNCT__ 1426 #define __FUNCT__ "DMComputeFunction_Matlab" 1427 /* 1428 DMComputeFunction_Matlab - Calls the function that has been set with 1429 DMSetFunctionMatlab(). 1430 1431 For linear problems x is null 1432 1433 .seealso: DMSetFunction(), DMGetFunction() 1434 */ 1435 PetscErrorCode DMComputeFunction_Matlab(DM dm,Vec x,Vec y) 1436 { 1437 PetscErrorCode ierr; 1438 DMMatlabContext *sctx; 1439 int nlhs = 1,nrhs = 4; 1440 mxArray *plhs[1],*prhs[4]; 1441 long long int lx = 0,ly = 0,ls = 0; 1442 1443 PetscFunctionBegin; 1444 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1445 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 1446 PetscCheckSameComm(dm,1,y,3); 1447 1448 /* call Matlab function in ctx with arguments x and y */ 1449 ierr = DMGetContext(dm,(void**)&sctx);CHKERRQ(ierr); 1450 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 1451 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 1452 ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr); 1453 prhs[0] = mxCreateDoubleScalar((double)ls); 1454 prhs[1] = mxCreateDoubleScalar((double)lx); 1455 prhs[2] = mxCreateDoubleScalar((double)ly); 1456 prhs[3] = mxCreateString(sctx->funcname); 1457 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr); 1458 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 1459 mxDestroyArray(prhs[0]); 1460 mxDestroyArray(prhs[1]); 1461 mxDestroyArray(prhs[2]); 1462 mxDestroyArray(prhs[3]); 1463 mxDestroyArray(plhs[0]); 1464 PetscFunctionReturn(0); 1465 } 1466 1467 1468 #undef __FUNCT__ 1469 #define __FUNCT__ "DMSetFunctionMatlab" 1470 /* 1471 DMSetFunctionMatlab - Sets the function evaluation routine 1472 1473 */ 1474 PetscErrorCode DMSetFunctionMatlab(DM dm,const char *func) 1475 { 1476 PetscErrorCode ierr; 1477 DMMatlabContext *sctx; 1478 1479 PetscFunctionBegin; 1480 /* currently sctx is memory bleed */ 1481 ierr = DMGetContext(dm,(void**)&sctx);CHKERRQ(ierr); 1482 if (!sctx) { 1483 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 1484 } 1485 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 1486 ierr = DMSetContext(dm,sctx);CHKERRQ(ierr); 1487 ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr); 1488 PetscFunctionReturn(0); 1489 } 1490 1491 #undef __FUNCT__ 1492 #define __FUNCT__ "DMComputeJacobian_Matlab" 1493 /* 1494 DMComputeJacobian_Matlab - Calls the function that has been set with 1495 DMSetJacobianMatlab(). 1496 1497 For linear problems x is null 1498 1499 .seealso: DMSetFunction(), DMGetFunction() 1500 */ 1501 PetscErrorCode DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str) 1502 { 1503 PetscErrorCode ierr; 1504 DMMatlabContext *sctx; 1505 int nlhs = 2,nrhs = 5; 1506 mxArray *plhs[2],*prhs[5]; 1507 long long int lx = 0,lA = 0,lB = 0,ls = 0; 1508 1509 PetscFunctionBegin; 1510 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1511 PetscValidHeaderSpecific(A,MAT_CLASSID,3); 1512 1513 /* call MATLAB function in ctx with arguments x, A, and B */ 1514 ierr = DMGetContext(dm,(void**)&sctx);CHKERRQ(ierr); 1515 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 1516 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 1517 ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr); 1518 ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr); 1519 prhs[0] = mxCreateDoubleScalar((double)ls); 1520 prhs[1] = mxCreateDoubleScalar((double)lx); 1521 prhs[2] = mxCreateDoubleScalar((double)lA); 1522 prhs[3] = mxCreateDoubleScalar((double)lB); 1523 prhs[4] = mxCreateString(sctx->jacname); 1524 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr); 1525 *str = (MatStructure) mxGetScalar(plhs[0]); 1526 ierr = (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr); 1527 mxDestroyArray(prhs[0]); 1528 mxDestroyArray(prhs[1]); 1529 mxDestroyArray(prhs[2]); 1530 mxDestroyArray(prhs[3]); 1531 mxDestroyArray(prhs[4]); 1532 mxDestroyArray(plhs[0]); 1533 mxDestroyArray(plhs[1]); 1534 PetscFunctionReturn(0); 1535 } 1536 1537 1538 #undef __FUNCT__ 1539 #define __FUNCT__ "DMSetJacobianMatlab" 1540 /* 1541 DMSetJacobianMatlab - Sets the Jacobian function evaluation routine 1542 1543 */ 1544 PetscErrorCode DMSetJacobianMatlab(DM dm,const char *func) 1545 { 1546 PetscErrorCode ierr; 1547 DMMatlabContext *sctx; 1548 1549 PetscFunctionBegin; 1550 /* currently sctx is memory bleed */ 1551 ierr = DMGetContext(dm,(void**)&sctx);CHKERRQ(ierr); 1552 if (!sctx) { 1553 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 1554 } 1555 ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr); 1556 ierr = DMSetContext(dm,sctx);CHKERRQ(ierr); 1557 ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr); 1558 PetscFunctionReturn(0); 1559 } 1560 #endif 1561