1 #include <petscsnes.h> 2 #include <petsc-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 45 v->workSize = 0; 46 v->workArray = PETSC_NULL; 47 v->ltogmap = PETSC_NULL; 48 v->ltogmapb = PETSC_NULL; 49 v->bs = 1; 50 v->coloringtype = IS_COLORING_GLOBAL; 51 52 *dm = v; 53 PetscFunctionReturn(0); 54 } 55 56 57 #undef __FUNCT__ 58 #define __FUNCT__ "DMSetVecType" 59 /*@C 60 DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() 61 62 Logically Collective on DMDA 63 64 Input Parameter: 65 + da - initial distributed array 66 . ctype - the vector type, currently either VECSTANDARD or VECCUSP 67 68 Options Database: 69 . -dm_vec_type ctype 70 71 Level: intermediate 72 73 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType 74 @*/ 75 PetscErrorCode DMSetVecType(DM da,const VecType ctype) 76 { 77 PetscErrorCode ierr; 78 79 PetscFunctionBegin; 80 PetscValidHeaderSpecific(da,DM_CLASSID,1); 81 ierr = PetscFree(da->vectype);CHKERRQ(ierr); 82 ierr = PetscStrallocpy(ctype,&da->vectype);CHKERRQ(ierr); 83 PetscFunctionReturn(0); 84 } 85 86 #undef __FUNCT__ 87 #define __FUNCT__ "DMSetOptionsPrefix" 88 /*@C 89 DMSetOptionsPrefix - Sets the prefix used for searching for all 90 DMDA options in the database. 91 92 Logically Collective on DMDA 93 94 Input Parameter: 95 + da - the DMDA context 96 - prefix - the prefix to prepend to all option names 97 98 Notes: 99 A hyphen (-) must NOT be given at the beginning of the prefix name. 100 The first character of all runtime options is AUTOMATICALLY the hyphen. 101 102 Level: advanced 103 104 .keywords: DMDA, set, options, prefix, database 105 106 .seealso: DMSetFromOptions() 107 @*/ 108 PetscErrorCode DMSetOptionsPrefix(DM dm,const char prefix[]) 109 { 110 PetscErrorCode ierr; 111 112 PetscFunctionBegin; 113 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 114 ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr); 115 PetscFunctionReturn(0); 116 } 117 118 #undef __FUNCT__ 119 #define __FUNCT__ "DMDestroy" 120 /*@ 121 DMDestroy - Destroys a vector packer or DMDA. 122 123 Collective on DM 124 125 Input Parameter: 126 . dm - the DM object to destroy 127 128 Level: developer 129 130 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 131 132 @*/ 133 PetscErrorCode DMDestroy(DM *dm) 134 { 135 PetscInt i, cnt = 0; 136 DMCoarsenHookLink link,next; 137 PetscErrorCode ierr; 138 139 PetscFunctionBegin; 140 if (!*dm) PetscFunctionReturn(0); 141 PetscValidHeaderSpecific((*dm),DM_CLASSID,1); 142 143 /* count all the circular references of DM and its contained Vecs */ 144 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 145 if ((*dm)->localin[i]) {cnt++;} 146 if ((*dm)->globalin[i]) {cnt++;} 147 } 148 if ((*dm)->x) { 149 PetscObject obj; 150 ierr = PetscObjectQuery((PetscObject)(*dm)->x,"DM",&obj);CHKERRQ(ierr); 151 if (obj == (PetscObject)*dm) cnt++; 152 } 153 154 if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);} 155 /* 156 Need this test because the dm references the vectors that 157 reference the dm, so destroying the dm calls destroy on the 158 vectors that cause another destroy on the dm 159 */ 160 if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0); 161 ((PetscObject) (*dm))->refct = 0; 162 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 163 if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()"); 164 ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr); 165 } 166 167 /* Destroy the list of hooks */ 168 for (link=(*dm)->coarsenhook; link; link=next) { 169 next = link->next; 170 ierr = PetscFree(link);CHKERRQ(ierr); 171 } 172 (*dm)->coarsenhook = PETSC_NULL; 173 174 if ((*dm)->ctx && (*dm)->ctxdestroy) { 175 ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr); 176 } 177 ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr); 178 ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr); 179 ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr); 180 ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr); 181 ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr); 182 ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr); 183 ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr); 184 ierr = PetscFree((*dm)->workArray);CHKERRQ(ierr); 185 /* if memory was published with AMS then destroy it */ 186 ierr = PetscObjectDepublish(*dm);CHKERRQ(ierr); 187 188 ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr); 189 ierr = PetscFree((*dm)->data);CHKERRQ(ierr); 190 ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr); 191 PetscFunctionReturn(0); 192 } 193 194 #undef __FUNCT__ 195 #define __FUNCT__ "DMSetUp" 196 /*@ 197 DMSetUp - sets up the data structures inside a DM object 198 199 Collective on DM 200 201 Input Parameter: 202 . dm - the DM object to setup 203 204 Level: developer 205 206 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 207 208 @*/ 209 PetscErrorCode DMSetUp(DM dm) 210 { 211 PetscErrorCode ierr; 212 213 PetscFunctionBegin; 214 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 215 if (dm->setupcalled) PetscFunctionReturn(0); 216 if (dm->ops->setup) { 217 ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr); 218 } 219 dm->setupcalled = PETSC_TRUE; 220 PetscFunctionReturn(0); 221 } 222 223 #undef __FUNCT__ 224 #define __FUNCT__ "DMSetFromOptions" 225 /*@ 226 DMSetFromOptions - sets parameters in a DM from the options database 227 228 Collective on DM 229 230 Input Parameter: 231 . dm - the DM object to set options for 232 233 Options Database: 234 + -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros 235 . -dm_vec_type <type> type of vector to create inside DM 236 . -dm_mat_type <type> type of matrix to create inside DM 237 - -dm_coloring_type <global or ghosted> 238 239 Level: developer 240 241 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 242 243 @*/ 244 PetscErrorCode DMSetFromOptions(DM dm) 245 { 246 PetscBool flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flg; 247 PetscErrorCode ierr; 248 char typeName[256] = MATAIJ; 249 250 PetscFunctionBegin; 251 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 252 ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr); 253 ierr = PetscOptionsBool("-dm_view", "Information on DM", "DMView", flg1, &flg1, PETSC_NULL);CHKERRQ(ierr); 254 ierr = PetscOptionsBool("-dm_view_detail", "Exhaustive mesh description", "DMView", flg2, &flg2, PETSC_NULL);CHKERRQ(ierr); 255 ierr = PetscOptionsBool("-dm_view_vtk", "Output mesh in VTK format", "DMView", flg3, &flg3, PETSC_NULL);CHKERRQ(ierr); 256 ierr = PetscOptionsBool("-dm_view_latex", "Output mesh in LaTeX TikZ format", "DMView", flg4, &flg4, PETSC_NULL);CHKERRQ(ierr); 257 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); 258 ierr = PetscOptionsList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr); 259 if (flg) { 260 ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr); 261 } 262 ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,typeName,typeName,sizeof typeName,&flg);CHKERRQ(ierr); 263 if (flg) { 264 ierr = PetscFree(dm->mattype);CHKERRQ(ierr); 265 ierr = PetscStrallocpy(typeName,&dm->mattype);CHKERRQ(ierr); 266 } 267 ierr = PetscOptionsEnum("-dm_is_coloring_type","Global or local coloring of Jacobian","ISColoringType",ISColoringTypes,(PetscEnum)dm->coloringtype,(PetscEnum*)&dm->coloringtype,PETSC_NULL);CHKERRQ(ierr); 268 if (dm->ops->setfromoptions) { 269 ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr); 270 } 271 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 272 ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr); 273 ierr = PetscOptionsEnd();CHKERRQ(ierr); 274 if (flg1) { 275 ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 276 } 277 if (flg2) { 278 PetscViewer viewer; 279 280 ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr); 281 ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 282 ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr); 283 ierr = DMView(dm, viewer);CHKERRQ(ierr); 284 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 285 } 286 if (flg3) { 287 PetscViewer viewer; 288 289 ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr); 290 ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 291 ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr); 292 ierr = PetscViewerFileSetName(viewer, "mesh.vtk");CHKERRQ(ierr); 293 ierr = DMView(dm, viewer);CHKERRQ(ierr); 294 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 295 } 296 if (flg4) { 297 PetscViewer viewer; 298 299 ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr); 300 ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 301 ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_LATEX);CHKERRQ(ierr); 302 ierr = PetscViewerFileSetName(viewer, "mesh.tex");CHKERRQ(ierr); 303 ierr = DMView(dm, viewer);CHKERRQ(ierr); 304 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 305 } 306 PetscFunctionReturn(0); 307 } 308 309 #undef __FUNCT__ 310 #define __FUNCT__ "DMView" 311 /*@C 312 DMView - Views a vector packer or DMDA. 313 314 Collective on DM 315 316 Input Parameter: 317 + dm - the DM object to view 318 - v - the viewer 319 320 Level: developer 321 322 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 323 324 @*/ 325 PetscErrorCode DMView(DM dm,PetscViewer v) 326 { 327 PetscErrorCode ierr; 328 329 PetscFunctionBegin; 330 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 331 if (!v) { 332 ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr); 333 } 334 if (dm->ops->view) { 335 ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr); 336 } 337 PetscFunctionReturn(0); 338 } 339 340 #undef __FUNCT__ 341 #define __FUNCT__ "DMCreateGlobalVector" 342 /*@ 343 DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object 344 345 Collective on DM 346 347 Input Parameter: 348 . dm - the DM object 349 350 Output Parameter: 351 . vec - the global vector 352 353 Level: beginner 354 355 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 356 357 @*/ 358 PetscErrorCode DMCreateGlobalVector(DM dm,Vec *vec) 359 { 360 PetscErrorCode ierr; 361 362 PetscFunctionBegin; 363 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 364 ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr); 365 PetscFunctionReturn(0); 366 } 367 368 #undef __FUNCT__ 369 #define __FUNCT__ "DMCreateLocalVector" 370 /*@ 371 DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object 372 373 Not Collective 374 375 Input Parameter: 376 . dm - the DM object 377 378 Output Parameter: 379 . vec - the local vector 380 381 Level: beginner 382 383 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 384 385 @*/ 386 PetscErrorCode DMCreateLocalVector(DM dm,Vec *vec) 387 { 388 PetscErrorCode ierr; 389 390 PetscFunctionBegin; 391 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 392 ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr); 393 PetscFunctionReturn(0); 394 } 395 396 #undef __FUNCT__ 397 #define __FUNCT__ "DMGetLocalToGlobalMapping" 398 /*@ 399 DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM. 400 401 Collective on DM 402 403 Input Parameter: 404 . dm - the DM that provides the mapping 405 406 Output Parameter: 407 . ltog - the mapping 408 409 Level: intermediate 410 411 Notes: 412 This mapping can then be used by VecSetLocalToGlobalMapping() or 413 MatSetLocalToGlobalMapping(). 414 415 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock() 416 @*/ 417 PetscErrorCode DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog) 418 { 419 PetscErrorCode ierr; 420 421 PetscFunctionBegin; 422 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 423 PetscValidPointer(ltog,2); 424 if (!dm->ltogmap) { 425 if (!dm->ops->createlocaltoglobalmapping) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping"); 426 ierr = (*dm->ops->createlocaltoglobalmapping)(dm);CHKERRQ(ierr); 427 } 428 *ltog = dm->ltogmap; 429 PetscFunctionReturn(0); 430 } 431 432 #undef __FUNCT__ 433 #define __FUNCT__ "DMGetLocalToGlobalMappingBlock" 434 /*@ 435 DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM. 436 437 Collective on DM 438 439 Input Parameter: 440 . da - the distributed array that provides the mapping 441 442 Output Parameter: 443 . ltog - the block mapping 444 445 Level: intermediate 446 447 Notes: 448 This mapping can then be used by VecSetLocalToGlobalMappingBlock() or 449 MatSetLocalToGlobalMappingBlock(). 450 451 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize() 452 @*/ 453 PetscErrorCode DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog) 454 { 455 PetscErrorCode ierr; 456 457 PetscFunctionBegin; 458 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 459 PetscValidPointer(ltog,2); 460 if (!dm->ltogmapb) { 461 PetscInt bs; 462 ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr); 463 if (bs > 1) { 464 if (!dm->ops->createlocaltoglobalmappingblock) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock"); 465 ierr = (*dm->ops->createlocaltoglobalmappingblock)(dm);CHKERRQ(ierr); 466 } else { 467 ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr); 468 ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr); 469 } 470 } 471 *ltog = dm->ltogmapb; 472 PetscFunctionReturn(0); 473 } 474 475 #undef __FUNCT__ 476 #define __FUNCT__ "DMGetBlockSize" 477 /*@ 478 DMGetBlockSize - Gets the inherent block size associated with a DM 479 480 Not Collective 481 482 Input Parameter: 483 . dm - the DM with block structure 484 485 Output Parameter: 486 . bs - the block size, 1 implies no exploitable block structure 487 488 Level: intermediate 489 490 .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock() 491 @*/ 492 PetscErrorCode DMGetBlockSize(DM dm,PetscInt *bs) 493 { 494 PetscFunctionBegin; 495 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 496 PetscValidPointer(bs,2); 497 if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet"); 498 *bs = dm->bs; 499 PetscFunctionReturn(0); 500 } 501 502 #undef __FUNCT__ 503 #define __FUNCT__ "DMCreateInterpolation" 504 /*@ 505 DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects 506 507 Collective on DM 508 509 Input Parameter: 510 + dm1 - the DM object 511 - dm2 - the second, finer DM object 512 513 Output Parameter: 514 + mat - the interpolation 515 - vec - the scaling (optional) 516 517 Level: developer 518 519 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 520 DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation. 521 522 For DMDA objects you can use this interpolation (more precisely the interpolation from the DMDAGetCoordinateDA()) to interpolate the mesh coordinate vectors 523 EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic. 524 525 526 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen() 527 528 @*/ 529 PetscErrorCode DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec) 530 { 531 PetscErrorCode ierr; 532 533 PetscFunctionBegin; 534 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 535 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 536 ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr); 537 PetscFunctionReturn(0); 538 } 539 540 #undef __FUNCT__ 541 #define __FUNCT__ "DMCreateInjection" 542 /*@ 543 DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects 544 545 Collective on DM 546 547 Input Parameter: 548 + dm1 - the DM object 549 - dm2 - the second, finer DM object 550 551 Output Parameter: 552 . ctx - the injection 553 554 Level: developer 555 556 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 557 DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection. 558 559 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation() 560 561 @*/ 562 PetscErrorCode DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx) 563 { 564 PetscErrorCode ierr; 565 566 PetscFunctionBegin; 567 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 568 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 569 ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr); 570 PetscFunctionReturn(0); 571 } 572 573 #undef __FUNCT__ 574 #define __FUNCT__ "DMCreateColoring" 575 /*@C 576 DMCreateColoring - Gets coloring for a DMDA or DMComposite 577 578 Collective on DM 579 580 Input Parameter: 581 + dm - the DM object 582 . ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL 583 - matype - either MATAIJ or MATBAIJ 584 585 Output Parameter: 586 . coloring - the coloring 587 588 Level: developer 589 590 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix() 591 592 @*/ 593 PetscErrorCode DMCreateColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring) 594 { 595 PetscErrorCode ierr; 596 597 PetscFunctionBegin; 598 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 599 if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet"); 600 ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr); 601 PetscFunctionReturn(0); 602 } 603 604 #undef __FUNCT__ 605 #define __FUNCT__ "DMCreateMatrix" 606 /*@C 607 DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite 608 609 Collective on DM 610 611 Input Parameter: 612 + dm - the DM object 613 - mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or 614 any type which inherits from one of these (such as MATAIJ) 615 616 Output Parameter: 617 . mat - the empty Jacobian 618 619 Level: beginner 620 621 Notes: This properly preallocates the number of nonzeros in the sparse matrix so you 622 do not need to do it yourself. 623 624 By default it also sets the nonzero structure and puts in the zero entries. To prevent setting 625 the nonzero pattern call DMDASetMatPreallocateOnly() 626 627 For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used 628 internally by PETSc. 629 630 For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires 631 the indices for the global numbering for DMDAs which is complicated. 632 633 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 634 635 @*/ 636 PetscErrorCode DMCreateMatrix(DM dm,const MatType mtype,Mat *mat) 637 { 638 PetscErrorCode ierr; 639 640 PetscFunctionBegin; 641 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 642 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 643 ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr); 644 #endif 645 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 646 PetscValidPointer(mat,3); 647 if (dm->mattype) { 648 ierr = (*dm->ops->creatematrix)(dm,dm->mattype,mat);CHKERRQ(ierr); 649 } else { 650 ierr = (*dm->ops->creatematrix)(dm,mtype,mat);CHKERRQ(ierr); 651 } 652 PetscFunctionReturn(0); 653 } 654 655 #undef __FUNCT__ 656 #define __FUNCT__ "DMSetMatrixPreallocateOnly" 657 /*@ 658 DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly 659 preallocated but the nonzero structure and zero values will not be set. 660 661 Logically Collective on DMDA 662 663 Input Parameter: 664 + dm - the DM 665 - only - PETSC_TRUE if only want preallocation 666 667 Level: developer 668 .seealso DMCreateMatrix() 669 @*/ 670 PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only) 671 { 672 PetscFunctionBegin; 673 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 674 dm->prealloc_only = only; 675 PetscFunctionReturn(0); 676 } 677 678 #undef __FUNCT__ 679 #define __FUNCT__ "DMGetWorkArray" 680 /*@C 681 DMGetWorkArray - Gets a work array guaranteed to be at least the input size 682 683 Not Collective 684 685 Input Parameters: 686 + dm - the DM object 687 - size - The minium size 688 689 Output Parameter: 690 . array - the work array 691 692 Level: developer 693 694 .seealso DMDestroy(), DMCreate() 695 @*/ 696 PetscErrorCode DMGetWorkArray(DM dm,PetscInt size,PetscScalar **array) 697 { 698 PetscErrorCode ierr; 699 700 PetscFunctionBegin; 701 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 702 PetscValidPointer(array,3); 703 if (size > dm->workSize) { 704 dm->workSize = size; 705 ierr = PetscFree(dm->workArray);CHKERRQ(ierr); 706 ierr = PetscMalloc(dm->workSize * sizeof(PetscScalar), &dm->workArray);CHKERRQ(ierr); 707 } 708 *array = dm->workArray; 709 PetscFunctionReturn(0); 710 } 711 712 #undef __FUNCT__ 713 #define __FUNCT__ "DMCreateDecompositionDM" 714 /*@C 715 DMCreateDecompositionDM - creates a DM that defines a decomposition of the original DM. 716 717 Not Collective 718 719 Input Parameters: 720 + dm - the DM object 721 - name - the name of the decomposition 722 723 Output Parameter: 724 . ddm - the decomposition DM (PETSC_NULL, if no such decomposition is known) 725 726 Level: advanced 727 728 .seealso DMDestroy(), DMCreate(), DMCreateDecomposition() 729 @*/ 730 PetscErrorCode DMCreateDecompositionDM(DM dm, const char* name, DM *ddm) 731 { 732 PetscErrorCode ierr; 733 734 PetscFunctionBegin; 735 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 736 PetscValidCharPointer(name,2); 737 PetscValidPointer(ddm,3); 738 if(!dm->ops->createdecompositiondm) { 739 *ddm = PETSC_NULL; 740 } 741 else { 742 ierr = (*dm->ops->createdecompositiondm)(dm,name,ddm); CHKERRQ(ierr); 743 } 744 PetscFunctionReturn(0); 745 } 746 747 #undef __FUNCT__ 748 #define __FUNCT__ "DMCreateFieldIS" 749 /*@C 750 DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field 751 752 Not collective 753 754 Input Parameter: 755 . dm - the DM object 756 757 Output Parameters: 758 + numFields - The number of fields (or PETSC_NULL if not requested) 759 . names - The name for each field (or PETSC_NULL if not requested) 760 - fields - The global indices for each field (or PETSC_NULL if not requested) 761 762 Level: intermediate 763 764 Notes: 765 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 766 PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with 767 PetscFree(). 768 769 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 770 @*/ 771 PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***names, IS **fields) 772 { 773 PetscErrorCode ierr; 774 775 PetscFunctionBegin; 776 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 777 if (numFields) {PetscValidPointer(numFields,2);} 778 if (names) {PetscValidPointer(names,3);} 779 if (fields) {PetscValidPointer(fields,4);} 780 ierr = (*dm->ops->createfieldis)(dm, numFields, names, fields);CHKERRQ(ierr); 781 PetscFunctionReturn(0); 782 } 783 784 785 #undef __FUNCT__ 786 #define __FUNCT__ "DMCreateDecomposition" 787 /*@C 788 DMCreateDecomposition - Returns a list of IS objects defining a decomposition of a problem into subproblems: 789 each IS contains the global indices of the dofs of the corresponding subproblem. 790 The optional list of DMs define the DM for each subproblem. 791 Generalizes DMCreateFieldIS(). 792 793 Not collective 794 795 Input Parameter: 796 . dm - the DM object 797 798 Output Parameters: 799 + len - The number of subproblems in the decomposition (or PETSC_NULL if not requested) 800 . namelist - The name for each subproblem (or PETSC_NULL if not requested) 801 . islist - The global indices for each subproblem (or PETSC_NULL if not requested) 802 - dmlist - The DMs for each subproblem (or PETSC_NULL, if not requested; if PETSC_NULL is returned, no DMs are defined) 803 804 Level: intermediate 805 806 Notes: 807 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 808 PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(), 809 and all of the arrays should be freed with PetscFree(). 810 811 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS() 812 @*/ 813 PetscErrorCode DMCreateDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 814 { 815 PetscErrorCode ierr; 816 817 PetscFunctionBegin; 818 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 819 if (len) {PetscValidPointer(len,2);} 820 if (namelist) {PetscValidPointer(namelist,3);} 821 if (islist) {PetscValidPointer(islist,4);} 822 if (dmlist) {PetscValidPointer(dmlist,5);} 823 if(!dm->ops->createdecomposition) { 824 ierr = (*dm->ops->createfieldis)(dm, len, namelist, islist);CHKERRQ(ierr); 825 /* By default there are no DMs associated with subproblems. */ 826 if(dmlist) *dmlist = PETSC_NULL; 827 } 828 else { 829 ierr = (*dm->ops->createdecomposition)(dm,len,namelist,islist,dmlist); CHKERRQ(ierr); 830 } 831 PetscFunctionReturn(0); 832 } 833 834 835 #undef __FUNCT__ 836 #define __FUNCT__ "DMRefine" 837 /*@ 838 DMRefine - Refines a DM object 839 840 Collective on DM 841 842 Input Parameter: 843 + dm - the DM object 844 - comm - the communicator to contain the new DM object (or PETSC_NULL) 845 846 Output Parameter: 847 . dmf - the refined DM 848 849 Level: developer 850 851 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 852 853 @*/ 854 PetscErrorCode DMRefine(DM dm,MPI_Comm comm,DM *dmf) 855 { 856 PetscErrorCode ierr; 857 858 PetscFunctionBegin; 859 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 860 ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr); 861 if (*dmf) { 862 (*dmf)->ops->creatematrix = dm->ops->creatematrix; 863 (*dmf)->ops->initialguess = dm->ops->initialguess; 864 (*dmf)->ops->function = dm->ops->function; 865 (*dmf)->ops->functionj = dm->ops->functionj; 866 if (dm->ops->jacobian != DMComputeJacobianDefault) { 867 (*dmf)->ops->jacobian = dm->ops->jacobian; 868 } 869 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr); 870 (*dmf)->ctx = dm->ctx; 871 (*dmf)->levelup = dm->levelup + 1; 872 } 873 PetscFunctionReturn(0); 874 } 875 876 #undef __FUNCT__ 877 #define __FUNCT__ "DMGetRefineLevel" 878 /*@ 879 DMGetRefineLevel - Get's the number of refinements that have generated this DM. 880 881 Not Collective 882 883 Input Parameter: 884 . dm - the DM object 885 886 Output Parameter: 887 . level - number of refinements 888 889 Level: developer 890 891 .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 892 893 @*/ 894 PetscErrorCode DMGetRefineLevel(DM dm,PetscInt *level) 895 { 896 PetscFunctionBegin; 897 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 898 *level = dm->levelup; 899 PetscFunctionReturn(0); 900 } 901 902 #undef __FUNCT__ 903 #define __FUNCT__ "DMGlobalToLocalBegin" 904 /*@ 905 DMGlobalToLocalBegin - Begins updating local vectors from global vector 906 907 Neighbor-wise Collective on DM 908 909 Input Parameters: 910 + dm - the DM object 911 . g - the global vector 912 . mode - INSERT_VALUES or ADD_VALUES 913 - l - the local vector 914 915 916 Level: beginner 917 918 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 919 920 @*/ 921 PetscErrorCode DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) 922 { 923 PetscErrorCode ierr; 924 925 PetscFunctionBegin; 926 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 927 ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 928 PetscFunctionReturn(0); 929 } 930 931 #undef __FUNCT__ 932 #define __FUNCT__ "DMGlobalToLocalEnd" 933 /*@ 934 DMGlobalToLocalEnd - Ends updating local vectors from global vector 935 936 Neighbor-wise Collective on DM 937 938 Input Parameters: 939 + dm - the DM object 940 . g - the global vector 941 . mode - INSERT_VALUES or ADD_VALUES 942 - l - the local vector 943 944 945 Level: beginner 946 947 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 948 949 @*/ 950 PetscErrorCode DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) 951 { 952 PetscErrorCode ierr; 953 954 PetscFunctionBegin; 955 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 956 ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 957 PetscFunctionReturn(0); 958 } 959 960 #undef __FUNCT__ 961 #define __FUNCT__ "DMLocalToGlobalBegin" 962 /*@ 963 DMLocalToGlobalBegin - updates global vectors from local vectors 964 965 Neighbor-wise Collective on DM 966 967 Input Parameters: 968 + dm - the DM object 969 . l - the local vector 970 . 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 971 base point. 972 - - the global vector 973 974 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 975 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 976 global array to the final global array with VecAXPY(). 977 978 Level: beginner 979 980 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin() 981 982 @*/ 983 PetscErrorCode DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g) 984 { 985 PetscErrorCode ierr; 986 987 PetscFunctionBegin; 988 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 989 ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 990 PetscFunctionReturn(0); 991 } 992 993 #undef __FUNCT__ 994 #define __FUNCT__ "DMLocalToGlobalEnd" 995 /*@ 996 DMLocalToGlobalEnd - updates global vectors from local vectors 997 998 Neighbor-wise Collective on DM 999 1000 Input Parameters: 1001 + dm - the DM object 1002 . l - the local vector 1003 . mode - INSERT_VALUES or ADD_VALUES 1004 - g - the global vector 1005 1006 1007 Level: beginner 1008 1009 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd() 1010 1011 @*/ 1012 PetscErrorCode DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g) 1013 { 1014 PetscErrorCode ierr; 1015 1016 PetscFunctionBegin; 1017 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1018 ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 1019 PetscFunctionReturn(0); 1020 } 1021 1022 #undef __FUNCT__ 1023 #define __FUNCT__ "DMComputeJacobianDefault" 1024 /*@ 1025 DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided 1026 1027 Collective on DM 1028 1029 Input Parameter: 1030 + dm - the DM object 1031 . x - location to compute Jacobian at; may be ignored for linear problems 1032 . A - matrix that defines the operator for the linear solve 1033 - B - the matrix used to construct the preconditioner 1034 1035 Level: developer 1036 1037 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1038 DMSetFunction() 1039 1040 @*/ 1041 PetscErrorCode DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 1042 { 1043 PetscErrorCode ierr; 1044 1045 PetscFunctionBegin; 1046 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1047 *stflag = SAME_NONZERO_PATTERN; 1048 ierr = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr); 1049 if (A != B) { 1050 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1051 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1052 } 1053 PetscFunctionReturn(0); 1054 } 1055 1056 #undef __FUNCT__ 1057 #define __FUNCT__ "DMCoarsen" 1058 /*@ 1059 DMCoarsen - Coarsens a DM object 1060 1061 Collective on DM 1062 1063 Input Parameter: 1064 + dm - the DM object 1065 - comm - the communicator to contain the new DM object (or PETSC_NULL) 1066 1067 Output Parameter: 1068 . dmc - the coarsened DM 1069 1070 Level: developer 1071 1072 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1073 1074 @*/ 1075 PetscErrorCode DMCoarsen(DM dm, MPI_Comm comm, DM *dmc) 1076 { 1077 PetscErrorCode ierr; 1078 DMCoarsenHookLink link; 1079 1080 PetscFunctionBegin; 1081 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1082 ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr); 1083 (*dmc)->ops->creatematrix = dm->ops->creatematrix; 1084 (*dmc)->ops->initialguess = dm->ops->initialguess; 1085 (*dmc)->ops->function = dm->ops->function; 1086 (*dmc)->ops->functionj = dm->ops->functionj; 1087 if (dm->ops->jacobian != DMComputeJacobianDefault) { 1088 (*dmc)->ops->jacobian = dm->ops->jacobian; 1089 } 1090 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr); 1091 (*dmc)->ctx = dm->ctx; 1092 (*dmc)->leveldown = dm->leveldown + 1; 1093 for (link=dm->coarsenhook; link; link=link->next) { 1094 if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);} 1095 } 1096 PetscFunctionReturn(0); 1097 } 1098 1099 #undef __FUNCT__ 1100 #define __FUNCT__ "DMCoarsenHookAdd" 1101 /*@ 1102 DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid 1103 1104 Logically Collective 1105 1106 Input Arguments: 1107 + fine - nonlinear solver context on which to run a hook when restricting to a coarser level 1108 . coarsenhook - function to run when setting up a coarser level 1109 . restricthook - function to run to update data on coarser levels (once per SNESSolve()) 1110 - ctx - [optional] user-defined context for provide data for the hooks (may be PETSC_NULL) 1111 1112 Calling sequence of coarsenhook: 1113 $ coarsenhook(DM fine,DM coarse,void *ctx); 1114 1115 + fine - fine level DM 1116 . coarse - coarse level DM to restrict problem to 1117 - ctx - optional user-defined function context 1118 1119 Calling sequence for restricthook: 1120 $ restricthook(DM fine,Mat mrestrict,Mat inject,DM coarse,void *ctx) 1121 1122 + fine - fine level DM 1123 . mrestrict - matrix restricting a fine-level solution to the coarse grid 1124 . inject - matrix restricting by applying the transpose of injection 1125 . coarse - coarse level DM to update 1126 - ctx - optional user-defined function context 1127 1128 Level: advanced 1129 1130 Notes: 1131 This function is only needed if auxiliary data needs to be set up on coarse grids. 1132 1133 If this function is called multiple times, the hooks will be run in the order they are added. 1134 1135 The 1136 1137 In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to 1138 extract the finest level information from its context (instead of from the SNES). 1139 1140 .seealso: SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 1141 @*/ 1142 PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx) 1143 { 1144 PetscErrorCode ierr; 1145 DMCoarsenHookLink link,*p; 1146 1147 PetscFunctionBegin; 1148 PetscValidHeaderSpecific(fine,DM_CLASSID,1); 1149 for (p=&fine->coarsenhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */ 1150 ierr = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr); 1151 link->coarsenhook = coarsenhook; 1152 link->restricthook = restricthook; 1153 link->ctx = ctx; 1154 link->next = PETSC_NULL; 1155 *p = link; 1156 PetscFunctionReturn(0); 1157 } 1158 1159 #undef __FUNCT__ 1160 #define __FUNCT__ "DMRestrict" 1161 /*@ 1162 DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd() 1163 1164 Collective if any hooks are 1165 1166 Input Arguments: 1167 + fine - finer DM to use as a base 1168 . restrct - restriction matrix, apply using MatRestrict() 1169 . inject - injection matrix, also use MatRestrict() 1170 - coarse - coarer DM to update 1171 1172 Level: developer 1173 1174 .seealso: DMCoarsenHookAdd(), MatRestrict() 1175 @*/ 1176 PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse) 1177 { 1178 PetscErrorCode ierr; 1179 DMCoarsenHookLink link; 1180 1181 PetscFunctionBegin; 1182 for (link=fine->coarsenhook; link; link=link->next) { 1183 if (link->restricthook) {ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);} 1184 } 1185 PetscFunctionReturn(0); 1186 } 1187 1188 #undef __FUNCT__ 1189 #define __FUNCT__ "DMGetCoarsenLevel" 1190 /*@ 1191 DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM. 1192 1193 Not Collective 1194 1195 Input Parameter: 1196 . dm - the DM object 1197 1198 Output Parameter: 1199 . level - number of coarsenings 1200 1201 Level: developer 1202 1203 .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1204 1205 @*/ 1206 PetscErrorCode DMGetCoarsenLevel(DM dm,PetscInt *level) 1207 { 1208 PetscFunctionBegin; 1209 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1210 *level = dm->leveldown; 1211 PetscFunctionReturn(0); 1212 } 1213 1214 1215 1216 #undef __FUNCT__ 1217 #define __FUNCT__ "DMRefineHierarchy" 1218 /*@C 1219 DMRefineHierarchy - Refines a DM object, all levels at once 1220 1221 Collective on DM 1222 1223 Input Parameter: 1224 + dm - the DM object 1225 - nlevels - the number of levels of refinement 1226 1227 Output Parameter: 1228 . dmf - the refined DM hierarchy 1229 1230 Level: developer 1231 1232 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1233 1234 @*/ 1235 PetscErrorCode DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[]) 1236 { 1237 PetscErrorCode ierr; 1238 1239 PetscFunctionBegin; 1240 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1241 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 1242 if (nlevels == 0) PetscFunctionReturn(0); 1243 if (dm->ops->refinehierarchy) { 1244 ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr); 1245 } else if (dm->ops->refine) { 1246 PetscInt i; 1247 1248 ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr); 1249 for (i=1; i<nlevels; i++) { 1250 ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr); 1251 } 1252 } else { 1253 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet"); 1254 } 1255 PetscFunctionReturn(0); 1256 } 1257 1258 #undef __FUNCT__ 1259 #define __FUNCT__ "DMCoarsenHierarchy" 1260 /*@C 1261 DMCoarsenHierarchy - Coarsens a DM object, all levels at once 1262 1263 Collective on DM 1264 1265 Input Parameter: 1266 + dm - the DM object 1267 - nlevels - the number of levels of coarsening 1268 1269 Output Parameter: 1270 . dmc - the coarsened DM hierarchy 1271 1272 Level: developer 1273 1274 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1275 1276 @*/ 1277 PetscErrorCode DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[]) 1278 { 1279 PetscErrorCode ierr; 1280 1281 PetscFunctionBegin; 1282 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1283 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 1284 if (nlevels == 0) PetscFunctionReturn(0); 1285 PetscValidPointer(dmc,3); 1286 if (dm->ops->coarsenhierarchy) { 1287 ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr); 1288 } else if (dm->ops->coarsen) { 1289 PetscInt i; 1290 1291 ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr); 1292 for (i=1; i<nlevels; i++) { 1293 ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr); 1294 } 1295 } else { 1296 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet"); 1297 } 1298 PetscFunctionReturn(0); 1299 } 1300 1301 #undef __FUNCT__ 1302 #define __FUNCT__ "DMCreateAggregates" 1303 /*@ 1304 DMCreateAggregates - Gets the aggregates that map between 1305 grids associated with two DMs. 1306 1307 Collective on DM 1308 1309 Input Parameters: 1310 + dmc - the coarse grid DM 1311 - dmf - the fine grid DM 1312 1313 Output Parameters: 1314 . rest - the restriction matrix (transpose of the projection matrix) 1315 1316 Level: intermediate 1317 1318 .keywords: interpolation, restriction, multigrid 1319 1320 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation() 1321 @*/ 1322 PetscErrorCode DMCreateAggregates(DM dmc, DM dmf, Mat *rest) 1323 { 1324 PetscErrorCode ierr; 1325 1326 PetscFunctionBegin; 1327 PetscValidHeaderSpecific(dmc,DM_CLASSID,1); 1328 PetscValidHeaderSpecific(dmf,DM_CLASSID,2); 1329 ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr); 1330 PetscFunctionReturn(0); 1331 } 1332 1333 #undef __FUNCT__ 1334 #define __FUNCT__ "DMSetApplicationContextDestroy" 1335 /*@C 1336 DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed 1337 1338 Not Collective 1339 1340 Input Parameters: 1341 + dm - the DM object 1342 - destroy - the destroy function 1343 1344 Level: intermediate 1345 1346 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1347 1348 @*/ 1349 PetscErrorCode DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**)) 1350 { 1351 PetscFunctionBegin; 1352 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1353 dm->ctxdestroy = destroy; 1354 PetscFunctionReturn(0); 1355 } 1356 1357 #undef __FUNCT__ 1358 #define __FUNCT__ "DMSetApplicationContext" 1359 /*@ 1360 DMSetApplicationContext - Set a user context into a DM object 1361 1362 Not Collective 1363 1364 Input Parameters: 1365 + dm - the DM object 1366 - ctx - the user context 1367 1368 Level: intermediate 1369 1370 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1371 1372 @*/ 1373 PetscErrorCode DMSetApplicationContext(DM dm,void *ctx) 1374 { 1375 PetscFunctionBegin; 1376 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1377 dm->ctx = ctx; 1378 PetscFunctionReturn(0); 1379 } 1380 1381 #undef __FUNCT__ 1382 #define __FUNCT__ "DMGetApplicationContext" 1383 /*@ 1384 DMGetApplicationContext - Gets a user context from a DM object 1385 1386 Not Collective 1387 1388 Input Parameter: 1389 . dm - the DM object 1390 1391 Output Parameter: 1392 . ctx - the user context 1393 1394 Level: intermediate 1395 1396 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1397 1398 @*/ 1399 PetscErrorCode DMGetApplicationContext(DM dm,void *ctx) 1400 { 1401 PetscFunctionBegin; 1402 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1403 *(void**)ctx = dm->ctx; 1404 PetscFunctionReturn(0); 1405 } 1406 1407 #undef __FUNCT__ 1408 #define __FUNCT__ "DMSetInitialGuess" 1409 /*@C 1410 DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers 1411 1412 Logically Collective on DM 1413 1414 Input Parameter: 1415 + dm - the DM object to destroy 1416 - f - the function to compute the initial guess 1417 1418 Level: intermediate 1419 1420 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1421 1422 @*/ 1423 PetscErrorCode DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec)) 1424 { 1425 PetscFunctionBegin; 1426 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1427 dm->ops->initialguess = f; 1428 PetscFunctionReturn(0); 1429 } 1430 1431 #undef __FUNCT__ 1432 #define __FUNCT__ "DMSetFunction" 1433 /*@C 1434 DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES 1435 1436 Logically Collective on DM 1437 1438 Input Parameter: 1439 + dm - the DM object 1440 - f - the function to compute (use PETSC_NULL to cancel a previous function that was set) 1441 1442 Level: intermediate 1443 1444 Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian 1445 computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian. 1446 1447 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1448 DMSetJacobian() 1449 1450 @*/ 1451 PetscErrorCode DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 1452 { 1453 PetscFunctionBegin; 1454 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1455 dm->ops->function = f; 1456 if (f) { 1457 dm->ops->functionj = f; 1458 } 1459 PetscFunctionReturn(0); 1460 } 1461 1462 #undef __FUNCT__ 1463 #define __FUNCT__ "DMSetJacobian" 1464 /*@C 1465 DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES 1466 1467 Logically Collective on DM 1468 1469 Input Parameter: 1470 + dm - the DM object to destroy 1471 - f - the function to compute the matrix entries 1472 1473 Level: intermediate 1474 1475 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1476 DMSetFunction() 1477 1478 @*/ 1479 PetscErrorCode DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*)) 1480 { 1481 PetscFunctionBegin; 1482 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1483 dm->ops->jacobian = f; 1484 PetscFunctionReturn(0); 1485 } 1486 1487 #undef __FUNCT__ 1488 #define __FUNCT__ "DMSetVariableBounds" 1489 /*@C 1490 DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI. 1491 1492 Logically Collective on DM 1493 1494 Input Parameter: 1495 + dm - the DM object 1496 - f - the function that computes variable bounds used by SNESVI (use PETSC_NULL to cancel a previous function that was set) 1497 1498 Level: intermediate 1499 1500 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1501 DMSetJacobian() 1502 1503 @*/ 1504 PetscErrorCode DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 1505 { 1506 PetscFunctionBegin; 1507 dm->ops->computevariablebounds = f; 1508 PetscFunctionReturn(0); 1509 } 1510 1511 #undef __FUNCT__ 1512 #define __FUNCT__ "DMHasVariableBounds" 1513 /*@ 1514 DMHasVariableBounds - does the DM object have a variable bounds function? 1515 1516 Not Collective 1517 1518 Input Parameter: 1519 . dm - the DM object to destroy 1520 1521 Output Parameter: 1522 . flg - PETSC_TRUE if the variable bounds function exists 1523 1524 Level: developer 1525 1526 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1527 1528 @*/ 1529 PetscErrorCode DMHasVariableBounds(DM dm,PetscBool *flg) 1530 { 1531 PetscFunctionBegin; 1532 *flg = (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE; 1533 PetscFunctionReturn(0); 1534 } 1535 1536 #undef __FUNCT__ 1537 #define __FUNCT__ "DMComputeVariableBounds" 1538 /*@C 1539 DMComputeVariableBounds - compute variable bounds used by SNESVI. 1540 1541 Logically Collective on DM 1542 1543 Input Parameters: 1544 + dm - the DM object to destroy 1545 - x - current solution at which the bounds are computed 1546 1547 Output parameters: 1548 + xl - lower bound 1549 - xu - upper bound 1550 1551 Level: intermediate 1552 1553 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1554 DMSetFunction(), DMSetVariableBounds() 1555 1556 @*/ 1557 PetscErrorCode DMComputeVariableBounds(DM dm, Vec xl, Vec xu) 1558 { 1559 PetscErrorCode ierr; 1560 PetscFunctionBegin; 1561 PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 1562 PetscValidHeaderSpecific(xu,VEC_CLASSID,2); 1563 if(dm->ops->computevariablebounds) { 1564 ierr = (*dm->ops->computevariablebounds)(dm, xl,xu); CHKERRQ(ierr); 1565 } 1566 else { 1567 ierr = VecSet(xl,SNES_VI_NINF); CHKERRQ(ierr); 1568 ierr = VecSet(xu,SNES_VI_INF); CHKERRQ(ierr); 1569 } 1570 PetscFunctionReturn(0); 1571 } 1572 1573 #undef __FUNCT__ 1574 #define __FUNCT__ "DMComputeInitialGuess" 1575 /*@ 1576 DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers 1577 1578 Collective on DM 1579 1580 Input Parameter: 1581 + dm - the DM object to destroy 1582 - x - the vector to hold the initial guess values 1583 1584 Level: developer 1585 1586 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat() 1587 1588 @*/ 1589 PetscErrorCode DMComputeInitialGuess(DM dm,Vec x) 1590 { 1591 PetscErrorCode ierr; 1592 PetscFunctionBegin; 1593 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1594 if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()"); 1595 ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr); 1596 PetscFunctionReturn(0); 1597 } 1598 1599 #undef __FUNCT__ 1600 #define __FUNCT__ "DMHasInitialGuess" 1601 /*@ 1602 DMHasInitialGuess - does the DM object have an initial guess function 1603 1604 Not Collective 1605 1606 Input Parameter: 1607 . dm - the DM object to destroy 1608 1609 Output Parameter: 1610 . flg - PETSC_TRUE if function exists 1611 1612 Level: developer 1613 1614 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1615 1616 @*/ 1617 PetscErrorCode DMHasInitialGuess(DM dm,PetscBool *flg) 1618 { 1619 PetscFunctionBegin; 1620 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1621 *flg = (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE; 1622 PetscFunctionReturn(0); 1623 } 1624 1625 #undef __FUNCT__ 1626 #define __FUNCT__ "DMHasFunction" 1627 /*@ 1628 DMHasFunction - does the DM object have a function 1629 1630 Not Collective 1631 1632 Input Parameter: 1633 . dm - the DM object to destroy 1634 1635 Output Parameter: 1636 . flg - PETSC_TRUE if function exists 1637 1638 Level: developer 1639 1640 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1641 1642 @*/ 1643 PetscErrorCode DMHasFunction(DM dm,PetscBool *flg) 1644 { 1645 PetscFunctionBegin; 1646 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1647 *flg = (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE; 1648 PetscFunctionReturn(0); 1649 } 1650 1651 #undef __FUNCT__ 1652 #define __FUNCT__ "DMHasJacobian" 1653 /*@ 1654 DMHasJacobian - does the DM object have a matrix function 1655 1656 Not Collective 1657 1658 Input Parameter: 1659 . dm - the DM object to destroy 1660 1661 Output Parameter: 1662 . flg - PETSC_TRUE if function exists 1663 1664 Level: developer 1665 1666 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1667 1668 @*/ 1669 PetscErrorCode DMHasJacobian(DM dm,PetscBool *flg) 1670 { 1671 PetscFunctionBegin; 1672 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1673 *flg = (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE; 1674 PetscFunctionReturn(0); 1675 } 1676 1677 #undef __FUNCT__ 1678 #define __FUNCT__ "DMSetVec" 1679 /*@C 1680 DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear. 1681 1682 Collective on DM 1683 1684 Input Parameter: 1685 + dm - the DM object 1686 - x - location to compute residual and Jacobian, if PETSC_NULL is passed to those routines; will be PETSC_NULL for linear problems. 1687 1688 Level: developer 1689 1690 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1691 DMSetFunction(), DMSetJacobian(), DMSetVariableBounds() 1692 1693 @*/ 1694 PetscErrorCode DMSetVec(DM dm,Vec x) 1695 { 1696 PetscErrorCode ierr; 1697 PetscFunctionBegin; 1698 if (x) { 1699 if (!dm->x) { 1700 ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 1701 } 1702 ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 1703 } 1704 else if(dm->x) { 1705 ierr = VecDestroy(&dm->x); CHKERRQ(ierr); 1706 } 1707 PetscFunctionReturn(0); 1708 } 1709 1710 1711 #undef __FUNCT__ 1712 #define __FUNCT__ "DMComputeFunction" 1713 /*@ 1714 DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES 1715 1716 Collective on DM 1717 1718 Input Parameter: 1719 + dm - the DM object to destroy 1720 . x - the location where the function is evaluationed, may be ignored for linear problems 1721 - b - the vector to hold the right hand side entries 1722 1723 Level: developer 1724 1725 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1726 DMSetJacobian() 1727 1728 @*/ 1729 PetscErrorCode DMComputeFunction(DM dm,Vec x,Vec b) 1730 { 1731 PetscErrorCode ierr; 1732 PetscFunctionBegin; 1733 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1734 if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()"); 1735 PetscStackPush("DM user function"); 1736 ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr); 1737 PetscStackPop; 1738 PetscFunctionReturn(0); 1739 } 1740 1741 1742 1743 #undef __FUNCT__ 1744 #define __FUNCT__ "DMComputeJacobian" 1745 /*@ 1746 DMComputeJacobian - compute the matrix entries for the solver 1747 1748 Collective on DM 1749 1750 Input Parameter: 1751 + dm - the DM object 1752 . x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM 1753 . A - matrix that defines the operator for the linear solve 1754 - B - the matrix used to construct the preconditioner 1755 1756 Level: developer 1757 1758 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1759 DMSetFunction() 1760 1761 @*/ 1762 PetscErrorCode DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 1763 { 1764 PetscErrorCode ierr; 1765 1766 PetscFunctionBegin; 1767 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1768 if (!dm->ops->jacobian) { 1769 ISColoring coloring; 1770 MatFDColoring fd; 1771 1772 ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr); 1773 ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr); 1774 ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 1775 ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr); 1776 ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr); 1777 ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr); 1778 1779 dm->fd = fd; 1780 dm->ops->jacobian = DMComputeJacobianDefault; 1781 1782 /* don't know why this is needed */ 1783 ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 1784 } 1785 if (!x) x = dm->x; 1786 ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr); 1787 1788 /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */ 1789 if (x) { 1790 if (!dm->x) { 1791 ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 1792 } 1793 ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 1794 } 1795 if (A != B) { 1796 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1797 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1798 } 1799 PetscFunctionReturn(0); 1800 } 1801 1802 1803 PetscFList DMList = PETSC_NULL; 1804 PetscBool DMRegisterAllCalled = PETSC_FALSE; 1805 1806 #undef __FUNCT__ 1807 #define __FUNCT__ "DMSetType" 1808 /*@C 1809 DMSetType - Builds a DM, for a particular DM implementation. 1810 1811 Collective on DM 1812 1813 Input Parameters: 1814 + dm - The DM object 1815 - method - The name of the DM type 1816 1817 Options Database Key: 1818 . -dm_type <type> - Sets the DM type; use -help for a list of available types 1819 1820 Notes: 1821 See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D). 1822 1823 Level: intermediate 1824 1825 .keywords: DM, set, type 1826 .seealso: DMGetType(), DMCreate() 1827 @*/ 1828 PetscErrorCode DMSetType(DM dm, const DMType method) 1829 { 1830 PetscErrorCode (*r)(DM); 1831 PetscBool match; 1832 PetscErrorCode ierr; 1833 1834 PetscFunctionBegin; 1835 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1836 ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr); 1837 if (match) PetscFunctionReturn(0); 1838 1839 if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1840 ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 1841 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method); 1842 1843 if (dm->ops->destroy) { 1844 ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr); 1845 dm->ops->destroy = PETSC_NULL; 1846 } 1847 ierr = (*r)(dm);CHKERRQ(ierr); 1848 ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr); 1849 PetscFunctionReturn(0); 1850 } 1851 1852 #undef __FUNCT__ 1853 #define __FUNCT__ "DMGetType" 1854 /*@C 1855 DMGetType - Gets the DM type name (as a string) from the DM. 1856 1857 Not Collective 1858 1859 Input Parameter: 1860 . dm - The DM 1861 1862 Output Parameter: 1863 . type - The DM type name 1864 1865 Level: intermediate 1866 1867 .keywords: DM, get, type, name 1868 .seealso: DMSetType(), DMCreate() 1869 @*/ 1870 PetscErrorCode DMGetType(DM dm, const DMType *type) 1871 { 1872 PetscErrorCode ierr; 1873 1874 PetscFunctionBegin; 1875 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1876 PetscValidCharPointer(type,2); 1877 if (!DMRegisterAllCalled) { 1878 ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1879 } 1880 *type = ((PetscObject)dm)->type_name; 1881 PetscFunctionReturn(0); 1882 } 1883 1884 #undef __FUNCT__ 1885 #define __FUNCT__ "DMConvert" 1886 /*@C 1887 DMConvert - Converts a DM to another DM, either of the same or different type. 1888 1889 Collective on DM 1890 1891 Input Parameters: 1892 + dm - the DM 1893 - newtype - new DM type (use "same" for the same type) 1894 1895 Output Parameter: 1896 . M - pointer to new DM 1897 1898 Notes: 1899 Cannot be used to convert a sequential DM to parallel or parallel to sequential, 1900 the MPI communicator of the generated DM is always the same as the communicator 1901 of the input DM. 1902 1903 Level: intermediate 1904 1905 .seealso: DMCreate() 1906 @*/ 1907 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M) 1908 { 1909 DM B; 1910 char convname[256]; 1911 PetscBool sametype, issame; 1912 PetscErrorCode ierr; 1913 1914 PetscFunctionBegin; 1915 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1916 PetscValidType(dm,1); 1917 PetscValidPointer(M,3); 1918 ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr); 1919 ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr); 1920 { 1921 PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL; 1922 1923 /* 1924 Order of precedence: 1925 1) See if a specialized converter is known to the current DM. 1926 2) See if a specialized converter is known to the desired DM class. 1927 3) See if a good general converter is registered for the desired class 1928 4) See if a good general converter is known for the current matrix. 1929 5) Use a really basic converter. 1930 */ 1931 1932 /* 1) See if a specialized converter is known to the current DM and the desired class */ 1933 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 1934 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 1935 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 1936 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 1937 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 1938 ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr); 1939 if (conv) goto foundconv; 1940 1941 /* 2) See if a specialized converter is known to the desired DM class. */ 1942 ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr); 1943 ierr = DMSetType(B, newtype);CHKERRQ(ierr); 1944 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 1945 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 1946 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 1947 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 1948 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 1949 ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr); 1950 if (conv) { 1951 ierr = DMDestroy(&B);CHKERRQ(ierr); 1952 goto foundconv; 1953 } 1954 1955 #if 0 1956 /* 3) See if a good general converter is registered for the desired class */ 1957 conv = B->ops->convertfrom; 1958 ierr = DMDestroy(&B);CHKERRQ(ierr); 1959 if (conv) goto foundconv; 1960 1961 /* 4) See if a good general converter is known for the current matrix */ 1962 if (dm->ops->convert) { 1963 conv = dm->ops->convert; 1964 } 1965 if (conv) goto foundconv; 1966 #endif 1967 1968 /* 5) Use a really basic converter. */ 1969 SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype); 1970 1971 foundconv: 1972 ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 1973 ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr); 1974 ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 1975 } 1976 ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr); 1977 PetscFunctionReturn(0); 1978 } 1979 1980 /*--------------------------------------------------------------------------------------------------------------------*/ 1981 1982 #undef __FUNCT__ 1983 #define __FUNCT__ "DMRegister" 1984 /*@C 1985 DMRegister - See DMRegisterDynamic() 1986 1987 Level: advanced 1988 @*/ 1989 PetscErrorCode DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM)) 1990 { 1991 char fullname[PETSC_MAX_PATH_LEN]; 1992 PetscErrorCode ierr; 1993 1994 PetscFunctionBegin; 1995 ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr); 1996 ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr); 1997 ierr = PetscStrcat(fullname, name);CHKERRQ(ierr); 1998 ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr); 1999 PetscFunctionReturn(0); 2000 } 2001 2002 2003 /*--------------------------------------------------------------------------------------------------------------------*/ 2004 #undef __FUNCT__ 2005 #define __FUNCT__ "DMRegisterDestroy" 2006 /*@C 2007 DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic(). 2008 2009 Not Collective 2010 2011 Level: advanced 2012 2013 .keywords: DM, register, destroy 2014 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic() 2015 @*/ 2016 PetscErrorCode DMRegisterDestroy(void) 2017 { 2018 PetscErrorCode ierr; 2019 2020 PetscFunctionBegin; 2021 ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr); 2022 DMRegisterAllCalled = PETSC_FALSE; 2023 PetscFunctionReturn(0); 2024 } 2025 2026 #if defined(PETSC_HAVE_MATLAB_ENGINE) 2027 #include <mex.h> 2028 2029 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext; 2030 2031 #undef __FUNCT__ 2032 #define __FUNCT__ "DMComputeFunction_Matlab" 2033 /* 2034 DMComputeFunction_Matlab - Calls the function that has been set with 2035 DMSetFunctionMatlab(). 2036 2037 For linear problems x is null 2038 2039 .seealso: DMSetFunction(), DMGetFunction() 2040 */ 2041 PetscErrorCode DMComputeFunction_Matlab(DM dm,Vec x,Vec y) 2042 { 2043 PetscErrorCode ierr; 2044 DMMatlabContext *sctx; 2045 int nlhs = 1,nrhs = 4; 2046 mxArray *plhs[1],*prhs[4]; 2047 long long int lx = 0,ly = 0,ls = 0; 2048 2049 PetscFunctionBegin; 2050 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2051 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 2052 PetscCheckSameComm(dm,1,y,3); 2053 2054 /* call Matlab function in ctx with arguments x and y */ 2055 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2056 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 2057 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2058 ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr); 2059 prhs[0] = mxCreateDoubleScalar((double)ls); 2060 prhs[1] = mxCreateDoubleScalar((double)lx); 2061 prhs[2] = mxCreateDoubleScalar((double)ly); 2062 prhs[3] = mxCreateString(sctx->funcname); 2063 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr); 2064 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2065 mxDestroyArray(prhs[0]); 2066 mxDestroyArray(prhs[1]); 2067 mxDestroyArray(prhs[2]); 2068 mxDestroyArray(prhs[3]); 2069 mxDestroyArray(plhs[0]); 2070 PetscFunctionReturn(0); 2071 } 2072 2073 2074 #undef __FUNCT__ 2075 #define __FUNCT__ "DMSetFunctionMatlab" 2076 /* 2077 DMSetFunctionMatlab - Sets the function evaluation routine 2078 2079 */ 2080 PetscErrorCode DMSetFunctionMatlab(DM dm,const char *func) 2081 { 2082 PetscErrorCode ierr; 2083 DMMatlabContext *sctx; 2084 2085 PetscFunctionBegin; 2086 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2087 /* currently sctx is memory bleed */ 2088 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2089 if (!sctx) { 2090 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 2091 } 2092 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2093 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 2094 ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr); 2095 PetscFunctionReturn(0); 2096 } 2097 2098 #undef __FUNCT__ 2099 #define __FUNCT__ "DMComputeJacobian_Matlab" 2100 /* 2101 DMComputeJacobian_Matlab - Calls the function that has been set with 2102 DMSetJacobianMatlab(). 2103 2104 For linear problems x is null 2105 2106 .seealso: DMSetFunction(), DMGetFunction() 2107 */ 2108 PetscErrorCode DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str) 2109 { 2110 PetscErrorCode ierr; 2111 DMMatlabContext *sctx; 2112 int nlhs = 2,nrhs = 5; 2113 mxArray *plhs[2],*prhs[5]; 2114 long long int lx = 0,lA = 0,lB = 0,ls = 0; 2115 2116 PetscFunctionBegin; 2117 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2118 PetscValidHeaderSpecific(A,MAT_CLASSID,3); 2119 2120 /* call MATLAB function in ctx with arguments x, A, and B */ 2121 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2122 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 2123 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2124 ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr); 2125 ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr); 2126 prhs[0] = mxCreateDoubleScalar((double)ls); 2127 prhs[1] = mxCreateDoubleScalar((double)lx); 2128 prhs[2] = mxCreateDoubleScalar((double)lA); 2129 prhs[3] = mxCreateDoubleScalar((double)lB); 2130 prhs[4] = mxCreateString(sctx->jacname); 2131 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr); 2132 *str = (MatStructure) mxGetScalar(plhs[0]); 2133 ierr = (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr); 2134 mxDestroyArray(prhs[0]); 2135 mxDestroyArray(prhs[1]); 2136 mxDestroyArray(prhs[2]); 2137 mxDestroyArray(prhs[3]); 2138 mxDestroyArray(prhs[4]); 2139 mxDestroyArray(plhs[0]); 2140 mxDestroyArray(plhs[1]); 2141 PetscFunctionReturn(0); 2142 } 2143 2144 2145 #undef __FUNCT__ 2146 #define __FUNCT__ "DMSetJacobianMatlab" 2147 /* 2148 DMSetJacobianMatlab - Sets the Jacobian function evaluation routine 2149 2150 */ 2151 PetscErrorCode DMSetJacobianMatlab(DM dm,const char *func) 2152 { 2153 PetscErrorCode ierr; 2154 DMMatlabContext *sctx; 2155 2156 PetscFunctionBegin; 2157 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2158 /* currently sctx is memory bleed */ 2159 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2160 if (!sctx) { 2161 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 2162 } 2163 ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr); 2164 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 2165 ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr); 2166 PetscFunctionReturn(0); 2167 } 2168 #endif 2169 2170 #undef __FUNCT__ 2171 #define __FUNCT__ "DMLoad" 2172 /*@C 2173 DMLoad - Loads a DM that has been stored in binary or HDF5 format 2174 with DMView(). 2175 2176 Collective on PetscViewer 2177 2178 Input Parameters: 2179 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or 2180 some related function before a call to DMLoad(). 2181 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 2182 HDF5 file viewer, obtained from PetscViewerHDF5Open() 2183 2184 Level: intermediate 2185 2186 Notes: 2187 Defaults to the DM DA. 2188 2189 Notes for advanced users: 2190 Most users should not need to know the details of the binary storage 2191 format, since DMLoad() and DMView() completely hide these details. 2192 But for anyone who's interested, the standard binary matrix storage 2193 format is 2194 .vb 2195 has not yet been determined 2196 .ve 2197 2198 In addition, PETSc automatically does the byte swapping for 2199 machines that store the bytes reversed, e.g. DEC alpha, freebsd, 2200 linux, Windows and the paragon; thus if you write your own binary 2201 read/write routines you have to swap the bytes; see PetscBinaryRead() 2202 and PetscBinaryWrite() to see how this may be done. 2203 2204 Concepts: vector^loading from file 2205 2206 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad() 2207 @*/ 2208 PetscErrorCode DMLoad(DM newdm, PetscViewer viewer) 2209 { 2210 PetscErrorCode ierr; 2211 2212 PetscFunctionBegin; 2213 PetscValidHeaderSpecific(newdm,DM_CLASSID,1); 2214 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2215 2216 if (!((PetscObject)newdm)->type_name) { 2217 ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr); 2218 } 2219 ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr); 2220 PetscFunctionReturn(0); 2221 } 2222 2223 /******************************** FEM Support **********************************/ 2224 2225 #undef __FUNCT__ 2226 #define __FUNCT__ "DMPrintCellVector" 2227 PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[]) { 2228 PetscInt f; 2229 PetscErrorCode ierr; 2230 2231 PetscFunctionBegin; 2232 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 2233 for(f = 0; f < len; ++f) { 2234 ierr = PetscPrintf(PETSC_COMM_SELF, " | %G |\n", PetscRealPart(x[f]));CHKERRQ(ierr); 2235 } 2236 PetscFunctionReturn(0); 2237 } 2238 2239 #undef __FUNCT__ 2240 #define __FUNCT__ "DMPrintCellMatrix" 2241 PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[]) { 2242 PetscInt f, g; 2243 PetscErrorCode ierr; 2244 2245 PetscFunctionBegin; 2246 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 2247 for(f = 0; f < rows; ++f) { 2248 ierr = PetscPrintf(PETSC_COMM_SELF, " |");CHKERRQ(ierr); 2249 for(g = 0; g < cols; ++g) { 2250 ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5G", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr); 2251 } 2252 ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr); 2253 } 2254 PetscFunctionReturn(0); 2255 } 2256 2257