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