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