1 #include <petscdmshell.h> /*I "petscdmshell.h" I*/ 2 #include <petscmat.h> 3 #include <petsc/private/dmimpl.h> 4 5 typedef struct { 6 Vec Xglobal; 7 Vec Xlocal; 8 Mat A; 9 VecScatter gtol; 10 VecScatter ltog; 11 VecScatter ltol; 12 void *ctx; 13 } DM_Shell; 14 15 /*@ 16 DMGlobalToLocalBeginDefaultShell - Uses the GlobalToLocal VecScatter context set by the user to begin a global to local scatter 17 Collective 18 19 Input Arguments: 20 + dm - shell DM 21 . g - global vector 22 . mode - InsertMode 23 - l - local vector 24 25 Level: advanced 26 27 Note: This is not normally called directly by user code, generally user code calls DMGlobalToLocalBegin() and DMGlobalToLocalEnd(). If the user provides their own custom routines to DMShellSetLocalToGlobal() then those routines might have reason to call this function. 28 29 .seealso: DMGlobalToLocalEndDefaultShell() 30 @*/ 31 PetscErrorCode DMGlobalToLocalBeginDefaultShell(DM dm,Vec g,InsertMode mode,Vec l) 32 { 33 PetscErrorCode ierr; 34 DM_Shell *shell = (DM_Shell*)dm->data; 35 36 PetscFunctionBegin; 37 if (!shell->gtol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()"); 38 ierr = VecScatterBegin(shell->gtol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 39 PetscFunctionReturn(0); 40 } 41 42 /*@ 43 DMGlobalToLocalEndDefaultShell - Uses the GlobalToLocal VecScatter context set by the user to end a global to local scatter 44 Collective 45 46 Input Arguments: 47 + dm - shell DM 48 . g - global vector 49 . mode - InsertMode 50 - l - local vector 51 52 Level: advanced 53 54 .seealso: DMGlobalToLocalBeginDefaultShell() 55 @*/ 56 PetscErrorCode DMGlobalToLocalEndDefaultShell(DM dm,Vec g,InsertMode mode,Vec l) 57 { 58 PetscErrorCode ierr; 59 DM_Shell *shell = (DM_Shell*)dm->data; 60 61 PetscFunctionBegin; 62 if (!shell->gtol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()"); 63 ierr = VecScatterEnd(shell->gtol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 64 PetscFunctionReturn(0); 65 } 66 67 /*@ 68 DMLocalToGlobalBeginDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to begin a local to global scatter 69 Collective 70 71 Input Arguments: 72 + dm - shell DM 73 . l - local vector 74 . mode - InsertMode 75 - g - global vector 76 77 Level: advanced 78 79 Note: This is not normally called directly by user code, generally user code calls DMLocalToGlobalBegin() and DMLocalToGlobalEnd(). If the user provides their own custom routines to DMShellSetLocalToGlobal() then those routines might have reason to call this function. 80 81 .seealso: DMLocalToGlobalEndDefaultShell() 82 @*/ 83 PetscErrorCode DMLocalToGlobalBeginDefaultShell(DM dm,Vec l,InsertMode mode,Vec g) 84 { 85 PetscErrorCode ierr; 86 DM_Shell *shell = (DM_Shell*)dm->data; 87 88 PetscFunctionBegin; 89 if (!shell->ltog) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()"); 90 ierr = VecScatterBegin(shell->ltog,l,g,mode,SCATTER_FORWARD);CHKERRQ(ierr); 91 PetscFunctionReturn(0); 92 } 93 94 /*@ 95 DMLocalToGlobalEndDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to end a local to global scatter 96 Collective 97 98 Input Arguments: 99 + dm - shell DM 100 . l - local vector 101 . mode - InsertMode 102 - g - global vector 103 104 Level: advanced 105 106 .seealso: DMLocalToGlobalBeginDefaultShell() 107 @*/ 108 PetscErrorCode DMLocalToGlobalEndDefaultShell(DM dm,Vec l,InsertMode mode,Vec g) 109 { 110 PetscErrorCode ierr; 111 DM_Shell *shell = (DM_Shell*)dm->data; 112 113 PetscFunctionBegin; 114 if (!shell->ltog) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()"); 115 ierr = VecScatterEnd(shell->ltog,l,g,mode,SCATTER_FORWARD);CHKERRQ(ierr); 116 PetscFunctionReturn(0); 117 } 118 119 /*@ 120 DMLocalToLocalBeginDefaultShell - Uses the LocalToLocal VecScatter context set by the user to begin a local to local scatter 121 Collective 122 123 Input Arguments: 124 + dm - shell DM 125 . g - the original local vector 126 - mode - InsertMode 127 128 Output Parameter: 129 . l - the local vector with correct ghost values 130 131 Level: advanced 132 133 Note: This is not normally called directly by user code, generally user code calls DMLocalToLocalBegin() and DMLocalToLocalEnd(). If the user provides their own custom routines to DMShellSetLocalToLocal() then those routines might have reason to call this function. 134 135 .seealso: DMLocalToLocalEndDefaultShell() 136 @*/ 137 PetscErrorCode DMLocalToLocalBeginDefaultShell(DM dm,Vec g,InsertMode mode,Vec l) 138 { 139 PetscErrorCode ierr; 140 DM_Shell *shell = (DM_Shell*)dm->data; 141 142 PetscFunctionBegin; 143 if (!shell->ltol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToLocalVecScatter()"); 144 ierr = VecScatterBegin(shell->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 145 PetscFunctionReturn(0); 146 } 147 148 /*@ 149 DMLocalToLocalEndDefaultShell - Uses the LocalToLocal VecScatter context set by the user to end a local to local scatter 150 Collective 151 152 Input Arguments: 153 + dm - shell DM 154 . g - the original local vector 155 - mode - InsertMode 156 157 Output Parameter: 158 . l - the local vector with correct ghost values 159 160 Level: advanced 161 162 .seealso: DMLocalToLocalBeginDefaultShell() 163 @*/ 164 PetscErrorCode DMLocalToLocalEndDefaultShell(DM dm,Vec g,InsertMode mode,Vec l) 165 { 166 PetscErrorCode ierr; 167 DM_Shell *shell = (DM_Shell*)dm->data; 168 169 PetscFunctionBegin; 170 if (!shell->ltol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()"); 171 ierr = VecScatterEnd(shell->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 172 PetscFunctionReturn(0); 173 } 174 175 static PetscErrorCode DMCreateMatrix_Shell(DM dm,Mat *J) 176 { 177 PetscErrorCode ierr; 178 DM_Shell *shell = (DM_Shell*)dm->data; 179 Mat A; 180 181 PetscFunctionBegin; 182 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 183 PetscValidPointer(J,2); 184 if (!shell->A) { 185 if (shell->Xglobal) { 186 PetscInt m,M; 187 ierr = PetscInfo(dm,"Naively creating matrix using global vector distribution without preallocation\n");CHKERRQ(ierr); 188 ierr = VecGetSize(shell->Xglobal,&M);CHKERRQ(ierr); 189 ierr = VecGetLocalSize(shell->Xglobal,&m);CHKERRQ(ierr); 190 ierr = MatCreate(PetscObjectComm((PetscObject)dm),&shell->A);CHKERRQ(ierr); 191 ierr = MatSetSizes(shell->A,m,m,M,M);CHKERRQ(ierr); 192 ierr = MatSetType(shell->A,dm->mattype);CHKERRQ(ierr); 193 ierr = MatSetUp(shell->A);CHKERRQ(ierr); 194 } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetMatrix(), DMShellSetCreateMatrix(), or provide a vector"); 195 } 196 A = shell->A; 197 ierr = MatDuplicate(A,MAT_SHARE_NONZERO_PATTERN,J);CHKERRQ(ierr); 198 ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 199 PetscFunctionReturn(0); 200 } 201 202 PetscErrorCode DMCreateGlobalVector_Shell(DM dm,Vec *gvec) 203 { 204 PetscErrorCode ierr; 205 DM_Shell *shell = (DM_Shell*)dm->data; 206 Vec X; 207 208 PetscFunctionBegin; 209 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 210 PetscValidPointer(gvec,2); 211 *gvec = NULL; 212 X = shell->Xglobal; 213 if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetGlobalVector() or DMShellSetCreateGlobalVector()"); 214 /* Need to create a copy in order to attach the DM to the vector */ 215 ierr = VecDuplicate(X,gvec);CHKERRQ(ierr); 216 ierr = VecZeroEntries(*gvec);CHKERRQ(ierr); 217 ierr = VecSetDM(*gvec,dm);CHKERRQ(ierr); 218 PetscFunctionReturn(0); 219 } 220 221 PetscErrorCode DMCreateLocalVector_Shell(DM dm,Vec *gvec) 222 { 223 PetscErrorCode ierr; 224 DM_Shell *shell = (DM_Shell*)dm->data; 225 Vec X; 226 227 PetscFunctionBegin; 228 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 229 PetscValidPointer(gvec,2); 230 *gvec = NULL; 231 X = shell->Xlocal; 232 if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetLocalVector() or DMShellSetCreateLocalVector()"); 233 /* Need to create a copy in order to attach the DM to the vector */ 234 ierr = VecDuplicate(X,gvec);CHKERRQ(ierr); 235 ierr = VecZeroEntries(*gvec);CHKERRQ(ierr); 236 ierr = VecSetDM(*gvec,dm);CHKERRQ(ierr); 237 PetscFunctionReturn(0); 238 } 239 240 /*@ 241 DMShellSetContext - set some data to be usable by this DM 242 243 Collective 244 245 Input Arguments: 246 + dm - shell DM 247 - ctx - the context 248 249 Level: advanced 250 251 .seealso: DMCreateMatrix(), DMShellGetContext() 252 @*/ 253 PetscErrorCode DMShellSetContext(DM dm,void *ctx) 254 { 255 DM_Shell *shell = (DM_Shell*)dm->data; 256 PetscErrorCode ierr; 257 PetscBool isshell; 258 259 PetscFunctionBegin; 260 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 261 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 262 if (!isshell) PetscFunctionReturn(0); 263 shell->ctx = ctx; 264 PetscFunctionReturn(0); 265 } 266 267 /*@ 268 DMShellGetContext - Returns the user-provided context associated to the DM 269 270 Collective 271 272 Input Argument: 273 . dm - shell DM 274 275 Output Argument: 276 . ctx - the context 277 278 Level: advanced 279 280 .seealso: DMCreateMatrix(), DMShellSetContext() 281 @*/ 282 PetscErrorCode DMShellGetContext(DM dm,void *ctx) 283 { 284 DM_Shell *shell = (DM_Shell*)dm->data; 285 PetscErrorCode ierr; 286 PetscBool isshell; 287 288 PetscFunctionBegin; 289 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 290 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 291 if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 292 *(void**)ctx = shell->ctx; 293 PetscFunctionReturn(0); 294 } 295 296 /*@ 297 DMShellSetMatrix - sets a template matrix associated with the DMShell 298 299 Collective 300 301 Input Arguments: 302 + dm - shell DM 303 - J - template matrix 304 305 Level: advanced 306 307 Developer Notes: 308 To avoid circular references, if J is already associated to the same DM, then MatDuplicate(SHARE_NONZERO_PATTERN) is called, followed by removing the DM reference from the private template. 309 310 .seealso: DMCreateMatrix(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext() 311 @*/ 312 PetscErrorCode DMShellSetMatrix(DM dm,Mat J) 313 { 314 DM_Shell *shell = (DM_Shell*)dm->data; 315 PetscErrorCode ierr; 316 PetscBool isshell; 317 DM mdm; 318 319 PetscFunctionBegin; 320 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 321 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 322 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 323 if (!isshell) PetscFunctionReturn(0); 324 if (J == shell->A) PetscFunctionReturn(0); 325 ierr = MatGetDM(J,&mdm);CHKERRQ(ierr); 326 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 327 ierr = MatDestroy(&shell->A);CHKERRQ(ierr); 328 if (mdm == dm) { 329 ierr = MatDuplicate(J,MAT_SHARE_NONZERO_PATTERN,&shell->A);CHKERRQ(ierr); 330 ierr = MatSetDM(shell->A,NULL);CHKERRQ(ierr); 331 } else shell->A = J; 332 PetscFunctionReturn(0); 333 } 334 335 /*@C 336 DMShellSetCreateMatrix - sets the routine to create a matrix associated with the shell DM 337 338 Logically Collective on dm 339 340 Input Arguments: 341 + dm - the shell DM 342 - func - the function to create a matrix 343 344 Level: advanced 345 346 .seealso: DMCreateMatrix(), DMShellSetMatrix(), DMShellSetContext(), DMShellGetContext() 347 @*/ 348 PetscErrorCode DMShellSetCreateMatrix(DM dm,PetscErrorCode (*func)(DM,Mat*)) 349 { 350 PetscFunctionBegin; 351 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 352 dm->ops->creatematrix = func; 353 PetscFunctionReturn(0); 354 } 355 356 /*@ 357 DMShellSetGlobalVector - sets a template global vector associated with the DMShell 358 359 Logically Collective on dm 360 361 Input Arguments: 362 + dm - shell DM 363 - X - template vector 364 365 Level: advanced 366 367 .seealso: DMCreateGlobalVector(), DMShellSetMatrix(), DMShellSetCreateGlobalVector() 368 @*/ 369 PetscErrorCode DMShellSetGlobalVector(DM dm,Vec X) 370 { 371 DM_Shell *shell = (DM_Shell*)dm->data; 372 PetscErrorCode ierr; 373 PetscBool isshell; 374 DM vdm; 375 376 PetscFunctionBegin; 377 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 378 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 379 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 380 if (!isshell) PetscFunctionReturn(0); 381 ierr = VecGetDM(X,&vdm);CHKERRQ(ierr); 382 /* 383 if the vector proposed as the new base global vector for the DM is a DM vector associated 384 with the same DM then the current base global vector for the DM is ok and if we replace it with the new one 385 we get a circular dependency that prevents the DM from being destroy when it should be. 386 This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided 387 DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries 388 to set its input vector (which is associated with the DM) as the base global vector. 389 Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien 390 for pointing out the problem. 391 */ 392 if (vdm == dm) PetscFunctionReturn(0); 393 ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 394 ierr = VecDestroy(&shell->Xglobal);CHKERRQ(ierr); 395 shell->Xglobal = X; 396 PetscFunctionReturn(0); 397 } 398 399 /*@ 400 DMShellGetGlobalVector - Returns the template global vector associated with the DMShell, or NULL if it was not set 401 402 Not collective 403 404 Input Arguments: 405 + dm - shell DM 406 - X - template vector 407 408 Level: advanced 409 410 .seealso: DMShellSetGlobalVector(), DMShellSetCreateGlobalVector(), DMCreateGlobalVector() 411 @*/ 412 PetscErrorCode DMShellGetGlobalVector(DM dm, Vec *X) 413 { 414 DM_Shell *shell = (DM_Shell *) dm->data; 415 PetscBool isshell; 416 PetscErrorCode ierr; 417 418 PetscFunctionBegin; 419 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 420 PetscValidPointer(X,2); 421 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 422 if (!isshell) PetscFunctionReturn(0); 423 *X = shell->Xglobal; 424 PetscFunctionReturn(0); 425 } 426 427 /*@C 428 DMShellSetCreateGlobalVector - sets the routine to create a global vector associated with the shell DM 429 430 Logically Collective 431 432 Input Arguments: 433 + dm - the shell DM 434 - func - the creation routine 435 436 Level: advanced 437 438 .seealso: DMShellSetGlobalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext() 439 @*/ 440 PetscErrorCode DMShellSetCreateGlobalVector(DM dm,PetscErrorCode (*func)(DM,Vec*)) 441 { 442 PetscFunctionBegin; 443 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 444 dm->ops->createglobalvector = func; 445 PetscFunctionReturn(0); 446 } 447 448 /*@ 449 DMShellSetLocalVector - sets a template local vector associated with the DMShell 450 451 Logically Collective on dm 452 453 Input Arguments: 454 + dm - shell DM 455 - X - template vector 456 457 Level: advanced 458 459 .seealso: DMCreateLocalVector(), DMShellSetMatrix(), DMShellSetCreateLocalVector() 460 @*/ 461 PetscErrorCode DMShellSetLocalVector(DM dm,Vec X) 462 { 463 DM_Shell *shell = (DM_Shell*)dm->data; 464 PetscErrorCode ierr; 465 PetscBool isshell; 466 DM vdm; 467 468 PetscFunctionBegin; 469 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 470 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 471 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 472 if (!isshell) PetscFunctionReturn(0); 473 ierr = VecGetDM(X,&vdm);CHKERRQ(ierr); 474 /* 475 if the vector proposed as the new base global vector for the DM is a DM vector associated 476 with the same DM then the current base global vector for the DM is ok and if we replace it with the new one 477 we get a circular dependency that prevents the DM from being destroy when it should be. 478 This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided 479 DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries 480 to set its input vector (which is associated with the DM) as the base global vector. 481 Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien 482 for pointing out the problem. 483 */ 484 if (vdm == dm) PetscFunctionReturn(0); 485 ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 486 ierr = VecDestroy(&shell->Xlocal);CHKERRQ(ierr); 487 shell->Xlocal = X; 488 PetscFunctionReturn(0); 489 } 490 491 /*@C 492 DMShellSetCreateLocalVector - sets the routine to create a local vector associated with the shell DM 493 494 Logically Collective 495 496 Input Arguments: 497 + dm - the shell DM 498 - func - the creation routine 499 500 Level: advanced 501 502 .seealso: DMShellSetLocalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext() 503 @*/ 504 PetscErrorCode DMShellSetCreateLocalVector(DM dm,PetscErrorCode (*func)(DM,Vec*)) 505 { 506 PetscFunctionBegin; 507 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 508 dm->ops->createlocalvector = func; 509 PetscFunctionReturn(0); 510 } 511 512 /*@C 513 DMShellSetGlobalToLocal - Sets the routines used to perform a global to local scatter 514 515 Logically Collective on dm 516 517 Input Arguments 518 + dm - the shell DM 519 . begin - the routine that begins the global to local scatter 520 - end - the routine that ends the global to local scatter 521 522 Notes: 523 If these functions are not provided but DMShellSetGlobalToLocalVecScatter() is called then 524 DMGlobalToLocalBeginDefaultShell()/DMGlobalToLocalEndDefaultShell() are used to to perform the transfers 525 526 Level: advanced 527 528 .seealso: DMShellSetLocalToGlobal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell() 529 @*/ 530 PetscErrorCode DMShellSetGlobalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) 531 { 532 PetscFunctionBegin; 533 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 534 dm->ops->globaltolocalbegin = begin; 535 dm->ops->globaltolocalend = end; 536 PetscFunctionReturn(0); 537 } 538 539 /*@C 540 DMShellSetLocalToGlobal - Sets the routines used to perform a local to global scatter 541 542 Logically Collective on dm 543 544 Input Arguments 545 + dm - the shell DM 546 . begin - the routine that begins the local to global scatter 547 - end - the routine that ends the local to global scatter 548 549 Notes: 550 If these functions are not provided but DMShellSetLocalToGlobalVecScatter() is called then 551 DMLocalToGlobalBeginDefaultShell()/DMLocalToGlobalEndDefaultShell() are used to to perform the transfers 552 553 Level: advanced 554 555 .seealso: DMShellSetGlobalToLocal() 556 @*/ 557 PetscErrorCode DMShellSetLocalToGlobal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) 558 { 559 PetscFunctionBegin; 560 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 561 dm->ops->localtoglobalbegin = begin; 562 dm->ops->localtoglobalend = end; 563 PetscFunctionReturn(0); 564 } 565 566 /*@C 567 DMShellSetLocalToLocal - Sets the routines used to perform a local to local scatter 568 569 Logically Collective on dm 570 571 Input Arguments 572 + dm - the shell DM 573 . begin - the routine that begins the local to local scatter 574 - end - the routine that ends the local to local scatter 575 576 Notes: 577 If these functions are not provided but DMShellSetLocalToLocalVecScatter() is called then 578 DMLocalToLocalBeginDefaultShell()/DMLocalToLocalEndDefaultShell() are used to to perform the transfers 579 580 Level: advanced 581 582 .seealso: DMShellSetGlobalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell() 583 @*/ 584 PetscErrorCode DMShellSetLocalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) 585 { 586 PetscFunctionBegin; 587 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 588 dm->ops->localtolocalbegin = begin; 589 dm->ops->localtolocalend = end; 590 PetscFunctionReturn(0); 591 } 592 593 /*@ 594 DMShellSetGlobalToLocalVecScatter - Sets a VecScatter context for global to local communication 595 596 Logically Collective on dm 597 598 Input Arguments 599 + dm - the shell DM 600 - gtol - the global to local VecScatter context 601 602 Level: advanced 603 604 .seealso: DMShellSetGlobalToLocal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell() 605 @*/ 606 PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol) 607 { 608 DM_Shell *shell = (DM_Shell*)dm->data; 609 PetscErrorCode ierr; 610 611 PetscFunctionBegin; 612 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 613 PetscValidHeaderSpecific(gtol,PETSCSF_CLASSID,2); 614 ierr = PetscObjectReference((PetscObject)gtol);CHKERRQ(ierr); 615 ierr = VecScatterDestroy(&shell->gtol);CHKERRQ(ierr); 616 shell->gtol = gtol; 617 PetscFunctionReturn(0); 618 } 619 620 /*@ 621 DMShellSetLocalToGlobalVecScatter - Sets a VecScatter context for local to global communication 622 623 Logically Collective on dm 624 625 Input Arguments 626 + dm - the shell DM 627 - ltog - the local to global VecScatter context 628 629 Level: advanced 630 631 .seealso: DMShellSetLocalToGlobal(), DMLocalToGlobalBeginDefaultShell(), DMLocalToGlobalEndDefaultShell() 632 @*/ 633 PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog) 634 { 635 DM_Shell *shell = (DM_Shell*)dm->data; 636 PetscErrorCode ierr; 637 638 PetscFunctionBegin; 639 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 640 PetscValidHeaderSpecific(ltog,PETSCSF_CLASSID,2); 641 ierr = PetscObjectReference((PetscObject)ltog);CHKERRQ(ierr); 642 ierr = VecScatterDestroy(&shell->ltog);CHKERRQ(ierr); 643 shell->ltog = ltog; 644 PetscFunctionReturn(0); 645 } 646 647 /*@ 648 DMShellSetLocalToLocalVecScatter - Sets a VecScatter context for local to local communication 649 650 Logically Collective on dm 651 652 Input Arguments 653 + dm - the shell DM 654 - ltol - the local to local VecScatter context 655 656 Level: advanced 657 658 .seealso: DMShellSetLocalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell() 659 @*/ 660 PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol) 661 { 662 DM_Shell *shell = (DM_Shell*)dm->data; 663 PetscErrorCode ierr; 664 665 PetscFunctionBegin; 666 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 667 PetscValidHeaderSpecific(ltol,PETSCSF_CLASSID,2); 668 ierr = PetscObjectReference((PetscObject)ltol);CHKERRQ(ierr); 669 ierr = VecScatterDestroy(&shell->ltol);CHKERRQ(ierr); 670 shell->ltol = ltol; 671 PetscFunctionReturn(0); 672 } 673 674 /*@C 675 DMShellSetCoarsen - Set the routine used to coarsen the shell DM 676 677 Logically Collective on dm 678 679 Input Arguments 680 + dm - the shell DM 681 - coarsen - the routine that coarsens the DM 682 683 Level: advanced 684 685 .seealso: DMShellSetRefine(), DMCoarsen(), DMShellGetCoarsen(), DMShellSetContext(), DMShellGetContext() 686 @*/ 687 PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*)) 688 { 689 PetscErrorCode ierr; 690 PetscBool isshell; 691 692 PetscFunctionBegin; 693 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 694 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 695 if (!isshell) PetscFunctionReturn(0); 696 dm->ops->coarsen = coarsen; 697 PetscFunctionReturn(0); 698 } 699 700 /*@C 701 DMShellGetCoarsen - Get the routine used to coarsen the shell DM 702 703 Logically Collective on dm 704 705 Input Argument: 706 . dm - the shell DM 707 708 Output Argument: 709 . coarsen - the routine that coarsens the DM 710 711 Level: advanced 712 713 .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine() 714 @*/ 715 PetscErrorCode DMShellGetCoarsen(DM dm, PetscErrorCode (**coarsen)(DM,MPI_Comm,DM*)) 716 { 717 PetscErrorCode ierr; 718 PetscBool isshell; 719 720 PetscFunctionBegin; 721 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 722 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 723 if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 724 *coarsen = dm->ops->coarsen; 725 PetscFunctionReturn(0); 726 } 727 728 /*@C 729 DMShellSetRefine - Set the routine used to refine the shell DM 730 731 Logically Collective on dm 732 733 Input Arguments 734 + dm - the shell DM 735 - refine - the routine that refines the DM 736 737 Level: advanced 738 739 .seealso: DMShellSetCoarsen(), DMRefine(), DMShellGetRefine(), DMShellSetContext(), DMShellGetContext() 740 @*/ 741 PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM,MPI_Comm,DM*)) 742 { 743 PetscErrorCode ierr; 744 PetscBool isshell; 745 746 PetscFunctionBegin; 747 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 748 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 749 if (!isshell) PetscFunctionReturn(0); 750 dm->ops->refine = refine; 751 PetscFunctionReturn(0); 752 } 753 754 /*@C 755 DMShellGetRefine - Get the routine used to refine the shell DM 756 757 Logically Collective on dm 758 759 Input Argument: 760 . dm - the shell DM 761 762 Output Argument: 763 . refine - the routine that refines the DM 764 765 Level: advanced 766 767 .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine() 768 @*/ 769 PetscErrorCode DMShellGetRefine(DM dm, PetscErrorCode (**refine)(DM,MPI_Comm,DM*)) 770 { 771 PetscErrorCode ierr; 772 PetscBool isshell; 773 774 PetscFunctionBegin; 775 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 776 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 777 if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 778 *refine = dm->ops->refine; 779 PetscFunctionReturn(0); 780 } 781 782 /*@C 783 DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator 784 785 Logically Collective on dm 786 787 Input Arguments 788 + dm - the shell DM 789 - interp - the routine to create the interpolation 790 791 Level: advanced 792 793 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateInterpolation(), DMShellSetCreateRestriction(), DMShellSetContext(), DMShellGetContext() 794 @*/ 795 PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM,DM,Mat*,Vec*)) 796 { 797 PetscErrorCode ierr; 798 PetscBool isshell; 799 800 PetscFunctionBegin; 801 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 802 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 803 if (!isshell) PetscFunctionReturn(0); 804 dm->ops->createinterpolation = interp; 805 PetscFunctionReturn(0); 806 } 807 808 /*@C 809 DMShellGetCreateInterpolation - Get the routine used to create the interpolation operator 810 811 Logically Collective on dm 812 813 Input Argument: 814 + dm - the shell DM 815 816 Output Argument: 817 - interp - the routine to create the interpolation 818 819 Level: advanced 820 821 .seealso: DMShellGetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext() 822 @*/ 823 PetscErrorCode DMShellGetCreateInterpolation(DM dm, PetscErrorCode (**interp)(DM,DM,Mat*,Vec*)) 824 { 825 PetscErrorCode ierr; 826 PetscBool isshell; 827 828 PetscFunctionBegin; 829 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 830 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 831 if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 832 *interp = dm->ops->createinterpolation; 833 PetscFunctionReturn(0); 834 } 835 836 /*@C 837 DMShellSetCreateRestriction - Set the routine used to create the restriction operator 838 839 Logically Collective on dm 840 841 Input Arguments 842 + dm - the shell DM 843 - striction- the routine to create the restriction 844 845 Level: advanced 846 847 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext() 848 @*/ 849 PetscErrorCode DMShellSetCreateRestriction(DM dm, PetscErrorCode (*restriction)(DM,DM,Mat*)) 850 { 851 PetscErrorCode ierr; 852 PetscBool isshell; 853 854 PetscFunctionBegin; 855 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 856 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 857 if (!isshell) PetscFunctionReturn(0); 858 dm->ops->createrestriction = restriction; 859 PetscFunctionReturn(0); 860 } 861 862 /*@C 863 DMShellGetCreateRestriction - Get the routine used to create the restriction operator 864 865 Logically Collective on dm 866 867 Input Argument: 868 + dm - the shell DM 869 870 Output Argument: 871 - restriction - the routine to create the restriction 872 873 Level: advanced 874 875 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetContext(), DMShellGetContext() 876 @*/ 877 PetscErrorCode DMShellGetCreateRestriction(DM dm, PetscErrorCode (**restriction)(DM,DM,Mat*)) 878 { 879 PetscErrorCode ierr; 880 PetscBool isshell; 881 882 PetscFunctionBegin; 883 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 884 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 885 if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 886 *restriction = dm->ops->createrestriction; 887 PetscFunctionReturn(0); 888 } 889 890 /*@C 891 DMShellSetCreateInjection - Set the routine used to create the injection operator 892 893 Logically Collective on dm 894 895 Input Arguments 896 + dm - the shell DM 897 - inject - the routine to create the injection 898 899 Level: advanced 900 901 .seealso: DMShellSetCreateInterpolation(), DMCreateInjection(), DMShellGetCreateInjection(), DMShellSetContext(), DMShellGetContext() 902 @*/ 903 PetscErrorCode DMShellSetCreateInjection(DM dm, PetscErrorCode (*inject)(DM,DM,Mat*)) 904 { 905 PetscErrorCode ierr; 906 PetscBool isshell; 907 908 PetscFunctionBegin; 909 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 910 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 911 if (!isshell) PetscFunctionReturn(0); 912 dm->ops->createinjection = inject; 913 PetscFunctionReturn(0); 914 } 915 916 /*@C 917 DMShellGetCreateInjection - Get the routine used to create the injection operator 918 919 Logically Collective on dm 920 921 Input Argument: 922 + dm - the shell DM 923 924 Output Argument: 925 - inject - the routine to create the injection 926 927 Level: advanced 928 929 .seealso: DMShellGetCreateInterpolation(), DMCreateInjection(), DMShellSetContext(), DMShellGetContext() 930 @*/ 931 PetscErrorCode DMShellGetCreateInjection(DM dm, PetscErrorCode (**inject)(DM,DM,Mat*)) 932 { 933 PetscErrorCode ierr; 934 PetscBool isshell; 935 936 PetscFunctionBegin; 937 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 938 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 939 if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 940 *inject = dm->ops->createinjection; 941 PetscFunctionReturn(0); 942 } 943 944 /*@C 945 DMShellSetCreateFieldDecomposition - Set the routine used to create a decomposition of fields for the shell DM 946 947 Logically Collective on dm 948 949 Input Arguments 950 + dm - the shell DM 951 - decomp - the routine to create the decomposition 952 953 Level: advanced 954 955 .seealso: DMCreateFieldDecomposition(), DMShellSetContext(), DMShellGetContext() 956 @*/ 957 PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,DM**)) 958 { 959 PetscErrorCode ierr; 960 PetscBool isshell; 961 962 PetscFunctionBegin; 963 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 964 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 965 if (!isshell) PetscFunctionReturn(0); 966 dm->ops->createfielddecomposition = decomp; 967 PetscFunctionReturn(0); 968 } 969 970 /*@C 971 DMShellSetCreateDomainDecomposition - Set the routine used to create a domain decomposition for the shell DM 972 973 Logically Collective on dm 974 975 Input Arguments 976 + dm - the shell DM 977 - decomp - the routine to create the decomposition 978 979 Level: advanced 980 981 .seealso: DMCreateDomainDecomposition(), DMShellSetContext(), DMShellGetContext() 982 @*/ 983 PetscErrorCode DMShellSetCreateDomainDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,IS**,DM**)) 984 { 985 PetscErrorCode ierr; 986 PetscBool isshell; 987 988 PetscFunctionBegin; 989 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 990 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 991 if (!isshell) PetscFunctionReturn(0); 992 dm->ops->createdomaindecomposition = decomp; 993 PetscFunctionReturn(0); 994 } 995 996 /*@C 997 DMShellSetCreateDomainDecompositionScatters - Set the routine used to create the scatter contexts for domain decomposition with a shell DM 998 999 Logically Collective on dm 1000 1001 Input Arguments 1002 + dm - the shell DM 1003 - scatter - the routine to create the scatters 1004 1005 Level: advanced 1006 1007 .seealso: DMCreateDomainDecompositionScatters(), DMShellSetContext(), DMShellGetContext() 1008 @*/ 1009 PetscErrorCode DMShellSetCreateDomainDecompositionScatters(DM dm, PetscErrorCode (*scatter)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**)) 1010 { 1011 PetscErrorCode ierr; 1012 PetscBool isshell; 1013 1014 PetscFunctionBegin; 1015 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1016 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 1017 if (!isshell) PetscFunctionReturn(0); 1018 dm->ops->createddscatters = scatter; 1019 PetscFunctionReturn(0); 1020 } 1021 1022 /*@C 1023 DMShellSetCreateSubDM - Set the routine used to create a sub DM from the shell DM 1024 1025 Logically Collective on dm 1026 1027 Input Arguments 1028 + dm - the shell DM 1029 - subdm - the routine to create the decomposition 1030 1031 Level: advanced 1032 1033 .seealso: DMCreateSubDM(), DMShellGetCreateSubDM(), DMShellSetContext(), DMShellGetContext() 1034 @*/ 1035 PetscErrorCode DMShellSetCreateSubDM(DM dm, PetscErrorCode (*subdm)(DM,PetscInt,const PetscInt[],IS*,DM*)) 1036 { 1037 PetscErrorCode ierr; 1038 PetscBool isshell; 1039 1040 PetscFunctionBegin; 1041 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1042 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 1043 if (!isshell) PetscFunctionReturn(0); 1044 dm->ops->createsubdm = subdm; 1045 PetscFunctionReturn(0); 1046 } 1047 1048 /*@C 1049 DMShellGetCreateSubDM - Get the routine used to create a sub DM from the shell DM 1050 1051 Logically Collective on dm 1052 1053 Input Argument: 1054 + dm - the shell DM 1055 1056 Output Argument: 1057 - subdm - the routine to create the decomposition 1058 1059 Level: advanced 1060 1061 .seealso: DMCreateSubDM(), DMShellSetCreateSubDM(), DMShellSetContext(), DMShellGetContext() 1062 @*/ 1063 PetscErrorCode DMShellGetCreateSubDM(DM dm, PetscErrorCode (**subdm)(DM,PetscInt,const PetscInt[],IS*,DM*)) 1064 { 1065 PetscErrorCode ierr; 1066 PetscBool isshell; 1067 1068 PetscFunctionBegin; 1069 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1070 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 1071 if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 1072 *subdm = dm->ops->createsubdm; 1073 PetscFunctionReturn(0); 1074 } 1075 1076 static PetscErrorCode DMDestroy_Shell(DM dm) 1077 { 1078 PetscErrorCode ierr; 1079 DM_Shell *shell = (DM_Shell*)dm->data; 1080 1081 PetscFunctionBegin; 1082 ierr = MatDestroy(&shell->A);CHKERRQ(ierr); 1083 ierr = VecDestroy(&shell->Xglobal);CHKERRQ(ierr); 1084 ierr = VecDestroy(&shell->Xlocal);CHKERRQ(ierr); 1085 ierr = VecScatterDestroy(&shell->gtol);CHKERRQ(ierr); 1086 ierr = VecScatterDestroy(&shell->ltog);CHKERRQ(ierr); 1087 ierr = VecScatterDestroy(&shell->ltol);CHKERRQ(ierr); 1088 /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 1089 ierr = PetscFree(shell);CHKERRQ(ierr); 1090 PetscFunctionReturn(0); 1091 } 1092 1093 static PetscErrorCode DMView_Shell(DM dm,PetscViewer v) 1094 { 1095 PetscErrorCode ierr; 1096 DM_Shell *shell = (DM_Shell*)dm->data; 1097 1098 PetscFunctionBegin; 1099 ierr = VecView(shell->Xglobal,v);CHKERRQ(ierr); 1100 PetscFunctionReturn(0); 1101 } 1102 1103 static PetscErrorCode DMLoad_Shell(DM dm,PetscViewer v) 1104 { 1105 PetscErrorCode ierr; 1106 DM_Shell *shell = (DM_Shell*)dm->data; 1107 1108 PetscFunctionBegin; 1109 ierr = VecCreate(PetscObjectComm((PetscObject)dm),&shell->Xglobal);CHKERRQ(ierr); 1110 ierr = VecLoad(shell->Xglobal,v);CHKERRQ(ierr); 1111 PetscFunctionReturn(0); 1112 } 1113 1114 PetscErrorCode DMCreateSubDM_Shell(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 1115 { 1116 PetscErrorCode ierr; 1117 1118 PetscFunctionBegin; 1119 if (subdm) {ierr = DMShellCreate(PetscObjectComm((PetscObject) dm), subdm);CHKERRQ(ierr);} 1120 ierr = DMCreateSectionSubDM(dm, numFields, fields, is, subdm);CHKERRQ(ierr); 1121 PetscFunctionReturn(0); 1122 } 1123 1124 PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm) 1125 { 1126 PetscErrorCode ierr; 1127 DM_Shell *shell; 1128 1129 PetscFunctionBegin; 1130 ierr = PetscNewLog(dm,&shell);CHKERRQ(ierr); 1131 dm->data = shell; 1132 1133 dm->ops->destroy = DMDestroy_Shell; 1134 dm->ops->createglobalvector = DMCreateGlobalVector_Shell; 1135 dm->ops->createlocalvector = DMCreateLocalVector_Shell; 1136 dm->ops->creatematrix = DMCreateMatrix_Shell; 1137 dm->ops->view = DMView_Shell; 1138 dm->ops->load = DMLoad_Shell; 1139 dm->ops->globaltolocalbegin = DMGlobalToLocalBeginDefaultShell; 1140 dm->ops->globaltolocalend = DMGlobalToLocalEndDefaultShell; 1141 dm->ops->localtoglobalbegin = DMLocalToGlobalBeginDefaultShell; 1142 dm->ops->localtoglobalend = DMLocalToGlobalEndDefaultShell; 1143 dm->ops->localtolocalbegin = DMLocalToLocalBeginDefaultShell; 1144 dm->ops->localtolocalend = DMLocalToLocalEndDefaultShell; 1145 dm->ops->createsubdm = DMCreateSubDM_Shell; 1146 ierr = DMSetMatType(dm,MATDENSE);CHKERRQ(ierr); 1147 PetscFunctionReturn(0); 1148 } 1149 1150 /*@ 1151 DMShellCreate - Creates a shell DM object, used to manage user-defined problem data 1152 1153 Collective 1154 1155 Input Parameter: 1156 . comm - the processors that will share the global vector 1157 1158 Output Parameters: 1159 . shell - the shell DM 1160 1161 Level: advanced 1162 1163 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateLocalVector(), DMShellSetContext(), DMShellGetContext() 1164 @*/ 1165 PetscErrorCode DMShellCreate(MPI_Comm comm,DM *dm) 1166 { 1167 PetscErrorCode ierr; 1168 1169 PetscFunctionBegin; 1170 PetscValidPointer(dm,2); 1171 ierr = DMCreate(comm,dm);CHKERRQ(ierr); 1172 ierr = DMSetType(*dm,DMSHELL);CHKERRQ(ierr); 1173 ierr = DMSetUp(*dm);CHKERRQ(ierr); 1174 PetscFunctionReturn(0); 1175 } 1176