1 2 #include <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 v->workSize = 0; 45 v->workArray = PETSC_NULL; 46 v->ltogmap = PETSC_NULL; 47 v->ltogmapb = PETSC_NULL; 48 v->bs = 1; 49 v->coloringtype = IS_COLORING_GLOBAL; 50 51 *dm = v; 52 PetscFunctionReturn(0); 53 } 54 55 56 #undef __FUNCT__ 57 #define __FUNCT__ "DMSetVecType" 58 /*@C 59 DMSetVecType - Sets the type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() 60 61 Logically Collective on DMDA 62 63 Input Parameter: 64 + da - initial distributed array 65 . ctype - the vector type, currently either VECSTANDARD or VECCUSP 66 67 Options Database: 68 . -dm_vec_type ctype 69 70 Level: intermediate 71 72 .seealso: DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDestroy(), DMDA, DMDAInterpolationType, VecType 73 @*/ 74 PetscErrorCode DMSetVecType(DM da,const VecType ctype) 75 { 76 PetscErrorCode ierr; 77 78 PetscFunctionBegin; 79 PetscValidHeaderSpecific(da,DM_CLASSID,1); 80 ierr = PetscFree(da->vectype);CHKERRQ(ierr); 81 ierr = PetscStrallocpy(ctype,&da->vectype);CHKERRQ(ierr); 82 PetscFunctionReturn(0); 83 } 84 85 #undef __FUNCT__ 86 #define __FUNCT__ "DMSetOptionsPrefix" 87 /*@C 88 DMSetOptionsPrefix - Sets the prefix used for searching for all 89 DMDA options in the database. 90 91 Logically Collective on DMDA 92 93 Input Parameter: 94 + da - the DMDA context 95 - prefix - the prefix to prepend to all option names 96 97 Notes: 98 A hyphen (-) must NOT be given at the beginning of the prefix name. 99 The first character of all runtime options is AUTOMATICALLY the hyphen. 100 101 Level: advanced 102 103 .keywords: DMDA, set, options, prefix, database 104 105 .seealso: DMSetFromOptions() 106 @*/ 107 PetscErrorCode DMSetOptionsPrefix(DM dm,const char prefix[]) 108 { 109 PetscErrorCode ierr; 110 111 PetscFunctionBegin; 112 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 113 ierr = PetscObjectSetOptionsPrefix((PetscObject)dm,prefix);CHKERRQ(ierr); 114 PetscFunctionReturn(0); 115 } 116 117 #undef __FUNCT__ 118 #define __FUNCT__ "DMDestroy" 119 /*@ 120 DMDestroy - Destroys a vector packer or DMDA. 121 122 Collective on DM 123 124 Input Parameter: 125 . dm - the DM object to destroy 126 127 Level: developer 128 129 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 130 131 @*/ 132 PetscErrorCode DMDestroy(DM *dm) 133 { 134 PetscInt i, cnt = 0; 135 PetscErrorCode ierr; 136 137 PetscFunctionBegin; 138 if (!*dm) PetscFunctionReturn(0); 139 PetscValidHeaderSpecific((*dm),DM_CLASSID,1); 140 141 /* count all the circular references of DM and its contained Vecs */ 142 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 143 if ((*dm)->localin[i]) {cnt++;} 144 if ((*dm)->globalin[i]) {cnt++;} 145 } 146 if ((*dm)->x) { 147 PetscObject obj; 148 ierr = PetscObjectQuery((PetscObject)(*dm)->x,"DM",&obj);CHKERRQ(ierr); 149 if (obj == (PetscObject)*dm) cnt++; 150 } 151 152 if (--((PetscObject)(*dm))->refct - cnt > 0) {*dm = 0; PetscFunctionReturn(0);} 153 /* 154 Need this test because the dm references the vectors that 155 reference the dm, so destroying the dm calls destroy on the 156 vectors that cause another destroy on the dm 157 */ 158 if (((PetscObject)(*dm))->refct < 0) PetscFunctionReturn(0); 159 ((PetscObject) (*dm))->refct = 0; 160 for (i=0; i<DM_MAX_WORK_VECTORS; i++) { 161 if ((*dm)->localout[i]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Destroying a DM that has a local vector obtained with DMGetLocalVector()"); 162 ierr = VecDestroy(&(*dm)->localin[i]);CHKERRQ(ierr); 163 } 164 165 if ((*dm)->ctx && (*dm)->ctxdestroy) { 166 ierr = (*(*dm)->ctxdestroy)(&(*dm)->ctx);CHKERRQ(ierr); 167 } 168 ierr = VecDestroy(&(*dm)->x);CHKERRQ(ierr); 169 ierr = MatFDColoringDestroy(&(*dm)->fd);CHKERRQ(ierr); 170 ierr = DMClearGlobalVectors(*dm);CHKERRQ(ierr); 171 ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmap);CHKERRQ(ierr); 172 ierr = ISLocalToGlobalMappingDestroy(&(*dm)->ltogmapb);CHKERRQ(ierr); 173 ierr = PetscFree((*dm)->vectype);CHKERRQ(ierr); 174 ierr = PetscFree((*dm)->mattype);CHKERRQ(ierr); 175 ierr = PetscFree((*dm)->workArray);CHKERRQ(ierr); 176 /* if memory was published with AMS then destroy it */ 177 ierr = PetscObjectDepublish(*dm);CHKERRQ(ierr); 178 179 ierr = (*(*dm)->ops->destroy)(*dm);CHKERRQ(ierr); 180 ierr = PetscFree((*dm)->data);CHKERRQ(ierr); 181 ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr); 182 PetscFunctionReturn(0); 183 } 184 185 #undef __FUNCT__ 186 #define __FUNCT__ "DMSetUp" 187 /*@ 188 DMSetUp - sets up the data structures inside a DM object 189 190 Collective on DM 191 192 Input Parameter: 193 . dm - the DM object to setup 194 195 Level: developer 196 197 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 198 199 @*/ 200 PetscErrorCode DMSetUp(DM dm) 201 { 202 PetscErrorCode ierr; 203 204 PetscFunctionBegin; 205 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 206 if (dm->setupcalled) PetscFunctionReturn(0); 207 if (dm->ops->setup) { 208 ierr = (*dm->ops->setup)(dm);CHKERRQ(ierr); 209 } 210 dm->setupcalled = PETSC_TRUE; 211 PetscFunctionReturn(0); 212 } 213 214 #undef __FUNCT__ 215 #define __FUNCT__ "DMSetFromOptions" 216 /*@ 217 DMSetFromOptions - sets parameters in a DM from the options database 218 219 Collective on DM 220 221 Input Parameter: 222 . dm - the DM object to set options for 223 224 Options Database: 225 + -dm_preallocate_only: Only preallocate the matrix for DMCreateMatrix(), but do not fill it with zeros 226 . -dm_vec_type <type> type of vector to create inside DM 227 . -dm_mat_type <type> type of matrix to create inside DM 228 - -dm_coloring_type <global or ghosted> 229 230 Level: developer 231 232 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 233 234 @*/ 235 PetscErrorCode DMSetFromOptions(DM dm) 236 { 237 PetscBool flg1 = PETSC_FALSE,flg; 238 PetscErrorCode ierr; 239 char typeName[256] = MATAIJ; 240 241 PetscFunctionBegin; 242 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 243 ierr = PetscObjectOptionsBegin((PetscObject)dm);CHKERRQ(ierr); 244 ierr = PetscOptionsBool("-dm_view", "Information on DM", "DMView", flg1, &flg1, PETSC_NULL);CHKERRQ(ierr); 245 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); 246 ierr = PetscOptionsList("-dm_vec_type","Vector type used for created vectors","DMSetVecType",VecList,dm->vectype,typeName,256,&flg);CHKERRQ(ierr); 247 if (flg) { 248 ierr = DMSetVecType(dm,typeName);CHKERRQ(ierr); 249 } 250 ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,typeName,typeName,sizeof typeName,&flg);CHKERRQ(ierr); 251 if (flg) { 252 ierr = PetscFree(dm->mattype);CHKERRQ(ierr); 253 ierr = PetscStrallocpy(typeName,&dm->mattype);CHKERRQ(ierr); 254 } 255 ierr = PetscOptionsEnum("-dm_is_coloring_type","Global or local coloring of Jacobian","ISColoringType",ISColoringTypes,(PetscEnum)dm->coloringtype,(PetscEnum*)&dm->coloringtype,PETSC_NULL);CHKERRQ(ierr); 256 if (dm->ops->setfromoptions) { 257 ierr = (*dm->ops->setfromoptions)(dm);CHKERRQ(ierr); 258 } 259 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 260 ierr = PetscObjectProcessOptionsHandlers((PetscObject) dm);CHKERRQ(ierr); 261 ierr = PetscOptionsEnd();CHKERRQ(ierr); 262 if (flg1) { 263 ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 264 } 265 PetscFunctionReturn(0); 266 } 267 268 #undef __FUNCT__ 269 #define __FUNCT__ "DMView" 270 /*@C 271 DMView - Views a vector packer or DMDA. 272 273 Collective on DM 274 275 Input Parameter: 276 + dm - the DM object to view 277 - v - the viewer 278 279 Level: developer 280 281 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 282 283 @*/ 284 PetscErrorCode DMView(DM dm,PetscViewer v) 285 { 286 PetscErrorCode ierr; 287 288 PetscFunctionBegin; 289 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 290 if (!v) { 291 ierr = PetscViewerASCIIGetStdout(((PetscObject)dm)->comm,&v);CHKERRQ(ierr); 292 } 293 if (dm->ops->view) { 294 ierr = (*dm->ops->view)(dm,v);CHKERRQ(ierr); 295 } 296 PetscFunctionReturn(0); 297 } 298 299 #undef __FUNCT__ 300 #define __FUNCT__ "DMCreateGlobalVector" 301 /*@ 302 DMCreateGlobalVector - Creates a global vector from a DMDA or DMComposite object 303 304 Collective on DM 305 306 Input Parameter: 307 . dm - the DM object 308 309 Output Parameter: 310 . vec - the global vector 311 312 Level: beginner 313 314 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 315 316 @*/ 317 PetscErrorCode DMCreateGlobalVector(DM dm,Vec *vec) 318 { 319 PetscErrorCode ierr; 320 321 PetscFunctionBegin; 322 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 323 ierr = (*dm->ops->createglobalvector)(dm,vec);CHKERRQ(ierr); 324 PetscFunctionReturn(0); 325 } 326 327 #undef __FUNCT__ 328 #define __FUNCT__ "DMCreateLocalVector" 329 /*@ 330 DMCreateLocalVector - Creates a local vector from a DMDA or DMComposite object 331 332 Not Collective 333 334 Input Parameter: 335 . dm - the DM object 336 337 Output Parameter: 338 . vec - the local vector 339 340 Level: beginner 341 342 .seealso DMDestroy(), DMView(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix() 343 344 @*/ 345 PetscErrorCode DMCreateLocalVector(DM dm,Vec *vec) 346 { 347 PetscErrorCode ierr; 348 349 PetscFunctionBegin; 350 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 351 ierr = (*dm->ops->createlocalvector)(dm,vec);CHKERRQ(ierr); 352 PetscFunctionReturn(0); 353 } 354 355 #undef __FUNCT__ 356 #define __FUNCT__ "DMGetLocalToGlobalMapping" 357 /*@ 358 DMGetLocalToGlobalMapping - Accesses the local-to-global mapping in a DM. 359 360 Collective on DM 361 362 Input Parameter: 363 . dm - the DM that provides the mapping 364 365 Output Parameter: 366 . ltog - the mapping 367 368 Level: intermediate 369 370 Notes: 371 This mapping can then be used by VecSetLocalToGlobalMapping() or 372 MatSetLocalToGlobalMapping(). 373 374 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMappingBlock() 375 @*/ 376 PetscErrorCode DMGetLocalToGlobalMapping(DM dm,ISLocalToGlobalMapping *ltog) 377 { 378 PetscErrorCode ierr; 379 380 PetscFunctionBegin; 381 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 382 PetscValidPointer(ltog,2); 383 if (!dm->ltogmap) { 384 if (!dm->ops->createlocaltoglobalmapping) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMapping"); 385 ierr = (*dm->ops->createlocaltoglobalmapping)(dm);CHKERRQ(ierr); 386 } 387 *ltog = dm->ltogmap; 388 PetscFunctionReturn(0); 389 } 390 391 #undef __FUNCT__ 392 #define __FUNCT__ "DMGetLocalToGlobalMappingBlock" 393 /*@ 394 DMGetLocalToGlobalMappingBlock - Accesses the blocked local-to-global mapping in a DM. 395 396 Collective on DM 397 398 Input Parameter: 399 . da - the distributed array that provides the mapping 400 401 Output Parameter: 402 . ltog - the block mapping 403 404 Level: intermediate 405 406 Notes: 407 This mapping can then be used by VecSetLocalToGlobalMappingBlock() or 408 MatSetLocalToGlobalMappingBlock(). 409 410 .seealso: DMCreateLocalVector(), DMGetLocalToGlobalMapping(), DMGetBlockSize(), VecSetBlockSize(), MatSetBlockSize() 411 @*/ 412 PetscErrorCode DMGetLocalToGlobalMappingBlock(DM dm,ISLocalToGlobalMapping *ltog) 413 { 414 PetscErrorCode ierr; 415 416 PetscFunctionBegin; 417 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 418 PetscValidPointer(ltog,2); 419 if (!dm->ltogmapb) { 420 PetscInt bs; 421 ierr = DMGetBlockSize(dm,&bs);CHKERRQ(ierr); 422 if (bs > 1) { 423 if (!dm->ops->createlocaltoglobalmappingblock) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"DM can not create LocalToGlobalMappingBlock"); 424 ierr = (*dm->ops->createlocaltoglobalmappingblock)(dm);CHKERRQ(ierr); 425 } else { 426 ierr = DMGetLocalToGlobalMapping(dm,&dm->ltogmapb);CHKERRQ(ierr); 427 ierr = PetscObjectReference((PetscObject)dm->ltogmapb);CHKERRQ(ierr); 428 } 429 } 430 *ltog = dm->ltogmapb; 431 PetscFunctionReturn(0); 432 } 433 434 #undef __FUNCT__ 435 #define __FUNCT__ "DMGetBlockSize" 436 /*@ 437 DMGetBlockSize - Gets the inherent block size associated with a DM 438 439 Not Collective 440 441 Input Parameter: 442 . dm - the DM with block structure 443 444 Output Parameter: 445 . bs - the block size, 1 implies no exploitable block structure 446 447 Level: intermediate 448 449 .seealso: ISCreateBlock(), VecSetBlockSize(), MatSetBlockSize(), DMGetLocalToGlobalMappingBlock() 450 @*/ 451 PetscErrorCode DMGetBlockSize(DM dm,PetscInt *bs) 452 { 453 PetscFunctionBegin; 454 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 455 PetscValidPointer(bs,2); 456 if (dm->bs < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DM does not have enough information to provide a block size yet"); 457 *bs = dm->bs; 458 PetscFunctionReturn(0); 459 } 460 461 #undef __FUNCT__ 462 #define __FUNCT__ "DMCreateInterpolation" 463 /*@ 464 DMCreateInterpolation - Gets interpolation matrix between two DMDA or DMComposite objects 465 466 Collective on DM 467 468 Input Parameter: 469 + dm1 - the DM object 470 - dm2 - the second, finer DM object 471 472 Output Parameter: 473 + mat - the interpolation 474 - vec - the scaling (optional) 475 476 Level: developer 477 478 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 479 DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the interpolation. 480 481 For DMDA objects you can use this interpolation (more precisely the interpolation from the DMDAGetCoordinateDA()) to interpolate the mesh coordinate vectors 482 EXCEPT in the periodic case where it does not make sense since the coordinate vectors are not periodic. 483 484 485 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMRefine(), DMCoarsen() 486 487 @*/ 488 PetscErrorCode DMCreateInterpolation(DM dm1,DM dm2,Mat *mat,Vec *vec) 489 { 490 PetscErrorCode ierr; 491 492 PetscFunctionBegin; 493 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 494 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 495 ierr = (*dm1->ops->createinterpolation)(dm1,dm2,mat,vec);CHKERRQ(ierr); 496 PetscFunctionReturn(0); 497 } 498 499 #undef __FUNCT__ 500 #define __FUNCT__ "DMCreateInjection" 501 /*@ 502 DMCreateInjection - Gets injection matrix between two DMDA or DMComposite objects 503 504 Collective on DM 505 506 Input Parameter: 507 + dm1 - the DM object 508 - dm2 - the second, finer DM object 509 510 Output Parameter: 511 . ctx - the injection 512 513 Level: developer 514 515 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 516 DMCoarsen(). The coordinates set into the DMDA are completely ignored in computing the injection. 517 518 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateColoring(), DMCreateMatrix(), DMCreateInterpolation() 519 520 @*/ 521 PetscErrorCode DMCreateInjection(DM dm1,DM dm2,VecScatter *ctx) 522 { 523 PetscErrorCode ierr; 524 525 PetscFunctionBegin; 526 PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 527 PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 528 ierr = (*dm1->ops->getinjection)(dm1,dm2,ctx);CHKERRQ(ierr); 529 PetscFunctionReturn(0); 530 } 531 532 #undef __FUNCT__ 533 #define __FUNCT__ "DMCreateColoring" 534 /*@C 535 DMCreateColoring - Gets coloring for a DMDA or DMComposite 536 537 Collective on DM 538 539 Input Parameter: 540 + dm - the DM object 541 . ctype - IS_COLORING_GHOSTED or IS_COLORING_GLOBAL 542 - matype - either MATAIJ or MATBAIJ 543 544 Output Parameter: 545 . coloring - the coloring 546 547 Level: developer 548 549 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateMatrix() 550 551 @*/ 552 PetscErrorCode DMCreateColoring(DM dm,ISColoringType ctype,const MatType mtype,ISColoring *coloring) 553 { 554 PetscErrorCode ierr; 555 556 PetscFunctionBegin; 557 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 558 if (!dm->ops->getcoloring) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No coloring for this type of DM yet"); 559 ierr = (*dm->ops->getcoloring)(dm,ctype,mtype,coloring);CHKERRQ(ierr); 560 PetscFunctionReturn(0); 561 } 562 563 #undef __FUNCT__ 564 #define __FUNCT__ "DMCreateMatrix" 565 /*@C 566 DMCreateMatrix - Gets empty Jacobian for a DMDA or DMComposite 567 568 Collective on DM 569 570 Input Parameter: 571 + dm - the DM object 572 - mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, or 573 any type which inherits from one of these (such as MATAIJ) 574 575 Output Parameter: 576 . mat - the empty Jacobian 577 578 Level: beginner 579 580 Notes: This properly preallocates the number of nonzeros in the sparse matrix so you 581 do not need to do it yourself. 582 583 By default it also sets the nonzero structure and puts in the zero entries. To prevent setting 584 the nonzero pattern call DMDASetMatPreallocateOnly() 585 586 For structured grid problems, when you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used 587 internally by PETSc. 588 589 For structured grid problems, in general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires 590 the indices for the global numbering for DMDAs which is complicated. 591 592 .seealso DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 593 594 @*/ 595 PetscErrorCode DMCreateMatrix(DM dm,const MatType mtype,Mat *mat) 596 { 597 PetscErrorCode ierr; 598 599 PetscFunctionBegin; 600 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 601 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 602 ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr); 603 #endif 604 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 605 PetscValidPointer(mat,3); 606 if (dm->mattype) { 607 ierr = (*dm->ops->creatematrix)(dm,dm->mattype,mat);CHKERRQ(ierr); 608 } else { 609 ierr = (*dm->ops->creatematrix)(dm,mtype,mat);CHKERRQ(ierr); 610 } 611 PetscFunctionReturn(0); 612 } 613 614 #undef __FUNCT__ 615 #define __FUNCT__ "DMSetMatrixPreallocateOnly" 616 /*@ 617 DMSetMatrixPreallocateOnly - When DMCreateMatrix() is called the matrix will be properly 618 preallocated but the nonzero structure and zero values will not be set. 619 620 Logically Collective on DMDA 621 622 Input Parameter: 623 + dm - the DM 624 - only - PETSC_TRUE if only want preallocation 625 626 Level: developer 627 .seealso DMCreateMatrix() 628 @*/ 629 PetscErrorCode DMSetMatrixPreallocateOnly(DM dm, PetscBool only) 630 { 631 PetscFunctionBegin; 632 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 633 dm->prealloc_only = only; 634 PetscFunctionReturn(0); 635 } 636 637 #undef __FUNCT__ 638 #define __FUNCT__ "DMGetWorkArray" 639 /*@C 640 DMGetWorkArray - Gets a work array guaranteed to be at least the input size 641 642 Not Collective 643 644 Input Parameters: 645 + dm - the DM object 646 - size - The minium size 647 648 Output Parameter: 649 . array - the work array 650 651 Level: developer 652 653 .seealso DMDestroy(), DMCreate() 654 @*/ 655 PetscErrorCode DMGetWorkArray(DM dm,PetscInt size,PetscScalar **array) 656 { 657 PetscErrorCode ierr; 658 659 PetscFunctionBegin; 660 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 661 PetscValidPointer(array,3); 662 if (size > dm->workSize) { 663 dm->workSize = size; 664 ierr = PetscFree(dm->workArray);CHKERRQ(ierr); 665 ierr = PetscMalloc(dm->workSize * sizeof(PetscScalar), &dm->workArray);CHKERRQ(ierr); 666 } 667 *array = dm->workArray; 668 PetscFunctionReturn(0); 669 } 670 671 #undef __FUNCT__ 672 #define __FUNCT__ "DMRefine" 673 /*@ 674 DMRefine - Refines a DM object 675 676 Collective on DM 677 678 Input Parameter: 679 + dm - the DM object 680 - comm - the communicator to contain the new DM object (or PETSC_NULL) 681 682 Output Parameter: 683 . dmf - the refined DM 684 685 Level: developer 686 687 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 688 689 @*/ 690 PetscErrorCode DMRefine(DM dm,MPI_Comm comm,DM *dmf) 691 { 692 PetscErrorCode ierr; 693 694 PetscFunctionBegin; 695 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 696 ierr = (*dm->ops->refine)(dm,comm,dmf);CHKERRQ(ierr); 697 (*dmf)->ops->initialguess = dm->ops->initialguess; 698 (*dmf)->ops->function = dm->ops->function; 699 (*dmf)->ops->functionj = dm->ops->functionj; 700 if (dm->ops->jacobian != DMComputeJacobianDefault) { 701 (*dmf)->ops->jacobian = dm->ops->jacobian; 702 } 703 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmf);CHKERRQ(ierr); 704 (*dmf)->ctx = dm->ctx; 705 (*dmf)->levelup = dm->levelup + 1; 706 PetscFunctionReturn(0); 707 } 708 709 #undef __FUNCT__ 710 #define __FUNCT__ "DMGetRefineLevel" 711 /*@ 712 DMGetRefineLevel - Get's the number of refinements that have generated this DM. 713 714 Not Collective 715 716 Input Parameter: 717 . dm - the DM object 718 719 Output Parameter: 720 . level - number of refinements 721 722 Level: developer 723 724 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 725 726 @*/ 727 PetscErrorCode DMGetRefineLevel(DM dm,PetscInt *level) 728 { 729 PetscFunctionBegin; 730 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 731 *level = dm->levelup; 732 PetscFunctionReturn(0); 733 } 734 735 #undef __FUNCT__ 736 #define __FUNCT__ "DMGlobalToLocalBegin" 737 /*@ 738 DMGlobalToLocalBegin - Begins updating local vectors from global vector 739 740 Neighbor-wise Collective on DM 741 742 Input Parameters: 743 + dm - the DM object 744 . g - the global vector 745 . mode - INSERT_VALUES or ADD_VALUES 746 - l - the local vector 747 748 749 Level: beginner 750 751 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 752 753 @*/ 754 PetscErrorCode DMGlobalToLocalBegin(DM dm,Vec g,InsertMode mode,Vec l) 755 { 756 PetscErrorCode ierr; 757 758 PetscFunctionBegin; 759 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 760 ierr = (*dm->ops->globaltolocalbegin)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 761 PetscFunctionReturn(0); 762 } 763 764 #undef __FUNCT__ 765 #define __FUNCT__ "DMGlobalToLocalEnd" 766 /*@ 767 DMGlobalToLocalEnd - Ends updating local vectors from global vector 768 769 Neighbor-wise Collective on DM 770 771 Input Parameters: 772 + dm - the DM object 773 . g - the global vector 774 . mode - INSERT_VALUES or ADD_VALUES 775 - l - the local vector 776 777 778 Level: beginner 779 780 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin() 781 782 @*/ 783 PetscErrorCode DMGlobalToLocalEnd(DM dm,Vec g,InsertMode mode,Vec l) 784 { 785 PetscErrorCode ierr; 786 787 PetscFunctionBegin; 788 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 789 ierr = (*dm->ops->globaltolocalend)(dm,g,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),l);CHKERRQ(ierr); 790 PetscFunctionReturn(0); 791 } 792 793 #undef __FUNCT__ 794 #define __FUNCT__ "DMLocalToGlobalBegin" 795 /*@ 796 DMLocalToGlobalBegin - updates global vectors from local vectors 797 798 Neighbor-wise Collective on DM 799 800 Input Parameters: 801 + dm - the DM object 802 . l - the local vector 803 . 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 804 base point. 805 - - the global vector 806 807 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 808 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 809 global array to the final global array with VecAXPY(). 810 811 Level: beginner 812 813 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalBegin() 814 815 @*/ 816 PetscErrorCode DMLocalToGlobalBegin(DM dm,Vec l,InsertMode mode,Vec g) 817 { 818 PetscErrorCode ierr; 819 820 PetscFunctionBegin; 821 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 822 ierr = (*dm->ops->localtoglobalbegin)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 823 PetscFunctionReturn(0); 824 } 825 826 #undef __FUNCT__ 827 #define __FUNCT__ "DMLocalToGlobalEnd" 828 /*@ 829 DMLocalToGlobalEnd - updates global vectors from local vectors 830 831 Neighbor-wise Collective on DM 832 833 Input Parameters: 834 + dm - the DM object 835 . l - the local vector 836 . mode - INSERT_VALUES or ADD_VALUES 837 - g - the global vector 838 839 840 Level: beginner 841 842 .seealso DMCoarsen(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMGlobalToLocalEnd(), DMGlobalToLocalEnd() 843 844 @*/ 845 PetscErrorCode DMLocalToGlobalEnd(DM dm,Vec l,InsertMode mode,Vec g) 846 { 847 PetscErrorCode ierr; 848 849 PetscFunctionBegin; 850 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 851 ierr = (*dm->ops->localtoglobalend)(dm,l,mode == INSERT_ALL_VALUES ? INSERT_VALUES : (mode == ADD_ALL_VALUES ? ADD_VALUES : mode),g);CHKERRQ(ierr); 852 PetscFunctionReturn(0); 853 } 854 855 #undef __FUNCT__ 856 #define __FUNCT__ "DMComputeJacobianDefault" 857 /*@ 858 DMComputeJacobianDefault - computes the Jacobian using the DMComputeFunction() if Jacobian computer is not provided 859 860 Collective on DM 861 862 Input Parameter: 863 + dm - the DM object 864 . x - location to compute Jacobian at; may be ignored for linear problems 865 . A - matrix that defines the operator for the linear solve 866 - B - the matrix used to construct the preconditioner 867 868 Level: developer 869 870 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 871 DMSetFunction() 872 873 @*/ 874 PetscErrorCode DMComputeJacobianDefault(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 875 { 876 PetscErrorCode ierr; 877 878 PetscFunctionBegin; 879 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 880 *stflag = SAME_NONZERO_PATTERN; 881 ierr = MatFDColoringApply(B,dm->fd,x,stflag,dm);CHKERRQ(ierr); 882 if (A != B) { 883 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 884 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 885 } 886 PetscFunctionReturn(0); 887 } 888 889 #undef __FUNCT__ 890 #define __FUNCT__ "DMCoarsen" 891 /*@ 892 DMCoarsen - Coarsens a DM object 893 894 Collective on DM 895 896 Input Parameter: 897 + dm - the DM object 898 - comm - the communicator to contain the new DM object (or PETSC_NULL) 899 900 Output Parameter: 901 . dmc - the coarsened DM 902 903 Level: developer 904 905 .seealso DMRefine(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 906 907 @*/ 908 PetscErrorCode DMCoarsen(DM dm, MPI_Comm comm, DM *dmc) 909 { 910 PetscErrorCode ierr; 911 912 PetscFunctionBegin; 913 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 914 ierr = (*dm->ops->coarsen)(dm, comm, dmc);CHKERRQ(ierr); 915 (*dmc)->ops->initialguess = dm->ops->initialguess; 916 (*dmc)->ops->function = dm->ops->function; 917 (*dmc)->ops->functionj = dm->ops->functionj; 918 if (dm->ops->jacobian != DMComputeJacobianDefault) { 919 (*dmc)->ops->jacobian = dm->ops->jacobian; 920 } 921 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)dm,(PetscObject)*dmc);CHKERRQ(ierr); 922 (*dmc)->ctx = dm->ctx; 923 (*dmc)->leveldown = dm->leveldown + 1; 924 PetscFunctionReturn(0); 925 } 926 927 #undef __FUNCT__ 928 #define __FUNCT__ "DMRefineHierarchy" 929 /*@C 930 DMRefineHierarchy - Refines a DM object, all levels at once 931 932 Collective on DM 933 934 Input Parameter: 935 + dm - the DM object 936 - nlevels - the number of levels of refinement 937 938 Output Parameter: 939 . dmf - the refined DM hierarchy 940 941 Level: developer 942 943 .seealso DMCoarsenHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 944 945 @*/ 946 PetscErrorCode DMRefineHierarchy(DM dm,PetscInt nlevels,DM dmf[]) 947 { 948 PetscErrorCode ierr; 949 950 PetscFunctionBegin; 951 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 952 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 953 if (nlevels == 0) PetscFunctionReturn(0); 954 if (dm->ops->refinehierarchy) { 955 ierr = (*dm->ops->refinehierarchy)(dm,nlevels,dmf);CHKERRQ(ierr); 956 } else if (dm->ops->refine) { 957 PetscInt i; 958 959 ierr = DMRefine(dm,((PetscObject)dm)->comm,&dmf[0]);CHKERRQ(ierr); 960 for (i=1; i<nlevels; i++) { 961 ierr = DMRefine(dmf[i-1],((PetscObject)dm)->comm,&dmf[i]);CHKERRQ(ierr); 962 } 963 } else { 964 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No RefineHierarchy for this DM yet"); 965 } 966 PetscFunctionReturn(0); 967 } 968 969 #undef __FUNCT__ 970 #define __FUNCT__ "DMCoarsenHierarchy" 971 /*@C 972 DMCoarsenHierarchy - Coarsens a DM object, all levels at once 973 974 Collective on DM 975 976 Input Parameter: 977 + dm - the DM object 978 - nlevels - the number of levels of coarsening 979 980 Output Parameter: 981 . dmc - the coarsened DM hierarchy 982 983 Level: developer 984 985 .seealso DMRefineHierarchy(), DMDestroy(), DMView(), DMCreateGlobalVector(), DMCreateInterpolation() 986 987 @*/ 988 PetscErrorCode DMCoarsenHierarchy(DM dm, PetscInt nlevels, DM dmc[]) 989 { 990 PetscErrorCode ierr; 991 992 PetscFunctionBegin; 993 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 994 if (nlevels < 0) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nlevels cannot be negative"); 995 if (nlevels == 0) PetscFunctionReturn(0); 996 PetscValidPointer(dmc,3); 997 if (dm->ops->coarsenhierarchy) { 998 ierr = (*dm->ops->coarsenhierarchy)(dm, nlevels, dmc);CHKERRQ(ierr); 999 } else if (dm->ops->coarsen) { 1000 PetscInt i; 1001 1002 ierr = DMCoarsen(dm,((PetscObject)dm)->comm,&dmc[0]);CHKERRQ(ierr); 1003 for (i=1; i<nlevels; i++) { 1004 ierr = DMCoarsen(dmc[i-1],((PetscObject)dm)->comm,&dmc[i]);CHKERRQ(ierr); 1005 } 1006 } else { 1007 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No CoarsenHierarchy for this DM yet"); 1008 } 1009 PetscFunctionReturn(0); 1010 } 1011 1012 #undef __FUNCT__ 1013 #define __FUNCT__ "DMCreateAggregates" 1014 /*@ 1015 DMCreateAggregates - Gets the aggregates that map between 1016 grids associated with two DMs. 1017 1018 Collective on DM 1019 1020 Input Parameters: 1021 + dmc - the coarse grid DM 1022 - dmf - the fine grid DM 1023 1024 Output Parameters: 1025 . rest - the restriction matrix (transpose of the projection matrix) 1026 1027 Level: intermediate 1028 1029 .keywords: interpolation, restriction, multigrid 1030 1031 .seealso: DMRefine(), DMCreateInjection(), DMCreateInterpolation() 1032 @*/ 1033 PetscErrorCode DMCreateAggregates(DM dmc, DM dmf, Mat *rest) 1034 { 1035 PetscErrorCode ierr; 1036 1037 PetscFunctionBegin; 1038 PetscValidHeaderSpecific(dmc,DM_CLASSID,1); 1039 PetscValidHeaderSpecific(dmf,DM_CLASSID,2); 1040 ierr = (*dmc->ops->getaggregates)(dmc, dmf, rest);CHKERRQ(ierr); 1041 PetscFunctionReturn(0); 1042 } 1043 1044 #undef __FUNCT__ 1045 #define __FUNCT__ "DMSetApplicationContextDestroy" 1046 /*@C 1047 DMSetApplicationContextDestroy - Sets a user function that will be called to destroy the application context when the DM is destroyed 1048 1049 Not Collective 1050 1051 Input Parameters: 1052 + dm - the DM object 1053 - destroy - the destroy function 1054 1055 Level: intermediate 1056 1057 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1058 1059 @*/ 1060 PetscErrorCode DMSetApplicationContextDestroy(DM dm,PetscErrorCode (*destroy)(void**)) 1061 { 1062 PetscFunctionBegin; 1063 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1064 dm->ctxdestroy = destroy; 1065 PetscFunctionReturn(0); 1066 } 1067 1068 #undef __FUNCT__ 1069 #define __FUNCT__ "DMSetApplicationContext" 1070 /*@ 1071 DMSetApplicationContext - Set a user context into a DM object 1072 1073 Not Collective 1074 1075 Input Parameters: 1076 + dm - the DM object 1077 - ctx - the user context 1078 1079 Level: intermediate 1080 1081 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1082 1083 @*/ 1084 PetscErrorCode DMSetApplicationContext(DM dm,void *ctx) 1085 { 1086 PetscFunctionBegin; 1087 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1088 dm->ctx = ctx; 1089 PetscFunctionReturn(0); 1090 } 1091 1092 #undef __FUNCT__ 1093 #define __FUNCT__ "DMGetApplicationContext" 1094 /*@ 1095 DMGetApplicationContext - Gets a user context from a DM object 1096 1097 Not Collective 1098 1099 Input Parameter: 1100 . dm - the DM object 1101 1102 Output Parameter: 1103 . ctx - the user context 1104 1105 Level: intermediate 1106 1107 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext() 1108 1109 @*/ 1110 PetscErrorCode DMGetApplicationContext(DM dm,void *ctx) 1111 { 1112 PetscFunctionBegin; 1113 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1114 *(void**)ctx = dm->ctx; 1115 PetscFunctionReturn(0); 1116 } 1117 1118 #undef __FUNCT__ 1119 #define __FUNCT__ "DMSetInitialGuess" 1120 /*@C 1121 DMSetInitialGuess - sets a function to compute an initial guess vector entries for the solvers 1122 1123 Logically Collective on DM 1124 1125 Input Parameter: 1126 + dm - the DM object to destroy 1127 - f - the function to compute the initial guess 1128 1129 Level: intermediate 1130 1131 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1132 1133 @*/ 1134 PetscErrorCode DMSetInitialGuess(DM dm,PetscErrorCode (*f)(DM,Vec)) 1135 { 1136 PetscFunctionBegin; 1137 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1138 dm->ops->initialguess = f; 1139 PetscFunctionReturn(0); 1140 } 1141 1142 #undef __FUNCT__ 1143 #define __FUNCT__ "DMSetFunction" 1144 /*@C 1145 DMSetFunction - sets a function to compute the right hand side vector entries for the KSP solver or nonlinear function for SNES 1146 1147 Logically Collective on DM 1148 1149 Input Parameter: 1150 + dm - the DM object 1151 - f - the function to compute (use PETSC_NULL to cancel a previous function that was set) 1152 1153 Level: intermediate 1154 1155 Notes: This sets both the function for function evaluations and the function used to compute Jacobians via finite differences if no Jacobian 1156 computer is provided with DMSetJacobian(). Canceling cancels the function, but not the function used to compute the Jacobian. 1157 1158 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1159 DMSetJacobian() 1160 1161 @*/ 1162 PetscErrorCode DMSetFunction(DM dm,PetscErrorCode (*f)(DM,Vec,Vec)) 1163 { 1164 PetscFunctionBegin; 1165 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1166 dm->ops->function = f; 1167 if (f) { 1168 dm->ops->functionj = f; 1169 } 1170 PetscFunctionReturn(0); 1171 } 1172 1173 #undef __FUNCT__ 1174 #define __FUNCT__ "DMSetJacobian" 1175 /*@C 1176 DMSetJacobian - sets a function to compute the matrix entries for the KSP solver or Jacobian for SNES 1177 1178 Logically Collective on DM 1179 1180 Input Parameter: 1181 + dm - the DM object to destroy 1182 - f - the function to compute the matrix entries 1183 1184 Level: intermediate 1185 1186 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1187 DMSetFunction() 1188 1189 @*/ 1190 PetscErrorCode DMSetJacobian(DM dm,PetscErrorCode (*f)(DM,Vec,Mat,Mat,MatStructure*)) 1191 { 1192 PetscFunctionBegin; 1193 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1194 dm->ops->jacobian = f; 1195 PetscFunctionReturn(0); 1196 } 1197 1198 #undef __FUNCT__ 1199 #define __FUNCT__ "DMComputeInitialGuess" 1200 /*@ 1201 DMComputeInitialGuess - computes an initial guess vector entries for the KSP solvers 1202 1203 Collective on DM 1204 1205 Input Parameter: 1206 + dm - the DM object to destroy 1207 - x - the vector to hold the initial guess values 1208 1209 Level: developer 1210 1211 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetRhs(), DMSetMat() 1212 1213 @*/ 1214 PetscErrorCode DMComputeInitialGuess(DM dm,Vec x) 1215 { 1216 PetscErrorCode ierr; 1217 PetscFunctionBegin; 1218 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1219 if (!dm->ops->initialguess) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetInitialGuess()"); 1220 ierr = (*dm->ops->initialguess)(dm,x);CHKERRQ(ierr); 1221 PetscFunctionReturn(0); 1222 } 1223 1224 #undef __FUNCT__ 1225 #define __FUNCT__ "DMHasInitialGuess" 1226 /*@ 1227 DMHasInitialGuess - does the DM object have an initial guess function 1228 1229 Not Collective 1230 1231 Input Parameter: 1232 . dm - the DM object to destroy 1233 1234 Output Parameter: 1235 . flg - PETSC_TRUE if function exists 1236 1237 Level: developer 1238 1239 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1240 1241 @*/ 1242 PetscErrorCode DMHasInitialGuess(DM dm,PetscBool *flg) 1243 { 1244 PetscFunctionBegin; 1245 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1246 *flg = (dm->ops->initialguess) ? PETSC_TRUE : PETSC_FALSE; 1247 PetscFunctionReturn(0); 1248 } 1249 1250 #undef __FUNCT__ 1251 #define __FUNCT__ "DMHasFunction" 1252 /*@ 1253 DMHasFunction - does the DM object have a function 1254 1255 Not Collective 1256 1257 Input Parameter: 1258 . dm - the DM object to destroy 1259 1260 Output Parameter: 1261 . flg - PETSC_TRUE if function exists 1262 1263 Level: developer 1264 1265 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1266 1267 @*/ 1268 PetscErrorCode DMHasFunction(DM dm,PetscBool *flg) 1269 { 1270 PetscFunctionBegin; 1271 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1272 *flg = (dm->ops->function) ? PETSC_TRUE : PETSC_FALSE; 1273 PetscFunctionReturn(0); 1274 } 1275 1276 #undef __FUNCT__ 1277 #define __FUNCT__ "DMHasJacobian" 1278 /*@ 1279 DMHasJacobian - does the DM object have a matrix function 1280 1281 Not Collective 1282 1283 Input Parameter: 1284 . dm - the DM object to destroy 1285 1286 Output Parameter: 1287 . flg - PETSC_TRUE if function exists 1288 1289 Level: developer 1290 1291 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetFunction(), DMSetJacobian() 1292 1293 @*/ 1294 PetscErrorCode DMHasJacobian(DM dm,PetscBool *flg) 1295 { 1296 PetscFunctionBegin; 1297 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1298 *flg = (dm->ops->jacobian) ? PETSC_TRUE : PETSC_FALSE; 1299 PetscFunctionReturn(0); 1300 } 1301 1302 #undef __FUNCT__ 1303 #define __FUNCT__ "DMComputeFunction" 1304 /*@ 1305 DMComputeFunction - computes the right hand side vector entries for the KSP solver or nonlinear function for SNES 1306 1307 Collective on DM 1308 1309 Input Parameter: 1310 + dm - the DM object to destroy 1311 . x - the location where the function is evaluationed, may be ignored for linear problems 1312 - b - the vector to hold the right hand side entries 1313 1314 Level: developer 1315 1316 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1317 DMSetJacobian() 1318 1319 @*/ 1320 PetscErrorCode DMComputeFunction(DM dm,Vec x,Vec b) 1321 { 1322 PetscErrorCode ierr; 1323 PetscFunctionBegin; 1324 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1325 if (!dm->ops->function) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to provide function with DMSetFunction()"); 1326 PetscStackPush("DM user function"); 1327 ierr = (*dm->ops->function)(dm,x,b);CHKERRQ(ierr); 1328 PetscStackPop; 1329 PetscFunctionReturn(0); 1330 } 1331 1332 1333 #undef __FUNCT__ 1334 #define __FUNCT__ "DMComputeJacobian" 1335 /*@ 1336 DMComputeJacobian - compute the matrix entries for the solver 1337 1338 Collective on DM 1339 1340 Input Parameter: 1341 + dm - the DM object 1342 . x - location to compute Jacobian at; will be PETSC_NULL for linear problems, for nonlinear problems if not provided then pulled from DM 1343 . A - matrix that defines the operator for the linear solve 1344 - B - the matrix used to construct the preconditioner 1345 1346 Level: developer 1347 1348 .seealso DMView(), DMCreateGlobalVector(), DMCreateInterpolation(), DMCreateColoring(), DMCreateMatrix(), DMGetApplicationContext(), DMSetInitialGuess(), 1349 DMSetFunction() 1350 1351 @*/ 1352 PetscErrorCode DMComputeJacobian(DM dm,Vec x,Mat A,Mat B,MatStructure *stflag) 1353 { 1354 PetscErrorCode ierr; 1355 1356 PetscFunctionBegin; 1357 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1358 if (!dm->ops->jacobian) { 1359 ISColoring coloring; 1360 MatFDColoring fd; 1361 1362 ierr = DMCreateColoring(dm,dm->coloringtype,MATAIJ,&coloring);CHKERRQ(ierr); 1363 ierr = MatFDColoringCreate(B,coloring,&fd);CHKERRQ(ierr); 1364 ierr = ISColoringDestroy(&coloring);CHKERRQ(ierr); 1365 ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))dm->ops->functionj,dm);CHKERRQ(ierr); 1366 ierr = PetscObjectSetOptionsPrefix((PetscObject)fd,((PetscObject)dm)->prefix);CHKERRQ(ierr); 1367 ierr = MatFDColoringSetFromOptions(fd);CHKERRQ(ierr); 1368 1369 dm->fd = fd; 1370 dm->ops->jacobian = DMComputeJacobianDefault; 1371 1372 /* don't know why this is needed */ 1373 ierr = PetscObjectDereference((PetscObject)dm);CHKERRQ(ierr); 1374 } 1375 if (!x) x = dm->x; 1376 ierr = (*dm->ops->jacobian)(dm,x,A,B,stflag);CHKERRQ(ierr); 1377 1378 /* if matrix depends on x; i.e. nonlinear problem, keep copy of input vector since needed by multigrid methods to generate coarse grid matrices */ 1379 if (x) { 1380 if (!dm->x) { 1381 ierr = DMCreateGlobalVector(dm,&dm->x);CHKERRQ(ierr); 1382 } 1383 ierr = VecCopy(x,dm->x);CHKERRQ(ierr); 1384 } 1385 if (A != B) { 1386 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1387 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1388 } 1389 PetscFunctionReturn(0); 1390 } 1391 1392 1393 PetscFList DMList = PETSC_NULL; 1394 PetscBool DMRegisterAllCalled = PETSC_FALSE; 1395 1396 #undef __FUNCT__ 1397 #define __FUNCT__ "DMSetType" 1398 /*@C 1399 DMSetType - Builds a DM, for a particular DM implementation. 1400 1401 Collective on DM 1402 1403 Input Parameters: 1404 + dm - The DM object 1405 - method - The name of the DM type 1406 1407 Options Database Key: 1408 . -dm_type <type> - Sets the DM type; use -help for a list of available types 1409 1410 Notes: 1411 See "petsc/include/petscdm.h" for available DM types (for instance, DM1D, DM2D, or DM3D). 1412 1413 Level: intermediate 1414 1415 .keywords: DM, set, type 1416 .seealso: DMGetType(), DMCreate() 1417 @*/ 1418 PetscErrorCode DMSetType(DM dm, const DMType method) 1419 { 1420 PetscErrorCode (*r)(DM); 1421 PetscBool match; 1422 PetscErrorCode ierr; 1423 1424 PetscFunctionBegin; 1425 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1426 ierr = PetscTypeCompare((PetscObject) dm, method, &match);CHKERRQ(ierr); 1427 if (match) PetscFunctionReturn(0); 1428 1429 if (!DMRegisterAllCalled) {ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1430 ierr = PetscFListFind(DMList, ((PetscObject)dm)->comm, method,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 1431 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown DM type: %s", method); 1432 1433 if (dm->ops->destroy) { 1434 ierr = (*dm->ops->destroy)(dm);CHKERRQ(ierr); 1435 } 1436 ierr = (*r)(dm);CHKERRQ(ierr); 1437 ierr = PetscObjectChangeTypeName((PetscObject)dm,method);CHKERRQ(ierr); 1438 PetscFunctionReturn(0); 1439 } 1440 1441 #undef __FUNCT__ 1442 #define __FUNCT__ "DMGetType" 1443 /*@C 1444 DMGetType - Gets the DM type name (as a string) from the DM. 1445 1446 Not Collective 1447 1448 Input Parameter: 1449 . dm - The DM 1450 1451 Output Parameter: 1452 . type - The DM type name 1453 1454 Level: intermediate 1455 1456 .keywords: DM, get, type, name 1457 .seealso: DMSetType(), DMCreate() 1458 @*/ 1459 PetscErrorCode DMGetType(DM dm, const DMType *type) 1460 { 1461 PetscErrorCode ierr; 1462 1463 PetscFunctionBegin; 1464 PetscValidHeaderSpecific(dm, DM_CLASSID,1); 1465 PetscValidCharPointer(type,2); 1466 if (!DMRegisterAllCalled) { 1467 ierr = DMRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1468 } 1469 *type = ((PetscObject)dm)->type_name; 1470 PetscFunctionReturn(0); 1471 } 1472 1473 #undef __FUNCT__ 1474 #define __FUNCT__ "DMConvert" 1475 /*@C 1476 DMConvert - Converts a DM to another DM, either of the same or different type. 1477 1478 Collective on DM 1479 1480 Input Parameters: 1481 + dm - the DM 1482 - newtype - new DM type (use "same" for the same type) 1483 1484 Output Parameter: 1485 . M - pointer to new DM 1486 1487 Notes: 1488 Cannot be used to convert a sequential DM to parallel or parallel to sequential, 1489 the MPI communicator of the generated DM is always the same as the communicator 1490 of the input DM. 1491 1492 Level: intermediate 1493 1494 .seealso: DMCreate() 1495 @*/ 1496 PetscErrorCode DMConvert(DM dm, const DMType newtype, DM *M) 1497 { 1498 DM B; 1499 char convname[256]; 1500 PetscBool sametype, issame; 1501 PetscErrorCode ierr; 1502 1503 PetscFunctionBegin; 1504 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1505 PetscValidType(dm,1); 1506 PetscValidPointer(M,3); 1507 ierr = PetscTypeCompare((PetscObject) dm, newtype, &sametype);CHKERRQ(ierr); 1508 ierr = PetscStrcmp(newtype, "same", &issame);CHKERRQ(ierr); 1509 { 1510 PetscErrorCode (*conv)(DM, const DMType, DM *) = PETSC_NULL; 1511 1512 /* 1513 Order of precedence: 1514 1) See if a specialized converter is known to the current DM. 1515 2) See if a specialized converter is known to the desired DM class. 1516 3) See if a good general converter is registered for the desired class 1517 4) See if a good general converter is known for the current matrix. 1518 5) Use a really basic converter. 1519 */ 1520 1521 /* 1) See if a specialized converter is known to the current DM and the desired class */ 1522 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 1523 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 1524 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 1525 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 1526 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 1527 ierr = PetscObjectQueryFunction((PetscObject)dm,convname,(void (**)(void))&conv);CHKERRQ(ierr); 1528 if (conv) goto foundconv; 1529 1530 /* 2) See if a specialized converter is known to the desired DM class. */ 1531 ierr = DMCreate(((PetscObject) dm)->comm, &B);CHKERRQ(ierr); 1532 ierr = DMSetType(B, newtype);CHKERRQ(ierr); 1533 ierr = PetscStrcpy(convname,"DMConvert_");CHKERRQ(ierr); 1534 ierr = PetscStrcat(convname,((PetscObject) dm)->type_name);CHKERRQ(ierr); 1535 ierr = PetscStrcat(convname,"_");CHKERRQ(ierr); 1536 ierr = PetscStrcat(convname,newtype);CHKERRQ(ierr); 1537 ierr = PetscStrcat(convname,"_C");CHKERRQ(ierr); 1538 ierr = PetscObjectQueryFunction((PetscObject)B,convname,(void (**)(void))&conv);CHKERRQ(ierr); 1539 if (conv) { 1540 ierr = DMDestroy(&B);CHKERRQ(ierr); 1541 goto foundconv; 1542 } 1543 1544 #if 0 1545 /* 3) See if a good general converter is registered for the desired class */ 1546 conv = B->ops->convertfrom; 1547 ierr = DMDestroy(&B);CHKERRQ(ierr); 1548 if (conv) goto foundconv; 1549 1550 /* 4) See if a good general converter is known for the current matrix */ 1551 if (dm->ops->convert) { 1552 conv = dm->ops->convert; 1553 } 1554 if (conv) goto foundconv; 1555 #endif 1556 1557 /* 5) Use a really basic converter. */ 1558 SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_SUP, "No conversion possible between DM types %s and %s", ((PetscObject) dm)->type_name, newtype); 1559 1560 foundconv: 1561 ierr = PetscLogEventBegin(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 1562 ierr = (*conv)(dm,newtype,M);CHKERRQ(ierr); 1563 ierr = PetscLogEventEnd(DM_Convert,dm,0,0,0);CHKERRQ(ierr); 1564 } 1565 ierr = PetscObjectStateIncrease((PetscObject) *M);CHKERRQ(ierr); 1566 PetscFunctionReturn(0); 1567 } 1568 1569 /*--------------------------------------------------------------------------------------------------------------------*/ 1570 1571 #undef __FUNCT__ 1572 #define __FUNCT__ "DMRegister" 1573 /*@C 1574 DMRegister - See DMRegisterDynamic() 1575 1576 Level: advanced 1577 @*/ 1578 PetscErrorCode DMRegister(const char sname[], const char path[], const char name[], PetscErrorCode (*function)(DM)) 1579 { 1580 char fullname[PETSC_MAX_PATH_LEN]; 1581 PetscErrorCode ierr; 1582 1583 PetscFunctionBegin; 1584 ierr = PetscStrcpy(fullname, path);CHKERRQ(ierr); 1585 ierr = PetscStrcat(fullname, ":");CHKERRQ(ierr); 1586 ierr = PetscStrcat(fullname, name);CHKERRQ(ierr); 1587 ierr = PetscFListAdd(&DMList, sname, fullname, (void (*)(void)) function);CHKERRQ(ierr); 1588 PetscFunctionReturn(0); 1589 } 1590 1591 1592 /*--------------------------------------------------------------------------------------------------------------------*/ 1593 #undef __FUNCT__ 1594 #define __FUNCT__ "DMRegisterDestroy" 1595 /*@C 1596 DMRegisterDestroy - Frees the list of DM methods that were registered by DMRegister()/DMRegisterDynamic(). 1597 1598 Not Collective 1599 1600 Level: advanced 1601 1602 .keywords: DM, register, destroy 1603 .seealso: DMRegister(), DMRegisterAll(), DMRegisterDynamic() 1604 @*/ 1605 PetscErrorCode DMRegisterDestroy(void) 1606 { 1607 PetscErrorCode ierr; 1608 1609 PetscFunctionBegin; 1610 ierr = PetscFListDestroy(&DMList);CHKERRQ(ierr); 1611 DMRegisterAllCalled = PETSC_FALSE; 1612 PetscFunctionReturn(0); 1613 } 1614 1615 #if defined(PETSC_HAVE_MATLAB_ENGINE) 1616 #include <mex.h> 1617 1618 typedef struct {char *funcname; char *jacname; mxArray *ctx;} DMMatlabContext; 1619 1620 #undef __FUNCT__ 1621 #define __FUNCT__ "DMComputeFunction_Matlab" 1622 /* 1623 DMComputeFunction_Matlab - Calls the function that has been set with 1624 DMSetFunctionMatlab(). 1625 1626 For linear problems x is null 1627 1628 .seealso: DMSetFunction(), DMGetFunction() 1629 */ 1630 PetscErrorCode DMComputeFunction_Matlab(DM dm,Vec x,Vec y) 1631 { 1632 PetscErrorCode ierr; 1633 DMMatlabContext *sctx; 1634 int nlhs = 1,nrhs = 4; 1635 mxArray *plhs[1],*prhs[4]; 1636 long long int lx = 0,ly = 0,ls = 0; 1637 1638 PetscFunctionBegin; 1639 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1640 PetscValidHeaderSpecific(y,VEC_CLASSID,3); 1641 PetscCheckSameComm(dm,1,y,3); 1642 1643 /* call Matlab function in ctx with arguments x and y */ 1644 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1645 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 1646 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 1647 ierr = PetscMemcpy(&ly,&y,sizeof(y));CHKERRQ(ierr); 1648 prhs[0] = mxCreateDoubleScalar((double)ls); 1649 prhs[1] = mxCreateDoubleScalar((double)lx); 1650 prhs[2] = mxCreateDoubleScalar((double)ly); 1651 prhs[3] = mxCreateString(sctx->funcname); 1652 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeFunctionInternal");CHKERRQ(ierr); 1653 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 1654 mxDestroyArray(prhs[0]); 1655 mxDestroyArray(prhs[1]); 1656 mxDestroyArray(prhs[2]); 1657 mxDestroyArray(prhs[3]); 1658 mxDestroyArray(plhs[0]); 1659 PetscFunctionReturn(0); 1660 } 1661 1662 1663 #undef __FUNCT__ 1664 #define __FUNCT__ "DMSetFunctionMatlab" 1665 /* 1666 DMSetFunctionMatlab - Sets the function evaluation routine 1667 1668 */ 1669 PetscErrorCode DMSetFunctionMatlab(DM dm,const char *func) 1670 { 1671 PetscErrorCode ierr; 1672 DMMatlabContext *sctx; 1673 1674 PetscFunctionBegin; 1675 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1676 /* currently sctx is memory bleed */ 1677 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1678 if (!sctx) { 1679 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 1680 } 1681 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 1682 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 1683 ierr = DMSetFunction(dm,DMComputeFunction_Matlab);CHKERRQ(ierr); 1684 PetscFunctionReturn(0); 1685 } 1686 1687 #undef __FUNCT__ 1688 #define __FUNCT__ "DMComputeJacobian_Matlab" 1689 /* 1690 DMComputeJacobian_Matlab - Calls the function that has been set with 1691 DMSetJacobianMatlab(). 1692 1693 For linear problems x is null 1694 1695 .seealso: DMSetFunction(), DMGetFunction() 1696 */ 1697 PetscErrorCode DMComputeJacobian_Matlab(DM dm,Vec x,Mat A,Mat B,MatStructure *str) 1698 { 1699 PetscErrorCode ierr; 1700 DMMatlabContext *sctx; 1701 int nlhs = 2,nrhs = 5; 1702 mxArray *plhs[2],*prhs[5]; 1703 long long int lx = 0,lA = 0,lB = 0,ls = 0; 1704 1705 PetscFunctionBegin; 1706 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1707 PetscValidHeaderSpecific(A,MAT_CLASSID,3); 1708 1709 /* call MATLAB function in ctx with arguments x, A, and B */ 1710 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1711 ierr = PetscMemcpy(&ls,&dm,sizeof(dm));CHKERRQ(ierr); 1712 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 1713 ierr = PetscMemcpy(&lA,&A,sizeof(A));CHKERRQ(ierr); 1714 ierr = PetscMemcpy(&lB,&B,sizeof(B));CHKERRQ(ierr); 1715 prhs[0] = mxCreateDoubleScalar((double)ls); 1716 prhs[1] = mxCreateDoubleScalar((double)lx); 1717 prhs[2] = mxCreateDoubleScalar((double)lA); 1718 prhs[3] = mxCreateDoubleScalar((double)lB); 1719 prhs[4] = mxCreateString(sctx->jacname); 1720 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscDMComputeJacobianInternal");CHKERRQ(ierr); 1721 *str = (MatStructure) mxGetScalar(plhs[0]); 1722 ierr = (PetscInt) mxGetScalar(plhs[1]);CHKERRQ(ierr); 1723 mxDestroyArray(prhs[0]); 1724 mxDestroyArray(prhs[1]); 1725 mxDestroyArray(prhs[2]); 1726 mxDestroyArray(prhs[3]); 1727 mxDestroyArray(prhs[4]); 1728 mxDestroyArray(plhs[0]); 1729 mxDestroyArray(plhs[1]); 1730 PetscFunctionReturn(0); 1731 } 1732 1733 1734 #undef __FUNCT__ 1735 #define __FUNCT__ "DMSetJacobianMatlab" 1736 /* 1737 DMSetJacobianMatlab - Sets the Jacobian function evaluation routine 1738 1739 */ 1740 PetscErrorCode DMSetJacobianMatlab(DM dm,const char *func) 1741 { 1742 PetscErrorCode ierr; 1743 DMMatlabContext *sctx; 1744 1745 PetscFunctionBegin; 1746 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1747 /* currently sctx is memory bleed */ 1748 ierr = DMGetApplicationContext(dm,&sctx);CHKERRQ(ierr); 1749 if (!sctx) { 1750 ierr = PetscMalloc(sizeof(DMMatlabContext),&sctx);CHKERRQ(ierr); 1751 } 1752 ierr = PetscStrallocpy(func,&sctx->jacname);CHKERRQ(ierr); 1753 ierr = DMSetApplicationContext(dm,sctx);CHKERRQ(ierr); 1754 ierr = DMSetJacobian(dm,DMComputeJacobian_Matlab);CHKERRQ(ierr); 1755 PetscFunctionReturn(0); 1756 } 1757 #endif 1758 1759 #undef __FUNCT__ 1760 #define __FUNCT__ "DMLoad" 1761 /*@C 1762 DMLoad - Loads a DM that has been stored in binary or HDF5 format 1763 with DMView(). 1764 1765 Collective on PetscViewer 1766 1767 Input Parameters: 1768 + newdm - the newly loaded DM, this needs to have been created with DMCreate() or 1769 some related function before a call to DMLoad(). 1770 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 1771 HDF5 file viewer, obtained from PetscViewerHDF5Open() 1772 1773 Level: intermediate 1774 1775 Notes: 1776 Defaults to the DM DA. 1777 1778 Notes for advanced users: 1779 Most users should not need to know the details of the binary storage 1780 format, since DMLoad() and DMView() completely hide these details. 1781 But for anyone who's interested, the standard binary matrix storage 1782 format is 1783 .vb 1784 has not yet been determined 1785 .ve 1786 1787 In addition, PETSc automatically does the byte swapping for 1788 machines that store the bytes reversed, e.g. DEC alpha, freebsd, 1789 linux, Windows and the paragon; thus if you write your own binary 1790 read/write routines you have to swap the bytes; see PetscBinaryRead() 1791 and PetscBinaryWrite() to see how this may be done. 1792 1793 Concepts: vector^loading from file 1794 1795 .seealso: PetscViewerBinaryOpen(), DMView(), MatLoad(), VecLoad() 1796 @*/ 1797 PetscErrorCode DMLoad(DM newdm, PetscViewer viewer) 1798 { 1799 PetscErrorCode ierr; 1800 1801 PetscFunctionBegin; 1802 PetscValidHeaderSpecific(newdm,DM_CLASSID,1); 1803 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 1804 1805 if (!((PetscObject)newdm)->type_name) { 1806 ierr = DMSetType(newdm, DMDA);CHKERRQ(ierr); 1807 } 1808 ierr = (*newdm->ops->load)(newdm,viewer);CHKERRQ(ierr); 1809 PetscFunctionReturn(0); 1810 } 1811 1812