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