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 /* the check below is tacky and incomplete */ 198 if (dm->mattype) { 199 PetscBool flg,aij,seqaij,mpiaij; 200 ierr = PetscObjectTypeCompare((PetscObject)A,dm->mattype,&flg);CHKERRQ(ierr); 201 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&seqaij);CHKERRQ(ierr); 202 ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&mpiaij);CHKERRQ(ierr); 203 ierr = PetscStrcmp(dm->mattype,MATAIJ,&aij);CHKERRQ(ierr); 204 if (!flg) { 205 if (!(aij && (seqaij || mpiaij))) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_NOTSAMETYPE,"Requested matrix of type %s, but only %s available",dm->mattype,((PetscObject)A)->type_name); 206 } 207 } 208 if (((PetscObject)A)->refct < 2) { /* We have an exclusive reference so we can give it out */ 209 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 210 ierr = MatZeroEntries(A);CHKERRQ(ierr); 211 *J = A; 212 } else { /* Need to create a copy, could use MAT_SHARE_NONZERO_PATTERN in most cases */ 213 ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,J);CHKERRQ(ierr); 214 ierr = MatZeroEntries(*J);CHKERRQ(ierr); 215 } 216 PetscFunctionReturn(0); 217 } 218 219 PetscErrorCode DMCreateGlobalVector_Shell(DM dm,Vec *gvec) 220 { 221 PetscErrorCode ierr; 222 DM_Shell *shell = (DM_Shell*)dm->data; 223 Vec X; 224 225 PetscFunctionBegin; 226 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 227 PetscValidPointer(gvec,2); 228 *gvec = 0; 229 X = shell->Xglobal; 230 if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetGlobalVector() or DMShellSetCreateGlobalVector()"); 231 if (((PetscObject)X)->refct < 2) { /* We have an exclusive reference so we can give it out */ 232 ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 233 ierr = VecZeroEntries(X);CHKERRQ(ierr); 234 *gvec = X; 235 } else { /* Need to create a copy, could use MAT_SHARE_NONZERO_PATTERN in most cases */ 236 ierr = VecDuplicate(X,gvec);CHKERRQ(ierr); 237 ierr = VecZeroEntries(*gvec);CHKERRQ(ierr); 238 } 239 ierr = VecSetDM(*gvec,dm);CHKERRQ(ierr); 240 PetscFunctionReturn(0); 241 } 242 243 PetscErrorCode DMCreateLocalVector_Shell(DM dm,Vec *gvec) 244 { 245 PetscErrorCode ierr; 246 DM_Shell *shell = (DM_Shell*)dm->data; 247 Vec X; 248 249 PetscFunctionBegin; 250 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 251 PetscValidPointer(gvec,2); 252 *gvec = 0; 253 X = shell->Xlocal; 254 if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetLocalVector() or DMShellSetCreateLocalVector()"); 255 if (((PetscObject)X)->refct < 2) { /* We have an exclusive reference so we can give it out */ 256 ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 257 ierr = VecZeroEntries(X);CHKERRQ(ierr); 258 *gvec = X; 259 } else { /* Need to create a copy, could use MAT_SHARE_NONZERO_PATTERN in most cases */ 260 ierr = VecDuplicate(X,gvec);CHKERRQ(ierr); 261 ierr = VecZeroEntries(*gvec);CHKERRQ(ierr); 262 } 263 ierr = VecSetDM(*gvec,dm);CHKERRQ(ierr); 264 PetscFunctionReturn(0); 265 } 266 267 /*@ 268 DMShellSetContext - set some data to be usable by this DM 269 270 Collective 271 272 Input Arguments: 273 + dm - shell DM 274 - ctx - the context 275 276 Level: advanced 277 278 .seealso: DMCreateMatrix(), DMShellGetContext() 279 @*/ 280 PetscErrorCode DMShellSetContext(DM dm,void *ctx) 281 { 282 DM_Shell *shell = (DM_Shell*)dm->data; 283 PetscErrorCode ierr; 284 PetscBool isshell; 285 286 PetscFunctionBegin; 287 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 288 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 289 if (!isshell) PetscFunctionReturn(0); 290 shell->ctx = ctx; 291 PetscFunctionReturn(0); 292 } 293 294 /*@ 295 DMShellGetContext - set some data to be usable by this DM 296 297 Collective 298 299 Input Argument: 300 . dm - shell DM 301 302 Output Argument: 303 . ctx - the context 304 305 Level: advanced 306 307 .seealso: DMCreateMatrix(), DMShellSetContext() 308 @*/ 309 PetscErrorCode DMShellGetContext(DM dm,void **ctx) 310 { 311 DM_Shell *shell = (DM_Shell*)dm->data; 312 PetscErrorCode ierr; 313 PetscBool isshell; 314 315 PetscFunctionBegin; 316 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 317 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 318 if (!isshell) PetscFunctionReturn(0); 319 *ctx = shell->ctx; 320 PetscFunctionReturn(0); 321 } 322 323 /*@ 324 DMShellSetMatrix - sets a template matrix associated with the DMShell 325 326 Collective 327 328 Input Arguments: 329 + dm - shell DM 330 - J - template matrix 331 332 Level: advanced 333 334 .seealso: DMCreateMatrix(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext() 335 @*/ 336 PetscErrorCode DMShellSetMatrix(DM dm,Mat J) 337 { 338 DM_Shell *shell = (DM_Shell*)dm->data; 339 PetscErrorCode ierr; 340 PetscBool isshell; 341 342 PetscFunctionBegin; 343 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 344 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 345 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 346 if (!isshell) PetscFunctionReturn(0); 347 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 348 ierr = MatDestroy(&shell->A);CHKERRQ(ierr); 349 shell->A = J; 350 PetscFunctionReturn(0); 351 } 352 353 /*@C 354 DMShellSetCreateMatrix - sets the routine to create a matrix associated with the shell DM 355 356 Logically Collective on DM 357 358 Input Arguments: 359 + dm - the shell DM 360 - func - the function to create a matrix 361 362 Level: advanced 363 364 .seealso: DMCreateMatrix(), DMShellSetMatrix(), DMShellSetContext(), DMShellGetContext() 365 @*/ 366 PetscErrorCode DMShellSetCreateMatrix(DM dm,PetscErrorCode (*func)(DM,Mat*)) 367 { 368 369 PetscFunctionBegin; 370 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 371 dm->ops->creatematrix = func; 372 PetscFunctionReturn(0); 373 } 374 375 /*@ 376 DMShellSetGlobalVector - sets a template global vector associated with the DMShell 377 378 Logically Collective on DM 379 380 Input Arguments: 381 + dm - shell DM 382 - X - template vector 383 384 Level: advanced 385 386 .seealso: DMCreateGlobalVector(), DMShellSetMatrix(), DMShellSetCreateGlobalVector() 387 @*/ 388 PetscErrorCode DMShellSetGlobalVector(DM dm,Vec X) 389 { 390 DM_Shell *shell = (DM_Shell*)dm->data; 391 PetscErrorCode ierr; 392 PetscBool isshell; 393 DM vdm; 394 395 PetscFunctionBegin; 396 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 397 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 398 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 399 if (!isshell) PetscFunctionReturn(0); 400 ierr = VecGetDM(X,&vdm);CHKERRQ(ierr); 401 /* 402 if the vector proposed as the new base global vector for the DM is a DM vector associated 403 with the same DM then the current base global vector for the DM is ok and if we replace it with the new one 404 we get a circular dependency that prevents the DM from being destroy when it should be. 405 This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided 406 DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries 407 to set its input vector (which is associated with the DM) as the base global vector. 408 Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien 409 for pointing out the problem. 410 */ 411 if (vdm == dm) PetscFunctionReturn(0); 412 ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 413 ierr = VecDestroy(&shell->Xglobal);CHKERRQ(ierr); 414 shell->Xglobal = X; 415 PetscFunctionReturn(0); 416 } 417 418 /*@C 419 DMShellSetCreateGlobalVector - sets the routine to create a global vector associated with the shell DM 420 421 Logically Collective 422 423 Input Arguments: 424 + dm - the shell DM 425 - func - the creation routine 426 427 Level: advanced 428 429 .seealso: DMShellSetGlobalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext() 430 @*/ 431 PetscErrorCode DMShellSetCreateGlobalVector(DM dm,PetscErrorCode (*func)(DM,Vec*)) 432 { 433 434 PetscFunctionBegin; 435 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 436 dm->ops->createglobalvector = func; 437 PetscFunctionReturn(0); 438 } 439 440 /*@ 441 DMShellSetLocalVector - sets a template local vector associated with the DMShell 442 443 Logically Collective on DM 444 445 Input Arguments: 446 + dm - shell DM 447 - X - template vector 448 449 Level: advanced 450 451 .seealso: DMCreateLocalVector(), DMShellSetMatrix(), DMShellSetCreateLocalVector() 452 @*/ 453 PetscErrorCode DMShellSetLocalVector(DM dm,Vec X) 454 { 455 DM_Shell *shell = (DM_Shell*)dm->data; 456 PetscErrorCode ierr; 457 PetscBool isshell; 458 DM vdm; 459 460 PetscFunctionBegin; 461 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 462 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 463 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 464 if (!isshell) PetscFunctionReturn(0); 465 ierr = VecGetDM(X,&vdm);CHKERRQ(ierr); 466 /* 467 if the vector proposed as the new base global vector for the DM is a DM vector associated 468 with the same DM then the current base global vector for the DM is ok and if we replace it with the new one 469 we get a circular dependency that prevents the DM from being destroy when it should be. 470 This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided 471 DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries 472 to set its input vector (which is associated with the DM) as the base global vector. 473 Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien 474 for pointing out the problem. 475 */ 476 if (vdm == dm) PetscFunctionReturn(0); 477 ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 478 ierr = VecDestroy(&shell->Xlocal);CHKERRQ(ierr); 479 shell->Xlocal = X; 480 PetscFunctionReturn(0); 481 } 482 483 /*@C 484 DMShellSetCreateLocalVector - sets the routine to create a local vector associated with the shell DM 485 486 Logically Collective 487 488 Input Arguments: 489 + dm - the shell DM 490 - func - the creation routine 491 492 Level: advanced 493 494 .seealso: DMShellSetLocalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext() 495 @*/ 496 PetscErrorCode DMShellSetCreateLocalVector(DM dm,PetscErrorCode (*func)(DM,Vec*)) 497 { 498 499 PetscFunctionBegin; 500 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 501 dm->ops->createlocalvector = func; 502 PetscFunctionReturn(0); 503 } 504 505 /*@C 506 DMShellSetGlobalToLocal - Sets the routines used to perform a global to local scatter 507 508 Logically Collective on DM 509 510 Input Arguments 511 + dm - the shell DM 512 . begin - the routine that begins the global to local scatter 513 - end - the routine that ends the global to local scatter 514 515 Notes: 516 If these functions are not provided but DMShellSetGlobalToLocalVecScatter() is called then 517 DMGlobalToLocalBeginDefaultShell()/DMGlobalToLocalEndDefaultShell() are used to to perform the transfers 518 519 Level: advanced 520 521 .seealso: DMShellSetLocalToGlobal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell() 522 @*/ 523 PetscErrorCode DMShellSetGlobalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) { 524 PetscFunctionBegin; 525 dm->ops->globaltolocalbegin = begin; 526 dm->ops->globaltolocalend = end; 527 PetscFunctionReturn(0); 528 } 529 530 /*@C 531 DMShellSetLocalToGlobal - Sets the routines used to perform a local to global scatter 532 533 Logically Collective on DM 534 535 Input Arguments 536 + dm - the shell DM 537 . begin - the routine that begins the local to global scatter 538 - end - the routine that ends the local to global scatter 539 540 Notes: 541 If these functions are not provided but DMShellSetLocalToGlobalVecScatter() is called then 542 DMLocalToGlobalBeginDefaultShell()/DMLocalToGlobalEndDefaultShell() are used to to perform the transfers 543 544 Level: advanced 545 546 .seealso: DMShellSetGlobalToLocal() 547 @*/ 548 PetscErrorCode DMShellSetLocalToGlobal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) { 549 PetscFunctionBegin; 550 dm->ops->localtoglobalbegin = begin; 551 dm->ops->localtoglobalend = end; 552 PetscFunctionReturn(0); 553 } 554 555 /*@C 556 DMShellSetLocalToLocal - Sets the routines used to perform a local to local scatter 557 558 Logically Collective on DM 559 560 Input Arguments 561 + dm - the shell DM 562 . begin - the routine that begins the local to local scatter 563 - end - the routine that ends the local to local scatter 564 565 Notes: 566 If these functions are not provided but DMShellSetLocalToLocalVecScatter() is called then 567 DMLocalToLocalBeginDefaultShell()/DMLocalToLocalEndDefaultShell() are used to to perform the transfers 568 569 Level: advanced 570 571 .seealso: DMShellSetGlobalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell() 572 @*/ 573 PetscErrorCode DMShellSetLocalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) { 574 PetscFunctionBegin; 575 dm->ops->localtolocalbegin = begin; 576 dm->ops->localtolocalend = end; 577 PetscFunctionReturn(0); 578 } 579 580 /*@ 581 DMShellSetGlobalToLocalVecScatter - Sets a VecScatter context for global to local communication 582 583 Logically Collective on DM 584 585 Input Arguments 586 + dm - the shell DM 587 - gtol - the global to local VecScatter context 588 589 Level: advanced 590 591 .seealso: DMShellSetGlobalToLocal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell() 592 @*/ 593 PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol) 594 { 595 DM_Shell *shell = (DM_Shell*)dm->data; 596 PetscErrorCode ierr; 597 598 PetscFunctionBegin; 599 ierr = PetscObjectReference((PetscObject)gtol);CHKERRQ(ierr); 600 /* Call VecScatterDestroy() to avoid a memory leak in case of re-setting. */ 601 ierr = VecScatterDestroy(&shell->gtol);CHKERRQ(ierr); 602 shell->gtol = gtol; 603 PetscFunctionReturn(0); 604 } 605 606 /*@ 607 DMShellSetLocalToGlobalVecScatter - Sets a VecScatter context for local to global communication 608 609 Logically Collective on DM 610 611 Input Arguments 612 + dm - the shell DM 613 - ltog - the local to global VecScatter context 614 615 Level: advanced 616 617 .seealso: DMShellSetLocalToGlobal(), DMLocalToGlobalBeginDefaultShell(), DMLocalToGlobalEndDefaultShell() 618 @*/ 619 PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog) 620 { 621 DM_Shell *shell = (DM_Shell*)dm->data; 622 PetscErrorCode ierr; 623 624 PetscFunctionBegin; 625 ierr = PetscObjectReference((PetscObject)ltog);CHKERRQ(ierr); 626 /* Call VecScatterDestroy() to avoid a memory leak in case of re-setting. */ 627 ierr = VecScatterDestroy(&shell->ltog);CHKERRQ(ierr); 628 shell->ltog = ltog; 629 PetscFunctionReturn(0); 630 } 631 632 /*@ 633 DMShellSetLocalToLocalVecScatter - Sets a VecScatter context for local to local communication 634 635 Logically Collective on DM 636 637 Input Arguments 638 + dm - the shell DM 639 - ltol - the local to local VecScatter context 640 641 Level: advanced 642 643 .seealso: DMShellSetLocalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell() 644 @*/ 645 PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol) 646 { 647 DM_Shell *shell = (DM_Shell*)dm->data; 648 PetscErrorCode ierr; 649 650 PetscFunctionBegin; 651 ierr = PetscObjectReference((PetscObject)ltol);CHKERRQ(ierr); 652 /* Call VecScatterDestroy() to avoid a memory leak in case of re-setting. */ 653 ierr = VecScatterDestroy(&shell->ltol);CHKERRQ(ierr); 654 shell->ltol = ltol; 655 PetscFunctionReturn(0); 656 } 657 658 /*@C 659 DMShellSetCoarsen - Set the routine used to coarsen the shell DM 660 661 Logically Collective on DM 662 663 Input Arguments 664 + dm - the shell DM 665 - coarsen - the routine that coarsens the DM 666 667 Level: advanced 668 669 .seealso: DMShellSetRefine(), DMCoarsen(), DMShellSetContext(), DMShellGetContext() 670 @*/ 671 PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*)) 672 { 673 PetscErrorCode ierr; 674 PetscBool isshell; 675 676 PetscFunctionBegin; 677 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 678 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 679 if (!isshell) PetscFunctionReturn(0); 680 dm->ops->coarsen = coarsen; 681 PetscFunctionReturn(0); 682 } 683 684 /*@C 685 DMShellGetCoarsen - Get the routine used to coarsen the shell DM 686 687 Logically Collective on DM 688 689 Input Arguments 690 . dm - the shell DM 691 692 Output Arguments 693 . coarsen - the routine that coarsens the DM 694 695 Level: advanced 696 697 .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine() 698 @*/ 699 PetscErrorCode DMShellGetCoarsen(DM dm, PetscErrorCode (**coarsen)(DM,MPI_Comm,DM*)) 700 { 701 PetscErrorCode ierr; 702 PetscBool isshell; 703 704 PetscFunctionBegin; 705 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 706 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 707 if (!isshell) PetscFunctionReturn(0); 708 *coarsen = dm->ops->coarsen; 709 PetscFunctionReturn(0); 710 } 711 712 /*@C 713 DMShellSetRefine - Set the routine used to refine the shell DM 714 715 Logically Collective on DM 716 717 Input Arguments 718 + dm - the shell DM 719 - refine - the routine that refines the DM 720 721 Level: advanced 722 723 .seealso: DMShellSetCoarsen(), DMRefine(), DMShellSetContext(), DMShellGetContext() 724 @*/ 725 PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM,MPI_Comm,DM*)) 726 { 727 PetscErrorCode ierr; 728 PetscBool isshell; 729 730 PetscFunctionBegin; 731 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 732 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 733 if (!isshell) PetscFunctionReturn(0); 734 dm->ops->refine = refine; 735 PetscFunctionReturn(0); 736 } 737 738 /*@C 739 DMShellGetRefine - Get the routine used to refine the shell DM 740 741 Logically Collective on DM 742 743 Input Arguments 744 . dm - the shell DM 745 746 Output Arguments 747 . refine - the routine that refines the DM 748 749 Level: advanced 750 751 .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine() 752 @*/ 753 PetscErrorCode DMShellGetRefine(DM dm, PetscErrorCode (**refine)(DM,MPI_Comm,DM*)) 754 { 755 PetscErrorCode ierr; 756 PetscBool isshell; 757 758 PetscFunctionBegin; 759 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 760 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 761 if (!isshell) PetscFunctionReturn(0); 762 *refine = dm->ops->refine; 763 PetscFunctionReturn(0); 764 } 765 766 /*@C 767 DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator 768 769 Logically Collective on DM 770 771 Input Arguments 772 + dm - the shell DM 773 - interp - the routine to create the interpolation 774 775 Level: advanced 776 777 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetCreateRestriction(), DMShellSetContext(), DMShellGetContext() 778 @*/ 779 PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM,DM,Mat*,Vec*)) 780 { 781 PetscErrorCode ierr; 782 PetscBool isshell; 783 784 PetscFunctionBegin; 785 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 786 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 787 if (!isshell) PetscFunctionReturn(0); 788 dm->ops->createinterpolation = interp; 789 PetscFunctionReturn(0); 790 } 791 792 /*@C 793 DMShellGetCreateInterpolation - Get the routine used to create the interpolation operator 794 795 Logically Collective on DM 796 797 Input Arguments 798 + dm - the shell DM 799 - interp - the routine to create the interpolation 800 801 Level: advanced 802 803 .seealso: DMShellGetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext() 804 @*/ 805 PetscErrorCode DMShellGetCreateInterpolation(DM dm, PetscErrorCode (**interp)(DM,DM,Mat*,Vec*)) 806 { 807 PetscErrorCode ierr; 808 PetscBool isshell; 809 810 PetscFunctionBegin; 811 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 812 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 813 if (!isshell) PetscFunctionReturn(0); 814 *interp = dm->ops->createinterpolation; 815 PetscFunctionReturn(0); 816 } 817 818 /*@C 819 DMShellSetCreateRestriction - Set the routine used to create the restriction operator 820 821 Logically Collective on DM 822 823 Input Arguments 824 + dm - the shell DM 825 - striction- the routine to create the restriction 826 827 Level: advanced 828 829 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetContext(), DMShellGetContext() 830 @*/ 831 PetscErrorCode DMShellSetCreateRestriction(DM dm, PetscErrorCode (*restriction)(DM,DM,Mat*)) 832 { 833 PetscErrorCode ierr; 834 PetscBool isshell; 835 836 PetscFunctionBegin; 837 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 838 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 839 if (!isshell) PetscFunctionReturn(0); 840 dm->ops->createrestriction = restriction; 841 PetscFunctionReturn(0); 842 } 843 844 /*@C 845 DMShellGetCreateRestriction - Get the routine used to create the restriction operator 846 847 Logically Collective on DM 848 849 Input Arguments 850 + dm - the shell DM 851 - striction- the routine to create the restriction 852 853 Level: advanced 854 855 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetContext(), DMShellGetContext() 856 @*/ 857 PetscErrorCode DMShellGetCreateRestriction(DM dm, PetscErrorCode (**restriction)(DM,DM,Mat*)) 858 { 859 PetscErrorCode ierr; 860 PetscBool isshell; 861 862 PetscFunctionBegin; 863 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 864 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 865 if (!isshell) PetscFunctionReturn(0); 866 *restriction = dm->ops->createrestriction; 867 PetscFunctionReturn(0); 868 } 869 870 /*@C 871 DMShellSetCreateInjection - Set the routine used to create the injection operator 872 873 Logically Collective on DM 874 875 Input Arguments 876 + dm - the shell DM 877 - inject - the routine to create the injection 878 879 Level: advanced 880 881 .seealso: DMShellSetCreateInterpolation(), DMCreateInjection(), DMShellSetContext(), DMShellGetContext() 882 @*/ 883 PetscErrorCode DMShellSetCreateInjection(DM dm, PetscErrorCode (*inject)(DM,DM,Mat*)) 884 { 885 PetscErrorCode ierr; 886 PetscBool isshell; 887 888 PetscFunctionBegin; 889 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 890 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 891 if (!isshell) PetscFunctionReturn(0); 892 dm->ops->getinjection = inject; 893 PetscFunctionReturn(0); 894 } 895 896 /*@C 897 DMShellGetCreateInjection - Get the routine used to create the injection operator 898 899 Logically Collective on DM 900 901 Input Arguments 902 + dm - the shell DM 903 - inject - the routine to create the injection 904 905 Level: advanced 906 907 .seealso: DMShellGetCreateInterpolation(), DMCreateInjection(), DMShellSetContext(), DMShellGetContext() 908 @*/ 909 PetscErrorCode DMShellGetCreateInjection(DM dm, PetscErrorCode (**inject)(DM,DM,Mat*)) 910 { 911 PetscErrorCode ierr; 912 PetscBool isshell; 913 914 PetscFunctionBegin; 915 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 916 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 917 if (!isshell) PetscFunctionReturn(0); 918 *inject = dm->ops->getinjection; 919 PetscFunctionReturn(0); 920 } 921 922 static PetscErrorCode DMHasCreateInjection_Shell(DM dm, PetscBool *flg) 923 { 924 PetscErrorCode ierr; 925 PetscBool isshell; 926 927 PetscFunctionBegin; 928 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 929 PetscValidPointer(flg,2); 930 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 931 if (!isshell) { 932 *flg = PETSC_FALSE; 933 } else { 934 *flg = dm->ops->getinjection ? PETSC_TRUE : PETSC_FALSE; 935 } 936 PetscFunctionReturn(0); 937 } 938 /*@C 939 DMShellSetCreateFieldDecomposition - Set the routine used to create a decomposition of fields for the shell DM 940 941 Logically Collective on DM 942 943 Input Arguments 944 + dm - the shell DM 945 - decomp - the routine to create the decomposition 946 947 Level: advanced 948 949 .seealso: DMCreateFieldDecomposition(), DMShellSetContext(), DMShellGetContext() 950 @*/ 951 PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,DM**)) 952 { 953 PetscErrorCode ierr; 954 PetscBool isshell; 955 956 PetscFunctionBegin; 957 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 958 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 959 if (!isshell) PetscFunctionReturn(0); 960 dm->ops->createfielddecomposition = decomp; 961 PetscFunctionReturn(0); 962 } 963 964 /*@C 965 DMShellSetCreateDomainDecomposition - Set the routine used to create a domain decomposition for the shell DM 966 967 Logically Collective on DM 968 969 Input Arguments 970 + dm - the shell DM 971 - decomp - the routine to create the decomposition 972 973 Level: advanced 974 975 .seealso: DMCreateDomainDecomposition(), DMShellSetContext(), DMShellGetContext() 976 @*/ 977 PetscErrorCode DMShellSetCreateDomainDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,IS**,DM**)) 978 { 979 PetscErrorCode ierr; 980 PetscBool isshell; 981 982 PetscFunctionBegin; 983 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 984 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 985 if (!isshell) PetscFunctionReturn(0); 986 dm->ops->createdomaindecomposition = decomp; 987 PetscFunctionReturn(0); 988 } 989 990 /*@C 991 DMShellSetCreateDomainDecompositionScatters - Set the routine used to create the scatter contexts for domain decomposition with a shell DM 992 993 Logically Collective on DM 994 995 Input Arguments 996 + dm - the shell DM 997 - scatter - the routine to create the scatters 998 999 Level: advanced 1000 1001 .seealso: DMCreateDomainDecompositionScatters(), DMShellSetContext(), DMShellGetContext() 1002 @*/ 1003 PetscErrorCode DMShellSetCreateDomainDecompositionScatters(DM dm, PetscErrorCode (*scatter)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**)) 1004 { 1005 PetscErrorCode ierr; 1006 PetscBool isshell; 1007 1008 PetscFunctionBegin; 1009 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1010 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 1011 if (!isshell) PetscFunctionReturn(0); 1012 dm->ops->createddscatters = scatter; 1013 PetscFunctionReturn(0); 1014 } 1015 1016 /*@C 1017 DMShellSetCreateSubDM - Set the routine used to create a sub DM from the shell DM 1018 1019 Logically Collective on DM 1020 1021 Input Arguments 1022 + dm - the shell DM 1023 - subdm - the routine to create the decomposition 1024 1025 Level: advanced 1026 1027 .seealso: DMCreateSubDM(), DMShellSetContext(), DMShellGetContext() 1028 @*/ 1029 PetscErrorCode DMShellSetCreateSubDM(DM dm, PetscErrorCode (*subdm)(DM,PetscInt,const PetscInt[],IS*,DM*)) 1030 { 1031 PetscErrorCode ierr; 1032 PetscBool isshell; 1033 1034 PetscFunctionBegin; 1035 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1036 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 1037 if (!isshell) PetscFunctionReturn(0); 1038 dm->ops->createsubdm = subdm; 1039 PetscFunctionReturn(0); 1040 } 1041 1042 /*@C 1043 DMShellGetCreateSubDM - Get the routine used to create a sub DM from the shell DM 1044 1045 Logically Collective on DM 1046 1047 Input Arguments 1048 + dm - the shell DM 1049 - subdm - the routine to create the decomposition 1050 1051 Level: advanced 1052 1053 .seealso: DMCreateSubDM(), DMShellSetCreateSubDM(), DMShellSetContext(), DMShellGetContext() 1054 @*/ 1055 PetscErrorCode DMShellGetCreateSubDM(DM dm, PetscErrorCode (**subdm)(DM,PetscInt,const PetscInt[],IS*,DM*)) 1056 { 1057 PetscErrorCode ierr; 1058 PetscBool isshell; 1059 1060 PetscFunctionBegin; 1061 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1062 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 1063 if (!isshell) PetscFunctionReturn(0); 1064 *subdm = dm->ops->createsubdm; 1065 PetscFunctionReturn(0); 1066 } 1067 1068 static PetscErrorCode DMDestroy_Shell(DM dm) 1069 { 1070 PetscErrorCode ierr; 1071 DM_Shell *shell = (DM_Shell*)dm->data; 1072 1073 PetscFunctionBegin; 1074 ierr = MatDestroy(&shell->A);CHKERRQ(ierr); 1075 ierr = VecDestroy(&shell->Xglobal);CHKERRQ(ierr); 1076 ierr = VecDestroy(&shell->Xlocal);CHKERRQ(ierr); 1077 ierr = VecScatterDestroy(&shell->gtol);CHKERRQ(ierr); 1078 ierr = VecScatterDestroy(&shell->ltog);CHKERRQ(ierr); 1079 ierr = VecScatterDestroy(&shell->ltol);CHKERRQ(ierr); 1080 /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 1081 ierr = PetscFree(shell);CHKERRQ(ierr); 1082 PetscFunctionReturn(0); 1083 } 1084 1085 static PetscErrorCode DMView_Shell(DM dm,PetscViewer v) 1086 { 1087 PetscErrorCode ierr; 1088 DM_Shell *shell = (DM_Shell*)dm->data; 1089 1090 PetscFunctionBegin; 1091 ierr = VecView(shell->Xglobal,v);CHKERRQ(ierr); 1092 PetscFunctionReturn(0); 1093 } 1094 1095 static PetscErrorCode DMLoad_Shell(DM dm,PetscViewer v) 1096 { 1097 PetscErrorCode ierr; 1098 DM_Shell *shell = (DM_Shell*)dm->data; 1099 1100 PetscFunctionBegin; 1101 ierr = VecCreate(PetscObjectComm((PetscObject)dm),&shell->Xglobal);CHKERRQ(ierr); 1102 ierr = VecLoad(shell->Xglobal,v);CHKERRQ(ierr); 1103 PetscFunctionReturn(0); 1104 } 1105 1106 PetscErrorCode DMCreateSubDM_Shell(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 1107 { 1108 PetscErrorCode ierr; 1109 1110 PetscFunctionBegin; 1111 if (subdm) {ierr = DMShellCreate(PetscObjectComm((PetscObject) dm), subdm);CHKERRQ(ierr);} 1112 ierr = DMCreateSectionSubDM(dm, numFields, fields, is, subdm);CHKERRQ(ierr); 1113 PetscFunctionReturn(0); 1114 } 1115 1116 PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm) 1117 { 1118 PetscErrorCode ierr; 1119 DM_Shell *shell; 1120 1121 PetscFunctionBegin; 1122 ierr = PetscNewLog(dm,&shell);CHKERRQ(ierr); 1123 dm->data = shell; 1124 1125 ierr = PetscObjectChangeTypeName((PetscObject)dm,DMSHELL);CHKERRQ(ierr); 1126 1127 dm->ops->destroy = DMDestroy_Shell; 1128 dm->ops->createglobalvector = DMCreateGlobalVector_Shell; 1129 dm->ops->createlocalvector = DMCreateLocalVector_Shell; 1130 dm->ops->creatematrix = DMCreateMatrix_Shell; 1131 dm->ops->view = DMView_Shell; 1132 dm->ops->load = DMLoad_Shell; 1133 dm->ops->globaltolocalbegin = DMGlobalToLocalBeginDefaultShell; 1134 dm->ops->globaltolocalend = DMGlobalToLocalEndDefaultShell; 1135 dm->ops->localtoglobalbegin = DMLocalToGlobalBeginDefaultShell; 1136 dm->ops->localtoglobalend = DMLocalToGlobalEndDefaultShell; 1137 dm->ops->localtolocalbegin = DMLocalToLocalBeginDefaultShell; 1138 dm->ops->localtolocalend = DMLocalToLocalEndDefaultShell; 1139 dm->ops->createsubdm = DMCreateSubDM_Shell; 1140 dm->ops->hascreateinjection = DMHasCreateInjection_Shell; 1141 PetscFunctionReturn(0); 1142 } 1143 1144 /*@ 1145 DMShellCreate - Creates a shell DM object, used to manage user-defined problem data 1146 1147 Collective on MPI_Comm 1148 1149 Input Parameter: 1150 . comm - the processors that will share the global vector 1151 1152 Output Parameters: 1153 . shell - the shell DM 1154 1155 Level: advanced 1156 1157 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateLocalVector(), DMShellSetContext(), DMShellGetContext() 1158 @*/ 1159 PetscErrorCode DMShellCreate(MPI_Comm comm,DM *dm) 1160 { 1161 PetscErrorCode ierr; 1162 1163 PetscFunctionBegin; 1164 PetscValidPointer(dm,2); 1165 ierr = DMCreate(comm,dm);CHKERRQ(ierr); 1166 ierr = DMSetType(*dm,DMSHELL);CHKERRQ(ierr); 1167 ierr = DMSetUp(*dm);CHKERRQ(ierr); 1168 PetscFunctionReturn(0); 1169 } 1170 1171