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