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