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