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