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