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