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