1 #include <petscsnes.h> 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 #undef __FUNCT__ 712 #define __FUNCT__ "DMCreateFieldIS" 713 /*@C 714 DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field 715 716 Not collective 717 718 Input Parameter: 719 . dm - the DM object 720 721 Output Parameters: 722 + numFields - The number of fields 723 . names - The name for each field 724 - fields - The global indices for each field 725 726 Level: intermediate 727 728 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 729 @*/ 730 PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, const char ***names, IS **fields) 731 { 732 PetscErrorCode ierr; 733 734 PetscFunctionBegin; 735 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 736 if (numFields) {PetscValidPointer(numFields,2);} 737 if (names) {PetscValidPointer(names,3);} 738 if (fields) {PetscValidPointer(fields,4);} 739 ierr = (*dm->ops->createfieldis)(dm, numFields, names, fields);CHKERRQ(ierr); 740 PetscFunctionReturn(0); 741 } 742 743 744 #undef __FUNCT__ 745 #define __FUNCT__ "DMRefine" 746 /*@ 747 DMRefine - Refines a DM object 748 749 Collective on DM 750 751 Input Parameter: 752 + dm - the DM object 753 - comm - the communicator to contain the new DM object (or PETSC_NULL) 754 755 Output Parameter: 756 . dmf - the refined DM 757 758 Level: developer 759 760 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 761 762 @*/ 763 PetscErrorCode DMRefine(DM dm,MPI_Comm comm,DM *dmf) 764 { 765 PetscErrorCode ierr; 766 767 PetscFunctionBegin; 768 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 769 ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr); 770 if (*dmf) { 771 (*dmf)->ops->initialguess = dm->ops->initialguess; 772 (*dmf)->ops->function = dm->ops->function; 773 (*dmf)->ops->functionj = dm->ops->functionj; 774 if (dm->ops->jacobian != DMComputeJacobianDefault) { 775 (*dmf)->ops->jacobian = dm->ops->jacobian; 776 } 777 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr); 778 (*dmf)->ctx = dm->ctx; 779 (*dmf)->levelup = dm->levelup + 1; 780 } 781 PetscFunctionReturn(0); 782 } 783 784 #undef __FUNCT__ 785 #define __FUNCT__ "DMGetRefineLevel" 786 /*@ 787 DMGetRefineLevel - Get's the number of refinements that have generated this DM. 788 789 Not Collective 790 791 Input Parameter: 792 . dm - the DM object 793 794 Output Parameter: 795 . level - number of refinements 796 797 Level: developer 798 799 .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 800 801 @*/ 802 PetscErrorCode DMGetRefineLevel(DM dm,PetscInt *level) 803 { 804 PetscFunctionBegin; 805 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 806 *level = dm->levelup; 807 PetscFunctionReturn(0); 808 } 809 810 #undef __FUNCT__ 811 #define __FUNCT__ "DMGlobalToLocalBegin" 812 /*@ 813 DMGlobalToLocalBegin - Begins updating local vectors from global vector 814 815 Neighbor-wise Collective on DM 816 817 Input Parameters: 818 + dm - the DM object 819 . g - the global vector 820 . mode - INSERT_VALUES or ADD_VALUES 821 - l - the local vector 822 823 824 Level: beginner 825 826 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 827 828 @*/ 829 PetscErrorCode DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) 830 { 831 PetscErrorCode ierr; 832 833 PetscFunctionBegin; 834 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 835 ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 836 PetscFunctionReturn(0); 837 } 838 839 #undef __FUNCT__ 840 #define __FUNCT__ "DMGlobalToLocalEnd" 841 /*@ 842 DMGlobalToLocalEnd - Ends updating local vectors from global vector 843 844 Neighbor-wise Collective on DM 845 846 Input Parameters: 847 + dm - the DM object 848 . g - the global vector 849 . mode - INSERT_VALUES or ADD_VALUES 850 - l - the local vector 851 852 853 Level: beginner 854 855 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 856 857 @*/ 858 PetscErrorCode DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) 859 { 860 PetscErrorCode ierr; 861 862 PetscFunctionBegin; 863 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 864 ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 865 PetscFunctionReturn(0); 866 } 867 868 #undef __FUNCT__ 869 #define __FUNCT__ "DMLocalToGlobalBegin" 870 /*@ 871 DMLocalToGlobalBegin - updates global vectors from local vectors 872 873 Neighbor-wise Collective on DM 874 875 Input Parameters: 876 + dm - the DM object 877 . l - the local vector 878 . 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 879 base point. 880 - - the global vector 881 882 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 883 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 884 global array to the final global array with VecAXPY(). 885 886 Level: beginner 887 888 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin() 889 890 @*/ 891 PetscErrorCode DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g) 892 { 893 PetscErrorCode ierr; 894 895 PetscFunctionBegin; 896 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 897 ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 898 PetscFunctionReturn(0); 899 } 900 901 #undef __FUNCT__ 902 #define __FUNCT__ "DMLocalToGlobalEnd" 903 /*@ 904 DMLocalToGlobalEnd - updates global vectors from local vectors 905 906 Neighbor-wise Collective on DM 907 908 Input Parameters: 909 + dm - the DM object 910 . l - the local vector 911 . mode - INSERT_VALUES or ADD_VALUES 912 - g - the global vector 913 914 915 Level: beginner 916 917 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd() 918 919 @*/ 920 PetscErrorCode DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g) 921 { 922 PetscErrorCode ierr; 923 924 PetscFunctionBegin; 925 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 926 ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 927 PetscFunctionReturn(0); 928 } 929 930 #undef __FUNCT__ 931 #define __FUNCT__ "DMComputeJacobianDefault" 932 /*@ 933 DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided 934 935 Collective on DM 936 937 Input Parameter: 938 + dm - the DM object 939 . x - location to compute Jacobian at; may be ignored for linear problems 940 . A - matrix that defines the operator for the linear solve 941 - B - the matrix used to construct the preconditioner 942 943 Level: developer 944 945 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 946 DMSetFunction() 947 948 @*/ 949 PetscErrorCode DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 950 { 951 PetscErrorCode ierr; 952 953 PetscFunctionBegin; 954 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 955 *stflag = SAME_NONZERO_PATTERN; 956 ierr = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr); 957 if (A != B) { 958 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 959 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 960 } 961 PetscFunctionReturn(0); 962 } 963 964 #undef __FUNCT__ 965 #define __FUNCT__ "DMCoarsen" 966 /*@ 967 DMCoarsen - Coarsens a DM object 968 969 Collective on DM 970 971 Input Parameter: 972 + dm - the DM object 973 - comm - the communicator to contain the new DM object (or PETSC_NULL) 974 975 Output Parameter: 976 . dmc - the coarsened DM 977 978 Level: developer 979 980 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 981 982 @*/ 983 PetscErrorCode DMCoarsen(DM dm, MPI_Comm comm, DM *dmc) 984 { 985 PetscErrorCode ierr; 986 DMCoarsenHookLink link; 987 988 PetscFunctionBegin; 989 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 990 ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr); 991 (*dmc)->ops->initialguess = dm->ops->initialguess; 992 (*dmc)->ops->function = dm->ops->function; 993 (*dmc)->ops->functionj = dm->ops->functionj; 994 if (dm->ops->jacobian != DMComputeJacobianDefault) { 995 (*dmc)->ops->jacobian = dm->ops->jacobian; 996 } 997 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr); 998 (*dmc)->ctx = dm->ctx; 999 (*dmc)->leveldown = dm->leveldown + 1; 1000 for (link=dm->coarsenhook; link; link=link->next) { 1001 if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);} 1002 } 1003 PetscFunctionReturn(0); 1004 } 1005 1006 #undef __FUNCT__ 1007 #define __FUNCT__ "DMCoarsenHookAdd" 1008 /*@ 1009 DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid 1010 1011 Logically Collective 1012 1013 Input Arguments: 1014 + fine - nonlinear solver context on which to run a hook when restricting to a coarser level 1015 . coarsenhook - function to run when setting up a coarser level 1016 . restricthook - function to run to update data on coarser levels (once per SNESSolve()) 1017 - ctx - [optional] user-defined context for provide data for the hooks (may be PETSC_NULL) 1018 1019 Calling sequence of coarsenhook: 1020 $ coarsenhook(DM fine,DM coarse,void *ctx); 1021 1022 + fine - fine level DM 1023 . coarse - coarse level DM to restrict problem to 1024 - ctx - optional user-defined function context 1025 1026 Calling sequence for restricthook: 1027 $ restricthook(DM fine,Mat mrestrict,Mat inject,DM coarse,void *ctx) 1028 1029 + fine - fine level DM 1030 . mrestrict - matrix restricting a fine-level solution to the coarse grid 1031 . inject - matrix restricting by applying the transpose of injection 1032 . coarse - coarse level DM to update 1033 - ctx - optional user-defined function context 1034 1035 Level: advanced 1036 1037 Notes: 1038 This function is only needed if auxiliary data needs to be set up on coarse grids. 1039 1040 If this function is called multiple times, the hooks will be run in the order they are added. 1041 1042 The 1043 1044 In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to 1045 extract the finest level information from its context (instead of from the SNES). 1046 1047 .seealso: SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 1048 @*/ 1049 PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx) 1050 { 1051 PetscErrorCode ierr; 1052 DMCoarsenHookLink link,*p; 1053 1054 PetscFunctionBegin; 1055 PetscValidHeaderSpecific(fine,DM_CLASSID,1); 1056 for (p=&fine->coarsenhook; *p; *p=(*p)->next) {} /* Scan to the end of the current list of hooks */ 1057 ierr = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr); 1058 link->coarsenhook = coarsenhook; 1059 link->restricthook = restricthook; 1060 link->ctx = ctx; 1061 *p = link; 1062 PetscFunctionReturn(0); 1063 } 1064 1065 #undef __FUNCT__ 1066 #define __FUNCT__ "DMRestrict" 1067 /*@ 1068 DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd() 1069 1070 Collective if any hooks are 1071 1072 Input Arguments: 1073 + fine - finer DM to use as a base 1074 . restrct - restriction matrix, apply using MatRestrict() 1075 . inject - injection matrix, also use MatRestrict() 1076 - coarse - coarer DM to update 1077 1078 Level: developer 1079 1080 .seealso: DMCoarsenHookAdd(), MatRestrict() 1081 @*/ 1082 PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse) 1083 { 1084 PetscErrorCode ierr; 1085 DMCoarsenHookLink link; 1086 1087 PetscFunctionBegin; 1088 for (link=fine->coarsenhook; link; link=link->next) { 1089 if (link->restricthook) {ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);} 1090 } 1091 PetscFunctionReturn(0); 1092 } 1093 1094 #undef __FUNCT__ 1095 #define __FUNCT__ "DMGetCoarsenLevel" 1096 /*@ 1097 DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM. 1098 1099 Not Collective 1100 1101 Input Parameter: 1102 . dm - the DM object 1103 1104 Output Parameter: 1105 . level - number of coarsenings 1106 1107 Level: developer 1108 1109 .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1110 1111 @*/ 1112 PetscErrorCode DMGetCoarsenLevel(DM dm,PetscInt *level) 1113 { 1114 PetscFunctionBegin; 1115 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1116 *level = dm->leveldown; 1117 PetscFunctionReturn(0); 1118 } 1119 1120 1121 1122 #undef __FUNCT__ 1123 #define __FUNCT__ "DMRefineHierarchy" 1124 /*@C 1125 DMRefineHierarchy - Refines a DM object, all levels at once 1126 1127 Collective on DM 1128 1129 Input Parameter: 1130 + dm - the DM object 1131 - nlevels - the number of levels of refinement 1132 1133 Output Parameter: 1134 . dmf - the refined DM hierarchy 1135 1136 Level: developer 1137 1138 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1139 1140 @*/ 1141 PetscErrorCode DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[]) 1142 { 1143 PetscErrorCode ierr; 1144 1145 PetscFunctionBegin; 1146 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1147 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 1148 if (nlevels == 0) PetscFunctionReturn(0); 1149 if (dm->ops->refinehierarchy) { 1150 ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr); 1151 } else if (dm->ops->refine) { 1152 PetscInt i; 1153 1154 ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr); 1155 for (i=1; i<nlevels; i++) { 1156 ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr); 1157 } 1158 } else { 1159 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet"); 1160 } 1161 PetscFunctionReturn(0); 1162 } 1163 1164 #undef __FUNCT__ 1165 #define __FUNCT__ "DMCoarsenHierarchy" 1166 /*@C 1167 DMCoarsenHierarchy - Coarsens a DM object, all levels at once 1168 1169 Collective on DM 1170 1171 Input Parameter: 1172 + dm - the DM object 1173 - nlevels - the number of levels of coarsening 1174 1175 Output Parameter: 1176 . dmc - the coarsened DM hierarchy 1177 1178 Level: developer 1179 1180 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1181 1182 @*/ 1183 PetscErrorCode DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[]) 1184 { 1185 PetscErrorCode ierr; 1186 1187 PetscFunctionBegin; 1188 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1189 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 1190 if (nlevels == 0) PetscFunctionReturn(0); 1191 PetscValidPointer(dmc,3); 1192 if (dm->ops->coarsenhierarchy) { 1193 ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr); 1194 } else if (dm->ops->coarsen) { 1195 PetscInt i; 1196 1197 ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr); 1198 for (i=1; i<nlevels; i++) { 1199 ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr); 1200 } 1201 } else { 1202 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet"); 1203 } 1204 PetscFunctionReturn(0); 1205 } 1206 1207 #undef __FUNCT__ 1208 #define __FUNCT__ "DMCreateAggregates" 1209 /*@ 1210 DMCreateAggregates - Gets the aggregates that map between 1211 grids associated with two DMs. 1212 1213 Collective on DM 1214 1215 Input Parameters: 1216 + dmc - the coarse grid DM 1217 - dmf - the fine grid DM 1218 1219 Output Parameters: 1220 . rest - the restriction matrix (transpose of the projection matrix) 1221 1222 Level: intermediate 1223 1224 .keywords: interpolation, restriction, multigrid 1225 1226 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation() 1227 @*/ 1228 PetscErrorCode DMCreateAggregates(DM dmc, DM dmf, Mat *rest) 1229 { 1230 PetscErrorCode ierr; 1231 1232 PetscFunctionBegin; 1233 PetscValidHeaderSpecific(dmc,DM_CLASSID,1); 1234 PetscValidHeaderSpecific(dmf,DM_CLASSID,2); 1235 ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr); 1236 PetscFunctionReturn(0); 1237 } 1238 1239 #undef __FUNCT__ 1240 #define __FUNCT__ "DMSetApplicationContextDestroy" 1241 /*@C 1242 DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed 1243 1244 Not Collective 1245 1246 Input Parameters: 1247 + dm - the DM object 1248 - destroy - the destroy function 1249 1250 Level: intermediate 1251 1252 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1253 1254 @*/ 1255 PetscErrorCode DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**)) 1256 { 1257 PetscFunctionBegin; 1258 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1259 dm->ctxdestroy = destroy; 1260 PetscFunctionReturn(0); 1261 } 1262 1263 #undef __FUNCT__ 1264 #define __FUNCT__ "DMSetApplicationContext" 1265 /*@ 1266 DMSetApplicationContext - Set a user context into a DM object 1267 1268 Not Collective 1269 1270 Input Parameters: 1271 + dm - the DM object 1272 - ctx - the user context 1273 1274 Level: intermediate 1275 1276 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1277 1278 @*/ 1279 PetscErrorCode DMSetApplicationContext(DM dm,void *ctx) 1280 { 1281 PetscFunctionBegin; 1282 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1283 dm->ctx = ctx; 1284 PetscFunctionReturn(0); 1285 } 1286 1287 #undef __FUNCT__ 1288 #define __FUNCT__ "DMGetApplicationContext" 1289 /*@ 1290 DMGetApplicationContext - Gets a user context from a DM object 1291 1292 Not Collective 1293 1294 Input Parameter: 1295 . dm - the DM object 1296 1297 Output Parameter: 1298 . ctx - the user context 1299 1300 Level: intermediate 1301 1302 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1303 1304 @*/ 1305 PetscErrorCode DMGetApplicationContext(DM dm,void *ctx) 1306 { 1307 PetscFunctionBegin; 1308 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1309 *(void**)ctx = dm->ctx; 1310 PetscFunctionReturn(0); 1311 } 1312 1313 #undef __FUNCT__ 1314 #define __FUNCT__ "DMSetInitialGuess" 1315 /*@C 1316 DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers 1317 1318 Logically Collective on DM 1319 1320 Input Parameter: 1321 + dm - the DM object to destroy 1322 - f - the function to compute the initial guess 1323 1324 Level: intermediate 1325 1326 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1327 1328 @*/ 1329 PetscErrorCode DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec)) 1330 { 1331 PetscFunctionBegin; 1332 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1333 dm->ops->initialguess = f; 1334 PetscFunctionReturn(0); 1335 } 1336 1337 #undef __FUNCT__ 1338 #define __FUNCT__ "DMSetFunction" 1339 /*@C 1340 DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES 1341 1342 Logically Collective on DM 1343 1344 Input Parameter: 1345 + dm - the DM object 1346 - f - the function to compute (use PETSC_NULL to cancel a previous function that was set) 1347 1348 Level: intermediate 1349 1350 Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian 1351 computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian. 1352 1353 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1354 DMSetJacobian() 1355 1356 @*/ 1357 PetscErrorCode DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 1358 { 1359 PetscFunctionBegin; 1360 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1361 dm->ops->function = f; 1362 if (f) { 1363 dm->ops->functionj = f; 1364 } 1365 PetscFunctionReturn(0); 1366 } 1367 1368 #undef __FUNCT__ 1369 #define __FUNCT__ "DMSetJacobian" 1370 /*@C 1371 DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES 1372 1373 Logically Collective on DM 1374 1375 Input Parameter: 1376 + dm - the DM object to destroy 1377 - f - the function to compute the matrix entries 1378 1379 Level: intermediate 1380 1381 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1382 DMSetFunction() 1383 1384 @*/ 1385 PetscErrorCode DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*)) 1386 { 1387 PetscFunctionBegin; 1388 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1389 dm->ops->jacobian = f; 1390 PetscFunctionReturn(0); 1391 } 1392 1393 #undef __FUNCT__ 1394 #define __FUNCT__ "DMSetVariableBounds" 1395 /*@C 1396 DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI. 1397 1398 Logically Collective on DM 1399 1400 Input Parameter: 1401 + dm - the DM object 1402 - f - the function that computes variable bounds used by SNESVI (use PETSC_NULL to cancel a previous function that was set) 1403 1404 Level: intermediate 1405 1406 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1407 DMSetJacobian() 1408 1409 @*/ 1410 PetscErrorCode DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 1411 { 1412 PetscFunctionBegin; 1413 dm->ops->computevariablebounds = f; 1414 PetscFunctionReturn(0); 1415 } 1416 1417 #undef __FUNCT__ 1418 #define __FUNCT__ "DMHasVariableBounds" 1419 /*@ 1420 DMHasVariableBounds - does the DM object have a variable bounds function? 1421 1422 Not Collective 1423 1424 Input Parameter: 1425 . dm - the DM object to destroy 1426 1427 Output Parameter: 1428 . flg - PETSC_TRUE if the variable bounds function exists 1429 1430 Level: developer 1431 1432 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1433 1434 @*/ 1435 PetscErrorCode DMHasVariableBounds(DM dm,PetscBool *flg) 1436 { 1437 PetscFunctionBegin; 1438 *flg = (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE; 1439 PetscFunctionReturn(0); 1440 } 1441 1442 #undef __FUNCT__ 1443 #define __FUNCT__ "DMComputeVariableBounds" 1444 /*@C 1445 DMComputeVariableBounds - compute variable bounds used by SNESVI. 1446 1447 Logically Collective on DM 1448 1449 Input Parameters: 1450 + dm - the DM object to destroy 1451 - x - current solution at which the bounds are computed 1452 1453 Output parameters: 1454 + xl - lower bound 1455 - xu - upper bound 1456 1457 Level: intermediate 1458 1459 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1460 DMSetFunction(), DMSetVariableBounds() 1461 1462 @*/ 1463 PetscErrorCode DMComputeVariableBounds(DM dm, Vec xl, Vec xu) 1464 { 1465 PetscErrorCode ierr; 1466 PetscFunctionBegin; 1467 PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 1468 PetscValidHeaderSpecific(xu,VEC_CLASSID,2); 1469 if(dm->ops->computevariablebounds) { 1470 ierr = (*dm->ops->computevariablebounds)(dm, xl,xu); CHKERRQ(ierr); 1471 } 1472 else { 1473 ierr = VecSet(xl,SNES_VI_NINF); CHKERRQ(ierr); 1474 ierr = VecSet(xu,SNES_VI_INF); CHKERRQ(ierr); 1475 } 1476 PetscFunctionReturn(0); 1477 } 1478 1479 #undef __FUNCT__ 1480 #define __FUNCT__ "DMComputeInitialGuess" 1481 /*@ 1482 DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers 1483 1484 Collective on DM 1485 1486 Input Parameter: 1487 + dm - the DM object to destroy 1488 - x - the vector to hold the initial guess values 1489 1490 Level: developer 1491 1492 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat() 1493 1494 @*/ 1495 PetscErrorCode DMComputeInitialGuess(DM dm,Vec x) 1496 { 1497 PetscErrorCode ierr; 1498 PetscFunctionBegin; 1499 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1500 if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()"); 1501 ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr); 1502 PetscFunctionReturn(0); 1503 } 1504 1505 #undef __FUNCT__ 1506 #define __FUNCT__ "DMHasInitialGuess" 1507 /*@ 1508 DMHasInitialGuess - does the DM object have an initial guess function 1509 1510 Not Collective 1511 1512 Input Parameter: 1513 . dm - the DM object to destroy 1514 1515 Output Parameter: 1516 . flg - PETSC_TRUE if function exists 1517 1518 Level: developer 1519 1520 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1521 1522 @*/ 1523 PetscErrorCode DMHasInitialGuess(DM dm,PetscBool *flg) 1524 { 1525 PetscFunctionBegin; 1526 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1527 *flg = (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE; 1528 PetscFunctionReturn(0); 1529 } 1530 1531 #undef __FUNCT__ 1532 #define __FUNCT__ "DMHasFunction" 1533 /*@ 1534 DMHasFunction - does the DM object have a function 1535 1536 Not Collective 1537 1538 Input Parameter: 1539 . dm - the DM object to destroy 1540 1541 Output Parameter: 1542 . flg - PETSC_TRUE if function exists 1543 1544 Level: developer 1545 1546 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1547 1548 @*/ 1549 PetscErrorCode DMHasFunction(DM dm,PetscBool *flg) 1550 { 1551 PetscFunctionBegin; 1552 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1553 *flg = (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE; 1554 PetscFunctionReturn(0); 1555 } 1556 1557 #undef __FUNCT__ 1558 #define __FUNCT__ "DMHasJacobian" 1559 /*@ 1560 DMHasJacobian - does the DM object have a matrix function 1561 1562 Not Collective 1563 1564 Input Parameter: 1565 . dm - the DM object to destroy 1566 1567 Output Parameter: 1568 . flg - PETSC_TRUE if function exists 1569 1570 Level: developer 1571 1572 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1573 1574 @*/ 1575 PetscErrorCode DMHasJacobian(DM dm,PetscBool *flg) 1576 { 1577 PetscFunctionBegin; 1578 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1579 *flg = (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE; 1580 PetscFunctionReturn(0); 1581 } 1582 1583 #undef __FUNCT__ 1584 #define __FUNCT__ "DMSetVec" 1585 /*@C 1586 DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear. 1587 1588 Collective on DM 1589 1590 Input Parameter: 1591 + dm - the DM object 1592 - x - location to compute residual and Jacobian, if PETSC_NULL is passed to those routines; will be PETSC_NULL for linear problems. 1593 1594 Level: developer 1595 1596 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1597 DMSetFunction(), DMSetJacobian(), DMSetVariableBounds() 1598 1599 @*/ 1600 PetscErrorCode DMSetVec(DM dm,Vec x) 1601 { 1602 PetscErrorCode ierr; 1603 PetscFunctionBegin; 1604 if (x) { 1605 if (!dm->x) { 1606 ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 1607 } 1608 ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 1609 } 1610 else if(dm->x) { 1611 ierr = VecDestroy(&dm->x); CHKERRQ(ierr); 1612 } 1613 PetscFunctionReturn(0); 1614 } 1615 1616 1617 #undef __FUNCT__ 1618 #define __FUNCT__ "DMComputeFunction" 1619 /*@ 1620 DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES 1621 1622 Collective on DM 1623 1624 Input Parameter: 1625 + dm - the DM object to destroy 1626 . x - the location where the function is evaluationed, may be ignored for linear problems 1627 - b - the vector to hold the right hand side entries 1628 1629 Level: developer 1630 1631 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1632 DMSetJacobian() 1633 1634 @*/ 1635 PetscErrorCode DMComputeFunction(DM dm,Vec x,Vec b) 1636 { 1637 PetscErrorCode ierr; 1638 PetscFunctionBegin; 1639 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1640 if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()"); 1641 PetscStackPush("DM user function"); 1642 ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr); 1643 PetscStackPop; 1644 PetscFunctionReturn(0); 1645 } 1646 1647 1648 1649 #undef __FUNCT__ 1650 #define __FUNCT__ "DMComputeJacobian" 1651 /*@ 1652 DMComputeJacobian - compute the matrix entries for the solver 1653 1654 Collective on DM 1655 1656 Input Parameter: 1657 + dm - the DM object 1658 . x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM 1659 . A - matrix that defines the operator for the linear solve 1660 - B - the matrix used to construct the preconditioner 1661 1662 Level: developer 1663 1664 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1665 DMSetFunction() 1666 1667 @*/ 1668 PetscErrorCode DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 1669 { 1670 PetscErrorCode ierr; 1671 1672 PetscFunctionBegin; 1673 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1674 if (!dm->ops->jacobian) { 1675 ISColoring coloring; 1676 MatFDColoring fd; 1677 1678 ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr); 1679 ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr); 1680 ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 1681 ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr); 1682 ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr); 1683 ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr); 1684 1685 dm->fd = fd; 1686 dm->ops->jacobian = DMComputeJacobianDefault; 1687 1688 /* don't know why this is needed */ 1689 ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 1690 } 1691 if (!x) x = dm->x; 1692 ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr); 1693 1694 /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */ 1695 if (x) { 1696 if (!dm->x) { 1697 ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 1698 } 1699 ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 1700 } 1701 if (A != B) { 1702 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1703 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1704 } 1705 PetscFunctionReturn(0); 1706 } 1707 1708 1709 PetscFList DMList = PETSC_NULL; 1710 PetscBool DMRegisterAllCalled = PETSC_FALSE; 1711 1712 #undef __FUNCT__ 1713 #define __FUNCT__ "DMSetType" 1714 /*@C 1715 DMSetType - Builds a DM, for a particular DM implementation. 1716 1717 Collective on DM 1718 1719 Input Parameters: 1720 + dm - The DM object 1721 - method - The name of the DM type 1722 1723 Options Database Key: 1724 . -dm_type <type> - Sets the DM type; use -help for a list of available types 1725 1726 Notes: 1727 See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D). 1728 1729 Level: intermediate 1730 1731 .keywords: DM, set, type 1732 .seealso: DMGetType(), DMCreate() 1733 @*/ 1734 PetscErrorCode DMSetType(DM dm, const DMType method) 1735 { 1736 PetscErrorCode (*r)(DM); 1737 PetscBool match; 1738 PetscErrorCode ierr; 1739 1740 PetscFunctionBegin; 1741 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1742 ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr); 1743 if (match) PetscFunctionReturn(0); 1744 1745 if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1746 ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 1747 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method); 1748 1749 if (dm->ops->destroy) { 1750 ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr); 1751 dm->ops->destroy = PETSC_NULL; 1752 } 1753 ierr = (*r)(dm);CHKERRQ(ierr); 1754 ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr); 1755 PetscFunctionReturn(0); 1756 } 1757 1758 #undef __FUNCT__ 1759 #define __FUNCT__ "DMGetType" 1760 /*@C 1761 DMGetType - Gets the DM type name (as a string) from the DM. 1762 1763 Not Collective 1764 1765 Input Parameter: 1766 . dm - The DM 1767 1768 Output Parameter: 1769 . type - The DM type name 1770 1771 Level: intermediate 1772 1773 .keywords: DM, get, type, name 1774 .seealso: DMSetType(), DMCreate() 1775 @*/ 1776 PetscErrorCode DMGetType(DM dm, const DMType *type) 1777 { 1778 PetscErrorCode ierr; 1779 1780 PetscFunctionBegin; 1781 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1782 PetscValidCharPointer(type,2); 1783 if (!DMRegisterAllCalled) { 1784 ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1785 } 1786 *type = ((PetscObject)dm)->type_name; 1787 PetscFunctionReturn(0); 1788 } 1789 1790 #undef __FUNCT__ 1791 #define __FUNCT__ "DMConvert" 1792 /*@C 1793 DMConvert - Converts a DM to another DM, either of the same or different type. 1794 1795 Collective on DM 1796 1797 Input Parameters: 1798 + dm - the DM 1799 - newtype - new DM type (use "same" for the same type) 1800 1801 Output Parameter: 1802 . M - pointer to new DM 1803 1804 Notes: 1805 Cannot be used to convert a sequential DM to parallel or parallel to sequential, 1806 the MPI communicator of the generated DM is always the same as the communicator 1807 of the input DM. 1808 1809 Level: intermediate 1810 1811 .seealso: DMCreate() 1812 @*/ 1813 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M) 1814 { 1815 DM B; 1816 char convname[256]; 1817 PetscBool sametype, issame; 1818 PetscErrorCode ierr; 1819 1820 PetscFunctionBegin; 1821 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1822 PetscValidType(dm,1); 1823 PetscValidPointer(M,3); 1824 ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr); 1825 ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr); 1826 { 1827 PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL; 1828 1829 /* 1830 Order of precedence: 1831 1) See if a specialized converter is known to the current DM. 1832 2) See if a specialized converter is known to the desired DM class. 1833 3) See if a good general converter is registered for the desired class 1834 4) See if a good general converter is known for the current matrix. 1835 5) Use a really basic converter. 1836 */ 1837 1838 /* 1) See if a specialized converter is known to the current DM and the desired class */ 1839 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 1840 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 1841 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 1842 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 1843 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 1844 ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr); 1845 if (conv) goto foundconv; 1846 1847 /* 2) See if a specialized converter is known to the desired DM class. */ 1848 ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr); 1849 ierr = DMSetType(B, newtype);CHKERRQ(ierr); 1850 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 1851 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 1852 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 1853 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 1854 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 1855 ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr); 1856 if (conv) { 1857 ierr = DMDestroy(&B);CHKERRQ(ierr); 1858 goto foundconv; 1859 } 1860 1861 #if 0 1862 /* 3) See if a good general converter is registered for the desired class */ 1863 conv = B->ops->convertfrom; 1864 ierr = DMDestroy(&B);CHKERRQ(ierr); 1865 if (conv) goto foundconv; 1866 1867 /* 4) See if a good general converter is known for the current matrix */ 1868 if (dm->ops->convert) { 1869 conv = dm->ops->convert; 1870 } 1871 if (conv) goto foundconv; 1872 #endif 1873 1874 /* 5) Use a really basic converter. */ 1875 SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype); 1876 1877 foundconv: 1878 ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 1879 ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr); 1880 ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 1881 } 1882 ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr); 1883 PetscFunctionReturn(0); 1884 } 1885 1886 /*--------------------------------------------------------------------------------------------------------------------*/ 1887 1888 #undef __FUNCT__ 1889 #define __FUNCT__ "DMRegister" 1890 /*@C 1891 DMRegister - See DMRegisterDynamic() 1892 1893 Level: advanced 1894 @*/ 1895 PetscErrorCode DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM)) 1896 { 1897 char fullname[PETSC_MAX_PATH_LEN]; 1898 PetscErrorCode ierr; 1899 1900 PetscFunctionBegin; 1901 ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr); 1902 ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr); 1903 ierr = PetscStrcat(fullname, name);CHKERRQ(ierr); 1904 ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr); 1905 PetscFunctionReturn(0); 1906 } 1907 1908 1909 /*--------------------------------------------------------------------------------------------------------------------*/ 1910 #undef __FUNCT__ 1911 #define __FUNCT__ "DMRegisterDestroy" 1912 /*@C 1913 DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic(). 1914 1915 Not Collective 1916 1917 Level: advanced 1918 1919 .keywords: DM, register, destroy 1920 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic() 1921 @*/ 1922 PetscErrorCode DMRegisterDestroy(void) 1923 { 1924 PetscErrorCode ierr; 1925 1926 PetscFunctionBegin; 1927 ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr); 1928 DMRegisterAllCalled = PETSC_FALSE; 1929 PetscFunctionReturn(0); 1930 } 1931 1932 #if defined(PETSC_HAVE_MATLAB_ENGINE) 1933 #include <mex.h> 1934 1935 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext; 1936 1937 #undef __FUNCT__ 1938 #define __FUNCT__ "DMComputeFunction_Matlab" 1939 /* 1940 DMComputeFunction_Matlab - Calls the function that has been set with 1941 DMSetFunctionMatlab(). 1942 1943 For linear problems x is null 1944 1945 .seealso: DMSetFunction(), DMGetFunction() 1946 */ 1947 PetscErrorCode DMComputeFunction_Matlab(DM dm,Vec x,Vec y) 1948 { 1949 PetscErrorCode ierr; 1950 DMMatlabContext *sctx; 1951 int nlhs = 1,nrhs = 4; 1952 mxArray *plhs[1],*prhs[4]; 1953 long long int lx = 0,ly = 0,ls = 0; 1954 1955 PetscFunctionBegin; 1956 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1957 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 1958 PetscCheckSameComm(dm,1,y,3); 1959 1960 /* call Matlab function in ctx with arguments x and y */ 1961 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1962 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 1963 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 1964 ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr); 1965 prhs[0] = mxCreateDoubleScalar((double)ls); 1966 prhs[1] = mxCreateDoubleScalar((double)lx); 1967 prhs[2] = mxCreateDoubleScalar((double)ly); 1968 prhs[3] = mxCreateString(sctx->funcname); 1969 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr); 1970 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 1971 mxDestroyArray(prhs[0]); 1972 mxDestroyArray(prhs[1]); 1973 mxDestroyArray(prhs[2]); 1974 mxDestroyArray(prhs[3]); 1975 mxDestroyArray(plhs[0]); 1976 PetscFunctionReturn(0); 1977 } 1978 1979 1980 #undef __FUNCT__ 1981 #define __FUNCT__ "DMSetFunctionMatlab" 1982 /* 1983 DMSetFunctionMatlab - Sets the function evaluation routine 1984 1985 */ 1986 PetscErrorCode DMSetFunctionMatlab(DM dm,const char *func) 1987 { 1988 PetscErrorCode ierr; 1989 DMMatlabContext *sctx; 1990 1991 PetscFunctionBegin; 1992 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1993 /* currently sctx is memory bleed */ 1994 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1995 if (!sctx) { 1996 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 1997 } 1998 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 1999 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 2000 ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr); 2001 PetscFunctionReturn(0); 2002 } 2003 2004 #undef __FUNCT__ 2005 #define __FUNCT__ "DMComputeJacobian_Matlab" 2006 /* 2007 DMComputeJacobian_Matlab - Calls the function that has been set with 2008 DMSetJacobianMatlab(). 2009 2010 For linear problems x is null 2011 2012 .seealso: DMSetFunction(), DMGetFunction() 2013 */ 2014 PetscErrorCode DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str) 2015 { 2016 PetscErrorCode ierr; 2017 DMMatlabContext *sctx; 2018 int nlhs = 2,nrhs = 5; 2019 mxArray *plhs[2],*prhs[5]; 2020 long long int lx = 0,lA = 0,lB = 0,ls = 0; 2021 2022 PetscFunctionBegin; 2023 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2024 PetscValidHeaderSpecific(A,MAT_CLASSID,3); 2025 2026 /* call MATLAB function in ctx with arguments x, A, and B */ 2027 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2028 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 2029 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2030 ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr); 2031 ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr); 2032 prhs[0] = mxCreateDoubleScalar((double)ls); 2033 prhs[1] = mxCreateDoubleScalar((double)lx); 2034 prhs[2] = mxCreateDoubleScalar((double)lA); 2035 prhs[3] = mxCreateDoubleScalar((double)lB); 2036 prhs[4] = mxCreateString(sctx->jacname); 2037 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr); 2038 *str = (MatStructure) mxGetScalar(plhs[0]); 2039 ierr = (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr); 2040 mxDestroyArray(prhs[0]); 2041 mxDestroyArray(prhs[1]); 2042 mxDestroyArray(prhs[2]); 2043 mxDestroyArray(prhs[3]); 2044 mxDestroyArray(prhs[4]); 2045 mxDestroyArray(plhs[0]); 2046 mxDestroyArray(plhs[1]); 2047 PetscFunctionReturn(0); 2048 } 2049 2050 2051 #undef __FUNCT__ 2052 #define __FUNCT__ "DMSetJacobianMatlab" 2053 /* 2054 DMSetJacobianMatlab - Sets the Jacobian function evaluation routine 2055 2056 */ 2057 PetscErrorCode DMSetJacobianMatlab(DM dm,const char *func) 2058 { 2059 PetscErrorCode ierr; 2060 DMMatlabContext *sctx; 2061 2062 PetscFunctionBegin; 2063 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2064 /* currently sctx is memory bleed */ 2065 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2066 if (!sctx) { 2067 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 2068 } 2069 ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr); 2070 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 2071 ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr); 2072 PetscFunctionReturn(0); 2073 } 2074 #endif 2075 2076 #undef __FUNCT__ 2077 #define __FUNCT__ "DMLoad" 2078 /*@C 2079 DMLoad - Loads a DM that has been stored in binary or HDF5 format 2080 with DMView(). 2081 2082 Collective on PetscViewer 2083 2084 Input Parameters: 2085 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or 2086 some related function before a call to DMLoad(). 2087 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 2088 HDF5 file viewer, obtained from PetscViewerHDF5Open() 2089 2090 Level: intermediate 2091 2092 Notes: 2093 Defaults to the DM DA. 2094 2095 Notes for advanced users: 2096 Most users should not need to know the details of the binary storage 2097 format, since DMLoad() and DMView() completely hide these details. 2098 But for anyone who's interested, the standard binary matrix storage 2099 format is 2100 .vb 2101 has not yet been determined 2102 .ve 2103 2104 In addition, PETSc automatically does the byte swapping for 2105 machines that store the bytes reversed, e.g. DEC alpha, freebsd, 2106 linux, Windows and the paragon; thus if you write your own binary 2107 read/write routines you have to swap the bytes; see PetscBinaryRead() 2108 and PetscBinaryWrite() to see how this may be done. 2109 2110 Concepts: vector^loading from file 2111 2112 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad() 2113 @*/ 2114 PetscErrorCode DMLoad(DM newdm, PetscViewer viewer) 2115 { 2116 PetscErrorCode ierr; 2117 2118 PetscFunctionBegin; 2119 PetscValidHeaderSpecific(newdm,DM_CLASSID,1); 2120 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2121 2122 if (!((PetscObject)newdm)->type_name) { 2123 ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr); 2124 } 2125 ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr); 2126 PetscFunctionReturn(0); 2127 } 2128 2129 /******************************** FEM Support **********************************/ 2130 2131 #undef __FUNCT__ 2132 #define __FUNCT__ "DMPrintCellVector" 2133 PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[]) { 2134 PetscInt f; 2135 PetscErrorCode ierr; 2136 2137 PetscFunctionBegin; 2138 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 2139 for(f = 0; f < len; ++f) { 2140 ierr = PetscPrintf(PETSC_COMM_SELF, " | %G |\n", PetscRealPart(x[f]));CHKERRQ(ierr); 2141 } 2142 PetscFunctionReturn(0); 2143 } 2144 2145 #undef __FUNCT__ 2146 #define __FUNCT__ "DMPrintCellMatrix" 2147 PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[]) { 2148 PetscInt f, g; 2149 PetscErrorCode ierr; 2150 2151 PetscFunctionBegin; 2152 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 2153 for(f = 0; f < rows; ++f) { 2154 ierr = PetscPrintf(PETSC_COMM_SELF, " |");CHKERRQ(ierr); 2155 for(g = 0; g < cols; ++g) { 2156 ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5G", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr); 2157 } 2158 ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr); 2159 } 2160 PetscFunctionReturn(0); 2161 } 2162