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