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