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