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 v->lf = PETSC_NULL; 52 v->lj = PETSC_NULL; 53 ierr = PetscSFCreate(comm, &v->sf);CHKERRQ(ierr); 54 ierr = PetscSFCreate(comm, &v->defaultSF);CHKERRQ(ierr); 55 v->defaultSection = PETSC_NULL; 56 v->defaultGlobalSection = PETSC_NULL; 57 58 *dm = v; 59 PetscFunctionReturn(0); 60 } 61 62 63 #undef __FUNCT__ 64 #define __FUNCT__ "DMSetVecType" 65 /*@C 66 DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() 67 68 Logically Collective on DMDA 69 70 Input Parameter: 71 + da - initial distributed array 72 . ctype - the vector type, currently either VECSTANDARD or VECCUSP 73 74 Options Database: 75 . -dm_vec_type ctype 76 77 Level: intermediate 78 79 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType 80 @*/ 81 PetscErrorCode DMSetVecType(DM da,const VecType ctype) 82 { 83 PetscErrorCode ierr; 84 85 PetscFunctionBegin; 86 PetscValidHeaderSpecific(da,DM_CLASSID,1); 87 ierr = PetscFree(da->vectype);CHKERRQ(ierr); 88 ierr = PetscStrallocpy(ctype,&da->vectype);CHKERRQ(ierr); 89 PetscFunctionReturn(0); 90 } 91 92 #undef __FUNCT__ 93 #define __FUNCT__ "DMSetMatType" 94 /*@C 95 DMSetMatType - Sets the type of matrix created with DMCreateMatrix() 96 97 Logically Collective on DM 98 99 Input Parameter: 100 + dm - the DM context 101 . ctype - the matrix type 102 103 Options Database: 104 . -dm_mat_type ctype 105 106 Level: intermediate 107 108 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMCreateMatrix(), DMSetMatrixPreallocateOnly(), MatType 109 @*/ 110 PetscErrorCode DMSetMatType(DM dm,const MatType ctype) 111 { 112 PetscErrorCode ierr; 113 PetscFunctionBegin; 114 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 115 ierr = PetscFree(dm->mattype);CHKERRQ(ierr); 116 ierr = PetscStrallocpy(ctype,&dm->mattype);CHKERRQ(ierr); 117 PetscFunctionReturn(0); 118 } 119 120 #undef __FUNCT__ 121 #define __FUNCT__ "DMSetOptionsPrefix" 122 /*@C 123 DMSetOptionsPrefix - Sets the prefix used for searching for all 124 DMDA options in the database. 125 126 Logically Collective on DMDA 127 128 Input Parameter: 129 + da - the DMDA context 130 - prefix - the prefix to prepend to all option names 131 132 Notes: 133 A hyphen (-) must NOT be given at the beginning of the prefix name. 134 The first character of all runtime options is AUTOMATICALLY the hyphen. 135 136 Level: advanced 137 138 .keywords: DMDA, set, options, prefix, database 139 140 .seealso: DMSetFromOptions() 141 @*/ 142 PetscErrorCode DMSetOptionsPrefix(DM dm,const char prefix[]) 143 { 144 PetscErrorCode ierr; 145 146 PetscFunctionBegin; 147 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 148 ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr); 149 PetscFunctionReturn(0); 150 } 151 152 #undef __FUNCT__ 153 #define __FUNCT__ "DMDestroy" 154 /*@ 155 DMDestroy - Destroys a vector packer or DMDA. 156 157 Collective on DM 158 159 Input Parameter: 160 . dm - the DM object to destroy 161 162 Level: developer 163 164 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 165 166 @*/ 167 PetscErrorCode DMDestroy(DM *dm) 168 { 169 PetscInt i, cnt = 0; 170 DMNamedVecLink nlink,nnext; 171 PetscErrorCode ierr; 172 173 PetscFunctionBegin; 174 if (!*dm) PetscFunctionReturn(0); 175 PetscValidHeaderSpecific((*dm),DM_CLASSID,1); 176 177 /* count all the circular references of DM and its contained Vecs */ 178 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 179 if ((*dm)->localin[i]) {cnt++;} 180 if ((*dm)->globalin[i]) {cnt++;} 181 } 182 for (nlink=(*dm)->namedglobal; nlink; nlink=nlink->next) cnt++; 183 if ((*dm)->x) { 184 PetscObject obj; 185 ierr = PetscObjectQuery((PetscObject)(*dm)->x,"DM",&obj);CHKERRQ(ierr); 186 if (obj == (PetscObject)*dm) cnt++; 187 } 188 189 if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);} 190 /* 191 Need this test because the dm references the vectors that 192 reference the dm, so destroying the dm calls destroy on the 193 vectors that cause another destroy on the dm 194 */ 195 if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0); 196 ((PetscObject) (*dm))->refct = 0; 197 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 198 if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()"); 199 ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr); 200 } 201 for (nlink=(*dm)->namedglobal; nlink; nlink=nnext) { /* Destroy the named vectors */ 202 nnext = nlink->next; 203 if (nlink->status != DMVEC_STATUS_IN) SETERRQ1(((PetscObject)*dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM still has Vec named '%s' checked out",nlink->name); 204 ierr = PetscFree(nlink->name);CHKERRQ(ierr); 205 ierr = VecDestroy(&nlink->X);CHKERRQ(ierr); 206 ierr = PetscFree(nlink);CHKERRQ(ierr); 207 } 208 (*dm)->namedglobal = PETSC_NULL; 209 210 /* Destroy the list of hooks */ 211 { 212 DMCoarsenHookLink link,next; 213 for (link=(*dm)->coarsenhook; link; link=next) { 214 next = link->next; 215 ierr = PetscFree(link);CHKERRQ(ierr); 216 } 217 (*dm)->coarsenhook = PETSC_NULL; 218 } 219 { 220 DMRefineHookLink link,next; 221 for (link=(*dm)->refinehook; link; link=next) { 222 next = link->next; 223 ierr = PetscFree(link);CHKERRQ(ierr); 224 } 225 (*dm)->refinehook = PETSC_NULL; 226 } 227 228 if ((*dm)->ctx && (*dm)->ctxdestroy) { 229 ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr); 230 } 231 ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr); 232 ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr); 233 ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr); 234 ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr); 235 ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr); 236 ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr); 237 ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr); 238 ierr = PetscFree((*dm)->workArray);CHKERRQ(ierr); 239 240 ierr = PetscSectionDestroy(&(*dm)->defaultSection);CHKERRQ(ierr); 241 ierr = PetscSectionDestroy(&(*dm)->defaultGlobalSection);CHKERRQ(ierr); 242 ierr = PetscSFDestroy(&(*dm)->sf);CHKERRQ(ierr); 243 ierr = PetscSFDestroy(&(*dm)->defaultSF);CHKERRQ(ierr); 244 /* if memory was published with AMS then destroy it */ 245 ierr = PetscObjectDepublish(*dm);CHKERRQ(ierr); 246 247 ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr); 248 /* We do not destroy (*dm)->data here so that we can reference count backend objects */ 249 ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr); 250 PetscFunctionReturn(0); 251 } 252 253 #undef __FUNCT__ 254 #define __FUNCT__ "DMSetUp" 255 /*@ 256 DMSetUp - sets up the data structures inside a DM object 257 258 Collective on DM 259 260 Input Parameter: 261 . dm - the DM object to setup 262 263 Level: developer 264 265 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 266 267 @*/ 268 PetscErrorCode DMSetUp(DM dm) 269 { 270 PetscErrorCode ierr; 271 272 PetscFunctionBegin; 273 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 274 if (dm->setupcalled) PetscFunctionReturn(0); 275 if (dm->ops->setup) { 276 ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr); 277 } 278 dm->setupcalled = PETSC_TRUE; 279 PetscFunctionReturn(0); 280 } 281 282 #undef __FUNCT__ 283 #define __FUNCT__ "DMSetFromOptions" 284 /*@ 285 DMSetFromOptions - sets parameters in a DM from the options database 286 287 Collective on DM 288 289 Input Parameter: 290 . dm - the DM object to set options for 291 292 Options Database: 293 + -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros 294 . -dm_vec_type <type> type of vector to create inside DM 295 . -dm_mat_type <type> type of matrix to create inside DM 296 - -dm_coloring_type <global or ghosted> 297 298 Level: developer 299 300 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 301 302 @*/ 303 PetscErrorCode DMSetFromOptions(DM dm) 304 { 305 PetscBool flg1 = PETSC_FALSE,flg2 = PETSC_FALSE,flg3 = PETSC_FALSE,flg4 = PETSC_FALSE,flg; 306 PetscErrorCode ierr; 307 char typeName[256] = MATAIJ; 308 309 PetscFunctionBegin; 310 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 311 ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr); 312 ierr = PetscOptionsBool("-dm_view", "Information on DM", "DMView", flg1, &flg1, PETSC_NULL);CHKERRQ(ierr); 313 ierr = PetscOptionsBool("-dm_view_detail", "Exhaustive mesh description", "DMView", flg2, &flg2, PETSC_NULL);CHKERRQ(ierr); 314 ierr = PetscOptionsBool("-dm_view_vtk", "Output mesh in VTK format", "DMView", flg3, &flg3, PETSC_NULL);CHKERRQ(ierr); 315 ierr = PetscOptionsBool("-dm_view_latex", "Output mesh in LaTeX TikZ format", "DMView", flg4, &flg4, PETSC_NULL);CHKERRQ(ierr); 316 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); 317 ierr = PetscOptionsList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr); 318 if (flg) { 319 ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr); 320 } 321 ierr = PetscOptionsList("-dm_mat_type","Matrix type used for created matrices","DMSetMatType",MatList,dm->mattype?dm->mattype:typeName,typeName,sizeof typeName,&flg);CHKERRQ(ierr); 322 if (flg) { 323 ierr = DMSetMatType(dm,typeName);CHKERRQ(ierr); 324 } 325 ierr = PetscOptionsEnum("-dm_is_coloring_type","Global or local coloring of Jacobian","ISColoringType",ISColoringTypes,(PetscEnum)dm->coloringtype,(PetscEnum*)&dm->coloringtype,PETSC_NULL);CHKERRQ(ierr); 326 if (dm->ops->setfromoptions) { 327 ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr); 328 } 329 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 330 ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr); 331 ierr = PetscOptionsEnd();CHKERRQ(ierr); 332 if (flg1) { 333 ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 334 } 335 if (flg2) { 336 PetscViewer viewer; 337 338 ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr); 339 ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 340 ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr); 341 ierr = DMView(dm, viewer);CHKERRQ(ierr); 342 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 343 } 344 if (flg3) { 345 PetscViewer viewer; 346 347 ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr); 348 ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 349 ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr); 350 ierr = PetscViewerFileSetName(viewer, "mesh.vtk");CHKERRQ(ierr); 351 ierr = DMView(dm, viewer);CHKERRQ(ierr); 352 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 353 } 354 if (flg4) { 355 PetscViewer viewer; 356 357 ierr = PetscViewerCreate(((PetscObject) dm)->comm, &viewer);CHKERRQ(ierr); 358 ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 359 ierr = PetscViewerSetFormat(viewer, PETSC_VIEWER_ASCII_LATEX);CHKERRQ(ierr); 360 ierr = PetscViewerFileSetName(viewer, "mesh.tex");CHKERRQ(ierr); 361 ierr = DMView(dm, viewer);CHKERRQ(ierr); 362 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 363 } 364 PetscFunctionReturn(0); 365 } 366 367 #undef __FUNCT__ 368 #define __FUNCT__ "DMView" 369 /*@C 370 DMView - Views a vector packer or DMDA. 371 372 Collective on DM 373 374 Input Parameter: 375 + dm - the DM object to view 376 - v - the viewer 377 378 Level: developer 379 380 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 381 382 @*/ 383 PetscErrorCode DMView(DM dm,PetscViewer v) 384 { 385 PetscErrorCode ierr; 386 387 PetscFunctionBegin; 388 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 389 if (!v) { 390 ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr); 391 } 392 if (dm->ops->view) { 393 ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr); 394 } 395 PetscFunctionReturn(0); 396 } 397 398 #undef __FUNCT__ 399 #define __FUNCT__ "DMCreateGlobalVector" 400 /*@ 401 DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object 402 403 Collective on DM 404 405 Input Parameter: 406 . dm - the DM object 407 408 Output Parameter: 409 . vec - the global vector 410 411 Level: beginner 412 413 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 414 415 @*/ 416 PetscErrorCode DMCreateGlobalVector(DM dm,Vec *vec) 417 { 418 PetscErrorCode ierr; 419 420 PetscFunctionBegin; 421 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 422 if (dm->defaultSection) { 423 PetscSection gSection; 424 PetscInt localSize; 425 426 ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr); 427 ierr = PetscSectionGetConstrainedStorageSize(dm->defaultGlobalSection, &localSize);CHKERRQ(ierr); 428 ierr = VecCreate(((PetscObject) dm)->comm, vec);CHKERRQ(ierr); 429 ierr = VecSetSizes(*vec, localSize, PETSC_DETERMINE);CHKERRQ(ierr); 430 /* ierr = VecSetType(*vec, dm->vectype);CHKERRQ(ierr); */ 431 ierr = VecSetFromOptions(*vec);CHKERRQ(ierr); 432 ierr = PetscObjectCompose((PetscObject) *vec, "DM", (PetscObject) dm);CHKERRQ(ierr); 433 /* ierr = VecSetLocalToGlobalMapping(*vec, dm->ltogmap);CHKERRQ(ierr); */ 434 /* ierr = VecSetLocalToGlobalMappingBlock(*vec, dm->ltogmapb);CHKERRQ(ierr); */ 435 /* ierr = VecSetOperation(*vec, VECOP_DUPLICATE, (void(*)(void)) VecDuplicate_MPI_DM);CHKERRQ(ierr); */ 436 } else { 437 ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr); 438 } 439 PetscFunctionReturn(0); 440 } 441 442 #undef __FUNCT__ 443 #define __FUNCT__ "DMCreateLocalVector" 444 /*@ 445 DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object 446 447 Not Collective 448 449 Input Parameter: 450 . dm - the DM object 451 452 Output Parameter: 453 . vec - the local vector 454 455 Level: beginner 456 457 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 458 459 @*/ 460 PetscErrorCode DMCreateLocalVector(DM dm,Vec *vec) 461 { 462 PetscErrorCode ierr; 463 464 PetscFunctionBegin; 465 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 466 if (dm->defaultSection) { 467 PetscInt localSize; 468 469 ierr = PetscSectionGetStorageSize(dm->defaultSection, &localSize);CHKERRQ(ierr); 470 ierr = VecCreate(PETSC_COMM_SELF, vec);CHKERRQ(ierr); 471 ierr = VecSetSizes(*vec, localSize, localSize);CHKERRQ(ierr); 472 ierr = VecSetFromOptions(*vec);CHKERRQ(ierr); 473 ierr = PetscObjectCompose((PetscObject) *vec, "DM", (PetscObject) dm);CHKERRQ(ierr); 474 } else { 475 ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr); 476 } 477 PetscFunctionReturn(0); 478 } 479 480 #undef __FUNCT__ 481 #define __FUNCT__ "DMGetLocalToGlobalMapping" 482 /*@ 483 DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM. 484 485 Collective on DM 486 487 Input Parameter: 488 . dm - the DM that provides the mapping 489 490 Output Parameter: 491 . ltog - the mapping 492 493 Level: intermediate 494 495 Notes: 496 This mapping can then be used by VecSetLocalToGlobalMapping() or 497 MatSetLocalToGlobalMapping(). 498 499 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock() 500 @*/ 501 PetscErrorCode DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog) 502 { 503 PetscErrorCode ierr; 504 505 PetscFunctionBegin; 506 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 507 PetscValidPointer(ltog,2); 508 if (!dm->ltogmap) { 509 PetscSection section, sectionGlobal; 510 511 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 512 if (section) { 513 PetscInt *ltog; 514 PetscInt pStart, pEnd, size, p, l; 515 516 ierr = DMGetDefaultGlobalSection(dm, §ionGlobal);CHKERRQ(ierr); 517 ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 518 ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr); 519 ierr = PetscMalloc(size * sizeof(PetscInt), <og);CHKERRQ(ierr); /* We want the local+overlap size */ 520 for(p = pStart, l = 0; p < pEnd; ++p) { 521 PetscInt dof, off, c; 522 523 /* Should probably use constrained dofs */ 524 ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 525 ierr = PetscSectionGetOffset(sectionGlobal, p, &off);CHKERRQ(ierr); 526 for(c = 0; c < dof; ++c, ++l) { 527 ltog[l] = off+c; 528 } 529 } 530 ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF, size, ltog, PETSC_OWN_POINTER, &dm->ltogmap);CHKERRQ(ierr); 531 ierr = PetscLogObjectParent(dm, dm->ltogmap);CHKERRQ(ierr); 532 } else { 533 if (!dm->ops->createlocaltoglobalmapping) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping"); 534 ierr = (*dm->ops->createlocaltoglobalmapping)(dm);CHKERRQ(ierr); 535 } 536 } 537 *ltog = dm->ltogmap; 538 PetscFunctionReturn(0); 539 } 540 541 #undef __FUNCT__ 542 #define __FUNCT__ "DMGetLocalToGlobalMappingBlock" 543 /*@ 544 DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM. 545 546 Collective on DM 547 548 Input Parameter: 549 . da - the distributed array that provides the mapping 550 551 Output Parameter: 552 . ltog - the block mapping 553 554 Level: intermediate 555 556 Notes: 557 This mapping can then be used by VecSetLocalToGlobalMappingBlock() or 558 MatSetLocalToGlobalMappingBlock(). 559 560 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize() 561 @*/ 562 PetscErrorCode DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog) 563 { 564 PetscErrorCode ierr; 565 566 PetscFunctionBegin; 567 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 568 PetscValidPointer(ltog,2); 569 if (!dm->ltogmapb) { 570 PetscInt bs; 571 ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr); 572 if (bs > 1) { 573 if (!dm->ops->createlocaltoglobalmappingblock) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock"); 574 ierr = (*dm->ops->createlocaltoglobalmappingblock)(dm);CHKERRQ(ierr); 575 } else { 576 ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr); 577 ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr); 578 } 579 } 580 *ltog = dm->ltogmapb; 581 PetscFunctionReturn(0); 582 } 583 584 #undef __FUNCT__ 585 #define __FUNCT__ "DMGetBlockSize" 586 /*@ 587 DMGetBlockSize - Gets the inherent block size associated with a DM 588 589 Not Collective 590 591 Input Parameter: 592 . dm - the DM with block structure 593 594 Output Parameter: 595 . bs - the block size, 1 implies no exploitable block structure 596 597 Level: intermediate 598 599 .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock() 600 @*/ 601 PetscErrorCode DMGetBlockSize(DM dm,PetscInt *bs) 602 { 603 PetscFunctionBegin; 604 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 605 PetscValidPointer(bs,2); 606 if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet"); 607 *bs = dm->bs; 608 PetscFunctionReturn(0); 609 } 610 611 #undef __FUNCT__ 612 #define __FUNCT__ "DMCreateInterpolation" 613 /*@ 614 DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects 615 616 Collective on DM 617 618 Input Parameter: 619 + dm1 - the DM object 620 - dm2 - the second, finer DM object 621 622 Output Parameter: 623 + mat - the interpolation 624 - vec - the scaling (optional) 625 626 Level: developer 627 628 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 629 DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation. 630 631 For DMDA objects you can use this interpolation (more precisely the interpolation from the DMDAGetCoordinateDA()) to interpolate the mesh coordinate vectors 632 EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic. 633 634 635 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen() 636 637 @*/ 638 PetscErrorCode DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec) 639 { 640 PetscErrorCode ierr; 641 642 PetscFunctionBegin; 643 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 644 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 645 ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr); 646 PetscFunctionReturn(0); 647 } 648 649 #undef __FUNCT__ 650 #define __FUNCT__ "DMCreateInjection" 651 /*@ 652 DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects 653 654 Collective on DM 655 656 Input Parameter: 657 + dm1 - the DM object 658 - dm2 - the second, finer DM object 659 660 Output Parameter: 661 . ctx - the injection 662 663 Level: developer 664 665 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 666 DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection. 667 668 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation() 669 670 @*/ 671 PetscErrorCode DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx) 672 { 673 PetscErrorCode ierr; 674 675 PetscFunctionBegin; 676 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 677 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 678 ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr); 679 PetscFunctionReturn(0); 680 } 681 682 #undef __FUNCT__ 683 #define __FUNCT__ "DMCreateColoring" 684 /*@C 685 DMCreateColoring - Gets coloring for a DMDA or DMComposite 686 687 Collective on DM 688 689 Input Parameter: 690 + dm - the DM object 691 . ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL 692 - matype - either MATAIJ or MATBAIJ 693 694 Output Parameter: 695 . coloring - the coloring 696 697 Level: developer 698 699 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix() 700 701 @*/ 702 PetscErrorCode DMCreateColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring) 703 { 704 PetscErrorCode ierr; 705 706 PetscFunctionBegin; 707 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 708 if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet"); 709 ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr); 710 PetscFunctionReturn(0); 711 } 712 713 #undef __FUNCT__ 714 #define __FUNCT__ "DMCreateMatrix" 715 /*@C 716 DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite 717 718 Collective on DM 719 720 Input Parameter: 721 + dm - the DM object 722 - mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or 723 any type which inherits from one of these (such as MATAIJ) 724 725 Output Parameter: 726 . mat - the empty Jacobian 727 728 Level: beginner 729 730 Notes: This properly preallocates the number of nonzeros in the sparse matrix so you 731 do not need to do it yourself. 732 733 By default it also sets the nonzero structure and puts in the zero entries. To prevent setting 734 the nonzero pattern call DMDASetMatPreallocateOnly() 735 736 For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used 737 internally by PETSc. 738 739 For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires 740 the indices for the global numbering for DMDAs which is complicated. 741 742 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 743 744 @*/ 745 PetscErrorCode DMCreateMatrix(DM dm,const MatType mtype,Mat *mat) 746 { 747 PetscErrorCode ierr; 748 749 PetscFunctionBegin; 750 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 751 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 752 ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr); 753 #endif 754 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 755 PetscValidPointer(mat,3); 756 if (dm->mattype) { 757 ierr = (*dm->ops->creatematrix)(dm,dm->mattype,mat);CHKERRQ(ierr); 758 } else { 759 ierr = (*dm->ops->creatematrix)(dm,mtype,mat);CHKERRQ(ierr); 760 } 761 PetscFunctionReturn(0); 762 } 763 764 #undef __FUNCT__ 765 #define __FUNCT__ "DMSetMatrixPreallocateOnly" 766 /*@ 767 DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly 768 preallocated but the nonzero structure and zero values will not be set. 769 770 Logically Collective on DMDA 771 772 Input Parameter: 773 + dm - the DM 774 - only - PETSC_TRUE if only want preallocation 775 776 Level: developer 777 .seealso DMCreateMatrix() 778 @*/ 779 PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only) 780 { 781 PetscFunctionBegin; 782 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 783 dm->prealloc_only = only; 784 PetscFunctionReturn(0); 785 } 786 787 #undef __FUNCT__ 788 #define __FUNCT__ "DMGetWorkArray" 789 /*@C 790 DMGetWorkArray - Gets a work array guaranteed to be at least the input size 791 792 Not Collective 793 794 Input Parameters: 795 + dm - the DM object 796 - size - The minium size 797 798 Output Parameter: 799 . array - the work array 800 801 Level: developer 802 803 .seealso DMDestroy(), DMCreate() 804 @*/ 805 PetscErrorCode DMGetWorkArray(DM dm,PetscInt size,PetscScalar **array) 806 { 807 PetscErrorCode ierr; 808 809 PetscFunctionBegin; 810 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 811 PetscValidPointer(array,3); 812 if (size > dm->workSize) { 813 dm->workSize = size; 814 ierr = PetscFree(dm->workArray);CHKERRQ(ierr); 815 ierr = PetscMalloc(dm->workSize * sizeof(PetscScalar), &dm->workArray);CHKERRQ(ierr); 816 } 817 *array = dm->workArray; 818 PetscFunctionReturn(0); 819 } 820 821 822 #undef __FUNCT__ 823 #define __FUNCT__ "DMSetNullSpaceConstructor" 824 PetscErrorCode DMSetNullSpaceConstructor(DM dm, PetscInt field, PetscErrorCode (*nullsp)(DM dm, PetscInt field, MatNullSpace *nullSpace)) 825 { 826 PetscFunctionBegin; 827 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 828 if (field >= 10) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %d >= 10 fields", field); 829 dm->nullspaceConstructors[field] = nullsp; 830 PetscFunctionReturn(0); 831 } 832 833 #undef __FUNCT__ 834 #define __FUNCT__ "DMCreateFieldIS" 835 /*@C 836 DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field 837 838 Not collective 839 840 Input Parameter: 841 . dm - the DM object 842 843 Output Parameters: 844 + numFields - The number of fields (or PETSC_NULL if not requested) 845 . fieldNames - The name for each field (or PETSC_NULL if not requested) 846 - fields - The global indices for each field (or PETSC_NULL if not requested) 847 848 Level: intermediate 849 850 Notes: 851 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 852 PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with 853 PetscFree(). 854 855 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 856 @*/ 857 PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields) 858 { 859 PetscSection section, sectionGlobal; 860 PetscErrorCode ierr; 861 862 PetscFunctionBegin; 863 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 864 if (numFields) { 865 PetscValidPointer(numFields,2); 866 *numFields = 0; 867 } 868 if (fieldNames) { 869 PetscValidPointer(fieldNames,3); 870 *fieldNames = PETSC_NULL; 871 } 872 if (fields) { 873 PetscValidPointer(fields,4); 874 *fields = PETSC_NULL; 875 } 876 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 877 if (section) { 878 PetscInt *fieldSizes, **fieldIndices; 879 PetscInt nF, f, pStart, pEnd, p; 880 881 ierr = DMGetDefaultGlobalSection(dm, §ionGlobal);CHKERRQ(ierr); 882 ierr = PetscSectionGetNumFields(section, &nF);CHKERRQ(ierr); 883 ierr = PetscMalloc2(nF,PetscInt,&fieldSizes,nF,PetscInt *,&fieldIndices);CHKERRQ(ierr); 884 ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr); 885 for(f = 0; f < nF; ++f) { 886 fieldSizes[f] = 0; 887 } 888 for(p = pStart; p < pEnd; ++p) { 889 PetscInt gdof; 890 891 ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr); 892 if (gdof > 0) { 893 for(f = 0; f < nF; ++f) { 894 PetscInt fdof, fcdof; 895 896 ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr); 897 ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr); 898 fieldSizes[f] += fdof-fcdof; 899 } 900 } 901 } 902 for(f = 0; f < nF; ++f) { 903 ierr = PetscMalloc(fieldSizes[f] * sizeof(PetscInt), &fieldIndices[f]);CHKERRQ(ierr); 904 fieldSizes[f] = 0; 905 } 906 for(p = pStart; p < pEnd; ++p) { 907 PetscInt gdof, goff; 908 909 ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr); 910 if (gdof > 0) { 911 ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 912 for(f = 0; f < nF; ++f) { 913 PetscInt fdof, fcdof, fc; 914 915 ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr); 916 ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr); 917 for(fc = 0; fc < fdof-fcdof; ++fc, ++fieldSizes[f]) { 918 fieldIndices[f][fieldSizes[f]] = goff++; 919 } 920 } 921 } 922 } 923 if (numFields) {*numFields = nF;} 924 if (fieldNames) { 925 ierr = PetscMalloc(nF * sizeof(char *), fieldNames);CHKERRQ(ierr); 926 for(f = 0; f < nF; ++f) { 927 const char *fieldName; 928 929 ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr); 930 ierr = PetscStrallocpy(fieldName, (char **) &(*fieldNames)[f]);CHKERRQ(ierr); 931 } 932 } 933 if (fields) { 934 ierr = PetscMalloc(nF * sizeof(IS), fields);CHKERRQ(ierr); 935 for(f = 0; f < nF; ++f) { 936 ierr = ISCreateGeneral(((PetscObject) dm)->comm, fieldSizes[f], fieldIndices[f], PETSC_OWN_POINTER, &(*fields)[f]);CHKERRQ(ierr); 937 } 938 } 939 ierr = PetscFree2(fieldSizes,fieldIndices);CHKERRQ(ierr); 940 } else { 941 if(dm->ops->createfieldis) {ierr = (*dm->ops->createfieldis)(dm, numFields, fieldNames, fields);CHKERRQ(ierr);} 942 } 943 PetscFunctionReturn(0); 944 } 945 946 #undef __FUNCT__ 947 #define __FUNCT__ "DMCreateFieldDecompositionDM" 948 /*@C 949 DMCreateFieldDecompositionDM - creates a DM that encapsulates a decomposition of the original DM into fields. 950 951 Not Collective 952 953 Input Parameters: 954 + dm - the DM object 955 - name - the name of the field decomposition 956 957 Output Parameter: 958 . ddm - the field decomposition DM (PETSC_NULL, if no such decomposition is known) 959 960 Level: advanced 961 962 .seealso DMDestroy(), DMCreate(), DMCreateFieldDecomposition(), DMCreateDomainDecompositionDM() 963 @*/ 964 PetscErrorCode DMCreateFieldDecompositionDM(DM dm, const char* name, DM *ddm) 965 { 966 PetscErrorCode ierr; 967 968 PetscFunctionBegin; 969 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 970 PetscValidCharPointer(name,2); 971 PetscValidPointer(ddm,3); 972 *ddm = PETSC_NULL; 973 if(dm->ops->createfielddecompositiondm) { 974 ierr = (*dm->ops->createfielddecompositiondm)(dm,name,ddm); CHKERRQ(ierr); 975 } 976 PetscFunctionReturn(0); 977 } 978 979 #undef __FUNCT__ 980 #define __FUNCT__ "DMCreateFieldDecomposition" 981 /*@C 982 DMCreateFieldDecomposition - Returns a list of IS objects defining a decomposition of a problem into subproblems 983 corresponding to different fields: each IS contains the global indices of the dofs of the 984 corresponding field. The optional list of DMs define the DM for each subproblem. 985 Generalizes DMCreateFieldIS(). 986 987 Not collective 988 989 Input Parameter: 990 . dm - the DM object 991 992 Output Parameters: 993 + len - The number of subproblems in the field decomposition (or PETSC_NULL if not requested) 994 . namelist - The name for each field (or PETSC_NULL if not requested) 995 . islist - The global indices for each field (or PETSC_NULL if not requested) 996 - dmlist - The DMs for each field subproblem (or PETSC_NULL, if not requested; if PETSC_NULL is returned, no DMs are defined) 997 998 Level: intermediate 999 1000 Notes: 1001 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 1002 PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(), 1003 and all of the arrays should be freed with PetscFree(). 1004 1005 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS() 1006 @*/ 1007 PetscErrorCode DMCreateFieldDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 1008 { 1009 PetscErrorCode ierr; 1010 1011 PetscFunctionBegin; 1012 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1013 if (len) {PetscValidPointer(len,2); *len = 0;} 1014 if (namelist) {PetscValidPointer(namelist,3); *namelist = 0;} 1015 if (islist) {PetscValidPointer(islist,4); *islist = 0;} 1016 if (dmlist) {PetscValidPointer(dmlist,5); *dmlist = 0;} 1017 if (!dm->ops->createfielddecomposition) { 1018 PetscSection section; 1019 PetscInt numFields, f; 1020 1021 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1022 if (section) {ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);} 1023 if (section && numFields && dm->ops->createsubdm) { 1024 *len = numFields; 1025 ierr = PetscMalloc3(numFields,char*,namelist,numFields,IS,islist,numFields,DM,dmlist);CHKERRQ(ierr); 1026 for(f = 0; f < numFields; ++f) { 1027 const char *fieldName; 1028 1029 ierr = DMCreateSubDM(dm, 1, &f, &(*islist)[f], &(*dmlist)[f]);CHKERRQ(ierr); 1030 ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr); 1031 ierr = PetscStrallocpy(fieldName, (char **) &(*namelist)[f]);CHKERRQ(ierr); 1032 } 1033 } else { 1034 ierr = DMCreateFieldIS(dm, len, namelist, islist);CHKERRQ(ierr); 1035 /* By default there are no DMs associated with subproblems. */ 1036 if(dmlist) *dmlist = PETSC_NULL; 1037 } 1038 } 1039 else { 1040 ierr = (*dm->ops->createfielddecomposition)(dm,len,namelist,islist,dmlist); CHKERRQ(ierr); 1041 } 1042 PetscFunctionReturn(0); 1043 } 1044 1045 #undef __FUNCT__ 1046 #define __FUNCT__ "DMCreateSubDM" 1047 /*@C 1048 DMCreateSubDM - Returns an IS and DM encapsulating a subproblem defined by the fields passed in. 1049 The fields are defined by DMCreateFieldIS(). 1050 1051 Not collective 1052 1053 Input Parameters: 1054 + dm - the DM object 1055 . numFields - number of fields in this subproblem 1056 - len - The number of subproblems in the decomposition (or PETSC_NULL if not requested) 1057 1058 Output Parameters: 1059 . is - The global indices for the subproblem 1060 - dm - The DM for the subproblem 1061 1062 Level: intermediate 1063 1064 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS() 1065 @*/ 1066 PetscErrorCode DMCreateSubDM(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm) 1067 { 1068 PetscErrorCode ierr; 1069 1070 PetscFunctionBegin; 1071 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1072 PetscValidPointer(fields,3); 1073 if (is) {PetscValidPointer(is,4);} 1074 if (subdm) {PetscValidPointer(subdm,5);} 1075 if (dm->ops->createsubdm) { 1076 ierr = (*dm->ops->createsubdm)(dm, numFields, fields, is, subdm); CHKERRQ(ierr); 1077 } else SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "This type has no DMCreateSubDM implementation defined"); 1078 PetscFunctionReturn(0); 1079 } 1080 1081 #undef __FUNCT__ 1082 #define __FUNCT__ "DMCreateDomainDecompositionDM" 1083 /*@C 1084 DMCreateDomainDecompositionDM - creates a DM that encapsulates a decomposition of the original DM into subdomains. 1085 1086 Not Collective 1087 1088 Input Parameters: 1089 + dm - the DM object 1090 - name - the name of the subdomain decomposition 1091 1092 Output Parameter: 1093 . ddm - the subdomain decomposition DM (PETSC_NULL, if no such decomposition is known) 1094 1095 Level: advanced 1096 1097 .seealso DMDestroy(), DMCreate(), DMCreateFieldDecomposition(), DMCreateDomainDecompositionDM() 1098 @*/ 1099 PetscErrorCode DMCreateDomainDecompositionDM(DM dm, const char* name, DM *ddm) 1100 { 1101 PetscErrorCode ierr; 1102 1103 PetscFunctionBegin; 1104 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1105 PetscValidCharPointer(name,2); 1106 PetscValidPointer(ddm,3); 1107 *ddm = PETSC_NULL; 1108 if(dm->ops->createdomaindecompositiondm) { 1109 ierr = (*dm->ops->createdomaindecompositiondm)(dm,name,ddm); CHKERRQ(ierr); 1110 } 1111 PetscFunctionReturn(0); 1112 } 1113 1114 #undef __FUNCT__ 1115 #define __FUNCT__ "DMCreateDomainDecomposition" 1116 /*@C 1117 DMCreateDomainDecomposition - Returns lists of IS objects defining a decomposition of a problem into subproblems 1118 corresponding to restrictions to pairs nested subdomains: each IS contains the global 1119 indices of the dofs of the corresponding subdomains. The inner subdomains conceptually 1120 define a nonoverlapping covering, while outer subdomains can overlap. 1121 The optional list of DMs define the DM for each subproblem. 1122 1123 Not collective 1124 1125 Input Parameter: 1126 . dm - the DM object 1127 1128 Output Parameters: 1129 + len - The number of subproblems in the domain decomposition (or PETSC_NULL if not requested) 1130 . namelist - The name for each subdomain (or PETSC_NULL if not requested) 1131 . innerislist - The global indices for each inner subdomain (or PETSC_NULL, if not requested) 1132 . outerislist - The global indices for each outer subdomain (or PETSC_NULL, if not requested) 1133 - dmlist - The DMs for each subdomain subproblem (or PETSC_NULL, if not requested; if PETSC_NULL is returned, no DMs are defined) 1134 1135 Level: intermediate 1136 1137 Notes: 1138 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 1139 PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(), 1140 and all of the arrays should be freed with PetscFree(). 1141 1142 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateDomainDecompositionDM(), DMCreateFieldDecomposition() 1143 @*/ 1144 PetscErrorCode DMCreateDomainDecomposition(DM dm, PetscInt *len, char ***namelist, IS **innerislist, IS **outerislist, DM **dmlist) 1145 { 1146 PetscErrorCode ierr; 1147 1148 PetscFunctionBegin; 1149 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1150 if (len) {PetscValidPointer(len,2); *len = PETSC_NULL;} 1151 if (namelist) {PetscValidPointer(namelist,3); *namelist = PETSC_NULL;} 1152 if (innerislist) {PetscValidPointer(innerislist,4); *innerislist = PETSC_NULL;} 1153 if (outerislist) {PetscValidPointer(outerislist,5); *outerislist = PETSC_NULL;} 1154 if (dmlist) {PetscValidPointer(dmlist,6); *dmlist = PETSC_NULL;} 1155 if(dm->ops->createdomaindecomposition) { 1156 ierr = (*dm->ops->createdomaindecomposition)(dm,len,namelist,innerislist,outerislist,dmlist); CHKERRQ(ierr); 1157 } 1158 PetscFunctionReturn(0); 1159 } 1160 1161 #undef __FUNCT__ 1162 #define __FUNCT__ "DMRefine" 1163 /*@ 1164 DMRefine - Refines a DM object 1165 1166 Collective on DM 1167 1168 Input Parameter: 1169 + dm - the DM object 1170 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 1171 1172 Output Parameter: 1173 . dmf - the refined DM, or PETSC_NULL 1174 1175 Note: If no refinement was done, the return value is PETSC_NULL 1176 1177 Level: developer 1178 1179 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1180 @*/ 1181 PetscErrorCode DMRefine(DM dm,MPI_Comm comm,DM *dmf) 1182 { 1183 PetscErrorCode ierr; 1184 DMRefineHookLink link; 1185 1186 PetscFunctionBegin; 1187 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1188 ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr); 1189 if (*dmf) { 1190 (*dmf)->ops->creatematrix = dm->ops->creatematrix; 1191 (*dmf)->ops->initialguess = dm->ops->initialguess; 1192 (*dmf)->ops->function = dm->ops->function; 1193 (*dmf)->ops->functionj = dm->ops->functionj; 1194 if (dm->ops->jacobian != DMComputeJacobianDefault) { 1195 (*dmf)->ops->jacobian = dm->ops->jacobian; 1196 } 1197 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr); 1198 (*dmf)->ctx = dm->ctx; 1199 (*dmf)->leveldown = dm->leveldown; 1200 (*dmf)->levelup = dm->levelup + 1; 1201 ierr = DMSetMatType(*dmf,dm->mattype);CHKERRQ(ierr); 1202 for (link=dm->refinehook; link; link=link->next) { 1203 if (link->refinehook) {ierr = (*link->refinehook)(dm,*dmf,link->ctx);CHKERRQ(ierr);} 1204 } 1205 } 1206 PetscFunctionReturn(0); 1207 } 1208 1209 #undef __FUNCT__ 1210 #define __FUNCT__ "DMRefineHookAdd" 1211 /*@ 1212 DMRefineHookAdd - adds a callback to be run when interpolating a nonlinear problem to a finer grid 1213 1214 Logically Collective 1215 1216 Input Arguments: 1217 + coarse - nonlinear solver context on which to run a hook when restricting to a coarser level 1218 . refinehook - function to run when setting up a coarser level 1219 . interphook - function to run to update data on finer levels (once per SNESSolve()) 1220 - ctx - [optional] user-defined context for provide data for the hooks (may be PETSC_NULL) 1221 1222 Calling sequence of refinehook: 1223 $ refinehook(DM coarse,DM fine,void *ctx); 1224 1225 + coarse - coarse level DM 1226 . fine - fine level DM to interpolate problem to 1227 - ctx - optional user-defined function context 1228 1229 Calling sequence for interphook: 1230 $ interphook(DM coarse,Mat interp,DM fine,void *ctx) 1231 1232 + coarse - coarse level DM 1233 . interp - matrix interpolating a coarse-level solution to the finer grid 1234 . fine - fine level DM to update 1235 - ctx - optional user-defined function context 1236 1237 Level: advanced 1238 1239 Notes: 1240 This function is only needed if auxiliary data needs to be passed to fine grids while grid sequencing 1241 1242 If this function is called multiple times, the hooks will be run in the order they are added. 1243 1244 .seealso: DMCoarsenHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 1245 @*/ 1246 PetscErrorCode DMRefineHookAdd(DM coarse,PetscErrorCode (*refinehook)(DM,DM,void*),PetscErrorCode (*interphook)(DM,Mat,DM,void*),void *ctx) 1247 { 1248 PetscErrorCode ierr; 1249 DMRefineHookLink link,*p; 1250 1251 PetscFunctionBegin; 1252 PetscValidHeaderSpecific(coarse,DM_CLASSID,1); 1253 for (p=&coarse->refinehook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */ 1254 ierr = PetscMalloc(sizeof(struct _DMRefineHookLink),&link);CHKERRQ(ierr); 1255 link->refinehook = refinehook; 1256 link->interphook = interphook; 1257 link->ctx = ctx; 1258 link->next = PETSC_NULL; 1259 *p = link; 1260 PetscFunctionReturn(0); 1261 } 1262 1263 #undef __FUNCT__ 1264 #define __FUNCT__ "DMInterpolate" 1265 /*@ 1266 DMInterpolate - interpolates user-defined problem data to a finer DM by running hooks registered by DMRefineHookAdd() 1267 1268 Collective if any hooks are 1269 1270 Input Arguments: 1271 + coarse - coarser DM to use as a base 1272 . restrct - interpolation matrix, apply using MatInterpolate() 1273 - fine - finer DM to update 1274 1275 Level: developer 1276 1277 .seealso: DMRefineHookAdd(), MatInterpolate() 1278 @*/ 1279 PetscErrorCode DMInterpolate(DM coarse,Mat interp,DM fine) 1280 { 1281 PetscErrorCode ierr; 1282 DMRefineHookLink link; 1283 1284 PetscFunctionBegin; 1285 for (link=fine->refinehook; link; link=link->next) { 1286 if (link->interphook) {ierr = (*link->interphook)(coarse,interp,fine,link->ctx);CHKERRQ(ierr);} 1287 } 1288 PetscFunctionReturn(0); 1289 } 1290 1291 #undef __FUNCT__ 1292 #define __FUNCT__ "DMGetRefineLevel" 1293 /*@ 1294 DMGetRefineLevel - Get's the number of refinements that have generated this DM. 1295 1296 Not Collective 1297 1298 Input Parameter: 1299 . dm - the DM object 1300 1301 Output Parameter: 1302 . level - number of refinements 1303 1304 Level: developer 1305 1306 .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1307 1308 @*/ 1309 PetscErrorCode DMGetRefineLevel(DM dm,PetscInt *level) 1310 { 1311 PetscFunctionBegin; 1312 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1313 *level = dm->levelup; 1314 PetscFunctionReturn(0); 1315 } 1316 1317 #undef __FUNCT__ 1318 #define __FUNCT__ "DMGlobalToLocalBegin" 1319 /*@ 1320 DMGlobalToLocalBegin - Begins updating local vectors from global vector 1321 1322 Neighbor-wise Collective on DM 1323 1324 Input Parameters: 1325 + dm - the DM object 1326 . g - the global vector 1327 . mode - INSERT_VALUES or ADD_VALUES 1328 - l - the local vector 1329 1330 1331 Level: beginner 1332 1333 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 1334 1335 @*/ 1336 PetscErrorCode DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) 1337 { 1338 PetscSF sf; 1339 PetscErrorCode ierr; 1340 1341 PetscFunctionBegin; 1342 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1343 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 1344 if (sf) { 1345 PetscScalar *lArray, *gArray; 1346 1347 if (mode == ADD_VALUES) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 1348 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 1349 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 1350 ierr = PetscSFBcastBegin(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr); 1351 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 1352 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 1353 } else { 1354 ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 1355 } 1356 PetscFunctionReturn(0); 1357 } 1358 1359 #undef __FUNCT__ 1360 #define __FUNCT__ "DMGlobalToLocalEnd" 1361 /*@ 1362 DMGlobalToLocalEnd - Ends updating local vectors from global vector 1363 1364 Neighbor-wise Collective on DM 1365 1366 Input Parameters: 1367 + dm - the DM object 1368 . g - the global vector 1369 . mode - INSERT_VALUES or ADD_VALUES 1370 - l - the local vector 1371 1372 1373 Level: beginner 1374 1375 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 1376 1377 @*/ 1378 PetscErrorCode DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) 1379 { 1380 PetscSF sf; 1381 PetscErrorCode ierr; 1382 PetscScalar *lArray, *gArray; 1383 1384 PetscFunctionBegin; 1385 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1386 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 1387 if (sf) { 1388 if (mode == ADD_VALUES) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 1389 1390 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 1391 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 1392 ierr = PetscSFBcastEnd(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr); 1393 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 1394 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 1395 } else { 1396 ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 1397 } 1398 PetscFunctionReturn(0); 1399 } 1400 1401 #undef __FUNCT__ 1402 #define __FUNCT__ "DMLocalToGlobalBegin" 1403 /*@ 1404 DMLocalToGlobalBegin - updates global vectors from local vectors 1405 1406 Neighbor-wise Collective on DM 1407 1408 Input Parameters: 1409 + dm - the DM object 1410 . l - the local vector 1411 . 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 1412 base point. 1413 - - the global vector 1414 1415 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 1416 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 1417 global array to the final global array with VecAXPY(). 1418 1419 Level: beginner 1420 1421 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin() 1422 1423 @*/ 1424 PetscErrorCode DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g) 1425 { 1426 PetscSF sf; 1427 PetscErrorCode ierr; 1428 1429 PetscFunctionBegin; 1430 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1431 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 1432 if (sf) { 1433 MPI_Op op; 1434 PetscScalar *lArray, *gArray; 1435 1436 switch(mode) { 1437 case INSERT_VALUES: 1438 case INSERT_ALL_VALUES: 1439 #if defined(PETSC_HAVE_MPI_REPLACE) 1440 op = MPI_REPLACE; break; 1441 #else 1442 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No support for INSERT_VALUES without an MPI-2 implementation"); 1443 #endif 1444 case ADD_VALUES: 1445 case ADD_ALL_VALUES: 1446 op = MPI_SUM; break; 1447 default: 1448 SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 1449 } 1450 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 1451 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 1452 ierr = PetscSFReduceBegin(sf, MPIU_SCALAR, lArray, gArray, op);CHKERRQ(ierr); 1453 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 1454 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 1455 } else { 1456 ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 1457 } 1458 PetscFunctionReturn(0); 1459 } 1460 1461 #undef __FUNCT__ 1462 #define __FUNCT__ "DMLocalToGlobalEnd" 1463 /*@ 1464 DMLocalToGlobalEnd - updates global vectors from local vectors 1465 1466 Neighbor-wise Collective on DM 1467 1468 Input Parameters: 1469 + dm - the DM object 1470 . l - the local vector 1471 . mode - INSERT_VALUES or ADD_VALUES 1472 - g - the global vector 1473 1474 1475 Level: beginner 1476 1477 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd() 1478 1479 @*/ 1480 PetscErrorCode DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g) 1481 { 1482 PetscSF sf; 1483 PetscErrorCode ierr; 1484 1485 PetscFunctionBegin; 1486 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1487 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 1488 if (sf) { 1489 MPI_Op op; 1490 PetscScalar *lArray, *gArray; 1491 1492 switch(mode) { 1493 case INSERT_VALUES: 1494 case INSERT_ALL_VALUES: 1495 #if defined(PETSC_HAVE_MPI_REPLACE) 1496 op = MPI_REPLACE; break; 1497 #else 1498 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No support for INSERT_VALUES without an MPI-2 implementation"); 1499 #endif 1500 case ADD_VALUES: 1501 case ADD_ALL_VALUES: 1502 op = MPI_SUM; break; 1503 default: 1504 SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 1505 } 1506 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 1507 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 1508 ierr = PetscSFReduceEnd(sf, MPIU_SCALAR, lArray, gArray, op);CHKERRQ(ierr); 1509 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 1510 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 1511 } else { 1512 ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 1513 } 1514 PetscFunctionReturn(0); 1515 } 1516 1517 #undef __FUNCT__ 1518 #define __FUNCT__ "DMComputeJacobianDefault" 1519 /*@ 1520 DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided 1521 1522 Collective on DM 1523 1524 Input Parameter: 1525 + dm - the DM object 1526 . x - location to compute Jacobian at; may be ignored for linear problems 1527 . A - matrix that defines the operator for the linear solve 1528 - B - the matrix used to construct the preconditioner 1529 1530 Level: developer 1531 1532 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1533 DMSetFunction() 1534 1535 @*/ 1536 PetscErrorCode DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 1537 { 1538 PetscErrorCode ierr; 1539 1540 PetscFunctionBegin; 1541 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1542 *stflag = SAME_NONZERO_PATTERN; 1543 ierr = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr); 1544 if (A != B) { 1545 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1546 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1547 } 1548 PetscFunctionReturn(0); 1549 } 1550 1551 #undef __FUNCT__ 1552 #define __FUNCT__ "DMCoarsen" 1553 /*@ 1554 DMCoarsen - Coarsens a DM object 1555 1556 Collective on DM 1557 1558 Input Parameter: 1559 + dm - the DM object 1560 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 1561 1562 Output Parameter: 1563 . dmc - the coarsened DM 1564 1565 Level: developer 1566 1567 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1568 1569 @*/ 1570 PetscErrorCode DMCoarsen(DM dm, MPI_Comm comm, DM *dmc) 1571 { 1572 PetscErrorCode ierr; 1573 DMCoarsenHookLink link; 1574 1575 PetscFunctionBegin; 1576 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1577 ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr); 1578 (*dmc)->ops->creatematrix = dm->ops->creatematrix; 1579 (*dmc)->ops->initialguess = dm->ops->initialguess; 1580 (*dmc)->ops->function = dm->ops->function; 1581 (*dmc)->ops->functionj = dm->ops->functionj; 1582 if (dm->ops->jacobian != DMComputeJacobianDefault) { 1583 (*dmc)->ops->jacobian = dm->ops->jacobian; 1584 } 1585 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr); 1586 (*dmc)->ctx = dm->ctx; 1587 (*dmc)->levelup = dm->levelup; 1588 (*dmc)->leveldown = dm->leveldown + 1; 1589 ierr = DMSetMatType(*dmc,dm->mattype);CHKERRQ(ierr); 1590 for (link=dm->coarsenhook; link; link=link->next) { 1591 if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);} 1592 } 1593 PetscFunctionReturn(0); 1594 } 1595 1596 #undef __FUNCT__ 1597 #define __FUNCT__ "DMCoarsenHookAdd" 1598 /*@ 1599 DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid 1600 1601 Logically Collective 1602 1603 Input Arguments: 1604 + fine - nonlinear solver context on which to run a hook when restricting to a coarser level 1605 . coarsenhook - function to run when setting up a coarser level 1606 . restricthook - function to run to update data on coarser levels (once per SNESSolve()) 1607 - ctx - [optional] user-defined context for provide data for the hooks (may be PETSC_NULL) 1608 1609 Calling sequence of coarsenhook: 1610 $ coarsenhook(DM fine,DM coarse,void *ctx); 1611 1612 + fine - fine level DM 1613 . coarse - coarse level DM to restrict problem to 1614 - ctx - optional user-defined function context 1615 1616 Calling sequence for restricthook: 1617 $ restricthook(DM fine,Mat mrestrict,Vec rscale,Mat inject,DM coarse,void *ctx) 1618 1619 + fine - fine level DM 1620 . mrestrict - matrix restricting a fine-level solution to the coarse grid 1621 . rscale - scaling vector for restriction 1622 . inject - matrix restricting by injection 1623 . coarse - coarse level DM to update 1624 - ctx - optional user-defined function context 1625 1626 Level: advanced 1627 1628 Notes: 1629 This function is only needed if auxiliary data needs to be set up on coarse grids. 1630 1631 If this function is called multiple times, the hooks will be run in the order they are added. 1632 1633 In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to 1634 extract the finest level information from its context (instead of from the SNES). 1635 1636 .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 1637 @*/ 1638 PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx) 1639 { 1640 PetscErrorCode ierr; 1641 DMCoarsenHookLink link,*p; 1642 1643 PetscFunctionBegin; 1644 PetscValidHeaderSpecific(fine,DM_CLASSID,1); 1645 for (p=&fine->coarsenhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */ 1646 ierr = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr); 1647 link->coarsenhook = coarsenhook; 1648 link->restricthook = restricthook; 1649 link->ctx = ctx; 1650 link->next = PETSC_NULL; 1651 *p = link; 1652 PetscFunctionReturn(0); 1653 } 1654 1655 #undef __FUNCT__ 1656 #define __FUNCT__ "DMRestrict" 1657 /*@ 1658 DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd() 1659 1660 Collective if any hooks are 1661 1662 Input Arguments: 1663 + fine - finer DM to use as a base 1664 . restrct - restriction matrix, apply using MatRestrict() 1665 . inject - injection matrix, also use MatRestrict() 1666 - coarse - coarer DM to update 1667 1668 Level: developer 1669 1670 .seealso: DMCoarsenHookAdd(), MatRestrict() 1671 @*/ 1672 PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse) 1673 { 1674 PetscErrorCode ierr; 1675 DMCoarsenHookLink link; 1676 1677 PetscFunctionBegin; 1678 for (link=fine->coarsenhook; link; link=link->next) { 1679 if (link->restricthook) {ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);} 1680 } 1681 PetscFunctionReturn(0); 1682 } 1683 1684 #undef __FUNCT__ 1685 #define __FUNCT__ "DMGetCoarsenLevel" 1686 /*@ 1687 DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM. 1688 1689 Not Collective 1690 1691 Input Parameter: 1692 . dm - the DM object 1693 1694 Output Parameter: 1695 . level - number of coarsenings 1696 1697 Level: developer 1698 1699 .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1700 1701 @*/ 1702 PetscErrorCode DMGetCoarsenLevel(DM dm,PetscInt *level) 1703 { 1704 PetscFunctionBegin; 1705 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1706 *level = dm->leveldown; 1707 PetscFunctionReturn(0); 1708 } 1709 1710 1711 1712 #undef __FUNCT__ 1713 #define __FUNCT__ "DMRefineHierarchy" 1714 /*@C 1715 DMRefineHierarchy - Refines a DM object, all levels at once 1716 1717 Collective on DM 1718 1719 Input Parameter: 1720 + dm - the DM object 1721 - nlevels - the number of levels of refinement 1722 1723 Output Parameter: 1724 . dmf - the refined DM hierarchy 1725 1726 Level: developer 1727 1728 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1729 1730 @*/ 1731 PetscErrorCode DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[]) 1732 { 1733 PetscErrorCode ierr; 1734 1735 PetscFunctionBegin; 1736 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1737 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 1738 if (nlevels == 0) PetscFunctionReturn(0); 1739 if (dm->ops->refinehierarchy) { 1740 ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr); 1741 } else if (dm->ops->refine) { 1742 PetscInt i; 1743 1744 ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr); 1745 for (i=1; i<nlevels; i++) { 1746 ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr); 1747 } 1748 } else { 1749 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet"); 1750 } 1751 PetscFunctionReturn(0); 1752 } 1753 1754 #undef __FUNCT__ 1755 #define __FUNCT__ "DMCoarsenHierarchy" 1756 /*@C 1757 DMCoarsenHierarchy - Coarsens a DM object, all levels at once 1758 1759 Collective on DM 1760 1761 Input Parameter: 1762 + dm - the DM object 1763 - nlevels - the number of levels of coarsening 1764 1765 Output Parameter: 1766 . dmc - the coarsened DM hierarchy 1767 1768 Level: developer 1769 1770 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1771 1772 @*/ 1773 PetscErrorCode DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[]) 1774 { 1775 PetscErrorCode ierr; 1776 1777 PetscFunctionBegin; 1778 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1779 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 1780 if (nlevels == 0) PetscFunctionReturn(0); 1781 PetscValidPointer(dmc,3); 1782 if (dm->ops->coarsenhierarchy) { 1783 ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr); 1784 } else if (dm->ops->coarsen) { 1785 PetscInt i; 1786 1787 ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr); 1788 for (i=1; i<nlevels; i++) { 1789 ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr); 1790 } 1791 } else { 1792 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet"); 1793 } 1794 PetscFunctionReturn(0); 1795 } 1796 1797 #undef __FUNCT__ 1798 #define __FUNCT__ "DMCreateAggregates" 1799 /*@ 1800 DMCreateAggregates - Gets the aggregates that map between 1801 grids associated with two DMs. 1802 1803 Collective on DM 1804 1805 Input Parameters: 1806 + dmc - the coarse grid DM 1807 - dmf - the fine grid DM 1808 1809 Output Parameters: 1810 . rest - the restriction matrix (transpose of the projection matrix) 1811 1812 Level: intermediate 1813 1814 .keywords: interpolation, restriction, multigrid 1815 1816 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation() 1817 @*/ 1818 PetscErrorCode DMCreateAggregates(DM dmc, DM dmf, Mat *rest) 1819 { 1820 PetscErrorCode ierr; 1821 1822 PetscFunctionBegin; 1823 PetscValidHeaderSpecific(dmc,DM_CLASSID,1); 1824 PetscValidHeaderSpecific(dmf,DM_CLASSID,2); 1825 ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr); 1826 PetscFunctionReturn(0); 1827 } 1828 1829 #undef __FUNCT__ 1830 #define __FUNCT__ "DMSetApplicationContextDestroy" 1831 /*@C 1832 DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed 1833 1834 Not Collective 1835 1836 Input Parameters: 1837 + dm - the DM object 1838 - destroy - the destroy function 1839 1840 Level: intermediate 1841 1842 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1843 1844 @*/ 1845 PetscErrorCode DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**)) 1846 { 1847 PetscFunctionBegin; 1848 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1849 dm->ctxdestroy = destroy; 1850 PetscFunctionReturn(0); 1851 } 1852 1853 #undef __FUNCT__ 1854 #define __FUNCT__ "DMSetApplicationContext" 1855 /*@ 1856 DMSetApplicationContext - Set a user context into a DM object 1857 1858 Not Collective 1859 1860 Input Parameters: 1861 + dm - the DM object 1862 - ctx - the user context 1863 1864 Level: intermediate 1865 1866 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1867 1868 @*/ 1869 PetscErrorCode DMSetApplicationContext(DM dm,void *ctx) 1870 { 1871 PetscFunctionBegin; 1872 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1873 dm->ctx = ctx; 1874 PetscFunctionReturn(0); 1875 } 1876 1877 #undef __FUNCT__ 1878 #define __FUNCT__ "DMGetApplicationContext" 1879 /*@ 1880 DMGetApplicationContext - Gets a user context from a DM object 1881 1882 Not Collective 1883 1884 Input Parameter: 1885 . dm - the DM object 1886 1887 Output Parameter: 1888 . ctx - the user context 1889 1890 Level: intermediate 1891 1892 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1893 1894 @*/ 1895 PetscErrorCode DMGetApplicationContext(DM dm,void *ctx) 1896 { 1897 PetscFunctionBegin; 1898 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1899 *(void**)ctx = dm->ctx; 1900 PetscFunctionReturn(0); 1901 } 1902 1903 #undef __FUNCT__ 1904 #define __FUNCT__ "DMSetInitialGuess" 1905 /*@C 1906 DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers 1907 1908 Logically Collective on DM 1909 1910 Input Parameter: 1911 + dm - the DM object to destroy 1912 - f - the function to compute the initial guess 1913 1914 Level: intermediate 1915 1916 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1917 1918 @*/ 1919 PetscErrorCode DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec)) 1920 { 1921 PetscFunctionBegin; 1922 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1923 dm->ops->initialguess = f; 1924 PetscFunctionReturn(0); 1925 } 1926 1927 #undef __FUNCT__ 1928 #define __FUNCT__ "DMSetFunction" 1929 /*@C 1930 DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES 1931 1932 Logically Collective on DM 1933 1934 Input Parameter: 1935 + dm - the DM object 1936 - f - the function to compute (use PETSC_NULL to cancel a previous function that was set) 1937 1938 Level: intermediate 1939 1940 Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian 1941 computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian. 1942 1943 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1944 DMSetJacobian() 1945 1946 @*/ 1947 PetscErrorCode DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 1948 { 1949 PetscFunctionBegin; 1950 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1951 dm->ops->function = f; 1952 if (f) { 1953 dm->ops->functionj = f; 1954 } 1955 PetscFunctionReturn(0); 1956 } 1957 1958 #undef __FUNCT__ 1959 #define __FUNCT__ "DMSetJacobian" 1960 /*@C 1961 DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES 1962 1963 Logically Collective on DM 1964 1965 Input Parameter: 1966 + dm - the DM object to destroy 1967 - f - the function to compute the matrix entries 1968 1969 Level: intermediate 1970 1971 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1972 DMSetFunction() 1973 1974 @*/ 1975 PetscErrorCode DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*)) 1976 { 1977 PetscFunctionBegin; 1978 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1979 dm->ops->jacobian = f; 1980 PetscFunctionReturn(0); 1981 } 1982 1983 #undef __FUNCT__ 1984 #define __FUNCT__ "DMSetVariableBounds" 1985 /*@C 1986 DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI. 1987 1988 Logically Collective on DM 1989 1990 Input Parameter: 1991 + dm - the DM object 1992 - f - the function that computes variable bounds used by SNESVI (use PETSC_NULL to cancel a previous function that was set) 1993 1994 Level: intermediate 1995 1996 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1997 DMSetJacobian() 1998 1999 @*/ 2000 PetscErrorCode DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 2001 { 2002 PetscFunctionBegin; 2003 dm->ops->computevariablebounds = f; 2004 PetscFunctionReturn(0); 2005 } 2006 2007 #undef __FUNCT__ 2008 #define __FUNCT__ "DMHasVariableBounds" 2009 /*@ 2010 DMHasVariableBounds - does the DM object have a variable bounds function? 2011 2012 Not Collective 2013 2014 Input Parameter: 2015 . dm - the DM object to destroy 2016 2017 Output Parameter: 2018 . flg - PETSC_TRUE if the variable bounds function exists 2019 2020 Level: developer 2021 2022 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 2023 2024 @*/ 2025 PetscErrorCode DMHasVariableBounds(DM dm,PetscBool *flg) 2026 { 2027 PetscFunctionBegin; 2028 *flg = (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE; 2029 PetscFunctionReturn(0); 2030 } 2031 2032 #undef __FUNCT__ 2033 #define __FUNCT__ "DMComputeVariableBounds" 2034 /*@C 2035 DMComputeVariableBounds - compute variable bounds used by SNESVI. 2036 2037 Logically Collective on DM 2038 2039 Input Parameters: 2040 + dm - the DM object to destroy 2041 - x - current solution at which the bounds are computed 2042 2043 Output parameters: 2044 + xl - lower bound 2045 - xu - upper bound 2046 2047 Level: intermediate 2048 2049 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 2050 DMSetFunction(), DMSetVariableBounds() 2051 2052 @*/ 2053 PetscErrorCode DMComputeVariableBounds(DM dm, Vec xl, Vec xu) 2054 { 2055 PetscErrorCode ierr; 2056 PetscFunctionBegin; 2057 PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 2058 PetscValidHeaderSpecific(xu,VEC_CLASSID,2); 2059 if(dm->ops->computevariablebounds) { 2060 ierr = (*dm->ops->computevariablebounds)(dm, xl,xu); CHKERRQ(ierr); 2061 } 2062 else { 2063 ierr = VecSet(xl,SNES_VI_NINF); CHKERRQ(ierr); 2064 ierr = VecSet(xu,SNES_VI_INF); CHKERRQ(ierr); 2065 } 2066 PetscFunctionReturn(0); 2067 } 2068 2069 #undef __FUNCT__ 2070 #define __FUNCT__ "DMComputeInitialGuess" 2071 /*@ 2072 DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers 2073 2074 Collective on DM 2075 2076 Input Parameter: 2077 + dm - the DM object to destroy 2078 - x - the vector to hold the initial guess values 2079 2080 Level: developer 2081 2082 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat() 2083 2084 @*/ 2085 PetscErrorCode DMComputeInitialGuess(DM dm,Vec x) 2086 { 2087 PetscErrorCode ierr; 2088 PetscFunctionBegin; 2089 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2090 if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()"); 2091 ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr); 2092 PetscFunctionReturn(0); 2093 } 2094 2095 #undef __FUNCT__ 2096 #define __FUNCT__ "DMHasInitialGuess" 2097 /*@ 2098 DMHasInitialGuess - does the DM object have an initial guess function 2099 2100 Not Collective 2101 2102 Input Parameter: 2103 . dm - the DM object to destroy 2104 2105 Output Parameter: 2106 . flg - PETSC_TRUE if function exists 2107 2108 Level: developer 2109 2110 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 2111 2112 @*/ 2113 PetscErrorCode DMHasInitialGuess(DM dm,PetscBool *flg) 2114 { 2115 PetscFunctionBegin; 2116 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2117 *flg = (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE; 2118 PetscFunctionReturn(0); 2119 } 2120 2121 #undef __FUNCT__ 2122 #define __FUNCT__ "DMHasFunction" 2123 /*@ 2124 DMHasFunction - does the DM object have a function 2125 2126 Not Collective 2127 2128 Input Parameter: 2129 . dm - the DM object to destroy 2130 2131 Output Parameter: 2132 . flg - PETSC_TRUE if function exists 2133 2134 Level: developer 2135 2136 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 2137 2138 @*/ 2139 PetscErrorCode DMHasFunction(DM dm,PetscBool *flg) 2140 { 2141 PetscFunctionBegin; 2142 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2143 *flg = (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE; 2144 PetscFunctionReturn(0); 2145 } 2146 2147 #undef __FUNCT__ 2148 #define __FUNCT__ "DMHasJacobian" 2149 /*@ 2150 DMHasJacobian - does the DM object have a matrix function 2151 2152 Not Collective 2153 2154 Input Parameter: 2155 . dm - the DM object to destroy 2156 2157 Output Parameter: 2158 . flg - PETSC_TRUE if function exists 2159 2160 Level: developer 2161 2162 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 2163 2164 @*/ 2165 PetscErrorCode DMHasJacobian(DM dm,PetscBool *flg) 2166 { 2167 PetscFunctionBegin; 2168 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2169 *flg = (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE; 2170 PetscFunctionReturn(0); 2171 } 2172 2173 #undef __FUNCT__ 2174 #define __FUNCT__ "DMSetVec" 2175 /*@C 2176 DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear. 2177 2178 Collective on DM 2179 2180 Input Parameter: 2181 + dm - the DM object 2182 - x - location to compute residual and Jacobian, if PETSC_NULL is passed to those routines; will be PETSC_NULL for linear problems. 2183 2184 Level: developer 2185 2186 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 2187 DMSetFunction(), DMSetJacobian(), DMSetVariableBounds() 2188 2189 @*/ 2190 PetscErrorCode DMSetVec(DM dm,Vec x) 2191 { 2192 PetscErrorCode ierr; 2193 PetscFunctionBegin; 2194 if (x) { 2195 if (!dm->x) { 2196 ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 2197 } 2198 ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 2199 } 2200 else if(dm->x) { 2201 ierr = VecDestroy(&dm->x); CHKERRQ(ierr); 2202 } 2203 PetscFunctionReturn(0); 2204 } 2205 2206 2207 #undef __FUNCT__ 2208 #define __FUNCT__ "DMComputeFunction" 2209 /*@ 2210 DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES 2211 2212 Collective on DM 2213 2214 Input Parameter: 2215 + dm - the DM object to destroy 2216 . x - the location where the function is evaluationed, may be ignored for linear problems 2217 - b - the vector to hold the right hand side entries 2218 2219 Level: developer 2220 2221 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 2222 DMSetJacobian() 2223 2224 @*/ 2225 PetscErrorCode DMComputeFunction(DM dm,Vec x,Vec b) 2226 { 2227 PetscErrorCode ierr; 2228 PetscFunctionBegin; 2229 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2230 if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()"); 2231 PetscStackPush("DM user function"); 2232 ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr); 2233 PetscStackPop; 2234 PetscFunctionReturn(0); 2235 } 2236 2237 2238 2239 #undef __FUNCT__ 2240 #define __FUNCT__ "DMComputeJacobian" 2241 /*@ 2242 DMComputeJacobian - compute the matrix entries for the solver 2243 2244 Collective on DM 2245 2246 Input Parameter: 2247 + dm - the DM object 2248 . x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM 2249 . A - matrix that defines the operator for the linear solve 2250 - B - the matrix used to construct the preconditioner 2251 2252 Level: developer 2253 2254 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 2255 DMSetFunction() 2256 2257 @*/ 2258 PetscErrorCode DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 2259 { 2260 PetscErrorCode ierr; 2261 2262 PetscFunctionBegin; 2263 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2264 if (!dm->ops->jacobian) { 2265 ISColoring coloring; 2266 MatFDColoring fd; 2267 const MatType mtype; 2268 2269 ierr = PetscObjectGetType((PetscObject)B,&mtype);CHKERRQ(ierr); 2270 ierr = DMCreateColoring(dm,dm->coloringtype,mtype,&coloring);CHKERRQ(ierr); 2271 ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr); 2272 ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 2273 ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr); 2274 ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr); 2275 ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr); 2276 2277 dm->fd = fd; 2278 dm->ops->jacobian = DMComputeJacobianDefault; 2279 2280 /* don't know why this is needed */ 2281 ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 2282 } 2283 if (!x) x = dm->x; 2284 ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr); 2285 2286 /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */ 2287 if (x) { 2288 if (!dm->x) { 2289 ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 2290 } 2291 ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 2292 } 2293 if (A != B) { 2294 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2295 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2296 } 2297 PetscFunctionReturn(0); 2298 } 2299 2300 2301 PetscFList DMList = PETSC_NULL; 2302 PetscBool DMRegisterAllCalled = PETSC_FALSE; 2303 2304 #undef __FUNCT__ 2305 #define __FUNCT__ "DMSetType" 2306 /*@C 2307 DMSetType - Builds a DM, for a particular DM implementation. 2308 2309 Collective on DM 2310 2311 Input Parameters: 2312 + dm - The DM object 2313 - method - The name of the DM type 2314 2315 Options Database Key: 2316 . -dm_type <type> - Sets the DM type; use -help for a list of available types 2317 2318 Notes: 2319 See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D). 2320 2321 Level: intermediate 2322 2323 .keywords: DM, set, type 2324 .seealso: DMGetType(), DMCreate() 2325 @*/ 2326 PetscErrorCode DMSetType(DM dm, const DMType method) 2327 { 2328 PetscErrorCode (*r)(DM); 2329 PetscBool match; 2330 PetscErrorCode ierr; 2331 2332 PetscFunctionBegin; 2333 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 2334 ierr = PetscObjectTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr); 2335 if (match) PetscFunctionReturn(0); 2336 2337 if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 2338 ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 2339 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method); 2340 2341 if (dm->ops->destroy) { 2342 ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr); 2343 dm->ops->destroy = PETSC_NULL; 2344 } 2345 ierr = (*r)(dm);CHKERRQ(ierr); 2346 ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr); 2347 PetscFunctionReturn(0); 2348 } 2349 2350 #undef __FUNCT__ 2351 #define __FUNCT__ "DMGetType" 2352 /*@C 2353 DMGetType - Gets the DM type name (as a string) from the DM. 2354 2355 Not Collective 2356 2357 Input Parameter: 2358 . dm - The DM 2359 2360 Output Parameter: 2361 . type - The DM type name 2362 2363 Level: intermediate 2364 2365 .keywords: DM, get, type, name 2366 .seealso: DMSetType(), DMCreate() 2367 @*/ 2368 PetscErrorCode DMGetType(DM dm, const DMType *type) 2369 { 2370 PetscErrorCode ierr; 2371 2372 PetscFunctionBegin; 2373 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 2374 PetscValidCharPointer(type,2); 2375 if (!DMRegisterAllCalled) { 2376 ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr); 2377 } 2378 *type = ((PetscObject)dm)->type_name; 2379 PetscFunctionReturn(0); 2380 } 2381 2382 #undef __FUNCT__ 2383 #define __FUNCT__ "DMConvert" 2384 /*@C 2385 DMConvert - Converts a DM to another DM, either of the same or different type. 2386 2387 Collective on DM 2388 2389 Input Parameters: 2390 + dm - the DM 2391 - newtype - new DM type (use "same" for the same type) 2392 2393 Output Parameter: 2394 . M - pointer to new DM 2395 2396 Notes: 2397 Cannot be used to convert a sequential DM to parallel or parallel to sequential, 2398 the MPI communicator of the generated DM is always the same as the communicator 2399 of the input DM. 2400 2401 Level: intermediate 2402 2403 .seealso: DMCreate() 2404 @*/ 2405 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M) 2406 { 2407 DM B; 2408 char convname[256]; 2409 PetscBool sametype, issame; 2410 PetscErrorCode ierr; 2411 2412 PetscFunctionBegin; 2413 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2414 PetscValidType(dm,1); 2415 PetscValidPointer(M,3); 2416 ierr = PetscObjectTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr); 2417 ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr); 2418 { 2419 PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL; 2420 2421 /* 2422 Order of precedence: 2423 1) See if a specialized converter is known to the current DM. 2424 2) See if a specialized converter is known to the desired DM class. 2425 3) See if a good general converter is registered for the desired class 2426 4) See if a good general converter is known for the current matrix. 2427 5) Use a really basic converter. 2428 */ 2429 2430 /* 1) See if a specialized converter is known to the current DM and the desired class */ 2431 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 2432 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 2433 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 2434 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 2435 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 2436 ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr); 2437 if (conv) goto foundconv; 2438 2439 /* 2) See if a specialized converter is known to the desired DM class. */ 2440 ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr); 2441 ierr = DMSetType(B, newtype);CHKERRQ(ierr); 2442 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 2443 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 2444 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 2445 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 2446 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 2447 ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr); 2448 if (conv) { 2449 ierr = DMDestroy(&B);CHKERRQ(ierr); 2450 goto foundconv; 2451 } 2452 2453 #if 0 2454 /* 3) See if a good general converter is registered for the desired class */ 2455 conv = B->ops->convertfrom; 2456 ierr = DMDestroy(&B);CHKERRQ(ierr); 2457 if (conv) goto foundconv; 2458 2459 /* 4) See if a good general converter is known for the current matrix */ 2460 if (dm->ops->convert) { 2461 conv = dm->ops->convert; 2462 } 2463 if (conv) goto foundconv; 2464 #endif 2465 2466 /* 5) Use a really basic converter. */ 2467 SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype); 2468 2469 foundconv: 2470 ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 2471 ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr); 2472 ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 2473 } 2474 ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr); 2475 PetscFunctionReturn(0); 2476 } 2477 2478 /*--------------------------------------------------------------------------------------------------------------------*/ 2479 2480 #undef __FUNCT__ 2481 #define __FUNCT__ "DMRegister" 2482 /*@C 2483 DMRegister - See DMRegisterDynamic() 2484 2485 Level: advanced 2486 @*/ 2487 PetscErrorCode DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM)) 2488 { 2489 char fullname[PETSC_MAX_PATH_LEN]; 2490 PetscErrorCode ierr; 2491 2492 PetscFunctionBegin; 2493 ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr); 2494 ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr); 2495 ierr = PetscStrcat(fullname, name);CHKERRQ(ierr); 2496 ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr); 2497 PetscFunctionReturn(0); 2498 } 2499 2500 2501 /*--------------------------------------------------------------------------------------------------------------------*/ 2502 #undef __FUNCT__ 2503 #define __FUNCT__ "DMRegisterDestroy" 2504 /*@C 2505 DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic(). 2506 2507 Not Collective 2508 2509 Level: advanced 2510 2511 .keywords: DM, register, destroy 2512 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic() 2513 @*/ 2514 PetscErrorCode DMRegisterDestroy(void) 2515 { 2516 PetscErrorCode ierr; 2517 2518 PetscFunctionBegin; 2519 ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr); 2520 DMRegisterAllCalled = PETSC_FALSE; 2521 PetscFunctionReturn(0); 2522 } 2523 2524 #if defined(PETSC_HAVE_MATLAB_ENGINE) 2525 #include <mex.h> 2526 2527 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext; 2528 2529 #undef __FUNCT__ 2530 #define __FUNCT__ "DMComputeFunction_Matlab" 2531 /* 2532 DMComputeFunction_Matlab - Calls the function that has been set with 2533 DMSetFunctionMatlab(). 2534 2535 For linear problems x is null 2536 2537 .seealso: DMSetFunction(), DMGetFunction() 2538 */ 2539 PetscErrorCode DMComputeFunction_Matlab(DM dm,Vec x,Vec y) 2540 { 2541 PetscErrorCode ierr; 2542 DMMatlabContext *sctx; 2543 int nlhs = 1,nrhs = 4; 2544 mxArray *plhs[1],*prhs[4]; 2545 long long int lx = 0,ly = 0,ls = 0; 2546 2547 PetscFunctionBegin; 2548 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2549 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 2550 PetscCheckSameComm(dm,1,y,3); 2551 2552 /* call Matlab function in ctx with arguments x and y */ 2553 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2554 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 2555 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2556 ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr); 2557 prhs[0] = mxCreateDoubleScalar((double)ls); 2558 prhs[1] = mxCreateDoubleScalar((double)lx); 2559 prhs[2] = mxCreateDoubleScalar((double)ly); 2560 prhs[3] = mxCreateString(sctx->funcname); 2561 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr); 2562 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2563 mxDestroyArray(prhs[0]); 2564 mxDestroyArray(prhs[1]); 2565 mxDestroyArray(prhs[2]); 2566 mxDestroyArray(prhs[3]); 2567 mxDestroyArray(plhs[0]); 2568 PetscFunctionReturn(0); 2569 } 2570 2571 2572 #undef __FUNCT__ 2573 #define __FUNCT__ "DMSetFunctionMatlab" 2574 /* 2575 DMSetFunctionMatlab - Sets the function evaluation routine 2576 2577 */ 2578 PetscErrorCode DMSetFunctionMatlab(DM dm,const char *func) 2579 { 2580 PetscErrorCode ierr; 2581 DMMatlabContext *sctx; 2582 2583 PetscFunctionBegin; 2584 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2585 /* currently sctx is memory bleed */ 2586 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2587 if (!sctx) { 2588 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 2589 } 2590 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2591 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 2592 ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr); 2593 PetscFunctionReturn(0); 2594 } 2595 2596 #undef __FUNCT__ 2597 #define __FUNCT__ "DMComputeJacobian_Matlab" 2598 /* 2599 DMComputeJacobian_Matlab - Calls the function that has been set with 2600 DMSetJacobianMatlab(). 2601 2602 For linear problems x is null 2603 2604 .seealso: DMSetFunction(), DMGetFunction() 2605 */ 2606 PetscErrorCode DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str) 2607 { 2608 PetscErrorCode ierr; 2609 DMMatlabContext *sctx; 2610 int nlhs = 2,nrhs = 5; 2611 mxArray *plhs[2],*prhs[5]; 2612 long long int lx = 0,lA = 0,lB = 0,ls = 0; 2613 2614 PetscFunctionBegin; 2615 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2616 PetscValidHeaderSpecific(A,MAT_CLASSID,3); 2617 2618 /* call MATLAB function in ctx with arguments x, A, and B */ 2619 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2620 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 2621 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2622 ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr); 2623 ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr); 2624 prhs[0] = mxCreateDoubleScalar((double)ls); 2625 prhs[1] = mxCreateDoubleScalar((double)lx); 2626 prhs[2] = mxCreateDoubleScalar((double)lA); 2627 prhs[3] = mxCreateDoubleScalar((double)lB); 2628 prhs[4] = mxCreateString(sctx->jacname); 2629 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr); 2630 *str = (MatStructure) mxGetScalar(plhs[0]); 2631 ierr = (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr); 2632 mxDestroyArray(prhs[0]); 2633 mxDestroyArray(prhs[1]); 2634 mxDestroyArray(prhs[2]); 2635 mxDestroyArray(prhs[3]); 2636 mxDestroyArray(prhs[4]); 2637 mxDestroyArray(plhs[0]); 2638 mxDestroyArray(plhs[1]); 2639 PetscFunctionReturn(0); 2640 } 2641 2642 2643 #undef __FUNCT__ 2644 #define __FUNCT__ "DMSetJacobianMatlab" 2645 /* 2646 DMSetJacobianMatlab - Sets the Jacobian function evaluation routine 2647 2648 */ 2649 PetscErrorCode DMSetJacobianMatlab(DM dm,const char *func) 2650 { 2651 PetscErrorCode ierr; 2652 DMMatlabContext *sctx; 2653 2654 PetscFunctionBegin; 2655 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2656 /* currently sctx is memory bleed */ 2657 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2658 if (!sctx) { 2659 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 2660 } 2661 ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr); 2662 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 2663 ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr); 2664 PetscFunctionReturn(0); 2665 } 2666 #endif 2667 2668 #undef __FUNCT__ 2669 #define __FUNCT__ "DMLoad" 2670 /*@C 2671 DMLoad - Loads a DM that has been stored in binary or HDF5 format 2672 with DMView(). 2673 2674 Collective on PetscViewer 2675 2676 Input Parameters: 2677 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or 2678 some related function before a call to DMLoad(). 2679 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 2680 HDF5 file viewer, obtained from PetscViewerHDF5Open() 2681 2682 Level: intermediate 2683 2684 Notes: 2685 Defaults to the DM DA. 2686 2687 Notes for advanced users: 2688 Most users should not need to know the details of the binary storage 2689 format, since DMLoad() and DMView() completely hide these details. 2690 But for anyone who's interested, the standard binary matrix storage 2691 format is 2692 .vb 2693 has not yet been determined 2694 .ve 2695 2696 In addition, PETSc automatically does the byte swapping for 2697 machines that store the bytes reversed, e.g. DEC alpha, freebsd, 2698 linux, Windows and the paragon; thus if you write your own binary 2699 read/write routines you have to swap the bytes; see PetscBinaryRead() 2700 and PetscBinaryWrite() to see how this may be done. 2701 2702 Concepts: vector^loading from file 2703 2704 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad() 2705 @*/ 2706 PetscErrorCode DMLoad(DM newdm, PetscViewer viewer) 2707 { 2708 PetscErrorCode ierr; 2709 2710 PetscFunctionBegin; 2711 PetscValidHeaderSpecific(newdm,DM_CLASSID,1); 2712 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2713 2714 if (!((PetscObject)newdm)->type_name) { 2715 ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr); 2716 } 2717 ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr); 2718 PetscFunctionReturn(0); 2719 } 2720 2721 /******************************** FEM Support **********************************/ 2722 2723 #undef __FUNCT__ 2724 #define __FUNCT__ "DMPrintCellVector" 2725 PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[]) { 2726 PetscInt f; 2727 PetscErrorCode ierr; 2728 2729 PetscFunctionBegin; 2730 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 2731 for(f = 0; f < len; ++f) { 2732 ierr = PetscPrintf(PETSC_COMM_SELF, " | %G |\n", PetscRealPart(x[f]));CHKERRQ(ierr); 2733 } 2734 PetscFunctionReturn(0); 2735 } 2736 2737 #undef __FUNCT__ 2738 #define __FUNCT__ "DMPrintCellMatrix" 2739 PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[]) { 2740 PetscInt f, g; 2741 PetscErrorCode ierr; 2742 2743 PetscFunctionBegin; 2744 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 2745 for(f = 0; f < rows; ++f) { 2746 ierr = PetscPrintf(PETSC_COMM_SELF, " |");CHKERRQ(ierr); 2747 for(g = 0; g < cols; ++g) { 2748 ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5G", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr); 2749 } 2750 ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr); 2751 } 2752 PetscFunctionReturn(0); 2753 } 2754 2755 #undef __FUNCT__ 2756 #define __FUNCT__ "DMGetLocalFunction" 2757 /*@C 2758 DMGetLocalFunction - Get the local residual function from this DM 2759 2760 Not collective 2761 2762 Input Parameter: 2763 . dm - The DM 2764 2765 Output Parameter: 2766 . lf - The local residual function 2767 2768 Calling sequence of lf: 2769 $ lf (SNES snes, Vec x, Vec f, void *ctx); 2770 2771 + snes - the SNES context 2772 . x - local vector with the state at which to evaluate residual 2773 . f - local vector to put residual in 2774 - ctx - optional user-defined function context 2775 2776 Level: intermediate 2777 2778 .seealso DMSetLocalFunction(), DMGetLocalJacobian(), DMSetLocalJacobian() 2779 @*/ 2780 PetscErrorCode DMGetLocalFunction(DM dm, PetscErrorCode (**lf)(DM, Vec, Vec, void *)) 2781 { 2782 PetscFunctionBegin; 2783 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2784 if (lf) *lf = dm->lf; 2785 PetscFunctionReturn(0); 2786 } 2787 2788 #undef __FUNCT__ 2789 #define __FUNCT__ "DMSetLocalFunction" 2790 /*@C 2791 DMSetLocalFunction - Set the local residual function from this DM 2792 2793 Not collective 2794 2795 Input Parameters: 2796 + dm - The DM 2797 - lf - The local residual function 2798 2799 Calling sequence of lf: 2800 $ lf (SNES snes, Vec x, Vec f, void *ctx); 2801 2802 + snes - the SNES context 2803 . x - local vector with the state at which to evaluate residual 2804 . f - local vector to put residual in 2805 - ctx - optional user-defined function context 2806 2807 Level: intermediate 2808 2809 .seealso DMGetLocalFunction(), DMGetLocalJacobian(), DMSetLocalJacobian() 2810 @*/ 2811 PetscErrorCode DMSetLocalFunction(DM dm, PetscErrorCode (*lf)(DM, Vec, Vec, void *)) 2812 { 2813 PetscFunctionBegin; 2814 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2815 dm->lf = lf; 2816 PetscFunctionReturn(0); 2817 } 2818 2819 #undef __FUNCT__ 2820 #define __FUNCT__ "DMGetLocalJacobian" 2821 /*@C 2822 DMGetLocalJacobian - Get the local Jacobian function from this DM 2823 2824 Not collective 2825 2826 Input Parameter: 2827 . dm - The DM 2828 2829 Output Parameter: 2830 . lj - The local Jacobian function 2831 2832 Calling sequence of lj: 2833 $ lj (SNES snes, Vec x, Mat J, Mat M, void *ctx); 2834 2835 + snes - the SNES context 2836 . x - local vector with the state at which to evaluate residual 2837 . J - matrix to put Jacobian in 2838 . M - matrix to use for defining Jacobian preconditioner 2839 - ctx - optional user-defined function context 2840 2841 Level: intermediate 2842 2843 .seealso DMSetLocalJacobian(), DMGetLocalFunction(), DMSetLocalFunction() 2844 @*/ 2845 PetscErrorCode DMGetLocalJacobian(DM dm, PetscErrorCode (**lj)(DM, Vec, Mat, Mat, void *)) 2846 { 2847 PetscFunctionBegin; 2848 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2849 if (lj) *lj = dm->lj; 2850 PetscFunctionReturn(0); 2851 } 2852 2853 #undef __FUNCT__ 2854 #define __FUNCT__ "DMSetLocalJacobian" 2855 /*@C 2856 DMSetLocalJacobian - Set the local Jacobian function from this DM 2857 2858 Not collective 2859 2860 Input Parameters: 2861 + dm - The DM 2862 - lj - The local Jacobian function 2863 2864 Calling sequence of lj: 2865 $ lj (SNES snes, Vec x, Mat J, Mat M, void *ctx); 2866 2867 + snes - the SNES context 2868 . x - local vector with the state at which to evaluate residual 2869 . J - matrix to put Jacobian in 2870 . M - matrix to use for defining Jacobian preconditioner 2871 - ctx - optional user-defined function context 2872 2873 Level: intermediate 2874 2875 .seealso DMGetLocalJacobian(), DMGetLocalFunction(), DMSetLocalFunction() 2876 @*/ 2877 PetscErrorCode DMSetLocalJacobian(DM dm, PetscErrorCode (*lj)(DM, Vec, Mat, Mat, void *)) 2878 { 2879 PetscFunctionBegin; 2880 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2881 dm->lj = lj; 2882 PetscFunctionReturn(0); 2883 } 2884 2885 #undef __FUNCT__ 2886 #define __FUNCT__ "DMGetDefaultSection" 2887 /*@ 2888 DMGetDefaultSection - Get the PetscSection encoding the local data layout for the DM. 2889 2890 Input Parameter: 2891 . dm - The DM 2892 2893 Output Parameter: 2894 . section - The PetscSection 2895 2896 Level: intermediate 2897 2898 Note: This gets a borrowed reference, so the user should not destroy this PetscSection. 2899 2900 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection() 2901 @*/ 2902 PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *section) { 2903 PetscFunctionBegin; 2904 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2905 PetscValidPointer(section, 2); 2906 *section = dm->defaultSection; 2907 PetscFunctionReturn(0); 2908 } 2909 2910 #undef __FUNCT__ 2911 #define __FUNCT__ "DMSetDefaultSection" 2912 /*@ 2913 DMSetDefaultSection - Set the PetscSection encoding the local data layout for the DM. 2914 2915 Input Parameters: 2916 + dm - The DM 2917 - section - The PetscSection 2918 2919 Level: intermediate 2920 2921 Note: Any existing Section will be destroyed 2922 2923 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection() 2924 @*/ 2925 PetscErrorCode DMSetDefaultSection(DM dm, PetscSection section) { 2926 PetscErrorCode ierr; 2927 2928 PetscFunctionBegin; 2929 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2930 ierr = PetscSectionDestroy(&dm->defaultSection);CHKERRQ(ierr); 2931 ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr); 2932 dm->defaultSection = section; 2933 PetscFunctionReturn(0); 2934 } 2935 2936 #undef __FUNCT__ 2937 #define __FUNCT__ "DMGetDefaultGlobalSection" 2938 /*@ 2939 DMGetDefaultGlobalSection - Get the PetscSection encoding the global data layout for the DM. 2940 2941 Input Parameter: 2942 . dm - The DM 2943 2944 Output Parameter: 2945 . section - The PetscSection 2946 2947 Level: intermediate 2948 2949 Note: This gets a borrowed reference, so the user should not destroy this PetscSection. 2950 2951 .seealso: DMSetDefaultSection(), DMGetDefaultSection() 2952 @*/ 2953 PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *section) { 2954 PetscErrorCode ierr; 2955 2956 PetscFunctionBegin; 2957 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2958 PetscValidPointer(section, 2); 2959 if (!dm->defaultGlobalSection) { 2960 ierr = PetscSectionCreateGlobalSection(dm->defaultSection, dm->sf, &dm->defaultGlobalSection);CHKERRQ(ierr); 2961 } 2962 *section = dm->defaultGlobalSection; 2963 PetscFunctionReturn(0); 2964 } 2965 2966 #undef __FUNCT__ 2967 #define __FUNCT__ "DMGetDefaultSF" 2968 /*@ 2969 DMGetDefaultSF - Get the PetscSF encoding the parallel dof overlap for the DM. If it has not been set, 2970 it is created from the default PetscSection layouts in the DM. 2971 2972 Input Parameter: 2973 . dm - The DM 2974 2975 Output Parameter: 2976 . sf - The PetscSF 2977 2978 Level: intermediate 2979 2980 Note: This gets a borrowed reference, so the user should not destroy this PetscSF. 2981 2982 .seealso: DMSetDefaultSF(), DMCreateDefaultSF() 2983 @*/ 2984 PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *sf) { 2985 PetscInt nroots; 2986 PetscErrorCode ierr; 2987 2988 PetscFunctionBegin; 2989 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2990 PetscValidPointer(sf, 2); 2991 ierr = PetscSFGetGraph(dm->defaultSF, &nroots, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 2992 if (nroots < 0) { 2993 PetscSection section, gSection; 2994 2995 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 2996 if (section) { 2997 ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr); 2998 ierr = DMCreateDefaultSF(dm, section, gSection);CHKERRQ(ierr); 2999 } else { 3000 *sf = PETSC_NULL; 3001 PetscFunctionReturn(0); 3002 } 3003 } 3004 *sf = dm->defaultSF; 3005 PetscFunctionReturn(0); 3006 } 3007 3008 #undef __FUNCT__ 3009 #define __FUNCT__ "DMSetDefaultSF" 3010 /*@ 3011 DMSetDefaultSF - Set the PetscSF encoding the parallel dof overlap for the DM 3012 3013 Input Parameters: 3014 + dm - The DM 3015 - sf - The PetscSF 3016 3017 Level: intermediate 3018 3019 Note: Any previous SF is destroyed 3020 3021 .seealso: DMGetDefaultSF(), DMCreateDefaultSF() 3022 @*/ 3023 PetscErrorCode DMSetDefaultSF(DM dm, PetscSF sf) { 3024 PetscErrorCode ierr; 3025 3026 PetscFunctionBegin; 3027 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3028 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 3029 ierr = PetscSFDestroy(&dm->defaultSF);CHKERRQ(ierr); 3030 dm->defaultSF = sf; 3031 PetscFunctionReturn(0); 3032 } 3033 3034 #undef __FUNCT__ 3035 #define __FUNCT__ "DMCreateDefaultSF" 3036 /*@C 3037 DMCreateDefaultSF - Create the PetscSF encoding the parallel dof overlap for the DM based upon the PetscSections 3038 describing the data layout. 3039 3040 Input Parameters: 3041 + dm - The DM 3042 . localSection - PetscSection describing the local data layout 3043 - globalSection - PetscSection describing the global data layout 3044 3045 Level: intermediate 3046 3047 .seealso: DMGetDefaultSF(), DMSetDefaultSF() 3048 @*/ 3049 PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection localSection, PetscSection globalSection) 3050 { 3051 MPI_Comm comm = ((PetscObject) dm)->comm; 3052 PetscLayout layout; 3053 const PetscInt *ranges; 3054 PetscInt *local; 3055 PetscSFNode *remote; 3056 PetscInt pStart, pEnd, p, nroots, nleaves, l; 3057 PetscMPIInt size, rank; 3058 PetscErrorCode ierr; 3059 3060 PetscFunctionBegin; 3061 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3062 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 3063 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 3064 ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr); 3065 ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr); 3066 ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 3067 ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr); 3068 ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr); 3069 ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 3070 ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr); 3071 for(p = pStart, nleaves = 0; p < pEnd; ++p) { 3072 PetscInt dof, cdof; 3073 3074 ierr = PetscSectionGetDof(globalSection, p, &dof);CHKERRQ(ierr); 3075 ierr = PetscSectionGetConstraintDof(globalSection, p, &cdof);CHKERRQ(ierr); 3076 nleaves += dof < 0 ? -(dof+1)-cdof : dof-cdof; 3077 } 3078 ierr = PetscMalloc(nleaves * sizeof(PetscInt), &local);CHKERRQ(ierr); 3079 ierr = PetscMalloc(nleaves * sizeof(PetscSFNode), &remote);CHKERRQ(ierr); 3080 for(p = pStart, l = 0; p < pEnd; ++p) { 3081 PetscInt *cind; 3082 PetscInt dof, gdof, cdof, dim, off, goff, d, c; 3083 3084 ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr); 3085 ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr); 3086 ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr); 3087 ierr = PetscSectionGetConstraintIndices(localSection, p, &cind);CHKERRQ(ierr); 3088 ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr); 3089 ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr); 3090 dim = dof-cdof; 3091 for(d = 0, c = 0; d < dof; ++d) { 3092 if ((c < cdof) && (cind[c] == d)) {++c; continue;} 3093 local[l+d-c] = off+d; 3094 } 3095 if (gdof < 0) { 3096 for(d = 0; d < dim; ++d, ++l) { 3097 PetscInt offset = -(goff+1) + d, r; 3098 3099 for(r = 0; r < size; ++r) { 3100 if ((offset >= ranges[r]) && (offset < ranges[r+1])) break; 3101 } 3102 remote[l].rank = r; 3103 remote[l].index = offset - ranges[r]; 3104 } 3105 } else { 3106 for(d = 0; d < dim; ++d, ++l) { 3107 remote[l].rank = rank; 3108 remote[l].index = goff+d - ranges[rank]; 3109 } 3110 } 3111 } 3112 ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 3113 ierr = PetscSFSetGraph(dm->defaultSF, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr); 3114 PetscFunctionReturn(0); 3115 } 3116