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 Notes: For DMDA objects this only works for "uniform refinement", that is the refined mesh was obtained DMRefine() or the coarse mesh was obtained by 441 DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation. 442 443 For DMDA objects you can use this interpolation (more precisely the interpolation from the DMDAGetCoordinateDA()) to interpolate the mesh coordinate vectors 444 EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic. 445 446 447 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix(), DMRefine(), DMCoarsen() 448 449 @*/ 450 PetscErrorCode DMGetInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec) 451 { 452 PetscErrorCode ierr; 453 454 PetscFunctionBegin; 455 ierr = (*dm1->ops->getinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr); 456 PetscFunctionReturn(0); 457 } 458 459 #undef __FUNCT__ 460 #define __FUNCT__ "DMGetInjection" 461 /*@ 462 DMGetInjection - Gets injection matrix between two DMDA or DMComposite objects 463 464 Collective on DM 465 466 Input Parameter: 467 + dm1 - the DM object 468 - dm2 - the second, finer DM object 469 470 Output Parameter: 471 . ctx - the injection 472 473 Level: developer 474 475 Notes: For DMDA objects this only works for "uniform refinement", that is the refined mesh was obtained DMRefine() or the coarse mesh was obtained by 476 DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection. 477 478 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetColoring(), DMGetMatrix(), DMGetInterpolation() 479 480 @*/ 481 PetscErrorCode DMGetInjection(DM dm1,DM dm2,VecScatter *ctx) 482 { 483 PetscErrorCode ierr; 484 485 PetscFunctionBegin; 486 ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr); 487 PetscFunctionReturn(0); 488 } 489 490 #undef __FUNCT__ 491 #define __FUNCT__ "DMGetColoring" 492 /*@C 493 DMGetColoring - Gets coloring for a DMDA or DMComposite 494 495 Collective on DM 496 497 Input Parameter: 498 + dm - the DM object 499 . ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL 500 - matype - either MATAIJ or MATBAIJ 501 502 Output Parameter: 503 . coloring - the coloring 504 505 Level: developer 506 507 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetMatrix() 508 509 @*/ 510 PetscErrorCode DMGetColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring) 511 { 512 PetscErrorCode ierr; 513 514 PetscFunctionBegin; 515 if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet"); 516 ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr); 517 PetscFunctionReturn(0); 518 } 519 520 #undef __FUNCT__ 521 #define __FUNCT__ "DMGetMatrix" 522 /*@C 523 DMGetMatrix - Gets empty Jacobian for a DMDA or DMComposite 524 525 Collective on DM 526 527 Input Parameter: 528 + dm - the DM object 529 - mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or 530 any type which inherits from one of these (such as MATAIJ) 531 532 Output Parameter: 533 . mat - the empty Jacobian 534 535 Level: developer 536 537 Notes: This properly preallocates the number of nonzeros in the sparse matrix so you 538 do not need to do it yourself. 539 540 By default it also sets the nonzero structure and puts in the zero entries. To prevent setting 541 the nonzero pattern call DMDASetMatPreallocateOnly() 542 543 For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used 544 internally by PETSc. 545 546 For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires 547 the indices for the global numbering for DMDAs which is complicated. 548 549 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 550 551 @*/ 552 PetscErrorCode DMGetMatrix(DM dm,const MatType mtype,Mat *mat) 553 { 554 PetscErrorCode ierr; 555 char ttype[256]; 556 PetscBool flg; 557 558 PetscFunctionBegin; 559 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 560 ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr); 561 #endif 562 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 563 PetscValidPointer(mat,3); 564 ierr = PetscStrncpy(ttype,mtype,sizeof(ttype));CHKERRQ(ierr); 565 ttype[sizeof(ttype)-1] = 0; 566 ierr = PetscOptionsBegin(((PetscObject)dm)->comm,((PetscObject)dm)->prefix,"DM options","Mat");CHKERRQ(ierr); 567 ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,ttype,ttype,sizeof(ttype),&flg);CHKERRQ(ierr); 568 ierr = PetscOptionsEnd(); 569 if (flg || mtype) { 570 ierr = (*dm->ops->getmatrix)(dm,ttype,mat);CHKERRQ(ierr); 571 } else { /* Let the implementation decide */ 572 ierr = (*dm->ops->getmatrix)(dm,PETSC_NULL,mat);CHKERRQ(ierr); 573 } 574 PetscFunctionReturn(0); 575 } 576 577 #undef __FUNCT__ 578 #define __FUNCT__ "DMSetMatrixPreallocateOnly" 579 /*@ 580 DMSetMatrixPreallocateOnly - When DMGetMatrix() is called the matrix will be properly 581 preallocated but the nonzero structure and zero values will not be set. 582 583 Logically Collective on DMDA 584 585 Input Parameter: 586 + dm - the DM 587 - only - PETSC_TRUE if only want preallocation 588 589 Level: developer 590 .seealso DMGetMatrix() 591 @*/ 592 PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only) 593 { 594 PetscFunctionBegin; 595 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 596 dm->prealloc_only = only; 597 PetscFunctionReturn(0); 598 } 599 600 #undef __FUNCT__ 601 #define __FUNCT__ "DMRefine" 602 /*@ 603 DMRefine - Refines a DM object 604 605 Collective on DM 606 607 Input Parameter: 608 + dm - the DM object 609 - comm - the communicator to contain the new DM object (or PETSC_NULL) 610 611 Output Parameter: 612 . dmf - the refined DM 613 614 Level: developer 615 616 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 617 618 @*/ 619 PetscErrorCode DMRefine(DM dm,MPI_Comm comm,DM *dmf) 620 { 621 PetscErrorCode ierr; 622 623 PetscFunctionBegin; 624 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 625 ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr); 626 (*dmf)->ops->initialguess = dm->ops->initialguess; 627 (*dmf)->ops->function = dm->ops->function; 628 (*dmf)->ops->functionj = dm->ops->functionj; 629 if (dm->ops->jacobian != DMComputeJacobianDefault) { 630 (*dmf)->ops->jacobian = dm->ops->jacobian; 631 } 632 (*dmf)->ctx = dm->ctx; 633 (*dmf)->levelup = dm->levelup + 1; 634 PetscFunctionReturn(0); 635 } 636 637 #undef __FUNCT__ 638 #define __FUNCT__ "DMGetRefineLevel" 639 /*@ 640 DMGetRefineLevel - Get's the number of refinements that have generated this DM. 641 642 Not Collective 643 644 Input Parameter: 645 . dm - the DM object 646 647 Output Parameter: 648 . level - number of refinements 649 650 Level: developer 651 652 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 653 654 @*/ 655 PetscErrorCode DMGetRefineLevel(DM dm,PetscInt *level) 656 { 657 PetscFunctionBegin; 658 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 659 *level = dm->levelup; 660 PetscFunctionReturn(0); 661 } 662 663 #undef __FUNCT__ 664 #define __FUNCT__ "DMGlobalToLocalBegin" 665 /*@ 666 DMGlobalToLocalBegin - Begins updating local vectors from global vector 667 668 Neighbor-wise Collective on DM 669 670 Input Parameters: 671 + dm - the DM object 672 . g - the global vector 673 . mode - INSERT_VALUES or ADD_VALUES 674 - l - the local vector 675 676 677 Level: beginner 678 679 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 680 681 @*/ 682 PetscErrorCode DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) 683 { 684 PetscErrorCode ierr; 685 686 PetscFunctionBegin; 687 ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode,l);CHKERRQ(ierr); 688 PetscFunctionReturn(0); 689 } 690 691 #undef __FUNCT__ 692 #define __FUNCT__ "DMGlobalToLocalEnd" 693 /*@ 694 DMGlobalToLocalEnd - Ends updating local vectors from global vector 695 696 Neighbor-wise Collective on DM 697 698 Input Parameters: 699 + dm - the DM object 700 . g - the global vector 701 . mode - INSERT_VALUES or ADD_VALUES 702 - l - the local vector 703 704 705 Level: beginner 706 707 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 708 709 @*/ 710 PetscErrorCode DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) 711 { 712 PetscErrorCode ierr; 713 714 PetscFunctionBegin; 715 ierr = (*dm->ops->globaltolocalend)(dm,g,mode,l);CHKERRQ(ierr); 716 PetscFunctionReturn(0); 717 } 718 719 #undef __FUNCT__ 720 #define __FUNCT__ "DMLocalToGlobalBegin" 721 /*@ 722 DMLocalToGlobalBegin - updates global vectors from local vectors 723 724 Neighbor-wise Collective on DM 725 726 Input Parameters: 727 + dm - the DM object 728 . l - the local vector 729 . 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 730 base point. 731 - - the global vector 732 733 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 734 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 735 global array to the final global array with VecAXPY(). 736 737 Level: beginner 738 739 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin() 740 741 @*/ 742 PetscErrorCode DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g) 743 { 744 PetscErrorCode ierr; 745 746 PetscFunctionBegin; 747 ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode,g);CHKERRQ(ierr); 748 PetscFunctionReturn(0); 749 } 750 751 #undef __FUNCT__ 752 #define __FUNCT__ "DMLocalToGlobalEnd" 753 /*@ 754 DMLocalToGlobalEnd - updates global vectors from local vectors 755 756 Neighbor-wise Collective on DM 757 758 Input Parameters: 759 + dm - the DM object 760 . l - the local vector 761 . mode - INSERT_VALUES or ADD_VALUES 762 - g - the global vector 763 764 765 Level: beginner 766 767 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd() 768 769 @*/ 770 PetscErrorCode DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g) 771 { 772 PetscErrorCode ierr; 773 774 PetscFunctionBegin; 775 ierr = (*dm->ops->localtoglobalend)(dm,l,mode,g);CHKERRQ(ierr); 776 PetscFunctionReturn(0); 777 } 778 779 #undef __FUNCT__ 780 #define __FUNCT__ "DMComputeJacobianDefault" 781 /*@ 782 DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided 783 784 Collective on DM 785 786 Input Parameter: 787 + dm - the DM object 788 . x - location to compute Jacobian at; may be ignored for linear problems 789 . A - matrix that defines the operator for the linear solve 790 - B - the matrix used to construct the preconditioner 791 792 Level: developer 793 794 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 795 DMSetFunction() 796 797 @*/ 798 PetscErrorCode DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 799 { 800 PetscErrorCode ierr; 801 PetscFunctionBegin; 802 *stflag = SAME_NONZERO_PATTERN; 803 ierr = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr); 804 if (A != B) { 805 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 806 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 807 } 808 PetscFunctionReturn(0); 809 } 810 811 #undef __FUNCT__ 812 #define __FUNCT__ "DMCoarsen" 813 /*@ 814 DMCoarsen - Coarsens a DM object 815 816 Collective on DM 817 818 Input Parameter: 819 + dm - the DM object 820 - comm - the communicator to contain the new DM object (or PETSC_NULL) 821 822 Output Parameter: 823 . dmc - the coarsened DM 824 825 Level: developer 826 827 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 828 829 @*/ 830 PetscErrorCode DMCoarsen(DM dm, MPI_Comm comm, DM *dmc) 831 { 832 PetscErrorCode ierr; 833 834 PetscFunctionBegin; 835 ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr); 836 (*dmc)->ops->initialguess = dm->ops->initialguess; 837 (*dmc)->ops->function = dm->ops->function; 838 (*dmc)->ops->functionj = dm->ops->functionj; 839 if (dm->ops->jacobian != DMComputeJacobianDefault) { 840 (*dmc)->ops->jacobian = dm->ops->jacobian; 841 } 842 (*dmc)->ctx = dm->ctx; 843 (*dmc)->leveldown = dm->leveldown + 1; 844 PetscFunctionReturn(0); 845 } 846 847 #undef __FUNCT__ 848 #define __FUNCT__ "DMRefineHierarchy" 849 /*@C 850 DMRefineHierarchy - Refines a DM object, all levels at once 851 852 Collective on DM 853 854 Input Parameter: 855 + dm - the DM object 856 - nlevels - the number of levels of refinement 857 858 Output Parameter: 859 . dmf - the refined DM hierarchy 860 861 Level: developer 862 863 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 864 865 @*/ 866 PetscErrorCode DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[]) 867 { 868 PetscErrorCode ierr; 869 870 PetscFunctionBegin; 871 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 872 if (nlevels == 0) PetscFunctionReturn(0); 873 if (dm->ops->refinehierarchy) { 874 ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr); 875 } else if (dm->ops->refine) { 876 PetscInt i; 877 878 ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr); 879 for (i=1; i<nlevels; i++) { 880 ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr); 881 } 882 } else { 883 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet"); 884 } 885 PetscFunctionReturn(0); 886 } 887 888 #undef __FUNCT__ 889 #define __FUNCT__ "DMCoarsenHierarchy" 890 /*@C 891 DMCoarsenHierarchy - Coarsens a DM object, all levels at once 892 893 Collective on DM 894 895 Input Parameter: 896 + dm - the DM object 897 - nlevels - the number of levels of coarsening 898 899 Output Parameter: 900 . dmc - the coarsened DM hierarchy 901 902 Level: developer 903 904 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMGetInterpolation() 905 906 @*/ 907 PetscErrorCode DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[]) 908 { 909 PetscErrorCode ierr; 910 911 PetscFunctionBegin; 912 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 913 if (nlevels == 0) PetscFunctionReturn(0); 914 PetscValidPointer(dmc,3); 915 if (dm->ops->coarsenhierarchy) { 916 ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr); 917 } else if (dm->ops->coarsen) { 918 PetscInt i; 919 920 ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr); 921 for (i=1; i<nlevels; i++) { 922 ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr); 923 } 924 } else { 925 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet"); 926 } 927 PetscFunctionReturn(0); 928 } 929 930 #undef __FUNCT__ 931 #define __FUNCT__ "DMGetAggregates" 932 /*@ 933 DMGetAggregates - Gets the aggregates that map between 934 grids associated with two DMs. 935 936 Collective on DM 937 938 Input Parameters: 939 + dmc - the coarse grid DM 940 - dmf - the fine grid DM 941 942 Output Parameters: 943 . rest - the restriction matrix (transpose of the projection matrix) 944 945 Level: intermediate 946 947 .keywords: interpolation, restriction, multigrid 948 949 .seealso: DMRefine(), DMGetInjection(), DMGetInterpolation() 950 @*/ 951 PetscErrorCode DMGetAggregates(DM dmc, DM dmf, Mat *rest) 952 { 953 PetscErrorCode ierr; 954 955 PetscFunctionBegin; 956 ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr); 957 PetscFunctionReturn(0); 958 } 959 960 #undef __FUNCT__ 961 #define __FUNCT__ "DMSetApplicationContext" 962 /*@ 963 DMSetApplicationContext - Set a user context into a DM object 964 965 Not Collective 966 967 Input Parameters: 968 + dm - the DM object 969 - ctx - the user context 970 971 Level: intermediate 972 973 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext() 974 975 @*/ 976 PetscErrorCode DMSetApplicationContext(DM dm,void *ctx) 977 { 978 PetscFunctionBegin; 979 dm->ctx = ctx; 980 PetscFunctionReturn(0); 981 } 982 983 #undef __FUNCT__ 984 #define __FUNCT__ "DMGetApplicationContext" 985 /*@ 986 DMGetApplicationContext - Gets a user context from a DM object 987 988 Not Collective 989 990 Input Parameter: 991 . dm - the DM object 992 993 Output Parameter: 994 . ctx - the user context 995 996 Level: intermediate 997 998 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext() 999 1000 @*/ 1001 PetscErrorCode DMGetApplicationContext(DM dm,void *ctx) 1002 { 1003 PetscFunctionBegin; 1004 *(void**)ctx = dm->ctx; 1005 PetscFunctionReturn(0); 1006 } 1007 1008 #undef __FUNCT__ 1009 #define __FUNCT__ "DMSetInitialGuess" 1010 /*@ 1011 DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers 1012 1013 Logically Collective on DM 1014 1015 Input Parameter: 1016 + dm - the DM object to destroy 1017 - f - the function to compute the initial guess 1018 1019 Level: intermediate 1020 1021 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1022 1023 @*/ 1024 PetscErrorCode DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec)) 1025 { 1026 PetscFunctionBegin; 1027 dm->ops->initialguess = f; 1028 PetscFunctionReturn(0); 1029 } 1030 1031 #undef __FUNCT__ 1032 #define __FUNCT__ "DMSetFunction" 1033 /*@ 1034 DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES 1035 1036 Logically Collective on DM 1037 1038 Input Parameter: 1039 + dm - the DM object 1040 - f - the function to compute (use PETSC_NULL to cancel a previous function that was set) 1041 1042 Level: intermediate 1043 1044 Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian 1045 computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian. 1046 1047 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1048 DMSetJacobian() 1049 1050 @*/ 1051 PetscErrorCode DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 1052 { 1053 PetscFunctionBegin; 1054 dm->ops->function = f; 1055 if (f) { 1056 dm->ops->functionj = f; 1057 } 1058 PetscFunctionReturn(0); 1059 } 1060 1061 #undef __FUNCT__ 1062 #define __FUNCT__ "DMSetJacobian" 1063 /*@ 1064 DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES 1065 1066 Logically Collective on DM 1067 1068 Input Parameter: 1069 + dm - the DM object to destroy 1070 - f - the function to compute the matrix entries 1071 1072 Level: intermediate 1073 1074 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1075 DMSetFunction() 1076 1077 @*/ 1078 PetscErrorCode DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*)) 1079 { 1080 PetscFunctionBegin; 1081 dm->ops->jacobian = f; 1082 PetscFunctionReturn(0); 1083 } 1084 1085 #undef __FUNCT__ 1086 #define __FUNCT__ "DMComputeInitialGuess" 1087 /*@ 1088 DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers 1089 1090 Collective on DM 1091 1092 Input Parameter: 1093 + dm - the DM object to destroy 1094 - x - the vector to hold the initial guess values 1095 1096 Level: developer 1097 1098 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat() 1099 1100 @*/ 1101 PetscErrorCode DMComputeInitialGuess(DM dm,Vec x) 1102 { 1103 PetscErrorCode ierr; 1104 PetscFunctionBegin; 1105 if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()"); 1106 ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr); 1107 PetscFunctionReturn(0); 1108 } 1109 1110 #undef __FUNCT__ 1111 #define __FUNCT__ "DMHasInitialGuess" 1112 /*@ 1113 DMHasInitialGuess - does the DM object have an initial guess function 1114 1115 Not Collective 1116 1117 Input Parameter: 1118 . dm - the DM object to destroy 1119 1120 Output Parameter: 1121 . flg - PETSC_TRUE if function exists 1122 1123 Level: developer 1124 1125 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1126 1127 @*/ 1128 PetscErrorCode DMHasInitialGuess(DM dm,PetscBool *flg) 1129 { 1130 PetscFunctionBegin; 1131 *flg = (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE; 1132 PetscFunctionReturn(0); 1133 } 1134 1135 #undef __FUNCT__ 1136 #define __FUNCT__ "DMHasFunction" 1137 /*@ 1138 DMHasFunction - does the DM object have a function 1139 1140 Not Collective 1141 1142 Input Parameter: 1143 . dm - the DM object to destroy 1144 1145 Output Parameter: 1146 . flg - PETSC_TRUE if function exists 1147 1148 Level: developer 1149 1150 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1151 1152 @*/ 1153 PetscErrorCode DMHasFunction(DM dm,PetscBool *flg) 1154 { 1155 PetscFunctionBegin; 1156 *flg = (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE; 1157 PetscFunctionReturn(0); 1158 } 1159 1160 #undef __FUNCT__ 1161 #define __FUNCT__ "DMHasJacobian" 1162 /*@ 1163 DMHasJacobian - does the DM object have a matrix function 1164 1165 Not Collective 1166 1167 Input Parameter: 1168 . dm - the DM object to destroy 1169 1170 Output Parameter: 1171 . flg - PETSC_TRUE if function exists 1172 1173 Level: developer 1174 1175 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1176 1177 @*/ 1178 PetscErrorCode DMHasJacobian(DM dm,PetscBool *flg) 1179 { 1180 PetscFunctionBegin; 1181 *flg = (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE; 1182 PetscFunctionReturn(0); 1183 } 1184 1185 #undef __FUNCT__ 1186 #define __FUNCT__ "DMComputeFunction" 1187 /*@ 1188 DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES 1189 1190 Collective on DM 1191 1192 Input Parameter: 1193 + dm - the DM object to destroy 1194 . x - the location where the function is evaluationed, may be ignored for linear problems 1195 - b - the vector to hold the right hand side entries 1196 1197 Level: developer 1198 1199 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1200 DMSetJacobian() 1201 1202 @*/ 1203 PetscErrorCode DMComputeFunction(DM dm,Vec x,Vec b) 1204 { 1205 PetscErrorCode ierr; 1206 PetscFunctionBegin; 1207 if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()"); 1208 PetscStackPush("DM user function"); 1209 ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr); 1210 PetscStackPop; 1211 PetscFunctionReturn(0); 1212 } 1213 1214 1215 #undef __FUNCT__ 1216 #define __FUNCT__ "DMComputeJacobian" 1217 /*@ 1218 DMComputeJacobian - compute the matrix entries for the solver 1219 1220 Collective on DM 1221 1222 Input Parameter: 1223 + dm - the DM object 1224 . x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM 1225 . A - matrix that defines the operator for the linear solve 1226 - B - the matrix used to construct the preconditioner 1227 1228 Level: developer 1229 1230 .seealso DMView(), DMCreateGlobalVector(), DMGetInterpolation(), DMGetColoring(), DMGetMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1231 DMSetFunction() 1232 1233 @*/ 1234 PetscErrorCode DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 1235 { 1236 PetscErrorCode ierr; 1237 1238 PetscFunctionBegin; 1239 if (!dm->ops->jacobian) { 1240 ISColoring coloring; 1241 MatFDColoring fd; 1242 1243 ierr = DMGetColoring(dm,IS_COLORING_GLOBAL,MATAIJ,&coloring);CHKERRQ(ierr); 1244 ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr); 1245 ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 1246 ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr); 1247 1248 dm->fd = fd; 1249 dm->ops->jacobian = DMComputeJacobianDefault; 1250 1251 /* don't know why this is needed */ 1252 ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 1253 } 1254 if (!x) x = dm->x; 1255 ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr); 1256 1257 /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */ 1258 if (x) { 1259 if (!dm->x) { 1260 ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 1261 } 1262 ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 1263 } 1264 if (A != B) { 1265 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1266 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1267 } 1268 PetscFunctionReturn(0); 1269 } 1270 1271 1272 PetscFList DMList = PETSC_NULL; 1273 PetscBool DMRegisterAllCalled = PETSC_FALSE; 1274 1275 #undef __FUNCT__ 1276 #define __FUNCT__ "DMSetType" 1277 /*@C 1278 DMSetType - Builds a DM, for a particular DM implementation. 1279 1280 Collective on DM 1281 1282 Input Parameters: 1283 + dm - The DM object 1284 - method - The name of the DM type 1285 1286 Options Database Key: 1287 . -dm_type <type> - Sets the DM type; use -help for a list of available types 1288 1289 Notes: 1290 See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D). 1291 1292 Level: intermediate 1293 1294 .keywords: DM, set, type 1295 .seealso: DMGetType(), DMCreate() 1296 @*/ 1297 PetscErrorCode DMSetType(DM dm, const DMType method) 1298 { 1299 PetscErrorCode (*r)(DM); 1300 PetscBool match; 1301 PetscErrorCode ierr; 1302 1303 PetscFunctionBegin; 1304 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1305 ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr); 1306 if (match) PetscFunctionReturn(0); 1307 1308 if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1309 ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 1310 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method); 1311 1312 if (dm->ops->destroy) { 1313 ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr); 1314 } 1315 ierr = (*r)(dm);CHKERRQ(ierr); 1316 ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr); 1317 PetscFunctionReturn(0); 1318 } 1319 1320 #undef __FUNCT__ 1321 #define __FUNCT__ "DMGetType" 1322 /*@C 1323 DMGetType - Gets the DM type name (as a string) from the DM. 1324 1325 Not Collective 1326 1327 Input Parameter: 1328 . dm - The DM 1329 1330 Output Parameter: 1331 . type - The DM type name 1332 1333 Level: intermediate 1334 1335 .keywords: DM, get, type, name 1336 .seealso: DMSetType(), DMCreate() 1337 @*/ 1338 PetscErrorCode DMGetType(DM dm, const DMType *type) 1339 { 1340 PetscErrorCode ierr; 1341 1342 PetscFunctionBegin; 1343 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1344 PetscValidCharPointer(type,2); 1345 if (!DMRegisterAllCalled) { 1346 ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1347 } 1348 *type = ((PetscObject)dm)->type_name; 1349 PetscFunctionReturn(0); 1350 } 1351 1352 #undef __FUNCT__ 1353 #define __FUNCT__ "DMConvert" 1354 /*@C 1355 DMConvert - Converts a DM to another DM, either of the same or different type. 1356 1357 Collective on DM 1358 1359 Input Parameters: 1360 + dm - the DM 1361 - newtype - new DM type (use "same" for the same type) 1362 1363 Output Parameter: 1364 . M - pointer to new DM 1365 1366 Notes: 1367 Cannot be used to convert a sequential DM to parallel or parallel to sequential, 1368 the MPI communicator of the generated DM is always the same as the communicator 1369 of the input DM. 1370 1371 Level: intermediate 1372 1373 .seealso: DMCreate() 1374 @*/ 1375 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M) 1376 { 1377 DM B; 1378 char convname[256]; 1379 PetscBool sametype, issame; 1380 PetscErrorCode ierr; 1381 1382 PetscFunctionBegin; 1383 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1384 PetscValidType(dm,1); 1385 PetscValidPointer(M,3); 1386 ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr); 1387 ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr); 1388 { 1389 PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL; 1390 1391 /* 1392 Order of precedence: 1393 1) See if a specialized converter is known to the current DM. 1394 2) See if a specialized converter is known to the desired DM class. 1395 3) See if a good general converter is registered for the desired class 1396 4) See if a good general converter is known for the current matrix. 1397 5) Use a really basic converter. 1398 */ 1399 1400 /* 1) See if a specialized converter is known to the current DM and the desired class */ 1401 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 1402 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 1403 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 1404 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 1405 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 1406 ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr); 1407 if (conv) goto foundconv; 1408 1409 /* 2) See if a specialized converter is known to the desired DM class. */ 1410 ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr); 1411 ierr = DMSetType(B, newtype);CHKERRQ(ierr); 1412 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 1413 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 1414 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 1415 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 1416 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 1417 ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr); 1418 if (conv) { 1419 ierr = DMDestroy(&B);CHKERRQ(ierr); 1420 goto foundconv; 1421 } 1422 1423 #if 0 1424 /* 3) See if a good general converter is registered for the desired class */ 1425 conv = B->ops->convertfrom; 1426 ierr = DMDestroy(&B);CHKERRQ(ierr); 1427 if (conv) goto foundconv; 1428 1429 /* 4) See if a good general converter is known for the current matrix */ 1430 if (dm->ops->convert) { 1431 conv = dm->ops->convert; 1432 } 1433 if (conv) goto foundconv; 1434 #endif 1435 1436 /* 5) Use a really basic converter. */ 1437 SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype); 1438 1439 foundconv: 1440 ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 1441 ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr); 1442 ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 1443 } 1444 ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr); 1445 PetscFunctionReturn(0); 1446 } 1447 1448 /*--------------------------------------------------------------------------------------------------------------------*/ 1449 1450 #undef __FUNCT__ 1451 #define __FUNCT__ "DMRegister" 1452 /*@C 1453 DMRegister - See DMRegisterDynamic() 1454 1455 Level: advanced 1456 @*/ 1457 PetscErrorCode DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM)) 1458 { 1459 char fullname[PETSC_MAX_PATH_LEN]; 1460 PetscErrorCode ierr; 1461 1462 PetscFunctionBegin; 1463 ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr); 1464 ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr); 1465 ierr = PetscStrcat(fullname, name);CHKERRQ(ierr); 1466 ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr); 1467 PetscFunctionReturn(0); 1468 } 1469 1470 1471 /*--------------------------------------------------------------------------------------------------------------------*/ 1472 #undef __FUNCT__ 1473 #define __FUNCT__ "DMRegisterDestroy" 1474 /*@C 1475 DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic(). 1476 1477 Not Collective 1478 1479 Level: advanced 1480 1481 .keywords: DM, register, destroy 1482 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic() 1483 @*/ 1484 PetscErrorCode DMRegisterDestroy(void) 1485 { 1486 PetscErrorCode ierr; 1487 1488 PetscFunctionBegin; 1489 ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr); 1490 DMRegisterAllCalled = PETSC_FALSE; 1491 PetscFunctionReturn(0); 1492 } 1493 1494 #if defined(PETSC_HAVE_MATLAB_ENGINE) 1495 #include <mex.h> 1496 1497 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext; 1498 1499 #undef __FUNCT__ 1500 #define __FUNCT__ "DMComputeFunction_Matlab" 1501 /* 1502 DMComputeFunction_Matlab - Calls the function that has been set with 1503 DMSetFunctionMatlab(). 1504 1505 For linear problems x is null 1506 1507 .seealso: DMSetFunction(), DMGetFunction() 1508 */ 1509 PetscErrorCode DMComputeFunction_Matlab(DM dm,Vec x,Vec y) 1510 { 1511 PetscErrorCode ierr; 1512 DMMatlabContext *sctx; 1513 int nlhs = 1,nrhs = 4; 1514 mxArray *plhs[1],*prhs[4]; 1515 long long int lx = 0,ly = 0,ls = 0; 1516 1517 PetscFunctionBegin; 1518 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1519 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 1520 PetscCheckSameComm(dm,1,y,3); 1521 1522 /* call Matlab function in ctx with arguments x and y */ 1523 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1524 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 1525 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 1526 ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr); 1527 prhs[0] = mxCreateDoubleScalar((double)ls); 1528 prhs[1] = mxCreateDoubleScalar((double)lx); 1529 prhs[2] = mxCreateDoubleScalar((double)ly); 1530 prhs[3] = mxCreateString(sctx->funcname); 1531 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr); 1532 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 1533 mxDestroyArray(prhs[0]); 1534 mxDestroyArray(prhs[1]); 1535 mxDestroyArray(prhs[2]); 1536 mxDestroyArray(prhs[3]); 1537 mxDestroyArray(plhs[0]); 1538 PetscFunctionReturn(0); 1539 } 1540 1541 1542 #undef __FUNCT__ 1543 #define __FUNCT__ "DMSetFunctionMatlab" 1544 /* 1545 DMSetFunctionMatlab - Sets the function evaluation routine 1546 1547 */ 1548 PetscErrorCode DMSetFunctionMatlab(DM dm,const char *func) 1549 { 1550 PetscErrorCode ierr; 1551 DMMatlabContext *sctx; 1552 1553 PetscFunctionBegin; 1554 /* currently sctx is memory bleed */ 1555 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1556 if (!sctx) { 1557 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 1558 } 1559 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 1560 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 1561 ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr); 1562 PetscFunctionReturn(0); 1563 } 1564 1565 #undef __FUNCT__ 1566 #define __FUNCT__ "DMComputeJacobian_Matlab" 1567 /* 1568 DMComputeJacobian_Matlab - Calls the function that has been set with 1569 DMSetJacobianMatlab(). 1570 1571 For linear problems x is null 1572 1573 .seealso: DMSetFunction(), DMGetFunction() 1574 */ 1575 PetscErrorCode DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str) 1576 { 1577 PetscErrorCode ierr; 1578 DMMatlabContext *sctx; 1579 int nlhs = 2,nrhs = 5; 1580 mxArray *plhs[2],*prhs[5]; 1581 long long int lx = 0,lA = 0,lB = 0,ls = 0; 1582 1583 PetscFunctionBegin; 1584 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1585 PetscValidHeaderSpecific(A,MAT_CLASSID,3); 1586 1587 /* call MATLAB function in ctx with arguments x, A, and B */ 1588 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1589 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 1590 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 1591 ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr); 1592 ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr); 1593 prhs[0] = mxCreateDoubleScalar((double)ls); 1594 prhs[1] = mxCreateDoubleScalar((double)lx); 1595 prhs[2] = mxCreateDoubleScalar((double)lA); 1596 prhs[3] = mxCreateDoubleScalar((double)lB); 1597 prhs[4] = mxCreateString(sctx->jacname); 1598 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr); 1599 *str = (MatStructure) mxGetScalar(plhs[0]); 1600 ierr = (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr); 1601 mxDestroyArray(prhs[0]); 1602 mxDestroyArray(prhs[1]); 1603 mxDestroyArray(prhs[2]); 1604 mxDestroyArray(prhs[3]); 1605 mxDestroyArray(prhs[4]); 1606 mxDestroyArray(plhs[0]); 1607 mxDestroyArray(plhs[1]); 1608 PetscFunctionReturn(0); 1609 } 1610 1611 1612 #undef __FUNCT__ 1613 #define __FUNCT__ "DMSetJacobianMatlab" 1614 /* 1615 DMSetJacobianMatlab - Sets the Jacobian function evaluation routine 1616 1617 */ 1618 PetscErrorCode DMSetJacobianMatlab(DM dm,const char *func) 1619 { 1620 PetscErrorCode ierr; 1621 DMMatlabContext *sctx; 1622 1623 PetscFunctionBegin; 1624 /* currently sctx is memory bleed */ 1625 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1626 if (!sctx) { 1627 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 1628 } 1629 ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr); 1630 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 1631 ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr); 1632 PetscFunctionReturn(0); 1633 } 1634 #endif 1635 1636 #undef __FUNCT__ 1637 #define __FUNCT__ "DMLoad" 1638 /*@C 1639 DMLoad - Loads a DM that has been stored in binary or HDF5 format 1640 with DMView(). 1641 1642 Collective on PetscViewer 1643 1644 Input Parameters: 1645 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or 1646 some related function before a call to DMLoad(). 1647 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 1648 HDF5 file viewer, obtained from PetscViewerHDF5Open() 1649 1650 Level: intermediate 1651 1652 Notes: 1653 Defaults to the DM DA. 1654 1655 Notes for advanced users: 1656 Most users should not need to know the details of the binary storage 1657 format, since DMLoad() and DMView() completely hide these details. 1658 But for anyone who's interested, the standard binary matrix storage 1659 format is 1660 .vb 1661 has not yet been determined 1662 .ve 1663 1664 In addition, PETSc automatically does the byte swapping for 1665 machines that store the bytes reversed, e.g. DEC alpha, freebsd, 1666 linux, Windows and the paragon; thus if you write your own binary 1667 read/write routines you have to swap the bytes; see PetscBinaryRead() 1668 and PetscBinaryWrite() to see how this may be done. 1669 1670 Concepts: vector^loading from file 1671 1672 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad() 1673 @*/ 1674 PetscErrorCode DMLoad(DM newdm, PetscViewer viewer) 1675 { 1676 PetscErrorCode ierr; 1677 1678 PetscFunctionBegin; 1679 PetscValidHeaderSpecific(newdm,DM_CLASSID,1); 1680 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 1681 1682 if (!((PetscObject)newdm)->type_name) { 1683 ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr); 1684 } 1685 ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr); 1686 PetscFunctionReturn(0); 1687 } 1688 1689