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,3); 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 *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 /*@C 400 DMShellSetCreateGlobalVector - sets the routine to create a global vector associated with the shell DM 401 402 Logically Collective 403 404 Input Arguments: 405 + dm - the shell DM 406 - func - the creation routine 407 408 Level: advanced 409 410 .seealso: DMShellSetGlobalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext() 411 @*/ 412 PetscErrorCode DMShellSetCreateGlobalVector(DM dm,PetscErrorCode (*func)(DM,Vec*)) 413 { 414 PetscFunctionBegin; 415 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 416 dm->ops->createglobalvector = func; 417 PetscFunctionReturn(0); 418 } 419 420 /*@ 421 DMShellSetLocalVector - sets a template local vector associated with the DMShell 422 423 Logically Collective on dm 424 425 Input Arguments: 426 + dm - shell DM 427 - X - template vector 428 429 Level: advanced 430 431 .seealso: DMCreateLocalVector(), DMShellSetMatrix(), DMShellSetCreateLocalVector() 432 @*/ 433 PetscErrorCode DMShellSetLocalVector(DM dm,Vec X) 434 { 435 DM_Shell *shell = (DM_Shell*)dm->data; 436 PetscErrorCode ierr; 437 PetscBool isshell; 438 DM vdm; 439 440 PetscFunctionBegin; 441 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 442 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 443 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 444 if (!isshell) PetscFunctionReturn(0); 445 ierr = VecGetDM(X,&vdm);CHKERRQ(ierr); 446 /* 447 if the vector proposed as the new base global vector for the DM is a DM vector associated 448 with the same DM then the current base global vector for the DM is ok and if we replace it with the new one 449 we get a circular dependency that prevents the DM from being destroy when it should be. 450 This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided 451 DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries 452 to set its input vector (which is associated with the DM) as the base global vector. 453 Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien 454 for pointing out the problem. 455 */ 456 if (vdm == dm) PetscFunctionReturn(0); 457 ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 458 ierr = VecDestroy(&shell->Xlocal);CHKERRQ(ierr); 459 shell->Xlocal = X; 460 PetscFunctionReturn(0); 461 } 462 463 /*@C 464 DMShellSetCreateLocalVector - sets the routine to create a local vector associated with the shell DM 465 466 Logically Collective 467 468 Input Arguments: 469 + dm - the shell DM 470 - func - the creation routine 471 472 Level: advanced 473 474 .seealso: DMShellSetLocalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext() 475 @*/ 476 PetscErrorCode DMShellSetCreateLocalVector(DM dm,PetscErrorCode (*func)(DM,Vec*)) 477 { 478 PetscFunctionBegin; 479 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 480 dm->ops->createlocalvector = func; 481 PetscFunctionReturn(0); 482 } 483 484 /*@C 485 DMShellSetGlobalToLocal - Sets the routines used to perform a global to local scatter 486 487 Logically Collective on dm 488 489 Input Arguments 490 + dm - the shell DM 491 . begin - the routine that begins the global to local scatter 492 - end - the routine that ends the global to local scatter 493 494 Notes: 495 If these functions are not provided but DMShellSetGlobalToLocalVecScatter() is called then 496 DMGlobalToLocalBeginDefaultShell()/DMGlobalToLocalEndDefaultShell() are used to to perform the transfers 497 498 Level: advanced 499 500 .seealso: DMShellSetLocalToGlobal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell() 501 @*/ 502 PetscErrorCode DMShellSetGlobalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) 503 { 504 PetscFunctionBegin; 505 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 506 dm->ops->globaltolocalbegin = begin; 507 dm->ops->globaltolocalend = end; 508 PetscFunctionReturn(0); 509 } 510 511 /*@C 512 DMShellSetLocalToGlobal - Sets the routines used to perform a local to global scatter 513 514 Logically Collective on dm 515 516 Input Arguments 517 + dm - the shell DM 518 . begin - the routine that begins the local to global scatter 519 - end - the routine that ends the local to global scatter 520 521 Notes: 522 If these functions are not provided but DMShellSetLocalToGlobalVecScatter() is called then 523 DMLocalToGlobalBeginDefaultShell()/DMLocalToGlobalEndDefaultShell() are used to to perform the transfers 524 525 Level: advanced 526 527 .seealso: DMShellSetGlobalToLocal() 528 @*/ 529 PetscErrorCode DMShellSetLocalToGlobal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) 530 { 531 PetscFunctionBegin; 532 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 533 dm->ops->localtoglobalbegin = begin; 534 dm->ops->localtoglobalend = end; 535 PetscFunctionReturn(0); 536 } 537 538 /*@C 539 DMShellSetLocalToLocal - Sets the routines used to perform a local to local scatter 540 541 Logically Collective on dm 542 543 Input Arguments 544 + dm - the shell DM 545 . begin - the routine that begins the local to local scatter 546 - end - the routine that ends the local to local scatter 547 548 Notes: 549 If these functions are not provided but DMShellSetLocalToLocalVecScatter() is called then 550 DMLocalToLocalBeginDefaultShell()/DMLocalToLocalEndDefaultShell() are used to to perform the transfers 551 552 Level: advanced 553 554 .seealso: DMShellSetGlobalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell() 555 @*/ 556 PetscErrorCode DMShellSetLocalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) 557 { 558 PetscFunctionBegin; 559 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 560 dm->ops->localtolocalbegin = begin; 561 dm->ops->localtolocalend = end; 562 PetscFunctionReturn(0); 563 } 564 565 /*@ 566 DMShellSetGlobalToLocalVecScatter - Sets a VecScatter context for global to local communication 567 568 Logically Collective on dm 569 570 Input Arguments 571 + dm - the shell DM 572 - gtol - the global to local VecScatter context 573 574 Level: advanced 575 576 .seealso: DMShellSetGlobalToLocal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell() 577 @*/ 578 PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol) 579 { 580 DM_Shell *shell = (DM_Shell*)dm->data; 581 PetscErrorCode ierr; 582 583 PetscFunctionBegin; 584 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 585 PetscValidHeaderSpecific(gtol,PETSCSF_CLASSID,2); 586 ierr = PetscObjectReference((PetscObject)gtol);CHKERRQ(ierr); 587 ierr = VecScatterDestroy(&shell->gtol);CHKERRQ(ierr); 588 shell->gtol = gtol; 589 PetscFunctionReturn(0); 590 } 591 592 /*@ 593 DMShellSetLocalToGlobalVecScatter - Sets a VecScatter context for local to global communication 594 595 Logically Collective on dm 596 597 Input Arguments 598 + dm - the shell DM 599 - ltog - the local to global VecScatter context 600 601 Level: advanced 602 603 .seealso: DMShellSetLocalToGlobal(), DMLocalToGlobalBeginDefaultShell(), DMLocalToGlobalEndDefaultShell() 604 @*/ 605 PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog) 606 { 607 DM_Shell *shell = (DM_Shell*)dm->data; 608 PetscErrorCode ierr; 609 610 PetscFunctionBegin; 611 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 612 PetscValidHeaderSpecific(ltog,PETSCSF_CLASSID,2); 613 ierr = PetscObjectReference((PetscObject)ltog);CHKERRQ(ierr); 614 ierr = VecScatterDestroy(&shell->ltog);CHKERRQ(ierr); 615 shell->ltog = ltog; 616 PetscFunctionReturn(0); 617 } 618 619 /*@ 620 DMShellSetLocalToLocalVecScatter - Sets a VecScatter context for local to local communication 621 622 Logically Collective on dm 623 624 Input Arguments 625 + dm - the shell DM 626 - ltol - the local to local VecScatter context 627 628 Level: advanced 629 630 .seealso: DMShellSetLocalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell() 631 @*/ 632 PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol) 633 { 634 DM_Shell *shell = (DM_Shell*)dm->data; 635 PetscErrorCode ierr; 636 637 PetscFunctionBegin; 638 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 639 PetscValidHeaderSpecific(ltol,PETSCSF_CLASSID,2); 640 ierr = PetscObjectReference((PetscObject)ltol);CHKERRQ(ierr); 641 ierr = VecScatterDestroy(&shell->ltol);CHKERRQ(ierr); 642 shell->ltol = ltol; 643 PetscFunctionReturn(0); 644 } 645 646 /*@C 647 DMShellSetCoarsen - Set the routine used to coarsen the shell DM 648 649 Logically Collective on dm 650 651 Input Arguments 652 + dm - the shell DM 653 - coarsen - the routine that coarsens the DM 654 655 Level: advanced 656 657 .seealso: DMShellSetRefine(), DMCoarsen(), DMShellGetCoarsen(), DMShellSetContext(), DMShellGetContext() 658 @*/ 659 PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*)) 660 { 661 PetscErrorCode ierr; 662 PetscBool isshell; 663 664 PetscFunctionBegin; 665 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 666 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 667 if (!isshell) PetscFunctionReturn(0); 668 dm->ops->coarsen = coarsen; 669 PetscFunctionReturn(0); 670 } 671 672 /*@C 673 DMShellGetCoarsen - Get the routine used to coarsen the shell DM 674 675 Logically Collective on dm 676 677 Input Argument: 678 . dm - the shell DM 679 680 Output Argument: 681 . coarsen - the routine that coarsens the DM 682 683 Level: advanced 684 685 .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine() 686 @*/ 687 PetscErrorCode DMShellGetCoarsen(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) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 696 *coarsen = dm->ops->coarsen; 697 PetscFunctionReturn(0); 698 } 699 700 /*@C 701 DMShellSetRefine - Set the routine used to refine the shell DM 702 703 Logically Collective on dm 704 705 Input Arguments 706 + dm - the shell DM 707 - refine - the routine that refines the DM 708 709 Level: advanced 710 711 .seealso: DMShellSetCoarsen(), DMRefine(), DMShellGetRefine(), DMShellSetContext(), DMShellGetContext() 712 @*/ 713 PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM,MPI_Comm,DM*)) 714 { 715 PetscErrorCode ierr; 716 PetscBool isshell; 717 718 PetscFunctionBegin; 719 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 720 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 721 if (!isshell) PetscFunctionReturn(0); 722 dm->ops->refine = refine; 723 PetscFunctionReturn(0); 724 } 725 726 /*@C 727 DMShellGetRefine - Get the routine used to refine the shell DM 728 729 Logically Collective on dm 730 731 Input Argument: 732 . dm - the shell DM 733 734 Output Argument: 735 . refine - the routine that refines the DM 736 737 Level: advanced 738 739 .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine() 740 @*/ 741 PetscErrorCode DMShellGetRefine(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) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 750 *refine = dm->ops->refine; 751 PetscFunctionReturn(0); 752 } 753 754 /*@C 755 DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator 756 757 Logically Collective on dm 758 759 Input Arguments 760 + dm - the shell DM 761 - interp - the routine to create the interpolation 762 763 Level: advanced 764 765 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateInterpolation(), DMShellSetCreateRestriction(), DMShellSetContext(), DMShellGetContext() 766 @*/ 767 PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM,DM,Mat*,Vec*)) 768 { 769 PetscErrorCode ierr; 770 PetscBool isshell; 771 772 PetscFunctionBegin; 773 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 774 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 775 if (!isshell) PetscFunctionReturn(0); 776 dm->ops->createinterpolation = interp; 777 PetscFunctionReturn(0); 778 } 779 780 /*@C 781 DMShellGetCreateInterpolation - Get the routine used to create the interpolation operator 782 783 Logically Collective on dm 784 785 Input Argument: 786 + dm - the shell DM 787 788 Output Argument: 789 - interp - the routine to create the interpolation 790 791 Level: advanced 792 793 .seealso: DMShellGetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext() 794 @*/ 795 PetscErrorCode DMShellGetCreateInterpolation(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) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 804 *interp = dm->ops->createinterpolation; 805 PetscFunctionReturn(0); 806 } 807 808 /*@C 809 DMShellSetCreateRestriction - Set the routine used to create the restriction operator 810 811 Logically Collective on dm 812 813 Input Arguments 814 + dm - the shell DM 815 - striction- the routine to create the restriction 816 817 Level: advanced 818 819 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext() 820 @*/ 821 PetscErrorCode DMShellSetCreateRestriction(DM dm, PetscErrorCode (*restriction)(DM,DM,Mat*)) 822 { 823 PetscErrorCode ierr; 824 PetscBool isshell; 825 826 PetscFunctionBegin; 827 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 828 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 829 if (!isshell) PetscFunctionReturn(0); 830 dm->ops->createrestriction = restriction; 831 PetscFunctionReturn(0); 832 } 833 834 /*@C 835 DMShellGetCreateRestriction - Get the routine used to create the restriction operator 836 837 Logically Collective on dm 838 839 Input Argument: 840 + dm - the shell DM 841 842 Output Argument: 843 - restriction - the routine to create the restriction 844 845 Level: advanced 846 847 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetContext(), DMShellGetContext() 848 @*/ 849 PetscErrorCode DMShellGetCreateRestriction(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) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 858 *restriction = dm->ops->createrestriction; 859 PetscFunctionReturn(0); 860 } 861 862 /*@C 863 DMShellSetCreateInjection - Set the routine used to create the injection operator 864 865 Logically Collective on dm 866 867 Input Arguments 868 + dm - the shell DM 869 - inject - the routine to create the injection 870 871 Level: advanced 872 873 .seealso: DMShellSetCreateInterpolation(), DMCreateInjection(), DMShellGetCreateInjection(), DMShellSetContext(), DMShellGetContext() 874 @*/ 875 PetscErrorCode DMShellSetCreateInjection(DM dm, PetscErrorCode (*inject)(DM,DM,Mat*)) 876 { 877 PetscErrorCode ierr; 878 PetscBool isshell; 879 880 PetscFunctionBegin; 881 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 882 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 883 if (!isshell) PetscFunctionReturn(0); 884 dm->ops->createinjection = inject; 885 PetscFunctionReturn(0); 886 } 887 888 /*@C 889 DMShellGetCreateInjection - Get the routine used to create the injection operator 890 891 Logically Collective on dm 892 893 Input Argument: 894 + dm - the shell DM 895 896 Output Argument: 897 - inject - the routine to create the injection 898 899 Level: advanced 900 901 .seealso: DMShellGetCreateInterpolation(), DMCreateInjection(), DMShellSetContext(), DMShellGetContext() 902 @*/ 903 PetscErrorCode DMShellGetCreateInjection(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) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 912 *inject = dm->ops->createinjection; 913 PetscFunctionReturn(0); 914 } 915 916 /*@C 917 DMShellSetCreateFieldDecomposition - Set the routine used to create a decomposition of fields for the shell DM 918 919 Logically Collective on dm 920 921 Input Arguments 922 + dm - the shell DM 923 - decomp - the routine to create the decomposition 924 925 Level: advanced 926 927 .seealso: DMCreateFieldDecomposition(), DMShellSetContext(), DMShellGetContext() 928 @*/ 929 PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,DM**)) 930 { 931 PetscErrorCode ierr; 932 PetscBool isshell; 933 934 PetscFunctionBegin; 935 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 936 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 937 if (!isshell) PetscFunctionReturn(0); 938 dm->ops->createfielddecomposition = decomp; 939 PetscFunctionReturn(0); 940 } 941 942 /*@C 943 DMShellSetCreateDomainDecomposition - Set the routine used to create a domain decomposition for the shell DM 944 945 Logically Collective on dm 946 947 Input Arguments 948 + dm - the shell DM 949 - decomp - the routine to create the decomposition 950 951 Level: advanced 952 953 .seealso: DMCreateDomainDecomposition(), DMShellSetContext(), DMShellGetContext() 954 @*/ 955 PetscErrorCode DMShellSetCreateDomainDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,IS**,DM**)) 956 { 957 PetscErrorCode ierr; 958 PetscBool isshell; 959 960 PetscFunctionBegin; 961 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 962 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 963 if (!isshell) PetscFunctionReturn(0); 964 dm->ops->createdomaindecomposition = decomp; 965 PetscFunctionReturn(0); 966 } 967 968 /*@C 969 DMShellSetCreateDomainDecompositionScatters - Set the routine used to create the scatter contexts for domain decomposition with a shell DM 970 971 Logically Collective on dm 972 973 Input Arguments 974 + dm - the shell DM 975 - scatter - the routine to create the scatters 976 977 Level: advanced 978 979 .seealso: DMCreateDomainDecompositionScatters(), DMShellSetContext(), DMShellGetContext() 980 @*/ 981 PetscErrorCode DMShellSetCreateDomainDecompositionScatters(DM dm, PetscErrorCode (*scatter)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**)) 982 { 983 PetscErrorCode ierr; 984 PetscBool isshell; 985 986 PetscFunctionBegin; 987 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 988 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 989 if (!isshell) PetscFunctionReturn(0); 990 dm->ops->createddscatters = scatter; 991 PetscFunctionReturn(0); 992 } 993 994 /*@C 995 DMShellSetCreateSubDM - Set the routine used to create a sub DM from the shell DM 996 997 Logically Collective on dm 998 999 Input Arguments 1000 + dm - the shell DM 1001 - subdm - the routine to create the decomposition 1002 1003 Level: advanced 1004 1005 .seealso: DMCreateSubDM(), DMShellGetCreateSubDM(), DMShellSetContext(), DMShellGetContext() 1006 @*/ 1007 PetscErrorCode DMShellSetCreateSubDM(DM dm, PetscErrorCode (*subdm)(DM,PetscInt,const PetscInt[],IS*,DM*)) 1008 { 1009 PetscErrorCode ierr; 1010 PetscBool isshell; 1011 1012 PetscFunctionBegin; 1013 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1014 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 1015 if (!isshell) PetscFunctionReturn(0); 1016 dm->ops->createsubdm = subdm; 1017 PetscFunctionReturn(0); 1018 } 1019 1020 /*@C 1021 DMShellGetCreateSubDM - Get the routine used to create a sub DM from the shell DM 1022 1023 Logically Collective on dm 1024 1025 Input Argument: 1026 + dm - the shell DM 1027 1028 Output Argument: 1029 - subdm - the routine to create the decomposition 1030 1031 Level: advanced 1032 1033 .seealso: DMCreateSubDM(), DMShellSetCreateSubDM(), DMShellSetContext(), DMShellGetContext() 1034 @*/ 1035 PetscErrorCode DMShellGetCreateSubDM(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) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 1044 *subdm = dm->ops->createsubdm; 1045 PetscFunctionReturn(0); 1046 } 1047 1048 static PetscErrorCode DMDestroy_Shell(DM dm) 1049 { 1050 PetscErrorCode ierr; 1051 DM_Shell *shell = (DM_Shell*)dm->data; 1052 1053 PetscFunctionBegin; 1054 ierr = MatDestroy(&shell->A);CHKERRQ(ierr); 1055 ierr = VecDestroy(&shell->Xglobal);CHKERRQ(ierr); 1056 ierr = VecDestroy(&shell->Xlocal);CHKERRQ(ierr); 1057 ierr = VecScatterDestroy(&shell->gtol);CHKERRQ(ierr); 1058 ierr = VecScatterDestroy(&shell->ltog);CHKERRQ(ierr); 1059 ierr = VecScatterDestroy(&shell->ltol);CHKERRQ(ierr); 1060 /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 1061 ierr = PetscFree(shell);CHKERRQ(ierr); 1062 PetscFunctionReturn(0); 1063 } 1064 1065 static PetscErrorCode DMView_Shell(DM dm,PetscViewer v) 1066 { 1067 PetscErrorCode ierr; 1068 DM_Shell *shell = (DM_Shell*)dm->data; 1069 1070 PetscFunctionBegin; 1071 ierr = VecView(shell->Xglobal,v);CHKERRQ(ierr); 1072 PetscFunctionReturn(0); 1073 } 1074 1075 static PetscErrorCode DMLoad_Shell(DM dm,PetscViewer v) 1076 { 1077 PetscErrorCode ierr; 1078 DM_Shell *shell = (DM_Shell*)dm->data; 1079 1080 PetscFunctionBegin; 1081 ierr = VecCreate(PetscObjectComm((PetscObject)dm),&shell->Xglobal);CHKERRQ(ierr); 1082 ierr = VecLoad(shell->Xglobal,v);CHKERRQ(ierr); 1083 PetscFunctionReturn(0); 1084 } 1085 1086 PetscErrorCode DMCreateSubDM_Shell(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 1087 { 1088 PetscErrorCode ierr; 1089 1090 PetscFunctionBegin; 1091 if (subdm) {ierr = DMShellCreate(PetscObjectComm((PetscObject) dm), subdm);CHKERRQ(ierr);} 1092 ierr = DMCreateSectionSubDM(dm, numFields, fields, is, subdm);CHKERRQ(ierr); 1093 PetscFunctionReturn(0); 1094 } 1095 1096 PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm) 1097 { 1098 PetscErrorCode ierr; 1099 DM_Shell *shell; 1100 1101 PetscFunctionBegin; 1102 ierr = PetscNewLog(dm,&shell);CHKERRQ(ierr); 1103 dm->data = shell; 1104 1105 dm->ops->destroy = DMDestroy_Shell; 1106 dm->ops->createglobalvector = DMCreateGlobalVector_Shell; 1107 dm->ops->createlocalvector = DMCreateLocalVector_Shell; 1108 dm->ops->creatematrix = DMCreateMatrix_Shell; 1109 dm->ops->view = DMView_Shell; 1110 dm->ops->load = DMLoad_Shell; 1111 dm->ops->globaltolocalbegin = DMGlobalToLocalBeginDefaultShell; 1112 dm->ops->globaltolocalend = DMGlobalToLocalEndDefaultShell; 1113 dm->ops->localtoglobalbegin = DMLocalToGlobalBeginDefaultShell; 1114 dm->ops->localtoglobalend = DMLocalToGlobalEndDefaultShell; 1115 dm->ops->localtolocalbegin = DMLocalToLocalBeginDefaultShell; 1116 dm->ops->localtolocalend = DMLocalToLocalEndDefaultShell; 1117 dm->ops->createsubdm = DMCreateSubDM_Shell; 1118 ierr = DMSetMatType(dm,MATDENSE);CHKERRQ(ierr); 1119 PetscFunctionReturn(0); 1120 } 1121 1122 /*@ 1123 DMShellCreate - Creates a shell DM object, used to manage user-defined problem data 1124 1125 Collective 1126 1127 Input Parameter: 1128 . comm - the processors that will share the global vector 1129 1130 Output Parameters: 1131 . shell - the shell DM 1132 1133 Level: advanced 1134 1135 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateLocalVector(), DMShellSetContext(), DMShellGetContext() 1136 @*/ 1137 PetscErrorCode DMShellCreate(MPI_Comm comm,DM *dm) 1138 { 1139 PetscErrorCode ierr; 1140 1141 PetscFunctionBegin; 1142 PetscValidPointer(dm,2); 1143 ierr = DMCreate(comm,dm);CHKERRQ(ierr); 1144 ierr = DMSetType(*dm,DMSHELL);CHKERRQ(ierr); 1145 ierr = DMSetUp(*dm);CHKERRQ(ierr); 1146 PetscFunctionReturn(0); 1147 } 1148