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