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 ierr = PetscFree((*dm)->data);CHKERRQ(ierr); 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__ "DMCreateFieldIS" 824 /*@C 825 DMCreateFieldIS - Creates a set of IS objects with the global indices of dofs for each field 826 827 Not collective 828 829 Input Parameter: 830 . dm - the DM object 831 832 Output Parameters: 833 + numFields - The number of fields (or PETSC_NULL if not requested) 834 . fieldNames - The name for each field (or PETSC_NULL if not requested) 835 - fields - The global indices for each field (or PETSC_NULL if not requested) 836 837 Level: intermediate 838 839 Notes: 840 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 841 PetscFree(), every entry of fields should be destroyed with ISDestroy(), and both arrays should be freed with 842 PetscFree(). 843 844 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 845 @*/ 846 PetscErrorCode DMCreateFieldIS(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields) 847 { 848 PetscSection section, sectionGlobal; 849 PetscErrorCode ierr; 850 851 PetscFunctionBegin; 852 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 853 if (numFields) { 854 PetscValidPointer(numFields,2); 855 *numFields = 0; 856 } 857 if (fieldNames) { 858 PetscValidPointer(fieldNames,3); 859 *fieldNames = PETSC_NULL; 860 } 861 if (fields) { 862 PetscValidPointer(fields,4); 863 *fields = PETSC_NULL; 864 } 865 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 866 if (section) { 867 PetscInt *fieldSizes, **fieldIndices; 868 PetscInt nF, f, pStart, pEnd, p; 869 870 ierr = DMGetDefaultGlobalSection(dm, §ionGlobal);CHKERRQ(ierr); 871 ierr = PetscSectionGetNumFields(section, &nF);CHKERRQ(ierr); 872 ierr = PetscMalloc2(nF,PetscInt,&fieldSizes,nF,PetscInt *,&fieldIndices);CHKERRQ(ierr); 873 ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr); 874 for(f = 0; f < nF; ++f) { 875 fieldSizes[f] = 0; 876 } 877 for(p = pStart; p < pEnd; ++p) { 878 PetscInt gdof; 879 880 ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr); 881 if (gdof > 0) { 882 for(f = 0; f < nF; ++f) { 883 PetscInt fdof, fcdof; 884 885 ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr); 886 ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr); 887 fieldSizes[f] += fdof-fcdof; 888 } 889 } 890 } 891 for(f = 0; f < nF; ++f) { 892 ierr = PetscMalloc(fieldSizes[f] * sizeof(PetscInt), &fieldIndices[f]);CHKERRQ(ierr); 893 fieldSizes[f] = 0; 894 } 895 for(p = pStart; p < pEnd; ++p) { 896 PetscInt gdof, goff; 897 898 ierr = PetscSectionGetDof(sectionGlobal, p, &gdof);CHKERRQ(ierr); 899 if (gdof > 0) { 900 ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 901 for(f = 0; f < nF; ++f) { 902 PetscInt fdof, fcdof, fc; 903 904 ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr); 905 ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr); 906 for(fc = 0; fc < fdof-fcdof; ++fc, ++fieldSizes[f]) { 907 fieldIndices[f][fieldSizes[f]] = goff++; 908 } 909 } 910 } 911 } 912 if (numFields) {*numFields = nF;} 913 if (fieldNames) { 914 ierr = PetscMalloc(nF * sizeof(char *), fieldNames);CHKERRQ(ierr); 915 for(f = 0; f < nF; ++f) { 916 const char *fieldName; 917 918 ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr); 919 ierr = PetscStrallocpy(fieldName, (char **) &(*fieldNames)[f]);CHKERRQ(ierr); 920 } 921 } 922 if (fields) { 923 ierr = PetscMalloc(nF * sizeof(IS), fields);CHKERRQ(ierr); 924 for(f = 0; f < nF; ++f) { 925 ierr = ISCreateGeneral(((PetscObject) dm)->comm, fieldSizes[f], fieldIndices[f], PETSC_OWN_POINTER, &(*fields)[f]);CHKERRQ(ierr); 926 } 927 } 928 ierr = PetscFree2(fieldSizes,fieldIndices);CHKERRQ(ierr); 929 } else { 930 if(dm->ops->createfieldis) {ierr = (*dm->ops->createfieldis)(dm, numFields, fieldNames, fields);CHKERRQ(ierr);} 931 } 932 PetscFunctionReturn(0); 933 } 934 935 936 #undef __FUNCT__ 937 #define __FUNCT__ "DMCreateFieldDecompositionDM" 938 /*@C 939 DMCreateFieldDecompositionDM - creates a DM that encapsulates a decomposition of the original DM into fields. 940 941 Not Collective 942 943 Input Parameters: 944 + dm - the DM object 945 - name - the name of the field decomposition 946 947 Output Parameter: 948 . ddm - the field decomposition DM (PETSC_NULL, if no such decomposition is known) 949 950 Level: advanced 951 952 .seealso DMDestroy(), DMCreate(), DMCreateFieldDecomposition(), DMCreateDomainDecompositionDM() 953 @*/ 954 PetscErrorCode DMCreateFieldDecompositionDM(DM dm, const char* name, DM *ddm) 955 { 956 PetscErrorCode ierr; 957 958 PetscFunctionBegin; 959 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 960 PetscValidCharPointer(name,2); 961 PetscValidPointer(ddm,3); 962 *ddm = PETSC_NULL; 963 if(dm->ops->createfielddecompositiondm) { 964 ierr = (*dm->ops->createfielddecompositiondm)(dm,name,ddm); CHKERRQ(ierr); 965 } 966 PetscFunctionReturn(0); 967 } 968 969 970 #undef __FUNCT__ 971 #define __FUNCT__ "DMCreateFieldDecomposition" 972 /*@C 973 DMCreateFieldDecomposition - Returns a list of IS objects defining a decomposition of a problem into subproblems 974 corresponding to different fields: each IS contains the global indices of the dofs of the 975 corresponding field. The optional list of DMs define the DM for each subproblem. 976 Generalizes DMCreateFieldIS(). 977 978 Not collective 979 980 Input Parameter: 981 . dm - the DM object 982 983 Output Parameters: 984 + len - The number of subproblems in the field decomposition (or PETSC_NULL if not requested) 985 . namelist - The name for each field (or PETSC_NULL if not requested) 986 . islist - The global indices for each field (or PETSC_NULL if not requested) 987 - dmlist - The DMs for each field subproblem (or PETSC_NULL, if not requested; if PETSC_NULL is returned, no DMs are defined) 988 989 Level: intermediate 990 991 Notes: 992 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 993 PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(), 994 and all of the arrays should be freed with PetscFree(). 995 996 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateFieldIS() 997 @*/ 998 PetscErrorCode DMCreateFieldDecomposition(DM dm, PetscInt *len, char ***namelist, IS **islist, DM **dmlist) 999 { 1000 PetscErrorCode ierr; 1001 1002 PetscFunctionBegin; 1003 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1004 if (len) {PetscValidPointer(len,2); *len = 0;} 1005 if (namelist) {PetscValidPointer(namelist,3); *namelist = 0;} 1006 if (islist) {PetscValidPointer(islist,4); *islist = 0;} 1007 if (dmlist) {PetscValidPointer(dmlist,5); *dmlist = 0;} 1008 if(!dm->ops->createfielddecomposition) { 1009 ierr = DMCreateFieldIS(dm, len, namelist, islist);CHKERRQ(ierr); 1010 /* By default there are no DMs associated with subproblems. */ 1011 if(dmlist) *dmlist = PETSC_NULL; 1012 } 1013 else { 1014 ierr = (*dm->ops->createfielddecomposition)(dm,len,namelist,islist,dmlist); CHKERRQ(ierr); 1015 } 1016 PetscFunctionReturn(0); 1017 } 1018 1019 #undef __FUNCT__ 1020 #define __FUNCT__ "DMCreateDomainDecompositionDM" 1021 /*@C 1022 DMCreateDomainDecompositionDM - creates a DM that encapsulates a decomposition of the original DM into subdomains. 1023 1024 Not Collective 1025 1026 Input Parameters: 1027 + dm - the DM object 1028 - name - the name of the subdomain decomposition 1029 1030 Output Parameter: 1031 . ddm - the subdomain decomposition DM (PETSC_NULL, if no such decomposition is known) 1032 1033 Level: advanced 1034 1035 .seealso DMDestroy(), DMCreate(), DMCreateFieldDecomposition(), DMCreateDomainDecompositionDM() 1036 @*/ 1037 PetscErrorCode DMCreateDomainDecompositionDM(DM dm, const char* name, DM *ddm) 1038 { 1039 PetscErrorCode ierr; 1040 1041 PetscFunctionBegin; 1042 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1043 PetscValidCharPointer(name,2); 1044 PetscValidPointer(ddm,3); 1045 *ddm = PETSC_NULL; 1046 if(dm->ops->createdomaindecompositiondm) { 1047 ierr = (*dm->ops->createdomaindecompositiondm)(dm,name,ddm); CHKERRQ(ierr); 1048 } 1049 PetscFunctionReturn(0); 1050 } 1051 1052 1053 #undef __FUNCT__ 1054 #define __FUNCT__ "DMCreateDomainDecomposition" 1055 /*@C 1056 DMCreateDomainDecomposition - Returns lists of IS objects defining a decomposition of a problem into subproblems 1057 corresponding to restrictions to pairs nested subdomains: each IS contains the global 1058 indices of the dofs of the corresponding subdomains. The inner subdomains conceptually 1059 define a nonoverlapping covering, while outer subdomains can overlap. 1060 The optional list of DMs define the DM for each subproblem. 1061 1062 Not collective 1063 1064 Input Parameter: 1065 . dm - the DM object 1066 1067 Output Parameters: 1068 + len - The number of subproblems in the domain decomposition (or PETSC_NULL if not requested) 1069 . namelist - The name for each subdomain (or PETSC_NULL if not requested) 1070 . innerislist - The global indices for each inner subdomain (or PETSC_NULL, if not requested) 1071 . outerislist - The global indices for each outer subdomain (or PETSC_NULL, if not requested) 1072 - dmlist - The DMs for each subdomain subproblem (or PETSC_NULL, if not requested; if PETSC_NULL is returned, no DMs are defined) 1073 1074 Level: intermediate 1075 1076 Notes: 1077 The user is responsible for freeing all requested arrays. In particular, every entry of names should be freed with 1078 PetscFree(), every entry of is should be destroyed with ISDestroy(), every entry of dm should be destroyed with DMDestroy(), 1079 and all of the arrays should be freed with PetscFree(). 1080 1081 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMCreateDomainDecompositionDM(), DMCreateFieldDecomposition() 1082 @*/ 1083 PetscErrorCode DMCreateDomainDecomposition(DM dm, PetscInt *len, char ***namelist, IS **innerislist, IS **outerislist, DM **dmlist) 1084 { 1085 PetscErrorCode ierr; 1086 1087 PetscFunctionBegin; 1088 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1089 if (len) {PetscValidPointer(len,2); *len = PETSC_NULL;} 1090 if (namelist) {PetscValidPointer(namelist,3); *namelist = PETSC_NULL;} 1091 if (innerislist) {PetscValidPointer(innerislist,4); *innerislist = PETSC_NULL;} 1092 if (outerislist) {PetscValidPointer(outerislist,5); *outerislist = PETSC_NULL;} 1093 if (dmlist) {PetscValidPointer(dmlist,6); *dmlist = PETSC_NULL;} 1094 if(dm->ops->createdomaindecomposition) { 1095 ierr = (*dm->ops->createdomaindecomposition)(dm,len,namelist,innerislist,outerislist,dmlist); CHKERRQ(ierr); 1096 } 1097 PetscFunctionReturn(0); 1098 } 1099 1100 #undef __FUNCT__ 1101 #define __FUNCT__ "DMRefine" 1102 /*@ 1103 DMRefine - Refines a DM object 1104 1105 Collective on DM 1106 1107 Input Parameter: 1108 + dm - the DM object 1109 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 1110 1111 Output Parameter: 1112 . dmf - the refined DM, or PETSC_NULL 1113 1114 Note: If no refinement was done, the return value is PETSC_NULL 1115 1116 Level: developer 1117 1118 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1119 @*/ 1120 PetscErrorCode DMRefine(DM dm,MPI_Comm comm,DM *dmf) 1121 { 1122 PetscErrorCode ierr; 1123 DMRefineHookLink link; 1124 1125 PetscFunctionBegin; 1126 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1127 ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr); 1128 if (*dmf) { 1129 (*dmf)->ops->creatematrix = dm->ops->creatematrix; 1130 (*dmf)->ops->initialguess = dm->ops->initialguess; 1131 (*dmf)->ops->function = dm->ops->function; 1132 (*dmf)->ops->functionj = dm->ops->functionj; 1133 if (dm->ops->jacobian != DMComputeJacobianDefault) { 1134 (*dmf)->ops->jacobian = dm->ops->jacobian; 1135 } 1136 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr); 1137 (*dmf)->ctx = dm->ctx; 1138 (*dmf)->leveldown = dm->leveldown; 1139 (*dmf)->levelup = dm->levelup + 1; 1140 ierr = DMSetMatType(*dmf,dm->mattype);CHKERRQ(ierr); 1141 for (link=dm->refinehook; link; link=link->next) { 1142 if (link->refinehook) {ierr = (*link->refinehook)(dm,*dmf,link->ctx);CHKERRQ(ierr);} 1143 } 1144 } 1145 PetscFunctionReturn(0); 1146 } 1147 1148 #undef __FUNCT__ 1149 #define __FUNCT__ "DMRefineHookAdd" 1150 /*@ 1151 DMRefineHookAdd - adds a callback to be run when interpolating a nonlinear problem to a finer grid 1152 1153 Logically Collective 1154 1155 Input Arguments: 1156 + coarse - nonlinear solver context on which to run a hook when restricting to a coarser level 1157 . refinehook - function to run when setting up a coarser level 1158 . interphook - function to run to update data on finer levels (once per SNESSolve()) 1159 - ctx - [optional] user-defined context for provide data for the hooks (may be PETSC_NULL) 1160 1161 Calling sequence of refinehook: 1162 $ refinehook(DM coarse,DM fine,void *ctx); 1163 1164 + coarse - coarse level DM 1165 . fine - fine level DM to interpolate problem to 1166 - ctx - optional user-defined function context 1167 1168 Calling sequence for interphook: 1169 $ interphook(DM coarse,Mat interp,DM fine,void *ctx) 1170 1171 + coarse - coarse level DM 1172 . interp - matrix interpolating a coarse-level solution to the finer grid 1173 . fine - fine level DM to update 1174 - ctx - optional user-defined function context 1175 1176 Level: advanced 1177 1178 Notes: 1179 This function is only needed if auxiliary data needs to be passed to fine grids while grid sequencing 1180 1181 If this function is called multiple times, the hooks will be run in the order they are added. 1182 1183 .seealso: DMCoarsenHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 1184 @*/ 1185 PetscErrorCode DMRefineHookAdd(DM coarse,PetscErrorCode (*refinehook)(DM,DM,void*),PetscErrorCode (*interphook)(DM,Mat,DM,void*),void *ctx) 1186 { 1187 PetscErrorCode ierr; 1188 DMRefineHookLink link,*p; 1189 1190 PetscFunctionBegin; 1191 PetscValidHeaderSpecific(coarse,DM_CLASSID,1); 1192 for (p=&coarse->refinehook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */ 1193 ierr = PetscMalloc(sizeof(struct _DMRefineHookLink),&link);CHKERRQ(ierr); 1194 link->refinehook = refinehook; 1195 link->interphook = interphook; 1196 link->ctx = ctx; 1197 link->next = PETSC_NULL; 1198 *p = link; 1199 PetscFunctionReturn(0); 1200 } 1201 1202 #undef __FUNCT__ 1203 #define __FUNCT__ "DMInterpolate" 1204 /*@ 1205 DMInterpolate - interpolates user-defined problem data to a finer DM by running hooks registered by DMRefineHookAdd() 1206 1207 Collective if any hooks are 1208 1209 Input Arguments: 1210 + coarse - coarser DM to use as a base 1211 . restrct - interpolation matrix, apply using MatInterpolate() 1212 - fine - finer DM to update 1213 1214 Level: developer 1215 1216 .seealso: DMRefineHookAdd(), MatInterpolate() 1217 @*/ 1218 PetscErrorCode DMInterpolate(DM coarse,Mat interp,DM fine) 1219 { 1220 PetscErrorCode ierr; 1221 DMRefineHookLink link; 1222 1223 PetscFunctionBegin; 1224 for (link=fine->refinehook; link; link=link->next) { 1225 if (link->interphook) {ierr = (*link->interphook)(coarse,interp,fine,link->ctx);CHKERRQ(ierr);} 1226 } 1227 PetscFunctionReturn(0); 1228 } 1229 1230 #undef __FUNCT__ 1231 #define __FUNCT__ "DMGetRefineLevel" 1232 /*@ 1233 DMGetRefineLevel - Get's the number of refinements that have generated this DM. 1234 1235 Not Collective 1236 1237 Input Parameter: 1238 . dm - the DM object 1239 1240 Output Parameter: 1241 . level - number of refinements 1242 1243 Level: developer 1244 1245 .seealso DMCoarsen(), DMGetCoarsenLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1246 1247 @*/ 1248 PetscErrorCode DMGetRefineLevel(DM dm,PetscInt *level) 1249 { 1250 PetscFunctionBegin; 1251 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1252 *level = dm->levelup; 1253 PetscFunctionReturn(0); 1254 } 1255 1256 #undef __FUNCT__ 1257 #define __FUNCT__ "DMGlobalToLocalBegin" 1258 /*@ 1259 DMGlobalToLocalBegin - Begins updating local vectors from global vector 1260 1261 Neighbor-wise Collective on DM 1262 1263 Input Parameters: 1264 + dm - the DM object 1265 . g - the global vector 1266 . mode - INSERT_VALUES or ADD_VALUES 1267 - l - the local vector 1268 1269 1270 Level: beginner 1271 1272 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 1273 1274 @*/ 1275 PetscErrorCode DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) 1276 { 1277 PetscSF sf; 1278 PetscErrorCode ierr; 1279 1280 PetscFunctionBegin; 1281 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1282 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 1283 if (sf) { 1284 PetscScalar *lArray, *gArray; 1285 1286 if (mode == ADD_VALUES) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 1287 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 1288 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 1289 ierr = PetscSFBcastBegin(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr); 1290 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 1291 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 1292 } else { 1293 ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 1294 } 1295 PetscFunctionReturn(0); 1296 } 1297 1298 #undef __FUNCT__ 1299 #define __FUNCT__ "DMGlobalToLocalEnd" 1300 /*@ 1301 DMGlobalToLocalEnd - Ends updating local vectors from global vector 1302 1303 Neighbor-wise Collective on DM 1304 1305 Input Parameters: 1306 + dm - the DM object 1307 . g - the global vector 1308 . mode - INSERT_VALUES or ADD_VALUES 1309 - l - the local vector 1310 1311 1312 Level: beginner 1313 1314 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 1315 1316 @*/ 1317 PetscErrorCode DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) 1318 { 1319 PetscSF sf; 1320 PetscErrorCode ierr; 1321 PetscScalar *lArray, *gArray; 1322 1323 PetscFunctionBegin; 1324 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1325 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 1326 if (sf) { 1327 if (mode == ADD_VALUES) SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 1328 1329 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 1330 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 1331 ierr = PetscSFBcastEnd(sf, MPIU_SCALAR, gArray, lArray);CHKERRQ(ierr); 1332 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 1333 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 1334 } else { 1335 ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 1336 } 1337 PetscFunctionReturn(0); 1338 } 1339 1340 #undef __FUNCT__ 1341 #define __FUNCT__ "DMLocalToGlobalBegin" 1342 /*@ 1343 DMLocalToGlobalBegin - updates global vectors from local vectors 1344 1345 Neighbor-wise Collective on DM 1346 1347 Input Parameters: 1348 + dm - the DM object 1349 . l - the local vector 1350 . 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 1351 base point. 1352 - - the global vector 1353 1354 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 1355 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 1356 global array to the final global array with VecAXPY(). 1357 1358 Level: beginner 1359 1360 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin() 1361 1362 @*/ 1363 PetscErrorCode DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g) 1364 { 1365 PetscSF sf; 1366 PetscErrorCode ierr; 1367 1368 PetscFunctionBegin; 1369 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1370 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 1371 if (sf) { 1372 MPI_Op op; 1373 PetscScalar *lArray, *gArray; 1374 1375 switch(mode) { 1376 case INSERT_VALUES: 1377 case INSERT_ALL_VALUES: 1378 #if defined(PETSC_HAVE_MPI_REPLACE) 1379 op = MPI_REPLACE; break; 1380 #else 1381 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No support for INSERT_VALUES without an MPI-2 implementation"); 1382 #endif 1383 case ADD_VALUES: 1384 case ADD_ALL_VALUES: 1385 op = MPI_SUM; break; 1386 default: 1387 SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 1388 } 1389 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 1390 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 1391 ierr = PetscSFReduceBegin(sf, MPIU_SCALAR, lArray, gArray, op);CHKERRQ(ierr); 1392 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 1393 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 1394 } else { 1395 ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 1396 } 1397 PetscFunctionReturn(0); 1398 } 1399 1400 #undef __FUNCT__ 1401 #define __FUNCT__ "DMLocalToGlobalEnd" 1402 /*@ 1403 DMLocalToGlobalEnd - updates global vectors from local vectors 1404 1405 Neighbor-wise Collective on DM 1406 1407 Input Parameters: 1408 + dm - the DM object 1409 . l - the local vector 1410 . mode - INSERT_VALUES or ADD_VALUES 1411 - g - the global vector 1412 1413 1414 Level: beginner 1415 1416 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd() 1417 1418 @*/ 1419 PetscErrorCode DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g) 1420 { 1421 PetscSF sf; 1422 PetscErrorCode ierr; 1423 1424 PetscFunctionBegin; 1425 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1426 ierr = DMGetDefaultSF(dm, &sf);CHKERRQ(ierr); 1427 if (sf) { 1428 MPI_Op op; 1429 PetscScalar *lArray, *gArray; 1430 1431 switch(mode) { 1432 case INSERT_VALUES: 1433 case INSERT_ALL_VALUES: 1434 #if defined(PETSC_HAVE_MPI_REPLACE) 1435 op = MPI_REPLACE; break; 1436 #else 1437 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No support for INSERT_VALUES without an MPI-2 implementation"); 1438 #endif 1439 case ADD_VALUES: 1440 case ADD_ALL_VALUES: 1441 op = MPI_SUM; break; 1442 default: 1443 SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid insertion mode %D", mode); 1444 } 1445 ierr = VecGetArray(l, &lArray);CHKERRQ(ierr); 1446 ierr = VecGetArray(g, &gArray);CHKERRQ(ierr); 1447 ierr = PetscSFReduceEnd(sf, MPIU_SCALAR, lArray, gArray, op);CHKERRQ(ierr); 1448 ierr = VecRestoreArray(l, &lArray);CHKERRQ(ierr); 1449 ierr = VecRestoreArray(g, &gArray);CHKERRQ(ierr); 1450 } else { 1451 ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 1452 } 1453 PetscFunctionReturn(0); 1454 } 1455 1456 #undef __FUNCT__ 1457 #define __FUNCT__ "DMComputeJacobianDefault" 1458 /*@ 1459 DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided 1460 1461 Collective on DM 1462 1463 Input Parameter: 1464 + dm - the DM object 1465 . x - location to compute Jacobian at; may be ignored for linear problems 1466 . A - matrix that defines the operator for the linear solve 1467 - B - the matrix used to construct the preconditioner 1468 1469 Level: developer 1470 1471 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1472 DMSetFunction() 1473 1474 @*/ 1475 PetscErrorCode DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 1476 { 1477 PetscErrorCode ierr; 1478 1479 PetscFunctionBegin; 1480 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1481 *stflag = SAME_NONZERO_PATTERN; 1482 ierr = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr); 1483 if (A != B) { 1484 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1485 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1486 } 1487 PetscFunctionReturn(0); 1488 } 1489 1490 #undef __FUNCT__ 1491 #define __FUNCT__ "DMCoarsen" 1492 /*@ 1493 DMCoarsen - Coarsens a DM object 1494 1495 Collective on DM 1496 1497 Input Parameter: 1498 + dm - the DM object 1499 - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 1500 1501 Output Parameter: 1502 . dmc - the coarsened DM 1503 1504 Level: developer 1505 1506 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1507 1508 @*/ 1509 PetscErrorCode DMCoarsen(DM dm, MPI_Comm comm, DM *dmc) 1510 { 1511 PetscErrorCode ierr; 1512 DMCoarsenHookLink link; 1513 1514 PetscFunctionBegin; 1515 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1516 ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr); 1517 (*dmc)->ops->creatematrix = dm->ops->creatematrix; 1518 (*dmc)->ops->initialguess = dm->ops->initialguess; 1519 (*dmc)->ops->function = dm->ops->function; 1520 (*dmc)->ops->functionj = dm->ops->functionj; 1521 if (dm->ops->jacobian != DMComputeJacobianDefault) { 1522 (*dmc)->ops->jacobian = dm->ops->jacobian; 1523 } 1524 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr); 1525 (*dmc)->ctx = dm->ctx; 1526 (*dmc)->levelup = dm->levelup; 1527 (*dmc)->leveldown = dm->leveldown + 1; 1528 ierr = DMSetMatType(*dmc,dm->mattype);CHKERRQ(ierr); 1529 for (link=dm->coarsenhook; link; link=link->next) { 1530 if (link->coarsenhook) {ierr = (*link->coarsenhook)(dm,*dmc,link->ctx);CHKERRQ(ierr);} 1531 } 1532 PetscFunctionReturn(0); 1533 } 1534 1535 #undef __FUNCT__ 1536 #define __FUNCT__ "DMCoarsenHookAdd" 1537 /*@ 1538 DMCoarsenHookAdd - adds a callback to be run when restricting a nonlinear problem to the coarse grid 1539 1540 Logically Collective 1541 1542 Input Arguments: 1543 + fine - nonlinear solver context on which to run a hook when restricting to a coarser level 1544 . coarsenhook - function to run when setting up a coarser level 1545 . restricthook - function to run to update data on coarser levels (once per SNESSolve()) 1546 - ctx - [optional] user-defined context for provide data for the hooks (may be PETSC_NULL) 1547 1548 Calling sequence of coarsenhook: 1549 $ coarsenhook(DM fine,DM coarse,void *ctx); 1550 1551 + fine - fine level DM 1552 . coarse - coarse level DM to restrict problem to 1553 - ctx - optional user-defined function context 1554 1555 Calling sequence for restricthook: 1556 $ restricthook(DM fine,Mat mrestrict,Vec rscale,Mat inject,DM coarse,void *ctx) 1557 1558 + fine - fine level DM 1559 . mrestrict - matrix restricting a fine-level solution to the coarse grid 1560 . rscale - scaling vector for restriction 1561 . inject - matrix restricting by injection 1562 . coarse - coarse level DM to update 1563 - ctx - optional user-defined function context 1564 1565 Level: advanced 1566 1567 Notes: 1568 This function is only needed if auxiliary data needs to be set up on coarse grids. 1569 1570 If this function is called multiple times, the hooks will be run in the order they are added. 1571 1572 In order to compose with nonlinear preconditioning without duplicating storage, the hook should be implemented to 1573 extract the finest level information from its context (instead of from the SNES). 1574 1575 .seealso: DMRefineHookAdd(), SNESFASGetInterpolation(), SNESFASGetInjection(), PetscObjectCompose(), PetscContainerCreate() 1576 @*/ 1577 PetscErrorCode DMCoarsenHookAdd(DM fine,PetscErrorCode (*coarsenhook)(DM,DM,void*),PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*),void *ctx) 1578 { 1579 PetscErrorCode ierr; 1580 DMCoarsenHookLink link,*p; 1581 1582 PetscFunctionBegin; 1583 PetscValidHeaderSpecific(fine,DM_CLASSID,1); 1584 for (p=&fine->coarsenhook; *p; p=&(*p)->next) {} /* Scan to the end of the current list of hooks */ 1585 ierr = PetscMalloc(sizeof(struct _DMCoarsenHookLink),&link);CHKERRQ(ierr); 1586 link->coarsenhook = coarsenhook; 1587 link->restricthook = restricthook; 1588 link->ctx = ctx; 1589 link->next = PETSC_NULL; 1590 *p = link; 1591 PetscFunctionReturn(0); 1592 } 1593 1594 #undef __FUNCT__ 1595 #define __FUNCT__ "DMRestrict" 1596 /*@ 1597 DMRestrict - restricts user-defined problem data to a coarser DM by running hooks registered by DMCoarsenHookAdd() 1598 1599 Collective if any hooks are 1600 1601 Input Arguments: 1602 + fine - finer DM to use as a base 1603 . restrct - restriction matrix, apply using MatRestrict() 1604 . inject - injection matrix, also use MatRestrict() 1605 - coarse - coarer DM to update 1606 1607 Level: developer 1608 1609 .seealso: DMCoarsenHookAdd(), MatRestrict() 1610 @*/ 1611 PetscErrorCode DMRestrict(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse) 1612 { 1613 PetscErrorCode ierr; 1614 DMCoarsenHookLink link; 1615 1616 PetscFunctionBegin; 1617 for (link=fine->coarsenhook; link; link=link->next) { 1618 if (link->restricthook) {ierr = (*link->restricthook)(fine,restrct,rscale,inject,coarse,link->ctx);CHKERRQ(ierr);} 1619 } 1620 PetscFunctionReturn(0); 1621 } 1622 1623 #undef __FUNCT__ 1624 #define __FUNCT__ "DMGetCoarsenLevel" 1625 /*@ 1626 DMGetCoarsenLevel - Get's the number of coarsenings that have generated this DM. 1627 1628 Not Collective 1629 1630 Input Parameter: 1631 . dm - the DM object 1632 1633 Output Parameter: 1634 . level - number of coarsenings 1635 1636 Level: developer 1637 1638 .seealso DMCoarsen(), DMGetRefineLevel(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1639 1640 @*/ 1641 PetscErrorCode DMGetCoarsenLevel(DM dm,PetscInt *level) 1642 { 1643 PetscFunctionBegin; 1644 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1645 *level = dm->leveldown; 1646 PetscFunctionReturn(0); 1647 } 1648 1649 1650 1651 #undef __FUNCT__ 1652 #define __FUNCT__ "DMRefineHierarchy" 1653 /*@C 1654 DMRefineHierarchy - Refines a DM object, all levels at once 1655 1656 Collective on DM 1657 1658 Input Parameter: 1659 + dm - the DM object 1660 - nlevels - the number of levels of refinement 1661 1662 Output Parameter: 1663 . dmf - the refined DM hierarchy 1664 1665 Level: developer 1666 1667 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1668 1669 @*/ 1670 PetscErrorCode DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[]) 1671 { 1672 PetscErrorCode ierr; 1673 1674 PetscFunctionBegin; 1675 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1676 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 1677 if (nlevels == 0) PetscFunctionReturn(0); 1678 if (dm->ops->refinehierarchy) { 1679 ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr); 1680 } else if (dm->ops->refine) { 1681 PetscInt i; 1682 1683 ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr); 1684 for (i=1; i<nlevels; i++) { 1685 ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr); 1686 } 1687 } else { 1688 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet"); 1689 } 1690 PetscFunctionReturn(0); 1691 } 1692 1693 #undef __FUNCT__ 1694 #define __FUNCT__ "DMCoarsenHierarchy" 1695 /*@C 1696 DMCoarsenHierarchy - Coarsens a DM object, all levels at once 1697 1698 Collective on DM 1699 1700 Input Parameter: 1701 + dm - the DM object 1702 - nlevels - the number of levels of coarsening 1703 1704 Output Parameter: 1705 . dmc - the coarsened DM hierarchy 1706 1707 Level: developer 1708 1709 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 1710 1711 @*/ 1712 PetscErrorCode DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[]) 1713 { 1714 PetscErrorCode ierr; 1715 1716 PetscFunctionBegin; 1717 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1718 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 1719 if (nlevels == 0) PetscFunctionReturn(0); 1720 PetscValidPointer(dmc,3); 1721 if (dm->ops->coarsenhierarchy) { 1722 ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr); 1723 } else if (dm->ops->coarsen) { 1724 PetscInt i; 1725 1726 ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr); 1727 for (i=1; i<nlevels; i++) { 1728 ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr); 1729 } 1730 } else { 1731 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet"); 1732 } 1733 PetscFunctionReturn(0); 1734 } 1735 1736 #undef __FUNCT__ 1737 #define __FUNCT__ "DMCreateAggregates" 1738 /*@ 1739 DMCreateAggregates - Gets the aggregates that map between 1740 grids associated with two DMs. 1741 1742 Collective on DM 1743 1744 Input Parameters: 1745 + dmc - the coarse grid DM 1746 - dmf - the fine grid DM 1747 1748 Output Parameters: 1749 . rest - the restriction matrix (transpose of the projection matrix) 1750 1751 Level: intermediate 1752 1753 .keywords: interpolation, restriction, multigrid 1754 1755 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation() 1756 @*/ 1757 PetscErrorCode DMCreateAggregates(DM dmc, DM dmf, Mat *rest) 1758 { 1759 PetscErrorCode ierr; 1760 1761 PetscFunctionBegin; 1762 PetscValidHeaderSpecific(dmc,DM_CLASSID,1); 1763 PetscValidHeaderSpecific(dmf,DM_CLASSID,2); 1764 ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr); 1765 PetscFunctionReturn(0); 1766 } 1767 1768 #undef __FUNCT__ 1769 #define __FUNCT__ "DMSetApplicationContextDestroy" 1770 /*@C 1771 DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed 1772 1773 Not Collective 1774 1775 Input Parameters: 1776 + dm - the DM object 1777 - destroy - the destroy function 1778 1779 Level: intermediate 1780 1781 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1782 1783 @*/ 1784 PetscErrorCode DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**)) 1785 { 1786 PetscFunctionBegin; 1787 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1788 dm->ctxdestroy = destroy; 1789 PetscFunctionReturn(0); 1790 } 1791 1792 #undef __FUNCT__ 1793 #define __FUNCT__ "DMSetApplicationContext" 1794 /*@ 1795 DMSetApplicationContext - Set a user context into a DM object 1796 1797 Not Collective 1798 1799 Input Parameters: 1800 + dm - the DM object 1801 - ctx - the user context 1802 1803 Level: intermediate 1804 1805 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1806 1807 @*/ 1808 PetscErrorCode DMSetApplicationContext(DM dm,void *ctx) 1809 { 1810 PetscFunctionBegin; 1811 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1812 dm->ctx = ctx; 1813 PetscFunctionReturn(0); 1814 } 1815 1816 #undef __FUNCT__ 1817 #define __FUNCT__ "DMGetApplicationContext" 1818 /*@ 1819 DMGetApplicationContext - Gets a user context from a DM object 1820 1821 Not Collective 1822 1823 Input Parameter: 1824 . dm - the DM object 1825 1826 Output Parameter: 1827 . ctx - the user context 1828 1829 Level: intermediate 1830 1831 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1832 1833 @*/ 1834 PetscErrorCode DMGetApplicationContext(DM dm,void *ctx) 1835 { 1836 PetscFunctionBegin; 1837 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1838 *(void**)ctx = dm->ctx; 1839 PetscFunctionReturn(0); 1840 } 1841 1842 #undef __FUNCT__ 1843 #define __FUNCT__ "DMSetInitialGuess" 1844 /*@C 1845 DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers 1846 1847 Logically Collective on DM 1848 1849 Input Parameter: 1850 + dm - the DM object to destroy 1851 - f - the function to compute the initial guess 1852 1853 Level: intermediate 1854 1855 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1856 1857 @*/ 1858 PetscErrorCode DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec)) 1859 { 1860 PetscFunctionBegin; 1861 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1862 dm->ops->initialguess = f; 1863 PetscFunctionReturn(0); 1864 } 1865 1866 #undef __FUNCT__ 1867 #define __FUNCT__ "DMSetFunction" 1868 /*@C 1869 DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES 1870 1871 Logically Collective on DM 1872 1873 Input Parameter: 1874 + dm - the DM object 1875 - f - the function to compute (use PETSC_NULL to cancel a previous function that was set) 1876 1877 Level: intermediate 1878 1879 Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian 1880 computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian. 1881 1882 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1883 DMSetJacobian() 1884 1885 @*/ 1886 PetscErrorCode DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 1887 { 1888 PetscFunctionBegin; 1889 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1890 dm->ops->function = f; 1891 if (f) { 1892 dm->ops->functionj = f; 1893 } 1894 PetscFunctionReturn(0); 1895 } 1896 1897 #undef __FUNCT__ 1898 #define __FUNCT__ "DMSetJacobian" 1899 /*@C 1900 DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES 1901 1902 Logically Collective on DM 1903 1904 Input Parameter: 1905 + dm - the DM object to destroy 1906 - f - the function to compute the matrix entries 1907 1908 Level: intermediate 1909 1910 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1911 DMSetFunction() 1912 1913 @*/ 1914 PetscErrorCode DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*)) 1915 { 1916 PetscFunctionBegin; 1917 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1918 dm->ops->jacobian = f; 1919 PetscFunctionReturn(0); 1920 } 1921 1922 #undef __FUNCT__ 1923 #define __FUNCT__ "DMSetVariableBounds" 1924 /*@C 1925 DMSetVariableBounds - sets a function to compute the the lower and upper bound vectors for SNESVI. 1926 1927 Logically Collective on DM 1928 1929 Input Parameter: 1930 + dm - the DM object 1931 - f - the function that computes variable bounds used by SNESVI (use PETSC_NULL to cancel a previous function that was set) 1932 1933 Level: intermediate 1934 1935 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1936 DMSetJacobian() 1937 1938 @*/ 1939 PetscErrorCode DMSetVariableBounds(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 1940 { 1941 PetscFunctionBegin; 1942 dm->ops->computevariablebounds = f; 1943 PetscFunctionReturn(0); 1944 } 1945 1946 #undef __FUNCT__ 1947 #define __FUNCT__ "DMHasVariableBounds" 1948 /*@ 1949 DMHasVariableBounds - does the DM object have a variable bounds function? 1950 1951 Not Collective 1952 1953 Input Parameter: 1954 . dm - the DM object to destroy 1955 1956 Output Parameter: 1957 . flg - PETSC_TRUE if the variable bounds function exists 1958 1959 Level: developer 1960 1961 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1962 1963 @*/ 1964 PetscErrorCode DMHasVariableBounds(DM dm,PetscBool *flg) 1965 { 1966 PetscFunctionBegin; 1967 *flg = (dm->ops->computevariablebounds) ? PETSC_TRUE : PETSC_FALSE; 1968 PetscFunctionReturn(0); 1969 } 1970 1971 #undef __FUNCT__ 1972 #define __FUNCT__ "DMComputeVariableBounds" 1973 /*@C 1974 DMComputeVariableBounds - compute variable bounds used by SNESVI. 1975 1976 Logically Collective on DM 1977 1978 Input Parameters: 1979 + dm - the DM object to destroy 1980 - x - current solution at which the bounds are computed 1981 1982 Output parameters: 1983 + xl - lower bound 1984 - xu - upper bound 1985 1986 Level: intermediate 1987 1988 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1989 DMSetFunction(), DMSetVariableBounds() 1990 1991 @*/ 1992 PetscErrorCode DMComputeVariableBounds(DM dm, Vec xl, Vec xu) 1993 { 1994 PetscErrorCode ierr; 1995 PetscFunctionBegin; 1996 PetscValidHeaderSpecific(xl,VEC_CLASSID,2); 1997 PetscValidHeaderSpecific(xu,VEC_CLASSID,2); 1998 if(dm->ops->computevariablebounds) { 1999 ierr = (*dm->ops->computevariablebounds)(dm, xl,xu); CHKERRQ(ierr); 2000 } 2001 else { 2002 ierr = VecSet(xl,SNES_VI_NINF); CHKERRQ(ierr); 2003 ierr = VecSet(xu,SNES_VI_INF); CHKERRQ(ierr); 2004 } 2005 PetscFunctionReturn(0); 2006 } 2007 2008 #undef __FUNCT__ 2009 #define __FUNCT__ "DMComputeInitialGuess" 2010 /*@ 2011 DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers 2012 2013 Collective on DM 2014 2015 Input Parameter: 2016 + dm - the DM object to destroy 2017 - x - the vector to hold the initial guess values 2018 2019 Level: developer 2020 2021 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat() 2022 2023 @*/ 2024 PetscErrorCode DMComputeInitialGuess(DM dm,Vec x) 2025 { 2026 PetscErrorCode ierr; 2027 PetscFunctionBegin; 2028 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2029 if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()"); 2030 ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr); 2031 PetscFunctionReturn(0); 2032 } 2033 2034 #undef __FUNCT__ 2035 #define __FUNCT__ "DMHasInitialGuess" 2036 /*@ 2037 DMHasInitialGuess - does the DM object have an initial guess function 2038 2039 Not Collective 2040 2041 Input Parameter: 2042 . dm - the DM object to destroy 2043 2044 Output Parameter: 2045 . flg - PETSC_TRUE if function exists 2046 2047 Level: developer 2048 2049 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 2050 2051 @*/ 2052 PetscErrorCode DMHasInitialGuess(DM dm,PetscBool *flg) 2053 { 2054 PetscFunctionBegin; 2055 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2056 *flg = (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE; 2057 PetscFunctionReturn(0); 2058 } 2059 2060 #undef __FUNCT__ 2061 #define __FUNCT__ "DMHasFunction" 2062 /*@ 2063 DMHasFunction - does the DM object have a function 2064 2065 Not Collective 2066 2067 Input Parameter: 2068 . dm - the DM object to destroy 2069 2070 Output Parameter: 2071 . flg - PETSC_TRUE if function exists 2072 2073 Level: developer 2074 2075 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 2076 2077 @*/ 2078 PetscErrorCode DMHasFunction(DM dm,PetscBool *flg) 2079 { 2080 PetscFunctionBegin; 2081 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2082 *flg = (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE; 2083 PetscFunctionReturn(0); 2084 } 2085 2086 #undef __FUNCT__ 2087 #define __FUNCT__ "DMHasJacobian" 2088 /*@ 2089 DMHasJacobian - does the DM object have a matrix function 2090 2091 Not Collective 2092 2093 Input Parameter: 2094 . dm - the DM object to destroy 2095 2096 Output Parameter: 2097 . flg - PETSC_TRUE if function exists 2098 2099 Level: developer 2100 2101 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 2102 2103 @*/ 2104 PetscErrorCode DMHasJacobian(DM dm,PetscBool *flg) 2105 { 2106 PetscFunctionBegin; 2107 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2108 *flg = (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE; 2109 PetscFunctionReturn(0); 2110 } 2111 2112 #undef __FUNCT__ 2113 #define __FUNCT__ "DMHasColoring" 2114 /*@ 2115 DMHasColoring - does the DM object have a method of providing a coloring? 2116 2117 Not Collective 2118 2119 Input Parameter: 2120 . dm - the DM object 2121 2122 Output Parameter: 2123 . flg - PETSC_TRUE if the DM has facilities for DMCreateColoring(). 2124 2125 Level: developer 2126 2127 .seealso DMHasFunction(), DMCreateColoring() 2128 2129 @*/ 2130 PetscErrorCode DMHasColoring(DM dm,PetscBool *flg) 2131 { 2132 PetscFunctionBegin; 2133 *flg = (dm->ops->getcoloring) ? PETSC_TRUE : PETSC_FALSE; 2134 PetscFunctionReturn(0); 2135 } 2136 2137 2138 #undef __FUNCT__ 2139 #define __FUNCT__ "DMSetVec" 2140 /*@C 2141 DMSetVec - set the vector at which to compute residual, Jacobian and VI bounds, if the problem is nonlinear. 2142 2143 Collective on DM 2144 2145 Input Parameter: 2146 + dm - the DM object 2147 - x - location to compute residual and Jacobian, if PETSC_NULL is passed to those routines; will be PETSC_NULL for linear problems. 2148 2149 Level: developer 2150 2151 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 2152 DMSetFunction(), DMSetJacobian(), DMSetVariableBounds() 2153 2154 @*/ 2155 PetscErrorCode DMSetVec(DM dm,Vec x) 2156 { 2157 PetscErrorCode ierr; 2158 PetscFunctionBegin; 2159 if (x) { 2160 if (!dm->x) { 2161 ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 2162 } 2163 ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 2164 } 2165 else if(dm->x) { 2166 ierr = VecDestroy(&dm->x); CHKERRQ(ierr); 2167 } 2168 PetscFunctionReturn(0); 2169 } 2170 2171 2172 #undef __FUNCT__ 2173 #define __FUNCT__ "DMComputeFunction" 2174 /*@ 2175 DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES 2176 2177 Collective on DM 2178 2179 Input Parameter: 2180 + dm - the DM object to destroy 2181 . x - the location where the function is evaluationed, may be ignored for linear problems 2182 - b - the vector to hold the right hand side entries 2183 2184 Level: developer 2185 2186 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 2187 DMSetJacobian() 2188 2189 @*/ 2190 PetscErrorCode DMComputeFunction(DM dm,Vec x,Vec b) 2191 { 2192 PetscErrorCode ierr; 2193 PetscFunctionBegin; 2194 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2195 if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()"); 2196 PetscStackPush("DM user function"); 2197 ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr); 2198 PetscStackPop; 2199 PetscFunctionReturn(0); 2200 } 2201 2202 2203 2204 #undef __FUNCT__ 2205 #define __FUNCT__ "DMComputeJacobian" 2206 /*@ 2207 DMComputeJacobian - compute the matrix entries for the solver 2208 2209 Collective on DM 2210 2211 Input Parameter: 2212 + dm - the DM object 2213 . x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM 2214 . A - matrix that defines the operator for the linear solve 2215 - B - the matrix used to construct the preconditioner 2216 2217 Level: developer 2218 2219 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 2220 DMSetFunction() 2221 2222 @*/ 2223 PetscErrorCode DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 2224 { 2225 PetscErrorCode ierr; 2226 2227 PetscFunctionBegin; 2228 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2229 if (!dm->ops->jacobian) { 2230 ISColoring coloring; 2231 MatFDColoring fd; 2232 const MatType mtype; 2233 2234 ierr = PetscObjectGetType((PetscObject)B,&mtype);CHKERRQ(ierr); 2235 ierr = DMCreateColoring(dm,dm->coloringtype,mtype,&coloring);CHKERRQ(ierr); 2236 ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr); 2237 ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 2238 ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr); 2239 ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr); 2240 ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr); 2241 2242 dm->fd = fd; 2243 dm->ops->jacobian = DMComputeJacobianDefault; 2244 2245 /* don't know why this is needed */ 2246 ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 2247 } 2248 if (!x) x = dm->x; 2249 ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr); 2250 2251 /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */ 2252 if (x) { 2253 if (!dm->x) { 2254 ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 2255 } 2256 ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 2257 } 2258 if (A != B) { 2259 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2260 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2261 } 2262 PetscFunctionReturn(0); 2263 } 2264 2265 2266 PetscFList DMList = PETSC_NULL; 2267 PetscBool DMRegisterAllCalled = PETSC_FALSE; 2268 2269 #undef __FUNCT__ 2270 #define __FUNCT__ "DMSetType" 2271 /*@C 2272 DMSetType - Builds a DM, for a particular DM implementation. 2273 2274 Collective on DM 2275 2276 Input Parameters: 2277 + dm - The DM object 2278 - method - The name of the DM type 2279 2280 Options Database Key: 2281 . -dm_type <type> - Sets the DM type; use -help for a list of available types 2282 2283 Notes: 2284 See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D). 2285 2286 Level: intermediate 2287 2288 .keywords: DM, set, type 2289 .seealso: DMGetType(), DMCreate() 2290 @*/ 2291 PetscErrorCode DMSetType(DM dm, const DMType method) 2292 { 2293 PetscErrorCode (*r)(DM); 2294 PetscBool match; 2295 PetscErrorCode ierr; 2296 2297 PetscFunctionBegin; 2298 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 2299 ierr = PetscObjectTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr); 2300 if (match) PetscFunctionReturn(0); 2301 2302 if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 2303 ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 2304 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method); 2305 2306 if (dm->ops->destroy) { 2307 ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr); 2308 dm->ops->destroy = PETSC_NULL; 2309 } 2310 ierr = (*r)(dm);CHKERRQ(ierr); 2311 ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr); 2312 PetscFunctionReturn(0); 2313 } 2314 2315 #undef __FUNCT__ 2316 #define __FUNCT__ "DMGetType" 2317 /*@C 2318 DMGetType - Gets the DM type name (as a string) from the DM. 2319 2320 Not Collective 2321 2322 Input Parameter: 2323 . dm - The DM 2324 2325 Output Parameter: 2326 . type - The DM type name 2327 2328 Level: intermediate 2329 2330 .keywords: DM, get, type, name 2331 .seealso: DMSetType(), DMCreate() 2332 @*/ 2333 PetscErrorCode DMGetType(DM dm, const DMType *type) 2334 { 2335 PetscErrorCode ierr; 2336 2337 PetscFunctionBegin; 2338 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 2339 PetscValidCharPointer(type,2); 2340 if (!DMRegisterAllCalled) { 2341 ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr); 2342 } 2343 *type = ((PetscObject)dm)->type_name; 2344 PetscFunctionReturn(0); 2345 } 2346 2347 #undef __FUNCT__ 2348 #define __FUNCT__ "DMConvert" 2349 /*@C 2350 DMConvert - Converts a DM to another DM, either of the same or different type. 2351 2352 Collective on DM 2353 2354 Input Parameters: 2355 + dm - the DM 2356 - newtype - new DM type (use "same" for the same type) 2357 2358 Output Parameter: 2359 . M - pointer to new DM 2360 2361 Notes: 2362 Cannot be used to convert a sequential DM to parallel or parallel to sequential, 2363 the MPI communicator of the generated DM is always the same as the communicator 2364 of the input DM. 2365 2366 Level: intermediate 2367 2368 .seealso: DMCreate() 2369 @*/ 2370 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M) 2371 { 2372 DM B; 2373 char convname[256]; 2374 PetscBool sametype, issame; 2375 PetscErrorCode ierr; 2376 2377 PetscFunctionBegin; 2378 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2379 PetscValidType(dm,1); 2380 PetscValidPointer(M,3); 2381 ierr = PetscObjectTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr); 2382 ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr); 2383 { 2384 PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL; 2385 2386 /* 2387 Order of precedence: 2388 1) See if a specialized converter is known to the current DM. 2389 2) See if a specialized converter is known to the desired DM class. 2390 3) See if a good general converter is registered for the desired class 2391 4) See if a good general converter is known for the current matrix. 2392 5) Use a really basic converter. 2393 */ 2394 2395 /* 1) See if a specialized converter is known to the current DM and the desired class */ 2396 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 2397 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 2398 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 2399 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 2400 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 2401 ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr); 2402 if (conv) goto foundconv; 2403 2404 /* 2) See if a specialized converter is known to the desired DM class. */ 2405 ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr); 2406 ierr = DMSetType(B, newtype);CHKERRQ(ierr); 2407 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 2408 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 2409 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 2410 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 2411 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 2412 ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr); 2413 if (conv) { 2414 ierr = DMDestroy(&B);CHKERRQ(ierr); 2415 goto foundconv; 2416 } 2417 2418 #if 0 2419 /* 3) See if a good general converter is registered for the desired class */ 2420 conv = B->ops->convertfrom; 2421 ierr = DMDestroy(&B);CHKERRQ(ierr); 2422 if (conv) goto foundconv; 2423 2424 /* 4) See if a good general converter is known for the current matrix */ 2425 if (dm->ops->convert) { 2426 conv = dm->ops->convert; 2427 } 2428 if (conv) goto foundconv; 2429 #endif 2430 2431 /* 5) Use a really basic converter. */ 2432 SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype); 2433 2434 foundconv: 2435 ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 2436 ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr); 2437 ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 2438 } 2439 ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr); 2440 PetscFunctionReturn(0); 2441 } 2442 2443 /*--------------------------------------------------------------------------------------------------------------------*/ 2444 2445 #undef __FUNCT__ 2446 #define __FUNCT__ "DMRegister" 2447 /*@C 2448 DMRegister - See DMRegisterDynamic() 2449 2450 Level: advanced 2451 @*/ 2452 PetscErrorCode DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM)) 2453 { 2454 char fullname[PETSC_MAX_PATH_LEN]; 2455 PetscErrorCode ierr; 2456 2457 PetscFunctionBegin; 2458 ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr); 2459 ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr); 2460 ierr = PetscStrcat(fullname, name);CHKERRQ(ierr); 2461 ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr); 2462 PetscFunctionReturn(0); 2463 } 2464 2465 2466 /*--------------------------------------------------------------------------------------------------------------------*/ 2467 #undef __FUNCT__ 2468 #define __FUNCT__ "DMRegisterDestroy" 2469 /*@C 2470 DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic(). 2471 2472 Not Collective 2473 2474 Level: advanced 2475 2476 .keywords: DM, register, destroy 2477 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic() 2478 @*/ 2479 PetscErrorCode DMRegisterDestroy(void) 2480 { 2481 PetscErrorCode ierr; 2482 2483 PetscFunctionBegin; 2484 ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr); 2485 DMRegisterAllCalled = PETSC_FALSE; 2486 PetscFunctionReturn(0); 2487 } 2488 2489 #if defined(PETSC_HAVE_MATLAB_ENGINE) 2490 #include <mex.h> 2491 2492 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext; 2493 2494 #undef __FUNCT__ 2495 #define __FUNCT__ "DMComputeFunction_Matlab" 2496 /* 2497 DMComputeFunction_Matlab - Calls the function that has been set with 2498 DMSetFunctionMatlab(). 2499 2500 For linear problems x is null 2501 2502 .seealso: DMSetFunction(), DMGetFunction() 2503 */ 2504 PetscErrorCode DMComputeFunction_Matlab(DM dm,Vec x,Vec y) 2505 { 2506 PetscErrorCode ierr; 2507 DMMatlabContext *sctx; 2508 int nlhs = 1,nrhs = 4; 2509 mxArray *plhs[1],*prhs[4]; 2510 long long int lx = 0,ly = 0,ls = 0; 2511 2512 PetscFunctionBegin; 2513 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2514 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 2515 PetscCheckSameComm(dm,1,y,3); 2516 2517 /* call Matlab function in ctx with arguments x and y */ 2518 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2519 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 2520 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2521 ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr); 2522 prhs[0] = mxCreateDoubleScalar((double)ls); 2523 prhs[1] = mxCreateDoubleScalar((double)lx); 2524 prhs[2] = mxCreateDoubleScalar((double)ly); 2525 prhs[3] = mxCreateString(sctx->funcname); 2526 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr); 2527 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2528 mxDestroyArray(prhs[0]); 2529 mxDestroyArray(prhs[1]); 2530 mxDestroyArray(prhs[2]); 2531 mxDestroyArray(prhs[3]); 2532 mxDestroyArray(plhs[0]); 2533 PetscFunctionReturn(0); 2534 } 2535 2536 2537 #undef __FUNCT__ 2538 #define __FUNCT__ "DMSetFunctionMatlab" 2539 /* 2540 DMSetFunctionMatlab - Sets the function evaluation routine 2541 2542 */ 2543 PetscErrorCode DMSetFunctionMatlab(DM dm,const char *func) 2544 { 2545 PetscErrorCode ierr; 2546 DMMatlabContext *sctx; 2547 2548 PetscFunctionBegin; 2549 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2550 /* currently sctx is memory bleed */ 2551 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2552 if (!sctx) { 2553 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 2554 } 2555 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2556 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 2557 ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr); 2558 PetscFunctionReturn(0); 2559 } 2560 2561 #undef __FUNCT__ 2562 #define __FUNCT__ "DMComputeJacobian_Matlab" 2563 /* 2564 DMComputeJacobian_Matlab - Calls the function that has been set with 2565 DMSetJacobianMatlab(). 2566 2567 For linear problems x is null 2568 2569 .seealso: DMSetFunction(), DMGetFunction() 2570 */ 2571 PetscErrorCode DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str) 2572 { 2573 PetscErrorCode ierr; 2574 DMMatlabContext *sctx; 2575 int nlhs = 2,nrhs = 5; 2576 mxArray *plhs[2],*prhs[5]; 2577 long long int lx = 0,lA = 0,lB = 0,ls = 0; 2578 2579 PetscFunctionBegin; 2580 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2581 PetscValidHeaderSpecific(A,MAT_CLASSID,3); 2582 2583 /* call MATLAB function in ctx with arguments x, A, and B */ 2584 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2585 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 2586 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2587 ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr); 2588 ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr); 2589 prhs[0] = mxCreateDoubleScalar((double)ls); 2590 prhs[1] = mxCreateDoubleScalar((double)lx); 2591 prhs[2] = mxCreateDoubleScalar((double)lA); 2592 prhs[3] = mxCreateDoubleScalar((double)lB); 2593 prhs[4] = mxCreateString(sctx->jacname); 2594 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr); 2595 *str = (MatStructure) mxGetScalar(plhs[0]); 2596 ierr = (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr); 2597 mxDestroyArray(prhs[0]); 2598 mxDestroyArray(prhs[1]); 2599 mxDestroyArray(prhs[2]); 2600 mxDestroyArray(prhs[3]); 2601 mxDestroyArray(prhs[4]); 2602 mxDestroyArray(plhs[0]); 2603 mxDestroyArray(plhs[1]); 2604 PetscFunctionReturn(0); 2605 } 2606 2607 2608 #undef __FUNCT__ 2609 #define __FUNCT__ "DMSetJacobianMatlab" 2610 /* 2611 DMSetJacobianMatlab - Sets the Jacobian function evaluation routine 2612 2613 */ 2614 PetscErrorCode DMSetJacobianMatlab(DM dm,const char *func) 2615 { 2616 PetscErrorCode ierr; 2617 DMMatlabContext *sctx; 2618 2619 PetscFunctionBegin; 2620 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2621 /* currently sctx is memory bleed */ 2622 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 2623 if (!sctx) { 2624 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 2625 } 2626 ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr); 2627 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 2628 ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr); 2629 PetscFunctionReturn(0); 2630 } 2631 #endif 2632 2633 #undef __FUNCT__ 2634 #define __FUNCT__ "DMLoad" 2635 /*@C 2636 DMLoad - Loads a DM that has been stored in binary or HDF5 format 2637 with DMView(). 2638 2639 Collective on PetscViewer 2640 2641 Input Parameters: 2642 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or 2643 some related function before a call to DMLoad(). 2644 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 2645 HDF5 file viewer, obtained from PetscViewerHDF5Open() 2646 2647 Level: intermediate 2648 2649 Notes: 2650 Defaults to the DM DA. 2651 2652 Notes for advanced users: 2653 Most users should not need to know the details of the binary storage 2654 format, since DMLoad() and DMView() completely hide these details. 2655 But for anyone who's interested, the standard binary matrix storage 2656 format is 2657 .vb 2658 has not yet been determined 2659 .ve 2660 2661 In addition, PETSc automatically does the byte swapping for 2662 machines that store the bytes reversed, e.g. DEC alpha, freebsd, 2663 linux, Windows and the paragon; thus if you write your own binary 2664 read/write routines you have to swap the bytes; see PetscBinaryRead() 2665 and PetscBinaryWrite() to see how this may be done. 2666 2667 Concepts: vector^loading from file 2668 2669 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad() 2670 @*/ 2671 PetscErrorCode DMLoad(DM newdm, PetscViewer viewer) 2672 { 2673 PetscErrorCode ierr; 2674 2675 PetscFunctionBegin; 2676 PetscValidHeaderSpecific(newdm,DM_CLASSID,1); 2677 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2678 2679 if (!((PetscObject)newdm)->type_name) { 2680 ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr); 2681 } 2682 ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr); 2683 PetscFunctionReturn(0); 2684 } 2685 2686 /******************************** FEM Support **********************************/ 2687 2688 #undef __FUNCT__ 2689 #define __FUNCT__ "DMPrintCellVector" 2690 PetscErrorCode DMPrintCellVector(PetscInt c, const char name[], PetscInt len, const PetscScalar x[]) { 2691 PetscInt f; 2692 PetscErrorCode ierr; 2693 2694 PetscFunctionBegin; 2695 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 2696 for(f = 0; f < len; ++f) { 2697 ierr = PetscPrintf(PETSC_COMM_SELF, " | %G |\n", PetscRealPart(x[f]));CHKERRQ(ierr); 2698 } 2699 PetscFunctionReturn(0); 2700 } 2701 2702 #undef __FUNCT__ 2703 #define __FUNCT__ "DMPrintCellMatrix" 2704 PetscErrorCode DMPrintCellMatrix(PetscInt c, const char name[], PetscInt rows, PetscInt cols, const PetscScalar A[]) { 2705 PetscInt f, g; 2706 PetscErrorCode ierr; 2707 2708 PetscFunctionBegin; 2709 ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D Element %s\n", c, name);CHKERRQ(ierr); 2710 for(f = 0; f < rows; ++f) { 2711 ierr = PetscPrintf(PETSC_COMM_SELF, " |");CHKERRQ(ierr); 2712 for(g = 0; g < cols; ++g) { 2713 ierr = PetscPrintf(PETSC_COMM_SELF, " % 9.5G", PetscRealPart(A[f*cols+g]));CHKERRQ(ierr); 2714 } 2715 ierr = PetscPrintf(PETSC_COMM_SELF, " |\n");CHKERRQ(ierr); 2716 } 2717 PetscFunctionReturn(0); 2718 } 2719 2720 #undef __FUNCT__ 2721 #define __FUNCT__ "DMGetLocalFunction" 2722 /*@C 2723 DMGetLocalFunction - Get the local residual function from this DM 2724 2725 Not collective 2726 2727 Input Parameter: 2728 . dm - The DM 2729 2730 Output Parameter: 2731 . lf - The local residual function 2732 2733 Calling sequence of lf: 2734 $ lf (SNES snes, Vec x, Vec f, void *ctx); 2735 2736 + snes - the SNES context 2737 . x - local vector with the state at which to evaluate residual 2738 . f - local vector to put residual in 2739 - ctx - optional user-defined function context 2740 2741 Level: intermediate 2742 2743 .seealso DMSetLocalFunction(), DMGetLocalJacobian(), DMSetLocalJacobian() 2744 @*/ 2745 PetscErrorCode DMGetLocalFunction(DM dm, PetscErrorCode (**lf)(DM, Vec, Vec, void *)) 2746 { 2747 PetscFunctionBegin; 2748 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2749 if (lf) *lf = dm->lf; 2750 PetscFunctionReturn(0); 2751 } 2752 2753 #undef __FUNCT__ 2754 #define __FUNCT__ "DMSetLocalFunction" 2755 /*@C 2756 DMSetLocalFunction - Set the local residual function from this DM 2757 2758 Not collective 2759 2760 Input Parameters: 2761 + dm - The DM 2762 - lf - The local residual function 2763 2764 Calling sequence of lf: 2765 $ lf (SNES snes, Vec x, Vec f, void *ctx); 2766 2767 + snes - the SNES context 2768 . x - local vector with the state at which to evaluate residual 2769 . f - local vector to put residual in 2770 - ctx - optional user-defined function context 2771 2772 Level: intermediate 2773 2774 .seealso DMGetLocalFunction(), DMGetLocalJacobian(), DMSetLocalJacobian() 2775 @*/ 2776 PetscErrorCode DMSetLocalFunction(DM dm, PetscErrorCode (*lf)(DM, Vec, Vec, void *)) 2777 { 2778 PetscFunctionBegin; 2779 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2780 dm->lf = lf; 2781 PetscFunctionReturn(0); 2782 } 2783 2784 #undef __FUNCT__ 2785 #define __FUNCT__ "DMGetLocalJacobian" 2786 /*@C 2787 DMGetLocalJacobian - Get the local Jacobian function from this DM 2788 2789 Not collective 2790 2791 Input Parameter: 2792 . dm - The DM 2793 2794 Output Parameter: 2795 . lj - The local Jacobian function 2796 2797 Calling sequence of lj: 2798 $ lj (SNES snes, Vec x, Mat J, Mat M, void *ctx); 2799 2800 + snes - the SNES context 2801 . x - local vector with the state at which to evaluate residual 2802 . J - matrix to put Jacobian in 2803 . M - matrix to use for defining Jacobian preconditioner 2804 - ctx - optional user-defined function context 2805 2806 Level: intermediate 2807 2808 .seealso DMSetLocalJacobian(), DMGetLocalFunction(), DMSetLocalFunction() 2809 @*/ 2810 PetscErrorCode DMGetLocalJacobian(DM dm, PetscErrorCode (**lj)(DM, Vec, Mat, Mat, void *)) 2811 { 2812 PetscFunctionBegin; 2813 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2814 if (lj) *lj = dm->lj; 2815 PetscFunctionReturn(0); 2816 } 2817 2818 #undef __FUNCT__ 2819 #define __FUNCT__ "DMSetLocalJacobian" 2820 /*@C 2821 DMSetLocalJacobian - Set the local Jacobian function from this DM 2822 2823 Not collective 2824 2825 Input Parameters: 2826 + dm - The DM 2827 - lj - The local Jacobian function 2828 2829 Calling sequence of lj: 2830 $ lj (SNES snes, Vec x, Mat J, Mat M, void *ctx); 2831 2832 + snes - the SNES context 2833 . x - local vector with the state at which to evaluate residual 2834 . J - matrix to put Jacobian in 2835 . M - matrix to use for defining Jacobian preconditioner 2836 - ctx - optional user-defined function context 2837 2838 Level: intermediate 2839 2840 .seealso DMGetLocalJacobian(), DMGetLocalFunction(), DMSetLocalFunction() 2841 @*/ 2842 PetscErrorCode DMSetLocalJacobian(DM dm, PetscErrorCode (*lj)(DM, Vec, Mat, Mat, void *)) 2843 { 2844 PetscFunctionBegin; 2845 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2846 dm->lj = lj; 2847 PetscFunctionReturn(0); 2848 } 2849 2850 #undef __FUNCT__ 2851 #define __FUNCT__ "DMGetDefaultSection" 2852 /*@ 2853 DMGetDefaultSection - Get the PetscSection encoding the local data layout for the DM. 2854 2855 Input Parameter: 2856 . dm - The DM 2857 2858 Output Parameter: 2859 . section - The PetscSection 2860 2861 Level: intermediate 2862 2863 Note: This gets a borrowed reference, so the user should not destroy this PetscSection. 2864 2865 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection() 2866 @*/ 2867 PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *section) { 2868 PetscFunctionBegin; 2869 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2870 PetscValidPointer(section, 2); 2871 *section = dm->defaultSection; 2872 PetscFunctionReturn(0); 2873 } 2874 2875 #undef __FUNCT__ 2876 #define __FUNCT__ "DMSetDefaultSection" 2877 /*@ 2878 DMSetDefaultSection - Set the PetscSection encoding the local data layout for the DM. 2879 2880 Input Parameters: 2881 + dm - The DM 2882 - section - The PetscSection 2883 2884 Level: intermediate 2885 2886 Note: Any existing Section will be destroyed 2887 2888 .seealso: DMSetDefaultSection(), DMGetDefaultGlobalSection() 2889 @*/ 2890 PetscErrorCode DMSetDefaultSection(DM dm, PetscSection section) { 2891 PetscErrorCode ierr; 2892 2893 PetscFunctionBegin; 2894 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2895 ierr = PetscSectionDestroy(&dm->defaultSection);CHKERRQ(ierr); 2896 ierr = PetscSectionDestroy(&dm->defaultGlobalSection);CHKERRQ(ierr); 2897 dm->defaultSection = section; 2898 PetscFunctionReturn(0); 2899 } 2900 2901 #undef __FUNCT__ 2902 #define __FUNCT__ "DMGetDefaultGlobalSection" 2903 /*@ 2904 DMGetDefaultGlobalSection - Get the PetscSection encoding the global data layout for the DM. 2905 2906 Input Parameter: 2907 . dm - The DM 2908 2909 Output Parameter: 2910 . section - The PetscSection 2911 2912 Level: intermediate 2913 2914 Note: This gets a borrowed reference, so the user should not destroy this PetscSection. 2915 2916 .seealso: DMSetDefaultSection(), DMGetDefaultSection() 2917 @*/ 2918 PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *section) { 2919 PetscErrorCode ierr; 2920 2921 PetscFunctionBegin; 2922 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2923 PetscValidPointer(section, 2); 2924 if (!dm->defaultGlobalSection) { 2925 ierr = PetscSectionCreateGlobalSection(dm->defaultSection, dm->sf, &dm->defaultGlobalSection);CHKERRQ(ierr); 2926 } 2927 *section = dm->defaultGlobalSection; 2928 PetscFunctionReturn(0); 2929 } 2930 2931 #undef __FUNCT__ 2932 #define __FUNCT__ "DMGetDefaultSF" 2933 /*@ 2934 DMGetDefaultSF - Get the PetscSF encoding the parallel dof overlap for the DM. If it has not been set, 2935 it is created from the default PetscSection layouts in the DM. 2936 2937 Input Parameter: 2938 . dm - The DM 2939 2940 Output Parameter: 2941 . sf - The PetscSF 2942 2943 Level: intermediate 2944 2945 Note: This gets a borrowed reference, so the user should not destroy this PetscSF. 2946 2947 .seealso: DMSetDefaultSF(), DMCreateDefaultSF() 2948 @*/ 2949 PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *sf) { 2950 PetscInt nroots; 2951 PetscErrorCode ierr; 2952 2953 PetscFunctionBegin; 2954 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2955 PetscValidPointer(sf, 2); 2956 ierr = PetscSFGetGraph(dm->defaultSF, &nroots, PETSC_NULL, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 2957 if (nroots < 0) { 2958 PetscSection section, gSection; 2959 2960 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 2961 if (section) { 2962 ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr); 2963 ierr = DMCreateDefaultSF(dm, section, gSection);CHKERRQ(ierr); 2964 } else { 2965 *sf = PETSC_NULL; 2966 PetscFunctionReturn(0); 2967 } 2968 } 2969 *sf = dm->defaultSF; 2970 PetscFunctionReturn(0); 2971 } 2972 2973 #undef __FUNCT__ 2974 #define __FUNCT__ "DMSetDefaultSF" 2975 /*@ 2976 DMSetDefaultSF - Set the PetscSF encoding the parallel dof overlap for the DM 2977 2978 Input Parameters: 2979 + dm - The DM 2980 - sf - The PetscSF 2981 2982 Level: intermediate 2983 2984 Note: Any previous SF is destroyed 2985 2986 .seealso: DMGetDefaultSF(), DMCreateDefaultSF() 2987 @*/ 2988 PetscErrorCode DMSetDefaultSF(DM dm, PetscSF sf) { 2989 PetscErrorCode ierr; 2990 2991 PetscFunctionBegin; 2992 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2993 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 2994 ierr = PetscSFDestroy(&dm->defaultSF);CHKERRQ(ierr); 2995 dm->defaultSF = sf; 2996 PetscFunctionReturn(0); 2997 } 2998 2999 #undef __FUNCT__ 3000 #define __FUNCT__ "DMCreateDefaultSF" 3001 /*@C 3002 DMCreateDefaultSF - Create the PetscSF encoding the parallel dof overlap for the DM based upon the PetscSections 3003 describing the data layout. 3004 3005 Input Parameters: 3006 + dm - The DM 3007 . localSection - PetscSection describing the local data layout 3008 - globalSection - PetscSection describing the global data layout 3009 3010 Level: intermediate 3011 3012 .seealso: DMGetDefaultSF(), DMSetDefaultSF() 3013 @*/ 3014 PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection localSection, PetscSection globalSection) 3015 { 3016 MPI_Comm comm = ((PetscObject) dm)->comm; 3017 PetscLayout layout; 3018 const PetscInt *ranges; 3019 PetscInt *local; 3020 PetscSFNode *remote; 3021 PetscInt pStart, pEnd, p, nroots, nleaves, l; 3022 PetscMPIInt size, rank; 3023 PetscErrorCode ierr; 3024 3025 PetscFunctionBegin; 3026 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3027 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 3028 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 3029 ierr = PetscSectionGetChart(globalSection, &pStart, &pEnd);CHKERRQ(ierr); 3030 ierr = PetscSectionGetConstrainedStorageSize(globalSection, &nroots);CHKERRQ(ierr); 3031 ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 3032 ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr); 3033 ierr = PetscLayoutSetLocalSize(layout, nroots);CHKERRQ(ierr); 3034 ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 3035 ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr); 3036 for(p = pStart, nleaves = 0; p < pEnd; ++p) { 3037 PetscInt dof, cdof; 3038 3039 ierr = PetscSectionGetDof(globalSection, p, &dof);CHKERRQ(ierr); 3040 ierr = PetscSectionGetConstraintDof(globalSection, p, &cdof);CHKERRQ(ierr); 3041 nleaves += dof < 0 ? -(dof+1)-cdof : dof-cdof; 3042 } 3043 ierr = PetscMalloc(nleaves * sizeof(PetscInt), &local);CHKERRQ(ierr); 3044 ierr = PetscMalloc(nleaves * sizeof(PetscSFNode), &remote);CHKERRQ(ierr); 3045 for(p = pStart, l = 0; p < pEnd; ++p) { 3046 PetscInt *cind; 3047 PetscInt dof, gdof, cdof, dim, off, goff, d, c; 3048 3049 ierr = PetscSectionGetDof(localSection, p, &dof);CHKERRQ(ierr); 3050 ierr = PetscSectionGetOffset(localSection, p, &off);CHKERRQ(ierr); 3051 ierr = PetscSectionGetConstraintDof(localSection, p, &cdof);CHKERRQ(ierr); 3052 ierr = PetscSectionGetConstraintIndices(localSection, p, &cind);CHKERRQ(ierr); 3053 ierr = PetscSectionGetDof(globalSection, p, &gdof);CHKERRQ(ierr); 3054 ierr = PetscSectionGetOffset(globalSection, p, &goff);CHKERRQ(ierr); 3055 dim = dof-cdof; 3056 for(d = 0, c = 0; d < dof; ++d) { 3057 if ((c < cdof) && (cind[c] == d)) {++c; continue;} 3058 local[l+d-c] = off+d; 3059 } 3060 if (gdof < 0) { 3061 for(d = 0; d < dim; ++d, ++l) { 3062 PetscInt offset = -(goff+1) + d, r; 3063 3064 for(r = 0; r < size; ++r) { 3065 if ((offset >= ranges[r]) && (offset < ranges[r+1])) break; 3066 } 3067 remote[l].rank = r; 3068 remote[l].index = offset - ranges[r]; 3069 } 3070 } else { 3071 for(d = 0; d < dim; ++d, ++l) { 3072 remote[l].rank = rank; 3073 remote[l].index = goff+d - ranges[rank]; 3074 } 3075 } 3076 } 3077 ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 3078 ierr = PetscSFSetGraph(dm->defaultSF, nroots, nleaves, local, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr); 3079 PetscFunctionReturn(0); 3080 } 3081