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