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