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