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