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