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 - 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 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 PetscFunctionBegin; 499 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 500 dm->ops->createlocalvector = func; 501 PetscFunctionReturn(0); 502 } 503 504 /*@C 505 DMShellSetGlobalToLocal - Sets the routines used to perform a global to local scatter 506 507 Logically Collective on dm 508 509 Input Arguments 510 + dm - the shell DM 511 . begin - the routine that begins the global to local scatter 512 - end - the routine that ends the global to local scatter 513 514 Notes: 515 If these functions are not provided but DMShellSetGlobalToLocalVecScatter() is called then 516 DMGlobalToLocalBeginDefaultShell()/DMGlobalToLocalEndDefaultShell() are used to to perform the transfers 517 518 Level: advanced 519 520 .seealso: DMShellSetLocalToGlobal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell() 521 @*/ 522 PetscErrorCode DMShellSetGlobalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) { 523 PetscFunctionBegin; 524 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 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 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 551 dm->ops->localtoglobalbegin = begin; 552 dm->ops->localtoglobalend = end; 553 PetscFunctionReturn(0); 554 } 555 556 /*@C 557 DMShellSetLocalToLocal - Sets the routines used to perform a local to local scatter 558 559 Logically Collective on dm 560 561 Input Arguments 562 + dm - the shell DM 563 . begin - the routine that begins the local to local scatter 564 - end - the routine that ends the local to local scatter 565 566 Notes: 567 If these functions are not provided but DMShellSetLocalToLocalVecScatter() is called then 568 DMLocalToLocalBeginDefaultShell()/DMLocalToLocalEndDefaultShell() are used to to perform the transfers 569 570 Level: advanced 571 572 .seealso: DMShellSetGlobalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell() 573 @*/ 574 PetscErrorCode DMShellSetLocalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) { 575 PetscFunctionBegin; 576 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 577 dm->ops->localtolocalbegin = begin; 578 dm->ops->localtolocalend = end; 579 PetscFunctionReturn(0); 580 } 581 582 /*@ 583 DMShellSetGlobalToLocalVecScatter - Sets a VecScatter context for global to local communication 584 585 Logically Collective on dm 586 587 Input Arguments 588 + dm - the shell DM 589 - gtol - the global to local VecScatter context 590 591 Level: advanced 592 593 .seealso: DMShellSetGlobalToLocal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell() 594 @*/ 595 PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol) 596 { 597 DM_Shell *shell = (DM_Shell*)dm->data; 598 PetscErrorCode ierr; 599 600 PetscFunctionBegin; 601 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 602 PetscValidHeaderSpecific(gtol,VEC_SCATTER_CLASSID,2); 603 ierr = PetscObjectReference((PetscObject)gtol);CHKERRQ(ierr); 604 ierr = VecScatterDestroy(&shell->gtol);CHKERRQ(ierr); 605 shell->gtol = gtol; 606 PetscFunctionReturn(0); 607 } 608 609 /*@ 610 DMShellSetLocalToGlobalVecScatter - Sets a VecScatter context for local to global communication 611 612 Logically Collective on dm 613 614 Input Arguments 615 + dm - the shell DM 616 - ltog - the local to global VecScatter context 617 618 Level: advanced 619 620 .seealso: DMShellSetLocalToGlobal(), DMLocalToGlobalBeginDefaultShell(), DMLocalToGlobalEndDefaultShell() 621 @*/ 622 PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog) 623 { 624 DM_Shell *shell = (DM_Shell*)dm->data; 625 PetscErrorCode ierr; 626 627 PetscFunctionBegin; 628 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 629 PetscValidHeaderSpecific(ltog,VEC_SCATTER_CLASSID,2); 630 ierr = PetscObjectReference((PetscObject)ltog);CHKERRQ(ierr); 631 ierr = VecScatterDestroy(&shell->ltog);CHKERRQ(ierr); 632 shell->ltog = ltog; 633 PetscFunctionReturn(0); 634 } 635 636 /*@ 637 DMShellSetLocalToLocalVecScatter - Sets a VecScatter context for local to local communication 638 639 Logically Collective on dm 640 641 Input Arguments 642 + dm - the shell DM 643 - ltol - the local to local VecScatter context 644 645 Level: advanced 646 647 .seealso: DMShellSetLocalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell() 648 @*/ 649 PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol) 650 { 651 DM_Shell *shell = (DM_Shell*)dm->data; 652 PetscErrorCode ierr; 653 654 PetscFunctionBegin; 655 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 656 PetscValidHeaderSpecific(ltol,VEC_SCATTER_CLASSID,2); 657 ierr = PetscObjectReference((PetscObject)ltol);CHKERRQ(ierr); 658 ierr = VecScatterDestroy(&shell->ltol);CHKERRQ(ierr); 659 shell->ltol = ltol; 660 PetscFunctionReturn(0); 661 } 662 663 /*@C 664 DMShellSetCoarsen - Set the routine used to coarsen the shell DM 665 666 Logically Collective on dm 667 668 Input Arguments 669 + dm - the shell DM 670 - coarsen - the routine that coarsens the DM 671 672 Level: advanced 673 674 .seealso: DMShellSetRefine(), DMCoarsen(), DMShellGetCoarsen(), DMShellSetContext(), DMShellGetContext() 675 @*/ 676 PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*)) 677 { 678 PetscErrorCode ierr; 679 PetscBool isshell; 680 681 PetscFunctionBegin; 682 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 683 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 684 if (!isshell) PetscFunctionReturn(0); 685 dm->ops->coarsen = coarsen; 686 PetscFunctionReturn(0); 687 } 688 689 /*@C 690 DMShellGetCoarsen - Get the routine used to coarsen the shell DM 691 692 Logically Collective on dm 693 694 Input Argument: 695 . dm - the shell DM 696 697 Output Argument: 698 . coarsen - the routine that coarsens the DM 699 700 Level: advanced 701 702 .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine() 703 @*/ 704 PetscErrorCode DMShellGetCoarsen(DM dm, PetscErrorCode (**coarsen)(DM,MPI_Comm,DM*)) 705 { 706 PetscErrorCode ierr; 707 PetscBool isshell; 708 709 PetscFunctionBegin; 710 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 711 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 712 if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 713 *coarsen = dm->ops->coarsen; 714 PetscFunctionReturn(0); 715 } 716 717 /*@C 718 DMShellSetRefine - Set the routine used to refine the shell DM 719 720 Logically Collective on dm 721 722 Input Arguments 723 + dm - the shell DM 724 - refine - the routine that refines the DM 725 726 Level: advanced 727 728 .seealso: DMShellSetCoarsen(), DMRefine(), DMShellGetRefine(), DMShellSetContext(), DMShellGetContext() 729 @*/ 730 PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM,MPI_Comm,DM*)) 731 { 732 PetscErrorCode ierr; 733 PetscBool isshell; 734 735 PetscFunctionBegin; 736 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 737 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 738 if (!isshell) PetscFunctionReturn(0); 739 dm->ops->refine = refine; 740 PetscFunctionReturn(0); 741 } 742 743 /*@C 744 DMShellGetRefine - Get the routine used to refine the shell DM 745 746 Logically Collective on dm 747 748 Input Argument: 749 . dm - the shell DM 750 751 Output Argument: 752 . refine - the routine that refines the DM 753 754 Level: advanced 755 756 .seealso: DMShellSetCoarsen(), DMCoarsen(), DMShellSetRefine(), DMRefine() 757 @*/ 758 PetscErrorCode DMShellGetRefine(DM dm, PetscErrorCode (**refine)(DM,MPI_Comm,DM*)) 759 { 760 PetscErrorCode ierr; 761 PetscBool isshell; 762 763 PetscFunctionBegin; 764 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 765 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 766 if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 767 *refine = dm->ops->refine; 768 PetscFunctionReturn(0); 769 } 770 771 /*@C 772 DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator 773 774 Logically Collective on dm 775 776 Input Arguments 777 + dm - the shell DM 778 - interp - the routine to create the interpolation 779 780 Level: advanced 781 782 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateInterpolation(), DMShellSetCreateRestriction(), DMShellSetContext(), DMShellGetContext() 783 @*/ 784 PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM,DM,Mat*,Vec*)) 785 { 786 PetscErrorCode ierr; 787 PetscBool isshell; 788 789 PetscFunctionBegin; 790 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 791 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 792 if (!isshell) PetscFunctionReturn(0); 793 dm->ops->createinterpolation = interp; 794 PetscFunctionReturn(0); 795 } 796 797 /*@C 798 DMShellGetCreateInterpolation - Get the routine used to create the interpolation operator 799 800 Logically Collective on dm 801 802 Input Argument: 803 + dm - the shell DM 804 805 Output Argument: 806 - interp - the routine to create the interpolation 807 808 Level: advanced 809 810 .seealso: DMShellGetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext() 811 @*/ 812 PetscErrorCode DMShellGetCreateInterpolation(DM dm, PetscErrorCode (**interp)(DM,DM,Mat*,Vec*)) 813 { 814 PetscErrorCode ierr; 815 PetscBool isshell; 816 817 PetscFunctionBegin; 818 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 819 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 820 if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 821 *interp = dm->ops->createinterpolation; 822 PetscFunctionReturn(0); 823 } 824 825 /*@C 826 DMShellSetCreateRestriction - Set the routine used to create the restriction operator 827 828 Logically Collective on dm 829 830 Input Arguments 831 + dm - the shell DM 832 - striction- the routine to create the restriction 833 834 Level: advanced 835 836 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellGetCreateRestriction(), DMShellSetContext(), DMShellGetContext() 837 @*/ 838 PetscErrorCode DMShellSetCreateRestriction(DM dm, PetscErrorCode (*restriction)(DM,DM,Mat*)) 839 { 840 PetscErrorCode ierr; 841 PetscBool isshell; 842 843 PetscFunctionBegin; 844 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 845 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 846 if (!isshell) PetscFunctionReturn(0); 847 dm->ops->createrestriction = restriction; 848 PetscFunctionReturn(0); 849 } 850 851 /*@C 852 DMShellGetCreateRestriction - Get the routine used to create the restriction operator 853 854 Logically Collective on dm 855 856 Input Argument: 857 + dm - the shell DM 858 859 Output Argument: 860 - restriction - the routine to create the restriction 861 862 Level: advanced 863 864 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetContext(), DMShellGetContext() 865 @*/ 866 PetscErrorCode DMShellGetCreateRestriction(DM dm, PetscErrorCode (**restriction)(DM,DM,Mat*)) 867 { 868 PetscErrorCode ierr; 869 PetscBool isshell; 870 871 PetscFunctionBegin; 872 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 873 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 874 if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 875 *restriction = dm->ops->createrestriction; 876 PetscFunctionReturn(0); 877 } 878 879 /*@C 880 DMShellSetCreateInjection - Set the routine used to create the injection operator 881 882 Logically Collective on dm 883 884 Input Arguments 885 + dm - the shell DM 886 - inject - the routine to create the injection 887 888 Level: advanced 889 890 .seealso: DMShellSetCreateInterpolation(), DMCreateInjection(), DMShellGetCreateInjection(), DMShellSetContext(), DMShellGetContext() 891 @*/ 892 PetscErrorCode DMShellSetCreateInjection(DM dm, PetscErrorCode (*inject)(DM,DM,Mat*)) 893 { 894 PetscErrorCode ierr; 895 PetscBool isshell; 896 897 PetscFunctionBegin; 898 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 899 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 900 if (!isshell) PetscFunctionReturn(0); 901 dm->ops->createinjection = inject; 902 PetscFunctionReturn(0); 903 } 904 905 /*@C 906 DMShellGetCreateInjection - Get the routine used to create the injection operator 907 908 Logically Collective on dm 909 910 Input Argument: 911 + dm - the shell DM 912 913 Output Argument: 914 - inject - the routine to create the injection 915 916 Level: advanced 917 918 .seealso: DMShellGetCreateInterpolation(), DMCreateInjection(), DMShellSetContext(), DMShellGetContext() 919 @*/ 920 PetscErrorCode DMShellGetCreateInjection(DM dm, PetscErrorCode (**inject)(DM,DM,Mat*)) 921 { 922 PetscErrorCode ierr; 923 PetscBool isshell; 924 925 PetscFunctionBegin; 926 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 927 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 928 if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 929 *inject = dm->ops->createinjection; 930 PetscFunctionReturn(0); 931 } 932 933 /*@C 934 DMShellSetCreateFieldDecomposition - Set the routine used to create a decomposition of fields for the shell DM 935 936 Logically Collective on dm 937 938 Input Arguments 939 + dm - the shell DM 940 - decomp - the routine to create the decomposition 941 942 Level: advanced 943 944 .seealso: DMCreateFieldDecomposition(), DMShellSetContext(), DMShellGetContext() 945 @*/ 946 PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,DM**)) 947 { 948 PetscErrorCode ierr; 949 PetscBool isshell; 950 951 PetscFunctionBegin; 952 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 953 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 954 if (!isshell) PetscFunctionReturn(0); 955 dm->ops->createfielddecomposition = decomp; 956 PetscFunctionReturn(0); 957 } 958 959 /*@C 960 DMShellSetCreateDomainDecomposition - Set the routine used to create a domain decomposition for the shell DM 961 962 Logically Collective on dm 963 964 Input Arguments 965 + dm - the shell DM 966 - decomp - the routine to create the decomposition 967 968 Level: advanced 969 970 .seealso: DMCreateDomainDecomposition(), DMShellSetContext(), DMShellGetContext() 971 @*/ 972 PetscErrorCode DMShellSetCreateDomainDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,IS**,DM**)) 973 { 974 PetscErrorCode ierr; 975 PetscBool isshell; 976 977 PetscFunctionBegin; 978 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 979 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 980 if (!isshell) PetscFunctionReturn(0); 981 dm->ops->createdomaindecomposition = decomp; 982 PetscFunctionReturn(0); 983 } 984 985 /*@C 986 DMShellSetCreateDomainDecompositionScatters - Set the routine used to create the scatter contexts for domain decomposition with a shell DM 987 988 Logically Collective on dm 989 990 Input Arguments 991 + dm - the shell DM 992 - scatter - the routine to create the scatters 993 994 Level: advanced 995 996 .seealso: DMCreateDomainDecompositionScatters(), DMShellSetContext(), DMShellGetContext() 997 @*/ 998 PetscErrorCode DMShellSetCreateDomainDecompositionScatters(DM dm, PetscErrorCode (*scatter)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**)) 999 { 1000 PetscErrorCode ierr; 1001 PetscBool isshell; 1002 1003 PetscFunctionBegin; 1004 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1005 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 1006 if (!isshell) PetscFunctionReturn(0); 1007 dm->ops->createddscatters = scatter; 1008 PetscFunctionReturn(0); 1009 } 1010 1011 /*@C 1012 DMShellSetCreateSubDM - Set the routine used to create a sub DM from the shell DM 1013 1014 Logically Collective on dm 1015 1016 Input Arguments 1017 + dm - the shell DM 1018 - subdm - the routine to create the decomposition 1019 1020 Level: advanced 1021 1022 .seealso: DMCreateSubDM(), DMShellGetCreateSubDM(), DMShellSetContext(), DMShellGetContext() 1023 @*/ 1024 PetscErrorCode DMShellSetCreateSubDM(DM dm, PetscErrorCode (*subdm)(DM,PetscInt,const PetscInt[],IS*,DM*)) 1025 { 1026 PetscErrorCode ierr; 1027 PetscBool isshell; 1028 1029 PetscFunctionBegin; 1030 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1031 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 1032 if (!isshell) PetscFunctionReturn(0); 1033 dm->ops->createsubdm = subdm; 1034 PetscFunctionReturn(0); 1035 } 1036 1037 /*@C 1038 DMShellGetCreateSubDM - Get the routine used to create a sub DM from the shell DM 1039 1040 Logically Collective on dm 1041 1042 Input Argument: 1043 + dm - the shell DM 1044 1045 Output Argument: 1046 - subdm - the routine to create the decomposition 1047 1048 Level: advanced 1049 1050 .seealso: DMCreateSubDM(), DMShellSetCreateSubDM(), DMShellSetContext(), DMShellGetContext() 1051 @*/ 1052 PetscErrorCode DMShellGetCreateSubDM(DM dm, PetscErrorCode (**subdm)(DM,PetscInt,const PetscInt[],IS*,DM*)) 1053 { 1054 PetscErrorCode ierr; 1055 PetscBool isshell; 1056 1057 PetscFunctionBegin; 1058 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1059 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 1060 if (!isshell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only use with DMSHELL type DMs"); 1061 *subdm = dm->ops->createsubdm; 1062 PetscFunctionReturn(0); 1063 } 1064 1065 static PetscErrorCode DMDestroy_Shell(DM dm) 1066 { 1067 PetscErrorCode ierr; 1068 DM_Shell *shell = (DM_Shell*)dm->data; 1069 1070 PetscFunctionBegin; 1071 ierr = MatDestroy(&shell->A);CHKERRQ(ierr); 1072 ierr = VecDestroy(&shell->Xglobal);CHKERRQ(ierr); 1073 ierr = VecDestroy(&shell->Xlocal);CHKERRQ(ierr); 1074 ierr = VecScatterDestroy(&shell->gtol);CHKERRQ(ierr); 1075 ierr = VecScatterDestroy(&shell->ltog);CHKERRQ(ierr); 1076 ierr = VecScatterDestroy(&shell->ltol);CHKERRQ(ierr); 1077 /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 1078 ierr = PetscFree(shell);CHKERRQ(ierr); 1079 PetscFunctionReturn(0); 1080 } 1081 1082 static PetscErrorCode DMView_Shell(DM dm,PetscViewer v) 1083 { 1084 PetscErrorCode ierr; 1085 DM_Shell *shell = (DM_Shell*)dm->data; 1086 1087 PetscFunctionBegin; 1088 ierr = VecView(shell->Xglobal,v);CHKERRQ(ierr); 1089 PetscFunctionReturn(0); 1090 } 1091 1092 static PetscErrorCode DMLoad_Shell(DM dm,PetscViewer v) 1093 { 1094 PetscErrorCode ierr; 1095 DM_Shell *shell = (DM_Shell*)dm->data; 1096 1097 PetscFunctionBegin; 1098 ierr = VecCreate(PetscObjectComm((PetscObject)dm),&shell->Xglobal);CHKERRQ(ierr); 1099 ierr = VecLoad(shell->Xglobal,v);CHKERRQ(ierr); 1100 PetscFunctionReturn(0); 1101 } 1102 1103 PetscErrorCode DMCreateSubDM_Shell(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 1104 { 1105 PetscErrorCode ierr; 1106 1107 PetscFunctionBegin; 1108 if (subdm) {ierr = DMShellCreate(PetscObjectComm((PetscObject) dm), subdm);CHKERRQ(ierr);} 1109 ierr = DMCreateSectionSubDM(dm, numFields, fields, is, subdm);CHKERRQ(ierr); 1110 PetscFunctionReturn(0); 1111 } 1112 1113 PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm) 1114 { 1115 PetscErrorCode ierr; 1116 DM_Shell *shell; 1117 1118 PetscFunctionBegin; 1119 ierr = PetscNewLog(dm,&shell);CHKERRQ(ierr); 1120 dm->data = shell; 1121 1122 dm->ops->destroy = DMDestroy_Shell; 1123 dm->ops->createglobalvector = DMCreateGlobalVector_Shell; 1124 dm->ops->createlocalvector = DMCreateLocalVector_Shell; 1125 dm->ops->creatematrix = DMCreateMatrix_Shell; 1126 dm->ops->view = DMView_Shell; 1127 dm->ops->load = DMLoad_Shell; 1128 dm->ops->globaltolocalbegin = DMGlobalToLocalBeginDefaultShell; 1129 dm->ops->globaltolocalend = DMGlobalToLocalEndDefaultShell; 1130 dm->ops->localtoglobalbegin = DMLocalToGlobalBeginDefaultShell; 1131 dm->ops->localtoglobalend = DMLocalToGlobalEndDefaultShell; 1132 dm->ops->localtolocalbegin = DMLocalToLocalBeginDefaultShell; 1133 dm->ops->localtolocalend = DMLocalToLocalEndDefaultShell; 1134 dm->ops->createsubdm = DMCreateSubDM_Shell; 1135 PetscFunctionReturn(0); 1136 } 1137 1138 /*@ 1139 DMShellCreate - Creates a shell DM object, used to manage user-defined problem data 1140 1141 Collective 1142 1143 Input Parameter: 1144 . comm - the processors that will share the global vector 1145 1146 Output Parameters: 1147 . shell - the shell DM 1148 1149 Level: advanced 1150 1151 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateLocalVector(), DMShellSetContext(), DMShellGetContext() 1152 @*/ 1153 PetscErrorCode DMShellCreate(MPI_Comm comm,DM *dm) 1154 { 1155 PetscErrorCode ierr; 1156 1157 PetscFunctionBegin; 1158 PetscValidPointer(dm,2); 1159 ierr = DMCreate(comm,dm);CHKERRQ(ierr); 1160 ierr = DMSetType(*dm,DMSHELL);CHKERRQ(ierr); 1161 ierr = DMSetUp(*dm);CHKERRQ(ierr); 1162 PetscFunctionReturn(0); 1163 } 1164