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