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