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